MFEM
v4.5.2
Finite element discretization library
|
#include <bilininteg.hpp>
Public Member Functions | |
VectorDiffusionIntegrator () | |
VectorDiffusionIntegrator (int vector_dimension) | |
Integrator with unit coefficient for caller-specified vector dimension. More... | |
VectorDiffusionIntegrator (Coefficient &q) | |
VectorDiffusionIntegrator (Coefficient &q, const IntegrationRule *ir) | |
VectorDiffusionIntegrator (Coefficient &q, int vector_dimension) | |
Integrator with scalar coefficient for caller-specified vector dimension. More... | |
VectorDiffusionIntegrator (VectorCoefficient &vq) | |
Integrator with VectorCoefficient . The vector dimension of the FiniteElementSpace is assumed to be the same as the dimension of the Vector . More... | |
VectorDiffusionIntegrator (MatrixCoefficient &mq) | |
Integrator with MatrixCoefficient . The vector dimension of the FiniteElementSpace is assumed to be the same as the dimension of the Matrix . More... | |
virtual void | AssembleElementMatrix (const FiniteElement &el, ElementTransformation &Trans, DenseMatrix &elmat) |
Given a particular Finite Element computes the element matrix elmat. More... | |
virtual void | AssembleElementVector (const FiniteElement &el, ElementTransformation &Tr, const Vector &elfun, Vector &elvect) |
Perform the local action of the BilinearFormIntegrator. Note that the default implementation in the base class is general but not efficient. More... | |
virtual void | AssemblePA (const FiniteElementSpace &fes) |
Method defining partial assembly. More... | |
virtual void | AssembleMF (const FiniteElementSpace &fes) |
Method defining matrix-free assembly. More... | |
virtual void | AssembleDiagonalPA (Vector &diag) |
Assemble diagonal and add it to Vector diag. More... | |
virtual void | AssembleDiagonalMF (Vector &diag) |
Assemble diagonal and add it to Vector diag. More... | |
virtual void | AddMultPA (const Vector &x, Vector &y) const |
Method for partially assembled action. More... | |
virtual void | AddMultMF (const Vector &x, Vector &y) const |
bool | SupportsCeed () const |
Indicates whether this integrator can use a Ceed backend. More... | |
virtual void | AssemblePA (const FiniteElementSpace &fes) |
Method defining partial assembly. More... | |
virtual void | AssemblePA (const FiniteElementSpace &trial_fes, const FiniteElementSpace &test_fes) |
virtual void | AssemblePA (const FiniteElementSpace &fes) |
Method defining partial assembly. More... | |
virtual void | AssemblePA (const FiniteElementSpace &trial_fes, const FiniteElementSpace &test_fes) |
Public Member Functions inherited from mfem::BilinearFormIntegrator | |
virtual void | AssemblePA (const FiniteElementSpace &trial_fes, const FiniteElementSpace &test_fes) |
virtual void | AssemblePAInteriorFaces (const FiniteElementSpace &fes) |
virtual void | AssemblePABoundaryFaces (const FiniteElementSpace &fes) |
virtual void | AssembleDiagonalPA_ADAt (const Vector &D, Vector &diag) |
Assemble diagonal of ADA^T (A is this integrator) and add it to diag. More... | |
virtual void | AddMultTransposePA (const Vector &x, Vector &y) const |
Method for partially assembled transposed action. More... | |
virtual void | AssembleEA (const FiniteElementSpace &fes, Vector &emat, const bool add=true) |
Method defining element assembly. More... | |
virtual void | AddMultTransposeMF (const Vector &x, Vector &y) const |
virtual void | AssembleEAInteriorFaces (const FiniteElementSpace &fes, Vector &ea_data_int, Vector &ea_data_ext, const bool add=true) |
virtual void | AssembleEABoundaryFaces (const FiniteElementSpace &fes, Vector &ea_data_bdr, const bool add=true) |
virtual void | AssembleElementMatrix2 (const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &elmat) |
virtual void | AssembleFaceMatrix (const FiniteElement &el1, const FiniteElement &el2, FaceElementTransformations &Trans, DenseMatrix &elmat) |
virtual void | AssembleFaceMatrix (const FiniteElement &trial_face_fe, const FiniteElement &test_fe1, const FiniteElement &test_fe2, FaceElementTransformations &Trans, DenseMatrix &elmat) |
virtual void | AssembleFaceVector (const FiniteElement &el1, const FiniteElement &el2, FaceElementTransformations &Tr, const Vector &elfun, Vector &elvect) |
Perform the local action of the BilinearFormIntegrator resulting from a face integral term. Note that the default implementation in the base class is general but not efficient. More... | |
virtual void | AssembleElementGrad (const FiniteElement &el, ElementTransformation &Tr, const Vector &elfun, DenseMatrix &elmat) |
Assemble the local gradient matrix. More... | |
virtual void | AssembleFaceGrad (const FiniteElement &el1, const FiniteElement &el2, FaceElementTransformations &Tr, const Vector &elfun, DenseMatrix &elmat) |
Assemble the local action of the gradient of the NonlinearFormIntegrator resulting from a face integral term. More... | |
virtual void | ComputeElementFlux (const FiniteElement &el, ElementTransformation &Trans, Vector &u, const FiniteElement &fluxelem, Vector &flux, bool with_coef=true, const IntegrationRule *ir=NULL) |
Virtual method required for Zienkiewicz-Zhu type error estimators. More... | |
virtual double | ComputeFluxEnergy (const FiniteElement &fluxelem, ElementTransformation &Trans, Vector &flux, Vector *d_energy=NULL) |
Virtual method required for Zienkiewicz-Zhu type error estimators. More... | |
virtual | ~BilinearFormIntegrator () |
virtual void | AssemblePA (const FiniteElementSpace &fes) |
Method defining partial assembly. More... | |
virtual void | AssemblePA (const FiniteElementSpace &trial_fes, const FiniteElementSpace &test_fes) |
Public Member Functions inherited from mfem::NonlinearFormIntegrator | |
virtual void | SetIntRule (const IntegrationRule *ir) |
Prescribe a fixed IntegrationRule to use (when ir != NULL) or let the integrator choose (when ir == NULL). More... | |
void | SetIntegrationRule (const IntegrationRule &ir) |
Prescribe a fixed IntegrationRule to use. More... | |
void | SetPAMemoryType (MemoryType mt) |
const IntegrationRule * | GetIntegrationRule () const |
Get the integration rule of the integrator (possibly NULL). More... | |
virtual double | GetElementEnergy (const FiniteElement &el, ElementTransformation &Tr, const Vector &elfun) |
Compute the local energy. More... | |
virtual void | AssembleGradPA (const Vector &x, const FiniteElementSpace &fes) |
Prepare the integrator for partial assembly (PA) gradient evaluations on the given FE space fes at the state x. More... | |
virtual double | GetLocalStateEnergyPA (const Vector &x) const |
Compute the local (to the MPI rank) energy with partial assembly. More... | |
virtual void | AddMultGradPA (const Vector &x, Vector &y) const |
Method for partially assembled gradient action. More... | |
virtual void | AssembleGradDiagonalPA (Vector &diag) const |
Method for computing the diagonal of the gradient with partial assembly. More... | |
ceed::Operator & | GetCeedOp () |
virtual | ~NonlinearFormIntegrator () |
Protected Attributes | |
Coefficient * | Q = NULL |
VectorCoefficient * | VQ = NULL |
MatrixCoefficient * | MQ = NULL |
const DofToQuad * | maps |
Not owned. More... | |
const GeometricFactors * | geom |
Not owned. More... | |
int | dim |
int | sdim |
int | ne |
int | dofs1D |
int | quad1D |
Vector | pa_data |
Protected Attributes inherited from mfem::NonlinearFormIntegrator | |
const IntegrationRule * | IntRule |
ceed::Operator * | ceedOp |
MemoryType | pa_mt = MemoryType::DEFAULT |
Additional Inherited Members | |
Protected Member Functions inherited from mfem::BilinearFormIntegrator | |
BilinearFormIntegrator (const IntegrationRule *ir=NULL) | |
Protected Member Functions inherited from mfem::NonlinearFormIntegrator | |
NonlinearFormIntegrator (const IntegrationRule *ir=NULL) | |
Integrator for
(Q grad u, grad v) = sum_i (Q grad u_i, grad v_i) e_i e_i^T
for vector FE spaces, where e_i is the unit vector in the i-th direction. The resulting local element matrix is square, of size vdim*dof
, where vdim
is the vector dimension space and dof
is the local degrees of freedom. The integrator is not aware of the true vector dimension and must use VectorCoefficient
, MatrixCoefficient
, or a caller-specified value to determine the vector space. For a scalar coefficient, the caller may manually specify the vector dimension or the vector dimension is assumed to be the spatial dimension (i.e. 2-dimension or 3-dimension).
Definition at line 2804 of file bilininteg.hpp.
|
inline |
Definition at line 2824 of file bilininteg.hpp.
|
inline |
Integrator with unit coefficient for caller-specified vector dimension.
If the vector dimension does not match the true dimension of the space, the resulting element matrix will be mathematically invalid.
Definition at line 2831 of file bilininteg.hpp.
|
inline |
Definition at line 2834 of file bilininteg.hpp.
|
inline |
Definition at line 2837 of file bilininteg.hpp.
|
inline |
Integrator with scalar coefficient for caller-specified vector dimension.
The element matrix is block-diagonal with vdim
copies of the element matrix integrated with the Coefficient
.
If the vector dimension does not match the true dimension of the space, the resulting element matrix will be mathematically invalid.
Definition at line 2848 of file bilininteg.hpp.
|
inline |
Integrator with VectorCoefficient
. The vector dimension of the FiniteElementSpace
is assumed to be the same as the dimension of the Vector
.
The element matrix is block-diagonal and each block is integrated with coefficient q_i.
If the vector dimension does not match the true dimension of the space, the resulting element matrix will be mathematically invalid.
Definition at line 2860 of file bilininteg.hpp.
|
inline |
Integrator with MatrixCoefficient
. The vector dimension of the FiniteElementSpace
is assumed to be the same as the dimension of the Matrix
.
The element matrix is populated in each block. Each block is integrated with coefficient q_ij.
If the vector dimension does not match the true dimension of the space, the resulting element matrix will be mathematically invalid.
Definition at line 2872 of file bilininteg.hpp.
Perform the action of integrator on the input x and add the result to the output y. Both x and y are E-vectors, i.e. they represent the element-wise discontinuous version of the FE space.
This method can be called only after the method AssembleMF() has been called.
Reimplemented from mfem::BilinearFormIntegrator.
Definition at line 52 of file bilininteg_vecdiffusion_mf.cpp.
Method for partially assembled action.
Perform the action of integrator on the input x and add the result to the output y. Both x and y are E-vectors, i.e. they represent the element-wise discontinuous version of the FE space.
This method can be called only after the method AssemblePA() has been called.
Reimplemented from mfem::BilinearFormIntegrator.
Definition at line 546 of file bilininteg_vecdiffusion.cpp.
|
virtual |
Assemble diagonal and add it to Vector diag.
Reimplemented from mfem::BilinearFormIntegrator.
Definition at line 65 of file bilininteg_vecdiffusion_mf.cpp.
|
virtual |
Assemble diagonal and add it to Vector diag.
Reimplemented from mfem::BilinearFormIntegrator.
Definition at line 773 of file bilininteg_vecdiffusion.cpp.
|
virtual |
Given a particular Finite Element computes the element matrix elmat.
Reimplemented from mfem::BilinearFormIntegrator.
Definition at line 2793 of file bilininteg.cpp.
|
virtual |
Perform the local action of the BilinearFormIntegrator. Note that the default implementation in the base class is general but not efficient.
Reimplemented from mfem::BilinearFormIntegrator.
Definition at line 2875 of file bilininteg.cpp.
|
virtual |
Method defining matrix-free assembly.
Used with BilinearFormIntegrators that have different spaces.The result of fully matrix-free assembly is stored internally so that it can be used later in the methods AddMultMF() and AddMultTransposeMF().
Reimplemented from mfem::BilinearFormIntegrator.
Definition at line 22 of file bilininteg_vecdiffusion_mf.cpp.
void mfem::NonlinearFormIntegrator::AssemblePA |
The result of the partial assembly is stored internally so that it can be used later in the methods AddMultPA(). Used with BilinearFormIntegrators that have different spaces.
Definition at line 31 of file nonlininteg.cpp.
void mfem::NonlinearFormIntegrator::AssemblePA |
Method defining partial assembly.
The result of the partial assembly is stored internally so that it can be used later in the methods AddMultPA().
Definition at line 25 of file nonlininteg.cpp.
void mfem::BilinearFormIntegrator::AssemblePA |
Method defining partial assembly.
The result of the partial assembly is stored internally so that it can be used later in the methods AddMultPA() and AddMultTransposePA().
Definition at line 23 of file bilininteg.cpp.
void mfem::BilinearFormIntegrator::AssemblePA |
Used with BilinearFormIntegrators that have different spaces.
Definition at line 29 of file bilininteg.cpp.
|
virtual |
Method defining partial assembly.
The result of the partial assembly is stored internally so that it can be used later in the methods AddMultPA() and AddMultTransposePA().
Reimplemented from mfem::BilinearFormIntegrator.
Definition at line 143 of file bilininteg_vecdiffusion.cpp.
|
inlinevirtual |
Indicates whether this integrator can use a Ceed backend.
Reimplemented from mfem::NonlinearFormIntegrator.
Definition at line 2888 of file bilininteg.hpp.
|
protected |
Definition at line 2814 of file bilininteg.hpp.
|
protected |
Definition at line 2814 of file bilininteg.hpp.
|
protected |
Not owned.
Definition at line 2813 of file bilininteg.hpp.
|
protected |
Not owned.
Definition at line 2812 of file bilininteg.hpp.
|
protected |
Definition at line 2809 of file bilininteg.hpp.
|
protected |
Definition at line 2814 of file bilininteg.hpp.
|
protected |
Definition at line 2815 of file bilininteg.hpp.
|
protected |
Definition at line 2807 of file bilininteg.hpp.
|
protected |
Definition at line 2814 of file bilininteg.hpp.
|
protected |
Definition at line 2814 of file bilininteg.hpp.
|
protected |
Definition at line 2808 of file bilininteg.hpp.