![]() |
MFEM v4.9.0
Finite element discretization library
|
Abstract boundary hyperbolic linear form integrator, assembling <ɑ/2 F(u,x) n - β |F(u,x) n|, v> terms for scalar finite elements. More...
#include <hyperbolic.hpp>
Public Member Functions | |
| BoundaryHyperbolicFlowIntegrator (const FluxFunction &flux, VectorCoefficient &u, real_t alpha=-1., int IntOrderOffset=0) | |
| Construct a new BoundaryHyperbolicFlowIntegrator object. | |
| BoundaryHyperbolicFlowIntegrator (const FluxFunction &flux, VectorCoefficient &u, real_t alpha, real_t beta, int IntOrderOffset=0) | |
| Construct a new BoundaryHyperbolicFlowIntegrator 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 | AssembleRHSElementVect (const FiniteElement &el, ElementTransformation &Tr, Vector &elvect) override |
| void | AssembleRHSElementVect (const FiniteElement &el, FaceElementTransformations &Tr, Vector &elvect) override |
| Implements <-F(u,x) n, v> with abstract F computed by FluxFunction::ComputeFluxDotN() of the flux function object. | |
| virtual void | AssembleRHSElementVect (const FiniteElement &el1, const FiniteElement &el2, FaceElementTransformations &Tr, Vector &elvect) |
Public Member Functions inherited from mfem::LinearFormIntegrator | |
| virtual bool | SupportsDevice () const |
| Method probing for assembly on device. | |
| virtual void | AssembleDevice (const FiniteElementSpace &fes, const Array< int > &markers, Vector &b) |
| Method defining assembly on device. | |
| 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::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::Integrator | |
| const IntegrationRule * | IntRule |
| NURBSMeshRules * | patchRules = nullptr |
Abstract boundary hyperbolic linear form integrator, assembling <ɑ/2 F(u,x) n - β |F(u,x) n|, v> terms for scalar finite elements.
This form integrator is coupled with a FluxFunction that evaluates the flux F at the boundary.
Note the upwinding is performed component-wise. For general boundary integration with a numerical flux, see BdrHyperbolicDirichletIntegrator.
Definition at line 515 of file hyperbolic.hpp.
|
inline |
Construct a new BoundaryHyperbolicFlowIntegrator object.
| [in] | flux | flux function |
| [in] | u | vector state coefficient |
| [in] | alpha | ɑ coefficient (β = ɑ/2) |
| [in] | IntOrderOffset | integration order offset |
Definition at line 542 of file hyperbolic.hpp.
| mfem::BoundaryHyperbolicFlowIntegrator::BoundaryHyperbolicFlowIntegrator | ( | const FluxFunction & | flux, |
| VectorCoefficient & | u, | ||
| real_t | alpha, | ||
| real_t | beta, | ||
| int | IntOrderOffset = 0 ) |
Construct a new BoundaryHyperbolicFlowIntegrator object.
| [in] | flux | flux function |
| [in] | u | vector state coefficient |
| [in] | alpha | ɑ coefficient |
| [in] | beta | β coefficient |
| [in] | IntOrderOffset | integration order offset |
Definition at line 605 of file hyperbolic.cpp.
|
overridevirtual |
Implements mfem::LinearFormIntegrator.
Definition at line 621 of file hyperbolic.cpp.
|
overridevirtual |
Implements <-F(u,x) n, v> with abstract F computed by FluxFunction::ComputeFluxDotN() of the flux function object.
| [in] | el | finite element |
| [in] | Tr | face element transformations |
| [out] | elvect | evaluated dual vector <F(u,x) n, v> |
Reimplemented from mfem::LinearFormIntegrator.
Definition at line 630 of file hyperbolic.cpp.
|
virtual |
Reimplemented from mfem::LinearFormIntegrator.
Definition at line 51 of file lininteg.cpp.
|
inline |
Get the associated flux function.
Definition at line 565 of file hyperbolic.hpp.
|
inline |
Get the maximum characteristic speed.
Definition at line 562 of file hyperbolic.hpp.
|
inline |
Reset the maximum characteristic speed to zero.
Definition at line 559 of file hyperbolic.hpp.