MaternKLExpansion#

class cuqipy_fenics.geometry.MaternKLExpansion(geometry, length_scale, num_terms, nu=0.5, boundary_conditions='Neumann', normalize=True)#

A geometry class that builds spectral representation of Matern covariance operator on the given input geometry. We create the representation using the stochastic partial differential operator, equation (15) in (Roininen, Huttunen and Lasanen, 2014). Zero Neumann boundary conditions are assumed for the stochastic partial differential equation (SPDE) and the default value of the smoothness parameter \(\nu\) is set to 0.5. To generate Matern field realizations, the method par2field() is used. The input p of this method need to be an n=dim i.i.d random variables that follow a normal distribution.

For more details about the formulation of the SPDE see: Roininen, L., Huttunen, J. M., & Lasanen, S. (2014). Whittle-Matérn priors for Bayesian statistical inversion with applications in electrical impedance tomography. Inverse Problems & Imaging, 8(2), 561.

Parameters:
  • geometry (cuqipy_fenics.geometry.Geometry) – An input geometry on which the Matern field representation is built (the geometry must have a mesh attribute)

  • length_scale (float) – Length scale parameter (controls correlation length)

  • num_terms (int) – Number of expansion terms to represent the Matern field realization

  • nu (float, default 0.5) – Smoothness parameter of the Matern field, must be greater then zero.

  • boundary_conditions (str) – Boundary conditions for the SPDE. Currently ‘Neumann’ for zero flux, and ‘zero’ for zero Dirichlet boundary conditions are supported.

  • normalize (bool, default True) – If True, the Matern field expansion modes are normalized to have a unit norm.

Example

import numpy as np
import matplotlib.pyplot as plt
from cuqipy_fenics.geometry import MaternKLExpansion, FEniCSContinuous
from cuqi.distribution import Gaussian
import dolfin as dl

mesh = dl.UnitSquareMesh(20,20)
V = dl.FunctionSpace(mesh, 'CG', 1)
geometry = FEniCSContinuous(V)
MaternGeometry = MaternKLExpansion(geometry,
                                length_scale = .2,
                                num_terms=128)

MaternField = Gaussian(mean=np.zeros(MaternGeometry.dim),
                cov=np.eye(MaternGeometry.dim),
                geometry= MaternGeometry)

samples = MaternField.sample()
samples.plot()
__init__(geometry, length_scale, num_terms, nu=0.5, boundary_conditions='Neumann', normalize=True)#

Methods

__init__(geometry, length_scale, num_terms)

fun2par(f)

The function to parameter map used to map function values back to parameters, if available.

fun2vec(fun)

Maps the function value (FEniCS object) to the corresponding vector representation of the function (ndarray of the function DOF values).

gradient(direction, wrt)

par2field(p)

Applies linear transformation of the parameters p to generate a realization of the Matern field (given that p is a sample of n=dim i.i.d random variables that follow a normal distribution)

par2fun(p)

The parameter to function map used to map parameters to function values in e.g. plotting.

plot(values[, is_par, plot_par])

Plots a function over the set defined by the geometry object.

plot_envelope(lo_values, hi_values[, ...])

Plots an envelope from lower and upper bounds over the set defined by the geometry object.

vec2fun(funvec)

Maps the vector representation of the function (ndarray of the function DOF values) to the function value (FEniCS object).

Attributes

axis_labels

boundary_conditions

eig_val

eig_vec

fun_dim

The dimension of the geometry (function space).

fun_is_array

Flag to indicate whether the function value is an array.

fun_shape

The shape of the geometry (function space).

function_space

funvec_dim

The dimension of the geometry (dimension of the vector representation of the function value).

funvec_shape

The shape of the geometry (shape of the vector representation of the function value).

grid

length_scale

mesh

normalize

nu

num_terms

par_dim

The dimension of the geometry (parameter space).

par_shape

The shape of the parameter space.

physical_dim

Returns the physical dimension of the geometry, e.g. 1, 2 or 3.

variables