Linear Maps

The construction of GMRFs involves various kinds of structured matrices. These structures may be leveraged in downstream computations to save compute and memory. But to make this possible, we need to actually keep track of these structures - which we achieve through diverse subtypes of LinearMap.

GaussianMarkovRandomFields.SymmetricBlockTridiagonalMapType
SymmetricBlockTridiagonalMap(
    diagonal_blocks::Tuple{LinearMap{T},Vararg{LinearMap{T},ND}},
    off_diagonal_blocks::Tuple{LinearMap{T},Vararg{LinearMap{T},NOD}},
)

A linear map representing a symmetric block tridiagonal matrix with diagonal blocks diagonal_blocks and lower off-diagonal blocks off_diagonal_blocks.

source
GaussianMarkovRandomFields.SSMBidiagonalMapType
SSMBidiagonalMap{T}(
    A::LinearMap{T},
    B::LinearMap{T},
    C::LinearMap{T},
    N_t::Int,
)

Represents the block-bidiagonal map given by the (Nt) x (Nt - 1) sized block structure:

\[\begin{bmatrix} A & 0 & \cdots & 0 \\ B & C & \cdots & 0 \\ 0 & B & C & 0 \\ \vdots & \vdots & \ddots & \vdots \\ 0 & 0 & \cdots & B \end{bmatrix}\]

which occurs as a square root in the discretization of GMRF-based state-space models. N_t is the total number of blocks along the rows.

source
GaussianMarkovRandomFields.LinearMapWithSqrtType
LinearMapWithSqrt{T}(
    A::LinearMap{T},
    A_sqrt::LinearMap{T},
)

A symmetric positive definite linear map A with known square root A_sqrt, i.e. A = A_sqrt * A_sqrt'. Behaves just like A, but taking the square root directly returns A_sqrt.

Arguments

  • A::LinearMap{T}: The linear map A.
  • A_sqrt::LinearMap{T}: The square root of A.
source
GaussianMarkovRandomFields.CholeskySqrtType
CholeskySqrt{T}(cho::Union{Cholesky{T},SparseArrays.CHOLMOD.Factor{T}})

A linear map representing the square root obtained from a Cholesky factorization, i.e. for A = L * L', this map represents L.

Arguments

  • cho::Union{Cholesky{T},SparseArrays.CHOLMOD.Factor{T}}: The Cholesky factorization of a symmetric positive definite matrix.
source
GaussianMarkovRandomFields.ADJacobianMapType
ADJacobianMap(f::Function, x₀::AbstractVector{T}, N_outputs::Int)

A linear map representing the Jacobian of f at x₀. Uses forward-mode AD in a matrix-free way, i.e. we do not actually store the Jacobian in memory and only compute JVPs.

Requires ForwardDiff.jl!

Arguments

  • f::Function: Function to differentiate.
  • x₀::AbstractVector{T}: Input vector at which to evaluate the Jacobian.
  • N_outputs::Int: Output dimension of f.
source
GaussianMarkovRandomFields.ADJacobianAdjointMapType
ADJacobianAdjointMap{T}(f::Function, x₀::AbstractVector{T}, N_outputs::Int)

A linear map representing the adjoint of the Jacobian of f at x₀. Uses reverse-mode AD in a matrix-free way, i.e. we do not actually store the Jacobian in memory and only compute VJPs.

Requires Zygote.jl!

Arguments

  • f::Function: Function to differentiate.
  • x₀::AbstractVector{T}: Input vector at which to evaluate the Jacobian.
  • N_outputs::Int: Output dimension of f.
source