![]() |
MFEM v4.8.0
Finite element discretization library
|
#include <lininteg.hpp>
Public Member Functions | |
| WhiteGaussianNoiseDomainLFIntegrator (int seed_=0) | |
| Sets the seed_ of the random number generator. A fixed seed allows for a reproducible sequence of white noise vectors. | |
| WhiteGaussianNoiseDomainLFIntegrator (MPI_Comm comm_, int seed_) | |
| Sets the MPI communicator comm_ and the seed_ of the random number generator. A fixed seed allows for a reproducible sequence of white noise vectors. | |
| WhiteGaussianNoiseDomainLFIntegrator (int seed_=0) | |
| Sets the seed_ of the random number generator. A fixed seed allows for a reproducible sequence of white noise vectors. | |
| void | SetSeed (int seed) |
| Sets/resets the seed of the random number generator. | |
| void | AssembleRHSElementVect (const FiniteElement &el, ElementTransformation &Tr, Vector &elvect) override |
| void | SaveFactors (int NE) |
| Saves the lower triangular matrices in the element-wise Cholesky decomposition. The parameter NE should be the number of elements in the mesh. | |
| void | ResetFactors (int NE) |
| Resets the array of saved lower triangular Cholesky decomposition matrices. The parameter NE should be the number of elements in the mesh. | |
| ~WhiteGaussianNoiseDomainLFIntegrator () | |
| virtual void | AssembleRHSElementVect (const FiniteElement &el, FaceElementTransformations &Tr, Vector &elvect) |
| virtual void | AssembleRHSElementVect (const FiniteElement &el1, const FiniteElement &el2, FaceElementTransformations &Tr, Vector &elvect) |
Public Member Functions inherited from mfem::LinearFormIntegrator | |
| virtual bool | SupportsDevice () const |
| Method probing for assembly on device. | |
| virtual void | AssembleDevice (const FiniteElementSpace &fes, const Array< int > &markers, Vector &b) |
| Method defining assembly on device. | |
| virtual | ~LinearFormIntegrator () |
Public Member Functions inherited from mfem::Integrator | |
| Integrator (const IntegrationRule *ir=NULL) | |
| Create a new Integrator, optionally providing a prescribed quadrature rule to use in assembly. | |
| virtual void | SetIntRule (const IntegrationRule *ir) |
| Prescribe a fixed IntegrationRule to use, or set to null to let the integrator choose an appropriate rule. | |
| void | SetIntegrationRule (const IntegrationRule &ir) |
| Prescribe a fixed IntegrationRule to use. Sets the NURBS patch integration rule to null. | |
| void | SetNURBSPatchIntRule (NURBSMeshRules *pr) |
| Sets an integration rule for use on NURBS patches. | |
| bool | HasNURBSPatchIntRule () const |
| Check if a NURBS patch integration rule has been set. | |
| const IntegrationRule * | GetIntRule () const |
| Directly return the IntRule pointer (possibly null) without checking for NURBS patch rules or falling back on a default. | |
| const IntegrationRule * | GetIntegrationRule () const |
| Equivalent to GetIntRule, but retained for backward compatibility with applications. | |
Additional Inherited Members | |
Protected Member Functions inherited from mfem::LinearFormIntegrator | |
| LinearFormIntegrator (const IntegrationRule *ir=NULL) | |
Protected Member Functions inherited from mfem::Integrator | |
| const IntegrationRule * | GetIntegrationRule (const FiniteElement &trial_fe, const FiniteElement &test_fe, const ElementTransformation &trans) const |
| Returns an integration rule based on the arguments and internal state of the Integrator object. | |
| const IntegrationRule * | GetIntegrationRule (const FiniteElement &el, const ElementTransformation &trans) const |
| Returns an integration rule based on the arguments and internal state. (Version for identical trial_fe and test_fe) | |
| virtual const IntegrationRule * | GetDefaultIntegrationRule (const FiniteElement &trial_fe, const FiniteElement &test_fe, const ElementTransformation &trans) const |
| Subclasses should override to choose a default integration rule. | |
Protected Attributes inherited from mfem::Integrator | |
| const IntegrationRule * | IntRule |
| NURBSMeshRules * | patchRules = nullptr |
Class for spatial white Gaussian noise integration.
The target problem is the linear SPDE \( a(u,v) = F(v)\) with \(F(v) := <\dot{W},v> \), where \(\dot{W}\) is spatial white Gaussian noise. When the Galerkin method is used to discretize this problem into a linear system of equations \(Ax = b\), the RHS is a Gaussian random vector \(b \sim N(0,M)\) whose covariance matrix is the same as the mass matrix \(M_{ij} = (v_i,v_j)\). This property can be ensured if \(b = H w\), where \(HH^{\mathrm{T}} = M\) and each component \(w_i\sim N(0,1)\).
There is much flexibility in how we may wish to define \(H\). In this PR, we define \(H = P^{\mathrm{T}} diag(L_e)\), where \(P\) is the local-to-global dof assembly matrix and \(\mathrm{diag}(L_e)\) is a block-diagonal matrix with \(L_e L_e^{\mathrm{T}} = M_e\), where \(M_e\) is the element mass matrix for element \(e\). A straightforward computation shows that \(HH^{\mathrm{T}} = P^{\mathrm{T}} diag(M_e) P = M\), as necessary.
Definition at line 643 of file lininteg.hpp.
|
inline |
Sets the seed_ of the random number generator. A fixed seed allows for a reproducible sequence of white noise vectors.
Definition at line 661 of file lininteg.hpp.
|
inline |
Sets the MPI communicator comm_ and the seed_ of the random number generator. A fixed seed allows for a reproducible sequence of white noise vectors.
Definition at line 670 of file lininteg.hpp.
|
inline |
Sets the seed_ of the random number generator. A fixed seed allows for a reproducible sequence of white noise vectors.
Definition at line 682 of file lininteg.hpp.
|
inline |
Definition at line 726 of file lininteg.hpp.
|
overridevirtual |
Given a particular Finite Element and a transformation (Tr) computes the element vector, elvect.
Implements mfem::LinearFormIntegrator.
Definition at line 1039 of file lininteg.cpp.
|
virtual |
Reimplemented from mfem::LinearFormIntegrator.
Definition at line 46 of file lininteg.cpp.
|
virtual |
Reimplemented from mfem::LinearFormIntegrator.
Definition at line 49 of file lininteg.cpp.
|
inline |
Resets the array of saved lower triangular Cholesky decomposition matrices. The parameter NE should be the number of elements in the mesh.
Definition at line 711 of file lininteg.hpp.
|
inline |
Saves the lower triangular matrices in the element-wise Cholesky decomposition. The parameter NE should be the number of elements in the mesh.
Definition at line 702 of file lininteg.hpp.
|
inline |
Sets/resets the seed of the random number generator.
Definition at line 689 of file lininteg.hpp.