GMRFs
GaussianMarkovRandomFields.AbstractGMRF
— TypeAbstractGMRF
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.
GaussianMarkovRandomFields.GMRF
— TypeGMRF(mean, precision, solver_blueprint=DefaultSolverBlueprint())
A Gaussian Markov Random Field with mean mean
and precision matrix precision
.
Arguments
mean::AbstractVector
: The mean vector of the GMRF.precision::LinearMap
: The precision matrix (inverse covariance) of the GMRF.solver_blueprint::AbstractSolverBlueprint
: Blueprint specifying how to construct the solver.
Type Parameters
T<:Real
: The numeric type (e.g., Float64).PrecisionMap<:LinearMap{T}
: The type of the precision matrix LinearMap.Solver<:AbstractSolver
: The type of the solver instance.
Fields
mean::Vector{T}
: The mean vector.precision::PrecisionMap
: The precision matrix as a LinearMap.solver::Solver
: The solver instance for computing GMRF quantities.
Notes
The solver is constructed automatically using construct_solver(solver_blueprint, mean, precision)
and is used to compute means, variances, samples, and other GMRF quantities efficiently.
GaussianMarkovRandomFields.precision_map
— Functionprecision_map(::AbstractGMRF)
Return the precision (inverse covariance) map of the GMRF.
GaussianMarkovRandomFields.precision_matrix
— Functionprecision_matrix(::AbstractGMRF)
Return the precision (inverse covariance) matrix of the GMRF.
Arithmetic
GaussianMarkovRandomFields.condition_on_observations
— Function" conditiononobservations( x::AbstractGMRF, A::Union{AbstractMatrix,LinearMap}, Qϵ::Union{AbstractMatrix,LinearMap,Real}, y::AbstractVector=spzeros(size(A)[1]), b::AbstractVector=spzeros(size(A)[1]); solverblueprint::AbstractSolverBlueprint=CholeskySolverBlueprint() )
Condition a GMRF x
on observations y = A * x + b + ϵ
where ϵ ~ N(0, Q_ϵ⁻¹)
.
Arguments
x::AbstractGMRF
: 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=spzeros(size(A)[1])
: The observationsy
; optional.b::AbstractVector=spzeros(size(A)[1])
: Offset vectorb
; optional.
Keyword arguments
solver_blueprint::AbstractSolverBlueprint=CholeskySolverBlueprint()
: The solver blueprint; optional.
Returns
A LinearConditionalGMRF
object representing the conditional GMRF x | (y = A * x + b + ϵ)
.
GaussianMarkovRandomFields.LinearConditionalGMRF
— TypeLinearConditionalGMRF{G}(
prior::G,
A::Union{AbstractMatrix,LinearMap},
Q_ϵ::Union{AbstractMatrix,LinearMap},
y::AbstractVector,
b::AbstractVector=spzeros(size(A, 1)),
solver_blueprint::AbstractSolverBlueprint=DefaultSolverBlueprint(),
) where {G<:AbstractGMRF}
A GMRF conditioned on observations y = A * x + b + ϵ
where ϵ ~ N(0, Q_ϵ)
.
Arguments
prior::G
: The prior GMRF.A::Union{AbstractMatrix,LinearMap}
: The matrixA
.Q_ϵ::Union{AbstractMatrix,LinearMap}
: The precision matrix of the noise termϵ
.y::AbstractVector=spzeros(size(A, 1))
: The observationsy
.b::AbstractVector=spzeros(size(A, 1))
: The offset vectorb
.solver_blueprint::AbstractSolverBlueprint=DefaultSolverBlueprint()
: The solver blueprint.
GaussianMarkovRandomFields.joint_gmrf
— Function" jointgmrf( 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
— TypeAbstractSpatiotemporalGMRF
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.
GaussianMarkovRandomFields.ConstantMeshSTGMRF
— TypeConstantMeshSTGMRF
A spatiotemporal GMRF with constant spatial discretization.
GaussianMarkovRandomFields.ImplicitEulerConstantMeshSTGMRF
— TypeImplicitEulerConstantMeshSTGMRF
A spatiotemporal GMRF with constant spatial discretization and an implicit Euler discretization of the temporal dynamics.
GaussianMarkovRandomFields.ConcreteConstantMeshSTGMRF
— TypeConcreteConstantMeshSTGMRF
A concrete implementation of a spatiotemporal GMRF with constant spatial discretization.
Quantities
GaussianMarkovRandomFields.N_t
— FunctionN_t(::AbstractSpatiotemporalGMRF)
Return the number of time points in the spatiotemporal GMRF.
GaussianMarkovRandomFields.time_means
— Functiontime_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
— Functiontime_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
— Functiontime_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
— Functiontime_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
— Functiondiscretization_at_time(::AbstractSpatiotemporalGMRF, t::Int)
Return the spatial discretization at time t
.
Utilities
GaussianMarkovRandomFields.spatial_to_spatiotemporal
— Functionspatial_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
— Functionkronecker_product_spatiotemporal_model(
Q_t::AbstractMatrix,
Q_s::AbstractMatrix,
spatial_disc::FEMDiscretization;
solver_blueprint = DefaultSolverBlueprint(),
)
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
solver_blueprint::AbstractSolverBlueprint=DefaultSolverBlueprint()
: The solver blueprint.