GMRFs
GaussianMarkovRandomFields.AbstractGMRF Type
AbstractGMRF
A Gaussian Markov Random Field (GMRF) is a special case of a multivariate normal distribution where the precision matrix is sparse. The zero entries in the precision correspond to conditional independencies.
sourceGaussianMarkovRandomFields.GMRF Type
GMRF(mean, precision, alg=LinearSolve.DefaultLinearSolver(); Q_sqrt=nothing, rbmc_strategy=RBMCStrategy(1000), linsolve_cache=nothing)
A Gaussian Markov Random Field with mean mean
and precision matrix precision
.
Arguments
mean::AbstractVector
: The mean vector of the GMRF.precision::Union{LinearMap, AbstractMatrix}
: The precision matrix (inverse covariance) of the GMRF.alg
: LinearSolve algorithm to use for linear system solving. Defaults toLinearSolve.DefaultLinearSolver()
.Q_sqrt::Union{Nothing, AbstractMatrix}
: Square root of precision matrix Q, used for sampling when algorithm doesn't support backward solve.rbmc_strategy
: RBMC strategy for marginal variance computation when selected inversion is unavailable. Defaults toRBMCStrategy(1000)
.linsolve_cache::Union{Nothing, LinearSolve.LinearCache}
: Existing LinearSolve cache to reuse. Ifnothing
, creates a new cache. Useful for iterative algorithms requiring factorization reuse.
Type Parameters
T<:Real
: The numeric type (e.g., Float64).PrecisionMap<:Union{LinearMap{T}, AbstractMatrix{T}}
: The type of the precision matrix.
Fields
mean::Vector{T}
: The mean vector.precision::PrecisionMap
: The precision matrix.Q_sqrt::Union{Nothing, AbstractMatrix{T}}
: Square root of precision matrix for sampling.linsolve_cache::LinearSolve.LinearCache
: The LinearSolve cache for efficient operations.rbmc_strategy
: RBMC strategy for variance computation fallback.
Notes
The LinearSolve cache is constructed automatically (if not provided) and is used to compute means, variances, samples, and other GMRF quantities efficiently. The algorithm choice determines which optimization strategies (selected inversion, backward solve) are available. When selected inversion is not supported, marginal variances are computed using the configured RBMC strategy. Providing an existing linsolve_cache
enables factorization reuse in iterative algorithms.
GaussianMarkovRandomFields.precision_map Function
precision_map(::AbstractGMRF)
Return the precision (inverse covariance) map of the GMRF.
sourceprecision_map(d::ConstrainedGMRF)
Return the precision map of the constrained GMRF. Note: This is singular due to the constraints, but we return it for interface compliance. In practice, this should rarely be used directly due to singularity.
sourceGaussianMarkovRandomFields.precision_matrix Function
precision_matrix(::AbstractGMRF)
Return the precision (inverse covariance) matrix of the GMRF.
sourceprecision_matrix(model::LatentModel; kwargs...)
Construct the precision matrix for the model given hyperparameter values.
Arguments
model
: The LatentModel instancekwargs...
: Hyperparameter values as keyword arguments
Returns
A precision matrix (AbstractMatrix or LinearMap) for use in GMRF construction.
sourceGaussianMarkovRandomFields.InformationVector Type
InformationVector(data::AbstractVector)
Wrapper type for information vectors (Q * μ) used in GMRF construction. This allows distinguishing between constructors that take mean vectors vs information vectors.
sourceGaussianMarkovRandomFields.information_vector Function
information_vector(d::GMRF)
Return the information vector (Q * μ) for the GMRF. If stored, returns the cached value; otherwise computes it.
sourceMetadata
GaussianMarkovRandomFields.MetaGMRF Type
MetaGMRF{M <: GMRFMetadata, T, P, G <: AbstractGMRF{T, P}} <: AbstractGMRF{T, P}
A wrapper that combines a core GMRF with metadata of type M. This allows for specialized behavior based on the metadata type while preserving the computational efficiency of the underlying GMRF.
Fields
gmrf::G
: The core computational GMRF (parametric type)metadata::M
: Domain-specific metadata
Usage
# Define metadata types
struct SpatialMetadata <: GMRFMetadata
coordinates::Matrix{Float64}
boundary_info::Vector{Int}
end
# Create wrapped GMRF
meta_gmrf = MetaGMRF(my_gmrf, SpatialMetadata(coords, boundary))
# Dispatch on metadata type for specialized behavior
function some_spatial_operation(mgmrf::MetaGMRF{SpatialMetadata})
# Access coordinates via mgmrf.metadata.coordinates
# Access GMRF via mgmrf.gmrf
end
GaussianMarkovRandomFields.GMRFMetadata Type
GMRFMetadata
Abstract base type for metadata that can be attached to GMRFs via MetaGMRF. Concrete subtypes should contain domain-specific information about the GMRF structure, coordinates, naming, etc.
sourceArithmetic
GaussianMarkovRandomFields.condition_on_observations Function
" condition_on_observations( x::GMRF, A::Union{AbstractMatrix,LinearMap}, Q_ϵ::Union{AbstractMatrix,LinearMap,Real}, y::AbstractVector=zeros(size(A)[1]), b::AbstractVector=zeros(size(A)[1]); # solver_blueprint parameter removed - no longer needed with LinearSolve )
Condition a GMRF x
on observations y = A * x + b + ϵ
where ϵ ~ N(0, Q_ϵ⁻¹)
.
Arguments
x::GMRF
: The GMRF to condition on.A::Union{AbstractMatrix,LinearMap}
: The matrixA
.Q_ϵ::Union{AbstractMatrix,LinearMap, Real}
: The precision matrix of the noise termϵ
. In case a real number is provided, it is interpreted as a scalar multiple of the identity matrix.y::AbstractVector=zeros(size(A)[1])
: The observationsy
; optional.b::AbstractVector=zeros(size(A)[1])
: Offset vectorb
; optional.
Keyword arguments
Note: solver_blueprint parameter removed - no longer needed with LinearSolve.jl
Returns
A GMRF
object representing the conditional GMRF x | (y = A * x + b + ϵ)
.
Notes
This function is deprecated. Use linear_condition
.
GaussianMarkovRandomFields.linear_condition Function
linear_condition(gmrf::GMRF; A, Q_ϵ, y, b=zeros(size(A, 1)), obs_precision_contrib=nothing)
Condition a GMRF on linear observations y = A * x + b + ϵ where ϵ ~ N(0, Q_ϵ^(-1)).
Arguments
gmrf::GMRF
: The prior GMRFA::Union{AbstractMatrix, LinearMap}
: Observation matrixQ_ϵ::Union{AbstractMatrix, LinearMap}
: Precision matrix of observation noisey::AbstractVector
: Observation valuesb::AbstractVector
: Offset vector (defaults to zeros)obs_precision_contrib
: Precomputed A' * Q_ϵ * A (optional optimization)
Returns
A new GMRF
representing the posterior distribution with updated mean and precision.
Performance Optimization
When A has special structure (e.g., index selection), users can precompute A' * Q_ϵ * A more efficiently and pass it to avoid redundant computation.
Notes
Uses information vector arithmetic for efficient conditioning without intermediate solves.
sourceGaussianMarkovRandomFields.joint_gmrf Function
" joint_gmrf( x1::AbstractGMRF, A::AbstractMatrix, Q_ϵ::AbstractMatrix, b::AbstractVector=spzeros(size(A)[1]) )
Return the joint GMRF of x1
and x2 = A * x1 + b + ϵ
where ϵ ~ N(0, Q_ϵ⁻¹)
.
Arguments
x1::AbstractGMRF
: The first GMRF.A::AbstractMatrix
: The matrixA
.Q_ϵ::AbstractMatrix
: The precision matrix of the noise termϵ
.b::AbstractVector=spzeros(size(A)[1])
: Offset vectorb
; optional.
Returns
A GMRF
object representing the joint GMRF of x1
and x2 = A * x1 + b + ϵ
.
Spatiotemporal setting
Types
GaussianMarkovRandomFields.AbstractSpatiotemporalGMRF Type
AbstractSpatiotemporalGMRF
A spatiotemporal GMRF is a GMRF that explicitly encodes the spatial and temporal structure of the underlying random field. All time points are modelled in one joint GMRF. It provides utilities to get statistics, draw samples and get the spatial discretization at a given time.
sourceGaussianMarkovRandomFields.ImplicitEulerConstantMeshSTGMRF Type
ImplicitEulerConstantMeshSTGMRF
A spatiotemporal GMRF with constant spatial discretization and an implicit Euler discretization of the temporal dynamics. Uses MetaGMRF for clean type structure.
sourceGaussianMarkovRandomFields.ConcreteConstantMeshSTGMRF Type
ConcreteConstantMeshSTGMRF
A concrete implementation of a spatiotemporal GMRF with constant spatial discretization. Uses MetaGMRF for clean type structure.
sourceQuantities
GaussianMarkovRandomFields.N_t Function
N_t(::AbstractSpatiotemporalGMRF)
Return the number of time points in the spatiotemporal GMRF.
sourceGaussianMarkovRandomFields.time_means Function
time_means(::AbstractSpatiotemporalGMRF)
Return the means of the spatiotemporal GMRF at each time point.
Returns
- A vector of means of length Nₜ, one for each time point.
GaussianMarkovRandomFields.time_vars Function
time_vars(::AbstractSpatiotemporalGMRF)
Return the marginal variances of the spatiotemporal GMRF at each time point.
Returns
- A vector of marginal variances of length Nₜ, one for each time point.
GaussianMarkovRandomFields.time_stds Function
time_stds(::AbstractSpatiotemporalGMRF)
Return the marginal standard deviations of the spatiotemporal GMRF at each time point.
Returns
- A vector of marginal standard deviations of length Nₜ, one for each time point.
GaussianMarkovRandomFields.time_rands Function
time_rands(::AbstractSpatiotemporalGMRF, rng::AbstractRNG)
Draw samples from the spatiotemporal GMRF at each time point.
Returns
- A vector of sample values of length Nₜ, one sample value vector for each time point.
GaussianMarkovRandomFields.discretization_at_time Function
discretization_at_time(::AbstractSpatiotemporalGMRF, t::Int)
Return the spatial discretization at time t
.
Utilities
GaussianMarkovRandomFields.spatial_to_spatiotemporal Function
spatial_to_spatiotemporal(
spatial_matrix::AbstractMatrix,
t_idx::Int,
N_t::Int,
)
Make a spatial matrix applicable to a spatiotemporal system at time index t_idx
. Results in a matrix that selects the spatial information exactly at time t_idx
.
Arguments
spatial_matrix::AbstractMatrix
: The spatial matrix.t_idx::Integer
: The time index.N_t::Integer
: The number of time points.
GaussianMarkovRandomFields.kronecker_product_spatiotemporal_model Function
kronecker_product_spatiotemporal_model(
Q_t::AbstractMatrix,
Q_s::AbstractMatrix,
spatial_disc::FEMDiscretization;
algorithm = nothing,
)
Create a spatiotemporal GMRF through a Kronecker product of the temporal and spatial precision matrices.
Arguments
Q_t::AbstractMatrix
: The temporal precision matrix.Q_s::AbstractMatrix
: The spatial precision matrix.spatial_disc::FEMDiscretization
: The spatial discretization.
Keyword arguments
algorithm
: The LinearSolve algorithm to use.