MFEM
v4.2.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 Member Functions | |
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 &irule) |
Prescribe a fixed IntegrationRule to use. More... | |
virtual void | AssembleElementVector (const FiniteElement &el, ElementTransformation &Tr, const Vector &elfun, Vector &elvect) |
Perform the local action of the NonlinearFormIntegrator. More... | |
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. 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 double | GetElementEnergy (const FiniteElement &el, ElementTransformation &Tr, const Vector &elfun) |
Compute the local energy. 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 | AddMultPA (const Vector &x, Vector &y) const |
Method for partially assembled action. More... | |
virtual | ~NonlinearFormIntegrator () |
Protected Member Functions | |
NonlinearFormIntegrator (const IntegrationRule *ir=NULL) | |
Protected Attributes | |
const IntegrationRule * | IntRule |
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 26 of file nonlininteg.hpp.
|
inlineprotected |
Definition at line 31 of file nonlininteg.hpp.
|
inlinevirtual |
Definition at line 91 of file nonlininteg.hpp.
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::DGTraceIntegrator, mfem::VectorDiffusionIntegrator, mfem::DivDivIntegrator, mfem::VectorDivergenceIntegrator, mfem::VectorFEMassIntegrator, mfem::CurlCurlIntegrator, mfem::VectorFEDivergenceIntegrator, mfem::VectorMassIntegrator, mfem::ConvectionIntegrator, mfem::MassIntegrator, mfem::DiffusionIntegrator, mfem::GradientIntegrator, mfem::MixedVectorWeakCurlIntegrator, mfem::MixedVectorCurlIntegrator, mfem::MixedVectorGradientIntegrator, mfem::VectorConvectionNLFIntegrator, mfem::TransposeIntegrator, and mfem::BilinearFormIntegrator.
Definition at line 31 of file nonlininteg.cpp.
|
virtual |
Assemble the local gradient matrix.
Reimplemented in mfem::TMOPComboIntegrator, mfem::TMOP_Integrator, mfem::VectorConvectionNLFIntegrator, mfem::HyperelasticNLFIntegrator, and mfem::BilinearFormIntegrator.
Definition at line 53 of file nonlininteg.cpp.
|
virtual |
Perform the local action of the NonlinearFormIntegrator.
Reimplemented in mfem::VectorDiffusionIntegrator, mfem::DiffusionIntegrator, mfem::TMOPComboIntegrator, mfem::TMOP_Integrator, mfem::VectorConvectionNLFIntegrator, mfem::HyperelasticNLFIntegrator, and mfem::BilinearFormIntegrator.
Definition at line 37 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 61 of file nonlininteg.cpp.
|
virtual |
Perform the local action of the NonlinearFormIntegrator resulting from a face integral term.
Reimplemented in mfem::BilinearFormIntegrator, and FaceIntegrator.
Definition at line 45 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::VectorDiffusionIntegrator, mfem::DivDivIntegrator, mfem::VectorFEMassIntegrator, mfem::CurlCurlIntegrator, mfem::VectorMassIntegrator, mfem::ConvectionIntegrator, mfem::MassIntegrator, mfem::DiffusionIntegrator, mfem::VectorConvectionNLFIntegrator, mfem::TransposeIntegrator, and mfem::BilinearFormIntegrator.
Definition at line 18 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::VectorDivergenceIntegrator, mfem::VectorFEMassIntegrator, mfem::VectorFEDivergenceIntegrator, mfem::GradientIntegrator, mfem::MixedVectorWeakCurlIntegrator, mfem::MixedVectorCurlIntegrator, mfem::MixedVectorGradientIntegrator, and mfem::BilinearFormIntegrator.
Definition at line 24 of file nonlininteg.cpp.
|
virtual |
Compute the local energy.
Reimplemented in mfem::VectorCurlCurlIntegrator, mfem::TMOPComboIntegrator, mfem::TMOP_Integrator, and mfem::HyperelasticNLFIntegrator.
Definition at line 70 of file nonlininteg.cpp.
|
inline |
Prescribe a fixed IntegrationRule to use.
Definition at line 40 of file nonlininteg.hpp.
|
inline |
Prescribe a fixed IntegrationRule to use (when ir != NULL) or let the integrator choose (when ir == NULL).
Definition at line 37 of file nonlininteg.hpp.
|
protected |
Definition at line 29 of file nonlininteg.hpp.