ParetoSmooth
Documentation for ParetoSmooth.
ParetoSmooth.ModelComparison
ParetoSmooth.Psis
ParetoSmooth.PsisLoo
ParetoSmooth.loo
ParetoSmooth.loo_compare
ParetoSmooth.loo_from_psis
ParetoSmooth.naive_lpd
ParetoSmooth.pointwise_log_likelihoods
ParetoSmooth.psis
ParetoSmooth.psis!
ParetoSmooth.psis_ess
ParetoSmooth.psis_loo
ParetoSmooth.relative_eff
ParetoSmooth.sup_ess
ParetoSmooth.ModelComparison
— TypeModelComparison
A struct containing the results of model comparison.
Fields
pointwise::KeyedArray
: AKeyedArray
of pointwise estimates. See [PsisLoo
]@ref.estimates::KeyedArray
: A table containing the results of model comparison, with the following columns –cv_elpd
: The difference in total leave-one-out cross validation scores between models.cv_avg
: The difference in average LOO-CV scores between models.weight
: A set of Akaike-like weights assigned to each model, which can be used in pseudo-Bayesian model averaging.
std_err::NamedTuple
: A named tuple containing the standard error ofcv_elpd
. Note that these estimators (incorrectly) assume all folds are independent, despite their substantial overlap, which creates a downward biased estimator. LOO-CV differences are not asymptotically normal, so these standard errors cannot be used to calculate a confidence interval.gmpd::NamedTuple
: The geometric mean of the predictive distribution. It equals the geometric mean of the probability assigned to each data point by the model, that is,exp(cv_avg)
. This measure is only meaningful for classifiers (variables with discrete outcomes). We can think of it as measuring how often the model was right: A model that always predicts incorrectly will have a GMPD of 0, while a model that always predicts correctly will have a GMPD of 1. However, the GMPD gives a model "Partial points" between 0 and 1 whenever the model assigns a probability other than 0 or 1 to the outcome that actually happened.
See also: PsisLoo
ParetoSmooth.Psis
— TypePsis{R<:Real, AT<:AbstractArray{R, 3}, VT<:AbstractVector{R}}
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 parameterk
of the generalized Pareto distribution.ess
: Estimated effective sample size for each LOO evaluation, based on the variance of the weights.sup_ess
: Estimated effective sample size for each LOO evaluation, based on the supremum norm, i.e. the size of the largest weight. More likely thaness
to warn when importance sampling has failed. However, it can have a high variance.r_eff
: The relative efficiency of the MCMC chain, i.e. ESS / posterior sample size.tail_len
: Vector indicating how large the "tail" is for each observation.posterior_sample_size
: How many draws from an MCMC chain were used for PSIS.data_size
: How many data points were used for PSIS.
ParetoSmooth.PsisLoo
— TypePsisLoo <: AbstractCV
A struct containing the results of leave-one-out cross validation computed with Pareto smoothed importance sampling.
Fields
estimates::KeyedArray
: A KeyedArray with columns:total, :se_total, :mean, :se_mean
, and rows:cv_elpd, :naive_lpd, :p_eff
. See# Extended help
for more.:cv_elpd
contains estimates for the out-of-sample prediction error, as estimated using leave-one-out cross validation.:naive_lpd
contains estimates of the in-sample prediction error.:p_eff
is the effective number of parameters – a model with ap_eff
of 2 is "about as overfit" as a model with 2 parameters and no regularization.
pointwise::KeyedArray
: AKeyedArray
of pointwise estimates with 5 columns –:cv_elpd
contains the estimated out-of-sample error for this point, as measured
:naive_lpd
contains the in-sample estimate of error for this point.:p_eff
is the difference in the two previous estimates.:ess
is the L2 effective sample size, which estimates the simulation error caused by using Monte Carlo estimates. It does not measure model performance.:inf_ess
is the supremum-based effective sample size, which estimates the simulation error caused by using Monte Carlo estimates. It is more robust than:ess
and should therefore be preferred. It does not measure model performance.: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
: APsis
object containing the results of Pareto-smoothed importance sampling.gmpd
: The geometric mean of the predictive density. It is defined as the geometric mean of the probability assigned to each data point by the model, i.e.exp(cv_avg)
. This measure is only interpretable for classifiers (variables with discrete outcomes). We can think of it as measuring how often the model was right: A model that always predicts incorrectly will have a GMPD of 0, while a model that always predicts correctly will have a GMPD of 1. However, the GMPD gives a model "Partial points" between 0 and 1 whenever the model assigns a probability other than 0 or 1 to the outcome that actually happened, making it a fully Bayesian measure of model quality.mcse
: A float containing the estimated Monte Carlo standard error for the total cross-validation estimate.
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 goodness-of-fit statistic (for this, see the average score).
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.loo
— Methodfunction loo(args...; kwargs...) -> PsisLoo
Compute an approximate leave-one-out cross-validation score.
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, so we recommend using psis_loo
instead if reproducibility is required.
ParetoSmooth.loo_compare
— Methodfunction loo_compare(
cv_results...;
sort_models::Bool=true,
best_to_worst::Bool=true,
[, model_names::Tuple{Symbol}]
) -> ModelComparison
Construct a model comparison table from several PsisLoo
objects.
Arguments
cv_results
: One or morePsisLoo
objects to be compared. Alternatively, a tuple or named tuple ofPsisLoo
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
: Sort models by total score.high_to_low
: Sort models from best to worst score. Iffalse
, reverse the order.
See also: ModelComparison
, PsisLoo
, psis_loo
ParetoSmooth.loo_from_psis
— Methodloo_from_psis(
log_likelihood::AbstractArray{<:Real}, psis_object::Psis;
chain_index::Vector{<:Integer}
)
Use a precalculated Psis
object to estimate the leave-one-out cross validation score.
Arguments
log_likelihood::Array
: A matrix or 3d array 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.
psis_object
: A precomputedPsis
object used to estimate the LOO-CV score.chain_index::Vector{Int}
: 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.
ParetoSmooth.pointwise_log_likelihoods
— Methodpointwise_log_likelihoods(
ll_fun::Function, samples::AbstractArray{<:Real,3}, data;
splat::Bool=true[, chain_index::Vector{<:Integer}]
)
Compute the pointwise log likelihoods.
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 step of the MCMC algorithm; the second dimension should indicate the parameter; and the third should indicate the chain.data
: An array of data points used to estimate the parameters of the model.splat
: Iftrue
(default),f
must be a function ofn
different parameters. Otherwise,f
is assumed to be a function of a single parameter vector.chain_index::Vector{Int}
: 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.
Returns
Array
: A three dimensional array of pointwise log-likelihoods.
ParetoSmooth.psis!
— Methodpsis!(
is_ratios::AbstractVector{<:Real};
tail_length::Integer, log_weights::Bool=true
) -> Real
Do PSIS on a single vector, smoothing its tail values in place before returning the estimated shape constant for the pareto_k
distribution. This does not normalize the log-weights.
Arguments
is_ratios::AbstractVector{<:Real}
: A vector of importance sampling ratios, scaled to have a maximum of 1.r_eff::AbstractVector{<:Real}
: The relative effective sample size, used to calculate the effective sample size. See [rel_eff]@ref for more information.log_weights::Bool
: A boolean indicating whether the input vector is a vector of log ratios, rather than raw importance sampling ratios.
Returns
Real
: ξ, the shape parameter for the GPD. Bigger numbers indicate thicker tails.
Notes
Unlike the methods for arrays, psis!
performs no checks to make sure the input values are valid.
ParetoSmooth.psis
— Methodpsis(
log_ratios::AbstractArray{T<:Real},
r_eff::AbstractVector{T};
source::String="mcmc"
) -> Psis
Implements Pareto-smoothed importance sampling (PSIS).
Arguments
Positional Arguments
log_ratios::AbstractArray
: A 2d or 3d array of (unnormalized) importance ratios on the log scale. Indices must be ordered as[data, step, chain]
. The chain index can be left off if there is only one chain, or if keyword argumentchain_index
is provided.r_eff::AbstractVector
: 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{Int}
: 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"].calc_ess::Bool=true
: Iffalse
, do not calculate ESS diagnostics. Attempting to access ESS diagnostics will return an empty array.checks::Bool=true
: Iftrue
, check inputs for possible errors. Disabling will improve performance slightly.
See also: [relative_eff
]@ref, [psis_loo
]@ref, [psis_ess
]@ref.
ParetoSmooth.psis_ess
— Methodfunction psis_ess(
weights::AbstractVector{T<:Real},
r_eff::AbstractVector{T}
) -> AbstractVector{T}
Calculate the (approximate) effective sample size of a PSIS sample, using the correction in Vehtari et al. 2019. This uses the entropy-based definition of ESS, measuring the K-L divergence of the proposal and target distributions.
Arguments
weights
: A set of normalized 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_loo
— Methodfunction psis_loo(
log_likelihood::AbstractArray{<:Real} [, args...];
[, chain_index::Vector{Int}, kwargs...]
) -> PsisLoo
Use Pareto-Smoothed Importance Sampling to calculate the leave-one-out cross validation score.
Arguments
log_likelihood::Array
: A matrix or 3d array 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 topsis
.chain_index::Vector{Int}
: 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.
kwargs...
: Keyword arguments to be passed topsis
.
ParetoSmooth.relative_eff
— Methodrelative_eff(
sample::AbstractArray{<:Real, 3};
source::Union{AbstractString, Symbol} = "default",
maxlag::Int = typemax(Int),
kwargs...,
)
Calculate the relative efficiency of an MCMC chain, i.e., the effective sample size divided by the nominal sample size.
If lowercase(String(source))
is "default"
or "mcmc"
, the relative effective sample size is computed with MCMCDiagnosticTools.ess
, using keyword arguments kind = :basic
, maxlag = maxlag
, and the remaining keyword arguments kwargs...
. Otherwise a vector of ones for each chain is returned.
Arguments
sample::AbstractArray{<:Real, 3}
: An array of log-likelihood values of the shape(parameters, draws, chains)
.
ParetoSmooth.sup_ess
— Methodfunction sup_ess(
weights::AbstractMatrix{T},
r_eff::AbstractVector{T}
) -> AbstractVector
Calculate the supremum-based effective sample size of a PSIS sample, i.e. the inverse of the maximum weight. This measure is more sensitive than the ess
from psis_ess
, but also much more variable. It uses the L-∞ norm.
Arguments
weights
: A set of importance sampling weights derived from PSIS.r_eff
: The relative efficiency of the MCMC chains; see also [relative_eff
]@ref.
ParetoSmooth.naive_lpd
— Functionnaive_lpd(log_likelihood::AbstractArray{<:Real}[, chain_index])
Calculate the naive (in-sample) estimate of the expected log probability density, otherwise known as the in-sample Bayes score. This method yields heavily biased results, and we advise against using it; it is included only for pedagogical purposes.
This method is unexported and can only be accessed by calling ParetoSmooth.naive_lpd
.
Arguments
log_likelihood::Array
: A matrix or 3d array 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.
chain_index::Vector{Int}
: 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.