Reference

All docstrings are copied here as a general reference.

Copulas.AMHGeneratorType
AMHGenerator{T}

Fields:

  • θ::Real - parameter

Constructor

AMHGenerator(θ)
AMHCopula(d,θ)

The AMH copula in dimension $d$ is parameterized by $\theta \in [-1,1)$. It is an Archimedean copula with generator :

\[\phi(t) = 1 - \frac{1-\theta}{e^{-t}-\theta}\]

It has a few special cases:

  • When θ = 0, it is the IndependentCopula

References:

  • [3] Nelsen, Roger B. An introduction to copulas. Springer, 2006.
source
Copulas.ArchimedeanCopulaType
ArchimedeanCopula{d, TG}

Fields: - G::TG : the generator <: Generator.

Constructor:

ArchimedeanCopula(d::Int,G::Generator)

For some Archimedean Generator G::Generator and some dimenson d, this class models the archimedean copula which has this generator. The constructor checks for validity by ensuring that max_monotony(G) ≥ d. The $d$-variate archimedean copula with generator $\phi$ writes:

\[C(\mathbf u) = \phi^{-1}\left(\sum_{i=1}^d \phi(u_i)\right)\]

The default sampling method is the Radial-simplex decomposition using the Williamson transformation of $\phi$.

There exists several known parametric generators that are implement in the package. For every NamedGenerator <: Generator implemented in the package, we provide a type alias `NamedCopula{d,...} = ArchimedeanCopula{d,NamedGenerator{...}} to be able to manipulate the classic archimedean copulas without too much hassle for known and usefull special cases.

A generic archimdean copula can be constructed as follows:

struct MyGenerator <: Generator end
ϕ(G::MyGenerator,t) = exp(-t) # your archimedean generator, can be any d-monotonous function.
max_monotony(G::MyGenerator) = Inf # could depend on generators parameters. 
C = ArchimedeanCopula(d,MyGenerator())

The obtained model can be used as follows:

spl = rand(C,1000)   # sampling
cdf(C,spl)           # cdf
pdf(C,spl)           # pdf
loglikelihood(C,spl) # llh

Bonus: If you know the Williamson d-transform of your generator and not your generator itself, you may take a look at WilliamsonGenerator that implements them. If you rather know the frailty distribution, take a look at WilliamsonFromFrailty.

References:

  • [21] Williamson, R. E. (1956). Multiply monotone functions and their Laplace transforms. Duke Math. J. 23 189–207. MR0077581
  • [20] McNeil, A. J., & Nešlehová, J. (2009). Multivariate Archimedean copulas, d-monotone functions and ℓ 1-norm symmetric distributions.
source
Copulas.ClaytonGeneratorType
ClaytonGenerator{T}

Fields:

  • θ::Real - parameter

Constructor

ClaytonGenerator(θ)
ClaytonCopula(d,θ)

The Clayton copula in dimension $d$ is parameterized by $\theta \in [-1/(d-1),\infty)$. It is an Archimedean copula with generator :

\[\phi(t) = \left(1+\mathrm{sign}(\theta)*t\right)^{-1\frac{1}{\theta}}\]

It has a few special cases:

  • When θ = -1/(d-1), it is the WCopula (Lower Frechet-Hoeffding bound)
  • When θ = 0, it is the IndependentCopula
  • When θ = ∞, is is the MCopula (Upper Frechet-Hoeffding bound)

References:

  • [3] Nelsen, Roger B. An introduction to copulas. Springer, 2006.
source
Copulas.EllipticalCopulaType
EllipticalCopula{d,MT}

This is an abstract type. It implements an interface for all Elliptical copulas. We construct internally elliptical copulas using the sklar's theorem, by considering the copula $C$ to be defined as :

\[C = F \circ (F_1^{-1},...,F_d^{-1}),\]

where $F$ and $F_1,...,F_d$ are respectively the multivariate distribution function of some elliptical random vector and the univariate distribution function of its marginals. For a type MyCop <: EllipitcalCopula, it is necessary to implement the following methods:

  • N(::Type{MyCOp}), returning the constructor of the elliptical random vector from its correlation matrix. For example, N(GaussianCopula) simply returns MvNormal from Distributions.jl.
  • U(::Type{MyCOp}), returning the constructor for the univariate marginal, usually in standardized form. For example, U(GaussianCopula) returns Normal from Distributions.jl.

From these two functions, the abstract type provides a fully functional copula.

Details

Recall the definition of spherical random vectors:

Definition (Spherical and elliptical random vectors): A random vector $\bm X$ is said to be spherical if for all orthogonal matrix $\bm A \in O_d(\mathbb R)$, $\bm A\bm X \sim \bm X$.

For every matrix $\bm B$ and vector $\bm c$, the random vector $\bm B \bm X + \bm c$ is then said to be elliptical.

Recall that spherical random vectors are random vectors which characteristic functions (c.f.) only depend on the norm of their arguments. Indeed, for any $\bm A \in O_d(\mathbb R)$,

\[\phi(\bm t) = \mathbb E\left(e^{\langle \bm t, \bm X \rangle}\right)= \mathbb E\left(e^{\langle \bm t, \bm A\bm X \rangle}\right) = \mathbb E\left(e^{\langle \bm A\bm t, \bm X \rangle}\right) = \phi(\bm A\bm t).\]

We can therefore express this characteristic function as $\phi(\bm t) = \psi(\lVert \bm t \rVert_2^2)$, where $\psi$ is a function that characterizes the spherical family, called the generator of the family. Any characteristic function that can be expressed as a function of the norm of its argument is the characteristic function of a spherical random vector, since $\lVert \bm A \bm t \rVert_2 = \lVert \bm t \rVert_2$ for any orthogonal matrix $\bm A$.

However, note that this is not how the underlying code is working, we do not check for validity of the proposed generator (we dont even use it). You can construct such an elliptical family using simply Sklar:

struct MyElliptical{d,T} <: EllipticalCopula{d,T}
    θ:T
end
U(::Type{MyElliptical{d,T}}) where {d,T} # Distribution of the univaraite marginals, Normal() for the Gaussian case. 
N(::Type{MyElliptical{d,T}}) where {d,T} # Distribution of the mutlivariate random vector, MvNormal(C.Σ) for the Gaussian case. 

These two functions are enough to implement the rest of the interface.

References:

  • [3] Nelsen, Roger B. An introduction to copulas. Springer, 2006.
source
Copulas.EmpiricalCopulaType
EmpiricalCopula{d,MT}

Fields:

  • u::MT - the matrix of observations.

Constructor

EmpiricalCopula(u;pseudos=true)

The EmpiricalCopula in dimension $d$ is parameterized by a pseudo-data matrix which should have shape (d,N). Its expression is given as :

\[C(\mathbf x) = \frac{1}{N}\sum_{i=1}^n \mathbf 1_{\mathbf u_i \le \mathbf x}\]

This function is very practical, be be aware that this is not a true copula (since $\mathbf u$ are only pseudo-observations). The constructor allows you to pass dirctly pseudo-observations (the default) or will compute them for you. You can then compute the cdf of the copula, and sample it through the standard interface.

References:

  • [3] Nelsen, Roger B. An introduction to copulas. Springer, 2006.
source
Copulas.FGMCopulaType
FGMCopula{d,T}

Fields:

  • θ::Real - parameter

Constructor

FGMCopula(d, θ)

The Multivariate Farlie-Gumbel-Morgenstern (FGM) copula of dimension d has $2^d-d-1$ parameters $\theta$ and function

\[C(\boldsymbol{u})=\prod_{i=1}^{d}u_i \left[1+ \sum_{k=2}^{d}\sum_{1 \leq j_1 < \cdots < j_k \leq d} \theta_{j_1 \cdots j_k} \bar{u}_{j_1}\cdots \bar{u}_{j_k} \right],\]

where $\bar{u}=1-u$.

More details about Farlie-Gumbel-Morgenstern (FGM) copula are found in :

Nelsen, Roger B. An introduction to copulas. Springer, 2006. Exercise 3.38.

We use the stochastic representation of the copula to obtain random samples.

Blier-Wong, C., Cossette, H., & Marceau, E. (2022). Stochastic representation of FGM copulas using multivariate Bernoulli random variables. Computational Statistics & Data Analysis, 173, 107506.

It has a few special cases:

  • When d=2 and θ = 0, it is the IndependentCopula.

References:

  • [3] Nelsen, Roger B. An introduction to copulas. Springer, 2006.
source
Copulas.FrankGeneratorType
FrankGenerator{T}

Fields:

  • θ::Real - parameter

Constructor

FrankGenerator(θ)
FrankCopula(d,θ)

The Frank copula in dimension $d$ is parameterized by $\theta \in [-\infty,\infty)$. It is an Archimedean copula with generator :

\[\phi(t) = -\frac{\log\left(1+e^{-t}(e^{-\theta-1})\right)}{ heta}\]

It has a few special cases:

  • When θ = -∞, it is the WCopula (Lower Frechet-Hoeffding bound)
  • When θ = 1, it is the IndependentCopula
  • When θ = ∞, is is the MCopula (Upper Frechet-Hoeffding bound)

References:

  • [3] Nelsen, Roger B. An introduction to copulas. Springer, 2006.
source
Copulas.GaussianCopulaType
GaussianCopula{d,MT}

Fields:

  • Σ::MT - covariance matrix

Constructor

GaussianCopula(Σ)

The Gaussian Copula is the copula of a Multivariate normal distribution. It is constructed as:

\[C(\mathbf{x}; \boldsymbol{\Sigma}) = F_{\Sigma}(F_{\Sigma,i}^{-1}(x_i),i\in 1,...d)\]

where $F_{\Sigma}$ is a cdf of a gaussian random vector and $F_{\Sigma,i}$ is the ith marginal cdf, while $\Sigma$ is the covariance matrix.

It can be constructed in Julia via:

C = GaussianCopula(Σ)

You can sample it, compute pdf and cdf, or even fit the distribution via:

u = rand(C,1000)
Random.rand!(C,u) # other calling syntax for rng.
pdf(C,u) # to get the density
cdf(C,u) # to get the distribution function 
Ĉ = fit(GaussianCopula,u) # to fit on the sampled data. 

GaussianCopulas have a special case:

  • When isdiag(Σ), the constructor returns an IndependentCopula(d)

References:

  • [3] Nelsen, Roger B. An introduction to copulas. Springer, 2006.
source
Copulas.GeneratorType
Generator

Abstract type. Implements the API for archimedean generators.

An Archimedean generator is simply a function $\phi :\mathbb R_+ \to [0,1]$ such that $\phi(0) = 1$ and $\phi(+\infty) = 0$.

To generate an archimedean copula in dimension $d$, the function also needs to be $d$-monotone, that is :

  • $\phi$ is $d-2$ times derivable.
  • $(-1)^k \phi^{(k)} \ge 0 \;\forall k \in \{1,..,d-2\},$ and if $(-1)^{d-2}\phi^{(d-2)}$ is a non-increasing and convex function.

The access to the function $\phi$ itself is done through the interface:

ϕ(G::Generator, t)

We do not check algorithmically that the proposed generators are d-monotonous. Instead, it is up to the person implementing the generator to tell the interface how big can $d$ be through the function

max_monotony(G::MyGenerator) = # some integer, the maximum d so that the generator is d-monotonous.

More methods can be implemented for performance, althouhg there are implement defaults in the package :

  • ϕ⁻¹( G::Generator, x) gives the inverse function of the generator.
  • ϕ⁽¹⁾(G::Generator, t) gives the first derivative.
  • ϕ⁽ᵏ⁾(G::Generator, k, t) gives the kth derivative.
  • williamson_dist(G::Generator, d) gives the Wiliamson d-transform of the generator, see WilliamsonTransforms.jl.

References:

  • [20] McNeil, A. J., & Nešlehová, J. (2009). Multivariate Archimedean copulas, d-monotone functions and ℓ 1-norm symmetric distributions.
source
Copulas.GumbelBarnettGeneratorType
GumbelBarnettGenerator{T}

Fields:

  • θ::Real - parameter

Constructor

GumbelBarnettGenerator(θ)
GumbelBarnettCopula(d,θ)

The Gumbel-Barnett copula is an archimdean copula with generator:

\[\phi(t) = \exp{θ^{-1}(1-e^{t})}, 0 \leq \theta \leq 1.\]

It has a few special cases:

  • When θ = 0, it is the IndependentCopula

References:

  • [4] Joe, H. (2014). Dependence modeling with copulas. CRC press, Page.437
  • [3] Nelsen, Roger B. An introduction to copulas. Springer, 2006.
source
Copulas.GumbelGeneratorType
GumbelGenerator{T}

Fields:

  • θ::Real - parameter

Constructor

GumbelGenerator(θ)
GumbelCopula(d,θ)

The Gumbel copula in dimension $d$ is parameterized by $\theta \in [1,\infty)$. It is an Archimedean copula with generator :

\[\phi(t) = \exp{-t^{\frac{1}{θ}}}\]

It has a few special cases:

  • When θ = 1, it is the IndependentCopula
  • When θ = ∞, is is the MCopula (Upper Frechet-Hoeffding bound)

References:

  • [3] Nelsen, Roger B. An introduction to copulas. Springer, 2006.
source
Copulas.IndependentGeneratorType
IndependentGenerator

Constructor

IndependentGenerator()
IndependentCopula(d)

The Independent Copula in dimension $d$ is the simplest copula, that has the form :

\[C(\mathbf{x}) = \prod_{i=1}^d x_i.\]

It happends to be an Archimedean Copula, with generator :

\[\phi(t) = \exp{-t}\]

References:

  • [3] Nelsen, Roger B. An introduction to copulas. Springer, 2006.
source
Copulas.InvGaussianGeneratorType
InvGaussianGenerator{T}

Fields:

  • θ::Real - parameter

Constructor

InvGaussianGenerator(θ)
InvGaussianCopula(d,θ)

The Inverse Gaussian copula in dimension $d$ is parameterized by $\theta \in [0,\infty)$. It is an Archimedean copula with generator :

\[\phi(t) = \exp{\frac{1-\sqrt{1+2θ^{2}t}}{θ}}.\]

More details about Inverse Gaussian Archimedean copula are found in :

Mai, Jan-Frederik, and Matthias Scherer. Simulating copulas: stochastic models, sampling algorithms, and applications. Vol. 6. # N/A, 2017. Page 74.

It has a few special cases:

  • When θ = 0, it is the IndependentCopula

References:

  • [3] Nelsen, Roger B. An introduction to copulas. Springer, 2006.
source
Copulas.JoeGeneratorType
JoeGenerator{T}

Fields:

  • θ::Real - parameter

Constructor

JoeGenerator(θ)
JoeCopula(d,θ)

The Joe copula in dimension $d$ is parameterized by $\theta \in [1,\infty)$. It is an Archimedean copula with generator :

\[\phi(t) = 1 - \left(1 - e^{-t}\right)^{\frac{1}{\theta}}\]

It has a few special cases:

  • When θ = 1, it is the IndependentCopula
  • When θ = ∞, is is the MCopula (Upper Frechet-Hoeffding bound)

References:

  • [3] Nelsen, Roger B. An introduction to copulas. Springer, 2006.
source
Copulas.MGeneratorType
MGenerator

Constructor

MGenerator()
MCopula(d)

The Upper Frechet-Hoeffding bound is the copula with the greatest value among all copulas. It correspond to comonotone random vectors.

For any copula $C$, if $W$ and $M$ are (respectively) the lower and uppder Frechet-Hoeffding bounds, we have that for all $\mathbf{u} \in [0,1]^d$,

\[W(\mathbf{u}) \le C(\mathbf{u}) \le M(\mathbf{u})\]

The two Frechet-Hoeffding bounds are also Archimedean copulas.

References:

  • [3] Nelsen, Roger B. An introduction to copulas. Springer, 2006.
source
Copulas.PlackettCopulaType
PlackettCopula{P}

Fields: - θ::Real - parameter

Constructor

PlackettCopula(θ)

Parameterized by $\theta > 0$ The Plackett copula is

\[C_{\theta}(u,v) = \frac{\left [1+(\theta-1)(u+v)\right]- \sqrt{[1+(\theta-1)(u+v)]^2-4uv\theta(\theta-1)}}{2(\theta-1)}\]

and for $\theta = 1$

\[C_{1}(u,v) = uv \]

It has a few special cases:

  • When θ = 0, is is the MCopula (Upper Frechet-Hoeffding bound)
  • When θ = 1, it is the IndependentCopula
  • When θ = ∞, is is the WCopula (Lower Frechet-Hoeffding bound)

References:

  • [4] Joe, H. (2014). Dependence modeling with copulas. CRC press, Page.164
  • [45] Johnson, Mark E. Multivariate statistical simulation: A guide to selecting and generating continuous multivariate distributions. Vol. 192. John Wiley & Sons, 1987. Page 193.
  • [3] Nelsen, Roger B. An introduction to copulas. Springer, 2006. Exercise 3.38.
source
Copulas.RafteryCopulaType
RafteryCopula{d, P}

Fields: - θ::Real - parameter

Constructor

RafteryCopula(d, θ)

The Multivariate Raftery Copula of dimension d is parameterized by $\theta \in [0,1]$

\[C_{\theta}(\mathbf{u}) = u_{(1)} + \frac{(1 - \theta)(1 - d)}{1 - \theta - d} \left(\prod_{j=1}^{d} u_j\right)^{\frac{1}{1-\theta}} - \sum_{i=2}^{d} \frac{\theta(1-\theta)}{(1-\theta-i)(2-\theta-i)} \left(\prod_{j=1}^{i-1}u_{(j)}\right)^{\frac{1}{1-\theta}}u_{(i)}^{\frac{2-\theta-i}{1-\theta}}\]

where $u_{(1)}, \ldots , u_{(d)}$ denote the order statistics of $u_1, \ldots ,u_d$. More details about Multivariate Raftery Copula are found in the references below.

It has a few special cases:

  • When θ = 0, it is the IndependentCopula.
  • When θ = 1, it is the the Fréchet upper bound

References:

  • [46] Saali, T., M. Mesfioui, and A. Shabri, 2023: Multivariate Extension of Raftery Copula. Mathematics, 11, 414, https://doi.org/10.3390/math11020414.
  • [3] Nelsen, Roger B. An introduction to copulas. Springer, 2006. Exercise 3.6.
source
Copulas.SklarDistType
SklarDist{CT,TplMargins}

Fields:

  • C::CT - The copula
  • m::TplMargins - a Tuple representing the marginal distributions

Constructor

SklarDist(C,m)

This function allows to construct a random vector specified, through the Sklar Theorem, by its marginals and its copula separately. See Sklar's theorem:

Theorem (Sklar 1959): For every random vector $\bm X$, there exists a copula $C$ such that

$\forall \bm x\in \mathbb R^d, F(\bm x) = C(F_{1}(x_{1}),...,F_{d}(x_{d})).$ The copula $C$ is uniquely determined on $\mathrm{Ran}(F_{1}) \times ... \times \mathrm{Ran}(F_{d})$, where $\mathrm{Ran}(F_i)$ denotes the range of the function $F_i$. In particular, if all marginals are absolutely continuous, $C$ is unique.

The obtain random vector follows Distributions.jl's API and can be sampled, pdf and cdf can be evaluated, etc... We even provide a fit function. See the folowing exemple code :

using Copulas, Distributions, Random
X₁ = Gamma(2,3)
X₂ = Pareto()
X₃ = LogNormal(0,1)
C = ClaytonCopula(3,0.7) # A 3-variate Clayton Copula with θ = 0.7
D = SklarDist(C,(X₁,X₂,X₃)) # The final distribution

simu = rand(D,1000) # Generate a dataset

# You may estimate a copula using the `fit` function:
D̂ = fit(SklarDist{ClaytonCopula,Tuple{Gamma,Normal,LogNormal}}, simu)

References:

  • [9] Sklar, M. (1959). Fonctions de répartition à n dimensions et leurs marges. In Annales de l'ISUP (Vol. 8, No. 3, pp. 229-231).
  • [3] Nelsen, Roger B. An introduction to copulas. Springer, 2006.
source
Copulas.SubsetCopulaType
SubsetCopula{d,CT}

Fields:

  • C::CT - The copula
  • dims::Tuple{Int64} - a Tuple representing which dimensions are used.

Constructor

SubsetCopula(C::Copula,dims)

This class allows to construct a random vector corresponding to a few dimensions of the starting copula. If $(X_1,...,X_n)$ is the random vector corresponding to the copula C, this returns the copula of ( $X_i$ for i in dims). The dependence structure is preserved. There are specialized methods for some copulas.

source
Copulas.SurvivalCopulaType
SurvivalCopula(C,indices)

Computes the survival version of any copula on given indices. From a copula $C$ in dimension $d$, and some indices $i_1,...i_k$ in ${1,...,d}$, the survival copula associated simply reverses its arguments on chosen indices. For exemple, for $d=4$ and indices $(2,3)$, we have:

\[S(u_1,...u_4) = C(u_1,1-u_2,1-u3,u_4)\]

This constructor allows to derive new "survival" families. For exemple, in bivariate cases, this allows to do "rotations". The obtained models can be treated as the starting one, i.e. as a random vector in [0,1]^d with uniforms marginals.

References:

  • [3] Nelsen, Roger B. An introduction to copulas. Springer, 2006.
source
Copulas.TCopulaType
TCopula{d,MT}

Fields:

  • df::Int - number of degree of freedom
  • Σ::MT - covariance matrix

Constructor

TCopula(df,Σ)

The Student's T Copula is the copula of a Multivariate Student distribution. It is constructed as :

\[C(\mathbf{x}; \boldsymbol{n,\Sigma}) = F_{n,\Sigma}(F_{n,\Sigma,i}^{-1}(x_i),i\in 1,...d)\]

where $F_{n,\Sigma}$ is a cdf of a multivariate student random vector with covariance matrix $\Sigma$ and $n$ degrees of freedom. and F_{n,\Sigma,i} is the ith marignal cdf.

It can be constructed in Julia via:

C = TCopula(2,Σ)

You can sample it, compute pdf and cdf, or even fit the distribution via:

u = rand(C,1000)
Random.rand!(C,u) # other calling syntax for rng.
pdf(C,u) # to get the density
cdf(C,u) # to get the distribution function 
Ĉ = fit(TCopula,u) # to fit on the sampled data. 

Except that currently it does not work since fit(Distributions.MvTDist,data) does not dispatch.

References:

  • [3] Nelsen, Roger B. An introduction to copulas. Springer, 2006.
source
Copulas.WGeneratorType
WGenerator

Constructor

WGenerator()
WCopula(d)

The Lower Frechet-Hoeffding bound is the copula with the lowest value among all copulas. Note that $W$ is only a proper copula when $d=2$, in greater dimensions it is still the (pointwise) lower bound, but not a copula anymore. For any copula $C$, if $W$ and $M$ are (respectively) the lower and uppder Frechet-Hoeffding bounds, we have that for all $\mathbf{u} \in [0,1]^d$,

\[W(\mathbf{u}) \le C(\mathbf{u}) \le M(\mathbf{u})\]

The two Frechet-Hoeffding bounds are also Archimedean copulas.

References:

  • [3] Nelsen, Roger B. An introduction to copulas. Springer, 2006.
source
Copulas.WilliamsonGeneratorType
WilliamsonGenerator{TX}
i𝒲{TX}

Fields:

  • X::TX – a random variable that represents its Williamson d-transform
  • d::Int – the dimension of the transformation.

Constructor

WilliamsonGenerator(X::Distributions.UnivariateDistribution, d)
i𝒲(X::Distributions.UnivariateDistribution,d)

The WilliamsonGenerator (alias i𝒲) allows to construct a d-monotonous archimedean generator from a positive random variable X::Distributions.UnivariateDistribution. The transformation, which is called the inverse Williamson transformation, is implemented in WilliamsonTransforms.jl.

For a univariate non-negative random variable $X$, with cumulative distribution function $F$ and an integer $d\ge 2$, the Williamson-d-transform of $X$ is the real function supported on $[0,\infty[$ given by:

\[\phi(t) = 𝒲_{d}(X)(t) = \int_{t}^{\infty} \left(1 - \frac{t}{x}\right)^{d-1} dF(x) = \mathbb E\left( (1 - \frac{t}{X})^{d-1}_+\right) \mathbb 1_{t > 0} + \left(1 - F(0)\right)\mathbb 1_{t <0}\]

This function has several properties:

  • We have that $\phi(0) = 1$ and $\phi(Inf) = 0$
  • $\phi$ is $d-2$ times derivable, and the signs of its derivatives alternates : $\forall k \in 0,...,d-2, (-1)^k \phi^{(k)} \ge 0$.
  • $\phi^{(d-2)}$ is convex.

These properties makes this function what is called a d-monotone archimedean generator, able to generate archimedean copulas in dimensions up to $d$. Our implementation provides this through the Generator interface: the function $\phi$ can be accessed by

G = WilliamsonGenerator(X, d)
ϕ(G,t)

Note that you'll always have:

max_monotony(WilliamsonGenerator(X,d)) === d

References:

  • [21] Williamson, R. E. (1956). Multiply monotone functions and their Laplace transforms. Duke Math. J. 23 189–207. MR0077581
  • [20] McNeil, Alexander J., and Johanna Nešlehová. "Multivariate Archimedean copulas, d-monotone functions and ℓ 1-norm symmetric distributions." (2009): 3059-3097.
source
Copulas.pseudosMethod
pseudos(sample)

Compute the pseudo-observations of a multivariate sample. Note that the sample has to be given in wide format (d,n), where d is the dimension and n the number of observations.

Warning: the order used is ordinal ranking like https://en.wikipedia.org/wiki/Ranking#Ordinalranking.28.221234.22_ranking.29, see StatsBase.ordinalrank for the ordering we use. If you want more flexibility, checkout NormalizeQuantiles.sampleranks.

source
Copulas.subsetdimsMethod
subsetdims(C::Copula,dims)
subsetdims(D::SklarDist, dims)

If $(X_1,...,X_n)$ is the random vector corresponding to the model C or D, this returns the distribution on ( $X_i$ for i in dims), preserving the dependence structure between the dimensions in dims. There are specialized methods for some copulas.

source
[3]
R. B. Nelsen. An Introduction to Copulas. 2nd ed Edition, Springer Series in Statistics (Springer, New York, 2006).
[4]
H. Joe. Dependence Modeling with Copulas (CRC press, 2014).
[9]
A. Sklar. Fonctions de Repartition à n Dimension et Leurs Marges. Université Paris 8, 1–3 (1959).
[20]
A. J. McNeil and J. Nešlehová. Multivariate Archimedean Copulas, d -Monotone Functions and L1 -Norm Symmetric Distributions. The Annals of Statistics 37, 3059–3097 (2009).
[21]
R. E. Williamson. On multiply monotone functions and their laplace transforms (Mathematics Division, Office of Scientific Research, US Air Force, 1955).
[45]
M. E. Johnson. Multivariate statistical simulation: A guide to selecting and generating continuous multivariate distributions. Vol. 192 (John Wiley & Sons, 1987).
[46]
T. Saali, M. Mesfioui and A. Shabri. Multivariate extension of Raftery copula. Mathematics 11, 414 (2023).