![]() |
MFEM v4.9.0
Finite element discretization library
|
Abstract boundary hyperbolic form integrator, assembling <F̂(u⁻,u_b,x) n, [v]> term for scalar finite elements at the boundary. More...
#include <hyperbolic.hpp>
Public Member Functions | |
| BdrHyperbolicDirichletIntegrator (const NumericalFlux &numFlux, VectorCoefficient &bdrState, const int IntOrderOffset=0, const real_t sign=1.) | |
| Construct a new BdrHyperbolicDirichletIntegrator object. | |
| void | ResetMaxCharSpeed () |
| Reset the maximum characteristic speed to zero. | |
| real_t | GetMaxCharSpeed () const |
| Get the maximum characteristic speed. | |
| const FluxFunction & | GetFluxFunction () const |
| Get the associated flux function. | |
| void | AssembleFaceVector (const FiniteElement &el1, const FiniteElement &el2, FaceElementTransformations &Tr, const Vector &elfun, Vector &elvect) override |
| Implements <-F̂(u⁻,u_b,x) n, [v]> with abstract F̂ computed by NumericalFlux::Eval() of the numerical flux object. | |
| void | AssembleFaceGrad (const FiniteElement &el1, const FiniteElement &el2, FaceElementTransformations &Tr, const Vector &elfun, DenseMatrix &elmat) override |
| Implements <-Ĵ(u⁻,u_b,x) n, [v]> with abstract Ĵ computed by NumericalFlux::Grad() of the numerical flux object. | |
Public Member Functions inherited from mfem::NonlinearFormIntegrator | |
| void | SetIntegrationMode (Mode m) |
| bool | Patchwise () const |
| void | SetPAMemoryType (MemoryType mt) |
| virtual void | AssembleElementVector (const FiniteElement &el, ElementTransformation &Tr, const Vector &elfun, Vector &elvect) |
| Perform the local action of the NonlinearFormIntegrator. | |
| virtual void | AssembleElementGrad (const FiniteElement &el, ElementTransformation &Tr, const Vector &elfun, DenseMatrix &elmat) |
| Assemble the local gradient matrix. | |
| 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 () |
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. | |
Public Attributes | |
| const int | num_equations |
Additional Inherited Members | |
Public Types inherited from mfem::NonlinearFormIntegrator | |
| enum | Mode { ELEMENTWISE = 0 , PATCHWISE = 1 , PATCHWISE_REDUCED = 2 } |
Protected Member Functions inherited from mfem::NonlinearFormIntegrator | |
| NonlinearFormIntegrator (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::NonlinearFormIntegrator | |
| Mode | integrationMode = Mode::ELEMENTWISE |
| ceed::Operator * | ceedOp |
| MemoryType | pa_mt = MemoryType::DEFAULT |
Protected Attributes inherited from mfem::Integrator | |
| const IntegrationRule * | IntRule |
| NURBSMeshRules * | patchRules = nullptr |
Abstract boundary hyperbolic form integrator, assembling <F̂(u⁻,u_b,x) n, [v]> term for scalar finite elements at the boundary.
This form integrator is coupled with a NumericalFlux that implements the numerical flux F̂ at the boundary faces. The flux F is obtained from the FluxFunction assigned to the aforementioned NumericalFlux with the given boundary coefficient for the state u_b.
Note the class can be used for imposing conditions on interior interfaces.
Definition at line 426 of file hyperbolic.hpp.
| mfem::BdrHyperbolicDirichletIntegrator::BdrHyperbolicDirichletIntegrator | ( | const NumericalFlux & | numFlux, |
| VectorCoefficient & | bdrState, | ||
| const int | IntOrderOffset = 0, | ||
| const real_t | sign = 1. ) |
Construct a new BdrHyperbolicDirichletIntegrator object.
| [in] | numFlux | numerical flux |
| [in] | bdrState | boundary state coefficient |
| [in] | IntOrderOffset | integration order offset |
| [in] | sign | sign of the convection term |
Definition at line 420 of file hyperbolic.cpp.
|
overridevirtual |
Implements <-Ĵ(u⁻,u_b,x) n, [v]> with abstract Ĵ computed by NumericalFlux::Grad() of the numerical flux object.
| [in] | el1 | finite element of the interior element |
| [in] | el2 | not used |
| [in] | Tr | face element transformations |
| [in] | elfun | local coefficient of basis for the interior element |
| [out] | elmat | evaluated Jacobian matrix <-Ĵ(u⁻,u_b,x) n, [v]> |
Reimplemented from mfem::NonlinearFormIntegrator.
Definition at line 524 of file hyperbolic.cpp.
|
overridevirtual |
Implements <-F̂(u⁻,u_b,x) n, [v]> with abstract F̂ computed by NumericalFlux::Eval() of the numerical flux object.
| [in] | el1 | finite element of the interior element |
| [in] | el2 | not used |
| [in] | Tr | face element transformations |
| [in] | elfun | local coefficient of basis for the interior element |
| [out] | elvect | evaluated dual vector <-F̂(u⁻,u_b,x) n, [v]> |
Reimplemented from mfem::NonlinearFormIntegrator.
Definition at line 445 of file hyperbolic.cpp.
|
inline |
Get the associated flux function.
Definition at line 472 of file hyperbolic.hpp.
|
inline |
Get the maximum characteristic speed.
Definition at line 469 of file hyperbolic.hpp.
|
inline |
Reset the maximum characteristic speed to zero.
Definition at line 466 of file hyperbolic.hpp.
| const int mfem::BdrHyperbolicDirichletIntegrator::num_equations |
Definition at line 449 of file hyperbolic.hpp.