MFEM v4.7.0
Finite element discretization library
|
#include <spde_solver.hpp>
Public Member Functions | |
SPDESolver (real_t nu, const Boundary &bc, ParFiniteElementSpace *fespace, real_t l1=0.1, real_t l2=0.1, real_t l3=0.1, real_t e1=0.0, real_t e2=0.0, real_t e3=0.0) | |
~SPDESolver () | |
Destructor. | |
void | Solve (ParLinearForm &b, ParGridFunction &x) |
void | SetupRandomFieldGenerator (int seed=0) |
void | GenerateRandomField (ParGridFunction &x) |
void | SetPrintLevel (int print_level) |
Set the print level. | |
Static Public Member Functions | |
static real_t | ConstructNormalizationCoefficient (real_t nu, real_t l1, real_t l2, real_t l3, int dim) |
static DenseMatrix | ConstructMatrixCoefficient (real_t l1, real_t l2, real_t l3, real_t e1, real_t e2, real_t e3, real_t nu, int dim) |
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.
mfem::spde::SPDESolver::SPDESolver | ( | real_t | nu, |
const Boundary & | bc, | ||
ParFiniteElementSpace * | fespace, | ||
real_t | l1 = 0.1, | ||
real_t | l2 = 0.1, | ||
real_t | l3 = 0.1, | ||
real_t | e1 = 0.0, | ||
real_t | e2 = 0.0, | ||
real_t | e3 = 0.0 ) |
Constructor.
nu | The coefficient nu, smoothness of the solution. |
bc | Boundary conditions. |
fespace | Finite element space. |
l1 | Correlation length in x |
l2 | Correlation length in y |
l3 | Correlation length in z |
e1 | Rotation angle in x |
e2 | Rotation angle in y |
e3 | Rotation angle in z |
Definition at line 351 of file spde_solver.cpp.
mfem::spde::SPDESolver::~SPDESolver | ( | ) |
Destructor.
Definition at line 814 of file spde_solver.cpp.
|
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 596 of file spde_solver.cpp.
|
static |
Construct the normalization coefficient eta of the white noise right hands side.
Definition at line 571 of file spde_solver.cpp.
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 555 of file spde_solver.cpp.
|
inline |
Set the print level.
Definition at line 137 of file spde_solver.hpp.
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 546 of file spde_solver.cpp.
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 472 of file spde_solver.cpp.