Title: | Multivariate Statistical Methods |
---|---|
Description: | Algorithms to build set partitions and commutator matrices and their use in the construction of multivariate d-Hermite polynomials; estimation and derivation of theoretical vector moments and vector cumulants of multivariate distributions; conversion formulae for multivariate moments and cumulants. Applications to estimation and derivation of multivariate measures of skewness and kurtosis; estimation and derivation of asymptotic covariances for d-variate Hermite polynomials, multivariate moments and cumulants and measures of skewness and kurtosis. The formulae implemented are discussed in Terdik (2021, ISBN:9783030813925), "Multivariate Statistical Methods". |
Authors: | Gyorgy Terdik [aut], Emanuele Taufer [aut, cre] |
Maintainer: | Emanuele Taufer <[email protected]> |
License: | GPL-3 |
Version: | 2.0.0 |
Built: | 2025-02-22 02:54:18 UTC |
Source: | https://github.com/cran/MultiStatM |
This function calculates the commutator index based on the specified type. The available types are "Kmn", "Kperm", "Mixing", and "Moment". Depending on the selected type, the corresponding specific function is called.
CommutatorIndx(Type, ...)
CommutatorIndx(Type, ...)
Type |
a string specifying the type of commutator index to be calculated. Must be one of "Kmn", "Kperm", "Mixing", or "Moment". |
... |
additional arguments passed to the specific commutator function. |
The function 'CommutatorIndx' acts as a wrapper to call specific commutator functions based on the input 'Type'.
Type "Kmn":
m
- Row-dimension.
n
- Col-dimension.
A vector of indexes to provide the commutation, transforming vec A to vec of the transposed A.
Type "Kperm":
perm
- Vector indicating the permutation of the order in the Kronecker product.
dims
- Vector indicating the dimensions of the vectors.
An index vector to produce the permutation of the Kronecker products of vectors of any length.
Type "Mixing":
x
- A vector of dimension prod(d1)*prod(d2)
.
d1
- Dimension of the first group of vectors.
d2
- Dimension of the second group of vectors.
A vector Kx representing the product of the moment commutator and the vector x.
Type "Moment":
x
- A vector of length d^n
where n is the length of el_rm
.
el_rm
- Type of a partition.
d
- Dimensionality of the underlying multivariate distribution.
A vector Kx representing the product of the moment commutator and the vector x.
A vector representing the commutator index.
Other Commutators:
CommutatorMatr()
# Kmn example A <- 1:6 CommutatorIndx(Type = "Kmn", m = 3, n = 2) # Kperm example a1 <- c(1, 2) a2 <- c(2, 3, 4) a3 <- c(1, 3) p1 <- a1 %x% a2 %x% a3 CommutatorIndx(Type = "Kperm", perm = c(3, 1, 2), dims = c(2, 3, 2)) # Mixing example d1 <- c(2, 3, 2) d2 <- c(3, 2, 2) x <- 1:(prod(d1) * prod(d2)) CommutatorIndx(Type = "Mixing", x = x, d1 = d1, d2 = d2) # Moment example n <- 4 r <- 2 m <- 1 d <- 2 PTA <- PartitionTypeAll(n) el_r <- PTA$eL_r[[r]][m, ] x <- 1:d^n CommutatorIndx(Type = "Moment", x = x, el_rm = el_r, d = d)
# Kmn example A <- 1:6 CommutatorIndx(Type = "Kmn", m = 3, n = 2) # Kperm example a1 <- c(1, 2) a2 <- c(2, 3, 4) a3 <- c(1, 3) p1 <- a1 %x% a2 %x% a3 CommutatorIndx(Type = "Kperm", perm = c(3, 1, 2), dims = c(2, 3, 2)) # Mixing example d1 <- c(2, 3, 2) d2 <- c(3, 2, 2) x <- 1:(prod(d1) * prod(d2)) CommutatorIndx(Type = "Mixing", x = x, d1 = d1, d2 = d2) # Moment example n <- 4 r <- 2 m <- 1 d <- 2 PTA <- PartitionTypeAll(n) el_r <- PTA$eL_r[[r]][m, ] x <- 1:d^n CommutatorIndx(Type = "Moment", x = x, el_rm = el_r, d = d)
This function generates various types of commutator matrices.
CommutatorMatr(Type, ...)
CommutatorMatr(Type, ...)
Type |
A string specifying the type of commutator matrix. Choices are "Kmn", "Kperm", "Mixing", or "Moment". |
... |
Additional arguments specific to the type of commutator matrix (see Details). |
The function CommutatorMatr
supports the following types of commutator matrices:
Description: Transforms vec(A)
to vec(A^T)
, where A^T
is the transpose of matrix A
. An option for sparse matrix is provided. By default, a non-sparse matrix is produced. Using sparse matrices increases computation times but requires far less memory.
Arguments:
m
(integer)Number of rows of the first matrix.
n
(integer)Number of columns of the first matrix.
useSparse
(logical, optional)If TRUE, returns a sparse matrix. Default is FALSE.
Description: Generates a commutation matrix for a specified permutation of matrix dimensions. An option for sparse matrix is provided. By default, a non-sparse matrix is produced. Using sparse matrices increases computation times but requires far less memory. Arguments:
perm
(integer vector)The permutation vector.
dims
(integer vector)The dimensions of the matrices involved.
useSparse
(logical, optional)If TRUE, returns a sparse matrix. Default is FALSE.
Description: Generates the Mixing commutation matrix used in linear algebra transformations involving tensor products. An option for sparse matrix is provided. By default, a non-sparse matrix is produced. Using sparse matrices increases computation times but requires far less memory. Arguments:
d1
(integer vector)Dimensions of the first set.
d2
(integer vector)Dimensions of the second set.
useSparse
(logical, optional)If TRUE, returns a sparse matrix. Default is FALSE.
Description: Generates the Moment commutation matrix based on partitioning of moments. An option for sparse matrix is provided. By default, a non-sparse matrix is produced. Using sparse matrices increases computation times but requires far less memory. Arguments:
el_rm
(integer vector)Elements of the partition.
d
(integer)Dimension of the partition.
useSparse
(logical, optional)If TRUE, returns a sparse matrix. Default is FALSE.
Depending on the type:
A commutation matrix of dimension . If
useSparse=TRUE
, an object of class "dgCMatrix" is produced.
A square permutation matrix of size prod(dims)
. If useSparse=TRUE
, an object of class "dgCMatrix" is produced.
A square matrix of dimension prod(d1) * prod(d2)
. If useSparse=TRUE
, an object of class "dgCMatrix" is produced.
A commutator matrix for moment formulae.
Other Commutators:
CommutatorIndx()
# Example for Kmn CommutatorMatr("Kmn", m = 3, n = 2) # Example for Kperm dims <- c(2, 3, 2) perm <- c(1, 3, 2) CommutatorMatr("Kperm", perm = perm, dims = dims) # Example for Mixing d1 <- c(2, 3, 2) d2 <- c(3, 2, 2) CommutatorMatr("Mixing", d1 = d1, d2 = d2) # Example for Moment n <- 4 r <- 2 m <- 1 d <- 2 PTA <- PartitionTypeAll(n) el_r <- PTA$eL_r[[r]][m,] CommutatorMatr("Moment", el_r = el_r, d = d)
# Example for Kmn CommutatorMatr("Kmn", m = 3, n = 2) # Example for Kperm dims <- c(2, 3, 2) perm <- c(1, 3, 2) CommutatorMatr("Kperm", perm = perm, dims = dims) # Example for Mixing d1 <- c(2, 3, 2) d2 <- c(3, 2, 2) CommutatorMatr("Mixing", d1 = d1, d2 = d2) # Example for Moment n <- 4 r <- 2 m <- 1 d <- 2 PTA <- PartitionTypeAll(n) el_r <- PTA$eL_r[[r]][m,] CommutatorMatr("Moment", el_r = el_r, d = d)
Obtains a vector of moments from a vector of cumulants for either univariate or multivariate data.
Cum2Mom(cumulants, Type = c("Univariate", "Multivariate"))
Cum2Mom(cumulants, Type = c("Univariate", "Multivariate"))
cumulants |
Either a vector of univariate cumulants or a list of vectors of multivariate cumulants. |
Type |
A character string specifying the type of cumulants provided. Use "Univariate" for univariate cumulants and "Multivariate" for multivariate cumulants. |
The vector of moments if Type
is "Univariate" or the list of vectors of moments if Type
is "Multivariate".
Gy. Terdik, Multivariate statistical methods - going beyond the linear, Springer 2021. Section 3.4.
Other Moments and cumulants:
EVSKSkewNorm()
,
EVSKUniS()
,
Mom2Cum()
,
MomCumCFUSN()
,
MomCumSkewNorm()
,
MomCumUniS()
,
MomCumZabs()
# Univariate example cum_x <- c(1, 2, 3, 4) Cum2Mom(cum_x, Type = "Univariate") # Multivariate example cum <- list(c(0,0), c(1,0,0,1), c(rep(0,8)), c(rep(0,16)), c(rep(0,32))) Cum2Mom(cum, Type = "Multivariate")
# Univariate example cum_x <- c(1, 2, 3, 4) Cum2Mom(cum_x, Type = "Univariate") # Multivariate example cum <- list(c(0,0), c(1,0,0,1), c(rep(0,8)), c(rep(0,16)), c(rep(0,32))) Cum2Mom(cum, Type = "Multivariate")
Eliminates the duplicated/q-plicated elements in a T-vector of multivariate moments and cumulants. Produces the same results as EliminMatr. Note EliminIndx does not provide the same results as unique()
EliminIndx(d, q)
EliminIndx(d, q)
d |
dimension of a vector x |
q |
power of the Kronecker product |
A vector of indexes of the distinct elements in the T-vector
Gy. Terdik, Multivariate statistical methods - going beyond the linear, Springer 2021. Section 1.3.2 Multi-Indexing, Elimination, and Duplication, p.21,(1.32)
Other Matrices and commutators:
EliminMatr()
,
QplicIndx()
,
QplicMatr()
,
SymIndx()
,
SymMatr()
,
UnivMomCum()
x<-c(1,0,3) y<-kronecker(x,kronecker(x,x)) y[EliminIndx(3,3)] ## Not the same results as unique(y)
x<-c(1,0,3) y<-kronecker(x,kronecker(x,x)) y[EliminIndx(3,3)] ## Not the same results as unique(y)
Eliminates the duplicated/q-plicated elements in a T-vector of multivariate moments and cumulants.
EliminMatr(d, q, useSparse = FALSE)
EliminMatr(d, q, useSparse = FALSE)
d |
dimension of a vector x |
q |
power of the Kronecker product |
useSparse |
TRUE or FALSE. |
Elimination matrix of order .
If
useSparse=TRUE
an object of the class "dgCMatrix" is produced.
Gy. Terdik, Multivariate statistical methods - going beyond the linear, Springer 2021. Section 1.3.2 Multi-Indexing, Elimination, and Duplication, p.21,(1.32)
Other Matrices and commutators:
EliminIndx()
,
QplicIndx()
,
QplicMatr()
,
SymIndx()
,
SymMatr()
,
UnivMomCum()
x<-c(1,2,3) y<-kronecker(kronecker(x,x),x) ## Distinct elements of y z<-as.matrix(EliminMatr(3,3))%*%y ## Restore eliminated elements in z as.vector(QplicMatr(3,3)%*%z)
x<-c(1,2,3) y<-kronecker(kronecker(x,x),x) ## Distinct elements of y z<-as.matrix(EliminMatr(3,3))%*%y ## Restore eliminated elements in z as.vector(QplicMatr(3,3)%*%z)
Computes the theoretical values of the mean vector, covariance, skewness vector, total skenwness, kurtosis vector and total kurtosis for the multivariate Skew Normal distribution
EVSKSkewNorm(omega, alpha)
EVSKSkewNorm(omega, alpha)
omega |
A |
alpha |
shape parameter d-vector |
A list of theoretical values for the mean vector, covariance, skewness vector, total skenwness, kurtosis vector and total kurtosis
Gy.Terdik, Multivariate statistical methods - Going beyond the linear, Springer 2021 (5.5) p.247
S. R. Jammalamadaka, E. Taufer, Gy. Terdik. On multivariate skewness and kurtosis. Sankhya A, 83(2), 607-644.
Other Moments and cumulants:
Cum2Mom()
,
EVSKUniS()
,
Mom2Cum()
,
MomCumCFUSN()
,
MomCumSkewNorm()
,
MomCumUniS()
,
MomCumZabs()
alpha<-c(10,5,0) omega<-diag(3) EVSKSkewNorm(omega,alpha)
alpha<-c(10,5,0) omega<-diag(3) EVSKSkewNorm(omega,alpha)
Cumulants (up to the 4th order), skewness, and kurtosis of the d-variate Uniform distribution on the sphere or the modulus of the d-variate Uniform distribution on the sphere.
EVSKUniS(d, nCum = TRUE, Type = c("Standard", "Modulus"))
EVSKUniS(d, nCum = TRUE, Type = c("Standard", "Modulus"))
d |
dimensions |
nCum |
if it is FALSE then moments (up to the 4th order) are calculated. |
Type |
specify the type of distribution: "Standard" for the Uniform distribution on the sphere, or "Modulus" for the modulus of the Uniform distribution on the sphere. |
A list of computed moments and cumulants.
When Type is "Standard":
EU1 |
Mean vector |
varU |
Covariance matrix |
Skew.U |
Skewness vector (always zero) |
Skew.tot |
Total skewness (always zero) |
Kurt.U |
Kurtosis vector |
Kurt.tot |
Total kurtosis |
When Type is "Modulus":
EU1 |
Mean vector |
varU |
Covariance matrix |
EU.k |
List of moments up to 4th order |
cumU.k |
List of cumulants up to 4th order |
skew.U |
Skewness vector |
kurt.U |
Kurtosis vector |
Gy. Terdik, Multivariate statistical methods - Going beyond the linear, Springer 2021 Proposition 5.3 p.297
S. R. Jammalamadaka, E. Taufer, Gy. Terdik. On multivariate skewness and kurtosis. Sankhya A, 83(2), 607-644.
Gy. Terdik, Multivariate statistical methods - Going beyond the linear, Springer 2021, Lemma 5.12 p.298
Other Moments and cumulants:
Cum2Mom()
,
EVSKSkewNorm()
,
Mom2Cum()
,
MomCumCFUSN()
,
MomCumSkewNorm()
,
MomCumUniS()
,
MomCumZabs()
# Example for Standard type EVSKUniS(d=3, Type="Standard") # Example for Modulus type EVSKUniS(d=3, Type="Modulus")
# Example for Standard type EVSKUniS(d=3, Type="Standard") # Example for Modulus type EVSKUniS(d=3, Type="Modulus")
Provides the coefficients of Hermite polynomials, either univariate or multivariate.
HermiteCoeff(Type, N, d = NULL)
HermiteCoeff(Type, N, d = NULL)
Type |
A character string specifying the type of Hermite polynomial. Must be either "Univariate" or "Multivariate". |
N |
The order of polynomial. Required for both types. |
d |
The dimension of the d-variate X. Required only for multivariate type. |
For 'Type = "Univariate"', returns a vector of coefficients of ,
, etc.
For 'Type = "Multivariate"', returns a list of matrices of coefficients for the d-variate polynomials from 1 to N.
Gy. Terdik, Multivariate statistical methods - going beyond the linear, Springer 2021. Sections 4.4 (4.24) and 4.6.2, p. 223, Remark 4.8
Other Hermite Polynomials:
HermiteCov12()
,
HermiteN()
,
HermiteN2X()
# Univariate example H_uni <- HermiteCoeff(Type = "Univariate", N = 5) # Multivariate example N <- 5; d <- 3 H_multi <- HermiteCoeff(Type = "Multivariate", N = N, d = d) X <- c(1:3) X3 <- kronecker(X, kronecker(X, X)) X5 <- kronecker(X3, kronecker(X, X)) Idv <- as.vector(diag(d)) # vector of variance matrix # value of H5 at X is vH5 <- H_multi[[1]] %*% X5 + H_multi[[2]] %*% kronecker(Idv, X3) + H_multi[[3]] %*% kronecker(kronecker(Idv, Idv), X)
# Univariate example H_uni <- HermiteCoeff(Type = "Univariate", N = 5) # Multivariate example N <- 5; d <- 3 H_multi <- HermiteCoeff(Type = "Multivariate", N = N, d = d) X <- c(1:3) X3 <- kronecker(X, kronecker(X, X)) X5 <- kronecker(X3, kronecker(X, X)) Idv <- as.vector(diag(d)) # vector of variance matrix # value of H5 at X is vH5 <- H_multi[[1]] %*% X5 + H_multi[[2]] %*% kronecker(Idv, X3) + H_multi[[3]] %*% kronecker(kronecker(Idv, Idv), X)
Computation of the covariance matrix between d-variate T-Hermite polynomials
and
.
HermiteCov12(SigX12, N)
HermiteCov12(SigX12, N)
SigX12 |
Covariance matrix of the Gaussian vectors X1 and X2 respectively of dimensions d1 and d2 |
N |
Common degree of the multivariate Hermite polynomials |
Covariance matrix of and
Gy.Terdik, Multivariate statistical methods - going beyond the linear, Springer 2021. (4.59), (4.66),
Other Hermite Polynomials:
HermiteCoeff()
,
HermiteN()
,
HermiteN2X()
Covmat<-matrix(c(1,0.8,0.8,1),2,2) Cov_X1_X2 <- HermiteCov12(Covmat,3)
Covmat<-matrix(c(1,0.8,0.8,1),2,2) Cov_X1_X2 <- HermiteCov12(Covmat,3)
Computes either univariate or multivariate Hermite polynomials up to a specified order.
HermiteN(x, N, Type, sigma2 = 1, Sig2 = diag(length(x)))
HermiteN(x, N, Type, sigma2 = 1, Sig2 = diag(length(x)))
x |
A scalar (for univariate) or a vector (for multivariate) at which to evaluate the Hermite polynomials. |
N |
The maximum order of the polynomials. |
Type |
A character string specifying the type of Hermite polynomials to compute. Can be either "Univariate" or "Multivariate". |
sigma2 |
The variance for univariate Hermite polynomials. Default is 1. (Only used if Type is "Univariate"). |
Sig2 |
The covariance matrix for multivariate Hermite polynomials. Default is the unit matrix diag(length(x)). (Only used if Type is "Multivariate"). |
Depending on the value of the 'Type' parameter, this function computes either the univariate or the multivariate Hermite polynomials.
Depending on the type, the function returns:
Univariate
: A vector of univariate Hermite polynomials with degrees from 1 to N evaluated at x.
Multivariate
: A list of multivariate polynomials of order from 1 to N evaluated at vector x.
Gy.Terdik, Multivariate statistical methods - going beyond the linear, Springer 2021. Section 4.1 (for univariate), Section 4.6.2, (4.73), p.223 (for multivariate).
Other Hermite Polynomials:
HermiteCoeff()
,
HermiteCov12()
,
HermiteN2X()
# Univariate example HermiteN(x = 1, N = 3, Type = "Univariate") # Multivariate example HermiteN(x = c(1, 3), N = 3, Type = "Multivariate", Sig2 = diag(2))
# Univariate example HermiteN(x = 1, N = 3, Type = "Univariate") # Multivariate example HermiteN(x = c(1, 3), N = 3, Type = "Multivariate", Sig2 = diag(2))
Compute the inverse of univariate or multivariate Hermite polynomials.
HermiteN2X(Type, H_N, N, Sig2 = NULL)
HermiteN2X(Type, H_N, N, Sig2 = NULL)
Type |
A string specifying the type of Hermite polynomial inversion. Must be either "Univariate" or "Multivariate". |
H_N |
Input Hermite polynomials. For univariate, it is a vector. For multivariate, it is a list. |
N |
The highest polynomial order. |
Sig2 |
The variance matrix of x for multivariate, or variance for univariate. Defaults to identity matrix for multivariate and 1 for univariate. |
This function computes the powers of x when Hermite polynomials are given. Depending on the type specified, it handles either univariate or multivariate Hermite polynomials.
A list of x powers: ,
, ... ,
for multivariate,
or a vector of x powers:
,
for univariate.
Gy.Terdik, Multivariate statistical methods - going beyond the linear, Springer 2021. Section 4.6.2, (4.72), p.223 and Section 4.4, (4.23), p.198.
Other Hermite Polynomials:
HermiteCoeff()
,
HermiteCov12()
,
HermiteN()
# Univariate example H_N_x <- c(1, 2, 3, 4) x_powers <- HermiteN2X(Type = "Univariate", H_N = H_N_x, N = 4, Sig2 = 1) # Multivariate example x <- c(1, 3) Sig2 <- diag(length(x)) N <- 4 H_N_X <- HermiteN(x, N, Type="Multivariate") x_ad_n <- HermiteN2X(Type = "Multivariate", H_N = H_N_X, N = N, Sig2 = Sig2)
# Univariate example H_N_x <- c(1, 2, 3, 4) x_powers <- HermiteN2X(Type = "Univariate", H_N = H_N_x, N = 4, Sig2 = 1) # Multivariate example x <- c(1, 3) Sig2 <- diag(length(x)) N <- 4 H_N_X <- HermiteN(x, N, Type="Multivariate") x_ad_n <- HermiteN2X(Type = "Multivariate", H_N = H_N_X, N = N, Sig2 = Sig2)
Obtains a vector of cumulants from a vector of moments for either univariate or multivariate data.
Mom2Cum(moments, Type = c("Univariate", "Multivariate"))
Mom2Cum(moments, Type = c("Univariate", "Multivariate"))
moments |
Either a vector of univariate moments or a list of vectors of multivariate moments. |
Type |
A character string specifying the type of moments provided. Use "Univariate" for univariate moments and "Multivariate" for multivariate moments. |
The vector of cumulants if Type
is "Univariate" or the list of vectors of cumulants if Type
is "Multivariate".
Gy. Terdik, Multivariate statistical methods - going beyond the linear, Springer 2021. Section 3.4.
Other Moments and cumulants:
Cum2Mom()
,
EVSKSkewNorm()
,
EVSKUniS()
,
MomCumCFUSN()
,
MomCumSkewNorm()
,
MomCumUniS()
,
MomCumZabs()
# Univariate example mu_x <- c(1, 2, 3, 4) Mom2Cum(mu_x, Type = "Univariate") # Multivariate example mu <- list(c(0,0), c(1,0,0,1), c(0,0,0,0,0,0,0,0), c(3,0,0,1,0,1,1,0,0,1,1,0,1,0,0,3), c(rep(0,32))) Mom2Cum(mu, Type = "Multivariate")
# Univariate example mu_x <- c(1, 2, 3, 4) Mom2Cum(mu_x, Type = "Univariate") # Multivariate example mu <- list(c(0,0), c(1,0,0,1), c(0,0,0,0,0,0,0,0), c(3,0,0,1,0,1,1,0,0,1,1,0,1,0,0,3), c(rep(0,32))) Mom2Cum(mu, Type = "Multivariate")
Provides the theoretical cumulants of the multivariate Canonical Fundamental Skew Normal distribution
MomCumCFUSN(r, d, p, Delta, nMu = FALSE)
MomCumCFUSN(r, d, p, Delta, nMu = FALSE)
r |
The highest cumulant order |
d |
The multivariate dimension and number of rows of the skewness matrix Delta |
p |
The number of cols of the skewness matrix Delta |
Delta |
The skewness matrix |
nMu |
If set to TRUE, the list of the first r d-variate moments is provided |
The list of theoretical cumulants in vector form
Gy.Terdik, Multivariate statistical methods - Going beyond the linear, Springer 2021, Lemma 5.3 p.251
Other Moments and cumulants:
Cum2Mom()
,
EVSKSkewNorm()
,
EVSKUniS()
,
Mom2Cum()
,
MomCumSkewNorm()
,
MomCumUniS()
,
MomCumZabs()
r <- 4; d <- 2; p <- 3 Lamd <- matrix(sample(1:50-25, d*p), nrow=d) ieg<- eigen(diag(p)+t(Lamd)%*%Lamd) V <- ieg$vectors Delta <-Lamd %*% V %*% diag(1/sqrt(ieg$values)) %*% t(V) MomCum <- MomCumCFUSN(r,d,p,Delta)
r <- 4; d <- 2; p <- 3 Lamd <- matrix(sample(1:50-25, d*p), nrow=d) ieg<- eigen(diag(p)+t(Lamd)%*%Lamd) V <- ieg$vectors Delta <-Lamd %*% V %*% diag(1/sqrt(ieg$values)) %*% t(V) MomCum <- MomCumCFUSN(r,d,p,Delta)
Computes the theoretical values of moments and cumulants up to the r-th order. Warning: if nMu = TRUE it can be very slow
MomCumSkewNorm(r = 4, omega, alpha, nMu = FALSE)
MomCumSkewNorm(r = 4, omega, alpha, nMu = FALSE)
r |
the highest moment and cumulant order |
omega |
A |
alpha |
shape parameter d-vector |
nMu |
if it is TRUE then moments are calculated as well |
A list of theoretical moments and cumulants
Gy.Terdik, Multivariate statistical methods - Going beyond the linear, Springer 2021 (5.5) p.247, Lemma 5.1 p. 246
S. R. Jammalamadaka, E. Taufer, Gy. Terdik. On multivariate skewness and kurtosis. Sankhya A, 83(2), 607-644.
Other Moments and cumulants:
Cum2Mom()
,
EVSKSkewNorm()
,
EVSKUniS()
,
Mom2Cum()
,
MomCumCFUSN()
,
MomCumUniS()
,
MomCumZabs()
alpha<-c(10,5,0) omega<-diag(3) MomCumSkewNorm(r=4,omega,alpha)
alpha<-c(10,5,0) omega<-diag(3) MomCumSkewNorm(r=4,omega,alpha)
By default, only moments are provided
MomCumUniS(r, d, nCum = FALSE)
MomCumUniS(r, d, nCum = FALSE)
r |
highest order of moments and cumulants |
d |
dimension |
nCum |
if it is TRUE then cumulants are calculated |
The list of moments and cumulants in vector form
Gy.Terdik, Multivariate statistical methods - Going beyond the linear, Springer 2021 Proposition 5.3 p.297
Other Moments and cumulants:
Cum2Mom()
,
EVSKSkewNorm()
,
EVSKUniS()
,
Mom2Cum()
,
MomCumCFUSN()
,
MomCumSkewNorm()
,
MomCumZabs()
# The first four moments for d=3 MomCumUniS(4,3,nCum=0) # The first four moments and cumulants for d=3 MomCumUniS(4,3,nCum=4)
# The first four moments for d=3 MomCumUniS(4,3,nCum=0) # The first four moments and cumulants for d=3 MomCumUniS(4,3,nCum=4)
Provides the theoretical moments and cumulants of the Central Folded Normal distribution. Depending on the choice of 'Type', either the univariate or d-variate distribution is used.
MomCumZabs(r, d, Type, nCum = FALSE)
MomCumZabs(r, d, Type, nCum = FALSE)
r |
The highest moment (cumulant) order. |
d |
Integer; the dimension of the distribution. Must be 1 when 'Type' is "Univariate" and greater than 1 when 'Type' is "Multivariate". |
Type |
Character; specifies the type of distribution. Must be either "Univariate" or "Multivariate". |
nCum |
Logical; if TRUE, then cumulants are calculated. |
A list containing moments and optionally cumulants.
For "Univariate" type:
MuZ
: The moments of the univariate Central Folded Normal distribution.
CumZ
: The cumulants of the univariate Central Folded Normal distribution.
For "Multivariate" type:
MuZ
: The moments of the d-variate Central Folded Normal distribution.
CumZ
: The cumulants of the d-variate Central Folded Normal distribution.
Gy.Terdik, Multivariate statistical methods - Going beyond the linear, Springer 2021, Proposition 5.1 p.242 and formula: p. 301
Other Moments and cumulants:
Cum2Mom()
,
EVSKSkewNorm()
,
EVSKUniS()
,
Mom2Cum()
,
MomCumCFUSN()
,
MomCumSkewNorm()
,
MomCumUniS()
# Univariate case: The first three moments MomCumZabs(3, 1, Type = "Univariate") # Univariate case: The first three moments and cumulants MomCumZabs(3, 1, Type = "Univariate",nCum = TRUE) # d-variate case: The first three moments MomCumZabs(3, 2, Type = "Multivariate" ) # d-variate case: The first three moments and cumulants MomCumZabs(3, d=2, Type = "Multivariate", nCum = TRUE)
# Univariate case: The first three moments MomCumZabs(3, 1, Type = "Univariate") # Univariate case: The first three moments and cumulants MomCumZabs(3, 1, Type = "Univariate",nCum = TRUE) # d-variate case: The first three moments MomCumZabs(3, 2, Type = "Multivariate" ) # d-variate case: The first three moments and cumulants MomCumZabs(3, d=2, Type = "Multivariate", nCum = TRUE)
For data formed by d-variate vectors x with sample covariance S and sample mean M,
it computes the values
MVStandardize(x)
MVStandardize(x)
x |
a multivariate data matrix, sample size is the number of rows |
a matrix of multivariate data with null mean vector and identity sample covariance matrix
x<-MASS::mvrnorm(1000,c(0,0,1,3),diag(4)) z<-MVStandardize(x) mu_z<- apply(z,2,mean) cov_z<- cov(z)
x<-MASS::mvrnorm(1000,c(0,0,1,3),diag(4)) z<-MVStandardize(x) mu_z<- apply(z,2,mean) cov_z<- cov(z)
A unified function to compute different types of partitions. Depending on the partition type specified, it calls the appropriate function: Partition_2Perm, Partition_DiagramsClosedNoLoops, Partition_Indecomposable, or Partition_Pairs.
Partitions(Type, ...)
Partitions(Type, ...)
Type |
A character string specifying the type of partion to compute. Choose from "2Perm", "Diagram", "Indecomp", "Pairs". |
... |
Additional arguments passed to the specific partition function:
|
Depending on the commutator type:
A vector with the elements 1 to N permuted according to L.
The list of partition matrices indecomposable with respect to L, representing diagrams without loops.
A list of partition matrices indecomposable with respect to L and a vector indicating the number of indecomposable partitions by sizes.
The list of partition matrices with blocks containing two elements. The list is empty if N is odd.
Other Partitions:
PartitionTypeAll()
,
PermutationInv()
# Example for 2Perm PA <- PartitionTypeAll(4) Partitions("2Perm", L = PA$Part.class[[3]]) # Example for Diagram L <- matrix(c(1,1,0,0,0,0,1,1),2,4,byrow=TRUE) Partitions("Diagram", L = L) # Example for Indecomp L <- matrix(c(1,1,0,0,0,0,1,1),2,4,byrow=TRUE) Partitions("Indecomp", L = L) # Example for Pairs Partitions("Pairs", N = 4)
# Example for 2Perm PA <- PartitionTypeAll(4) Partitions("2Perm", L = PA$Part.class[[3]]) # Example for Diagram L <- matrix(c(1,1,0,0,0,0,1,1),2,4,byrow=TRUE) Partitions("Diagram", L = L) # Example for Indecomp L <- matrix(c(1,1,0,0,0,0,1,1),2,4,byrow=TRUE) Partitions("Indecomp", L = L) # Example for Pairs Partitions("Pairs", N = 4)
Generates all partitions of N
numbers and classify them by type
PartitionTypeAll(N)
PartitionTypeAll(N)
N |
The (integer) number of elements to be partitioned |
Part.class
The list of all possible partitions given as partition matrices
S_N_r
A vector with the number of partitions of size r=1, r=2, etc. (Stirling numbers of second kind )
eL_r
A list of partition types with respect to partitions of size r=1, r=2, etc.
S_r_j
Vectors of number of partitions with given types grouped by partitions of size r=1, r=2, etc.
Gy. Terdik, Multivariate statistical methods - going beyond the linear, Springer 2021. Case 1.4, p.31 and Example 1.18, p.32.
Other Partitions:
Partitions()
,
PermutationInv()
# See Example 1.18, p. 32, reference below PTA<-PartitionTypeAll(4) # Partitions generated PTA$Part.class # Partitions of size 2 includes two types PTA$eL_r[[2]] # Number of partitions with r=1 blocks, r=2 blocks, etc- PTA$S_N_r # Number of different types collected by partitions of size r=1, r=2, etc. PTA$S_r_j # Partitions with size r=2, includes two types (above) each with number PTA$S_r_j[[2]]
# See Example 1.18, p. 32, reference below PTA<-PartitionTypeAll(4) # Partitions generated PTA$Part.class # Partitions of size 2 includes two types PTA$eL_r[[2]] # Number of partitions with r=1 blocks, r=2 blocks, etc- PTA$S_N_r # Number of different types collected by partitions of size r=1, r=2, etc. PTA$S_r_j # Partitions with size r=2, includes two types (above) each with number PTA$S_r_j[[2]]
Inverse of a Permutation
PermutationInv(permutation0)
PermutationInv(permutation0)
permutation0 |
A permutation of numbers 1:n |
A vector containing the inverse permutation of permutation0
Gy. Terdik, Multivariate statistical methods - going beyond the linear, Springer 2021, Remark 1.1, p.2
Other Partitions:
PartitionTypeAll()
,
Partitions()
Restores the duplicated/q-plicated elements which are eliminated by EliminMatr or EliminIndx in a T-product of vectors of dimension d. It produces the same results as QplicMatr.
QplicIndx(d, q)
QplicIndx(d, q)
d |
dimension of the vectors in the T-product |
q |
power of the Kronecker product |
A vector (T-vector) with all elements previously eliminated by EliminIndx
Gy. Terdik, Multivariate statistical methods - going beyond the linear, Springer 2021, p.21, (1.31)
Other Matrices and commutators:
EliminIndx()
,
EliminMatr()
,
QplicMatr()
,
SymIndx()
,
SymMatr()
,
UnivMomCum()
x<-c(1,2,3) y<-kronecker(kronecker(x,x),x) ## Distinct elements of y z<-y[EliminIndx(3,3)] ## Restore eliminated elements in z z[QplicIndx(3,3)]
x<-c(1,2,3) y<-kronecker(kronecker(x,x),x) ## Distinct elements of y z<-y[EliminIndx(3,3)] ## Restore eliminated elements in z z[QplicIndx(3,3)]
Restores the duplicated/q-plicated elements which are eliminated by EliminMatr in a T-product of vectors of dimension d.
QplicMatr(d, q, useSparse = FALSE)
QplicMatr(d, q, useSparse = FALSE)
d |
dimension of a vector x |
q |
power of the Kronecker product |
useSparse |
TRUE or FALSE. |
Note: since the algorithm of elimination is not unique, q-plication works together with the function EliminMatr only.
Qplication matrix of order , see (1.30), p.15.
If
useSparse=TRUE
an object of the class "dgCMatrix" is produced.
Gy. Terdik, Multivariate statistical methods - going beyond the linear, Springer 2021, p.21, (1.31)
Other Matrices and commutators:
EliminIndx()
,
EliminMatr()
,
QplicIndx()
,
SymIndx()
,
SymMatr()
,
UnivMomCum()
x<-c(1,2,3) y<-kronecker(kronecker(x,x),x) ## Distinct elements of y z<-as.matrix(EliminMatr(3,3))%*%y ## Restore eliminated elements in z as.vector(QplicMatr(3,3)%*%z)
x<-c(1,2,3) y<-kronecker(kronecker(x,x),x) ## Distinct elements of y z<-as.matrix(EliminMatr(3,3))%*%y ## Restore eliminated elements in z as.vector(QplicMatr(3,3)%*%z)
Generate random d-vectors from the multivariate Canonical Fundamental Skew-Normal (CFUSN) distribution
rCFUSN(n, Delta)
rCFUSN(n, Delta)
n |
The number of variates to be generated |
Delta |
Correlation matrix, the skewness matrix Delta |
A random matrix
Gy.Terdik, Multivariate statistical methods - Going beyond the linear, Springer 2021 (5.5) p.247
S. R. Jammalamadaka, E. Taufer, Gy. Terdik. On multivariate skewness and kurtosis. Sankhya A, 83(2), 607-644.
Other Random generation:
rCFUSSD()
,
rSkewNorm()
,
rUniS()
d <- 2; p <- 3 Lamd <- matrix(sample(1:50-25, d*p), nrow=d) ieg<- eigen(diag(p)+t(Lamd)%*%Lamd) V <- ieg$vectors Delta <-Lamd %*% V %*% diag(1/sqrt(ieg$values)) %*% t(V) x<-rCFUSN(20,Delta)
d <- 2; p <- 3 Lamd <- matrix(sample(1:50-25, d*p), nrow=d) ieg<- eigen(diag(p)+t(Lamd)%*%Lamd) V <- ieg$vectors Delta <-Lamd %*% V %*% diag(1/sqrt(ieg$values)) %*% t(V) x<-rCFUSN(20,Delta)
Generate random d-vectors from the multivariate Canonical Fundamental Skew-Spherical distribution (CFUSSD) with Gamma generator
rCFUSSD(n, d, p, a, b, Delta)
rCFUSSD(n, d, p, a, b, Delta)
n |
sample size |
d |
dimension |
p |
dimension of the first term of (5.5) |
a |
shape parameter of the Gamma generator |
b |
scale parameter of the Gamma generator |
Delta |
skewness matrix |
A matrix of random numbers
Gy.Terdik, Multivariate statistical methods - Going beyond the linear, Springer 2021, (5.36) p. 266, (see p.247 for Delta)
Other Random generation:
rCFUSN()
,
rSkewNorm()
,
rUniS()
n <- 10^3; d <- 2; p <- 3 ; a <- 1; b <- 1 Lamd <- matrix(sample(1:50-25, d*p), nrow=d) ieg<- eigen(diag(p)+t(Lamd)%*%Lamd) V <- ieg$vectors Delta <-Lamd %*% V %*% diag(1/sqrt(ieg$values)) %*% t(V) rCFUSSD(20,d,p,1,1,Delta)
n <- 10^3; d <- 2; p <- 3 ; a <- 1; b <- 1 Lamd <- matrix(sample(1:50-25, d*p), nrow=d) ieg<- eigen(diag(p)+t(Lamd)%*%Lamd) V <- ieg$vectors Delta <-Lamd %*% V %*% diag(1/sqrt(ieg$values)) %*% t(V) rCFUSSD(20,d,p,1,1,Delta)
Generate random d-vectors from the multivariate Skew Normal distribution
rSkewNorm(n, omega, alpha)
rSkewNorm(n, omega, alpha)
n |
sample size |
omega |
correlation matrix with d dimension |
alpha |
shape parameter vector of dimension d |
A random matrix
Azzalini, A. with the collaboration of Capitanio, A. (2014). The Skew-Normal and Related Families. Cambridge University Press, IMS Monographs series.
Gy.H.Terdik, Multivariate statistical methods - Going beyond the linear, Springer 2021, Section 5.1.2
Other Random generation:
rCFUSN()
,
rCFUSSD()
,
rUniS()
alpha<-c(10,5,0) omega<-diag(3) x<-rSkewNorm(20,omega,alpha)
alpha<-c(10,5,0) omega<-diag(3) x<-rSkewNorm(20,omega,alpha)
Generate random d-vectors from the Uniform distribution on the sphere
rUniS(n, d)
rUniS(n, d)
n |
sample size |
d |
dimension |
A random matrix
Gy.Terdik, Multivariate statistical methods - Going beyond the linear, Springer 2021
S. R. Jammalamadaka, E. Taufer, Gy. Terdik. On multivariate skewness and kurtosis. Sankhya A, 83(2), 607-644.
Other Random generation:
rCFUSN()
,
rCFUSSD()
,
rSkewNorm()
Provides estimates of mean, variance, skewness and kurtosis vectors for d-variate data
SampleEVSK(X)
SampleEVSK(X)
X |
d-variate data vector |
The list of the estimated mean, variance, skewness and kurtosis vectors
Gy.Terdik, Multivariate statistical methods - going beyond the linear, Springer 2021, Sections 6.4.1 and 6.5.1
Other Estimation:
SampleKurt()
,
SampleMomCum()
,
SampleSkew()
,
VarianceKurt()
,
VarianceSkew()
x<- MASS::mvrnorm(100,rep(0,3), 3*diag(rep(1,3))) EVSK<-SampleEVSK(x) names(EVSK) EVSK$estSkew
x<- MASS::mvrnorm(100,rep(0,3), 3*diag(rep(1,3))) EVSK<-SampleEVSK(x) names(EVSK) EVSK$estSkew
Provides the truncated Gram-Charlier approximation to a multivariate density. Approximation can be up to the first k=8 cumulants.
SampleGC(X, k = 4, cum = NULL)
SampleGC(X, k = 4, cum = NULL)
X |
A matrix of d-variate data |
k |
the order of the approximation, by default set to 4; (k must not be smaller than 3 or greater than 8) |
cum |
if NULL (default) the cumulant vector is estimated from X.
If |
The vector of the Gram-Charlier density evaluated at X
Gy.Terdik, Multivariate statistical methods - Going beyond the linear, Springer 2021. Section 4.7.
# Gram-Charlier density approximation (k=4) of data generated from # a bivariate skew-gaussian distribution n<-50 alpha<-c(10,0) omega<-diag(2) X<-rSkewNorm(n,omega,alpha) EC<-SampleEVSK(X) fy4<-SampleGC(X[1:5,],cum=EC)
# Gram-Charlier density approximation (k=4) of data generated from # a bivariate skew-gaussian distribution n<-50 alpha<-c(10,0) omega<-diag(2) X<-rSkewNorm(n,omega,alpha) EC<-SampleEVSK(X) fy4<-SampleGC(X[1:5,],cum=EC)
The vector x is standardized and the N-th d-variate polynomial is computed
SampleHermiteN(x, N)
SampleHermiteN(x, N)
x |
a d-variate data vector |
N |
the order of the d-variate Hermite polynomial |
The vector of the N-th d-variate polynomial
x<-MASS::mvrnorm(100,rep(0,3),diag(3)) H3<-SampleHermiteN(x,3)
x<-MASS::mvrnorm(100,rep(0,3),diag(3)) H3<-SampleHermiteN(x,3)
Estimates the sample kurtosis index based on the specified method: Mardia, MRSz, or Total.
SampleKurt(x, Type = c("Mardia", "MRSz", "Total"))
SampleKurt(x, Type = c("Mardia", "MRSz", "Total"))
x |
A matrix of multivariate data. |
Type |
A character string specifying the type of kurtosis index to estimate. Use "Mardia" for Mardia's kurtosis index, "MRSz" for the Mori-Rohatgi-Szekely kurtosis matrix, or "Total" for the total kurtosis index. |
A list containing the estimated kurtosis index or matrix and the associated p-value under the Gaussian hypothesis.
Mardia.Kurtosis |
The kurtosis index when |
MRSz.Kurtosis |
The kurtosis matrix when |
Total.Kurtosis |
The total kurtosis index when |
p.value |
The p-value under the Gaussian hypothesis for the estimated kurtosis. |
Gy.Terdik, Multivariate statistical methods - going beyond the linear, Springer 2021. Example 6.1 and 6.9.
Other Estimation:
SampleEVSK()
,
SampleMomCum()
,
SampleSkew()
,
VarianceKurt()
,
VarianceSkew()
# Mardia's kurtosis example x <- matrix(rnorm(100*5), ncol=5) SampleKurt(x, Type = "Mardia") # MRSz's kurtosis example SampleKurt(x, Type = "MRSz") # Total kurtosis example SampleKurt(x, Type = "Total")
# Mardia's kurtosis example x <- matrix(rnorm(100*5), ncol=5) SampleKurt(x, Type = "Mardia") # MRSz's kurtosis example SampleKurt(x, Type = "MRSz") # Total kurtosis example SampleKurt(x, Type = "Total")
Provides estimates of univariate and multivariate moments and cumulants up to order r. By default data are standardized; using only demeaned or raw data is also possible.
SampleMomCum(X, r, centering = FALSE, scaling = TRUE)
SampleMomCum(X, r, centering = FALSE, scaling = TRUE)
X |
d-vector data |
r |
The highest moment order (r >2) |
centering |
set to T (and scaling = F) if only centering is needed |
scaling |
set to T (and centering=F) if standardization of multivariate data is needed |
estMu.r
: the list of the multivariate moments up to order r
estCum.r
: the list of the multivariate cumulants up to order r
Gy.Terdik, Multivariate statistical methods - going beyond the linear, Springer 2021.
Other Estimation:
SampleEVSK()
,
SampleKurt()
,
SampleSkew()
,
VarianceKurt()
,
VarianceSkew()
## generate random data from a 3-variate skew normal distribution alpha<-c(10,5,0) omega<-diag(3) x<-rSkewNorm(50,omega,alpha) ## estimate the first three moments and cumulants from raw (uncentered and unstandardized) data SampleMomCum(x,3,centering=FALSE,scaling=FALSE) ## estimate the first three moments and cumulants from standardized data SampleMomCum(x,3,centering=FALSE,scaling=TRUE)
## generate random data from a 3-variate skew normal distribution alpha<-c(10,5,0) omega<-diag(3) x<-rSkewNorm(50,omega,alpha) ## estimate the first three moments and cumulants from raw (uncentered and unstandardized) data SampleMomCum(x,3,centering=FALSE,scaling=FALSE) ## estimate the first three moments and cumulants from standardized data SampleMomCum(x,3,centering=FALSE,scaling=TRUE)
Estimates the sample skewness index based on the specified method: Mardia or MRSz.
SampleSkew(x, Type = c("Mardia", "MRSz"))
SampleSkew(x, Type = c("Mardia", "MRSz"))
x |
A matrix of multivariate data. |
Type |
A character string specifying the type of skewness index to estimate. Use "Mardia" for Mardia's skewness index or "MRSz" for the Mori-Rohatgi-Szekely skewness vector and index. |
A list containing the estimated skewness index or vector and the associated p-value under the Gaussian hypothesis.
Mardia.Skewness |
The skewness index when |
MRSz.Skewness.Vector |
The skewness vector when |
MRSz.Skewness.Index |
The skewness index when |
p.value |
The p-value under the Gaussian hypothesis for the estimated skewness. |
Gy.Terdik, Multivariate statistical methods - going beyond the linear, Springer 2021. Example 6.1 and 6.2.
S. R. Jammalamadaka, E. Taufer, Gy. Terdik. On multivariate skewness and kurtosis. Sankhya A, 83(2), 607-644.
Other Estimation:
SampleEVSK()
,
SampleKurt()
,
SampleMomCum()
,
VarianceKurt()
,
VarianceSkew()
# Mardia's skewness example x <- matrix(rnorm(100*5), ncol=5) SampleSkew(x, Type = "Mardia") # MRSz's skewness example SampleSkew(x, Type = "MRSz")
# Mardia's skewness example x <- matrix(rnorm(100*5), ncol=5) SampleSkew(x, Type = "Mardia") # MRSz's skewness example SampleSkew(x, Type = "MRSz")
Provides the estimated covariance matrices of the data-estimated skewness and kurtosis vectors.
SampleVarianceSkewKurt(X)
SampleVarianceSkewKurt(X)
X |
A matrix of d-variate data |
The list of covariance matrices of the skewness and kurtosis vectors
Gy.Terdik, Multivariate statistical methods - going beyond the linear, Springer 2021.
Vector symmetrizing a T-product of vectors of the same dimension d. Produces the same results as SymMatr
SymIndx(x, d, n)
SymIndx(x, d, n)
x |
the vector to be symmetrized of dimension d^n |
d |
size of the single vectors in the product |
n |
power of the T-product |
A vector with the symmetrized version of x of dimension d^n
Gy. Terdik, Multivariate statistical methods - going beyond the linear, Springer 2021.Section 1.3.1 Symmetrization, p.14. (1.29)
Other Matrices and commutators:
EliminIndx()
,
EliminMatr()
,
QplicIndx()
,
QplicMatr()
,
SymMatr()
,
UnivMomCum()
a<-c(1,2) b<-c(2,3) c<-kronecker(kronecker(a,a),b) ## The symmetrized version of c is SymIndx(c,2,3)
a<-c(1,2) b<-c(2,3) c<-kronecker(kronecker(a,a),b) ## The symmetrized version of c is SymIndx(c,2,3)
Based on Chacon and Duong (2015) efficient recursive algorithms for functionals based on higher order derivatives. An option for sparse matrix is provided. By using sparse matrices far less memory is required and faster computation times are obtained
SymMatr(d, n, useSparse = FALSE)
SymMatr(d, n, useSparse = FALSE)
d |
dimension of a vector x |
n |
power of the Kronecker product |
useSparse |
TRUE or FALSE. If TRUE an object of the class "dgCMatrix" is produced. |
A Symmetrizer matrix with order . If
useSparse=TRUE
an object of the class "dgCMatrix" is produced.
Chacon, J. E., and Duong, T. (2015). Efficient recursive algorithms for functionals based on higher order derivatives of the multivariate Gaussian density. Statistics and Computing, 25(5), 959-974.
Gy. Terdik, Multivariate statistical methods - going beyond the linear, Springer 2021.Section 1.3.1 Symmetrization, p.14. (1.29)
Other Matrices and commutators:
EliminIndx()
,
EliminMatr()
,
QplicIndx()
,
QplicMatr()
,
SymIndx()
,
UnivMomCum()
a<-c(1,2) b<-c(2,3) c<-kronecker(kronecker(a,a),b) ## The symmetrized version of c is as.vector(SymMatr(2,3)%*%c)
a<-c(1,2) b<-c(2,3) c<-kronecker(kronecker(a,a),b) ## The symmetrized version of c is as.vector(SymMatr(2,3)%*%c)
A vector of indexes to select the moments and cumulants of the single components of the random vector X for which a T-vector of moments and cumulants is available
UnivMomCum(d, q)
UnivMomCum(d, q)
d |
dimension of a vector X |
q |
power of the Kronecker product |
A vector of indexes
Other Matrices and commutators:
EliminIndx()
,
EliminMatr()
,
QplicIndx()
,
QplicMatr()
,
SymIndx()
,
SymMatr()
## For a 3-variate skewness and kurtosis vectors estimated from data, extract ## the skewness and kurtosis estimates for each of the single components of the vector alpha<-c(10,5,0) omega<-diag(rep(1,3)) X<-rSkewNorm(200, omega, alpha) EVSK<-SampleEVSK(X) ## Get the univariate skewness and kurtosis for X1,X2,X3 EVSK$estSkew[UnivMomCum(3,3)] EVSK$estKurt[UnivMomCum(3,4)]
## For a 3-variate skewness and kurtosis vectors estimated from data, extract ## the skewness and kurtosis estimates for each of the single components of the vector alpha<-c(10,5,0) omega<-diag(rep(1,3)) X<-rSkewNorm(200, omega, alpha) EVSK<-SampleEVSK(X) ## Get the univariate skewness and kurtosis for X1,X2,X3 EVSK$estSkew[UnivMomCum(3,3)] EVSK$estKurt[UnivMomCum(3,4)]
Warning: the function requires 8! computations, for d>3, the timing required maybe large.
VarianceKurt(cum)
VarianceKurt(cum)
cum |
The theoretical/estimated cumulants up to the 8th order in vector form |
The matrix of theoretical/estimated variance
Gy.Terdik, Multivariate statistical methods - going beyond the linear, Springer 2021. Ch. 6, formula (6.26)
Other Estimation:
SampleEVSK()
,
SampleKurt()
,
SampleMomCum()
,
SampleSkew()
,
VarianceSkew()
Asymptotic Variance of the estimated skewness vector
VarianceSkew(cum)
VarianceSkew(cum)
cum |
The theoretical/estimated cumulants up to order 6 in vector form |
The matrix of theoretical/estimated variance
Gy.Terdik, Multivariate statistical methods - Going beyond the linear, Springer 2021. Ch.6, formula (6.13)
Other Estimation:
SampleEVSK()
,
SampleKurt()
,
SampleMomCum()
,
SampleSkew()
,
VarianceKurt()
alpha<-c(10,5) omega<-diag(rep(1,2)) MC <- MomCumSkewNorm(r = 6,omega,alpha) cum <- MC$CumX VS <- VarianceSkew(cum)
alpha<-c(10,5) omega<-diag(rep(1,2)) MC <- MomCumSkewNorm(r = 6,omega,alpha) cum <- MC$CumX VS <- VarianceSkew(cum)