Public API
This page lists all public docstrings exposed by the package.
Copulas.AMHCopula Type
AMHGenerator{T}, AMHCopula{d, T}Fields:
- θ::Real - parameter
Constructors:
AMHGenerator(θ) # Constructs the generator.
AMHCopula(d,θ) # Construct the copulaThe AMH Copula in dimension d is parameterized by θ ∈ [-1,1). It is an Archimedean copula with generator:
Special cases:
- When θ = 0, it collapses to independence.
References:
- [3] Nelsen, Roger B. An introduction to copulas. Springer, 2006.
Copulas.ArchimaxCopula Type
ArchimaxCopula{d, TG, TT}Fields
gen::TGArchimedean generator(implements ϕ,ϕ⁻¹, derivatives)tail::TTExtreme-value tail (implements PickandsA/ STDFℓ)
Constructor
ArchimaxCopula(d, gen::Generator, tail::Tail)Definition (bivariate shown). Let
Reductions
If
(i.e., tail = NoTail()), this is the Archimedean copula with generatorgen.If
(i.e., gen = IndependentGenerator()), this is the extreme-value copula with tailtail.
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.
Copulas.ArchimedeanCopula Type
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
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:
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) # llhBonus: 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.
Copulas.AsymGalambosCopula Type
AsymGalambosTail{T}, AsymGalambosCopula{T}Fields:
α::Real — dependence parameter
θ₁::Real — asymmetry weight in [0,1]
θ₂::Real — asymmetry weight in [0,1]
Constructor
AsymGalambosCopula(α, θ₁, θ₂)
ExtremeValueCopula(2, AsymGalambosTail(α, θ₁, θ₂))The (bivariate) asymmetric Galambos extreme-value copula is parameterized by
Special cases:
α = 0 ⇒ IndependentCopula
θ₁ = θ₂ = 0 ⇒ IndependentCopula
θ₁ = θ₂ = 1 ⇒ GalambosCopula
References:
- [48] Families of min-stable multivariate exponential and multivariate extreme value distributions. Statist. Probab, 1990.
Copulas.AsymLogCopula Type
AsymLogTail{T}, AsymLogCopula{T}Fields:
α::Real — dependence parameter (α ≥ 1)
θ₁::Real — asymmetry weight in [0,1]
θ₂::Real — asymmetry weight in [0,1]
Constructor
AsymLogCopula(α, θ₁, θ₂)
ExtremeValueCopula(2, AsymLogTail(α, θ₁, θ₂))The (bivariate) asymmetric logistic extreme–value copula is parameterized by α ∈ [1, ∞) and θ₁, θ₂ ∈ [0,1]. Its Pickands dependence function is
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.
Copulas.AsymMixedCopula Type
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
θ₁ ≥ 0
θ₁ + θ₂ ≤ 1
θ₁ + 2θ₂ ≤ 1
θ₁ + 3θ₂ ≥ 0
Its Pickands dependence function is
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.
Copulas.BB10Copula Type
BB10Generator{T}, BB10Copula{d, T}Fields:
θ::Real - parameter
δ::Real - parameter
Constructor
BB10Generator(θ, δ)
BB10Copula(d, θ, δ)The BB10 copula has parameters
References:
- [4] Joe, H. (2014). Dependence modeling with copulas. CRC press, Page.206-207
Copulas.BB1Copula Type
BB1Generator{T}, BB1Copula{d, T}Fields:
θ::Real - parameter
δ::Real - parameter
Constructor
BB1Generator(θ, δ)
BB1Copula(d, θ, δ)The BB1 copula is parameterized by
Special cases:
- When δ = 1, it is the ClaytonCopula with parameter
.
References:
- [4] Joe, H. (2014). Dependence modeling with copulas. CRC press, Page.190-192
Copulas.BB2Copula Type
BB2Generator{T}, BB2Copula{d, T}Fields:
θ::Real - parameter
δ::Real - parameter
Constructor
BB2Generator(θ, δ)
BB2Copula(d, θ, δ)The BB2 copula has parameters
References:
- [4] Joe, H. (2014). Dependence modeling with copulas. CRC press, Page.193-194
Copulas.BB3Copula Type
BB3Generator{T}, BB3Copula{d, T}Fields:
θ::Real - parameter
δ::Real - parameter
Constructor
BB3Generator(θ, δ)
BB3Copula(d, θ, δ)The BB3 copula has parameters
References:
- [4] Joe, H. (2014). Dependence modeling with copulas. CRC press, Page.195-196
Copulas.BB4Copula Type
BB4Copula{T}Fields: - θ::Real - dependence parameter (θ ≥ 0) - δ::Real - shape parameter (δ > 0)
Constructor
BB4Copula(θ, δ)The BB4 copula is a two-parameter Archimax copula constructed from the Galambos tail and the Clayton generator. Its distribution function is
for
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
Copulas.BB5Copula Type
BB5Copula{T}Fields: - θ::Real - dependence parameter (θ ≥ 1) - δ::Real - shape parameter (δ > 0)
Constructor
BB5Copula(θ, δ)The BB5 copula is a two-parameter Archimax copula, constructed from the Galambos tail and the Gumbel generator. Its distribution function is
where
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
Copulas.BB6Copula Type
BB6Generator{T}, BB6Copula{d, T}Fields:
θ::Real - parameter
δ::Real - parameter
Constructor
BB6Generator(θ, δ)
BB6Copula(d, θ, δ)The BB6 copula has parameters
References:
- [4] Joe, H. (2014). Dependence modeling with copulas. CRC press, Page.200-201
Copulas.BB7Copula Type
BB7Generator{T}, BB7Copula{d, T}Fields:
θ::Real - parameter
δ::Real - parameter
Constructor
BB7Generator(θ, δ)
BB7Copula(d, θ, δ)The BB7 copula is parameterized by
References:
- [4] Joe, H. (2014). Dependence modeling with copulas. CRC press, Page.202-203
Copulas.BB8Copula Type
BB8Generator{T}, BB8Copula{d, T}Fields:
ϑ::Real - parameter
δ::Real - parameter
Constructor
BB8Generator(ϑ, δ)
BB8Copula(d, ϑ, δ)The BB8 copula has parameters
where
References:
- [4] Joe, H. (2014). Dependence modeling with copulas. CRC press, Page.204-205
Copulas.BB9Copula Type
BB9Generator{T}, BB9Copula{d, T}Fields:
ϑ::Real - parameter
δ::Real - parameter
Constructor
BB9Generator(ϑ, δ)
BB9Copula(d, ϑ, δ)The BB9 copula has parameters
References:
- [4] Joe, H. (2014). Dependence modeling with copulas. CRC press, Page.205-206
Copulas.BC2Copula Type
BC2Tail{T}, BC2Copula{T}Fields:
a::Real — parameter (a ∈ [0,1])
b::Real — parameter (b ∈ [0,1])
Constructor
BC2Copula(a, b)
ExtremeValueCopula(2, BC2Tail(a, b))The bivariate BC2 extreme-value copula is parameterized by two parameters
References:
- [50] Mai, J. F., & Scherer, M. (2011). Bivariate extreme-value copulas with discrete Pickands dependence measure. Extremes, 14, 311-324. Springer, 2011.
Copulas.BernsteinCopula Type
BernsteinCopula{d}Fields:
m::NTuple{d,Int}- polynomial degrees (smoothing parameters)weights::Array{Float64, d}- precomputed grid of box measures
Constructor
BernsteinCopula(C; m=10)
BernsteinCopula(data; m=10)The Bernstein copula in dimension
It is a polynomial approximation of the base copula
Implementation notes:
The grid of box measures (weights) is fully precomputed and stored as an
-dimensional array at construction. This enables fast evaluation of the copula and its density, but can be memory-intensive for large or . The choice of
mcontrols the smoothness of the approximation: largermyields finer approximation but exponentially increases memory and computation cost (boxes). For high dimensions or large
, memory usage may become prohibitive; see documentation for scaling behavior. If
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.
Copulas.BetaCopula Type
BetaCopula{d, MT}Fields:
ranks::MT- ranks matrix (d × n), each row contains integers 1..n
Constructor
BetaCopula(u)The empirical beta copula in dimension
where
Notes:
This is always a valid copula for any finite sample size
n.Supports
cdf,logpdfat 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.
Copulas.CheckerboardCopula Type
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. TypicallyDict{NTuple{d,Int}, Float64}built withStatsBase.proportionmap.
Constructor:
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
c(u) = w_k × ∏_i m[i] when u ∈ box k, and 0 otherwise.The CDF admits the multilinear overlap form
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
misnothing, we usem = fill(n, d)wheren = size(X, 2).When
pseudo_values=true(default),Xmust already be pseudo-observations in [0,1]. Otherwise pass raw data and setpseudo_values=falseto convert viapseudos(X).Each
m[i]must dividento 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.
Copulas.ClaytonCopula Type
ClaytonGenerator{T}, ClaytonCopula{d, T}Fields:
- θ::Real - parameter
Constructor
ClaytonGenerator(θ)
ClaytonCopula(d, θ)The Clayton copula in dimension
with the continuous extension
Special cases (for the copula in dimension
When
, it is the WCopula (Lower Fréchet–Hoeffding bound) When
, 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.
Copulas.CopulaModel Type
CopulaModel{CT, TM, TD} <: StatsBase.StatisticalModelA 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 (orSklarDist).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.
Copulas.CuadrasAugeCopula Type
CuadrasAugeTail{T}, CuadrasAugeCopula{T}Fields:
- θ::Real — dependence parameter, θ ∈ [0,1]
Constructor
CuadrasAugeCopula(θ)
ExtremeValueCopula(2, CuadrasAugeTail(θ))The (bivariate) Cuadras-Augé extreme-value copula is parameterized by
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.
Copulas.EmpiricalCopula Type
EmpiricalCopula{d, MT}Fields:
u::MT— pseudo-observation matrix of size(d, N).
Constructor
EmpiricalCopula(u; pseudo_values=true)The empirical copula in dimension
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
Notes:
This is an empirical object based on pseudo-observations; it is not necessarily a true copula for finite
but is widely used for nonparametric inference. Supports
cdf,logpdfat observed points, random sampling, and subsetting.
References:
- [3] Nelsen, Roger B. An introduction to copulas. Springer, 2006.
Copulas.EmpiricalEVCopula Type
EmpiricalEVTail, EmpiricalEVCopulaFields:
tgrid::Vector{Float64}— evaluation grid in (0,1)Ahat::Vector{Float64}— estimated Pickands function values ontgridslope::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
Â(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
Copulas.FGMCopula Type
FGMCopula{d,T}Fields:
- θ::Real - parameter
Constructor
FGMCopula(d, θ)The multivariate Farlie–Gumbel–Morgenstern (FGM) copula of dimension d has
where
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.
Copulas.FrankCopula Type
FrankGenerator{T}, FrankCopula{d, T}Fields:
- θ::Real - parameter
Constructor
FrankGenerator(θ)
FrankCopula(d,θ)The Frank copula in dimension
Special cases:
When
, it is the WCopula (Lower Fréchet–Hoeffding bound) When
, 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.
Copulas.GalambosCopula Type
GalambosTail{T}, GalambosCopula{T}Fields:
- θ::Real — dependence parameter, θ ≥ 0
Constructor
GalambosCopula(θ)
ExtremeValueCopula(2, GalambosTail(θ))The (bivariate) Galambos extreme-value copula is parameterized by
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.
Copulas.GaussianCopula Type
GaussianCopula{d, MT}Fields:
Σ::MT— correlation matrix (the constructor coerces the input to a correlation matrix).
Constructors
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 ρ:
Σ = 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
where
Example usage:
C = GaussianCopula(Σ)
u = rand(C, 1000)
pdf(C, u); cdf(C, u)
Ĉ = fit(GaussianCopula, u)Special case:
- If
isdiag(Σ), the constructor returnsIndependentCopula(d).
References:
- [3] Nelsen, Roger B. An introduction to copulas. Springer, 2006.
Copulas.GumbelBarnettCopula Type
GumbelBarnettGenerator{T}, GumbelBarnettCopula{d, T}Fields:
- θ::Real - parameter
Constructor
GumbelBarnettGenerator(θ)
GumbelBarnettCopula(d,θ)The Gumbel-Barnett copula is an archimdean copula with generator:
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.
Copulas.GumbelCopula Type
GumbelGenerator{T}, GumbelCopula{d, T}Fields:
- θ::Real - parameter
Constructor
GumbelGenerator(θ)
GumbelCopula(d,θ)The Gumbel copula in dimension
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.
Copulas.HuslerReissCopula Type
HuslerReissTail{T}, HuslerReissCopula{T}Fields:
- θ::Real — dependence parameter, θ ≥ 0
Constructor
HuslerReissCopula(θ)
ExtremeValueCopula(2, HuslerReissTail(θ))The (bivariate) Hüsler-Reiss extreme-value copula is parameterized by
where
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.
Copulas.IndependentCopula Type
IndependentCopula(d)The independent copula in dimension
It is Archimedean with generator
References:
- [3] Nelsen, Roger B. An introduction to copulas. Springer, 2006.
Copulas.InvGaussianCopula Type
InvGaussianGenerator{T}, InvGaussianCopula{d, T}Fields:
- θ::Real - parameter
Constructor
InvGaussianGenerator(θ)
InvGaussianCopula(d,θ)The Inverse Gaussian copula in dimension
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.Special cases:
- When θ = 0, it is the IndependentCopula
References:
- [3] Nelsen, Roger B. An introduction to copulas. Springer, 2006.
Copulas.JoeCopula Type
JoeGenerator{T}, JoeCopula{d, T}Fields:
- θ::Real - parameter
Constructor
JoeGenerator(θ)
JoeCopula(d,θ)The Joe copula in dimension
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.
Copulas.LogCopula Type
LogTail{T}, LogCopula{T}Fields:
- θ::Real — dependence parameter, θ ∈ [0,1]
Constructor
LogCopula(θ)
ExtremeValueCopula(2, LogTail(θ))The (bivariate) Mixed extreme-value copula is parameterized by
Special cases:
- θ = 0 ⇒ IndependentCopula
References:
- [49] : Tawn, Jonathan A. "Bivariate extreme value theory: models and estimation." Biometrika 75.3 (1988): 397-415.
Copulas.MCopula Type
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
Both Fréchet–Hoeffding bounds are Archimedean copulas.
References:
- [3] Nelsen, Roger B. An introduction to copulas. Springer, 2006.
Copulas.MOCopula Type
MOTail{T}, MOCopula{T}Fields:
λ₁::Real — parameter ≥ 0
λ₂::Real — parameter ≥ 0
λ₁₂::Real — parameter ≥ 0
Constructor
MOCopula(λ₁, λ₂, λ₁₂)
ExtremeValueCopula(2, MOTail(λ₁, λ₂, λ₁₂))The (bivariate) Marshall-Olkin extreme-value copula is parameterized by
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.
Copulas.MixedCopula Type
MixedTail{T}, MixedCopula{T}Fields:
- θ::Real — dependence parameter, θ ∈ [0,1]
Constructor
MixedCopula(θ)
ExtremeValueCopula(2, MixedTail(θ))The (bivariate) Mixed extreme-value copula is parameterized by
Special cases:
- θ = 0 ⇒ IndependentCopula
References:
- [49] : Tawn, Jonathan A. "Bivariate extreme value theory: models and estimation." Biometrika 75.3 (1988): 397-415.
Copulas.PlackettCopula Type
PlackettCopula{P}Fields: - θ::Real - parameter
Constructor
PlackettCopula(θ)The Plackett copula is parameterized by
and for
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.
Copulas.RafteryCopula Type
RafteryCopula{d, P}Fields: - θ::Real - parameter
Constructor
RafteryCopula(d, θ)The multivariate Raftery copula of dimension d is parameterized by
where
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.
Copulas.SklarDist Type
SklarDist{CT,TplMargins}Fields:
C::CT- The copulam::TplMargins- a Tuple representing the marginal distributions
Constructor
SklarDist(C,m)Construct a joint distribution via Sklar's theorem from marginals and a copula. See Sklar's theorem:
For every random vector
The resulting random vector follows the Distributions.jl API (rand/cdf/pdf/logpdf). A fit method is also provided. Example:
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.
Copulas.SurvivalCopula Type
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:
SurvivalCopula(ClaytonCopula(4, θ), (2,3))which is equivalent to the explicit type form:
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:
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
.
References:
- [3] Nelsen (2006), An introduction to copulas.
Copulas.TCopula Type
TCopula{d, df, MT}Fields:
df::Int— degrees of freedomΣ::MT— correlation matrix
Constructor
TCopula(df, Σ)The Student t copula is the copula of a multivariate Student t distribution. It is defined by
where
Example usage:
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.
Copulas.TiltedGenerator Type
TiltedGenerator(G, p, sJ)Archimedean generator tilted by conditioning on p components fixed at values with cumulative generator sum sJ = ∑ ϕ⁻¹(u_j). It defines
ϕ_tilt(t) = ϕ^{(p)}(sJ + t) / ϕ^{(p)}(sJ)and higher derivatives accordingly:
ϕ_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.
sourceCopulas.WCopula Type
WCopulaThe lower Fréchet–Hoeffding bound is the copula with the smallest value among all copulas. Note that
Both Fréchet–Hoeffding bounds are Archimedean copulas.
References:
- [3] Nelsen, Roger B. An introduction to copulas. Springer, 2006.
Copulas.WilliamsonGenerator Type
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
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
This function has several properties:
We have that
and is times derivable, and the signs of its derivatives alternates : . is convex.
These properties makes this function what is called a d-monotone archimedean generator, able to generate archimedean copulas in dimensions up to Generator interface: the function
G = WilliamsonGenerator(X, d)
ϕ(G,t)Note that you'll always have:
max_monotony(WilliamsonGenerator(X,d)) === dSpecial case (finite-support discrete X)
If
X isa Distributions.DiscreteUnivariateDistributionandsupport(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
WilliamsonGeneratoris 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.
Copulas.tEVCopula Type
tEVTail{Tdf,Tρ}, tEVCopula{T}Fields:
ν::Real — degrees of freedom (ν > 0)
ρ::Real — correlation parameter (ρ ∈ (-1,1])
Constructor
tEVCopula(ν, ρ)
ExtremeValueCopula(2, tEVTail(ν, ρ))The (bivariate) extreme-t copula is parameterized by
Where
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.
Copulas.EmpiricalGenerator Method
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
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
EmpiricalGeneratortype should instead treat the result as aGenerator.
References
[41]
[42]
[36] Genest, Neslehova and Ziegel (2011), Inference in Multivariate Archimedean Copula Models
Copulas.condition Method
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 copulaX::SklarDist: joint distribution with copulaX.Cand marginalsX.mjs: indices of conditioned coordinates (tuple, NTuple, or vector)u_js: values in [0,1] forU_js(when conditioning a copula)x_js: values on original scale forX_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 aDistortionon [0,1] describingU_i | U_js = u_js.condition(X, js, x_js)returns an unconditional univariate distribution forX_i | X_js = x_js, computed as the push-forwardD(X.m[i])whereD = condition(C, js, u_js)andu_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 aSklarDist(ConditionalCopula, distortions).condition(X, js, x_js)returns the conditional joint distribution on the original scale as aSklarDistwith copulaConditionalCopula(C, js, u_js)and appropriately distorted marginalsD_k(X.m[i_k]).
Notes
For best performance, pass
jsandu_jsas NTuple to keepp = length(js)known at compile time. The specialized methodcondition(::Copula{2}, j, u_j)exploits this for the commonD = 2, d = 1case.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” isConditionalCopula(C, js, u_js), i.e., the copula of that conditional distribution.
Copulas.inverse_rosenblatt Method
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.
Copulas.pseudos Method
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.
Copulas.rosenblatt Method
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.
Copulas.subsetdims Method
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
SubsetCopulaor a newSklarDistobject corresponding to the selected dimensions. Ifp == 1, returns aUniformdistribution 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.
sourceReferences
R. B. Nelsen. An Introduction to Copulas. 2nd ed Edition, Springer Series in Statistics (Springer, New York, 2006).
H. Joe. Dependence Modeling with Copulas (CRC press, 2014).
A. Sklar. Fonctions de Repartition à n Dimension et Leurs Marges. Université Paris 8, 1–3 (1959).
M. Rosenblatt. Remarks on a multivariate transformation. Annals of Mathematical Statistics 23, 470–472 (1952).
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).
R. E. Williamson. On multiply monotone functions and their laplace transforms (Mathematics Division, Office of Scientific Research, US Air Force, 1955).
C. Genest, J. Nešlehová and J. Ziegel. Inference in Multivariate Archimedean Copula Models. TEST 20, 223–256 (2011).
A. J. McNeil and J. Nešlehová. Multivariate Archimedean copulas,
-monotone functions and -norm symmetric distributions. Annals of Statistics 37, 3059–3097 (2009). R. E. Williamson. Multiply monotone functions and their Laplace transforms. Duke Mathematical Journal 23, 189–207 (1956).
H. Joe. Families of min-stable multivariate exponential and multivariate extreme value distributions. Statistics & probability letters 9, 75–81 (1990).
J. A. Tawn. Bivariate extreme value theory: models and estimation. Biometrika 75, 397–415 (1988).
J.-F. Mai and M. Scherer. Bivariate extreme-value copulas with discrete Pickands dependence measure. Extremes 14, 311–324 (2011).
J.-F. Mai and M. Scherer. Simulating copulas: stochastic models, sampling algorithms, and applications. Vol. 4 (World Scientific, 2012).
J. Galambos. Order statistics of samples from multivariate distributions. Journal of the American Statistical Association 70, 674–680 (1975).
J. Hüsler and R.-D. Reiss. Maxima of normal random vectors: between independence and complete dependence. Statistics & Probability Letters 7, 283–286 (1989).
A. K. Nikoloulopoulos, H. Joe and H. Li. Extreme value properties of multivariate t copulas. Extremes 12, 129–148 (2009).
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).
A. Charpentier, A.-L. Fougères, C. Genest and J. Nešlehová. Multivariate archimax copulas. Journal of Multivariate Analysis 126, 118–136 (2014).
J. Segers, M. Sibuya and H. Tsukahara. The empirical beta copula. Journal of Multivariate Analysis 155, 35–51 (2017).
A. Sancetta and S. Satchell. The Bernstein copula and its applications to modeling and approximations of multivariate distributions. Econometric theory 20, 535–562 (2004).
M. E. Johnson. Multivariate statistical simulation: A guide to selecting and generating continuous multivariate distributions. Vol. 192 (John Wiley & Sons, 1987).
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).
T. Saali, M. Mesfioui and A. Shabri. Multivariate extension of Raftery copula. Mathematics 11, 414 (2023).