Title: | Statistical Depth Functions for Multivariate Analysis |
---|---|
Description: | Data depth concept offers a variety of powerful and user friendly tools for robust exploration and inference for multivariate data. The offered techniques may be successfully used in cases of lack of our knowledge on parametric models generating data due to their nature. The package consist of among others implementations of several data depth techniques involving multivariate quantile-quantile plots, multivariate scatter estimators, multivariate Wilcoxon tests and robust regressions. |
Authors: | Zygmunt Zawadzki [aut, cre], Daniel Kosiorowski [aut], Krzysztof Slomczynski [ctb], Mateusz Bocian [ctb], Anna Wegrzynkiewicz [ctb] |
Maintainer: | Zygmunt Zawadzki <[email protected]> |
License: | GPL-2 |
Version: | 2.1.5 |
Built: | 2024-11-11 07:26:46 UTC |
Source: | CRAN |
Add fitted line to a plot. This is overloaded function for robust regression methods from package depthproc.
## S4 method for signature 'RobReg' abline( a = NULL, b = NULL, h = NULL, v = NULL, reg = NULL, coef = NULL, untf = FALSE, ... )
## S4 method for signature 'RobReg' abline( a = NULL, b = NULL, h = NULL, v = NULL, reg = NULL, coef = NULL, untf = FALSE, ... )
a |
an object of class RobReg |
b |
not used. |
h |
not supported. |
v |
not supported. |
reg |
not supported. |
coef |
not supported. |
untf |
not supported. |
... |
Arguments to be passed to methods, such as graphical parameters (see par). |
Create a matrix from DepthCurve and DepthCurveList.
as.matrix(x, ...) ## S4 method for signature 'DepthCurveList' as.matrix(x)
as.matrix(x, ...) ## S4 method for signature 'DepthCurveList' as.matrix(x)
x |
an object of class that inherits from DepthCurveList (ScaleCurveList or AsymmetryCurveList). |
... |
other arguments passed to standard as.matrix function. |
Produces an asymmetry curve estimated from given data.
asymmetryCurve( x, y = NULL, alpha = seq(0, 1, 0.01), movingmedian = FALSE, name = "X", name_y = "Y", depth_params = list(method = "Projection") )
asymmetryCurve( x, y = NULL, alpha = seq(0, 1, 0.01), movingmedian = FALSE, name = "X", name_y = "Y", depth_params = list(method = "Projection") )
x |
The data as a matrix or data frame. If it is a matrix or data frame, then each row is viewed as one multivariate observation. |
y |
Additional matrix of multivariate data. |
alpha |
An ordered vector containing indices of central regins used for asymmetry curve calculation. |
movingmedian |
Logical. For default FALSE only one depth median is used to compute asymmetry norm. If TRUE — for every central area, a new depth median will be used — this approach needs much more time. |
name |
Name of set X — used in plot legend |
name_y |
Name of set Y — used in plot legend |
depth_params |
list of parameters for function depth (method, threads, ndir, la, lb, pdim, mean, cov, exact). |
method |
Character string which determines the depth function used. The method can be "Projection" (the default), "Mahalanobis", "Euclidean", "Tukey" or "LP". For details see |
For sample depth function ,
,
,
,
denoting
— central region, we can define the asymmetry curve
, for
being nonparametric scale and asymmetry functional correspondingly, where
— denotes constant,
— denotes mean vector, denotes multivariate median induced by depth function and
— denotes a volume.
Asymmetry curve takes uses function convhulln from package geometry for computing a volume of convex hull containing central region.
Daniel Kosiorowski, Mateusz Bocian, Anna Wegrzynkiewicz and Zygmunt Zawadzki from Cracow University of Economics.
Serfling R. J. Multivariate Symmetry and Asymmetry, Encyclopedia of Statistical Science, S Kotz, C.B. Read, N. Balakrishnan, B. Vidakovic (eds), 2nd, ed., John Wiley.
Liu, R.Y., Parelius, J.M. and Singh, K. (1999), Multivariate analysis by data depth: Descriptive statistics, graphics and inference (with discussion), Ann. Statist., 27, 783–858.
Chaudhuri, P. (1996), On a Geometric Notion of Quantiles for Multivariate Data, Journal of the American Statistical Association, 862–872.
Dyckerhoff, R. (2004), Data Depths Satisfying the Projection Property, Allgemeines Statistisches Archiv., 88, 163–190.
# EXAMPLE 1 library(sn) xi <- c(0, 0) alpha <- c(2, -5) Omega <- diag(2) * 5 n <- 500 X <- mvrnorm(n, xi, Omega) # normal distribution Y <- rmst(n, xi, Omega, alpha, nu = 1) asymmetryCurve(X, Y, name = "NORM", name_y = "S_T(2, -5, 10)") # EXAMPLE 2 data(under5.mort) data(inf.mort) data(maesles.imm) data1990 <- cbind(under5.mort[, 1], inf.mort[, 1], maesles.imm[, 1]) data2011 <- cbind(under5.mort[, 22], inf.mort[, 22], maesles.imm[, 22]) as1990 <- asymmetryCurve(data1990, name = "scale curve 1990") as2011 <- asymmetryCurve(data2011, name = "scale curve 2011") figure <- getPlot(combineDepthCurves(as1990, as2011)) + ggtitle("Scale curves") figure
# EXAMPLE 1 library(sn) xi <- c(0, 0) alpha <- c(2, -5) Omega <- diag(2) * 5 n <- 500 X <- mvrnorm(n, xi, Omega) # normal distribution Y <- rmst(n, xi, Omega, alpha, nu = 1) asymmetryCurve(X, Y, name = "NORM", name_y = "S_T(2, -5, 10)") # EXAMPLE 2 data(under5.mort) data(inf.mort) data(maesles.imm) data1990 <- cbind(under5.mort[, 1], inf.mort[, 1], maesles.imm[, 1]) data2011 <- cbind(under5.mort[, 22], inf.mort[, 22], maesles.imm[, 22]) as1990 <- asymmetryCurve(data1990, name = "scale curve 1990") as2011 <- asymmetryCurve(data2011, name = "scale curve 2011") figure <- getPlot(combineDepthCurves(as1990, as2011)) + ggtitle("Scale curves") figure
AsymmetryCurve is a class that stores results of asymmetryCurve function.
The mechanism of creating plots with multiple curves is shown in DepthCurve-class (same mechanism is applied for ScaleCurve).
Class that stores result of function binningDepth2D(...)
freq
Matrix with number of elements in certain bin.
mid_x
Middle values on x-axis.
mid_y
Middle values on y-axis.
breaks_x
Boundaries of bins.
breaks_y
Boundaries of bins.
input_data
Binned data.
max_depth_x
Point with maximum depth on x-axis.
max_depth_y
Point with maximum depth on y-axis.
A robust method of decreasing a sample size and therefore a complexity of a statistical procedure. The method may be used within a kernel density or a predictive distribution estimation.
binningDepth2D( x, binmethod = "LocDepth", nbins = 8, k = 1, remove_borders = FALSE, depth_params = list(method = "LP") )
binningDepth2D( x, binmethod = "LocDepth", nbins = 8, k = 1, remove_borders = FALSE, depth_params = list(method = "LP") )
x |
bivariate matrix containing data. Each row is viewed as one two-dimensional observation. |
binmethod |
A method for calculation center and dispersion measures. "LocDepth" uses location-scale depth, MAD uses median and MAD in each dimension. |
nbins |
number of bins in each dimension |
k |
responsible for tightness of bins. |
remove_borders |
Logical, include or not marginal bins |
depth_params |
other arguments passed to depthMedian |
Let us recall, that binning is a popular method of decreasing a sample size. To bin a window of points
to a grid
we simply assign each sample point
to the nearest grid point
. When binning is completed, each grid point
has an associated number
, which is the sum of all the points that have been assigned to
. This procedure replaces the data
with the smaller set
. Although simple binning can speed up the computation, it is criticized for a lack of precise approximate control over the accuracy of the approximation. Robust binning however stresses properties of the majority of the data and decreases the computational complexity of the DSA at the same time.
For a 1D window , let
denote a 2D window created basing on
and consisted of
pairs of observations and the
lagged observations
. Robust 2D binning of the
is a very useful technique in a context of robust estimation of the predictive distribution of a time series (see Kosiorowski:2013b).
Assume we analyze a data stream using a moving window of a fixed length
, i.e.,
and the derivative window
. In a first step we calculate the weighted sample
depth for
. Next we choose equally spaced grid of points
in this way that
covers fraction of the
central points of
w.r.t. the calculated
depth, i.e., it covers
for certain prefixed threshold
. For both
and
we perform a simple binning using following bins:
,
, ...,
. For robust binning we reject "border" classes and further use only midpoints and binned frequencies for classes
,
, ...,
.
freq: a matrix containing the binned frequencies
mid_x: mid points for x
mid_y: mid points for y
breaks_x: breaks for x
breaks_y: breaks for y
input_data: max_depth_x and max_depth_y:
Daniel Kosiorowski and Zygmunt Zawadzki from Cracow University of Economics.
Hall, P., Wand, M. P. (1996) On the Accuracy of Binned Kernel Density Estimators, Journal of Multivariate Analysis archive, Volume 56 Issue 2, 165–184
Holmstrom, L. (2000) The Accuracy and the Computational Complexity of a Multivariate Binned Kernel Density Estimator, Journal of Multivariate Analysis, Volume 72, Issue 2, 264–309, doi:10.1006/jmva.1999.1863. (https://www.sciencedirect.com/science/article/pii/S0047259X99918638)
# EXAMPLE 1 Sigma1 <- matrix(c(10, 3, 3, 2), 2, 2) X1 <- mvrnorm(n = 8500, mu = c(0, 0), Sigma1) Sigma2 <- matrix(c(10, 0, 0, 2), 2, 2) X2 <- mvrnorm(n = 1500, mu = c(-10, 6), Sigma2) BALLOT <- rbind(X1, X2) train <- sample(1:10000, 500) data <- BALLOT[train, ] plot(data) b1 <- binningDepth2D(data, remove_borders = FALSE, nbins = 12, k = 1) b2 <- binningDepth2D(data, nbins = 12, k = 1, remove_borders = TRUE) plot(b1) plot(b2) # EXAMPLE 2 data(under5.mort) data(maesles.imm) data2011 <- cbind(under5.mort[, 22], maesles.imm[, 22]) plot(binningDepth2D(data2011, nbins = 8, k = 0.5, remove_borders = TRUE))
# EXAMPLE 1 Sigma1 <- matrix(c(10, 3, 3, 2), 2, 2) X1 <- mvrnorm(n = 8500, mu = c(0, 0), Sigma1) Sigma2 <- matrix(c(10, 0, 0, 2), 2, 2) X2 <- mvrnorm(n = 1500, mu = c(-10, 6), Sigma2) BALLOT <- rbind(X1, X2) train <- sample(1:10000, 500) data <- BALLOT[train, ] plot(data) b1 <- binningDepth2D(data, remove_borders = FALSE, nbins = 12, k = 1) b2 <- binningDepth2D(data, nbins = 12, k = 1, remove_borders = TRUE) plot(b1) plot(b2) # EXAMPLE 2 data(under5.mort) data(maesles.imm) data2011 <- cbind(under5.mort[, 22], maesles.imm[, 22]) plot(binningDepth2D(data2011, nbins = 8, k = 0.5, remove_borders = TRUE))
Adds plots
combineDepthCurves(x, y, .list = NULL) ## S4 method for signature 'ANY,ANY,list' combineDepthCurves(x, y, .list = NULL) ## S4 method for signature 'DepthCurveList,DepthCurve,ANY' combineDepthCurves(x, y, .list = NULL) ## S4 method for signature 'DepthCurve,DepthCurveList,ANY' combineDepthCurves(x, y, .list = NULL) ## S4 method for signature 'DepthCurve,DepthCurve,ANY' combineDepthCurves(x, y, .list = NULL)
combineDepthCurves(x, y, .list = NULL) ## S4 method for signature 'ANY,ANY,list' combineDepthCurves(x, y, .list = NULL) ## S4 method for signature 'DepthCurveList,DepthCurve,ANY' combineDepthCurves(x, y, .list = NULL) ## S4 method for signature 'DepthCurve,DepthCurveList,ANY' combineDepthCurves(x, y, .list = NULL) ## S4 method for signature 'DepthCurve,DepthCurve,ANY' combineDepthCurves(x, y, .list = NULL)
x |
object |
y |
object |
.list |
list of plots to combine. |
See DepthCurve-class
for description.
Weighted by depth (outlyingness) multivariate location and scatter estimators.
CovLP(x, pdim = 2, la = 1, lb = 1)
CovLP(x, pdim = 2, la = 1, lb = 1)
x |
The data as a matrix or data frame. If it is a matrix or data frame, then each row is viewed as one multivariate observation. |
pdim |
The parameter of the weighted |
la |
parameter of a simple weight function |
lb |
parameter of a simple weight function |
Using depth function one can define a depth-weighted location and scatter estimators. In case of location estimator we have
Subsequently, a depth-weighted scatter estimator is defined as
where is a suitable weight function that can be different from
.
The DepthProc package offers these estimators for weighted depth. Note that
and
include multivariate versions of trimmed means and covariance matrices. Their sample counterparts take the form
where are sample depth weights,
.
loc: Robust Estimate of Location:
cov: Robust Estimate of Covariance:
Returns depth weighted covariance matrix.
Daniel Kosiorowski and Zygmunt Zawadzki from Cracow University of Economics.
depthContour
and depthPersp
for depth graphics.
# EXAMPLE 1 x <- mvrnorm(n = 100, mu = c(0, 0), Sigma = 3 * diag(2)) cov_x <- CovLP(x, 2, 1, 1) # EXAMPLE 2 data(under5.mort, inf.mort, maesles.imm) data1990 <- na.omit(cbind(under5.mort[, 1], inf.mort[, 1], maesles.imm[, 1])) CovLP(data1990)
# EXAMPLE 1 x <- mvrnorm(n = 100, mu = c(0, 0), Sigma = 3 * diag(2)) cov_x <- CovLP(x, 2, 1, 1) # EXAMPLE 2 data(under5.mort, inf.mort, maesles.imm) data1990 <- na.omit(cbind(under5.mort[, 1], inf.mort[, 1], maesles.imm[, 1])) CovLP(data1990)
Air pollution with PM10 in Cracow within day and night in December 2016
data("cracow.airpollution")
data("cracow.airpollution")
data frame containing 744 rows.
1. Kosiorowski D, Rydlewski J P, Zawadzki Z (2017). 'Functional Outliers Detection By The Example Of Air Quality Monitoring.' submitted.
2. Kosiorowski D, Szlachtowska E (2017). 'K- local Median Algorithm for Functional Data in Empirical Analysis of Air Pollution Data.' Proceedings from the 11th Professor A. Zelias International Conference, pp. 153-162. Cracow University of Economics.
Produces a normal DD plot of a multivariate dataset.
ddMvnorm( x, size = nrow(x), robust = FALSE, alpha = 0.05, title = "ddMvnorm", depth_params = list() )
ddMvnorm( x, size = nrow(x), robust = FALSE, alpha = 0.05, title = "ddMvnorm", depth_params = list() )
x |
The data sample for DD plot. |
size |
size of theoretical set |
robust |
Logical. Default |
alpha |
cutoff point for robust measure of covariance. |
title |
title of a plot. |
depth_params |
list of parameters for function depth (method, threads, ndir, la, lb, pdim, mean, cov, exact). |
In the first step the location and scale of x are estimated and theoretical sample from normal distribution with those parameters is generated. The plot presents the depth of empirical points with respect to dataset x and with respect to the theoretical sample.
Returns the normal depth versus depth plot of multivariate dataset x
.
Daniel Kosiorowski, Mateusz Bocian, Anna Wegrzynkiewicz and Zygmunt Zawadzki from Cracow University of Economics.
Liu, R.Y., Parelius, J.M. and Singh, K. (1999), Multivariate analysis by data depth: Descriptive statistics, graphics and inference (with discussion), Ann. Statist., 27, 783–858.
Liu, R.Y., Singh K. (1993), A Quality Index Based on Data Depth and Multivariate Rank Test, Journal of the American Statistical Association vol. 88.
ddPlot
to generate ddPlot to compare to datasets or to compare a dataset with other distributions.
# EXAMPLE 1 norm <- mvrnorm(1000, c(0, 0, 0), diag(3)) con <- mvrnorm(100, c(1, 2, 5), 3 * diag(3)) sample <- rbind(norm, con) ddMvnorm(sample, robust = TRUE) # EXAMPLE 2 data(under5.mort, inf.mort, maesles.imm) data1990 <- na.omit(cbind(under5.mort[, 1], inf.mort[, 1], maesles.imm[, 1])) ddMvnorm(data1990, robust = FALSE)
# EXAMPLE 1 norm <- mvrnorm(1000, c(0, 0, 0), diag(3)) con <- mvrnorm(100, c(1, 2, 5), 3 * diag(3)) sample <- rbind(norm, con) ddMvnorm(sample, robust = TRUE) # EXAMPLE 2 data(under5.mort, inf.mort, maesles.imm) data1990 <- na.omit(cbind(under5.mort[, 1], inf.mort[, 1], maesles.imm[, 1])) ddMvnorm(data1990, robust = FALSE)
Produces a DD plot which allows to compare two multivariate datasets or to compare a subject dataset with theoretical distribution.
ddPlot( x, y, scale = FALSE, location = FALSE, name = "X", name_y = "Y", title = "Depth vs. depth plot", depth_params = list() )
ddPlot( x, y, scale = FALSE, location = FALSE, name = "X", name_y = "Y", title = "Depth vs. depth plot", depth_params = list() )
x |
The first or only data sample for ddPlot. |
y |
The second data sample. |
scale |
logical. determines whether the dispersion is to be aligned. |
location |
determines whether the location is to be aligned to 0 vector with depth median. |
name |
name for data set x. It will be passed to drawing function. |
name_y |
as above for y |
title |
title of the plot. |
depth_params |
list of parameters for function depth (method, threads, ndir, la, lb, pdim, mean, cov, exact). |
For two probability distributions and
, both in
, we can define
depth vs. depth
plot being very useful generalization of the one dimensional quantile-quantile plot:
Its sample counterpart calculated for two samples from
, and
from
is defined as
Daniel Kosiorowski, Mateusz Bocian, Anna Wegrzynkiewicz and Zygmunt Zawadzki from Cracow University of Economics.
Liu, R.Y., Parelius, J.M. and Singh, K. (1999), Multivariate analysis by data depth: Descriptive statistics, graphics and inference (with discussion), Ann. Statist., 27, 822–831.
Liu, R.Y., Singh K. (1993), A Quality Index Based on Data Depth and Multivariate Rank Test, Journal of the American Statistical Association vol. 88.
library(sn) library(mvtnorm) # EXAMPLE 1: Location difference standard <- mvrnorm(1000, c(0, 0), diag(2)) shift <- mvrnorm(1000, c(0.5, 0), diag(2)) ddPlot(x = standard, y = shift, title = "Difference in position") ddPlot(x = standard, y = shift, location = TRUE, title = "Location aligned") # EXAMPLE 2: Scale difference standard <- mvrnorm(1000, c(0, 0), diag(2)) scale <- mvrnorm(1000, c(0, 0), 4 * diag(2)) ddPlot(x = standard, y = scale) ddPlot(x = standard, y = scale, scale = TRUE)
library(sn) library(mvtnorm) # EXAMPLE 1: Location difference standard <- mvrnorm(1000, c(0, 0), diag(2)) shift <- mvrnorm(1000, c(0.5, 0), diag(2)) ddPlot(x = standard, y = shift, title = "Difference in position") ddPlot(x = standard, y = shift, location = TRUE, title = "Location aligned") # EXAMPLE 2: Scale difference standard <- mvrnorm(1000, c(0, 0), diag(2)) scale <- mvrnorm(1000, c(0, 0), 4 * diag(2)) ddPlot(x = standard, y = scale) ddPlot(x = standard, y = scale, scale = TRUE)
Class fro DDPlot
X
Object of class Depth-class.
Y
Object of class Depth-class.
title
title of a plot.
This function calculates deepest regression estimator for simple regression.
deepReg2d(x, y)
deepReg2d(x, y)
x |
Independent variable. |
y |
Dependent variable. |
Function originates from an original algorithm proposed by Rousseeuw and Hubert. Let denotes a sample considered from a following semiparametric model:
, we calculate a depth of a fit
as
, where
denotes the regression residual,
,
.
The deepest regression estimator
is defined as
Daniel Kosiorowski, Mateusz Bocian, Anna Wegrzynkiewicz and Zygmunt Zawadzki from Cracow University of Economics.
Rousseeuw J.P., Hubert M. (1998), Regression Depth, Journal of The American Statistical Association, vol.94.
# EXAMPLE 1 data(pension) plot(pension) abline( lm(Reserves ~ Income, data = pension), lty = 3, lwd = 2) # lm abline( deepReg2d(pension[, 1], pension[, 2]), lwd = 2) # deepreg2d # EXAMPLE 2 data(under5.mort) data(inf.mort) data(maesles.imm) data2011 <- na.omit( cbind(under5.mort[, 22], inf.mort[, 22], maesles.imm[, 22])) x <- data2011[, 3] y <- data2011[, 2] plot( x, y, cex = 1.2, ylab = "infant mortality rate per 1000 live birth", xlab = "against masles immunized percentage", main = "Projection Depth Trimmed vs. LS regressions" ) abline(lm(x ~ y), lwd = 2, col = "black") # lm abline( deepReg2d (x, y), lwd = 2, col = "red" ) # trimmed reg legend( "bottomleft", c("LS", "DeepReg"), fill = c("black", "red"), cex = 1.4, bty = "n" )
# EXAMPLE 1 data(pension) plot(pension) abline( lm(Reserves ~ Income, data = pension), lty = 3, lwd = 2) # lm abline( deepReg2d(pension[, 1], pension[, 2]), lwd = 2) # deepreg2d # EXAMPLE 2 data(under5.mort) data(inf.mort) data(maesles.imm) data2011 <- na.omit( cbind(under5.mort[, 22], inf.mort[, 22], maesles.imm[, 22])) x <- data2011[, 3] y <- data2011[, 2] plot( x, y, cex = 1.2, ylab = "infant mortality rate per 1000 live birth", xlab = "against masles immunized percentage", main = "Projection Depth Trimmed vs. LS regressions" ) abline(lm(x ~ y), lwd = 2, col = "black") # lm abline( deepReg2d (x, y), lwd = 2, col = "red" ) # trimmed reg legend( "bottomleft", c("LS", "DeepReg"), fill = c("black", "red"), cex = 1.4, bty = "n" )
Class for robust regression methods from depthproc package
coef
coefficients of fitted model
depth
regression depth of the fitted values
Calculate depth functions.
depth(u, X, method = "Projection", threads = -1, ...)
depth(u, X, method = "Projection", threads = -1, ...)
u |
Numerical vector or matrix whose depth is to be calculated. Dimension has to be the same as that of the observations. |
X |
The data as a matrix, data frame or list. If it is a matrix or data frame, then each row is viewed as one multivariate observation. If it is a list, all components must be numerical vectors of equal length (coordinates of observations). |
method |
Character string which determines the depth function. |
threads |
number of threads used in parallel computations. Default value -1 means that all possible cores will be used. |
... |
parameters specific to method — see |
The Mahalanobis depth
where denotes the sample covariance matrix
.
A symmetric projection depth of a point
,
is defined as
where Med denotes the univariate median, =
. Its sample version denoted by
or
is obtained by replacing
by its empirical counterpart
calculated from the sample
.
Next interesting depth is the weighted depth. The weighted
depth
of a point
,
generated by
dimensional random vector
with distribution
, is defined as
where
is a suitable weight function on
, and
stands for the
norm (when p = 2 we have usual Euclidean norm). We assume that
is non-decreasing and continuous on
with
, and for
satisfying
. Examples of the weight functions are:
,
or
. The empirical version of the weighted
depth is obtained by replacing distribution
of
in
by its empirical counterpart calculated from the sample
...
The Projection and Tukey's depths are calculated using an approximate algorithm. Calculations of Mahalanobis, Euclidean and depths are exact. Returns the depth of multivariate point u with respect to data set X.
Daniel Kosiorowski, Mateusz Bocian, Anna Wegrzynkiewicz and Zygmunt Zawadzki from Cracow University of Economics.
Liu, R.Y., Parelius, J.M. and Singh, K. (1999), Multivariate analysis by data depth: Descriptive statistics, graphics and inference (with discussion), Ann. Statist., 27, 783–858.
Mosler K (2013). Depth statistics. In C Becker, R Fried, K S (eds.), Robustness and Complex Data Structures, Festschrift in Honour of Ursula Gather, pp. 17–34. Springer.
Rousseeuw, P.J. and Struyf, A. (1998), Computing location depth and regression depth in higher dimensions, Stat. Comput., 8, 193–203.
Zuo, Y. and Serfling, R. (2000), General Notions of Statistical Depth Functions, Ann. Statist., 28, no. 2, 461–482.
depthContour
and depthPersp
for depth graphics.
library(robustbase) # Calculation of Projection depth data(starsCYG, package = "robustbase") depth(t(colMeans(starsCYG)), starsCYG) # Also for matrices depth(starsCYG, starsCYG) # Projection depth applied to a large bivariate data set x <- matrix(rnorm(9999), nc = 3) depth(x, x)
library(robustbase) # Calculation of Projection depth data(starsCYG, package = "robustbase") depth(t(colMeans(starsCYG)), starsCYG) # Also for matrices depth(starsCYG, starsCYG) # Projection depth applied to a large bivariate data set x <- matrix(rnorm(9999), nc = 3) depth(x, x)
Virtual class with structure for every depth class from depthproc package.
u
data set.
X
reference set.
method
depth type.
Draws an approximate contours of depth for bivariate data.
depthContour( x, xlim = extendrange(x[, 1], f = 0.1), ylim = extendrange(x[, 2], f = 0.1), n = 50, pmean = TRUE, mcol = "blue", pdmedian = TRUE, mecol = "brown", legend = TRUE, points = FALSE, colors = heat_hcl, levels = 10, depth_params = list(), graph_params = list(), contour_method = c("auto", "convexhull", "contour") )
depthContour( x, xlim = extendrange(x[, 1], f = 0.1), ylim = extendrange(x[, 2], f = 0.1), n = 50, pmean = TRUE, mcol = "blue", pdmedian = TRUE, mecol = "brown", legend = TRUE, points = FALSE, colors = heat_hcl, levels = 10, depth_params = list(), graph_params = list(), contour_method = c("auto", "convexhull", "contour") )
x |
Bivariate data |
xlim |
Determines the width of x-axis. |
ylim |
Determines the width of y-axis. |
n |
Number of points in each coordinate direction to be used in contour plot. |
pmean |
Logical. If TRUE mean will be marked. |
mcol |
Determines the color of lines describing the mean. |
pdmedian |
Logical. If TRUE depth median will be marked. |
mecol |
Determines the color of lines describing the depth median. |
legend |
Logical. If TRUE legend for mean and depth median will be drawn. |
points |
Logical. If TRUE points from matrix x will be drawn. |
colors |
function for colors pallete (e.g. gray.colors). |
levels |
number of levels for color scale. |
depth_params |
list of parameters for function depth (method, threads, ndir, la, lb, pdim, mean, cov, exact). |
graph_params |
list of graphical parameters for functions filled.contour and contour (e.g. lwd, lty, main). |
contour_method |
determines the method used to draw the contour lines. The default value ("auto") tries
to determine the best method for given depth function.
"convexhull" uses a convex hull algorithm to determine boundaries.
"contour" uses the algorithm from |
The set of all points that have depth at least is called
-trimmed region. The
-trimmed region w.r.t.
is denoted by
, i.e.,
Daniel Kosiorowski, Mateusz Bocian, Anna Wegrzynkiewicz and Zygmunt Zawadzki from Cracow University of Economics.
# EXAMPLE 1 set.seed(123) x <- mvrnorm(1000, c(0, 0), diag(2)) depthContour(x, colors = gray.colors) # with points depthContour(x, points = TRUE) depthContour(x, points = FALSE, levels = 10) # EXAMPLE 2 data(inf.mort, maesles.imm) data1990 <- na.omit(cbind(inf.mort[, 1], maesles.imm[, 1])) depthContour(data1990, n = 50, pmean = TRUE, mcol = "blue", pdmedian = TRUE, mecol = "brown", legend = TRUE, points = TRUE, depth_params = list(method = "LP"), graph_params = list( xlab = "infant mortality rate per 1000 live birth", ylab = "against masles immunized percentage", main = "L2 depth, UN Fourth Goal 2011 year")) #EXAMPLE 3 data("france") depthContour(france, depth_params = list(method = "Tukey"), points = TRUE )
# EXAMPLE 1 set.seed(123) x <- mvrnorm(1000, c(0, 0), diag(2)) depthContour(x, colors = gray.colors) # with points depthContour(x, points = TRUE) depthContour(x, points = FALSE, levels = 10) # EXAMPLE 2 data(inf.mort, maesles.imm) data1990 <- na.omit(cbind(inf.mort[, 1], maesles.imm[, 1])) depthContour(data1990, n = 50, pmean = TRUE, mcol = "blue", pdmedian = TRUE, mecol = "brown", legend = TRUE, points = TRUE, depth_params = list(method = "LP"), graph_params = list( xlab = "infant mortality rate per 1000 live birth", ylab = "against masles immunized percentage", main = "L2 depth, UN Fourth Goal 2011 year")) #EXAMPLE 3 data("france") depthContour(france, depth_params = list(method = "Tukey"), points = TRUE )
This page describes mechanism behavior of ScaleCurve and AsymmetryCurve
DepthCurve is a virtual class that contains methods (getPlot(...) and plot(...)) for rendering single curve such as ScaleCurve or AsymmetryCurve. Such object can be combined by overloaded operator '
depth
object of Depth-class
name
name of dataset used on plot
title
title of a plot
alpha
central area values
library(mvtnorm) x <- mvrnorm(n = 100, mu = c(0, 0), Sigma = 2 * diag(2)) y <- rmvt(n = 100, sigma = diag(2), df = 4) s1 <- scaleCurve(x, depth_params = list(method = "Projection")) s2 <- scaleCurve(y, depth_params = list(method = "Projection"), name = "Set2") sc_list <- combineDepthCurves(s1, s2) # Add one curve to another plot(sc_list) # Draw plot with two curves z <- mvrnorm(n = 100, mu = c(0, 0), Sigma = 1 * diag(2)) s3 <- scaleCurve(z, depth_params = list(method = "Projection")) plot(combineDepthCurves(sc_list, s3)) # Add third curve and draw a plot
library(mvtnorm) x <- mvrnorm(n = 100, mu = c(0, 0), Sigma = 2 * diag(2)) y <- rmvt(n = 100, sigma = diag(2), df = 4) s1 <- scaleCurve(x, depth_params = list(method = "Projection")) s2 <- scaleCurve(y, depth_params = list(method = "Projection"), name = "Set2") sc_list <- combineDepthCurves(s1, s2) # Add one curve to another plot(sc_list) # Draw plot with two curves z <- mvrnorm(n = 100, mu = c(0, 0), Sigma = 1 * diag(2)) s3 <- scaleCurve(z, depth_params = list(method = "Projection")) plot(combineDepthCurves(sc_list, s3)) # Add third curve and draw a plot
DepthCurveList is a special container for DepthCurve objects. See DepthCurve-class
Experimental function used to fit depth weighted density estimator.
depthDensity(x, y, nx = 5, ny = 32, xg = NULL, yg = NULL, ...)
depthDensity(x, y, nx = 5, ny = 32, xg = NULL, yg = NULL, ...)
x |
numeric vector |
y |
numeric vector |
nx |
the number of equally spaced points at which the density is to be estimated in x-dimension. |
ny |
the number of equally spaced points at which the density is to be estimated in x-dimension. |
xg |
vector of point at which the density is to be estimated. |
yg |
vector of point at which the density is to be estimated. |
... |
arguments passed to depthLocal. |
Kosiorowski D. and Zawadzki Z. (2014) Notes on optimality of predictive distribution pseudo-estimators in the CHARME models and automatic trading strategies, FindEcon2014, submitted
## Not run: # .sampleData is special function for creating # data for testing conditional denisty estimators data <- DepthProc:::.sampleData(1:5, 100) x <- data[, 1] y <- data[, 2] plot(x, y) dep <- depthDensity(x, y) plot(dep, type = "raw") plot(dep, type = "depth") ## End(Not run)
## Not run: # .sampleData is special function for creating # data for testing conditional denisty estimators data <- DepthProc:::.sampleData(1:5, 100) x <- data[, 1] y <- data[, 2] plot(x, y) dep <- depthDensity(x, y) plot(dep, type = "raw") plot(dep, type = "depth") ## End(Not run)
Computes the euclidean depth of a point or vectors of points with respect to a multivariate data set.
depthEuclid(u, X)
depthEuclid(u, X)
u |
Numerical vector or matrix whose depth is to be calculated. Dimension has to be the same as that of the observations. |
X |
The data as a matrix, data frame or list. If it is a matrix or data frame, then each row is viewed as one multivariate observation. If it is a list, all components must be numerical vectors of equal length (coordinates of observations). |
Calculation of Euclidean depth is exact.
Returns the depth of multivariate point u
with respect to data set X
.
Daniel Kosiorowski, Mateusz Bocian, Anna Wegrzynkiewicz and Zygmunt Zawadzki from Cracow University of Economics.
x <- matrix(rnorm(9999), nc = 3) depthEuclid(x, x)
x <- matrix(rnorm(9999), nc = 3) depthEuclid(x, x)
Computes local version of depth according to proposals of Paindaveine and Van Bever — see referencess.
depthLocal( u, X, beta = 0.5, depth_params1 = list(method = "Projection"), depth_params2 = depth_params1 )
depthLocal( u, X, beta = 0.5, depth_params1 = list(method = "Projection"), depth_params2 = depth_params1 )
u |
Numerical vector or matrix whose depth is to be calculated. Dimension has to be the same as that of the observations. |
X |
The data as a matrix, data frame. If it is a matrix or data frame, then each row is viewed as one multivariate observation. |
beta |
cutoff value for neighbourhood |
depth_params1 |
list of parameters for function depth (method, threads, ndir, la, lb, pdim, mean, cov, exact). |
depth_params2 |
as above — default is depth_params1. |
A successful concept of local depth was proposed by Paindaveine and Van Bever (2012). For defining a neighbourhood of a point authors proposed using idea of symmetrisation of a distribution (a sample) with respect to a point in which depth is calculated. In their approach instead of a distribution , a distribution
is used. For any
, let us introduce the smallest depth region bigger or equal to
,
where . Then for a locality parameter
we can take a neighbourhood of a point
as
.
Formally, let be a depth function. Then the local depth with the locality parameter
and w.r.t. a point
is defined as
where is cond. distr. of
conditioned on
.
Paindaveine, D., Van Bever, G. (2013) From depth to local depth : a focus on centrality. Journal of the American Statistical Association 105, 1105–1119.
## Not run: # EXAMPLE 1 data <- mvrnorm(100, c(0, 5), diag(2) * 5) # By default depth_params2 = depth_params1 depthLocal(data, data, depth_params1 = list(method = "LP")) depthLocal(data, data, depth_params1 = list(method = "LP"), depth_params2 = list(method = "Projection")) # Depth contour depthContour(data, depth_params = list(method = "Local", depth_params1 = list(method = "LP"))) # EXAMPLE 2 data(inf.mort, maesles.imm) data1990 <- na.omit(cbind(inf.mort[, 1], maesles.imm[, 1])) depthContour(data1990, depth_params = list( method = "Local", depth_params1 = list(method = "LP"), beta = 0.3 )) # EXAMPLE 3 Sigma1 <- matrix(c(10, 3, 3, 2), 2, 2) X1 <- mvrnorm(n = 8500, mu = c(0, 0), Sigma1) Sigma2 <- matrix(c(10, 0, 0, 2), 2, 2) X2 <- mvrnorm(n = 1500, mu = c(-10, 6), Sigma2) BALLOT <- rbind(X1, X2) train <- sample(1:10000, 100) data <- BALLOT[train, ] depthContour(data, depth_params = list( method = "Local", beta = 0.3, depth_params1 = list(method = "Projection") )) ## End(Not run)
## Not run: # EXAMPLE 1 data <- mvrnorm(100, c(0, 5), diag(2) * 5) # By default depth_params2 = depth_params1 depthLocal(data, data, depth_params1 = list(method = "LP")) depthLocal(data, data, depth_params1 = list(method = "LP"), depth_params2 = list(method = "Projection")) # Depth contour depthContour(data, depth_params = list(method = "Local", depth_params1 = list(method = "LP"))) # EXAMPLE 2 data(inf.mort, maesles.imm) data1990 <- na.omit(cbind(inf.mort[, 1], maesles.imm[, 1])) depthContour(data1990, depth_params = list( method = "Local", depth_params1 = list(method = "LP"), beta = 0.3 )) # EXAMPLE 3 Sigma1 <- matrix(c(10, 3, 3, 2), 2, 2) X1 <- mvrnorm(n = 8500, mu = c(0, 0), Sigma1) Sigma2 <- matrix(c(10, 0, 0, 2), 2, 2) X2 <- mvrnorm(n = 1500, mu = c(-10, 6), Sigma2) BALLOT <- rbind(X1, X2) train <- sample(1:10000, 100) data <- BALLOT[train, ] depthContour(data, depth_params = list( method = "Local", beta = 0.3, depth_params1 = list(method = "Projection") )) ## End(Not run)
Computes the LP depth of a point or vectors of points with respect to a multivariate data set.
depthLP(u, X, pdim = 2, la = 1, lb = 1, threads = -1, func = NULL)
depthLP(u, X, pdim = 2, la = 1, lb = 1, threads = -1, func = NULL)
u |
Numerical vector or matrix whose depth is to be calculated. Dimension has to be the same as that of the observations. |
X |
The data as a matrix, data frame or list. If it is a matrix or data frame, then each row is viewed as one multivariate observation. If it is a list, all components must be numerical vectors of equal length (coordinates of observations). |
pdim |
dimension used in calculating depth function. |
la |
slope the weighing function. |
lb |
intercept in the weighing function. |
threads |
number of threads used in parallel computations. Default value -1 means that all possible cores will be used. |
func |
the weighing function. Currently it is not supported. |
Returns the depth of multivariate point u
with respect to data set X
.
Daniel Kosiorowski, Mateusz Bocian, Anna Wegrzynkiewicz and Zygmunt Zawadzki from Cracow University of Economics.
x <- matrix(rnorm(3000), ncol = 3) # Same results depthLP(x, x, pdim = 2)
x <- matrix(rnorm(3000), ncol = 3) # Same results depthLP(x, x, pdim = 2)
Computes the mahalanobis depth of a point or vectors of points with respect to a multivariate data set.
depthMah(u, X, cov = NULL, mean = NULL, threads = -1)
depthMah(u, X, cov = NULL, mean = NULL, threads = -1)
u |
Numerical vector or matrix whose depth is to be calculated. Dimension has to be the same as that of the observations. |
X |
The data as a matrix, data frame or list. If it is a matrix or data frame, then each row is viewed as one multivariate observation. If it is a list, all components must be numerical vectors of equal length (coordinates of observations). |
cov |
custom covariance matrix passed. If NULL standard calculations will be based on standard covariance estimator. |
mean |
custom mean vector. If null — mean average will be used. |
threads |
number of threads used in parallel computations. Default value -1 means that all possible cores will be used. |
Calculation of Mahalanobis depth is exact.
Returns the depth of multivariate point u
with respect to data set X
.
Daniel Kosiorowski, Mateusz Bocian, Anna Wegrzynkiewicz and Zygmunt Zawadzki from Cracow University of Economics.
x <- matrix(rnorm(9999), nc = 3) depthMah(x, x)
x <- matrix(rnorm(9999), nc = 3) depthMah(x, x)
Return point with maximum depth function value. If multiple points have the same value, mean average of them will be returned.
depthMedian(x, depth_params = list(), convex = FALSE) ## S4 method for signature 'matrix' depthMedian(x, depth_params = list(), convex = FALSE) ## S4 method for signature 'data.frame' depthMedian(x, depth_params = list(), convex = FALSE) ## S4 method for signature 'Depth' depthMedian(x, convex = FALSE)
depthMedian(x, depth_params = list(), convex = FALSE) ## S4 method for signature 'matrix' depthMedian(x, depth_params = list(), convex = FALSE) ## S4 method for signature 'data.frame' depthMedian(x, depth_params = list(), convex = FALSE) ## S4 method for signature 'Depth' depthMedian(x, convex = FALSE)
x |
object of class Depth or matrix. |
depth_params |
list of parameters for function depth (method, threads, ndir, la, lb, pdim, mean, cov, exact). |
convex |
logical. If true, than centroid of the convex hull created from deepest points is returned. |
# depthMedian for matrix x <- matrix(rnorm(600), nc = 3) depthMedian(x) # depthMedian works with object of class Depth dp <- depth(x) depthMedian(dp)
# depthMedian for matrix x <- matrix(rnorm(600), nc = 3) depthMedian(x) # depthMedian works with object of class Depth dp <- depth(x) depthMedian(dp)
Draws a perspective plot of depth function over x-y plane.
depthPersp( x, plot_method = "lattice", xlim = extendrange(x[, 1], f = 0.1), ylim = extendrange(x[, 2], f = 0.1), n = 50, xlab = "x", ylab = "y", plot_title = NULL, colors = heat_hcl, depth_params = list(), graph_params = list() )
depthPersp( x, plot_method = "lattice", xlim = extendrange(x[, 1], f = 0.1), ylim = extendrange(x[, 2], f = 0.1), n = 50, xlab = "x", ylab = "y", plot_title = NULL, colors = heat_hcl, depth_params = list(), graph_params = list() )
x |
bivariate data |
plot_method |
there are two options "lattice", and "rgl" — see details |
xlim |
limits for x-axis |
ylim |
limits for y-axis |
n |
number of points that will be used to create plot ( |
xlab |
description of x-axis |
ylab |
description of y-axis |
plot_title |
plot title (default NULL means paste(depth_params$method, "depth")) |
colors |
function for colors pallete (e.g. gray.colors). |
depth_params |
list of parameters for function depth ("method", "threads", "ndir", "la", "lb", "pdim", "mean", "cov", "exact"). |
graph_params |
list of graphical parameters for functions rgl::persp3d and lattice::wireframe. |
plot_method — rgl package is not in depends list beacuse it may cause problems when OpenGL is not supported. To use plot_method = "rgl" you must load this package on your own.
Daniel Kosiorowski, Mateusz Bocian, Anna Wegrzynkiewicz and Zygmunt Zawadzki from Cracow University of Economics.
# EXAMPLE 1 x <- mvrnorm(100, c(0, 0), diag(2)) depthPersp(x, depth_params = list(method = "Euclidean")) # EXAMPLE 2 data(inf.mort, maesles.imm) data1990 <- na.omit(cbind(inf.mort[, 1], maesles.imm[, 1])) ## Not run: library(rgl) depthPersp(data1990, plot_method = "rgl", depth_params = list(method = "Projection")) ## End(Not run)
# EXAMPLE 1 x <- mvrnorm(100, c(0, 0), diag(2)) depthPersp(x, depth_params = list(method = "Euclidean")) # EXAMPLE 2 data(inf.mort, maesles.imm) data1990 <- na.omit(cbind(inf.mort[, 1], maesles.imm[, 1])) ## Not run: library(rgl) depthPersp(data1990, plot_method = "rgl", depth_params = list(method = "Projection")) ## End(Not run)
Computes the Projection depth of a point or vectors of points with respect to a multivariate data set.
depthProjection(u, X, ndir = 1000, threads = -1)
depthProjection(u, X, ndir = 1000, threads = -1)
u |
Numerical vector or matrix whose depth is to be calculated. Dimension has to be the same as that of the observations. |
X |
The data as a matrix, data frame or list. If it is a matrix or data frame, then each row is viewed as one multivariate observation. If it is a list, all components must be numerical vectors of equal length (coordinates of observations). |
ndir |
number of directions used in computations |
threads |
number of threads used in parallel computations. Default value -1 means that all possible cores will be used. |
Irrespective of dimension, Projection and Tukey's depth is obtained by approximate calculation.
Returns the depth of multivariate point u
with respect to data set X
.
Daniel Kosiorowski, Mateusz Bocian, Anna Wegrzynkiewicz and Zygmunt Zawadzki from Cracow University of Economics.
x <- matrix(rnorm(3000), nc = 3) a <- depthProjection(x, x, ndir = 2000)
x <- matrix(rnorm(3000), nc = 3) a <- depthProjection(x, x, ndir = 2000)
Computes the Tukey depth of a point or vectors of points with respect to a multivariate data set.
depthTukey(u, X, ndir = 1000, threads = -1, exact = FALSE)
depthTukey(u, X, ndir = 1000, threads = -1, exact = FALSE)
u |
Numerical vector or matrix whose depth is to be calculated. Dimension has to be the same as that of the observations. |
X |
The data as a matrix, data frame or list. If it is a matrix or data frame, then each row is viewed as one multivariate observation. If it is a list, all components must be numerical vectors of equal length (coordinates of observations). |
ndir |
number of directions used in computations |
threads |
number of threads used in parallel computations. Default value -1 means that all possible cores will be used. |
exact |
if TRUE exact alhorithm will be used . Currently it works only for 2 dimensional data set. |
Irrespective of dimension, Projection and Tukey's depth is obtained by approximate calculation.
Returns the depth of multivariate point u
with respect to data set X
.
Daniel Kosiorowski, Mateusz Bocian, Anna Wegrzynkiewicz and Zygmunt Zawadzki from Cracow University of Economics.
## Not run: x <- matrix(rnorm(3000), nc = 3) depthTukey(x, ndir = 2000) ## End(Not run) # Exact algorithm in 2d x <- matrix(rnorm(2000), nc = 2) depthTukey(x, exact = TRUE)
## Not run: x <- matrix(rnorm(3000), nc = 3) depthTukey(x, ndir = 2000) ## End(Not run) # Exact algorithm in 2d x <- matrix(rnorm(2000), nc = 2) depthTukey(x, exact = TRUE)
Functional boxplot based on Modified Band Depth
fncBoxPlot(u, X = NULL, bands = c(0, 0.5), method = "MBD", byrow = NULL, ...)
fncBoxPlot(u, X = NULL, bands = c(0, 0.5), method = "MBD", byrow = NULL, ...)
u |
data matrix |
X |
reference set. If null u will be used as reference. |
bands |
limits for bands |
method |
depth method |
byrow |
byrow |
... |
other arguments passed to fncDepth |
# some data: x <- matrix(rnorm(200), ncol = 10) fncBoxPlot(x, bands = c(0, 0.5, 1), method = "FM") fncBoxPlot(x, bands = c(0, 0.5, 1), method = "FM", byrow = FALSE) colnames(x) <- paste0("f", 1:ncol(x)) fncBoxPlot(x, bands = c(0, 0.5, 1), method = "FM") # fncBoxPlot handles zoo and xts objects library(xts) x <- matrix(rnorm(200), ncol = 10) time <- as.POSIXct(1:ncol(x) * 86400, origin = "1970-01-01") x_xts <- xts(t(x), order.by = time) fncBoxPlot(x_xts, bands = c(0, 0.5, 1), method = "FM") data("katowice.airpollution") pl <- fncBoxPlot(katowice.airpollution, bands = c(0, 0.5, 1), method = "MBD") pl + ggtitle("Air pollution in Katowice") + labs(y= "pollination ", x = "hour ")
# some data: x <- matrix(rnorm(200), ncol = 10) fncBoxPlot(x, bands = c(0, 0.5, 1), method = "FM") fncBoxPlot(x, bands = c(0, 0.5, 1), method = "FM", byrow = FALSE) colnames(x) <- paste0("f", 1:ncol(x)) fncBoxPlot(x, bands = c(0, 0.5, 1), method = "FM") # fncBoxPlot handles zoo and xts objects library(xts) x <- matrix(rnorm(200), ncol = 10) time <- as.POSIXct(1:ncol(x) * 86400, origin = "1970-01-01") x_xts <- xts(t(x), order.by = time) fncBoxPlot(x_xts, bands = c(0, 0.5, 1), method = "FM") data("katowice.airpollution") pl <- fncBoxPlot(katowice.airpollution, bands = c(0, 0.5, 1), method = "MBD") pl + ggtitle("Air pollution in Katowice") + labs(y= "pollination ", x = "hour ")
Calculates depth functions.
fncDepth(u, X = NULL, method = "MBD", byrow = NULL, ...) ## S3 method for class 'matrix' fncDepth(u, X = NULL, method = "MBD", byrow = NULL, ...) ## S3 method for class 'zoo' fncDepth(u, X = NULL, method = "MBD", byrow = NULL, ...)
fncDepth(u, X = NULL, method = "MBD", byrow = NULL, ...) ## S3 method for class 'matrix' fncDepth(u, X = NULL, method = "MBD", byrow = NULL, ...) ## S3 method for class 'zoo' fncDepth(u, X = NULL, method = "MBD", byrow = NULL, ...)
u |
data |
X |
reference set. If null u will be used as reference. |
method |
depth method - "MBD" (default), or "FM" (Frainman-Muniz depth) |
byrow |
logical or character. |
... |
additional arguments passed to fncDepthFM. |
x <- matrix(rnorm(60), ncol = 20) fncDepth(x, method = "FM", dep1d = "Mahalanobis") fncDepth(x, byrow = FALSE) # zoo and xts library(xts) data(sample_matrix) sample.xts <- as.xts(sample_matrix, descr = "my new xts object") fncDepth(sample.xts)
x <- matrix(rnorm(60), ncol = 20) fncDepth(x, method = "FM", dep1d = "Mahalanobis") fncDepth(x, byrow = FALSE) # zoo and xts library(xts) data(sample_matrix) sample.xts <- as.xts(sample_matrix, descr = "my new xts object") fncDepth(sample.xts)
Computes Frainman-Muniz depth for functional data.
fncDepthFM(u, X, dep1d_params = list(method = "Projection"))
fncDepthFM(u, X, dep1d_params = list(method = "Projection"))
u |
Numerical vector or matrix whose depth is to be calculated. Dimension has to be the same as that of the observations. |
X |
The data as a matrix. If it is a matrix or data frame, then each row is viewed as one multivariate observation. |
dep1d_params |
parameters passed to depth function used in one dimension. |
x <- matrix(rnorm(60), nc = 20) fncDepthFM(x)
x <- matrix(rnorm(60), nc = 20) fncDepthFM(x)
Computes the modified band depth.
fncDepthMBD(u, X)
fncDepthMBD(u, X)
u |
Numerical vector or matrix whose depth is to be calculated. Dimension has to be the same as that of the observations. |
X |
The data as a matrix. If it is a matrix or data frame, then each row is viewed as one multivariate observation. |
x <- matrix(rnorm(60), nc = 20) fncDepthMBD(x) fncDepthMBD(x, x)
x <- matrix(rnorm(60), nc = 20) fncDepthMBD(x) fncDepthMBD(x, x)
Calculate functional median based on data depth.
fncDepthMedian(u, X = NULL, method = "MBD", byrow = NULL, unique = TRUE, ...)
fncDepthMedian(u, X = NULL, method = "MBD", byrow = NULL, unique = TRUE, ...)
u |
data matrix |
X |
reference set. If null u will be used as reference. |
method |
depth method |
byrow |
byrow |
unique |
if true |
... |
other arguments passed to fncDepth |
x <- matrix(rnorm(600), nc = 20) md <- fncDepthMedian(x, method = "FM", dep1d = "Mahalanobis")
x <- matrix(rnorm(600), nc = 20) md <- fncDepthMedian(x, method = "FM", dep1d = "Mahalanobis")
Extract bands from functional depth object.
fncGetBand(obj, band = 0.5)
fncGetBand(obj, band = 0.5)
obj |
object that inherits from FunctionalDepth. |
band |
single numeric value. |
x <- matrix(rnorm(600), nc = 20) obj <- fncDepth(x, method = "FM", dep1d = "Mahalanobis") fncGetBand(obj)
x <- matrix(rnorm(600), nc = 20) obj <- fncDepth(x, method = "FM", dep1d = "Mahalanobis") fncGetBand(obj)
Relation between minimum wage (MW) and unemployment rate (UR) in France.
data(france)
data(france)
data frame containing 17 rows and two column. MW is a minimum wage, and UR is an unemployment rate.
Virtual class with structure for every functional depth class from depthproc package. Inherits from Depth-class
.
index
numeric, or time-based object.
Create an object of class ggplot from DepthCurve and DepthCurveList.
getPlot(object) ## S4 method for signature 'AsymmetryCurveList' getPlot(object) ## S4 method for signature 'DDPlot' getPlot(object) ## S4 method for signature 'ScaleCurveList' getPlot(object)
getPlot(object) ## S4 method for signature 'AsymmetryCurveList' getPlot(object) ## S4 method for signature 'DDPlot' getPlot(object) ## S4 method for signature 'ScaleCurveList' getPlot(object)
object |
a DDPlot ScaleCurve or AsymmetryCurve object class. |
Infant mortality rate (0–1 year) per 1,000 live births
data(inf.mort)
data(inf.mort)
A data frame with 654 rows and 4 variables
http://mdgs.un.org/unsd/mdg/Data.aspx
Internet view data
data(internet.users)
data(internet.users)
data frame containing 17518 rows and 6 columns — 17518 working days of the Internet service considered with respect to variables: service, month, day, hour, unique users and page views.
Kosiorowski, Rydlewski, Snarska (2016), Detecting a Structural Change in Functional Time Series Using Local Wilcoxon Statistic, arXiv: 1604.03776v2
Air pollution in Katowice city by hour.
data("katowice.airpollution")
data("katowice.airpollution")
data frame containing 181 rows (days) and 24 columns. Each column is an air pollination for given hour.
This function add one location-scale contour to the existing plot.
lsdAddContour(x, cont = NULL, ...) ## S4 method for signature 'LSDepthContour' lsdAddContour(x, cont = NULL, ...)
lsdAddContour(x, cont = NULL, ...) ## S4 method for signature 'LSDepthContour' lsdAddContour(x, cont = NULL, ...)
x |
object of class LSDepthContour |
cont |
depth of contour to plot |
... |
other arguments passed to polygon function |
smp <- rf(100, 5, 10) x <- lsdSampleDepthContours(smp) plot(x) lsdAddContour(x, 0.1, col = "grey50") lsdAddContour(x, 0.3, col = "grey10", border = "red", lwd = 4)
smp <- rf(100, 5, 10) x <- lsdSampleDepthContours(smp) plot(x) lsdAddContour(x, 0.1, col = "grey50") lsdAddContour(x, 0.3, col = "grey10", border = "red", lwd = 4)
Class used to store maximum location-scale depth results.
max_depth
maximum Student depth value.
mu
location estimate in the deepest point.
sigma
scale estimate in the deepest point.
Class used to store result of location-scale depth contours.
cont_depth
depth values used to calculate contours.
sample
original sample used to calculate depth contours.
.Data
list with estimated values of scale-depth contours.
Get numeric values of the location-scale depth contour from existing object of LSDepthContour class.
lsdGetContour(x, cont) ## S4 method for signature 'LSDepthContour' lsdGetContour(x, cont)
lsdGetContour(x, cont) ## S4 method for signature 'LSDepthContour' lsdGetContour(x, cont)
x |
object of class LSDepthContour |
cont |
single numeric — depth of contour to return |
Calculations are based on lsdepth algorithm written by Ch. Muller.
dcont <- lsdSampleDepthContours(rf(200, 4, 7), depth = c(0.1, 0.2)) # get contour that is present in dcont object lsdGetContour(dcont, 0.1) # get contour that is not present in dcont # it will be automatically calculated lsdGetContour(dcont, 0.3)
dcont <- lsdSampleDepthContours(rf(200, 4, 7), depth = c(0.1, 0.2)) # get contour that is present in dcont object lsdGetContour(dcont, 0.1) # get contour that is not present in dcont # it will be automatically calculated lsdGetContour(dcont, 0.3)
Calculate sample one-dimensional Mizera and Muller Student depth contours.
lsdSampleDepthContours(x, depth = c(0.1, 0.2, 0.3, 0.4), lengthmu = 1000)
lsdSampleDepthContours(x, depth = c(0.1, 0.2, 0.3, 0.4), lengthmu = 1000)
x |
one dimensional vector with sample |
depth |
depth level for contours |
lengthmu |
number of points to evalute depth |
Calculations are based on lsdepth algorithm written by Ch. Muller.
Mizera, I., Muller, C. H., 2004. Location-scale depth (with discussion). Journal of the American Statistical Association 99, 949–966.
# EXAMPLE 1 # F-distribution dcont <- lsdSampleDepthContours(rf(200, 4, 7)) plot(dcont) # EXAMPLE 2 # normal distribution - more contours calculated dcont_norm <- lsdSampleDepthContours(rnorm(100), seq(0.05, 0.4, 0.05)) plot(dcont_norm)
# EXAMPLE 1 # F-distribution dcont <- lsdSampleDepthContours(rf(200, 4, 7)) plot(dcont) # EXAMPLE 2 # normal distribution - more contours calculated dcont_norm <- lsdSampleDepthContours(rnorm(100), seq(0.05, 0.4, 0.05)) plot(dcont_norm)
Calculates the maximum Student depth estimator of location and scale for one dimensional data (an alternative for MED and MAD or for the mean and standard deviation).
lsdSampleMaxDepth(x, iter = 100, eps = 1e-04, p_length = 10)
lsdSampleMaxDepth(x, iter = 100, eps = 1e-04, p_length = 10)
x |
one dimensional vector with sample |
iter |
maximum number of iterations in algorith for calculation Location-Scale Depth |
eps |
tolerance level |
p_length |
is the maximum length of the precision step at the end |
Calculations are based on lsdepth algorithm written by Ch. Muller.
Mizera, I., Muller, C. H., 2004. Location-scale depth (with discussion). Journal of the American Statistical Association 99, 949–966.
x <- rnorm(100) lsdSampleMaxDepth(x) y <- rf(100, 4, 10) lsdSampleMaxDepth(y)
x <- rnorm(100) lsdSampleMaxDepth(x) y <- rf(100, 4, 10) lsdSampleMaxDepth(y)
Children 1 year old immunized against measles, percentage
data(maesles.imm)
data(maesles.imm)
A data frame with 654 rows and 4 variables
http://mdgs.un.org/unsd/mdg/Data.aspx
Depth based multivariate Wilcoxon test for a scale difference.
mWilcoxonTest(x, y, alternative = "two.sided", depth_params = list())
mWilcoxonTest(x, y, alternative = "two.sided", depth_params = list())
x |
data matrix |
y |
data matrix |
alternative |
a character string specifying the alternative hypothesis, must be one of "two.sided" (default), "greater" or "less". |
depth_params |
list of parameters for function depth (method, threads, ndir, la, lb, pdim, mean, cov, exact). |
Having two samples and
using any depth function, we can compute depth values in a combined sample
, assuming the empirical distribution calculated basing on all observations, or only on observations belonging to one of the samples
or
.
For example if we observe depths are more likely to cluster tightly around the center of the combined sample, while
depths are more likely to scatter outlying positions, then we conclude
was drawn from a distribution with larger scale.
Properties of the DD plot based statistics in the i.i.d setting were studied by Li & Liu (2004). Authors proposed several DD-plot based statistics and presented bootstrap arguments for their consistency and good effectiveness in comparison to Hotelling and multivariate analogues of Ansari-Bradley and Tukey-Siegel statistics. Asymptotic distributions of depth based multivariate Wilcoxon rank-sum test statistic under the null and general alternative hypotheses were obtained by Zuo & He (2006). Several properties of the depth based rang test involving its unbiasedness was critically discussed by Jureckova & Kalina (2012). Basing on DD-plot object, which is available within the DepthProc it is possible to define several multivariate generalizations of one-dimensional rank and order statistics in an easy way. These generalizations cover well known Wilcoxon rang-sum statistic.
The depth based multivariate Wilcoxon rang sum test is especially useful for the multivariate scale changes detection and was introduced among other by Liu & Singh (2003) and intensively studied by Jureckowa & Kalina (2012) and Zuo & He (2006) in the i.i.d. setting.
For the samples ,
, their
,
, depths w.r.t. a combined sample
the Wilcoxon statistic is defined as
, where
denotes the rang of the i-th observation,
in the combined sample
.
The distribution of is symmetric about
, its variance is
.
Jureckova J, Kalina J (2012). Nonparametric multivariate rank tests and their unbiasedness. Bernoulli, 18(1), 229–251. Li J, Liu RY (2004). New nonparametric tests of multivariate locations and scales using data depth. Statistical Science, 19(4), 686–696. Liu RY, Singh K (1995). A quality index based on data depth and multivariate rank tests. Journal of American Statistical Association, 88, 252–260. Zuo Y, He X (2006). On the limiting distributions of multivariate depth-based rank sum statistics and related tests. The Annals of Statistics, 34, 2879–2896.
# EXAMPLE 1 x <- mvrnorm(100, c(0, 0), diag(2)) y <- mvrnorm(100, c(0, 0), diag(2) * 1.4) mWilcoxonTest(x, y) mWilcoxonTest(x, y, depth_params = list(method = "LP")) # EXAMPLE 2 data(under5.mort) data(inf.mort) data(maesles.imm) data2011 <- na.omit(cbind(under5.mort[, 22], inf.mort[, 22], maesles.imm[, 22])) data1990 <- na.omit(cbind(under5.mort[, 1], inf.mort[, 1], maesles.imm[, 1])) mWilcoxonTest(data2011, data1990)
# EXAMPLE 1 x <- mvrnorm(100, c(0, 0), diag(2)) y <- mvrnorm(100, c(0, 0), diag(2) * 1.4) mWilcoxonTest(x, y) mWilcoxonTest(x, y, depth_params = list(method = "LP")) # EXAMPLE 2 data(under5.mort) data(inf.mort) data(maesles.imm) data2011 <- na.omit(cbind(under5.mort[, 22], inf.mort[, 22], maesles.imm[, 22])) data1990 <- na.omit(cbind(under5.mort[, 1], inf.mort[, 1], maesles.imm[, 1])) mWilcoxonTest(data2011, data1990)
Plot Depth curve
plot(x, y, ...) ## S4 method for signature 'DDPlot,ANY' plot(x) ## S4 method for signature 'DepthCurve,ANY' plot(x) ## S4 method for signature 'DepthCurveList,ANY' plot(x)
plot(x, y, ...) ## S4 method for signature 'DDPlot,ANY' plot(x) ## S4 method for signature 'DepthCurve,ANY' plot(x) ## S4 method for signature 'DepthCurveList,ANY' plot(x)
x |
object that inherits from DepthCurve class (ScaleCurve or AsymmetryCurve), or DDPlot class. |
y |
not supported. |
... |
not supported. |
x <- mvrnorm(n = 100, mu = c(0, 0), Sigma = 3 * diag(2)) sc <- scaleCurve(x) plot(sc)
x <- mvrnorm(n = 100, mu = c(0, 0), Sigma = 3 * diag(2)) sc <- scaleCurve(x) plot(sc)
Binning 2d
## S4 method for signature 'BinnDepth2d,ANY' plot(x, ..., alpha = 0.1, bg_col = "red", add_mid = TRUE)
## S4 method for signature 'BinnDepth2d,ANY' plot(x, ..., alpha = 0.1, bg_col = "red", add_mid = TRUE)
x |
object of class BinnDepth2d |
... |
graphical parameters passed to plot |
alpha |
alpha value for rgb function |
bg_col |
backgroud color |
add_mid |
logical. If TRUE centers of binns will be marked. |
tmp <- binningDepth2D(x = mvrnorm(100, rep(0, 2), diag(2))) plot(tmp)
tmp <- binningDepth2D(x = mvrnorm(100, rep(0, 2), diag(2))) plot(tmp)
Create plot for DepthDensity. See depthDensity
for more information.
## S4 method for signature 'DepthDensity,ANY' plot(x, type = "depth", ...)
## S4 method for signature 'DepthDensity,ANY' plot(x, type = "depth", ...)
x |
object of class DepthDensity |
type |
type of density that will be plotted. "depth" is a depth scaled density, and "raw" is denisty without scaling. |
... |
graphical arguments. |
Create location-scale depth plot. See lsdSampleDepthContours
for more information.
## S4 method for signature 'LSDepthContour,ANY' plot( x, cont = NULL, ratio = 1, mu_min = NULL, mu_max = NULL, col = NULL, border = NULL, ... )
## S4 method for signature 'LSDepthContour,ANY' plot( x, cont = NULL, ratio = 1, mu_min = NULL, mu_max = NULL, col = NULL, border = NULL, ... )
x |
object of class LSDepthContour |
cont |
plotted contours. Default NULL means that all contours stored in x will be plotted. |
ratio |
ratio |
mu_min |
mu_min |
mu_max |
mu_max |
col |
vectors with area colors passed to polygon function |
border |
vector with colors for borders |
... |
other parameters passed to polygon |
smp <- rf(100, 5, 10) x <- lsdSampleDepthContours(smp) plot(x, col = paste0("grey", col = rev(seq(10, 40, 10))))
smp <- rf(100, 5, 10) x <- lsdSampleDepthContours(smp) plot(x, col = paste0("grey", col = rev(seq(10, 40, 10))))
Virtual class for robust regression methods from depthproc package
coef
coefficients of fitted model
This function generates random numbers from p-dimensional unit sphere.
runifsphere(n, p = 2)
runifsphere(n, p = 2)
n |
number of random samples. |
p |
dimension of the unit sphere. |
Daniel Kosiorowski, Mateusz Bocian, Anna Wegrzynkiewicz and Zygmunt Zawadzki from Cracow University of Economics.
x <- runifsphere(n = 100) plot(x)
x <- runifsphere(n = 100) plot(x)
Draws a scale curve: measure of dispersion.
scaleCurve( x, y = NULL, alpha = seq(0, 1, 0.01), name = "X", name_y = "Y", title = "Scale Curve", depth_params = list(method = "Projection") )
scaleCurve( x, y = NULL, alpha = seq(0, 1, 0.01), name = "X", name_y = "Y", title = "Scale Curve", depth_params = list(method = "Projection") )
x |
Multivariate data as a matrix. |
y |
Additional matrix with multivariate data. |
alpha |
Vector with values of central area to be used in computation. |
name |
Name of matrix X used in legend. |
name_y |
Name of matrix Y used in legend. |
title |
title of the plot. |
depth_params |
list of parameters for function depth (method, threads, ndir, la, lb, pdim, mean, cov, exact). |
For sample depth function ,
,
,
,
denoting
— central region, we can define the scale curve
, for
The scale curve is a two-dimensional method of describing the dispersion of random vector around the depth induced median.
Function scalecurve for determining the volumes of the convex hull containing points from alpha central regions, uses function convhulln from geometry package.
The minimal dimension of data in X or Y is 2.
ggplot2 package is used to draw a plot.
Returns the volume of the convex hull containing subsequent central points of X
.
Daniel Kosiorowski, Mateusz Bocian, Anna Wegrzynkiewicz and Zygmunt Zawadzki from Cracow University of Economics.
Liu, R.Y., Parelius, J.M. and Singh, K. (1999), Multivariate analysis by data depth: Descriptive statistics, graphics and inference (with discussion), Ann. Statist., 27, 783–858.
Chaudhuri, P. (1996), On a Geometric Notion of Quantiles for Multivariate Data, Journal of the American Statistical Association, 862–872.
Dyckerhoff, R. (2004), Data Depths Satisfying the Projection Property, Allgemeines Statistisches Archiv., 88, 163–190.
depthContour
and depthPersp
for depth graphics.
library(mvtnorm) x <- mvrnorm(n = 100, mu = c(0, 0), Sigma = 3 * diag(2)) y <- rmvt(n = 100, sigma = diag(2), df = 2) scaleCurve(x, y, depth_params = list(method = "Projection")) # Comparing two scale curves # normal distribution and mixture of normal distributions x <- mvrnorm(100, c(0, 0), diag(2)) y <- mvrnorm(80, c(0, 0), diag(2)) z <- mvrnorm(20, c(5, 5), diag(2)) scaleCurve(x, rbind(y, z), name = "N", name_y = "Mixture of N", depth_params = list(method = "Projection"))
library(mvtnorm) x <- mvrnorm(n = 100, mu = c(0, 0), Sigma = 3 * diag(2)) y <- rmvt(n = 100, sigma = diag(2), df = 2) scaleCurve(x, y, depth_params = list(method = "Projection")) # Comparing two scale curves # normal distribution and mixture of normal distributions x <- mvrnorm(100, c(0, 0), diag(2)) y <- mvrnorm(80, c(0, 0), diag(2)) z <- mvrnorm(20, c(5, 5), diag(2)) scaleCurve(x, rbind(y, z), name = "N", name_y = "Mixture of N", depth_params = list(method = "Projection"))
ScaleCurve is a class that stores results of scaleCurve function.
ScaleCurve intherits behviour from numeric vector, so raw values of ScaleCurve can be accessed via as.numeric(...).
The mechanism of creating plots with multiple curves is shown in DepthCurve-class (same mechanism is applied for AsymmetryCurve).
library(mvtnorm) x <- mvrnorm(n = 100, mu = c(0, 0), Sigma = 2 * diag(2)) y <- rmvt(n = 100, sigma = diag(2), df = 4) s1 <- scaleCurve(x, depth_params = list(method = "Projection")) s2 <- scaleCurve(y, depth_params = list(method = "Projection"), name = "Set2") sc_list <- combineDepthCurves(s1, s2) # Add one curve to another plot(sc_list) # Draw plot with two curves z <- mvrnorm(n = 100, mu = c(0, 0), Sigma = 1 * diag(2)) s3 <- scaleCurve(z, depth_params = list(method = "Projection")) plot(combineDepthCurves(sc_list, s3)) # Add third curve and draw a plot
library(mvtnorm) x <- mvrnorm(n = 100, mu = c(0, 0), Sigma = 2 * diag(2)) y <- rmvt(n = 100, sigma = diag(2), df = 4) s1 <- scaleCurve(x, depth_params = list(method = "Projection")) s2 <- scaleCurve(y, depth_params = list(method = "Projection"), name = "Set2") sc_list <- combineDepthCurves(s1, s2) # Add one curve to another plot(sc_list) # Draw plot with two curves z <- mvrnorm(n = 100, mu = c(0, 0), Sigma = 1 * diag(2)) s3 <- scaleCurve(z, depth_params = list(method = "Projection")) plot(combineDepthCurves(sc_list, s3)) # Add third curve and draw a plot
Computes projection trimmed regression in 2 dimensions.
trimProjReg2d(x, y, alpha = 0.1)
trimProjReg2d(x, y, alpha = 0.1)
x |
Independent variable |
y |
Dependent variable |
alpha |
Percentage of trimmed observations |
Zygmunt Zawadzki from Cracow University of Economics.
# EXAMPLE 1 data(pension) plot(pension) abline(lm(Reserves ~ Income, data = pension), lty = 3, lwd = 2) # lm abline(trimProjReg2d(pension[, 1], pension[, 2]), lwd = 2) # trimprojreg2d legend("bottomright", c("OLS", "TrimLS"), lty = 1:2) # EXAMPLE 2 data(under5.mort) data(inf.mort) data(maesles.imm) data2011 <- na.omit(cbind(under5.mort[, 22], inf.mort[, 22], maesles.imm[, 22])) x <- data2011[, 3] y <- data2011[, 2] plot(x, y, cex = 1.2, ylab = "infant mortality rate per 1000 live birth", xlab = "against masles immunized percentage", main = "Projection Depth Trimmed vs. LS regressions") abline(lm(x ~ y), lwd = 2, col = "black") # lm abline(trimProjReg2d(x, y), lwd = 2, col = "red") # trimmed reg legend("bottomleft", c("LS", "TrimReg"), fill = c("black", "red"), cex = 1.4, bty = "n") ##### Comparsion of a few regression methods ##### library(DepthProc) library(MASS) data("france") plot(UR ~ MW, pch = 19, data = france) # linear regression lm.fit <- lm(UR ~ MW, data = france) abline(lm.fit, lwd=2, cex=3, col='red') # M-estimator rlm.fit <- rlm(UR ~ MW, data = france) abline(rlm.fit, lwd = 2,col = "blue") # LMS lqs.lms <- lqs(UR ~ MW, method = "lms", data = france) #least median of squares# lqs.lts <- lqs(UR ~ MW, method = "lts", data = france) #least trimmed squares# abline(lqs.lms, lwd = 2, col="green") abline(lqs.lts, lwd = 2, col="pink") # Lowess lines(lowess(france$MW, france$UR, f = 0.5, iter = 0), lwd = 2) # loess # Depth trimmed regression trim.reg <- trimProjReg2d(france$MW, france$UR) #trimprojreg2d abline(trim.reg, lwd = 4, col = 'orange')
# EXAMPLE 1 data(pension) plot(pension) abline(lm(Reserves ~ Income, data = pension), lty = 3, lwd = 2) # lm abline(trimProjReg2d(pension[, 1], pension[, 2]), lwd = 2) # trimprojreg2d legend("bottomright", c("OLS", "TrimLS"), lty = 1:2) # EXAMPLE 2 data(under5.mort) data(inf.mort) data(maesles.imm) data2011 <- na.omit(cbind(under5.mort[, 22], inf.mort[, 22], maesles.imm[, 22])) x <- data2011[, 3] y <- data2011[, 2] plot(x, y, cex = 1.2, ylab = "infant mortality rate per 1000 live birth", xlab = "against masles immunized percentage", main = "Projection Depth Trimmed vs. LS regressions") abline(lm(x ~ y), lwd = 2, col = "black") # lm abline(trimProjReg2d(x, y), lwd = 2, col = "red") # trimmed reg legend("bottomleft", c("LS", "TrimReg"), fill = c("black", "red"), cex = 1.4, bty = "n") ##### Comparsion of a few regression methods ##### library(DepthProc) library(MASS) data("france") plot(UR ~ MW, pch = 19, data = france) # linear regression lm.fit <- lm(UR ~ MW, data = france) abline(lm.fit, lwd=2, cex=3, col='red') # M-estimator rlm.fit <- rlm(UR ~ MW, data = france) abline(rlm.fit, lwd = 2,col = "blue") # LMS lqs.lms <- lqs(UR ~ MW, method = "lms", data = france) #least median of squares# lqs.lts <- lqs(UR ~ MW, method = "lts", data = france) #least trimmed squares# abline(lqs.lms, lwd = 2, col="green") abline(lqs.lts, lwd = 2, col="pink") # Lowess lines(lowess(france$MW, france$UR, f = 0.5, iter = 0), lwd = 2) # loess # Depth trimmed regression trim.reg <- trimProjReg2d(france$MW, france$UR) #trimprojreg2d abline(trim.reg, lwd = 4, col = 'orange')
Children under 5 months mortality rate per 1,000 live births
data(under5.mort)
data(under5.mort)
A data frame with 654 rows and 4 variables
http://mdgs.un.org/unsd/mdg/Data.aspx
US Labour dataset
data(USLABOUR)
data(USLABOUR)
A data frame with 654 rows and 4 variables
U.S.Department of Labor — Bureau of Labour Statistics FRED