MFEM v4.7.0
Finite element discretization library
|
This class is used to express the local action of a general nonlinear finite element operator. In addition it may provide the capability to assemble the local gradient operator and to compute the local energy. More...
#include <nonlininteg.hpp>
Public Types | |
enum | Mode { ELEMENTWISE = 0 , PATCHWISE = 1 , PATCHWISE_REDUCED = 2 } |
Public Member Functions | |
virtual void | SetIntRule (const IntegrationRule *ir) |
Prescribe a fixed IntegrationRule to use (when ir != NULL) or let the integrator choose (when ir == NULL). | |
void | SetIntegrationMode (Mode m) |
void | SetNURBSPatchIntRule (NURBSMeshRules *pr) |
For patchwise integration, SetNURBSPatchIntRule must be called. | |
bool | HasNURBSPatchIntRule () const |
bool | Patchwise () const |
void | SetIntegrationRule (const IntegrationRule &ir) |
Prescribe a fixed IntegrationRule to use. | |
void | SetPAMemoryType (MemoryType mt) |
const IntegrationRule * | GetIntegrationRule () const |
Get the integration rule of the integrator (possibly NULL). | |
virtual void | AssembleElementVector (const FiniteElement &el, ElementTransformation &Tr, const Vector &elfun, Vector &elvect) |
Perform the local action of the NonlinearFormIntegrator. | |
virtual void | AssembleFaceVector (const FiniteElement &el1, const FiniteElement &el2, FaceElementTransformations &Tr, const Vector &elfun, Vector &elvect) |
Perform the local action of the NonlinearFormIntegrator resulting from a face integral term. | |
virtual void | AssembleElementGrad (const FiniteElement &el, ElementTransformation &Tr, const Vector &elfun, DenseMatrix &elmat) |
Assemble the local gradient matrix. | |
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. | |
virtual real_t | GetElementEnergy (const FiniteElement &el, ElementTransformation &Tr, const Vector &elfun) |
Compute the local energy. | |
virtual void | AssemblePA (const FiniteElementSpace &fes) |
Method defining partial assembly. | |
virtual void | AssemblePA (const FiniteElementSpace &trial_fes, const FiniteElementSpace &test_fes) |
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. | |
virtual real_t | GetLocalStateEnergyPA (const Vector &x) const |
Compute the local (to the MPI rank) energy with partial assembly. | |
virtual void | AddMultPA (const Vector &x, Vector &y) const |
Method for partially assembled action. | |
virtual void | AddMultGradPA (const Vector &x, Vector &y) const |
Method for partially assembled gradient action. | |
virtual void | AssembleGradDiagonalPA (Vector &diag) const |
Method for computing the diagonal of the gradient with partial assembly. | |
virtual bool | SupportsCeed () const |
Indicates whether this integrator can use a Ceed backend. | |
virtual void | AssembleMF (const FiniteElementSpace &fes) |
Method defining fully unassembled operator. | |
virtual void | AddMultMF (const Vector &x, Vector &y) const |
ceed::Operator & | GetCeedOp () |
virtual | ~NonlinearFormIntegrator () |
Protected Member Functions | |
NonlinearFormIntegrator (const IntegrationRule *ir=NULL) | |
Protected Attributes | |
const IntegrationRule * | IntRule |
Mode | integrationMode = Mode::ELEMENTWISE |
NURBSMeshRules * | patchRules = nullptr |
ceed::Operator * | ceedOp |
MemoryType | pa_mt = MemoryType::DEFAULT |
This class is used to express the local action of a general nonlinear finite element operator. In addition it may provide the capability to assemble the local gradient operator and to compute the local energy.
Definition at line 27 of file nonlininteg.hpp.
Enumerator | |
---|---|
ELEMENTWISE | Element-wise integration (default) |
PATCHWISE | Patch-wise integration (NURBS meshes) |
PATCHWISE_REDUCED | Patch-wise integration (NURBS meshes) with reduced integration rules. |
Definition at line 30 of file nonlininteg.hpp.
|
inlineprotected |
Definition at line 51 of file nonlininteg.hpp.
|
inlinevirtual |
Definition at line 169 of file nonlininteg.hpp.
Method for partially assembled gradient action.
All arguments are E-vectors. This method can be called only after the method AssembleGradPA() has been called.
[in] | x | The gradient Operator is applied to the Vector x. |
[in,out] | y | The result Vector: \( y += G x \). |
Reimplemented in mfem::TMOP_Integrator, and mfem::TMOPComboIntegrator.
Definition at line 51 of file nonlininteg.cpp.
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 in mfem::BilinearFormIntegrator, mfem::ConvectionIntegrator, mfem::DiffusionIntegrator, mfem::MassIntegrator, mfem::SumIntegrator, mfem::VectorConvectionNLFIntegrator, mfem::VectorDiffusionIntegrator, and mfem::VectorMassIntegrator.
Definition at line 69 of file nonlininteg.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 in mfem::BilinearFormIntegrator, mfem::ConvectionIntegrator, mfem::CurlCurlIntegrator, mfem::DGTraceIntegrator, mfem::DiffusionIntegrator, mfem::DivDivIntegrator, mfem::ElasticityComponentIntegrator, mfem::ElasticityIntegrator, mfem::GradientIntegrator, mfem::GradientInterpolator, mfem::IdentityInterpolator, mfem::MassIntegrator, mfem::MixedScalarCurlIntegrator, mfem::MixedVectorCurlIntegrator, mfem::MixedVectorGradientIntegrator, mfem::MixedVectorWeakCurlIntegrator, mfem::SumIntegrator, mfem::TMOP_Integrator, mfem::TMOPComboIntegrator, mfem::TransposeIntegrator, mfem::VectorConvectionNLFIntegrator, mfem::VectorDiffusionIntegrator, mfem::VectorDivergenceIntegrator, mfem::VectorFEDivergenceIntegrator, mfem::VectorFEMassIntegrator, and mfem::VectorMassIntegrator.
Definition at line 45 of file nonlininteg.cpp.
|
virtual |
Assemble the local gradient matrix.
Reimplemented in mfem::BilinearFormIntegrator, mfem::common::PUMPLaplacian, mfem::common::ScreenedPoisson, mfem::ConvectiveVectorConvectionNLFIntegrator, mfem::HyperelasticNLFIntegrator, mfem::pLaplace, mfem::pLaplaceAD< CQVectAutoDiff >, mfem::pLaplaceSL< sizeres >, mfem::SkewSymmetricVectorConvectionNLFIntegrator, mfem::TMOP_Integrator, mfem::TMOPComboIntegrator, and mfem::VectorConvectionNLFIntegrator.
Definition at line 91 of file nonlininteg.cpp.
|
virtual |
Perform the local action of the NonlinearFormIntegrator.
Reimplemented in mfem::BilinearFormIntegrator, mfem::common::PUMPLaplacian, mfem::common::ScreenedPoisson, mfem::DiffusionIntegrator, mfem::HyperbolicFormIntegrator, mfem::HyperelasticNLFIntegrator, mfem::pLaplace, mfem::pLaplaceAD< CQVectAutoDiff >, mfem::pLaplaceSL< sizeres >, mfem::TMOP_Integrator, mfem::TMOPComboIntegrator, mfem::VectorConvectionNLFIntegrator, and mfem::VectorDiffusionIntegrator.
Definition at line 75 of file nonlininteg.cpp.
|
virtual |
Assemble the local action of the gradient of the NonlinearFormIntegrator resulting from a face integral term.
Reimplemented in mfem::BilinearFormIntegrator.
Definition at line 99 of file nonlininteg.cpp.
|
virtual |
Perform the local action of the NonlinearFormIntegrator resulting from a face integral term.
Reimplemented in mfem::BilinearFormIntegrator, and mfem::HyperbolicFormIntegrator.
Definition at line 83 of file nonlininteg.cpp.
|
virtual |
Method for computing the diagonal of the gradient with partial assembly.
The result Vector diag is an E-Vector. This method can be called only after the method AssembleGradPA() has been called.
[in,out] | diag | The result Vector: \( diag += diag(G) \). |
Reimplemented in mfem::TMOP_Integrator, and mfem::TMOPComboIntegrator.
Definition at line 57 of file nonlininteg.cpp.
|
virtual |
Prepare the integrator for partial assembly (PA) gradient evaluations on the given FE space fes at the state x.
The result of the partial assembly is stored internally so that it can be used later in the methods AddMultGradPA() and AssembleGradDiagonalPA(). The state Vector x is an E-vector.
Reimplemented in mfem::TMOP_Integrator, and mfem::TMOPComboIntegrator.
Definition at line 38 of file nonlininteg.cpp.
|
virtual |
Method defining fully unassembled operator.
Reimplemented in mfem::BilinearFormIntegrator, mfem::ConvectionIntegrator, mfem::DiffusionIntegrator, mfem::MassIntegrator, mfem::SumIntegrator, mfem::VectorConvectionNLFIntegrator, mfem::VectorDiffusionIntegrator, and mfem::VectorMassIntegrator.
Definition at line 63 of file nonlininteg.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().
Reimplemented in mfem::BilinearFormIntegrator, mfem::ConvectionIntegrator, mfem::CurlCurlIntegrator, mfem::DGDiffusionIntegrator, mfem::DiffusionIntegrator, mfem::DivDivIntegrator, mfem::ElasticityComponentIntegrator, mfem::ElasticityIntegrator, mfem::GradientIntegrator, mfem::GradientInterpolator, mfem::IdentityInterpolator, mfem::MassIntegrator, mfem::MixedScalarCurlIntegrator, mfem::MixedVectorCurlIntegrator, mfem::MixedVectorGradientIntegrator, mfem::MixedVectorWeakCurlIntegrator, mfem::SumIntegrator, mfem::TMOP_Integrator, mfem::TMOPComboIntegrator, mfem::TransposeIntegrator, mfem::VectorConvectionNLFIntegrator, mfem::VectorDiffusionIntegrator, mfem::VectorDivergenceIntegrator, mfem::VectorFEDivergenceIntegrator, mfem::VectorFEMassIntegrator, and mfem::VectorMassIntegrator.
Definition at line 25 of file nonlininteg.cpp.
|
virtual |
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.
Reimplemented in mfem::BilinearFormIntegrator, mfem::ConvectionIntegrator, mfem::CurlCurlIntegrator, mfem::DGDiffusionIntegrator, mfem::DiffusionIntegrator, mfem::DivDivIntegrator, mfem::ElasticityComponentIntegrator, mfem::ElasticityIntegrator, mfem::GradientIntegrator, mfem::GradientInterpolator, mfem::IdentityInterpolator, mfem::MassIntegrator, mfem::MixedScalarCurlIntegrator, mfem::MixedVectorCurlIntegrator, mfem::MixedVectorGradientIntegrator, mfem::MixedVectorWeakCurlIntegrator, mfem::SumIntegrator, mfem::TMOP_Integrator, mfem::TMOPComboIntegrator, mfem::TransposeIntegrator, mfem::VectorConvectionNLFIntegrator, mfem::VectorDiffusionIntegrator, mfem::VectorDivergenceIntegrator, mfem::VectorFEDivergenceIntegrator, mfem::VectorFEMassIntegrator, and mfem::VectorMassIntegrator.
Definition at line 31 of file nonlininteg.cpp.
|
inline |
Definition at line 167 of file nonlininteg.hpp.
|
virtual |
Compute the local energy.
Reimplemented in mfem::common::PUMPLaplacian, mfem::common::ScreenedPoisson, mfem::HyperelasticNLFIntegrator, mfem::pLaplace, mfem::pLaplaceAD< CQVectAutoDiff >, mfem::pLaplaceSL< sizeres >, mfem::TMOP_Integrator, mfem::TMOPComboIntegrator, and mfem::VectorCurlCurlIntegrator.
Definition at line 108 of file nonlininteg.cpp.
|
inline |
Get the integration rule of the integrator (possibly NULL).
Definition at line 75 of file nonlininteg.hpp.
Compute the local (to the MPI rank) energy with partial assembly.
Here the state x is an E-vector. This method can be called only after the method AssemblePA() has been called.
Reimplemented in mfem::TMOP_Integrator, and mfem::TMOPComboIntegrator.
Definition at line 18 of file nonlininteg.cpp.
|
inline |
Definition at line 63 of file nonlininteg.hpp.
|
inline |
Definition at line 65 of file nonlininteg.hpp.
|
inline |
Definition at line 59 of file nonlininteg.hpp.
|
inline |
Prescribe a fixed IntegrationRule to use.
Definition at line 68 of file nonlininteg.hpp.
|
inlinevirtual |
Prescribe a fixed IntegrationRule to use (when ir != NULL) or let the integrator choose (when ir == NULL).
Reimplemented in mfem::InverseIntegrator, mfem::LumpedIntegrator, mfem::SumIntegrator, and mfem::TransposeIntegrator.
Definition at line 57 of file nonlininteg.hpp.
|
inline |
For patchwise integration, SetNURBSPatchIntRule must be called.
Definition at line 62 of file nonlininteg.hpp.
|
inline |
Set the memory type used for GeometricFactors and other large allocations in PA extensions.
Definition at line 72 of file nonlininteg.hpp.
|
inlinevirtual |
Indicates whether this integrator can use a Ceed backend.
Reimplemented in mfem::ConvectionIntegrator, mfem::DiffusionIntegrator, mfem::MassIntegrator, mfem::VectorDiffusionIntegrator, and mfem::VectorMassIntegrator.
Definition at line 154 of file nonlininteg.hpp.
|
protected |
Definition at line 47 of file nonlininteg.hpp.
|
protected |
Definition at line 41 of file nonlininteg.hpp.
|
protected |
Definition at line 39 of file nonlininteg.hpp.
|
protected |
Definition at line 49 of file nonlininteg.hpp.
|
protected |
Definition at line 44 of file nonlininteg.hpp.