List of useful methods
The following lists docstrings for various important functions.
Once a DataModel
object has been defined, it can subsequently be used to compute various quantities as follows:
StatsAPI.loglikelihood
— Methodloglikelihood(DM::DataModel, θ::AbstractVector) -> Real
Calculates the logarithm of the likelihood $L$, i.e. $\ell(\mathrm{data} \, | \, \theta) \coloneqq \mathrm{ln} \big( L(\mathrm{data} \, | \, \theta) \big)$ given a DataModel
and a parameter configuration $\theta$.
InformationGeometry.MLE
— MethodMLE(DM::DataModel) -> Vector
Returns the parameter configuration $\theta_\text{MLE} \in \mathcal{M}$ which is estimated to have the highest likelihood of producing the observed data (under the assumption that the specified model captures the true relationship present in the data). For performance reasons, the maximum likelihood estimate is stored as a part of the DataModel
type.
InformationGeometry.MLEuncert
— MethodMLEuncert(DM::AbstractDataModel, mle::AbstractVector=MLE(DM), F::AbstractMatrix=FisherMetric(DM, mle))
Returns vector of type Measurements.Measurement
where the parameter uncertainties are approximated via the diagonal of the inverse Fisher metric. That is, the stated uncertainties are a linearized symmetric approximation of the true parameter uncertainties around the MLE.
InformationGeometry.LogLikeMLE
— MethodLogLikeMLE(DM::DataModel) -> Real
Returns the value of the log-likelihood $\ell$ when evaluated at the maximum likelihood estimate, i.e. $\ell(\mathrm{data} \, | \, \theta_\text{MLE})$. For performance reasons, this value is stored as a part of the DataModel
type.
Various geometric quantities which are intrinsic to the parameter manifold $\mathcal{M}$ can be computed as a result of the Fisher metric $g$ (and subsequent choice of the Levi-Civita connection) such as the Riemann and Ricci tensors and the Ricci scalar $R$.
InformationGeometry.Score
— MethodScore(DM::DataModel, θ::AbstractVector{<:Number}; ADmode::Val=Val(:ForwardDiff))
Calculates the gradient of the log-likelihood $\ell$ with respect to a set of parameters $\theta$. ADmode=Val(false)
computes the Score using the Jacobian dmodel
provided in DM
, i.e. by having to separately evaluate both the model
as well as dmodel
. Other choices of ADmode
directly compute the Score by differentiating the formula the log-likelihood, i.e. only one evaluation on a dual variable is performed.
InformationGeometry.FisherMetric
— MethodFisherMetric(DM::DataModel, θ::AbstractVector{<:Number})
Computes the Fisher metric $g$ given a DataModel
and a parameter configuration $\theta$ under the assumption that the likelihood $L(\mathrm{data} \, | \, \theta)$ is a multivariate normal distribution.
\[g_{ab}(\theta) \coloneqq -\int_{\mathcal{D}} \mathrm{d}^m y_{\mathrm{data}} \, L(y_{\mathrm{data}} \,|\, \theta) \, \frac{\partial^2 \, \mathrm{ln}(L)}{\partial \theta^a \, \partial \theta^b} = -\mathbb{E} \bigg( \frac{\partial^2 \, \mathrm{ln}(L)}{\partial \theta^a \, \partial \theta^b} \bigg)\]
InformationGeometry.GeometricDensity
— MethodGeometricDensity(DM::AbstractDataModel, θ::AbstractVector) -> Real
Computes the square root of the determinant of the Fisher metric $\sqrt{\mathrm{det}\big(g(\theta)\big)}$ at the point $\theta$.
InformationGeometry.ChristoffelSymbol
— MethodChristoffelSymbol(DM::DataModel, θ::AbstractVector; BigCalc::Bool=false)
ChristoffelSymbol(Metric::Function, θ::AbstractVector; BigCalc::Bool=false)
Calculates the components of the $(1,2)$ Christoffel symbol $\Gamma$ at a point $\theta$ (i.e. the Christoffel symbol "of the second kind") through finite differencing of the Metric
. Accurate to ≈ 3e-11. BigCalc=true
increases accuracy through BigFloat
calculation.
InformationGeometry.Riemann
— MethodRiemann(DM::DataModel, θ::AbstractVector; BigCalc::Bool=false)
Riemann(Metric::Function, θ::AbstractVector; BigCalc::Bool=false)
Calculates the components of the $(1,3)$ Riemann tensor by finite differencing of the Metric
. BigCalc=true
increases accuracy through BigFloat calculation.
InformationGeometry.Ricci
— MethodRicci(DM::DataModel, θ::AbstractVector; BigCalc::Bool=false)
Ricci(Metric::Function, θ::AbstractVector; BigCalc::Bool=false)
Calculates the components of the $(0,2)$ Ricci tensor by finite differencing of the Metric
. BigCalc=true
increases accuracy through BigFloat
calculation.
InformationGeometry.RicciScalar
— MethodRicciScalar(DM::DataModel, θ::AbstractVector; BigCalc::Bool=false) -> Real
RicciScalar(Metric::Function, θ::AbstractVector; BigCalc::Bool=false) -> Real
Calculates the Ricci scalar by finite differencing of the Metric
. BigCalc=true
increases accuracy through BigFloat
calculation.
Further, studying the geodesics associated with a metric manifold can yield insights into its geometry.
InformationGeometry.GeodesicDistance
— MethodGeodesicDistance(DM::DataModel, P::AbstractVector{<:Number}, Q::AbstractVector{<:Number}; tol::Real=1e-10)
GeodesicDistance(Metric::Function, P::AbstractVector{<:Number}, Q::AbstractVector{<:Number}; tol::Real=1e-10)
Computes the length of a geodesic connecting the points P
and Q
.
InformationGeometry.AIC
— MethodAIC(DM::DataModel, θ::AbstractVector) -> Real
Calculates the Akaike Information Criterion given a parameter configuration $\theta$ defined by $\mathrm{AIC} = 2 \, \mathrm{length}(\theta) -2 \, \ell(\mathrm{data} \, | \, \theta)$. Lower values for the AIC indicate that the associated model function is more likely to be correct. For linearly parametrized models and small sample sizes, it is advisable to instead use the AICc which is more accurate.
InformationGeometry.AICc
— MethodAICc(DM::DataModel, θ::AbstractVector) -> Real
Computes Akaike Information Criterion with an added correction term that prevents the AIC from selecting models with too many parameters (i.e. overfitting) in the case of small sample sizes. $\mathrm{AICc} = \mathrm{AIC} + \frac{2\mathrm{length}(\theta)^2 + 2 \mathrm{length}(\theta)}{N - \mathrm{length}(\theta) - 1}$ where $N$ is the number of data points. Whereas AIC constitutes a first order estimate of the information loss, the AICc constitutes a second order estimate. However, this particular correction term assumes that the model is linearly parametrized.
InformationGeometry.BIC
— MethodBIC(DM::DataModel, θ::AbstractVector) -> Real
Calculates the Bayesian Information Criterion given a parameter configuration $\theta$ defined by $\mathrm{BIC} = \mathrm{ln}(N) \cdot \mathrm{length}(\theta) -2 \, \ell(\mathrm{data} \, | \, \theta)$ where $N$ is the number of data points.
InformationGeometry.IsLinearParameter
— FunctionIsLinearParameter(DM::DataModel, MLE::AbstractVector) -> BitVector
Checks with respect to which parameters the model function model(x,θ)
is linear and returns vector of booleans where true
indicates linearity. This test is performed by comparing the Jacobians of the model for two random configurations $\theta_1, \theta_2 \in \mathcal{M}$ column by column.
InformationGeometry.ConfidenceRegionVolume
— FunctionConfidenceRegionVolume(DM::AbstractDataModel, Confnum::Real; N::Int=Int(1e5), WE::Bool=true, Approx::Bool=false, kwargs...) -> Real
Computes coordinate-invariant volume of confidence region associated with level Confnum
via Monte Carlo by integrating the geometric density factor. For likelihoods which are particularly expensive to evaluate, Approx=true
can improve the performance by approximating the confidence region via polygons.
InformationGeometry.Pullback
— FunctionPullback(DM::AbstractDataModel, ω::AbstractVector{<:Number}, θ::AbstractVector) -> Vector
Pull-back of a covector to the parameter manifold $T^*\mathcal{M} \longleftarrow T^*\mathcal{D}$.
Pullback(DM::DataModel, G::AbstractArray{<:Number,2}, θ::AbstractVector) -> Matrix
Pull-back of a (0,2)-tensor G
to the parameter manifold.
InformationGeometry.Pushforward
— FunctionPushforward(DM::DataModel, X::AbstractVector, θ::AbstractVector) -> Vector
Calculates the push-forward of a vector X
from the parameter manifold to the data space $T\mathcal{M} \longrightarrow T\mathcal{D}$.
In many applied settings, one often does not have a dataset of sufficient size for all parameters in the model to be "practically identifiable", which means that bounded confidence regions may only exist for very low confidence levels (e.g. up to $0.1\sigma$). In such cases, it is still possible to compute radial geodesics emanating from the MLE to study the geometry of the parameter space.
A slightly more robust alternative to using geodesics is given by the so-called profile likelihood method. Essentially, it consists of pinning one of the parameters at particular values on a grid, while optimizing the remaining parameters to maximize the likelihood function at every step. Ultimately, one ends up with one-dimensional slices of the parameter manifold along which the likelihood decays most slowly.
Missing docstring for ParameterProfiles(::DataModel,::Float64)
. Check Documenter's build log for details.
InformationGeometry.ProfileBox
— FunctionProfileBox(DM::AbstractDataModel, Fs::AbstractVector{<:AbstractInterpolation}, Confnum::Real=2.) -> HyperCube
Constructs HyperCube
which bounds the confidence region associated with the confidence level Confnum
from the interpolated likelihood profiles.