MFEM v4.9.0
Finite element discretization library
Loading...
Searching...
No Matches
mfem::BoundaryHyperbolicFlowIntegrator Class Reference

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>

Inheritance diagram for mfem::BoundaryHyperbolicFlowIntegrator:
[legend]
Collaboration diagram for mfem::BoundaryHyperbolicFlowIntegrator:
[legend]

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 FluxFunctionGetFluxFunction () 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 IntegrationRuleGetIntRule () const
 Directly return the IntRule pointer (possibly null) without checking for NURBS patch rules or falling back on a default.
 
const IntegrationRuleGetIntegrationRule () 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 IntegrationRuleGetIntegrationRule (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 IntegrationRuleGetIntegrationRule (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 IntegrationRuleGetDefaultIntegrationRule (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 IntegrationRuleIntRule
 
NURBSMeshRulespatchRules = nullptr
 

Detailed Description

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.

Constructor & Destructor Documentation

◆ BoundaryHyperbolicFlowIntegrator() [1/2]

mfem::BoundaryHyperbolicFlowIntegrator::BoundaryHyperbolicFlowIntegrator ( const FluxFunction & flux,
VectorCoefficient & u,
real_t alpha = -1.,
int IntOrderOffset = 0 )
inline

Construct a new BoundaryHyperbolicFlowIntegrator object.

Parameters
[in]fluxflux function
[in]uvector state coefficient
[in]alphaɑ coefficient (β = ɑ/2)
[in]IntOrderOffsetintegration order offset

Definition at line 542 of file hyperbolic.hpp.

◆ BoundaryHyperbolicFlowIntegrator() [2/2]

mfem::BoundaryHyperbolicFlowIntegrator::BoundaryHyperbolicFlowIntegrator ( const FluxFunction & flux,
VectorCoefficient & u,
real_t alpha,
real_t beta,
int IntOrderOffset = 0 )

Construct a new BoundaryHyperbolicFlowIntegrator object.

Parameters
[in]fluxflux function
[in]uvector state coefficient
[in]alphaɑ coefficient
[in]betaβ coefficient
[in]IntOrderOffsetintegration order offset

Definition at line 605 of file hyperbolic.cpp.

Member Function Documentation

◆ AssembleRHSElementVect() [1/3]

void mfem::BoundaryHyperbolicFlowIntegrator::AssembleRHSElementVect ( const FiniteElement & el,
ElementTransformation & Tr,
Vector & elvect )
overridevirtual
Warning
Boundary element integration not implemented, use AssembleRHSElementVect(const FiniteElement&, FaceElementTransformations &, Vector &) instead

Implements mfem::LinearFormIntegrator.

Definition at line 621 of file hyperbolic.cpp.

◆ AssembleRHSElementVect() [2/3]

void mfem::BoundaryHyperbolicFlowIntegrator::AssembleRHSElementVect ( const FiniteElement & el,
FaceElementTransformations & Tr,
Vector & elvect )
overridevirtual

Implements <-F(u,x) n, v> with abstract F computed by FluxFunction::ComputeFluxDotN() of the flux function object.

Parameters
[in]elfinite element
[in]Trface element transformations
[out]elvectevaluated dual vector <F(u,x) n, v>

Reimplemented from mfem::LinearFormIntegrator.

Definition at line 630 of file hyperbolic.cpp.

◆ AssembleRHSElementVect() [3/3]

void mfem::LinearFormIntegrator::AssembleRHSElementVect ( const FiniteElement & el1,
const FiniteElement & el2,
FaceElementTransformations & Tr,
Vector & elvect )
virtual

Reimplemented from mfem::LinearFormIntegrator.

Definition at line 51 of file lininteg.cpp.

◆ GetFluxFunction()

const FluxFunction & mfem::BoundaryHyperbolicFlowIntegrator::GetFluxFunction ( ) const
inline

Get the associated flux function.

Definition at line 565 of file hyperbolic.hpp.

◆ GetMaxCharSpeed()

real_t mfem::BoundaryHyperbolicFlowIntegrator::GetMaxCharSpeed ( ) const
inline

Get the maximum characteristic speed.

Definition at line 562 of file hyperbolic.hpp.

◆ ResetMaxCharSpeed()

void mfem::BoundaryHyperbolicFlowIntegrator::ResetMaxCharSpeed ( )
inline

Reset the maximum characteristic speed to zero.

Definition at line 559 of file hyperbolic.hpp.


The documentation for this class was generated from the following files: