Extremes.AbstractFittedExtremeValueModelType
AbstractFittedExtremeValueModel{T<:AbstractExtremeValueModel}

Abstract type containing the fitted extreme value model types.

  • BayesianAbstractExtremeValueModel
  • MaximumLikelihoodAbstractExtremeValueModel
  • pwmAbstractExtremeValueModel
Extremes.BlockMaximaMethod
BlockMaxima{GeneralizedExtremeValue}(data::Vector{<:Real};
    locationcov::Vector{Variable} = Vector{Variable}(),
    logscalecov::Vector{Variable} = Vector{Variable}(),
    shapecov::Vector{Variable} = Vector{Variable}())::BlockMaxima

Creates a BlockMaxima structure.

Extremes.BlockMaximaMethod
BlockMaxima{Gumbel}(data::Vector{<:Real};
    locationcov::Vector{Variable} = Vector{Variable}(),
    logscalecov::Vector{Variable} = Vector{Variable}())::BlockMaxima

Creates a BlockMaxima{Gumbel} structure.

Extremes.ClusterType
Cluster(u₁::Real,u₂::Real,position::Vector{<:Int},value::Vector{<:Real})

Cluster type.

Extremes.FlatType
Flat()

Construct a Flat <: ContinuousUnivariateDistribution object representing an improper uniform distribution on the real line.

Extremes.ThresholdExceedanceMethod
ThresholdExceedance(exceedances::Vector{<:Real};
    logscalecov::Vector{<:DataItem} = Vector{Variable}(),
    shapecov::Vector{<:DataItem} = Vector{Variable}())::ThresholdExceedance

Creates a ThresholdExceedance structure.

Extremes.VariableType
Variable(name::String, value :: Vector{<:Real})

Construct a Variable type

Extremes.VariableStdMethod
VariableStd(name::String, z::Vector{<:Real})::VariableStd

Construct a VariableStd type from the standardized vector z with the name name.

Base.lengthMethod
length(c::Cluster)

Compute the cluster length.

Base.maximumMethod
max(c::Cluster)

Compute the cluster maximum.

Base.showMethod
Base.show(io::IO, obj::AbstractExtremeValueModel)

Override of the show function for the objects of type AbstractExtremeValueModel.

Base.showMethod
Base.show(io::IO, obj::AbstractFittedExtremeValueModel)

Override of the show function for the objects of type AbstractFittedExtremeValueModel.

Base.showMethod
Base.show(io::IO, obj::Cluster)

Override of the show function for the objects of type AbstractExtremeValueModel.

Base.showMethod
Base.show(io::IO, obj::ReturnLevel)

Override of the show function for the objects of type ReturnLevel.

Base.sumMethod
sum(c::Cluster)

Compute the cluster sum.

Distributions.locationMethod
location(fm::AbstractFittedExtremeValueModel)

Return the location parameters of the fitted model.

Distributions.scaleMethod
scale(fm::AbstractFittedExtremeValueModel)

Return the scale parameters of the fitted model.

Distributions.shapeMethod
shape(fm::AbstractFittedExtremeValueModel)

Return the shape parameters of the fitted model.

Distributions.shapeMethod
shape(pd::Gumbel)

Return the Gumbel distribution shape parameter value, i.e. 0.

Extremes.aicMethod
aic(fm:::MaximumLikelihoodAbstractExtremeValueModel)

Compute the Akaike information criterion (AIC) of the fitted model by maximum likelihood method.

Details

The AIC is defined as follows:

$AIC = 2 k - 2 \log \hat{L};$

where $k$ is the number of estimated parameters and $\hat{L}$ is the maximized value of the likelihood function for the model.

Extremes.bicMethod
bic(fm:::MaximumLikelihoodAbstractExtremeValueModel)

Compute the Bayesian information criterion (BIC) of the fitted model by maximum likelihood method.

Details

The BIC is defined as follows:

$BIC = k \log n - 2 \log \hat{L};$

where $k$ is the number of estimated parameters, $n$ is the number of data and $\hat{L}$ is the maximized value of the likelihood function for the model.

Extremes.buildVariablesMethod
buildVariables(df::DataFrame, ids::Vector{Symbol})::Vector{Variable}

Build the Variable type from the columns ids of the DataFrame df.

Example

julia> df = Extremes.dataset("fremantle")
julia> Extremes.buildVariables(df, [:Year, :SOI])
Extremes.checknonstationarityMethod
checknonstationarity(model::AbstractExtremeValueModel)

Check if the extreme value model model is nonstationary.

Extremes.checkstationarityMethod
checkstationarity(model::AbstractExtremeValueModel)

Check if the extreme value model model is stationary.

Extremes.cintFunction
cint(..., confidencelevel::Real=.95)

Compute confidence interval or credible interval in the case of Bayesian estimation.

The function can be applied on any AbstractFittedExtremeValueModel subtype to obtain a confidence interval on the model parameters. It can also be applied on ReturnLevel type to obtain a confidence interval on the return level.

Implementation

The method used for computing the interval depends on the estimation method. In the case of maximum likelihood estimation, the confidence intervals are computed using the Wald approximation based on the approximate parameter estimates covariance matrix. In the case of Bayesian estimation, the return interval is the highest posterior density estimate based on the MCMC sample. In the case of probability weighted moment estimation, the intervals are computed using a boostrap procedure.

Extremes.computeparamfunctionMethod
computeparamfunction(covariates::Vector{Variable})

Establish the parameter as function of the corresponding covariates.

Extremes.datasetMethod
dataset(name::String)::DataFrame

Load the dataset associated with name.

Some datasets used by Coles (2001) are available using the following names:

  • portpirie: annual maximum sea-levels in Port Pirie,
  • glass: breaking strengths of glass fibers
  • fremantle: annual maximum sea-levels in Fremantle
  • rain: daily rainfall accumulations in south-west England
  • wooster: daily minimum temperatures recorded in Wooster
  • dowjones: daily closing prices of the Dow Jones Index

These datasets have been retrieved using the R package ismev.

Examples

julia> Extremes.dataset("portpirie")
Extremes.deltaMethod
delta(g::Function, θ̂::AbstractVector{<:Real}, H::AbstractPDMat)

Compute the variance of the function g of estimated paramters θ̂ with negative observed information matrix H.

Detail

Hcorresponds to the Hessian matrix of the negative log likelihood.

Extremes.ecdfMethod
ecdf(y::Vector{<:Real})::Tuple{Vector{<:Real}, Vector{<:Real}}

Compute the empirical cumulative distribution function using the Gumbel formula.

The empirical quantiles are computed using the Gumbel plotting positions as as recommended by Makkonen (2006).

Example

julia> (x, F̂) = Extremes.ecdf(y)

Reference

Makkonen, L. (2006). Plotting positions in extreme value analysis. Journal of Applied Meteorology and Climatology, 45(2), 334-340.

Extremes.findposteriormodeMethod
findposteriormode(fm::BayesianAbstractExtremeValueModel)::Vector{<:Real}

Find the maximum a posteriori probability (MAP) estimate.

Extremes.fitMethod
fit(model::AbstractExtremeValueModel; initialvalues::Vector{<:Real})::MaximumLikelihoodAbstractExtremeValueModel

Fit the extreme value model by maximum likelihood.

Extremes.fitMethod
fit(model::AbstractExtremeValueModel)::MaximumLikelihoodAbstractExtremeValueModel

Fit the extreme value model by maximum likelihood.

Extremes.fitbayesMethod
fitbayes(model::AbstractExtremeValueModel; niter::Int=5000, warmup::Int=2000)::BayesianAbstractExtremeValueModel

Fit the extreme value model under the Bayesian paradigm.

Extremes.fitpwmfunctionMethod
fitpwmfunction(fm::pwmAbstractExtremeValueModel{BlockMaxima{GeneralizedExtremeValue})::Function

Returns the corresponding fitpwm function.

Extremes.fitpwmfunctionMethod
fitpwmfunction(fm::pwmAbstractExtremeValueModel{BlockMaxima{Gumbel}})::Function

Returns the corresponding fitpwm function.

Extremes.fitpwmfunctionMethod
fitpwmfunction(fm::pwmAbstractExtremeValueModel{ThresholdExceedance})::Function

Returns the corresponding fitpwm function.

Extremes.getclusterMethod
getcluster(y::Vector{<:Real}, u₁::Real, u₂::Real)

Extract the clusters from vector y.

A cluster is defined as a sequence of values higher than threshold u₂ with at least a value higher than threshold u₁.

See also Cluster.

Extremes.getclusterMethod
getcluster(y::Vector{<:Real}, u::Real; runlength::Int=1)::Vector{Cluster}

Extract the clusters from vector y.

Threshold exceedances separated by fewer than r non-exceedances belong to the same cluster. The value r is corresponds to the runlength parameter. This approach is referred to as the runs declustering scheme (see Coles, 2001 sec. 5.3.2).

See also Cluster.

Extremes.getdistributionFunction
getdistribution(model::AbstractExtremeValueModel, θ::Vector{<:Real})
getdistribution(fm::AbstractFittedExtremeValueModel)

Return the distributions corresponding to the model or the fitted model.

If an extreme value model is provided, the distributions corresponding to the parameter vector θ are returned. If a fitted extreme value model is provident, the distributions corresponding to the parameter estimates are returned.

Implementation

In the stationary case, a single extreme value distribution is returned.

In the non-stationary case, a vector of extreme value distributions is returned, one for each data value.

In the Bayesian fitted model case, a array of distributions is returned where each column corresponds to a MCMC iteration.

Extremes.getdistributionMethod
getdistribution(fittedmodel::MaximumLikelihoodAbstractExtremeValueModel)::Vector{<:Distribution}

Return the fitted distribution in case of stationarity or the vector of fitted distribution in case of non-stationarity.

Extremes.getdistributionMethod
getdistribution(fittedmodel::pwmAbstractExtremeValueModel)::Vector{<:Distribution}

Return the fitted distribution for the model fitted with the probability weigthed moments.

Extremes.getinitialvalueFunction
getinitialvalue(model::AbstractExtremeValueModel)

Get an initial estimates of the model parameters.

Extremes.getinitialvalueMethod
getinitialvalue(::Type{GeneralizedExtremeValue},y::Vector{<:Real})::Vector{<:Real}

Compute the initial values of the GEV parameters given the data y.

If the probability weighted moment estimations are valid, then those values are returned. Otherwise, the probability weighted moment estimations of the Gumbel distribution are returned.

Example

julia> Extremes.getinitialvalue(GeneralizedExtremeValue, y)
Extremes.getinitialvalueMethod
getinitialvalue(::Type{GeneralizedPareto},y::Vector{<:Real})::Vector{<:Real}

Compute the initial values of the GP distribution parameters given the data y.

If the probability weighted moment estimations are valid, then those values are returned. Otherwise, the moment estimation of the exponential distribution is returned.

Example

julia> Extremes.getinitialvalue(GeneralizedPareto, y)
Extremes.gevfitFunction
gevfit()

Estimate the GEV parameters by maximum likelihood.

Implementation

The function uses Nelder-Mead solver implemented in the Optim.jl package to find the point where the log-likelihood is maximal.

The GEV parameters can be modeled as function of covariates as follows:

\[μ = X₁ × β₁,\]

\[ϕ = X₂ × β₂,\]

\[ξ = X₃ × β₃.\]

The covariates are standardized before estimating the parameters to help fit the model. They are transformed back on their original scales before returning the fitted model.

See also gevfit for the other methods, gevfitpwm, gevfitbayes and BlockMaxima.

Extremes.gevfitMethod
gevfit(model::{BlockMaxima{GeneralizedExtremeValue}}, initialvalues::Vector{<:Real})

Estimate the parameters of the BlockMaxima model using the given initialvalues.

Extremes.gevfitMethod
gevfit(df::DataFrame,
    datacol::Symbol,
    locationcovid = Vector{Symbol}(),
    logscalecovid = Vector{Symbol}(),
    shapecovid = Vector{Symbol}()
    )

Estimate the GEV parameters.

Arguments

  • df::DataFrame: The dataframe containing the data.
  • datacol::Symbol: The symbol of the column of df containing the block maxima data.
  • initialvalues::Vector{<:Real}: Vector of parameters initial values.
  • locationcovid::Vector{Symbol} = Vector{Symbol}(): The symbols of the columns of df containing the covariates of the location parameter.
  • logscalecovid::Vector{Symbol} = Vector{Symbol}(): The symbols of the columns of df containing the covariates of the log-scale parameter.
  • shapecovid::Vector{Symbol} = Vector{Symbol}(): The symbols of the columns of df containing the covariates of the shape parameter.
Extremes.gevfitMethod
gevfit(df::DataFrame,
    datacol::Symbol,
    locationcovid = Vector{Symbol}(),
    logscalecovid = Vector{Symbol}(),
    shapecovid = Vector{Symbol}()
    )

Estimate the GEV parameters.

Arguments

  • df::DataFrame: The dataframe containing the data.
  • datacol::Symbol: The symbol of the column of df containing the block maxima data.
  • locationcovid::Vector{Symbol} = Vector{Symbol}(): The symbols of the columns of df containing the covariates of the location parameter.
  • logscalecovid::Vector{Symbol} = Vector{Symbol}(): The symbols of the columns of df containing the covariates of the log-scale parameter.
  • shapecovid::Vector{Symbol} = Vector{Symbol}(): The symbols of the columns of df containing the covariates of the shape parameter.
Extremes.gevfitMethod
gevfit(y,
    initialvalues,
    locationcov = Vector{Variable}(),
    logscalecov = Vector{Variable}(),
    shapecov = Vector{Variable}()
    )

Estimate the GEV parameters.

Arguments

  • y::Vector{<:Real}: the vector of block maxima.
  • initialvalues::Vector{<:Real}: Vector of parameters initial values.
  • locationcov::Vector{<:DataItem} = Vector{Variable}(): The covariates of the location parameter.
  • logscalecov::Vector{<:DataItem} = Vector{Variable}(): The covariates of the log-scale parameter.
  • shapecov::Vector{<:DataItem} = Vector{Variable}(): The covariates of the shape parameter.
Extremes.gevfitMethod
gevfit(y,
    locationcov = Vector{Variable}(),
    logscalecov = Vector{Variable}(),
    shapecov = Vector{Variable}()
    )

Estimate the GEV parameters.

Arguments

  • y::Vector{<:Real}: The vector of block maxima.
  • locationcov::Vector{<:DataItem} = Vector{Variable}(): The covariates of the location parameter.
  • logscalecov::Vector{<:DataItem} = Vector{Variable}(): The covariates of the log-scale parameter.
  • shapecov::Vector{<:DataItem} = Vector{Variable}(): The covariates of the shape parameter.
Extremes.gevfitbayesFunction
gevfitbayes(..., niter::Int=5000, warmup::Int=2000)

Generate a sample from the GEV parameters' posterior distribution.

Arguments

  • niter::Int = 5000: The total number of MCMC iterations.
  • warmup::Int = 2000: The number of warmup iterations (burn-in).

Implementation

The function uses the No-U-Turn Sampler (NUTS; Hoffman and Gelman, 2014) implemented in the Mamba.jl package to generate a random sample from the posterior distribution.

Currently, only the improper uniform prior is implemented, i.e.

\[f_{(β₁,β₂,β₃)}(β₁,β₂,β₃) ∝ 1,\]

where

\[μ = X₁ × β₁,\]

\[ϕ = X₂ × β₂,\]

\[ξ = X₃ × β₃.\]

In the stationary case, this improper prior yields to a proper posterior if the sample size is larger than 3 (Northrop and Attalides, 2016).

The covariates are standardized before estimating the parameters to help fit the model. They are transformed back on their original scales before returning the fitted model.

See also gevfitbayes for the other methods, gevfitpwm, gevfit and BlockMaxima.

References

Hoffman M. D. and Gelman A. (2014). The No-U-Turn sampler: adaptively setting path lengths in Hamiltonian Monte Carlo. Journal of Machine Learning Research, 15:1593–1623.

Paul J. Northrop P. J. and Attalides N. (2016). Posterior propriety in Bayesian extreme value analyses using reference priors. Statistica Sinica, 26:721-743.

Extremes.gevfitbayesMethod
gevfitbayes(model::BlockMaxima{GeneralizedExtremeValue};
    niter::Int=5000,
    warmup::Int=2000)

Generate a sample from the BlockMaxima model parameters' posterior distribution.

Extremes.gevfitbayesMethod
gevfitbayes(df::DataFrame,
    datacol::Symbol,
    locationcovid = Vector{Symbol}(),
    logscalecovid = Vector{Symbol}(),
    shapecovid = Vector{Symbol}(),
    niter::Int=5000,
    warmup::Int=2000)

Generate a sample from the GEV parameters' posterior distribution.

Arguments

  • df::DataFrame: The dataframe containing the data.
  • datacol::Symbol: The symbol of the column of df containing the block maxima data.
  • locationcovid::Vector{Symbol} = Vector{Symbol}(): The symbols of the columns of df containing the covariates of the location parameter.
  • logscalecovid::Vector{Symbol} = Vector{Symbol}(): The symbols of the columns of df containing the covariates of the log-scale parameter.
  • shapecovid::Vector{Symbol} = Vector{Symbol}(): The symbols of the columns of df containing the covariates of the shape parameter.
Extremes.gevfitbayesMethod
gevfitbayes(y,
    locationcov = Vector{Variable}(),
    logscalecov = Vector{Variable}(),
    shapecov = Vector{Variable}(),
    niter::Int=5000,
    warmup::Int=2000
    )

Generate a sample from the GEV parameters' posterior distribution.

Arguments

  • y::Vector{<:Real}: The vector of block maxima.
  • locationcov::Vector{<:DataItem} = Vector{Variable}(): The covariates of the location parameter.
  • logscalecov::Vector{<:DataItem} = Vector{Variable}(): The covariates of the log-scale parameter.
  • shapecov::Vector{<:DataItem} = Vector{Variable}(): The covariates of the shape parameter.
Extremes.gevfitpwmFunction
gevfitpwm(...)

Estimate the GEV parameters with the probability weighted moments.

Implementation

Estimation with the probability weighted moments, as described by Hosking et al. (1985), is only possible in the stationary case.

See also gevfitpwm for the other methods, gevfit, gevfitbayes and BlockMaxima.

Reference

Hosking, J. R. M., Wallis, J. R. and Wood, E. F. (1985). Estimation of the generalized extreme-value distribution by the method of probability-weighted moments. Technometrics, 27:251-261.

Extremes.gevfitpwmMethod
gevfitpwm(model::BlockMaxima{GeneralizedExtremeValue})

Estimate the GEV parameters with the probability weighted moments.

Extremes.gevfitpwmMethod
gevfitpwm(df::DataFrame, datacol::Symbol)

Estimate the GEV parameters with the probability weighted moments.

Block maxima data are in the column datacol of the dataframe df.

Extremes.gevfitpwmMethod
gevfitpwm(y::Vector{<:Real})

Estimate the GEV parameters with the probability weighted moments.

Extremes.gpfitFunction
gpfit(...)

Estimate the GP parameters by maximum likelihood.

Data provided must be the exceedances above the threshold, i.e. the data above the threshold minus the threshold.

Implementation

The function uses Nelder-Mead solver implemented in the Optim.jl package to find the point where the log-likelihood is maximal.

The GP parameters can be modeled as function of covariates as follows:

\[ϕ = X₂ × β₂,\]

\[ξ = X₃ × β₃.\]

The covariates are standardized before estimating the parameters to help fit the model. They are transformed back on their original scales before returning the fitted model.

See also gpfit for the other methods, gpfitpwm, gpfitbayes and ThresholdExceedance.

Extremes.gpfitMethod
gpfit(df::DataFrame,
    datacol::Symbol,
    logscalecovid = Vector{Symbol}(),
    shapecovid = Vector{Symbol}()
    )

Estimate the GP parameters

Arguments

  • df::DataFrame: The dataframe containing the data.
  • datacol::Symbol: The symbol of the column of df containing the exceedances.
  • initialvalues::Vector{<:Real}: Vector of parameters initial values.
  • logscalecovid::Vector{Symbol} = Vector{Symbol}(): The symbols of the columns of df containing the covariates of the log-scale parameter.
  • shapecovid::Vector{Symbol} = Vector{Symbol}(): The symbols of the columns of df containing the covariates of the shape parameter.
Extremes.gpfitMethod
gpfit(df::DataFrame,
    datacol::Symbol,
    logscalecovid = Vector{Symbol}(),
    shapecovid = Vector{Symbol}()
    )

Estimate the GP parameters

Arguments

  • df::DataFrame: The dataframe containing the data.
  • datacol::Symbol: The symbol of the column of df containing the exceedances.
  • logscalecovid::Vector{Symbol} = Vector{Symbol}(): The symbols of the columns of df containing the covariates of the log-scale parameter.
  • shapecovid::Vector{Symbol} = Vector{Symbol}(): The symbols of the columns of df containing the covariates of the shape parameter.
Extremes.gpfitMethod
gpfit(model::ThresholdExceedance, initialvalues::Vector{<:Real})

Estimate the parameters of the ThresholdExceedance model using the given initialvalues.

Extremes.gpfitMethod
gpfit(y,
    initialvalues;
    logscalecov = Vector{Variable}(),
    shapecov = Vector{Variable}()
    )

Estimate the GP parameters

Arguments

  • y::Vector{<:Real}: The vector of exceedances.
  • initialvalues::Vector{<:Real}: The vector of parameters initial values.
  • logscalecov::Vector{<:DataItem} = Vector{Variable}(): The covariates of the log-scale parameter.
  • shapecov::Vector{<:DataItem} = Vector{Variable}(): The covariates of the shape parameter.
Extremes.gpfitMethod
gpfit(y,
    logscalecov = Vector{Variable}(),
    shapecov = Vector{Variable}()
    )

Estimate the GP parameters

Arguments

  • y::Vector{<:Real}: The vector of exceedances.
  • logscalecov::Vector{<:DataItem} = Vector{Variable}(): The covariates of the log-scale parameter.
  • shapecov::Vector{<:DataItem} = Vector{Variable}(): The covariates of the shape parameter.
Extremes.gpfitbayesFunction
gpfitbayes(..., niter::Int=5000, warmup::Int=2000)

Generate a sample from the GP parameters' posterior distribution.

Data provided must be the exceedances above the threshold, i.e. the data above the threshold minus the threshold.

Arguments

  • niter::Int = 5000: The total number of MCMC iterations
  • warmup::Int = 2000: The number of warmup iterations (burn-in).

Implementation

The function uses the No-U-Turn Sampler (NUTS; Hoffman and Gelman, 2014) implemented in the Mamba.jl package to generate a random sample from the posterior distribution.

Currently, only the improper uniform prior is implemented, i.e.

\[f_{(β₂,β₃)}(β₂,β₃) ∝ 1,\]

where

\[ϕ = X₂ × β₂,\]

\[ξ = X₃ × β₃.\]

In the stationary case, this improper prior yields to a proper posterior if the sample size is larger than 2 (Northrop and Attalides, 2016).

The covariates are standardized before estimating the parameters to help fit the model. They are transformed back on their original scales before returning the fitted model.

See also gpfitbayes for the other methods, gpfitpwm, gpfit and ThresholdExceedance.

Reference

Hoffman M. D. and Gelman A. (2014). The No-U-Turn sampler: adaptively setting path lengths in Hamiltonian Monte Carlo. Journal of Machine Learning Research, 15:1593–1623.

Paul J. Northrop P. J. and Attalides N. (2016). Posterior propriety in Bayesian extreme value analyses using reference priors. Statistica Sinica, 26:721-743.

Extremes.gpfitbayesMethod
gpfitbayes(df::DataFrame,
    datacol::Symbol,
    logscalecovid = Vector{Symbol}(),
    shapecovid = Vector{Symbol}(),
    niter::Int=5000,
    warmup::Int=2000)

Generate a sample from the GP parameters' posterior distribution.

Arguments

  • df::DataFrame: The dataframe containing the data.
  • datacol::Symbol: The symbol of the column of df containing the exceedances.
  • logscalecovid::Vector{Symbol} = Vector{Symbol}(): The symbols of the columns of df containing the covariates of the log-scale parameter.
  • shapecovid::Vector{Symbol} = Vector{Symbol}(): The symbols of the columns of df containing the covariates of the shape parameter.
Extremes.gpfitbayesMethod
gpfitbayes(model::ThresholdExceedance;
    niter::Int=5000,
    warmup::Int=2000)

Generate a sample from the GP parameters' posterior distribution.

Extremes.gpfitbayesMethod
gpfitbayes(y,
    logscalecov = Vector{Variable}(),
    shapecov = Vector{Variable}(),
    niter::Int=5000,
    warmup::Int=2000
    )

Generate a sample from the GP parameters' posterior distribution.

Arguments

  • y::Vector{<:Real}: The vector of exceedances.
  • logscalecov::Vector{<:DataItem} = Vector{Variable}(): The covariates of the log-scale parameter.
  • shapecov::Vector{<:DataItem} = Vector{Variable}(): The covariates of the shape parameter.
Extremes.gpfitpwmFunction
gpfitpwm(...)

Estimate the GP parameters with the probability weighted moments.

Implementation

Estimation with the probability weighted moments, as described by Hosking and Wallis (1987), is only possible in the stationary case.

See also gpfitpwm for the other methods, gpfit, gpfitbayes and ThresholdExceedance.

Reference

Hosking, J. R. M. and Wallis, J. R. (1987). Parameter and quantile estimation for the Generalized Pareto distribution, Technometrics, 29:339-349.

Extremes.gpfitpwmMethod
gpfitpwm(df::DataFrame, datacol::Symbol)

Estimate the GP parameters with the probability weighted moments.

Block maxima data are in the column datacol of the dataframe df.

Extremes.gpfitpwmMethod
gpfitpwm(model::ThresholdExceedance)

Estimate the GP parameters with the probability weighted moments.

Extremes.gpfitpwmMethod
gpfitpwm(y::Vector{<:Real})

Estimate the GP parameters with the probability weighted moments.

Extremes.gumbelfitFunction
gumbelfit()

Estimate the Gumbel parameters by maximum likelihood.

Details

The Gumbel distribution is a particular case of the generalized extreme value distribution when the shape parameter is zero.

Extreme value theory

In extreme value theory, it is best to avoid the choice of a sub-family of extreme value familu as the Gumbel. This is because the choice of family is made with the data at hand, and when extrapolating to large quantiles, i.e. larger than the range of the data, the uncertainty associated with this choice is not taken into account. If the data suggest that the Gumbel family is the best one, this does not imply that the other families are not plausible. In applications, the confidence intervals on the shape parameter are often wide, representing the difficulty of discriminating the tail behavior using only the limited number of data. Therefore, the use of the GEV distribution for the block maxima model makes more sense. As Coles (2001) also argued in Page 64, this is "...the safest option is to accept there is uncertainty about the value of the shape parameter ... and to prefer the inference based on the GEV model. The larger measures of uncertainty generated by the GEV model then provide a more realistic quantification of the genuine uncertainties involved in model extrapolation."

Implementation

The function uses Nelder-Mead solver implemented in the Optim.jl package to find the point where the log-likelihood is maximal.

The Gumbel parameters can be modeled as function of covariates as follows:

\[μ = X₁ × β₁,\]

\[ϕ = X₂ × β₂,\]

The covariates are standardized before estimating the parameters to help fit the model. They are transformed back on their original scales before returning the fitted model.

See also gumbelfit for the other methods and gevfit for estiming the parameters of the GEV distribution.

Extremes.gumbelfitMethod
gumbelfit(model::{BlockMaxima{Gumbel}}, initialvalues::Vector{<:Real})

Estimate the parameters of the BlockMaxima model using the given initialvalues.

Extremes.gumbelfitMethod
gumbelfit(df::DataFrame,
    datacol::Symbol,
    locationcovid = Vector{Symbol}(),
    logscalecovid = Vector{Symbol}()
    )

Estimate the Gumbel parameters.

Arguments

  • df::DataFrame: The dataframe containing the data.
  • datacol::Symbol: The symbol of the column of df containing the block maxima data.
  • initialvalues::Vector{<:Real}: Vector of parameters initial values.
  • locationcovid::Vector{Symbol} = Vector{Symbol}(): The symbols of the columns of df containing the covariates of the location parameter.
  • logscalecovid::Vector{Symbol} = Vector{Symbol}(): The symbols of the columns of df containing the covariates of the log-scale parameter.
Extremes.gumbelfitMethod
gumbelfit(df::DataFrame,
    datacol::Symbol,
    locationcovid = Vector{Symbol}(),
    logscalecovid = Vector{Symbol}(),
    shapecovid = Vector{Symbol}()
    )

Estimate the Gumbel parameters.

Arguments

  • df::DataFrame: The dataframe containing the data.
  • datacol::Symbol: The symbol of the column of df containing the block maxima data.
  • locationcovid::Vector{Symbol} = Vector{Symbol}(): The symbols of the columns of df containing the covariates of the location parameter.
  • logscalecovid::Vector{Symbol} = Vector{Symbol}(): The symbols of the columns of df containing the covariates of the log-scale parameter.
Extremes.gumbelfitMethod
gumbelfit(y,
    initialvalues,
    locationcov = Vector{Variable}(),
    logscalecov = Vector{Variable}()
    )

Estimate the Gumbel parameters.

Arguments

  • y::Vector{<:Real}: the vector of block maxima.
  • initialvalues::Vector{<:Real}: Vector of parameters initial values.
  • locationcov::Vector{<:DataItem} = Vector{Variable}(): The covariates of the location parameter.
  • logscalecov::Vector{<:DataItem} = Vector{Variable}(): The covariates of the log-scale parameter.
Extremes.gumbelfitMethod
gumbelfit(y,
    locationcov = Vector{Variable}(),
    logscalecov = Vector{Variable}()
    )

Estimate the Gumbel parameters.

Arguments

  • y::Vector{<:Real}: The vector of block maxima.
  • locationcov::Vector{<:DataItem} = Vector{Variable}(): The covariates of the location parameter.
  • logscalecov::Vector{<:DataItem} = Vector{Variable}(): The covariates of the log-scale parameter.
Extremes.gumbelfitbayesFunction
gumbelfitbayes(..., niter::Int=5000, warmup::Int=2000)

Generate a sample from the Gumbel parameters' posterior distribution.

Arguments

  • niter::Int = 5000: The total number of MCMC iterations.
  • warmup::Int = 2000: The number of warmup iterations (burn-in).

Implementation

The function uses the No-U-Turn Sampler (NUTS; Hoffman and Gelman, 2014) implemented in the Mamba.jl package to generate a random sample from the posterior distribution.

Currently, only the improper uniform prior is implemented, i.e.

\[f_{(β₁,β₂,β₃)}(β₁,β₂,β₃) ∝ 1,\]

where

\[μ = X₁ × β₁,\]

\[ϕ = X₂ × β₂,\]

In the stationary case, this improper prior yields to a proper posterior if the sample size is larger than 3 (Northrop and Attalides, 2016).

The covariates are standardized before estimating the parameters to help fit the model. They are transformed back on their original scales before returning the fitted model.

See also gevfitbayes for the other methods and BlockMaxima.

References

Hoffman M. D. and Gelman A. (2014). The No-U-Turn sampler: adaptively setting path lengths in Hamiltonian Monte Carlo. Journal of Machine Learning Research, 15:1593–1623.

Paul J. Northrop P. J. and Attalides N. (2016). Posterior propriety in Bayesian extreme value analyses using reference priors. Statistica Sinica, 26:721-743.

Extremes.gumbelfitbayesMethod
gumbelfitbayes(model::BlockMaxima{Gumbel};
    niter::Int=5000,
    warmup::Int=2000)

Generate a sample from the BlockMaxima model parameters' posterior distribution.

Extremes.gumbelfitbayesMethod
gumbelfitbayes(df::DataFrame,
    datacol::Symbol,
    locationcovid = Vector{Symbol}(),
    logscalecovid = Vector{Symbol}(),
    niter::Int=5000,
    warmup::Int=2000)

Generate a sample from the Gumbel parameters' posterior distribution.

Arguments

  • df::DataFrame: The dataframe containing the data.
  • datacol::Symbol: The symbol of the column of df containing the block maxima data.
  • locationcovid::Vector{Symbol} = Vector{Symbol}(): The symbols of the columns of df containing the covariates of the location parameter.
  • logscalecovid::Vector{Symbol} = Vector{Symbol}(): The symbols of the columns of df containing the covariates of the log-scale parameter.
Extremes.gumbelfitbayesMethod
gumbelfitbayes(y,
    locationcov = Vector{Variable}(),
    logscalecov = Vector{Variable}(),
    niter::Int=5000,
    warmup::Int=2000
    )

Generate a sample from the Gumbel parameters' posterior distribution.

Arguments

  • y::Vector{<:Real}: The vector of block maxima.
  • locationcov::Vector{<:DataItem} = Vector{Variable}(): The covariates of the location parameter.
  • logscalecov::Vector{<:DataItem} = Vector{Variable}(): The covariates of the log-scale parameter.
Extremes.gumbelfitpwmFunction
gumbelfitpwm(...)

Estimate the Gumbel parameters with the probability weighted moments.

Implementation

Estimation with the probability weighted moments, as described by Landwehr *et al. (1979), is only possible in the stationary case.

See also gumbelfitpwm for the other methods, gevfit, gevfitbayes and BlockMaxima.

Reference

Landwehr, J. M., Matalas, N. C. and Wallis, J. R. (1979). Probability weighted moments compared with some traditional techniques in estimating Gumbel parameters and quantiles. Water Resources Research, 15:1055–1064.

Extremes.gumbelfitpwmMethod
gumbelfitpwm(model::BlockMaxima{Gumbel})

Estimate the Gumbel parameters with the probability weighted moments.

Extremes.gumbelfitpwmMethod
gumbelfitpwm(df::DataFrame, datacol::Symbol)::pwmAbstractExtremeValueModel

Estimate the Gumbel parameters with the probability weighted moments.

Block maxima data are in the column datacol of the dataframe df.

Extremes.gumbelfitpwmMethod
gumbelfitpwm(y::Vector{<:Real})

Estimate the Gumbel parameters with the probability weighted moments.

Extremes.hessianMethod
hessian(model::MaximumLikelihoodAbstractExtremeValueModel)::PDMat{Float64, Matrix{Float64}}

Calculates the Hessian matrix associated with the MaximumLikelihoodAbstractExtremeValueModel model.

Extremes.hisplot_dataFunction
histplot_data(fm::fittedModel)

Return the histogram plot data in a Dictionary.

Extremes.histplotMethod
histplot(fm::AbstractFittedExtremeValueModel)

Histogram plot

Extremes.loglikeMethod
loglike(model::AbstractExtremeValueModel, θ::Vector{<:Real})

Compute the model loglikelihood AbstractExtremeValueModelluated at θ.

Extremes.loglikeMethod
loglike(fd::MaximumLikelihoodAbstractExtremeValueModel)::Real

Compute the model loglikelihood AbstractExtremeValueModelluated at θ̂ if the maximum likelihood method has been used.

Extremes.mergeMethod
merge(c₁::Cluster, c₂::Cluster)

Merge cluster c₁ and c₂ into a single cluster.

Extremes.mrlplot_dataFunction
mrlplot_data(y::Vector{<:Real}, steps::Int = 100)::DataFrame

Compute the mean residual life from vector y.

The set of thresholds ranges from minimum(y) to the second-to-last larger value in steps number of steps.

Extremes.parametervarFunction
parametervar(fm::pwmAbstractExtremeValueModel, nboot::Int=1000)

Estimate the parameter estimates covariance matrix by bootstrap.

Extremes.parametervarMethod
parametervar(fm::BayesianAbstractExtremeValueModel)::Array{Float64, 2}

Compute the covariance parameters estimate of the fitted model fm.

Extremes.parametervarMethod
parametervar(fm::MaximumLikelihoodAbstractExtremeValueModel)::Array{Float64, 2}

Compute the covariance parameters estimate of the fitted model fm.

Extremes.paramindexFunction
paramindex(model::AbstractExtremeValueModel)

Return the postitions corresponding to the location, scale and shape parameter.

Extremes.probplotMethod
probplot(fm::AbstractFittedExtremeValueModel)

Probability plot

Extremes.probplot_dataFunction
probplot_data(fm::fittedModel)

Return the probability plot data in a DataFrame.

Extremes.pwmMethod
pwm(x::Vector{<:Real},p::Int,r::Int,s::Int)

Compute the empirical probability weighted moments.

The probability weighted moments are defined by Landwehr *et al. (1979) as follows:

\[M_{p,r,s} = \mathbb{E}\left[ X^p F^r(X) \left\{ 1-F(X) \right\}^s \right].\]

The unbiased empirical estimates is computed using the formula given by Landwehr *et al. (1979).

Reference

Landwehr, J. M., Matalas, N. C. and Wallis, J. R. (1979). Probability weighted moments compared with some traditional techniques in estimating Gumbel parameters and quantiles. Water Resources Research, 15:1055–1064.

Extremes.qqplot_dataFunction
qqplot_data(fm::fittedModel)

Return the quantile-quantile plot data in a DataFrame.

Extremes.qqplotciFunction
qqplotci(fm::AbstractFittedExtremeValueModel, α::Real=.05)

Quantile-Quantile plot along with the confidence/credible interval of level 1-α.

Note

This function is currently only available for stationary models.

See also returnlevelplotci and qqplot.

Example

using Distributions, Extremes

pd = GeneralizedExtremeValue(0,1,0)
y = rand(pd, 300)
fm = gevfit(y)

qqplotci(fm)
Extremes.quantilAbstractExtremeValueModelrFunction
quantilAbstractExtremeValueModelr(fm::pwmAbstractExtremeValueModel, level::Real, nboot::Int=1000)::Vector{<:Real}

Compute the approximate variance of the quantile of level level from the fitted model fm by bootstrap.

Extremes.quantilevarMethod
quantilevar(fd::MaximumLikelihoodAbstractExtremeValueModel, level::Real)::Vector{<:Real}

Compute the variance of the quantile of level level from the fitted model fm.

Extremes.returnlevelMethod
returnlevel(fm::BayesianAbstractExtremeValueModel{ThresholdExceedance}, threshold::Real, nobservation::Int,
    nobsperblock::Int, returnPeriod::Real)::ReturnLevel

Compute the return level corresponding to the return period returnPeriod from the fitted model fm.

The threshold should be a scalar. A varying threshold is not yet implemented.

Extremes.returnlevelMethod
returnlevel(fm::MaximumLikelihoodAbstractExtremeValueModel{ThresholdExceedance}, threshold::Real, nobservation::Int,
    nobsperblock::Int, returnPeriod::Real)::ReturnLevel

Compute the return level corresponding to the return period returnPeriod from the fitted model fm.

The threshold should be a scalar. A varying threshold is not yet implemented.

Extremes.returnlevelMethod
returnlevel(fm::pwmAbstractExtremeValueModel{ThresholdExceedance, T} where T<:Distribution, threshold::Real, nobservation::Int,
    nobsperblock::Int, returnPeriod::Real)::ReturnLevel

Compute the return level corresponding to the return period returnPeriod from the fitted model fm.

The threshold should be a scalar. A varying threshold is not yet implemented.

Extremes.returnlevelMethod
returnlevel(fm::BayesianAbstractExtremeValueModel{BlockMaxima}, returnPeriod::Real)::ReturnLevel

Compute the return level corresponding to the return period returnPeriod from the fitted model fm.

Extremes.returnlevelMethod
returnlevel(fm::MaximumLikelihoodAbstractExtremeValueModel{BlockMaxima}, returnPeriod::Real)::ReturnLevel

Compute the return level corresponding to the return period returnPeriod from the fitted model fm.

Extremes.returnlevelMethod
returnlevel(fm::pwmAbstractExtremeValueModel{BlockMaxima, T} where T<:Distribution, returnPeriod::Real)::ReturnLevel

Compute the confidence intervel for the return level corresponding to the return period returnPeriod from the fitted model fm with confidence level confidencelevel.

Extremes.returnlevelplotciFunction
returnlevelplotci(fm::AbstractFittedExtremeValueModel, α::Real=.05)

Return level plot along with the confidence/credible interval of level 1-α.

Note

This function is currently only available for stationary models.

See also returnlevelplotci and qqplot.

Example

using Distributions, Extremes

pd = GeneralizedExtremeValue(0,1,0)
y = rand(pd, 300)
fm = gevfit(y)

returnlevelplotci(fm)
Extremes.showAbstractExtremeValueModelMethod
showAbstractExtremeValueModel(io::IO, obj::BlockMaxima{GeneralizedExtremeValue}; prefix::String = "")

Displays a BlockMaxima{GeneralizedExtremeValue} with the prefix prefix before every line.

Extremes.showAbstractExtremeValueModelMethod
showAbstractExtremeValueModel(io::IO, obj::BlockMaxima{Gumbel}; prefix::String = "")

Displays a BlockMaxima{Gumbel} with the prefix prefix before every line.

Extremes.showAbstractExtremeValueModelMethod
showAbstractExtremeValueModel(io::IO, obj::ThresholdExceedance; prefix::String = "")

Displays a ThresholdExceedance with the prefix prefix before every line.

Extremes.showAbstractFittedExtremeValueModelMethod
showAbstractFittedExtremeValueModel(io::IO, obj::BayesianAbstractExtremeValueModel; prefix::String = "")

Displays a BayesianAbstractExtremeValueModel with the prefix prefix before every line.

Extremes.showAbstractFittedExtremeValueModelMethod
showAbstractFittedExtremeValueModel(io::IO, obj::MaximumLikelihoodAbstractExtremeValueModel; prefix::String = "")

Displays a MaximumLikelihoodAbstractExtremeValueModel with the prefix prefix before every line.

Extremes.showAbstractFittedExtremeValueModelMethod
showAbstractFittedExtremeValueModel(io::IO, obj::pwmAbstractExtremeValueModel; prefix::String = "")

Displays a pwmAbstractExtremeValueModel with the prefix prefix before every line.

Extremes.showChainMethod
showChain(io::IO, obj::MambaLite.Chains; prefix::String = "")

Displays a MambaLite.Chains with the prefix prefix before every line.

Extremes.showClusterMethod
showCluster(io::IO, obj::Cluster; prefix::String = " ")

Displays a Cluster with the prefix prefix before every line.

Extremes.showparamfunMethod
showparamfun(name::String, param::paramfun)::String

Constructs a string describing a parameter param with name name.

Extremes.slicematrixMethod
slicematrix(A::AbstractMatrix{T}; dims::Int=1)::Array{Array{T,1},1} where T

Convert a Matrix in a Vector of Vectors.

The slicing dimension is defined with the optional dims keyword.

Examples

julia> A = rand(1:10,5,4)
julia> V₁ = Extremes.slicematrix(A)
julia> V₂ = Extremes.slicematrix(A, dims=2)
Extremes.standardizeMethod
standardize(v::Variable)

Standardize the values of the Variable.

Implementation

The Variable values are standardized by substracting the empirical mean and dividing by the empirical standard deviation. A VariableStd type is returned.

See also Variable, VariableStd and reconstruct.

Extremes.transformMethod
transform(fm::BayesianAbstractExtremeValueModel{BlockMaxima{GeneralizedExtremeValue}})::BayesianAbstractExtremeValueModel

Transform the fitted model for the original covariate scales.

Extremes.transformMethod
transform(fm::BayesianAbstractExtremeValueModel{BlockMaxima{Gumbel}})::BayesianAbstractExtremeValueModel

Transform the fitted model for the original covariate scales.

Extremes.transformMethod
transform(fm::BayesianAbstractExtremeValueModel{ThresholdExceedance})

Transform the fitted model for the original covariate scales.

Extremes.transformMethod
transform(fm::MaximumLikelihoodAbstractExtremeValueModel{BlockMaxima})::MaximumLikelihoodAbstractExtremeValueModel

Transform the fitted model for the original covariate scales.

Extremes.transformMethod
transform(fm::MaximumLikelihoodAbstractExtremeValueModel{BlockMaxima})::MaximumLikelihoodAbstractExtremeValueModel

Transform the fitted model for the original covariate scales.

Extremes.transformMethod
transform(fm::MaximumLikelihoodAbstractExtremeValueModel{ThresholdExceedance})

Transform the fitted model for the original covariate scales.

Extremes.unslicematrixMethod
unslicematrix(B::Array{Array{T,1},1}; dims::Int=1)::AbstractMatrix{T} where T

Convert a Vector of Vectors in a Matrix if all vectors share the same length.

The concatenation dimension is defined with the optional dims keyword.

Examples

julia> V = Extremes.slicematrix(rand(1:10, 5, 4))
julia> A₁ = Extremes.unslicematrix(V)
julia> A₂ = Extremes.unslicematrix(V, dims=2)
Extremes.validatelengthMethod
validatelength(length::Real, explvariables::Vector{Variable})

Validate that the explanatory variables are of length length.

Extremes.validatestationarityMethod
validatestationarity(model::T)::T where T<:AbstractExtremeValueModel

Throw warning if the model is nonstationary.

Statistics.quantileMethod
quantile(model::AbstractExtremeValueModel, θ::Vector{<:Real}, p::Real)

Compute the quantile of level p from the model AbstractExtremeValueModelluated at `θ"".

If the model is non-stationary, then the effective quantiles are returned, i.e. one for each covariate value.

Statistics.quantileMethod
quantile(fm::BayesianAbstractExtremeValueModel,p::Real)::Real

Compute the quantile of level p from the fitted Bayesian model fm.

If the model is stationary, then a quantile is returned for each MCMC steps.

If the model is non-stationary, a matrix of quantiles is returned, where each row corresponds to a MCMC step and each column to a covariate.

Statistics.quantileMethod
quantile(fm::MaximumLikelihoodAbstractExtremeValueModel, p::Real)::Vector{<:Real}

Compute the quantile of level p from the fitted model by maximum likelihood. In the case of non-stationarity, the effective quantiles are returned.

Statistics.quantileMethod
quantile(fm::pwmAbstractExtremeValueModel, p::Real)::Vector{<:Real}

Compute the quantile of level p from the fitted model. For the probability weighted moment method, the model has to be stationary.

StatsAPI.paramsMethod
params(fm::AbstractFittedExtremeValueModel)

Return the parameters of the fitted model.