statFEM_analysis package
statFEM_analysis.oneDim module
- statFEM_analysis.oneDim.BigPhiMat(J, grid)[source]
Function to compute the \(\Phi\) matrix.
- Parameters
J (int) – controls the FE mesh size (\(h=1/J\))
grid (array) – grid in the definition of \(\Phi\)
- Returns
\(\Phi\) matrix returned as a sparse csr_matrix.
- Return type
csr_matrix
- class statFEM_analysis.oneDim.MyExpression(*args, **kwargs)[source]
Bases:
dolfin.function.expression.UserExpression
Class to allow users to user their own functions to create a FEniCS UserExpression.
- statFEM_analysis.oneDim.c_post(x, y, c, Y, B)[source]
This function evaluates the posterior covariance at \((x,y)\)
- Parameters
x (float) – one of two points to evaluate the covariance at
y (float) – one of two points to evaluate the covariance at
c (function) – function which returns the prior covariance at any given pair \((x,y)\)
Y (array) – vector of sensor locations
B (array) – the matrix \(\epsilon^{2}I+C_{Y}\) to be inverted in order to obtain the posterior
- Returns
posterior covariance evaluated at \((x,y)\)
- Return type
float
- statFEM_analysis.oneDim.c_post_fem_assembler(J, k_f, grid, Y, ε, par, trans)[source]
Function to assemble the statFEM posterior cov mat on a reference grid specified by grid.
- Parameters
J (int) – controls the FE mesh size (\(h=1/J\))
k_f (function) – the covariance function for the forcing \(f\)
grid (array) – the fixed reference grid \(\{x_{i}\}_{i=1}^{N}\) on which to assemble the posterior cov mat
Y (array) – vector of sensor locations
par (bool) – argument indicating whether the computation of the forcing cov mat should be done in parallel
trans (bool) – argument indicating whether the computation of the forcing cov mat should be computed assuming k_f is translation invariant or not
- Returns
statFEM posterior cov matrix evaluated on grid
- Return type
array
- statFEM_analysis.oneDim.cov_assembler(J, k_f, grid, parallel, translation_inv)[source]
Function to assemble the approximate FEM covariance matrix on the reference grid.
- Parameters
J (int) – controls the FE mesh size (\(h=1/J\))
k_f (function) – the covariance function for the forcing \(f\)
grid (array) – the reference grid where the FEM cov matrix should be computed on
parallel (bool) – argument indicating whether the intermediate computation of \(C_f\) should be done in parallel
translation_inv (bool) – argument indicating whether the intermediate computation of \(C_f\) should be computed assuming k_f is translation invariant or not
- Returns
approximate FEM covariance matrix on the reference grid
- Return type
array
- statFEM_analysis.oneDim.fem_cov_assembler_post(J, k_f, Y, parallel, translation_inv)[source]
Function to create the matrix \(C_{Y,h}\) and the vector function \(c^{(h)}\) required for the statFEM posterior mean.
- Parameters
J (int) – controls the FE mesh size (\(h=1/J\))
k_f (function) – the covariance function for the forcing \(f\)
Y (array) – vector of sensor locations
parallel (bool) – argument indicating whether the computation of the forcing cov mat should be done in parallel
translation_inv (bool) – argument indicating whether the computation of the forcing cov mat should be computed assuming k_f is translation invariant or not
- Returns
returns tuple of matrix \(C_{Y,h}\) and the vector function \(c^{(h)}\)
- Return type
tuple
- statFEM_analysis.oneDim.gen_sensor(ε, m, k, Y, u_quad, grid, par=False, trans=True, tol=1e-09, maxiter=50, require_f=False)[source]
Function to generate noisy sensor observations of the solution u on a sensor grid Y.
- Parameters
ϵ (float) – controls the amount of sensor noise
m (function) – mean function for the forcing \(f\)
k (function) – cov function for the forcing \(f\)
Y (array) – vector of sensor locations
u_quad (function) – function to accurately compute the solution \(u\) given a realisation of the forcing \(f\)
grid (array) – grid where forcing \(f\) is sampled on
par (bool, optional) – argument indicating whether the computation of the forcing cov matrix should be done in parallel, by default False
trans (bool, optional) – argument indicating whether the computation of the forcing cov matrix should be computed assuming k is translation invariant or not, by default True
tol (_type_, optional) – controls the size of the tiny diagonal perturbation added to forcing cov matrix to ensure it is strictly positive definite, by default 1e-9
maxiter (int, optional) – parameter which controls the accuracy of the quadrature used in u_quad, by default 50
require_f (bool, optional) – argument indicating whether or not to also return the realisation of the forcing \(f\), by default False
- Returns
If require_f is set to True returns a tuple of two arrays containing the generated noisy sensor observations and the realisation of the forcing. If require_f is set to False returns just the generated noisy sensor observations.
- Return type
tuple or array
- statFEM_analysis.oneDim.kernMat(k, grid, parallel=True, translation_inv=False)[source]
Function to compute the covariance matrix \(K\) corresponding to the covariance kernel \(k\) on a grid.
- Parameters
k (function) – forcing covariance function
grid (array) – grid on which to evaluate the covariance function
parallel (bool, optional) – specifies whether or not the cov kernel is translation invariant, by default True
translation_inv (bool, optional) – specifies whether or not the cov kernel is translation invariant, by default False
- Returns
Covariance matrix K with \(ij-\) where \(x_i\) is the \(i\)-th point of the grid.
- Return type
array
- statFEM_analysis.oneDim.m_post(x, m, c, v, Y, B)[source]
This function evalutes the posterior mean at the point \(x\).
- Parameters
x (float) – point where the posterior mean will be evaluated
m (function) – function which computes the prior mean at a given point y
c (function) – function which returns the vector (c(x,y)) for y in Y (note: c is the prior covariance function)
v (array) – vector of noisy sensor readings
Y (array) – vector of sensor locations
B (array) – the matrix \(\epsilon^{2}I+C_Y\) to be inverted in order to obtain the posterior
- Returns
posterior mean at the point \(x\)
- Return type
float
- statFEM_analysis.oneDim.m_post_fem_assembler(J, f_bar, k_f, ε, Y, v_dat, par=False, trans=True)[source]
Function to assemble the statFEM posterior mean function.
- Parameters
J (int) – controls the FE mesh size (\(h=1/J\))
f_bar (function) – the mean function for the forcing \(f\)
k_f (function) – the covariance function for the forcing \(f\)
ϵ (float) – controls the amount of sensor noise
Y (array) – vector of sensor locations
v_dat (array) – vector of noisy sensor observations
par (bool, optional) – argument passed to fem_cov_assembler_post’s argument parallel, by default False
trans (bool, optional) – argument passed to fem_cov_assembler_post’s argument translation_inv, by default True
- Returns
statFEM posterior mean function
- Return type
function
- statFEM_analysis.oneDim.mean_assembler(h, f_bar)[source]
This function assembles the mean for the statFEM prior for our 1-D problem.
- Parameters
h (float) – mesh size h
f_bar (FEniCS function/expression/constant) – mean function for the forcing
- Returns
mean of the approximate statFEM prior
- Return type
FEniCS function
- statFEM_analysis.oneDim.post_fem_cov_assembler(J, k_f, grid, Y, parallel, translation_inv)[source]
Function which assembles the matrices \(\Sigma_X,\Sigma_{XY}\), and \(\Sigma_Y\) required for the statFEM posterior covariance.
- Parameters
J (int) – controls the FE mesh size (\(h=1/J\))
k_f (function) – the covariance function for the forcing \(f\)
grid (array) – the fixed reference grid \(\{x_{i}\}_{i=1}^{N}\) on which to assemble the posterior cov mat
Y (array) – vector of sensor locations
parallel (bool) – argument indicating whether the computation of the forcing cov mat should be done in parallel
translation_inv (bool) – argument indicating whether the computation of the forcing cov mat should be computed assuming k_f is translation invariant or not
- Returns
tuple of the matrices \(\Sigma_X,\Sigma_{XY}\), and \(\Sigma_Y\)
- Return type
tuple
- statFEM_analysis.oneDim.sample_gp(n_sim, m, k, grid, par=False, trans=True, tol=1e-09)[source]
Function to sample a GP with mean \(m\) and cov \(k\) on a grid.
- Parameters
n_sim (int) – number of trajectories to be sampled
m (function) – mean function for the GP
k (function) – cov function for the GP
grid (array) – grid of points on which to sample the GP
par (bool, optional) – argument indicating whether the computation of the cov matrix should be done in parallel, by default False
trans (bool, optional) – argument indicating whether the computation of the cov matrix should be computed assuming k is translation invariant or not, by default True
tol (_type_, optional) – controls the size of the tiny diagonal perturbation added to cov matrix to ensure it is strictly positive definite, by default 1e-9
- Returns
array containing n_sim samples from a GP with mean \(m\) and cov \(k\) evaluated on the grid grid
- Return type
array
statFEM_analysis.twoDim module
- statFEM_analysis.twoDim.BigPhiMat(J, grid)[source]
Function to compute the \(\Phi\) matrix.
- Parameters
J (int) – controls the FE mesh size (\(h=1/J^{2}\))
grid (array) – grid in the definition of \(\Phi\)
- Returns
\(\Phi\) matrix returned as a sparse csr_matrix.
- Return type
csr_matrix
- class statFEM_analysis.twoDim.MyExpression(*args, **kwargs)[source]
Bases:
dolfin.function.expression.UserExpression
Class to allow users to user their own functions to create a FEniCS UserExpression.
- statFEM_analysis.twoDim.c_post(x, y, c, Y, B)[source]
This function evaluates the posterior covariance at \((x,y)\)
- Parameters
x (float) – one of two points to evaluate the covariance at
y (float) – one of two points to evaluate the covariance at
c (function) – function which returns the prior covariance at any given pair \((x,y)\)
Y (array) – vector of sensor locations
B (array) – the matrix \(\epsilon^{2}I+C_{Y}\) to be inverted in order to obtain the posterior
- Returns
posterior covariance evaluated at \((x,y)\)
- Return type
float
- statFEM_analysis.twoDim.c_post_fem_assembler(J, k_f, grid, Y, ε, par, trans)[source]
Function to assemble the statFEM posterior cov mat on a reference grid specified by grid.
- Parameters
J (int) – controls the FE mesh size (\(h=1/J^{2}\))
k_f (function) – the covariance function for the forcing \(f\)
grid (array) – the fixed reference grid \(\{x_{i}\}_{i=1}^{N}\) on which to assemble the posterior cov mat
Y (array) – vector of sensor locations
par (bool) – argument indicating whether the computation of the forcing cov mat should be done in parallel
trans (bool) – argument indicating whether the computation of the forcing cov mat should be computed assuming k_f is translation invariant or not
- Returns
statFEM posterior cov matrix evaluated on grid
- Return type
array
- statFEM_analysis.twoDim.cov_assembler(J, k_f, grid, parallel, translation_inv)[source]
Function to assemble the approximate FEM covariance matrix on the reference grid.
- Parameters
J (int) – controls the FE mesh size (\(h=1/J^{2}\))
k_f (function) – the covariance function for the forcing \(f\)
grid (array) – the reference grid where the FEM cov matrix should be computed on
parallel (bool) – argument indicating whether the intermediate computation of \(C_f\) should be done in parallel
translation_inv (bool) – argument indicating whether the intermediate computation of \(C_f\) should be computed assuming k_f is translation invariant or not
- Returns
approximate FEM covariance matrix on the reference grid
- Return type
array
- statFEM_analysis.twoDim.fem_cov_assembler_post(J, k_f, Y, parallel, translation_inv)[source]
Function to create the matrix \(C_{Y,h}\) and the vector function \(c^{(h)}\) required for the statFEM posterior mean.
- Parameters
J (int) – controls the FE mesh size (\(h=1/J^{2}\))
k_f (function) – the covariance function for the forcing \(f\)
Y (array) – vector of sensor locations
parallel (bool) – argument indicating whether the computation of the forcing cov mat should be done in parallel
translation_inv (bool) – argument indicating whether the computation of the forcing cov mat should be computed assuming k_f is translation invariant or not
- Returns
returns tuple of matrix \(C_{Y,h}\) and the vector function \(c^{(h)}\)
- Return type
tuple
- statFEM_analysis.twoDim.gen_sensor(ε, m, k, Y, J, par, trans, tol=1e-09, require=False)[source]
Function to generate noisy sensor observations of the solution u on a sensor grid Y.
- Parameters
ϵ (float) – controls the amount of sensor noise
m (function) – mean function for the forcing \(f\)
k (function) – cov function for the forcing \(f\)
Y (array) – vector of sensor locations
J (int) – controls the FE mesh size (\(h=1/J^{2}\))
par (bool, optional) – argument indicating whether the computation of the forcing cov matrix should be done in parallel, by default False
trans (bool, optional) – argument indicating whether the computation of the forcing cov matrix should be computed assuming k is translation invariant or not, by default True
tol (_type_, optional) – controls the size of the tiny diagonal perturbation added to forcing cov matrix to ensure it is strictly positive definite, by default 1e-9
require_f (bool, optional) – argument indicating whether or not to also return the realisation of the forcing \(f\), by default False
- Returns
If require_f is set to True returns a tuple of two arrays containing the generated noisy sensor observations and the realisation of the forcing. If require_f is set to False returns just the generated noisy sensor observations.
- Return type
tuple or array
- statFEM_analysis.twoDim.kernMat(k, grid, parallel=True, translation_inv=False)[source]
Function to compute the covariance matrix \(K\) corresponding to the covariance kernel \(k\) on a grid.
- Parameters
k (function) – forcing covariance function
grid (array) – grid on which to evaluate the covariance function
parallel (bool, optional) – specifies whether or not the cov kernel is translation invariant, by default True
translation_inv (bool, optional) – specifies whether or not the cov kernel is translation invariant, by default False
- Returns
Covariance matrix K with \(ij-\) where \(x_i\) is the \(i\)-th point of the grid.
- Return type
array
- statFEM_analysis.twoDim.m_post(x, m, c, v, Y, B)[source]
This function evalutes the posterior mean at the point \(x\).
- Parameters
x (float) – point where the posterior mean will be evaluated
m (function) – function which computes the prior mean at a given point y
c (function) – function which returns the vector (c(x,y)) for y in Y (note: c is the prior covariance function)
v (array) – vector of noisy sensor readings
Y (array) – vector of sensor locations
B (array) – the matrix \(\epsilon^{2}I+C_Y\) to be inverted in order to obtain the posterior
- Returns
posterior mean at the point \(x\)
- Return type
float
- statFEM_analysis.twoDim.m_post_fem_assembler(J, f_bar, k_f, ε, Y, v_dat, par=False, trans=True)[source]
Function to assemble the statFEM posterior mean function.
- Parameters
J (int) – controls the FE mesh size (\(h=1/J^{2}\))
f_bar (function) – the mean function for the forcing \(f\)
k_f (function) – the covariance function for the forcing \(f\)
ϵ (float) – controls the amount of sensor noise
Y (array) – vector of sensor locations
v_dat (array) – vector of noisy sensor observations
par (bool, optional) – argument passed to fem_cov_assembler_post’s argument parallel, by default False
trans (bool, optional) – argument passed to fem_cov_assembler_post’s argument translation_inv, by default True
- Returns
statFEM posterior mean function
- Return type
function
- statFEM_analysis.twoDim.mean_assembler(h, f_bar)[source]
This function assembles the mean for the statFEM prior for our 1-D problem.
- Parameters
h (float) – mesh size h
f_bar (FEniCS function/expression/constant) – mean function for the forcing
- Returns
mean of the approximate statFEM prior
- Return type
FEniCS function
- statFEM_analysis.twoDim.post_fem_cov_assembler(J, k_f, grid, Y, parallel, translation_inv)[source]
Function which assembles the matrices \(\Sigma_X,\Sigma_{XY}\), and \(\Sigma_Y\) required for the statFEM posterior covariance.
- Parameters
J (int) – controls the FE mesh size (\(h=1/J^{2}\))
k_f (function) – the covariance function for the forcing \(f\)
grid (array) – the fixed reference grid \(\{x_{i}\}_{i=1}^{N}\) on which to assemble the posterior cov mat
Y (array) – vector of sensor locations
parallel (bool) – argument indicating whether the computation of the forcing cov mat should be done in parallel
translation_inv (bool) – argument indicating whether the computation of the forcing cov mat should be computed assuming k_f is translation invariant or not
- Returns
tuple of the matrices \(\Sigma_X,\Sigma_{XY}\), and \(\Sigma_Y\)
- Return type
tuple
- statFEM_analysis.twoDim.sample_gp(n_sim, m, k, grid, par, trans, tol=1e-09)[source]
Function to sample a GP with mean \(m\) and cov \(k\) on a grid.
- Parameters
n_sim (int) – number of trajectories to be sampled
m (function) – mean function for the GP
k (function) – cov function for the GP
grid (array) – grid of points on which to sample the GP
par (bool, optional) – argument indicating whether the computation of the cov matrix should be done in parallel, by default False
trans (bool, optional) – argument indicating whether the computation of the cov matrix should be computed assuming k is translation invariant or not, by default True
tol (_type_, optional) – controls the size of the tiny diagonal perturbation added to cov matrix to ensure it is strictly positive definite, by default 1e-9
- Returns
array containing n_sim samples from a GP with mean \(m\) and cov \(k\) evaluated on the grid grid
- Return type
array
statFEM_analysis.maxDist module
- statFEM_analysis.maxDist.wass(a, b, n_bins)[source]
This function computes an approximation of the 2-Wasserstein distance between two datasets.
- Parameters
a (array) – dataset
b (array) – dataset
n_bins (int) – controls number of bins used to create the histograms for the dataset.
- Returns
estimate of the 2-Wasserstein distance between the two data-sets a and b
- Return type
float