![]() |
MFEM v4.8.0
Finite element discretization library
|
#include <lininteg.hpp>
Public Member Functions | |
| VectorDomainLFIntegrator (VectorCoefficient &QF) | |
| Constructs a domain integrator with a given VectorCoefficient. | |
| bool | SupportsDevice () const override |
| Method probing for assembly on device. | |
| void | AssembleDevice (const FiniteElementSpace &fes, const Array< int > &markers, Vector &b) override |
| Method defining assembly on device. | |
| void | AssembleRHSElementVect (const FiniteElement &el, ElementTransformation &Tr, Vector &elvect) override |
| void | AssembleDeltaElementVect (const FiniteElement &fe, ElementTransformation &Trans, Vector &elvect) override |
| Assemble the delta coefficient at the IntegrationPoint set in Trans which is assumed to map to the delta coefficient center. | |
| 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::DeltaLFIntegrator | |
| bool | IsDelta () const |
| Returns true if the derived class instance uses a delta coefficient. | |
| void | GetDeltaCenter (Vector ¢er) |
| Returns the center of the delta coefficient. | |
Public Member Functions inherited from mfem::LinearFormIntegrator | |
| virtual | ~LinearFormIntegrator () |
Public Member Functions inherited from mfem::Integrator | |
| Integrator (const IntegrationRule *ir=NULL) | |
| Create a new Integrator, optionally providing a prescribed quadrature rule to use in assembly. | |
| virtual void | SetIntRule (const IntegrationRule *ir) |
| Prescribe a fixed IntegrationRule to use, or set to null to let the integrator choose an appropriate rule. | |
| void | SetIntegrationRule (const IntegrationRule &ir) |
| Prescribe a fixed IntegrationRule to use. Sets the NURBS patch integration rule to null. | |
| void | SetNURBSPatchIntRule (NURBSMeshRules *pr) |
| Sets an integration rule for use on NURBS patches. | |
| bool | HasNURBSPatchIntRule () const |
| Check if a NURBS patch integration rule has been set. | |
| const IntegrationRule * | GetIntRule () const |
| Directly return the IntRule pointer (possibly null) without checking for NURBS patch rules or falling back on a default. | |
| const IntegrationRule * | GetIntegrationRule () const |
| Equivalent to GetIntRule, but retained for backward compatibility with applications. | |
Additional Inherited Members | |
Protected Member Functions inherited from mfem::DeltaLFIntegrator | |
| DeltaLFIntegrator (Coefficient &q, const IntegrationRule *ir=NULL) | |
| This constructor should be used by derived classes that use a scalar DeltaCoefficient. | |
| DeltaLFIntegrator (VectorCoefficient &vq, const IntegrationRule *ir=NULL) | |
| This constructor should be used by derived classes that use a VectorDeltaCoefficient. | |
Protected Member Functions inherited from mfem::LinearFormIntegrator | |
| LinearFormIntegrator (const IntegrationRule *ir=NULL) | |
Protected Member Functions inherited from mfem::Integrator | |
| const IntegrationRule * | GetIntegrationRule (const FiniteElement &trial_fe, const FiniteElement &test_fe, const ElementTransformation &trans) const |
| Returns an integration rule based on the arguments and internal state of the Integrator object. | |
| const IntegrationRule * | GetIntegrationRule (const FiniteElement &el, const ElementTransformation &trans) const |
| Returns an integration rule based on the arguments and internal state. (Version for identical trial_fe and test_fe) | |
| virtual const IntegrationRule * | GetDefaultIntegrationRule (const FiniteElement &trial_fe, const FiniteElement &test_fe, const ElementTransformation &trans) const |
| Subclasses should override to choose a default integration rule. | |
Protected Attributes inherited from mfem::DeltaLFIntegrator | |
| DeltaCoefficient * | delta |
| VectorDeltaCoefficient * | vec_delta |
Protected Attributes inherited from mfem::Integrator | |
| const IntegrationRule * | IntRule |
| NURBSMeshRules * | patchRules = nullptr |
Class for domain integration of \( L(v) := (f, v) \), where \( f = (f_1,\dots,f_n)\) and \( v = (v_1,\dots,v_n) \).
Definition at line 251 of file lininteg.hpp.
|
inline |
Constructs a domain integrator with a given VectorCoefficient.
Definition at line 259 of file lininteg.hpp.
|
overridevirtual |
Assemble the delta coefficient at the IntegrationPoint set in Trans which is assumed to map to the delta coefficient center.
Implements mfem::DeltaLFIntegrator.
Definition at line 311 of file lininteg.cpp.
|
overridevirtual |
Method defining assembly on device.
Reimplemented from mfem::LinearFormIntegrator.
Definition at line 255 of file lininteg_domain.cpp.
|
overridevirtual |
Given a particular Finite Element and a transformation (Tr) computes the element right hand side element vector, elvect.
Implements mfem::LinearFormIntegrator.
Definition at line 269 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.
|
inlineoverridevirtual |
Method probing for assembly on device.
Reimplemented from mfem::LinearFormIntegrator.
Definition at line 262 of file lininteg.hpp.