Skip to content

Extreme Value family

Extreme value copulas are fundamental in the study of rare and extreme events due to their ability to model dependency in situations of extreme risk. This package provides a wide selection of bivariate extreme value copulas; multivariate cases are not yet implemented. Feel free to open an issue or propose a pull request if you want to contribute a multivariate case.

Info: Only Bivariate

The implementation here only deals with bivariate extreme value copulas. Multivariate cases are more tedious to implement, but not impossible: if you want to propose an implementation, we can provide guidance on how to merge it here. Do not hesitate to reach us on GitHub.

A bivariate extreme value copula [44] C has the following characteristic property:

C(u1t,u2t)=C(u1,u2)t,t>0.

It can be represented through its stable tail dependence function ():

C(u1,u2)=exp{(log(u1),log(u2))},

or through a convex function A:[0,1][1/2,1] satisfying max(t,t1)A(t)1, called its Pickands dependence function:

C(u1,u2)=exp{log(u1u2)A(log(u1)log(u1u2))},

In the context of bivariate extreme value copulas, the functions and A are related as follows:

(u1,u2)=(u1+u2)A(u1u1+u2).
Tip: Only `A` is needed

In our implementation, it is sufficient to provide the Pickands dependence function A to construct the extreme value copula and have it work correctly. Providing the other functions would, of course, improve performance.

In this package, there is an abstract type ExtremeValueCopula that provides a foundation for defining bivariate extreme value copulas. Many extreme value copulas are already implemented for you! See this list to get an overview.

If you do not find the one you need, you may define it yourself by subtyping ExtremeValueCopula. The API requires only a method for the Pickands function A(C::ExtremeValueCopula) = .... By providing this function, you can easily create a new extreme value copula that fits your specific needs:

julia
struct MyExtremeValueCopula{P} <: ExtremeValueCopula{P}
    θ::P
end

A(C::ExtremeValueCopula, t) = (t^C.θ + (1 - t)^C.θ)^(1/C.θ) # This is the Pickands function of the Logistic (Gumbel) Copula
Info: Empirical EV estimator available

When you have bivariate data and want a nonparametric Extreme Value copula, you can estimate the Pickands function from pseudo-observations using EmpiricalEVTail and plug it into an ExtremeValueCopula. For convenience, EmpiricalEVCopula(u) builds the EV copula in one step. See the bestiary entry for EmpiricalEVTail and ExtremeValueCopula docs for details.

Advanced Concepts

Here, we present some important concepts from the theory of extreme value copulas that are useful for the development of this package.

Let (X,Y)C where C is a bivariate extreme value copula. We have the following result from [45]:

Property: Ghoudi 1998

Let (X,Y)C, where C is an extreme value copula. The joint distribution of X and Z=log(X)log(XY) is given by:

P(Zz,Xx)=G(z,x)=(z+z(1z)A(z)A(z))xA(z)/z,0x,z1

where A(z) denotes the derivative of A(z) at point z.

Since A is a convex function defined on [0,1] and satisfies 1A(z)1, by extension, we define A(1) as the supremum of A(z) over (0,1). By setting x=1 in the previous result, we obtain the marginal distribution of Z: P(Zz)=GZ(z)=z+z(1z)A(z)A(z),0z1.

This result was demonstrated by Deheuvels (1991) [46] in the case where A admits a second derivative.

Simulation of Bivariate Extreme Value Distributions

To simulate a bivariate extreme value distribution C(x,y), note that if F1 and F2 are univariate extreme value distributions, then the pair (F11(X),F21(Y)) is distributed according to a bivariate extreme value distribution. The proposed algorithm in Ghoudi, 1998 [45] allows simulating such a distribution.

Assume A has a second derivative, making the distribution absolutely continuous. In this case, Z is also absolutely continuous and has a density gZ(z) given by:

gZ(z)=ddzGZ(z)=1+(1z)1(A(z)zA(z))

The conditional distribution of W given Z is:

F(w|z)=1gZ(z)ddzF(z,w),

which simplifies to:

F(w|z)=wz(1z)A(z)A(z)gZ(z)+(wwlogw)(1z(1z)A(z)A(z)gZ(z))

Given Z, the distribution of W is uniform on (0,1) with probability p(Z) and equals the product of two independent uniforms on (0,1) with probability 1p(Z), where:

p(z)=z(1z)A(z)A(z)gZ(z)

Since gZ(z) is the derivative of the cumulative distribution function of Z, it holds that 0p(z)1.

For the class of Extreme Value Copulas, We follow the methodology proposed by Ghoudi,1998. page 191. [45]. Here, is a detailed algorithm for sampling from bivariate Extreme Value Copulas:

::: algorithm Bivariate Extreme Value Copulas sampling

  • Simulate U1,U2U[0,1]

  • Simulate ZGZ(z)

  • Select W=U1 with probability p(Z) and W=U1U2 with probability 1p(Z)

  • Return X=WZ/A(Z) and Y=W(1Z)/A(Z)

:::

Note that all functions present in the algorithm were previously defined to ensure that the implemented methodology has a solid theoretical basis.

Copulas.Tail Type
julia
Tail

Abstract type. Implements the API for stable tail dependence functions (STDFs) of extreme-value copulas in dimension d.

A STDF is a function :R+d[0,) that is 1-homogeneous ((t·x)=t·(x) for all t0), convex, and satisfies the bounds max(x1,,xd)(x)x1++xd (in particular (ei)=1).

Pickands representation. By homogeneity, for x0 let x1=x1++xd and ω=x/x1Δd1. There exists a Pickands dependence function A:Δd1[0,1] (convex, max(ωi)A(ω)1) such that (x)=x1·A(ω). For d=2, A reduces to a convex function on [0,1] with max(t,1t)A(t)1 and A(0)=A(1)=1.

Interface.

  • A(tail::Tail, ω::NTuple{d,Real}) — Pickands function on the simplex \Delta_{d-1}. (For d=2, a convenience A(tail::Tail{2}, t::Real) may be provided.)

  • ℓ(tail::Tail, x::NTuple{d,Real}) — STDF. By default the package defines ℓ(tail, x) = ‖x‖₁ * A(tail, x/‖x‖₁) when A is available.

We do not algorithmically verify convexity/bounds; implementers are responsible for validity.

Additional helpers (with defaults).

  • For d=2: dA, d²A via AD; stable logpdf/rand (Ghoudi sampler).

  • In any d: cdf(u) = exp(-ℓ(-log.(u))).

References:

  • Pickands (1981); Gudendorf & Segers (2010); Ghoudi, Khoudraji & Rivest (1998); de Haan & Ferreira (2006).

  • Rasell

source
Copulas.ExtremeValueCopula Type
julia
ExtremeValueCopula{d, TT}

Constructor

julia
ExtremeValueCopula(d, tail::Tail)

Extreme-value copulas model tail dependence via a stable tail dependence function (STDF) or, equivalently, via a Pickands dependence function A. In any dimension d, the copula cdf is

C(u)=\xp(\ll(logu1,,logud)).

For d=2, write x=logu, y=logv, s=x+y, and t=x/s. The relation between and A is

(x,y)=sA(t),A:[0,1][1/2,1],A(0)=A(1)=1, A convex.

Usage

  • Provide any valid tail tail::Tail (which implements A and/or ) to construct the copula.

  • Sampling, cdf, and logpdf follow the standard Distributions.jl API.

Example

julia
C = ExtremeValueCopula(2, GalambosTail(θ))
U = rand(C, 1000)
logpdf.(Ref(C), eachcol(U))

References:

  • [44] G., & Segers, J. (2010). Extreme-value copulas. In Copula Theory and Its Applications (pp. 127-145). Springer.

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

  • [47] Mai, J. F., & Scherer, M. (2014). Financial engineering with copulas explained (p. 168). London: Palgrave Macmillan.

source

Conditionals and distortions

For any copula C, the conditional copula and the univariate conditional distortions are given by partial-derivative ratios. If we condition on a set J with m=|J| components and write I={1,,d}J, then

CIJ(uIuJ)=muJC(uI,uJ)/muJC(1I,uJ),

and, for a single coordinate iI (setting the other coordinates in I to 1), the conditional distortion is

HiJ(uuJ)=muJC(,ui=u,,uJ)/muJC(1,uJ).

For a bivariate extreme value copula with Pickands function A, the copula is

C(u1,u2)=exp{log(u1u2)A(logu1log(u1u2))},

so the above derivatives can be written explicitly in terms of A (and A when it exists) by the chain rule. In the implementation, these derivatives are obtained directly from this representation, using analytic formulas when available and automatic differentiation otherwise.

Visual illustrations

Pickands dependence functions A(t)

julia
using Copulas, Plots, Distributions
ts = range(0.0, 1.0; length=401)
Cs = (
    GalambosCopula(2, 0.8),    # upper tail dep.
    HuslerReissCopula(2, 1.0), # intermediate
    LogCopula(2, 1.6),         # asymmetric
)
labels = ("Galambos(0.8)", "Hüsler–Reiss(1.0)", "Log(1.6)")
plot(size=(700, 300))
for (i, C) in enumerate(Cs)
    plot!(ts, Copulas.A.(C.tail, ts); label=labels[i])
end
plot!(ts, max.(ts, 1 .- ts); label="bounds", ls=:dash, color=:black)
plot!(ts, ones(length(ts)); label="1", ls=:dot, color=:gray)

Sample scatter (uniform scale)

julia
C = GalambosCopula(2, 1.0)
plot(C, title="Galambos copula sample")

Conditional distortion (EV example)

julia
C = HuslerReissCopula(2, 1.2)
u2 = 0.4
D = condition(C, 2, u2)
ts = range(0.0, 1.0; length=401)
plot(ts, cdf.(Ref(D), ts); xlabel="u", ylabel="H_{1|2}(u|u₂=0.4)",
     title="Conditional distortion for Hüsler–Reiss")

Rosenblatt sanity check (EV)

julia
using StatsBase
U = rand(C, 2000)
S = reduce(hcat, (rosenblatt(C, U[:, i]) for i in 1:size(U,2)))
ts = range(0.0, 1.0; length=401)
EC = [ecdf(S[k, :]) for k in 1:2]
plot(ts, ts; label="Uniform", color=:blue, alpha=0.6, size=(650,300))
plot!(ts, EC[1].(ts); seriestype=:steppost, label="s₁", color=:black)
plot!(ts, EC[2].(ts); seriestype=:steppost, label="s₂", color=:gray)

Available models

MTail

Copulas.MTail Type
julia
MTail

Corresponds to the MCopula viewed as an etreme value copula.

source

NoTail

Copulas.NoTail Type
julia
NoTail

Corresponds to the case where the pickads function is identically One, which means no particular tail behavior.

source

AsymGalambosTail

Copulas.AsymGalambosTail 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

AsymLogTail

Copulas.AsymLogTail 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

AsymMixedTail

Copulas.AsymMixedTail 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

BC2Tail

Copulas.BC2Tail 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

CuadrasAugeTail

Copulas.CuadrasAugeTail 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

GalambosTail

Copulas.GalambosTail 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

HuslerReissTail

Copulas.HuslerReissTail 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

LogTail

Copulas.LogTail 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

MixedTail

Copulas.MixedTail 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

MOTail

Copulas.MOTail 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

tEVTail

Copulas.tEVTail 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

EmpiricalEVTail

Copulas.EmpiricalEVTail 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

References

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

  2. G. Gudendorf and J. Segers. Extreme-value copulas. In: Copula Theory and Its Applications: Proceedings of the Workshop Held in Warsaw, 25-26 September 2009 (Springer, 2010); pp. 127–145.

  3. K. Ghoudi, A. Khoudraji and E. L.-P. Rivest. Propriétés statistiques des copules de valeurs extrêmes bidimensionnelles. Canadian Journal of Statistics 26, 187–197 (1998).

  4. P. Deheuvels. On the limiting behavior of the Pickands estimator for bivariate extreme-value distributions. Statistics & Probability Letters 12, 429–439 (1991).

  5. J.-F. Mai and M. Scherer. Financial engineering with copulas explained (Springer, 2014).

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

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

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

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

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

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

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