![]() |
MFEM v4.7.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. | |
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. | |
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) |
![]() | |
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 void | SetIntRule (const IntegrationRule *ir) |
const IntegrationRule * | GetIntRule () |
virtual | ~LinearFormIntegrator () |
Additional Inherited Members | |
![]() | |
LinearFormIntegrator (const IntegrationRule *ir=NULL) | |
![]() | |
const IntegrationRule * | IntRule |
Class for spatial white Gaussian noise integration.
The target problem is the linear SPDE
There is much flexibility in how we may wish to define
Definition at line 627 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 645 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 654 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 666 of file lininteg.hpp.
|
inline |
Definition at line 710 of file lininteg.hpp.
|
virtual |
Given a particular Finite Element and a transformation (Tr) computes the element vector, elvect.
Implements mfem::LinearFormIntegrator.
Definition at line 1005 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 695 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 686 of file lininteg.hpp.
|
inline |
Sets/resets the seed of the random number generator.
Definition at line 673 of file lininteg.hpp.