MFEM
v4.5.2
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. More... | |
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. More... | |
WhiteGaussianNoiseDomainLFIntegrator (int seed_=0) | |
Sets the seed_ of the random number generator. A fixed seed allows for a reproducible sequence of white noise vectors. More... | |
void | SetSeed (int seed) |
Sets/resets the seed of the random number generator. More... | |
virtual void | AssembleRHSElementVect (const FiniteElement &el, ElementTransformation &Tr, Vector &elvect) |
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. More... | |
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. More... | |
~WhiteGaussianNoiseDomainLFIntegrator () | |
virtual void | AssembleRHSElementVect (const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)=0 |
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 () |
Method probing for assembly on device. More... | |
virtual void | AssembleDevice (const FiniteElementSpace &fes, const Array< int > &markers, Vector &b) |
Method defining assembly on device. More... | |
virtual void | AssembleRHSElementVect (const FiniteElement &el, FaceElementTransformations &Tr, Vector &elvect) |
virtual void | AssembleRHSElementVect (const FiniteElement &el1, const FiniteElement &el2, FaceElementTransformations &Tr, Vector &elvect) |
virtual void | SetIntRule (const IntegrationRule *ir) |
const IntegrationRule * | GetIntRule () |
virtual | ~LinearFormIntegrator () |
Additional Inherited Members | |
Protected Member Functions inherited from mfem::LinearFormIntegrator | |
LinearFormIntegrator (const IntegrationRule *ir=NULL) | |
Protected Attributes inherited from mfem::LinearFormIntegrator | |
const IntegrationRule * | IntRule |
Class for spatial white Gaussian noise integration.
The target problem is the linear SPDE a(u,v) = F(v) with F(v) := <Ẇ,v>, where Ẇ 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~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ᵀ = M and each component w_i~N(0,1).
There is much flexibility in how we may wish to define H. In this PR, we define H = Pᵀ diag(L_e), where P is the local-to-global dof assembly matrix and diag(L_e) is a block-diagonal matrix with L_e L_eᵀ = M_e, where M_e is the element mass matrix for element e. A straightforward computation shows that HHᵀ = Pᵀ diag(M_e) P = M, as necessary.
Definition at line 615 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 633 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 642 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 654 of file lininteg.hpp.
|
inline |
Definition at line 698 of file lininteg.hpp.
virtual void mfem::LinearFormIntegrator::AssembleRHSElementVect |
Given a particular Finite Element and a transformation (Tr) computes the element vector, elvect.
void mfem::LinearFormIntegrator::AssembleRHSElementVect |
Definition at line 32 of file lininteg.cpp.
void mfem::LinearFormIntegrator::AssembleRHSElementVect |
Definition at line 26 of file lininteg.cpp.
|
virtual |
Given a particular Finite Element and a transformation (Tr) computes the element vector, elvect.
Implements mfem::LinearFormIntegrator.
Definition at line 1007 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 683 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 674 of file lininteg.hpp.
|
inline |
Sets/resets the seed of the random number generator.
Definition at line 661 of file lininteg.hpp.