ParetoSmooth

Documentation for ParetoSmooth.

ParetoSmooth.ModelComparisonType
ModelComparison

A struct containing the results of model comparison.

Fields

  • pointwise::KeyedArray:
  • estimates::KeyedArray: A table containing the results of model comparison, with the following columns –
    • cv_est: The difference in total leave-one-out cross validation scores between models.
    • se_loo_diff: The standard error for each row of cv_est.
    • weight: A set of Akaike-like weights assigned to each model, which can be used in pseudo-Bayesian model averaging.

Example

┌───────┬──────────┬───────┬────────┐
│       │   cv_est │ se_cv │ weight │
├───────┼──────────┼───────┼────────┤
│ m5_1t │     0.00 │  0.00 │   0.67 │
│ m5_3t │    -0.69 │  0.42 │   0.33 │
│ m5_2t │    -6.68 │  4.74 │   0.00 │
└───────┴──────────┴───────┴────────┘

See also: PsisLoo

ParetoSmooth.PsisType
Psis{V<:AbstractVector{F},I<:Integer} where {F<:Real}

A struct containing the results of Pareto-smoothed importance sampling.

Fields

  • weights: A vector of smoothed, truncated, and normalized importance sampling weights.
  • pareto_k: Estimates of the shape parameter k of the generalized Pareto distribution.
  • ess: Estimated effective sample size for each LOO evaluation.
  • tail_len: Vector indicating how large the "tail" is for each observation.
  • dims: Named tuple of length 2 containing s (posterior sample size) and n (number of observations).
ParetoSmooth.PsisLooType
PsisLoo{
    F <: Real,
    AF <: AbstractArray{F},
    VF <: AbstractVector{F},
    I <: Integer,
    VI <: AbstractVector{I},
} <: AbstractCV

A struct containing the results of jackknife (leave-one-out) cross validation using Pareto smoothed importance sampling.

Fields

  • estimates::KeyedArray: A KeyedArray with columns :total, :se_total, :mean, :se_mean, and rows :cv_est, :naive_est, :overfit. See # Extended help for more.
    • :cv_est contains estimates for the out-of-sample prediction error, as predicted using the jackknife (LOO-CV).
    • :naive_est contains estimates of the in-sample prediction error.
    • :overfit is the difference between the previous two estimators, and estimates the amount of overfitting. When using the log probability score, it is equal to the effective number of parameters – a model with an overfit of 2 is "about as overfit" as a model with 2 independent parameters that have a flat prior.
  • pointwise::KeyedArray: A KeyedArray of pointwise estimates with 5 columns –
    • :cv_est contains the estimated out-of-sample error for this point, as measured
    using leave-one-out cross validation.
    • :naive_est contains the in-sample estimate of error for this point.
    • :overfit is the difference in the two previous estimates.
    • :ess is the effective sample size, which estimates the simulation error caused by using Monte Carlo estimates. It does not measure the accuracy of predictions.
    • :pareto_k is the estimated value for the parameter ξ of the generalized Pareto distribution. Values above .7 indicate that PSIS has failed to approximate the true distribution.
  • psis_object::Psis: A Psis object containing the results of Pareto-smoothed importance sampling.

Extended help

The total score depends on the sample size, and summarizes the weight of evidence for or against a model. Total scores are on an interval scale, meaning that only differences of scores are meaningful. It is not possible to interpret a total score by looking at it. The total score is not a relative goodness-of-fit statistic (for this, see the average score).

The overfit is equal to the difference between the in-sample and out-of-sample predictive accuracy. When using the log probability score, it is equal to the "effective number of parameters" – a model with an overfit of 2 is "about as overfit" as a model with 2 free parameters and flat priors.

The average score is the total score, divided by the sample size. It estimates the expected log score, i.e. the expectation of the log probability density of observing the next point. The average score is a relative goodness-of-fit statistic which does not depend on sample size.

Unlike for chi-square goodness of fit tests, models do not have to be nested for model comparison using cross-validation methods.

See also: [loo]@ref, [bayes_cv]@ref, [psis_loo]@ref, [Psis]@ref

ParetoSmooth.PsisLooMethodType
PsisLooMethod

Use Pareto-smoothed importance sampling together with leave-one-out cross validation to estimate the out-of-sample predictive accuracy.

ParetoSmooth._def_tail_lengthMethod
_def_tail_length(log_ratios::AbstractVector, r_eff::Real) -> tail_len::Integer

Define the tail length as in Vehtari et al. (2019).

ParetoSmooth._do_psis_i!Method
_do_psis_i!(is_ratios::AbstractVector{Real}, tail_length::Integer) -> T

Do PSIS on a single vector, smoothing its tail values.

Arguments

  • is_ratios::AbstractVector{Real}: A vector of importance sampling ratios, scaled to have a maximum of 1.

Returns

  • T<:Real: ξ, the shape parameter for the GPD; big numbers indicate thick tails.

Extended help

Additional information can be found in the LOO package from R.

ParetoSmooth._psis_smooth_tail!Method
_psis_smooth_tail!(tail::AbstractVector{T}, cutoff::T) where {T<:Real} -> ξ::T

Takes an already sorted vector of observations from the tail and smooths it in place with PSIS before returning shape parameter ξ.

ParetoSmooth.gpd_quantileMethod
gpd_quantile(p::T, k::T, sigma::T) where {T<:Real} -> T

Compute the p quantile of the Generalized Pareto Distribution (GPD).

Arguments

  • p: A scalar between 0 and 1.
  • ξ: A scalar shape parameter.
  • σ: A scalar scale parameter.

Returns

A quantile of the Generalized Pareto Distribution.

ParetoSmooth.gpdfitMethod
gpdfit(
    sample::AbstractVector{T<:Real}; 
    wip::Bool=true, 
    min_grid_pts::Integer=30, 
    sort_sample::Bool=false
) -> (ξ::T, σ::T)

Return a named list of estimates for the parameters ξ (shape) and σ (scale) of the generalized Pareto distribution (GPD), assuming the location parameter is 0.

Arguments

  • sample::AbstractVector: A numeric vector. The sample from which to estimate the parameters.
  • wip::Bool = true: Logical indicating whether to adjust ξ based on a weakly informative Gaussian prior centered on 0.5. Defaults to true.
  • min_grid_pts::Integer = 30: The minimum number of grid points used in the fitting algorithm. The actual number used is min_grid_pts + ⌊sqrt(length(sample))⌋.

Note

Estimation method taken from Zhang, J. and Stephens, M.A. (2009). The parameter ξ is the negative of k.

ParetoSmooth.looMethod
function loo(args...; method=PsisLooMethod(), kwargs...) -> PsisLoo

Compute the approximate leave-one-out cross-validation score using the specified method.

Currently, this function only serves to call psis_loo, but this could change in the future. The default methods or return type may change without warning; thus, we recommend using psis_loo instead if reproducibility is required.

See also: psis_loo, PsisLoo.

ParetoSmooth.loo_compareMethod
function loo_compare(
    cv_results::PsisLoo...;
    sort_models::Bool=true,
    [, model_names::Tuple{Symbol}]
) -> ModelComparison

Construct a model comparison table from several PsisLoo objects.

Arguments

- `cv_results`: One or more [`PsisLoo`](@ref) objects to be compared. Alternatively,
a tuple or named tuple of `PsisLoo` objects can be passed. If a named tuple is passed,
these names will be used to label each model. 
- `model_names`: A vector or tuple of strings or symbols used to identify models. If
none, models are numbered using the order of the arguments.
- `sort_models=true`: Sort models by total score (best first).

See also: ModelComparison, PsisLoo, psis_loo

ParetoSmooth.pointwise_log_likelihoodsMethod
pointwise_log_likelihoods(
    ll_fun::Function, 
    samples::AbstractArray{<:Real,3}, 
    data;
    splat::Bool=true
)

Compute the pointwise log likelihood.

Arguments

  • ll_fun::Function: A function taking a single data point and returning the log-likelihood

of that point. This function must take the form f(θ[1], ..., θ[n], data), where θ is the parameter vector. See also the splat keyword argument.

  • samples::AbstractArray: A three dimensional array of MCMC samples. Here, the first dimension should indicate the iteration of the MCMC ; the second dimension should indicate the parameter ; and the third dimension represents the chains.
  • data: A vector of data used to estimate the parameters of the model.
  • splat: If true (default), f must be a function of n different parameters. Otherwise, f is assumed to be a function of a single parameter vector.

Returns

  • Array: A three dimensional array of pointwise log-likelihoods.
ParetoSmooth.psisMethod
psis(
    log_ratios::AbstractArray{T<:Real}, 
    r_eff::AbstractVector; 
    source::String="mcmc", 
    log_weights::Bool=false
) -> Psis

Implements Pareto-smoothed importance sampling (PSIS).

Arguments

Positional Arguments

  • log_ratios::AbstractArray: A 2d or 3d array of importance ratios on the log scale (for PSIS-LOO these are negative log-likelihood values). Indices must be ordered as [data, step, chain]: log_ratios[1, 2, 3] should be the log-likelihood of the first data point, evaluated at the second step in the third chain. Chain indices can be left off if there is only one chain, or if keyword argument chain_index is provided.
  • r_eff::AbstractArray{T}: An (optional) vector of relative effective sample sizes used

in ESS calculations. If left empty, calculated automatically using the FFTESS method from InferenceDiagnostics.jl. See relative_eff to calculate these values.

Keyword Arguments

  • chain_index::Vector: An optional vector of integers specifying which chain each step belongs to. For instance, chain_index[step] should return 2 if log_likelihood[:, step] belongs to the second chain.

  • source::String="mcmc": A string or symbol describing the source of the sample being used. If "mcmc", adjusts ESS for autocorrelation. Otherwise, samples are assumed to be independent. Currently permitted values are ["mcmc", "vi", "other"].

  • log_weights::Bool=false: Return the log weights, rather than the PSIS weights.

See also: [relative_eff]@ref, [psis_loo]@ref, [psis_ess]@ref.

ParetoSmooth.psis_essMethod
function psis_ess(
    weights::AbstractVector{T},
    r_eff::AbstractVector{T}
) -> AbstractVector{T}

Calculate the (approximate) effective sample size of a PSIS sample, using the correction in Vehtari et al. 2019.

Arguments

  • weights: A set of importance sampling weights derived from PSIS.
  • r_eff: The relative efficiency of the MCMC chains from which PSIS samples were derived.

See ?relative_eff to calculate r_eff.

ParetoSmooth.psis_looMethod
function psis_loo(
    log_likelihood::Array{Real} [, args...];
    [, chain_index::Vector{Integer}, kwargs...]
) -> PsisLoo

Use Pareto-Smoothed Importance Sampling to calculate the leave-one-out cross validation score.

Arguments

  • log_likelihood::Array: An array or matrix of log-likelihood values indexed as [data, step, chain]. The chain argument can be left off if chain_index is provided or if all posterior samples were drawn from a single chain.
  • args...: Positional arguments to be passed to psis.
  • chain_index::Vector: An optional vector of integers specifying which chain each step belongs to. For instance, chain_index[3] should return 2 if log_likelihood[:, 3] belongs to the second chain.
  • kwargs...: Keyword arguments to be passed to psis.

See also: psis, loo, PsisLoo.

ParetoSmooth.relative_effMethod
relative_eff(
    sample::AbstractArray{Real, 3}; 
    method=MCMCDiagnosticTools.FFTESSMethod()
)

Calculate the relative efficiency of an MCMC chain, i.e. the effective sample size divided by the nominal sample size.