Estimation Methods
These types and functions define the methods for covariance smoothing.
CovIterSolvers.AbstractEstimateMethod
— TypeAbstractEstimateMethod{T}
Supertype for all estimation methods in this package.
Subtypes
BSplineMethod{T}
: Methods based on B-splines.RBFKernelMethod{T}
: Methods based on Radial Basis Functions (RBF) kernels.
CovIterSolvers.BSplineMethod
— TypeBSplineMethod(order, knots)
An estimation method that uses a B-spline basis.
Arguments
order::Int
: The order of the B-spline (e.g., 4 for cubic B-splines).knots::AbstractVector
: The knot vector defining the B-spline basis.
Examples
julia> knots = 0.0:0.1:1.0;
julia> method = BSplineMethod(4, knots);
julia> method.order
4
CovIterSolvers.RBFKernelMethod
— TypeRBFKernelMethod{T} <: AbstractEstimateMethod{T}
Abstract supertype for methods based on a Radial Basis Function (RBF) kernel.
All RBF kernels define a function K(x, y)
that depends only on the distance between x
and y
, i.e., K(x, y) = f(||x - y||)
.
See also GaussianKernel
, LaplacianKernel
, MaternKernel
, CustomKernel
for concrete implementations.
CovIterSolvers.GaussianKernel
— TypeGaussianKernel(γ, [trunc_dist=nothing])
A Gaussian Radial Basis Function (RBF) kernel defined by K(x, y) = exp(-γ * ||x - y||^2)
.
Arguments
γ::Real
: The positive kernel parameter.trunc_dist::Union{Nothing, Real}
: An optional truncation distance. If provided, the kernel will return zero for||x - y|| > trunc_dist
.
Examples
julia> gk = GaussianKernel(0.5)
GaussianKernel{Float64}(0.5, nothing)
julia> gk_trunc = GaussianKernel(0.5, 10.0)
GaussianKernel{Float64}(0.5, 10.0)
CovIterSolvers.LaplacianKernel
— TypeLaplacianKernel(γ, [trunc_dist=nothing])
A Laplacian Radial Basis Function (RBF) kernel defined by K(x, y) = exp(-γ * ||x - y||)
.
Arguments
γ::Real
: The positive kernel parameter.trunc_dist::Union{Nothing, Real}
: An optional truncation distance. If provided, the kernel will return zero for||x - y|| > trunc_dist
.
Examples
julia> lk = LaplacianKernel(0.5)
LaplacianKernel{Float64}(0.5, nothing)
julia> lk_trunc = LaplacianKernel(0.5, 10.0)
LaplacianKernel{Float64}(0.5, 10.0)
CovIterSolvers.MaternKernel
— TypeMaternKernel(ν, γ, [trunc_dist=nothing])
A Matern Radial Basis Function (RBF) kernel defined by K(x, y) = (2^(1-ν) / Γ(ν)) * (γ * ||x - y||)^ν * K_{ν}(γ * ||x - y||)
, where K_{ν}
is the modified Bessel function of the second kind.
Arguments
ν::Real
: The smoothness parameter (ν = 1/2: Laplacian, ν → ∞: Gaussian).γ::Real
: The positive kernel parameter.trunc_dist::Union{Nothing, Real}
: An optional truncation distance. If provided, the kernel will return zero for||x - y|| > trunc_dist
.
Examples
julia> mk = MaternKernel(1.5, 0.5)
MaternKernel{Float64}(1.5, 0.5, nothing)
julia> mk_trunc = MaternKernel(1.5, 0.5, 10.0)
MaternKernel{Float64}(1.5, 0.5, 10.0)
CovIterSolvers.CustomKernel
— TypeCustomKernel(f, [trunc_dist=nothing])
A custom Radial Basis Function (RBF) kernel defined by a user-provided function f
.
The kernel is defined as K(x, y) = f(||x - y||)
.
Arguments
f::Function
: A function that takes a distance of typeT
and returns a real number.trunc_dist::Union{Nothing, T}
: An optional truncation distance. If provided, the kernel will return zero for||x - y|| > trunc_dist
.
Examples
julia> f(x) = exp(-x^2); # Example custom kernel function
julia> ck = CustomKernel(f)
CustomKernel{Float64,typeof(f)}(f, nothing)
julia> ck_trunc = CustomKernel(f, 10.0)
CustomKernel{Float64,typeof(f)}(f, 10.0)
CovIterSolvers.mean_fwd
— Functionmean_fwd(loc, myspline)
Construct a forward mapping matrix Φ
for a B-spline basis.
Each row of Φ
corresponds to a location in loc
, and each column corresponds to a B-spline basis function. The entry Φ[i, j]
is the value of the j-th basis function evaluated at the i-th location.
Arguments
loc::BlockVector{T}
: A block vector of locations.myspline::BSplineMethod{T}
: ABSplineMethod
object containing the B-spline order and knot vector.
Returns
BlockMatrix{T}
: A sparse block matrix representing the B-spline evaluation.
Examples
julia> knots = 0.0:0.2:1.0; # Simplified knots for a clear example
julia> myspline = BSplineMethod(4, knots);
julia> loc = loc_grid([2, 1]; seed=1);
julia> Φ = mean_fwd(loc, myspline);
julia> size(Φ)
(3, 8)
julia> Φ[1, 3] ≈ 0.1565989814693208
true
mean_fwd(loc, kernel)
Construct a Gram matrix K
from an RBF kernel and a set of locations.
The entry K[i, j]
is the value of the kernel evaluated at the distance between the i-th and j-th locations, i.e., K[i, j] = kernel(||loc[i] - loc[j]||)
.
Arguments
loc::BlockVector{T}
: A block vector of evaluation locations.kernel::RBFKernelMethod{T}
: An RBF kernel object. Ifkernel.trunc_dist
is set, the resulting matrix will be sparse.
Examples
julia> loc = loc_grid([2, 1]; seed=1);
julia> kernel = GaussianKernel(50.0); # A kernel with high decay
julia> K = mean_fwd(loc, kernel);
julia> size(K)
(3, 3)
julia> K[1, 2] ≈ 0.022251307501408822
true
CovIterSolvers.eval_fwd
— Functioneval_fwd(eval_points, myspline)
eval_fwd(g, myspline)
eval_fwd(eval_points, loc, kernel)
eval_fwd(g, loc, kernel)
Construct a forward mapping matrix E
by evaluating a basis or kernel at specific points.
This function has two main modes of operation:
- B-spline Basis Evaluation: When given a
BSplineMethod
, it computes the matrix of B-spline basis functions evaluated at each point ineval_points
. - RBF Kernel Evaluation: When given an
RBFKernelMethod
, it computes the Gram matrix of kernel evaluations between each point ineval_points
and each point inloc
.
Arguments
eval_points::AbstractVector
: A vector of points at which to evaluate the basis/kernel.g::Int
: An integer to create a regular grid ofg
evaluation points from 0 to 1.myspline::BSplineMethod
: A B-spline object containing the order and knots.loc::BlockVector
: A block vector of source locations for the RBF kernel.kernel::RBFKernelMethod
: An RBF kernel object.
Returns
BlockMatrix
: A sparse block matrix representing the evaluation of the basis or kernel.
Examples
julia> myeval = range(0, 1; length=10);
julia> knots = 0.0:0.2:1.0;
julia> myspline = BSplineMethod(4, knots);
julia> basis = BSplineBasis(BSplineOrder(order), knots);
julia> E1 = eval_fwd(myeval, myspline);
julia> E2 = eval_fwd(10, myspline);
julia> E1 == E2
true
julia> E1[1, 1] ≈ 1.0
true
CovIterSolvers.eval_covariance
— Functioneval_covariance(E, A)
Computes the covariance matrix Σ = E * A * E'
.
Arguments
E::BlockMatrix{T}
: The forward evaluation matrix.A::BlockDiagonal{T}
: A block-diagonal matrix, typically representing the prior covariance of the coefficients.
Returns
Matrix{T}
: The resulting dense covariance matrixΣ
.
Examples
julia> E = BlockMatrix(ones(3, 2), [3], [2]);
julia> A = BlockDiagonal([[1.0 0.5; 0.5 1.0]]);
julia> Σ = eval_covariance(E, A)
3×3 Matrix{Float64}:
3.0 3.0 3.0
3.0 3.0 3.0
3.0 3.0 3.0
CovIterSolvers.compute_kernel
— Functioncompute_kernel(dist, kernel)
Compute the value of an RBF kernel for a given distance.
This function uses multiple dispatch to select the correct kernel formula based on the type of the kernel
object.
Arguments
dist::Real
: The distance between two points.kernel::RBFKernelMethod
: The kernel object containing its parameters.
Returns
- The scalar value of the kernel evaluation.
Examples
julia> dist = 0.4;
julia> kernel = GaussianKernel(1.0);
julia> compute_kernel(dist, kernel) ≈ 0.8521437889662113
true