Package 'CholWishart'

Title: Cholesky Decomposition of the Wishart Distribution
Description: Sampling from the Cholesky factorization of a Wishart random variable, sampling from the inverse Wishart distribution, sampling from the Cholesky factorization of an inverse Wishart random variable, sampling from the pseudo Wishart distribution, sampling from the generalized inverse Wishart distribution, computing densities for the Wishart and inverse Wishart distributions, and computing the multivariate gamma and digamma functions. Provides a header file so the C functions can be called directly from other programs.
Authors: Geoffrey Thompson [aut, cre] , R Core Team [ctb]
Maintainer: Geoffrey Thompson <[email protected]>
License: GPL (>= 3)
Version: 1.1.4
Built: 2024-11-18 06:27:28 UTC
Source: CRAN

Help Index


Density for Random Wishart Distributed Matrices

Description

Compute the density of an observation of a random Wishart distributed matrix (dWishart) or an observation from the inverse Wishart distribution (dInvWishart).

Usage

dWishart(x, df, Sigma, log = TRUE)

dInvWishart(x, df, Sigma, log = TRUE)

Arguments

x

positive definite p×pp \times p observations for density estimation - either one matrix or a 3-D array.

df

numeric parameter, "degrees of freedom".

Sigma

positive definite p×pp \times p "scale" matrix, the matrix parameter of the distribution.

log

logical, whether to return value on the log scale.

Details

Note there are different ways of parameterizing the Inverse Wishart distribution, so check which one you need. Here, If XIWp(Σ,ν)X \sim IW_p(\Sigma, \nu) then X1Wp(Σ1,ν)X^{-1} \sim W_p(\Sigma^{-1}, \nu). Dawid (1981) has a different definition: if XWp(Σ1,ν)X \sim W_p(\Sigma^{-1}, \nu) and ν>p1\nu > p - 1, then X1=YIW(Σ,δ)X^{-1} = Y \sim IW(\Sigma, \delta), where δ=νp+1\delta = \nu - p + 1.

Value

Density or log of density

Functions

  • dInvWishart(): density for the inverse Wishart distribution.

References

Dawid, A. (1981). Some Matrix-Variate Distribution Theory: Notational Considerations and a Bayesian Application. Biometrika, 68(1), 265-274. doi:10.2307/2335827

Gupta, A. K. and D. K. Nagar (1999). Matrix variate distributions. Chapman and Hall.

Mardia, K. V., J. T. Kent, and J. M. Bibby (1979) Multivariate Analysis, London: Academic Press.

Examples

set.seed(20180222)
A <- rWishart(1, 10, diag(4))[, , 1]
A
dWishart(x = A, df = 10, Sigma = diag(4L), log = TRUE)
dInvWishart(x = solve(A), df = 10, Sigma = diag(4L), log = TRUE)

Multivariate Gamma Function

Description

A special mathematical function related to the gamma function, generalized for multivariate gammas. lmvgamma is the log of the multivariate gamma, mvgamma.

The multivariate gamma function for a dimension p is defined as:

Γp(a)=πp(p1)/4j=1pΓ[a+(1j)/2]\Gamma_{p}(a)=\pi^{p(p-1)/4}\prod_{j=1}^{p} \Gamma [a+(1-j)/2]

For p=1p = 1, this is the same as the usual gamma function.

Usage

lmvgamma(x, p)

mvgamma(x, p)

Arguments

x

non-negative numeric vector, matrix, or array

p

positive integer, dimension of a square matrix

Value

For lmvgamma log of multivariate gamma of dimension p for each entry of x. For non-log variant, use mvgamma.

Functions

  • mvgamma(): Multivariate gamma function.

References

A. K. Gupta and D. K. Nagar 1999. Matrix variate distributions. Chapman and Hall.

Multivariate gamma function. In Wikipedia, The Free Encyclopedia,from https://en.wikipedia.org/w/index.php?title=Multivariate_gamma_function

See Also

gamma and lgamma

Examples

lgamma(1:12)
lmvgamma(1:12, 1L)
mvgamma(1:12, 1L)
gamma(1:12)

Multivariate Digamma Function

Description

A special mathematical function related to the gamma function, generalized for multivariate distributions. The multivariate digamma function is the derivative of the log of the multivariate gamma function; for p=1p = 1 it is the same as the univariate digamma function.

ψp(a)=i=1pψ(a+(1i)/2)\psi_{p}(a)=\sum _{i=1}^{p}\psi(a+(1-i)/2)

where ψ\psi is the univariate digamma function (the derivative of the log-gamma function).

Usage

mvdigamma(x, p)

Arguments

x

non-negative numeric vector, matrix, or array

p

positive integer, dimension of a square matrix

Value

vector of values of multivariate digamma function.

References

A. K. Gupta and D. K. Nagar 1999. Matrix variate distributions. Chapman and Hall.

Multivariate gamma function. In Wikipedia, The Free Encyclopedia,from https://en.wikipedia.org/w/index.php?title=Multivariate_gamma_function

See Also

gamma, lgamma, digamma, and mvgamma

Examples

digamma(1:10)
mvdigamma(1:10, 1L)

Cholesky Factor of Random Wishart Distributed Matrices

Description

Generate n random matrices, distributed according to the Cholesky factorization of a Wishart distribution with parameters Sigma and df, Wp(Sigma,df)W_p(Sigma, df) (known as the Bartlett decomposition in the context of Wishart random matrices).

Usage

rCholWishart(n, df, Sigma)

Arguments

n

integer sample size.

df

numeric parameter, "degrees of freedom".

Sigma

positive definite p×pp \times p "scale" matrix, the matrix parameter of the distribution.

Value

a numeric array, say R, of dimension p×p×np \times p \times n, where each R[,,i] is a Cholesky decomposition of a sample from the Wishart distribution Wp(Sigma,df)W_p(Sigma, df). Based on a modification of the existing code for the rWishart function.

References

Anderson, T. W. (2003). An Introduction to Multivariate Statistical Analysis (3rd ed.). Hoboken, N. J.: Wiley Interscience.

Mardia, K. V., J. T. Kent, and J. M. Bibby (1979) Multivariate Analysis, London: Academic Press.

A. K. Gupta and D. K. Nagar 1999. Matrix variate distributions. Chapman and Hall.

See Also

rWishart, rInvCholWishart

Examples

# How it is parameterized:
set.seed(20180211)
A <- rCholWishart(1L, 10, 3 * diag(5L))[, , 1]
A
set.seed(20180211)
B <- rInvCholWishart(1L, 10, 1 / 3 * diag(5L))[, , 1]
B
crossprod(A) %*% crossprod(B)

set.seed(20180211)
C <- chol(stats::rWishart(1L, 10, 3 * diag(5L))[, , 1])
C

Random Generalized Inverse Wishart Distributed Matrices

Description

Generate n random matrices, distributed according to the generalized inverse Wishart distribution with parameters Sigma and df, Wp(Σ,df)W_p(\Sigma, df), with sample size df less than the dimension p.

Let XiX_i, i=1,2,...,dfi = 1, 2, ..., df be df observations of a multivariate normal distribution with mean 0 and covariance Sigma. Then XiXi\sum X_i X_i' is distributed as a pseudo Wishart Wp(Σ,df)W_p(\Sigma, df). Sometimes this is called a singular Wishart distribution, however, that can be confused with the case where Σ\Sigma itself is singular. Then the generalized inverse Wishart distribution is the natural extension of the inverse Wishart using the Moore-Penrose pseudo-inverse. This can generate samples for positive semi-definite Σ\Sigma however, a function dedicated to generating singular normal random distributions or singular pseudo Wishart distributions should be used if that is desired.

Note there are different ways of parameterizing the Inverse Wishart distribution, so check which one you need. Here, if XIWp(Σ,ν)X \sim IW_p(\Sigma, \nu) then X1Wp(Σ1,ν)X^{-1} \sim W_p(\Sigma^{-1}, \nu). Dawid (1981) has a different definition: if XWp(Σ1,ν)X \sim W_p(\Sigma^{-1}, \nu) and ν>p1\nu > p - 1, then X1=YIW(Σ,δ)X^{-1} = Y \sim IW(\Sigma, \delta), where δ=νp+1\delta = \nu - p + 1.

Usage

rGenInvWishart(n, df, Sigma)

Arguments

n

integer sample size.

df

integer parameter, "degrees of freedom", should be less than the dimension of p

Sigma

positive semi-definite p×pp \times p "scale" matrix, the matrix parameter of the distribution.

Value

a numeric array, say R, of dimension p×p×np \times p \times n, where each R[,,i] is a realization of the pseudo Wishart distribution Wp(Sigma,df)W_p(Sigma, df).

References

Diaz-Garcia, Jose A, Ramon Gutierrez Jaimez, and Kanti V Mardia. 1997. “Wishart and Pseudo-Wishart Distributions and Some Applications to Shape Theory.” Journal of Multivariate Analysis 63 (1): 73–87. doi:10.1006/jmva.1997.1689.

Bodnar, T., Mazur, S., Podgórski, K. "Singular inverse Wishart distribution and its application to portfolio theory", Journal of Multivariate Analysis, Volume 143, 2016, Pages 314-326, ISSN 0047-259X, doi:10.1016/j.jmva.2015.09.021.

Bodnar, T., Okhrin, Y., "Properties of the singular, inverse and generalized inverse partitioned Wishart distributions", Journal of Multivariate Analysis, Volume 99, Issue 10, 2008, Pages 2389-2405, ISSN 0047-259X, doi:10.1016/j.jmva.2008.02.024.

Uhlig, Harald. "On Singular Wishart and Singular Multivariate Beta Distributions." Ann. Statist. 22 (1994), no. 1, 395–405. doi:10.1214/aos/1176325375.

See Also

rWishart, rInvWishart, and rPseudoWishart

Examples

set.seed(20181228)
A <- rGenInvWishart(1L, 4L, 5.0 * diag(5L))[, , 1]
A
# A should be singular
eigen(A)$values
set.seed(20181228)
B <- rPseudoWishart(1L, 4L, 5.0 * diag(5L))[, , 1]

# A should be a Moore-Penrose pseudo-inverse of B
B
# this should be equal to B
B %*% A %*% B
# this should be equal to A
A %*% B %*% A

Cholesky Factor of Random Inverse Wishart Distributed Matrices

Description

Generate n random matrices, distributed according to the Cholesky factor of an inverse Wishart distribution with parameters Sigma and df, Wp(Sigma,df)W_p(Sigma, df).

Note there are different ways of parameterizing the Inverse Wishart distribution, so check which one you need. Here, if XIWp(Σ,ν)X \sim IW_p(\Sigma, \nu) then X1Wp(Σ1,ν)X^{-1} \sim W_p(\Sigma^{-1}, \nu). Dawid (1981) has a different definition: if XWp(Σ1,ν)X \sim W_p(\Sigma^{-1}, \nu) and ν>p1\nu > p - 1, then X1=YIW(Σ,δ)X^{-1} = Y \sim IW(\Sigma, \delta), where δ=νp+1\delta = \nu - p + 1.

Usage

rInvCholWishart(n, df, Sigma)

Arguments

n

integer sample size.

df

numeric parameter, "degrees of freedom".

Sigma

positive definite p×pp \times p "scale" matrix, the matrix parameter of the distribution.

Value

a numeric array, say R, of dimension p×p×np \times p \times n, where each R[,,i] is a Cholesky decomposition of a realization of the Wishart distribution Wp(Sigma,df)W_p(Sigma, df). Based on a modification of the existing code for the rWishart function

References

Anderson, T. W. (2003). An Introduction to Multivariate Statistical Analysis (3rd ed.). Hoboken, N. J.: Wiley Interscience.

Dawid, A. (1981). Some Matrix-Variate Distribution Theory: Notational Considerations and a Bayesian Application. Biometrika, 68(1), 265-274. doi:10.2307/2335827

Gupta, A. K. and D. K. Nagar (1999). Matrix variate distributions. Chapman and Hall.

Mardia, K. V., J. T. Kent, and J. M. Bibby (1979) Multivariate Analysis, London: Academic Press.

See Also

rWishart and rCholWishart

Examples

# How it is parameterized:
set.seed(20180211)
A <- rCholWishart(1L, 10, 3 * diag(5L))[, , 1]
A
set.seed(20180211)
B <- rInvCholWishart(1L, 10, 1 / 3 * diag(5L))[, , 1]
B
crossprod(A) %*% crossprod(B)

set.seed(20180211)
C <- chol(stats::rWishart(1L, 10, 3 * diag(5L))[, , 1])
C

Random Inverse Wishart Distributed Matrices

Description

Generate n random matrices, distributed according to the inverse Wishart distribution with parameters Sigma and df, Wp(Sigma,df)W_p(Sigma, df).

Note there are different ways of parameterizing the Inverse Wishart distribution, so check which one you need. Here, if XIWp(Σ,ν)X \sim IW_p(\Sigma, \nu) then X1Wp(Σ1,ν)X^{-1} \sim W_p(\Sigma^{-1}, \nu). Dawid (1981) has a different definition: if XWp(Σ1,ν)X \sim W_p(\Sigma^{-1}, \nu) and ν>p1\nu > p - 1, then X1=YIW(Σ,δ)X^{-1} = Y \sim IW(\Sigma, \delta), where δ=νp+1\delta = \nu - p + 1.

Usage

rInvWishart(n, df, Sigma)

Arguments

n

integer sample size.

df

numeric parameter, "degrees of freedom".

Sigma

positive definite p×pp \times p "scale" matrix, the matrix parameter of the distribution.

Value

a numeric array, say R, of dimension p×p×np \times p \times n, where each R[,,i] is a realization of the inverse Wishart distribution IWp(Sigma,df)IW_p(Sigma, df). Based on a modification of the existing code for the rWishart function.

References

Dawid, A. (1981). Some Matrix-Variate Distribution Theory: Notational Considerations and a Bayesian Application. Biometrika, 68(1), 265-274. doi:10.2307/2335827

Gupta, A. K. and D. K. Nagar (1999). Matrix variate distributions. Chapman and Hall.

Mardia, K. V., J. T. Kent, and J. M. Bibby (1979) Multivariate Analysis, London: Academic Press.

See Also

rWishart, rCholWishart, and rInvCholWishart

Examples

set.seed(20180221)
A <- rInvWishart(1L, 10, 5 * diag(5L))[, , 1]
set.seed(20180221)
B <- stats::rWishart(1L, 10, .2 * diag(5L))[, , 1]

A %*% B

Random Pseudo Wishart Distributed Matrices

Description

Generate n random matrices, distributed according to the pseudo Wishart distribution with parameters Sigma and df, Wp(Σ,df)W_p(\Sigma, df), with sample size df less than the dimension p.

Let XiX_i, i=1,2,...,dfi = 1, 2, ..., df be df observations of a multivariate normal distribution with mean 0 and covariance Sigma. Then XiXi\sum X_i X_i' is distributed as a pseudo Wishart Wp(Σ,df)W_p(\Sigma, df). Sometimes this is called a singular Wishart distribution, however, that can be confused with the case where Σ\Sigma itself is singular. If cases with a singular Σ\Sigma are desired, this function cannot provide them.

Usage

rPseudoWishart(n, df, Sigma)

Arguments

n

integer sample size.

df

integer parameter, "degrees of freedom", should be less than the dimension of p

Sigma

positive definite p×pp \times p "scale" matrix, the matrix parameter of the distribution.

Value

a numeric array, say R, of dimension p×p×np \times p \times n, where each R[,,i] is a realization of the pseudo Wishart distribution Wp(Sigma,df)W_p(Sigma, df).

References

Diaz-Garcia, Jose A, Ramon Gutierrez Jaimez, and Kanti V Mardia. 1997. “Wishart and Pseudo-Wishart Distributions and Some Applications to Shape Theory.” Journal of Multivariate Analysis 63 (1): 73–87. doi:10.1006/jmva.1997.1689.

Uhlig, Harald. "On Singular Wishart and Singular Multivariate Beta Distributions." Ann. Statist. 22 (1994), no. 1, 395–405. doi:10.1214/aos/1176325375.

See Also

rWishart, rInvWishart, and rGenInvWishart

Examples

set.seed(20181227)
A <- rPseudoWishart(1L, 4L, 5.0 * diag(5L))[, , 1]
# A should be singular
eigen(A)$values