Skip to content

Public API

This page lists all public docstrings exposed by the package.

Copulas.AMHCopula Type
julia
AMHGenerator{T}, AMHCopula{d, T}

Fields:

  • θ::Real - parameter

Constructors:

julia
AMHGenerator(θ)  # Constructs the generator. 
AMHCopula(d,θ)   # Construct the copula

The AMH Copula in dimension d is parameterized by θ ∈ [-1,1). It is an Archimedean copula with generator:

ϕ(t)=11θetθ.

Special cases:

  • When θ = 0, it collapses to independence.

References:

  • [3] Nelsen, Roger B. An introduction to copulas. Springer, 2006.
source
Copulas.ArchimaxCopula Type
julia
ArchimaxCopula{d, TG, TT}

Fields

  • gen::TG Archimedean generator ϕ (implements ϕ, ϕ⁻¹, derivatives)

  • tail::TT Extreme-value tail (implements Pickands A / STDF )

Constructor

julia
ArchimaxCopula(d, gen::Generator, tail::Tail)

Definition (bivariate shown). Let xi=ϕ1(ui) and denote the STDF by . The cdf is

C(u1,u2)=ϕ((x1,x2)).

Reductions

  • If (x)=x1+x2 (i.e., tail = NoTail()), this is the Archimedean copula with generator gen.

  • If ϕ(s)=es (i.e., gen = IndependentGenerator()), this is the extreme-value copula with tail tail.

params(::ArchimaxCopula) concatenates the parameter tuples of gen and tail.

References:

  • [55] Capéraà, Fougères & Genest (2000), Bivariate Distributions with Given Extreme Value Attractor.

  • [56] Charpentier, Fougères & Genest (2014), Multivariate Archimax Copulas.

  • [51] Mai, J. F., & Scherer, M. (2012). Simulating copulas: stochastic models, sampling algorithms, and applications.

source
Copulas.ArchimedeanCopula Type
julia
ArchimedeanCopula{d, TG}

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

Constructor:

julia
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 ϕ writes:

C(u)=ϕ1(i=1dϕ(ui))

The default sampling method is the Radial-simplex decomposition using the Williamson transformation of ϕ.

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 archimedean copula can be constructed as follows:

julia
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:

julia
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:

  • [28] Williamson, R. E. (1956). Multiply monotone functions and their Laplace transforms. Duke Math. J. 23 189–207. MR0077581

  • [14] McNeil, A. J., & Nešlehová, J. (2009). Multivariate Archimedean copulas, d-monotone functions and ℓ 1-norm symmetric distributions.

source
Copulas.AsymGalambosCopula Type
julia
AsymGalambosTail{T}, AsymGalambosCopula{T}

Fields:

  • α::Real — dependence parameter

  • θ₁::Real — asymmetry weight in [0,1]

  • θ₂::Real — asymmetry weight in [0,1]

Constructor

julia
AsymGalambosCopula(α, θ₁, θ₂)
ExtremeValueCopula(2, AsymGalambosTail(α, θ₁, θ₂))

The (bivariate) asymmetric Galambos extreme-value copula is parameterized by α[0,) and θ1,θ2[0,1]. It is an EV copula with Pickands function

A(t)=1((θ1t)α+(θ2(1t))α)1/α,t[0,1].

Special cases:

  • α = 0 ⇒ IndependentCopula

  • θ₁ = θ₂ = 0 ⇒ IndependentCopula

  • θ₁ = θ₂ = 1 ⇒ GalambosCopula

References:

  • [48] Families of min-stable multivariate exponential and multivariate extreme value distributions. Statist. Probab, 1990.
source
Copulas.AsymLogCopula Type
julia
AsymLogTail{T}, AsymLogCopula{T}

Fields:

  • α::Real — dependence parameter (α ≥ 1)

  • θ₁::Real — asymmetry weight in [0,1]

  • θ₂::Real — asymmetry weight in [0,1]

Constructor

julia
AsymLogCopula(α, θ₁, θ₂)
ExtremeValueCopula(2, AsymLogTail(α, θ₁, θ₂))

The (bivariate) asymmetric logistic extreme–value copula is parameterized by α ∈ [1, ∞) and θ₁, θ₂ ∈ [0,1]. Its Pickands dependence function is

A(t)=(θ1α(1t)α+θ2αtα)1/α+(θ1θ2)t+1θ1,t[0,1].

Special cases:

  • θ₁ = θ₂ = 1 ⇒ symmetric Logistic (Gumbel) copula

  • θ₁ = θ₂ = 0 ⇒ independence (A(t) ≡ 1)

References:

  • [49] : Tawn, Jonathan A. "Bivariate extreme value theory: models and estimation." Biometrika 75.3 (1988): 397-415.
source
Copulas.AsymMixedCopula Type
julia
AsymMixedTail{T}, AsymMixedCopula{T}

Fields:

  • θ₁::Real — parameter

  • θ₂::Real — parameter

Constructor

AsymMixedCopula(θ₁, θ₂) ExtremeValueCopula(2, AsymMixedTail(θ₁, θ₂))

The (bivariate) asymmetric Mixed extreme-value copula is parameterized by two parameters θ1, θ2 subject to the following constraints:

  • θ₁ ≥ 0

  • θ₁ + θ₂ ≤ 1

  • θ₁ + 2θ₂ ≤ 1

  • θ₁ + 3θ₂ ≥ 0

Its Pickands dependence function is

A(t)=θ2t3+θ1t2(θ1+θ2)t+1,t[0,1].

Special cases:

  • θ₁ = θ₂ = 0 ⇒ IndependentCopula

  • θ₂ = 0 ⇒ symmetric Mixed copula

References:

  • [49] : Tawn, Jonathan A. "Bivariate extreme value theory: models and estimation." Biometrika 75.3 (1988): 397-415.
source
Copulas.BB10Copula Type
julia
BB10Generator{T}, BB10Copula{d, T}

Fields:

  • θ::Real - parameter

  • δ::Real - parameter

Constructor

julia
BB10Generator(θ, δ)
BB10Copula(d, θ, δ)

The BB10 copula has parameters θ(0,) and δ[0,1]. It is an Archimedean copula with generator:

ϕ(t)=(1δetδ)1/θ,

References:

  • [4] Joe, H. (2014). Dependence modeling with copulas. CRC press, Page.206-207
source
Copulas.BB1Copula Type
julia
BB1Generator{T}, BB1Copula{d, T}

Fields:

  • θ::Real - parameter

  • δ::Real - parameter

Constructor

julia
BB1Generator(θ, δ)
BB1Copula(d, θ, δ)

The BB1 copula is parameterized by θ(0,) and δ[1,). It is an Archimedean copula with generator:

ϕ(t)=(1+t1/δ)1/θ,

Special cases:

  • When δ = 1, it is the ClaytonCopula with parameter θ.

References:

  • [4] Joe, H. (2014). Dependence modeling with copulas. CRC press, Page.190-192
source
Copulas.BB2Copula Type
julia
BB2Generator{T}, BB2Copula{d, T}

Fields:

  • θ::Real - parameter

  • δ::Real - parameter

Constructor

julia
BB2Generator(θ, δ)
BB2Copula(d, θ, δ)

The BB2 copula has parameters θ,δ(0,). It is an Archimedean copula with generator:

ϕ(t)=[1+δ1log(1+t)]1θ,

References:

  • [4] Joe, H. (2014). Dependence modeling with copulas. CRC press, Page.193-194
source
Copulas.BB3Copula Type
julia
BB3Generator{T}, BB3Copula{d, T}

Fields:

  • θ::Real - parameter

  • δ::Real - parameter

Constructor

julia
BB3Generator(θ, δ)
BB3Copula(d, θ, δ)

The BB3 copula has parameters θ[1,) and δ(0,). It is an Archimedean copula with generator:

ϕ(t)=exp([δ1log(1+t)]1θ),

References:

  • [4] Joe, H. (2014). Dependence modeling with copulas. CRC press, Page.195-196
source
Copulas.BB4Copula Type
julia
BB4Copula{T}

Fields: - θ::Real - dependence parameter (θ ≥ 0) - δ::Real - shape parameter (δ > 0)

Constructor

julia
BB4Copula(θ, δ)

The BB4 copula is a two-parameter Archimax copula constructed from the Galambos tail and the Clayton generator. Its distribution function is

C(u,v;θ,δ)=(uθ+vθ1[(uθ1)δ+(vθ1)δ]1/δ)1/θ,θ0,δ>0.

for 0u,v1.

Special cases:

  • As δ → 0+, reduces to the Clayton Copula (Archimedean).

  • As θ → 0+, reduces to the Galambos copula (extreme-value).

  • As θ → ∞ or δ → ∞, approaches the M Copula.

References:

  • [4] Joe, H. (2014). Dependence modeling with copulas. CRC press, Page.197-198
source
Copulas.BB5Copula Type
julia
BB5Copula{T}

Fields: - θ::Real - dependence parameter (θ ≥ 1) - δ::Real - shape parameter (δ > 0)

Constructor

julia
BB5Copula(θ, δ)

The BB5 copula is a two-parameter Archimax copula, constructed from the Galambos tail and the Gumbel generator. Its distribution function is

C(u,v;θ,δ)=exp{[xθ+yθ(xθδ+yθδ)1/δ]1/θ},θ1,δ>0,

where x=log(u) and y=log(v).

Special cases:

  • As δ → 0⁺, reduces to the Gumbel copula (extreme-value and Archimedean).

  • As θ = 1, reduces to the Galambos copula.

  • As θ → ∞ or δ → ∞, converges to the M copula.

References:

  • [4] Joe, H. (2014). Dependence modeling with copulas. CRC press, Page.197-198
source
Copulas.BB6Copula Type
julia
BB6Generator{T}, BB6Copula{d, T}

Fields:

  • θ::Real - parameter

  • δ::Real - parameter

Constructor

julia
BB6Generator(θ, δ)
BB6Copula(d, θ, δ)

The BB6 copula has parameters θ,δ[1,). It is an Archimedean copula with generator:

ϕ(t)=1[1exp(t1δ)]1θ

References:

  • [4] Joe, H. (2014). Dependence modeling with copulas. CRC press, Page.200-201
source
Copulas.BB7Copula Type
julia
BB7Generator{T}, BB7Copula{d, T}

Fields:

  • θ::Real - parameter

  • δ::Real - parameter

Constructor

julia
BB7Generator(θ, δ)
BB7Copula(d, θ, δ)

The BB7 copula is parameterized by θ[1,) and δ(0,). It is an Archimedean copula with generator:

ϕ(t)=1[1(1+t)1/δ]1/θ.

References:

  • [4] Joe, H. (2014). Dependence modeling with copulas. CRC press, Page.202-203
source
Copulas.BB8Copula Type
julia
BB8Generator{T}, BB8Copula{d, T}

Fields:

  • ϑ::Real - parameter

  • δ::Real - parameter

Constructor

julia
BB8Generator(ϑ, δ)
BB8Copula(d, ϑ, δ)

The BB8 copula has parameters ϑ[1,) and δ(0,1]. It is an Archimedean copula with generator:

ϕ(t)=δ1[1(1ηexp(t))1ϑ],

where η=1(1δ)ϑ.

References:

  • [4] Joe, H. (2014). Dependence modeling with copulas. CRC press, Page.204-205
source
Copulas.BB9Copula Type
julia
BB9Generator{T}, BB9Copula{d, T}

Fields:

  • ϑ::Real - parameter

  • δ::Real - parameter

Constructor

julia
BB9Generator(ϑ, δ)
BB9Copula(d, ϑ, δ)

The BB9 copula has parameters ϑ[1,) and δ(0,). It is an Archimedean copula with generator:

ϕ(t)=exp((δϑ+t)1ϑ+δ1),

References:

  • [4] Joe, H. (2014). Dependence modeling with copulas. CRC press, Page.205-206
source
Copulas.BC2Copula Type
julia
BC2Tail{T}, BC2Copula{T}

Fields:

  • a::Real — parameter (a ∈ [0,1])

  • b::Real — parameter (b ∈ [0,1])

Constructor

julia
BC2Copula(a, b)
ExtremeValueCopula(2, BC2Tail(a, b))

The bivariate BC2 extreme-value copula is parameterized by two parameters a,b[0,1]. Its Pickands dependence function is

A(t)=max{at,b(1t)}+max{(1a)t,(1b)(1t)},t[0,1].

References:

  • [50] Mai, J. F., & Scherer, M. (2011). Bivariate extreme-value copulas with discrete Pickands dependence measure. Extremes, 14, 311-324. Springer, 2011.
source
Copulas.BernsteinCopula Type
julia
BernsteinCopula{d}

Fields:

  • m::NTuple{d,Int} - polynomial degrees (smoothing parameters)

  • weights::Array{Float64, d} - precomputed grid of box measures

Constructor

julia
BernsteinCopula(C; m=10)
BernsteinCopula(data; m=10)

The Bernstein copula in dimension d is defined as

Bm(C)(u)=s1=0m1sd=0mdC(s1m1,,sdmd)j=1d(mjsj)ujsj(1uj)mjsj.

It is a polynomial approximation of the base copula C using the multivariate Bernstein operator.

Implementation notes:

  • The grid of box measures (weights) is fully precomputed and stored as an n-dimensional array at construction. This enables fast evaluation of the copula and its density, but can be memory-intensive for large d or m.

  • The choice of m controls the smoothness of the approximation: larger m yields finer approximation but exponentially increases memory and computation cost (jmj boxes).

  • For high dimensions or large m, memory usage may become prohibitive; see documentation for scaling behavior.

  • If C is an EmpiricalCopula, the constructor produces the empirical Bernstein copula, a smoothed version of the empirical copula.

  • Supports cdf, logpdf, and random generation via mixtures of beta distributions.

References:

  • [60] Sancetta, A., & Satchell, S. (2004). The Bernstein copula and its applications to modeling and approximations of multivariate distributions. Econometric Theory, 20(3), 535-562.

  • [59] Segers, J., Sibuya, M., & Tsukahara, H. (2017). The empirical beta copula. Journal of Multivariate Analysis, 155, 35-51.

source
Copulas.BetaCopula Type
julia
BetaCopula{d, MT}

Fields:

  • ranks::MT - ranks matrix (d × n), each row contains integers 1..n

Constructor

julia
BetaCopula(u)

The empirical beta copula in dimension d is defined as

Cnβ(u)=1ni=1nj=1dFn,Rij(uj),

where Rij is the rank of observation i in margin j, and UBeta(r,n+1r).

Notes:

  • This is always a valid copula for any finite sample size n.

  • Supports cdf, logpdf at observed points and random sampling.

References:

  • [59] Segers, J., Sibuya, M., & Tsukahara, H. (2017). The empirical beta copula. Journal of Multivariate Analysis, 155, 35-51.
source
Copulas.CheckerboardCopula Type
julia
CheckerboardCopula{d, T}

Fields:

  • m::Vector{Int} — length d; number of partitions per dimension (grid resolution).

  • boxes::Dict{NTuple{d,Int}, T} — dictionary-like mapping from grid box indices to empirical weights. Typically Dict{NTuple{d,Int}, Float64} built with StatsBase.proportionmap.

Constructor:

julia
CheckerboardCopula(X; m=nothing, pseudo_values=true)

Builds a piecewise-constant (histogram) copula on a regular grid. The unit cube in each dimension i is partitioned into m[i] equal bins. Each observation is assigned to a box k ∈ ∏_i {0, …, m[i]-1}; the empirical box weights w_k sum to 1. The copula density is constant inside each box, with

julia
c(u) = w_k × ∏_i m[i]  when u  box k,  and  0 otherwise.

The CDF admits the multilinear overlap form

julia
C(u) = ∑_k w_k × ∏_i clamp(m[i]·u_i  k_i, 0, 1),

which this type evaluates directly without storing all grid corners.

Notes:

  • If m is nothing, we use m = fill(n, d) where n = size(X, 2).

  • When pseudo_values=true (default), X must already be pseudo-observations in [0,1]. Otherwise pass raw data and set pseudo_values=false to convert via pseudos(X).

  • Each m[i] must divide n to produce a valid checkerboard on the sample grid; this is enforced by the constructor.

References

  • Neslehova (2007). On rank correlation measures for non-continuous random variables.

  • Durante, Sanchez & Sempi (2013) Multivariate patchwork copulas: a unified approach with applications to partial comonotonicity.

  • Segers, Sibuya & Tsukahara (2017). The empirical beta copula. J. Multivariate Analysis, 155, 35-51.

  • Genest, Neslehova & Rémillard (2017) Asymptotic behavior of the empirical multilinear copula process under broad conditions.

  • Cuberos, Masiello & Maume-Deschamps (2019) Copulas checker-type approximations: application to quantiles estimation of aggregated variables.

  • Fredricks & Hofert (2025). On the checkerboard copula and maximum entropy.

source
Copulas.ClaytonCopula Type
julia
ClaytonGenerator{T}, ClaytonCopula{d, T}

Fields:

  • θ::Real - parameter

Constructor

julia
ClaytonGenerator(θ)
ClaytonCopula(d, θ)

The Clayton copula in dimension d is parameterized by θ[1/(d1),) (with the independence case as the limit θ0). It is an Archimedean copula with generator

ϕ(t)=(1+θt)1/θ

with the continuous extension ϕ(t)=et at θ=0.

Special cases (for the copula in dimension d):

  • When θ=1/(d1), it is the WCopula (Lower Fréchet–Hoeffding bound)

  • When θ0, it is the IndependentCopula

  • When θ, it is the MCopula (Upper Fréchet–Hoeffding bound)

References:

  • [3] Nelsen, Roger B. An introduction to copulas. Springer, 2006.
source
Copulas.CopulaModel Type
julia
CopulaModel{CT, TM, TD} <: StatsBase.StatisticalModel

A fitted copula model.

This type stores the result of fitting a copula (or a Sklar distribution) to pseudo-observations or raw data, together with auxiliary information useful for statistical inference and model comparison.

Fields

  • result::CT — the fitted copula (or SklarDist).

  • n::Int — number of observations used in the fit.

  • ll::Float64 — log-likelihood at the optimum.

  • method::Symbol — fitting method used (e.g. :mle, :itau, :deheuvels).

  • vcov::Union{Nothing, AbstractMatrix} — estimated covariance of the parameters, if available.

  • converged::Bool — whether the optimizer reported convergence.

  • iterations::Int — number of iterations used in optimization.

  • elapsed_sec::Float64 — time spent in fitting.

  • method_details::NamedTuple — additional method-specific metadata (grid size, pseudo-values, etc.).

CopulaModel implements the standard StatsBase.StatisticalModel interface: StatsBase.nobs, StatsBase.coef, StatsBase.coefnames, StatsBase.vcov, StatsBase.aic, StatsBase.bic, StatsBase.deviance, etc.

See also Distributions.fit.

source
Copulas.CuadrasAugeCopula Type
julia
CuadrasAugeTail{T}, CuadrasAugeCopula{T}

Fields:

  • θ::Real — dependence parameter, θ ∈ [0,1]

Constructor

julia
CuadrasAugeCopula(θ)
ExtremeValueCopula(2, CuadrasAugeTail(θ))

The (bivariate) Cuadras-Augé extreme-value copula is parameterized by θ[0,1]. Its Pickands dependence function is

A(t)=max{t,1t}+(1θ)min{t,1t},t[0,1].

Special cases:

  • θ = 0 ⇒ IndependentCopula

  • θ = 1 ⇒ MCopula (comonotone copula)

References:

  • [51] Mai, J. F., & Scherer, M. (2012). Simulating copulas: stochastic models, sampling algorithms, and applications (Vol. 4). World Scientific.
source
Copulas.EmpiricalCopula Type
julia
EmpiricalCopula{d, MT}

Fields:

  • u::MT — pseudo-observation matrix of size (d, N).

Constructor

julia
EmpiricalCopula(u; pseudo_values=true)

The empirical copula in dimension d is defined from a matrix of pseudo-observations u=(ui,j)1id, 1jN with entries in [0,1]. Its distribution function is

C(x)=1Nj=1N1{u,jx},

where the inequality is componentwise. If pseudo_values=false, the constructor first ranks the raw data into pseudo-observations; otherwise it assumes u already contains pseudo-observations in [0,1].

Notes:

  • This is an empirical object based on pseudo-observations; it is not necessarily a true copula for finite N but is widely used for nonparametric inference.

  • Supports cdf, logpdf at observed points, random sampling, and subsetting.

References:

  • [3] Nelsen, Roger B. An introduction to copulas. Springer, 2006.
source
Copulas.EmpiricalEVCopula Type
julia
EmpiricalEVTail, EmpiricalEVCopula

Fields:

  • tgrid::Vector{Float64} — evaluation grid in (0,1)

  • Ahat::Vector{Float64} — estimated Pickands function values on tgrid

  • slope::Vector{Float64} — per-segment slopes for linear interpolation

Constructor

EmpiricalEVTail(u; method=:ols, grid=401, eps=1e-3, pseudo_values=true) ExtremeValueCopula(2, EmpiricalEVTail(u; ...))

The empirical extreme-value (EV) copula (bivariate) is defined from pseudo-observations u = (U₁, U₂) and a nonparametric estimator of the Pickands dependence function. Supported estimators are:

  • :pickands — classical Pickands estimator

  • :cfg — Capéraà–Fougères–Genest (CFG) estimator

  • :ols — OLS-intercept estimator

For stability, the estimated function is always projected onto the class of valid Pickands functions (convex, bounded between max(t,1-t) and 1, with endpoints fixed at 1).

Its Pickands function is

julia
Â(t),  t  (0,1),

evaluated via piecewise linear interpolation on the grid tgrid.

References

  • [caperaa1997nonparametric] Capéraà, Fougères, Genest (1997) Biometrika

  • [gudendorf2011nonparametric] Gudendorf, Segers (2011) Journal of Multivariate Analysis

source
Copulas.FGMCopula Type
julia
FGMCopula{d,T}

Fields:

  • θ::Real - parameter

Constructor

julia
FGMCopula(d, θ)

The multivariate Farlie–Gumbel–Morgenstern (FGM) copula of dimension d has 2dd1 parameters θ and

C(u)=i=1dui[1+k=2d1j1<<jkdθj1jku¯j1u¯jk],

where u¯=1u.

Special cases:

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

More details about Farlie-Gumbel-Morgenstern (FGM) copula are found in [3]. We use the stochastic representation from [73] to obtain random samples.

References:

  • [3] Nelsen, Roger B. An introduction to copulas. Springer, 2006.

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

source
Copulas.FrankCopula Type
julia
FrankGenerator{T}, FrankCopula{d, T}

Fields:

  • θ::Real - parameter

Constructor

julia
FrankGenerator(θ)
FrankCopula(d,θ)

The Frank copula in dimension d is parameterized by θ(,) (with independence as the limit θ0). It is an Archimedean copula with generator

ϕ(t)=1θlog(1(1eθ)et).

Special cases:

  • When θ, it is the WCopula (Lower Fréchet–Hoeffding bound)

  • When θ0, it is the IndependentCopula

  • When θ, it is the MCopula (Upper Fréchet–Hoeffding bound)

References:

  • [3] Nelsen, Roger B. An introduction to copulas. Springer, 2006.
source
Copulas.GalambosCopula Type
julia
GalambosTail{T}, GalambosCopula{T}

Fields:

  • θ::Real — dependence parameter, θ ≥ 0

Constructor

julia
GalambosCopula(θ)
ExtremeValueCopula(2, GalambosTail(θ))

The (bivariate) Galambos extreme-value copula is parameterized by θ[0,). Its Pickands dependence function is

A(t)=1(tθ+(1t)θ)1/θ,t(0,1).

Special cases:

  • θ = 0 ⇒ IndependentCopula

  • θ = ∞ ⇒ MCopula (upper Fréchet-Hoeffding bound)

References:

  • [52] Galambos, J. (1975). Order statistics of samples from multivariate distributions. Journal of the American Statistical Association, 70(351a), 674-680.
source
Copulas.GaussianCopula Type
julia
GaussianCopula{d, MT}

Fields:

  • Σ::MT — correlation matrix (the constructor coerces the input to a correlation matrix).

Constructors

julia
GaussianCopula(Σ)
GaussianCopula(d, ρ)
GaussianCopula(d::Integer, ρ::Real)

Where Σ is a (symmetric) covariance or correlation matrix. The two-argument form with (d, ρ) builds the equicorrelation matrix with ones on the diagonal and constant off-diagonal correlation ρ:

julia
Σ = fill(ρ, d, d); Σ[diagind(Σ)] .= 1
C = GaussianCopula(d, ρ)            # == GaussianCopula(Σ)

Validity domain (equicorrelated PD matrix): -1/(d-1) < ρ < 1. The boundary ρ = -1/(d-1) is singular and rejected. If ρ == 0, this returns IndependentCopula(d) (same fast-path as when passing a diagonal matrix).

The Gaussian copula is the copula of a multivariate normal distribution. It is defined by

C(x;Σ)=FΣ(FΣ,11(x1),,FΣ,d1(xd)),

where FΣ is the cdf of a centered multivariate normal with covariance/correlation Σ and FΣ,i its i-th marginal cdf.

Example usage:

julia
C = GaussianCopula(Σ)
u = rand(C, 1000)
pdf(C, u); cdf(C, u)
Ĉ = fit(GaussianCopula, u)

Special case:

  • If isdiag(Σ), the constructor returns IndependentCopula(d).

References:

  • [3] Nelsen, Roger B. An introduction to copulas. Springer, 2006.
source
Copulas.GumbelBarnettCopula Type
julia
GumbelBarnettGenerator{T}, GumbelBarnettCopula{d, T}

Fields:

  • θ::Real - parameter

Constructor

julia
GumbelBarnettGenerator(θ)
GumbelBarnettCopula(d,θ)

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

ϕ(t)=exp(θ1(1et)),0θ1.

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.GumbelCopula Type
julia
GumbelGenerator{T}, GumbelCopula{d, T}

Fields:

  • θ::Real - parameter

Constructor

julia
GumbelGenerator(θ)
GumbelCopula(d,θ)

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

ϕ(t)=exp(t1/θ).

It has a few special cases:

  • When θ = 1, it is the IndependentCopula

  • When θ → ∞, it is the MCopula (Upper Fréchet–Hoeffding bound)

References:

  • [3] Nelsen, Roger B. An introduction to copulas. Springer, 2006.
source
Copulas.HuslerReissCopula Type
julia
HuslerReissTail{T}, HuslerReissCopula{T}

Fields:

  • θ::Real — dependence parameter, θ ≥ 0

Constructor

julia
HuslerReissCopula(θ)
ExtremeValueCopula(2, HuslerReissTail(θ))

The (bivariate) Hüsler-Reiss extreme-value copula is parameterized by θ[0,). Its Pickands dependence function is

A(t)=tΦ(θ1+12θlog(t1t))+(1t)Φ(θ1+12θlog(1tt))

where Φ is the standard normal cdf.

Special cases:

  • θ = 0 ⇒ IndependentCopula

  • θ = ∞ ⇒ MCopula (upper Fréchet-Hoeffding bound)

References:

  • [53] Hüsler, J., & Reiss, R. D. (1989). Maxima of normal random vectors: between independence and complete dependence. Statistics & Probability Letters, 7(4), 283-286.
source
Copulas.IndependentCopula Type
julia
IndependentCopula(d)

The independent copula in dimension d has distribution function

C(x)=i=1dxi.

It is Archimedean with generator ψ(s)=es.

References:

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

Fields:

  • θ::Real - parameter

Constructor

julia
InvGaussianGenerator(θ)
InvGaussianCopula(d,θ)

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

ϕ(t)=exp(11+2θ2tθ).

More details about Inverse Gaussian Archimedean copula are found in :

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

Special cases:

  • When θ = 0, it is the IndependentCopula

References:

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

Fields:

  • θ::Real - parameter

Constructor

julia
JoeGenerator(θ)
JoeCopula(d,θ)

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

ϕ(t)=1(1et)1/θ.

It has a few special cases:

  • When θ = 1, it is the IndependentCopula

  • When θ = ∞, it is the MCopula (Upper Fréchet–Hoeffding bound)

References:

  • [3] Nelsen, Roger B. An introduction to copulas. Springer, 2006.
source
Copulas.LogCopula Type
julia
LogTail{T}, LogCopula{T}

Fields:

  • θ::Real — dependence parameter, θ ∈ [0,1]

Constructor

julia
LogCopula(θ)
ExtremeValueCopula(2, LogTail(θ))

The (bivariate) Mixed extreme-value copula is parameterized by θ[0,1]. Its Pickands dependence function is

A(t)=θt2θt+1,t[0,1].

Special cases:

  • θ = 0 ⇒ IndependentCopula

References:

  • [49] : Tawn, Jonathan A. "Bivariate extreme value theory: models and estimation." Biometrika 75.3 (1988): 397-415.
source
Copulas.MCopula Type
julia
MCopula(d)

The upper Fréchet–Hoeffding bound is the copula with the largest value among all copulas; it corresponds to comonotone random vectors. For any copula C and all u[0,1]d,

W(u)C(u)M(u).

Both Fréchet–Hoeffding bounds are Archimedean copulas.

References:

  • [3] Nelsen, Roger B. An introduction to copulas. Springer, 2006.
source
Copulas.MOCopula Type
julia
MOTail{T}, MOCopula{T}

Fields:

  • λ₁::Real — parameter ≥ 0

  • λ₂::Real — parameter ≥ 0

  • λ₁₂::Real — parameter ≥ 0

Constructor

julia
MOCopula(λ₁, λ₂, λ₁₂)
ExtremeValueCopula(2, MOTail(λ₁, λ₂, λ₁₂))

The (bivariate) Marshall-Olkin extreme-value copula is parameterized by λi[0,),i=1,2,{1,2} Its Pickands dependence function is

A(t)=λ1(1t)λ1+λ1,2+λ2tλ2+λ1,2+λ1,2max{1tλ1+λ1,2,tλ2+λ1,2}

Special cases:

  • If λ₁₂ = 0, reduces to an asymmetric independence-like form.

  • If λ₁ = λ₂ = 0, degenerates to complete dependence.

References:

  • [51] Mai, J. F., & Scherer, M. (2012). Simulating copulas: stochastic models, sampling algorithms, and applications (Vol. 4). World Scientific.
source
Copulas.MixedCopula Type
julia
MixedTail{T}, MixedCopula{T}

Fields:

  • θ::Real — dependence parameter, θ ∈ [0,1]

Constructor

julia
MixedCopula(θ)
ExtremeValueCopula(2, MixedTail(θ))

The (bivariate) Mixed extreme-value copula is parameterized by θ[0,1]. Its Pickands dependence function is

A(t)=θt2θt+1,t[0,1].

Special cases:

  • θ = 0 ⇒ IndependentCopula

References:

  • [49] : Tawn, Jonathan A. "Bivariate extreme value theory: models and estimation." Biometrika 75.3 (1988): 397-415.
source
Copulas.PlackettCopula Type
julia
PlackettCopula{P}

Fields: - θ::Real - parameter

Constructor

julia
PlackettCopula(θ)

The Plackett copula is parameterized by θ>0 and is defined by

Cθ(u,v)=[1+(θ1)(u+v)][1+(θ1)(u+v)]24uvθ(θ1)2(θ1)

and for θ=1 we have C1(u,v)=uv.

Special cases:

  • θ = 0: MCopula (upper Fréchet–Hoeffding bound)

  • θ = 1: IndependentCopula

  • θ = ∞: WCopula (lower Fréchet–Hoeffding bound)

References:

  • [4] Joe, H. (2014). Dependence modeling with copulas. CRC press, Page.164

  • [72] 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.RafteryCopula Type
julia
RafteryCopula{d, P}

Fields: - θ::Real - parameter

Constructor

julia
RafteryCopula(d, θ)

The multivariate Raftery copula of dimension d is parameterized by θ[0,1].

Cθ(u)=u(1)+(1θ)(1d)1θd(j=1duj)11θi=2dθ(1θ)(1θi)(2θi)(j=1i1u(j))11θu(i)2θi1θ

where u(1),,u(d) denote the order statistics of u1,,ud. More details about Multivariate Raftery Copula are found in the references below.

Special cases:

  • When θ = 0, it is the IndependentCopula.

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

References:

  • [74] 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.SklarDist Type
julia
SklarDist{CT,TplMargins}

Fields:

  • C::CT - The copula

  • m::TplMargins - a Tuple representing the marginal distributions

Constructor

julia
SklarDist(C,m)

Construct a joint distribution via Sklar's theorem from marginals and a copula. See Sklar's theorem:

Theorem: Sklar 1959

For every random vector X, there exists a copula C such that

xRd,F(x)=C(F1(x1),...,Fd(xd)). The copula C is uniquely determined on Ran(F1)×...×Ran(Fd), where Ran(Fi) denotes the range of the function Fi. In particular, if all marginals are absolutely continuous, C is unique.

The resulting random vector follows the Distributions.jl API (rand/cdf/pdf/logpdf). A fit method is also provided. Example:

julia
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:
= 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.SurvivalCopula Type
julia
SurvivalCopula(C, flips)
SurvivalCopula{d,CT,flips}

Construct the survival (flipped) version of a copula by flipping the arguments at the given indices.

Type-level encoding: The indices to flip are encoded at the type level as a tuple of integers, e.g. SurvivalCopula{4,ClaytonCopula,(2,3)}. This enables compile-time specialization and dispatch, and ensures that the flipping pattern is part of the type.

Sugar constructor: The ergonomic constructor SurvivalCopula(C, flips::Tuple) infers the type parameters from the arguments, so you can write:

julia
SurvivalCopula(ClaytonCopula(4, θ), (2,3))

which is equivalent to the explicit type form:

julia
SurvivalCopula{4,ClaytonCopula,(2,3)}(ClaytonCopula(4, θ))

For a copula C in dimension d and indices i₁, ..., iₖ ∈ 1:d, the survival copula flips the corresponding arguments:

S(u1,,ud)=C(v1,,vd),vj={1ujjflipsujotherwise

Notes:

  • In the bivariate case, this includes the usual 90/180/270-degree "rotations" of a copula family.

  • The resulting object is handled like the base copula: same API (cdf, pdf/logpdf, rand, fit) and uniform marginals in [0,1]d.

References:

  • [3] Nelsen (2006), An introduction to copulas.
source
Copulas.TCopula Type
julia
TCopula{d, df, MT}

Fields:

  • df::Int — degrees of freedom

  • Σ::MT — correlation matrix

Constructor

julia
TCopula(df, Σ)

The Student t copula is the copula of a multivariate Student t distribution. It is defined by

C(x;ν,Σ)=Fν,Σ(Fν,Σ,11(x1),,Fν,Σ,d1(xd)),

where Fν,Σ is the cdf of a centered multivariate t with correlation Σ and ν degrees of freedom.

Example usage:

julia
C = TCopula(2, Σ)
u = rand(C, 1000)
pdf(C, u); cdf(C, u)
Ĉ = fit(TCopula, u)

References:

  • [3] Nelsen, Roger B. An introduction to copulas. Springer, 2006.
source
Copulas.TiltedGenerator Type
julia
TiltedGenerator(G, p, sJ)

Archimedean generator tilted by conditioning on p components fixed at values with cumulative generator sum sJ = ∑ ϕ⁻¹(u_j). It defines

julia
ϕ_tilt(t) = ϕ^{(p)}(sJ + t) / ϕ^{(p)}(sJ)

and higher derivatives accordingly:

julia
ϕ_tilt^{(k)}(t) = ϕ^{(k+p)}(sJ + t) / ϕ^{(p)}(sJ)

which yields the conditional copula within the Archimedean family for the remaining d-p variables. You will get a TiltedGenerator if you condition() an archimedean copula.

source
Copulas.WCopula Type
julia
WCopula

The lower Fréchet–Hoeffding bound is the copula with the smallest value among all copulas. Note that W is a proper copula only when d=2; for d>2 it remains the pointwise lower bound but is not itself a copula. For any copula C and all u[0,1]d,

W(u)C(u)M(u).

Both Fréchet–Hoeffding bounds are Archimedean copulas.

References:

  • [3] Nelsen, Roger B. An introduction to copulas. Springer, 2006.
source
Copulas.WilliamsonGenerator Type
julia
WilliamsonGenerator{TX, d} (alias 𝒲{TX, d})

Fields:

  • X::TX – a random variable that represents its Williamson d-transform

The type parameter d::Int is the dimension of the transformation.

Constructor

julia
WilliamsonGenerator(X::Distributions.UnivariateDistribution, d)
𝒲(X::Distributions.UnivariateDistribution,d)
WilliamsonGenerator(atoms::AbstractVector, weights::AbstractVector, d)
𝒲(atoms::AbstractVector, weights::AbstractVector, d)

The WilliamsonGenerator (alias 𝒲) 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 fully generically in the package.

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

ϕ(t)=𝒲d(X)(t)=t(1tx)d1dF(x)=E((1tX)+d1)1t>0+(1F(0))1t<0

This function has several properties:

  • We have that ϕ(0)=1 and ϕ(Inf)=0

  • ϕ is d2 times derivable, and the signs of its derivatives alternates : k0,...,d2,(1)kϕ(k)0.

  • ϕ(d2) 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 ϕ can be accessed by

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

Note that you'll always have:

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

Special case (finite-support discrete X)

  • If X isa Distributions.DiscreteUnivariateDistribution and support(X) is finite, or if you pass directly atoms and weights to the constructor, the produced generator is piecewise-polynomial ϕ(t) = ∑_j w_j · (1 − t/r_j)_+^(d−1) matching the Williamson transform of a discrete radial law. It has specialized methods.

  • For infinite-support discrete distributions or when the support is not accessible as a finite iterable, the standard WilliamsonGenerator is constructed.

References:

  • [28] Williamson, R. E. (1956). Multiply monotone functions and their Laplace transforms. Duke Math. J. 23 189–207. MR0077581

  • [14] McNeil, Alexander J., and Johanna Nešlehová. "Multivariate Archimedean copulas, d-monotone functions and ℓ 1-norm symmetric distributions." (2009): 3059-3097.

source
Copulas.tEVCopula Type
julia
tEVTail{Tdf,Tρ}, tEVCopula{T}

Fields:

  • ν::Real — degrees of freedom (ν > 0)

  • ρ::Real — correlation parameter (ρ ∈ (-1,1])

Constructor

julia
tEVCopula(ν, ρ)
ExtremeValueCopula(2, tEVTail(ν, ρ))

The (bivariate) extreme-t copula is parameterized by ν>0 and \rho \in (-1,1]``. Its Pickands dependence function is

A(x)=xtν+1(Zx)+(1x)tν+1(Z1x)

Where tν+1 is the cumulative distribution function (CDF) of the standard t distribution with ν+1 degrees of freedom and

Zx=1+ν1ρ2((x1x)1/νρ)

Special cases:

  • ρ = 0 ⇒ IndependentCopula

  • ρ = 1 ⇒ M Copula (upper Fréchet-Hoeffding bound)

References:

  • [54] Nikoloulopoulos, A. K., Joe, H., & Li, H. (2009). Extreme value properties of multivariate t copulas. Extremes, 12, 129-148.
source
Copulas.EmpiricalGenerator Method
julia
EmpiricalGenerator(u::AbstractMatrix)

Nonparametric Archimedean generator fit via inversion of the empirical Kendall distribution.

This function returns a WilliamsonGenerator{TX, d} whose underlying distribution TX is a Distributions.DiscreteNonParametric, rather than a separate struct. The returned object still implements all optimized methods (ϕ, derivatives, inverses) via specialized dispatch on WilliamsonGenerator{<:DiscreteNonParametric}.

Usage

julia
G = EmpiricalGenerator(u)

where u::AbstractMatrix is a d×n matrix of observations (already on copula or pseudo scale).

Notes

  • The recovered discrete radial support is rescaled so its largest atom equals 1 (scale is not identifiable).

  • We keep the old documentation entry point for backward compatibility; existing code that relied on the EmpiricalGenerator type should instead treat the result as a Generator.

References

  • [41]

  • [42]

  • [36] Genest, Neslehova and Ziegel (2011), Inference in Multivariate Archimedean Copula Models

source
Copulas.condition Method
julia
    condition(C::Copula{D}, js, u_js)
    condition(X::SklarDist, js, x_js)

Construct conditional distributions with respect to a copula, either on the uniform scale (when passing a Copula) or on the original data scale (when passing a SklarDist).

Arguments

  • C::Copula{D}: D-variate copula

  • X::SklarDist: joint distribution with copula X.C and marginals X.m

  • js: indices of conditioned coordinates (tuple, NTuple, or vector)

  • u_js: values in [0,1] for U_js (when conditioning a copula)

  • x_js: values on original scale for X_js (when conditioning a SklarDist)

  • j, u_j, x_j: 1D convenience overloads for the common p = 1 case

Returns

  • If the number of remaining coordinates d = D - length(js) is 1:

    • condition(C, js, u_js) returns a Distortion on [0,1] describing U_i | U_js = u_js.

    • condition(X, js, x_js) returns an unconditional univariate distribution for X_i | X_js = x_js, computed as the push-forward D(X.m[i]) where D = condition(C, js, u_js) and u_js = cdf.(X.m[js], x_js).

  • If d > 1:

    • condition(C, js, u_js) returns the conditional joint distribution on the uniform scale as a SklarDist(ConditionalCopula, distortions).

    • condition(X, js, x_js) returns the conditional joint distribution on the original scale as a SklarDist with copula ConditionalCopula(C, js, u_js) and appropriately distorted marginals D_k(X.m[i_k]).

Notes

  • For best performance, pass js and u_js as NTuple to keep p = length(js) known at compile time. The specialized method condition(::Copula{2}, j, u_j) exploits this for the common D = 2, d = 1 case.

  • Specializations are provided for many copula families (Independent, Gaussian, t, Archimedean, several bivariate families). Others fall back to an automatic differentiation based construction.

  • This function returns the conditional joint distribution H_{I|J}(· | u_J). The “conditional copula” is ConditionalCopula(C, js, u_js), i.e., the copula of that conditional distribution.

source
Copulas.inverse_rosenblatt Method
julia
inverse_rosenblatt(C::Copula, u)

Computes the inverse rosenblatt transform associated to the copula C on the vector u. Formally, assuming that U ∼ Π, the independence copula, the result should be distributed as C. Also look at rosenblatt(C, u) for the inverse transformation. The interface proposes faster versions for matrix inputs u.

Generic inverse Rosenblatt using conditional distortions: U₁ = S₁, U_k = H_{k|1:(k-1)}^{-1}(S_k | U₁:U_{k-1}). Specialized families may provide faster overrides.

References:

  • [13] Rosenblatt, M. (1952). Remarks on a multivariate transformation. Annals of Mathematical Statistics, 23(3), 470-472.

  • [4] Joe, H. (2014). Dependence Modeling with Copulas. CRC Press. (Section 2.10)

  • [14] McNeil, A. J., & Nešlehová, J. (2009). Multivariate Archimedean copulas, d-monotone functions and ℓ 1-norm symmetric distributions.

source
Copulas.pseudos Method
julia
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#Ordinal_ranking_.28.221234.22_ranking.29, see StatsBase.ordinalrank for the ordering we use. If you want more flexibility, checkout NormalizeQuantiles.sampleranks.

source
Copulas.rosenblatt Method
julia
rosenblatt(C::Copula, u)

Computes the rosenblatt transform associated to the copula C on the vector u. Formally, assuming that U ∼ C, the result should be uniformely distributed on the unit hypercube. The importance of this transofrmation comes from its bijectivity: inverse_rosenblatt(C, rand(d)) is equivalent to rand(C). The interface proposes faster versions for matrix inputs u.

Generic Rosenblatt transform using conditional distortions: S₁ = U₁, S_k = H_{k|1:(k-1)}(U_k | U₁:U_{k-1}). Specialized families may provide faster overrides.

  • [13] Rosenblatt, M. (1952). Remarks on a multivariate transformation. Annals of Mathematical Statistics, 23(3), 470-472.

  • [4] Joe, H. (2014). Dependence Modeling with Copulas. CRC Press. (Section 2.10)

  • [14] McNeil, A. J., & Nešlehová, J. (2009). Multivariate Archimedean copulas, d-monotone functions and ℓ 1-norm symmetric distributions.

source
Copulas.subsetdims Method
julia
subsetdims(C::Copula, dims::NTuple{p, Int})
subsetdims(D::SklarDist, dims)

Return a new copula or Sklar distribution corresponding to the subset of dimensions specified by dims.

Arguments

  • C::Copula: The original copula object.

  • D::SklarDist: The original Sklar distribution.

  • dims::NTuple{p, Int}: Tuple of indices representing the dimensions to keep.

Returns

  • A SubsetCopula or a new SklarDist object corresponding to the selected dimensions. If p == 1, returns a Uniform distribution or the corresponding marginal.

Details

This function extracts the dependence structure among the specified dimensions from the original copula or Sklar distribution. Specialized methods exist for some copula types to ensure efficiency and correctness.

source

References

  1. R. B. Nelsen. An Introduction to Copulas. 2nd ed Edition, Springer Series in Statistics (Springer, New York, 2006).

  2. H. Joe. Dependence Modeling with Copulas (CRC press, 2014).

  3. A. Sklar. Fonctions de Repartition à n Dimension et Leurs Marges. Université Paris 8, 1–3 (1959).

  4. M. Rosenblatt. Remarks on a multivariate transformation. Annals of Mathematical Statistics 23, 470–472 (1952).

  5. 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).

  6. R. E. Williamson. On multiply monotone functions and their laplace transforms (Mathematics Division, Office of Scientific Research, US Air Force, 1955).

  7. C. Genest, J. Nešlehová and J. Ziegel. Inference in Multivariate Archimedean Copula Models. TEST 20, 223–256 (2011).

  8. A. J. McNeil and J. Nešlehová. Multivariate Archimedean copulas, d-monotone functions and 1-norm symmetric distributions. Annals of Statistics 37, 3059–3097 (2009).

  9. R. E. Williamson. Multiply monotone functions and their Laplace transforms. Duke Mathematical Journal 23, 189–207 (1956).

  10. H. Joe. Families of min-stable multivariate exponential and multivariate extreme value distributions. Statistics & probability letters 9, 75–81 (1990).

  11. J. A. Tawn. Bivariate extreme value theory: models and estimation. Biometrika 75, 397–415 (1988).

  12. J.-F. Mai and M. Scherer. Bivariate extreme-value copulas with discrete Pickands dependence measure. Extremes 14, 311–324 (2011).

  13. J.-F. Mai and M. Scherer. Simulating copulas: stochastic models, sampling algorithms, and applications. Vol. 4 (World Scientific, 2012).

  14. J. Galambos. Order statistics of samples from multivariate distributions. Journal of the American Statistical Association 70, 674–680 (1975).

  15. J. Hüsler and R.-D. Reiss. Maxima of normal random vectors: between independence and complete dependence. Statistics & Probability Letters 7, 283–286 (1989).

  16. A. K. Nikoloulopoulos, H. Joe and H. Li. Extreme value properties of multivariate t copulas. Extremes 12, 129–148 (2009).

  17. P. Capéraà, A.-L. Fougères and C. Genest. Bivariate distributions with given extreme value attractor. Journal of Multivariate Analysis 72, 30–49 (2000).

  18. A. Charpentier, A.-L. Fougères, C. Genest and J. Nešlehová. Multivariate archimax copulas. Journal of Multivariate Analysis 126, 118–136 (2014).

  19. J. Segers, M. Sibuya and H. Tsukahara. The empirical beta copula. Journal of Multivariate Analysis 155, 35–51 (2017).

  20. A. Sancetta and S. Satchell. The Bernstein copula and its applications to modeling and approximations of multivariate distributions. Econometric theory 20, 535–562 (2004).

  21. M. E. Johnson. Multivariate statistical simulation: A guide to selecting and generating continuous multivariate distributions. Vol. 192 (John Wiley & Sons, 1987).

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

  23. T. Saali, M. Mesfioui and A. Shabri. Multivariate extension of Raftery copula. Mathematics 11, 414 (2023).