MFEM  v4.6.0
Finite element discretization library
Public Member Functions | Static Public Member Functions | List of all members
mfem::spde::SPDESolver Class Reference

#include <spde_solver.hpp>

Public Member Functions

 SPDESolver (double nu, const Boundary &bc, ParFiniteElementSpace *fespace, double l1=0.1, double l2=0.1, double l3=0.1, double e1=0.0, double e2=0.0, double e3=0.0)
 
 ~SPDESolver ()
 Destructor. More...
 
void Solve (ParLinearForm &b, ParGridFunction &x)
 
void SetupRandomFieldGenerator (int seed=0)
 
void GenerateRandomField (ParGridFunction &x)
 
void SetPrintLevel (int print_level)
 Set the print level. More...
 

Static Public Member Functions

static double ConstructNormalizationCoefficient (double nu, double l1, double l2, double l3, int dim)
 
static DenseMatrix ConstructMatrixCoefficient (double l1, double l2, double l3, double e1, double e2, double e3, double nu, int dim)
 

Detailed Description

Solver for the SPDE method based on a rational approximation with the AAA algorithm. The SPDE method is described in the paper Lindgren, F., Rue, H., Lindström, J. (2011). An explicit link between Gaussian fields and Gaussian Markov random fields: the stochastic partial differential equation approach. Journal of the Royal Statistical Society: Series B (Statistical Methodology), 73(4), 423–498. https://doi.org/10.1111/j.1467-9868.2011.00777.x

The solver solves the SPDE problem defined as (A)^-alpha u = b where A is A = div ( Theta(x) grad + Id ) u(x) and alpha is given as alpha = (2 nu + dim) / 2. Theta (anisotropy tensor) and nu (smoothness) can be specified in the constructor. Traditionally, the SPDE method requires the specification of a white noise right hands side. SPDESolver accepts arbitrary right hand sides but the solver has only been tested with white noise.

Definition at line 91 of file spde_solver.hpp.

Constructor & Destructor Documentation

◆ SPDESolver()

mfem::spde::SPDESolver::SPDESolver ( double  nu,
const Boundary bc,
ParFiniteElementSpace fespace,
double  l1 = 0.1,
double  l2 = 0.1,
double  l3 = 0.1,
double  e1 = 0.0,
double  e2 = 0.0,
double  e3 = 0.0 
)

Constructor.

Parameters
nuThe coefficient nu, smoothness of the solution.
bcBoundary conditions.
fespaceFinite element space.
l1Correlation length in x
l2Correlation length in y
l3Correlation length in z
e1Rotation angle in x
e2Rotation angle in y
e3Rotation angle in z

Definition at line 350 of file spde_solver.cpp.

◆ ~SPDESolver()

mfem::spde::SPDESolver::~SPDESolver ( )

Destructor.

Definition at line 813 of file spde_solver.cpp.

Member Function Documentation

◆ ConstructMatrixCoefficient()

DenseMatrix mfem::spde::SPDESolver::ConstructMatrixCoefficient ( double  l1,
double  l2,
double  l3,
double  e1,
double  e2,
double  e3,
double  nu,
int  dim 
)
static

Construct the second order tensor (matrix coefficient) Theta from the equation R^T(e1,e2,e3) diag(l1,l2,l3) R (e1,e2,e3).

Definition at line 595 of file spde_solver.cpp.

◆ ConstructNormalizationCoefficient()

double mfem::spde::SPDESolver::ConstructNormalizationCoefficient ( double  nu,
double  l1,
double  l2,
double  l3,
int  dim 
)
static

Construct the normalization coefficient eta of the white noise right hands side.

Definition at line 570 of file spde_solver.cpp.

◆ GenerateRandomField()

void mfem::spde::SPDESolver::GenerateRandomField ( ParGridFunction x)

Generate a random field. Calls back to solve but generates the stochastic load b internally.

Definition at line 554 of file spde_solver.cpp.

◆ SetPrintLevel()

void mfem::spde::SPDESolver::SetPrintLevel ( int  print_level)
inline

Set the print level.

Definition at line 137 of file spde_solver.hpp.

◆ SetupRandomFieldGenerator()

void mfem::spde::SPDESolver::SetupRandomFieldGenerator ( int  seed = 0)

Set up the random field generator. If called more than once it resets the generator.

Definition at line 545 of file spde_solver.cpp.

◆ Solve()

void mfem::spde::SPDESolver::Solve ( ParLinearForm b,
ParGridFunction x 
)

Solve the SPDE for a given right hand side b. May alter b if the exponent (alpha) is larger than 1. We avoid copying by default. If you need b later on, make a copy of it before calling this function.

Definition at line 471 of file spde_solver.cpp.


The documentation for this class was generated from the following files: