MFEM  v4.6.0
Finite element discretization library
Public Member Functions | Static Public Member Functions | Protected Attributes | List of all members
mfem::DGTraceIntegrator Class Reference

#include <bilininteg.hpp>

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

Public Member Functions

 DGTraceIntegrator (VectorCoefficient &u_, double a)
 Construct integrator with rho = 1, b = 0.5*a. More...
 
 DGTraceIntegrator (VectorCoefficient &u_, double a, double b)
 Construct integrator with rho = 1. More...
 
 DGTraceIntegrator (Coefficient &rho_, VectorCoefficient &u_, double a, double b)
 
virtual void AssembleFaceMatrix (const FiniteElement &el1, const FiniteElement &el2, FaceElementTransformations &Trans, DenseMatrix &elmat)
 
virtual void AssemblePAInteriorFaces (const FiniteElementSpace &fes)
 
virtual void AssemblePABoundaryFaces (const FiniteElementSpace &fes)
 
virtual void AddMultTransposePA (const Vector &x, Vector &y) const
 Method for partially assembled transposed action. More...
 
virtual void AddMultPA (const Vector &, Vector &) const
 Method for partially assembled action. More...
 
virtual void AssembleEAInteriorFaces (const FiniteElementSpace &fes, Vector &ea_data_int, Vector &ea_data_ext, const bool add)
 
virtual void AssembleEABoundaryFaces (const FiniteElementSpace &fes, Vector &ea_data_bdr, const bool add)
 
virtual void AssembleFaceMatrix (const FiniteElement &el1, const FiniteElement &el2, FaceElementTransformations &Trans, DenseMatrix &elmat)
 
virtual void AssembleFaceMatrix (const FiniteElement &trial_face_fe, const FiniteElement &test_fe1, const FiniteElement &test_fe2, FaceElementTransformations &Trans, DenseMatrix &elmat)
 
virtual void AssemblePA (const FiniteElementSpace &fes)
 Method defining partial assembly. More...
 
virtual void AssemblePA (const FiniteElementSpace &trial_fes, const FiniteElementSpace &test_fes)
 
virtual void AssemblePA (const FiniteElementSpace &fes)
 Method defining partial assembly. More...
 
virtual void AssemblePA (const FiniteElementSpace &trial_fes, const FiniteElementSpace &test_fes)
 
- Public Member Functions inherited from mfem::BilinearFormIntegrator
virtual void AssemblePA (const FiniteElementSpace &fes)
 Method defining partial assembly. More...
 
virtual void AssemblePA (const FiniteElementSpace &trial_fes, const FiniteElementSpace &test_fes)
 
virtual void AssembleNURBSPA (const FiniteElementSpace &fes)
 Method defining partial assembly on NURBS patches. More...
 
virtual void AssemblePABoundary (const FiniteElementSpace &fes)
 
virtual void AssembleDiagonalPA (Vector &diag)
 Assemble diagonal and add it to Vector diag. More...
 
virtual void AssembleDiagonalPA_ADAt (const Vector &D, Vector &diag)
 Assemble diagonal of ADA^T (A is this integrator) and add it to diag. More...
 
virtual void AddMultNURBSPA (const Vector &x, Vector &y) const
 Method for partially assembled action on NURBS patches. More...
 
virtual void AssembleEA (const FiniteElementSpace &fes, Vector &emat, const bool add=true)
 Method defining element assembly. More...
 
virtual void AssembleMF (const FiniteElementSpace &fes)
 Method defining matrix-free assembly. More...
 
virtual void AddMultMF (const Vector &x, Vector &y) const
 
virtual void AddMultTransposeMF (const Vector &x, Vector &y) const
 
virtual void AssembleDiagonalMF (Vector &diag)
 Assemble diagonal and add it to Vector diag. More...
 
virtual void AssembleElementMatrix (const FiniteElement &el, ElementTransformation &Trans, DenseMatrix &elmat)
 Given a particular Finite Element computes the element matrix elmat. More...
 
virtual void AssembleElementMatrix2 (const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &elmat)
 
virtual void AssemblePatchMatrix (const int patch, const FiniteElementSpace &fes, SparseMatrix *&smat)
 
virtual void AssembleFaceMatrix (const FiniteElement &trial_face_fe, const FiniteElement &test_fe1, const FiniteElement &test_fe2, FaceElementTransformations &Trans, DenseMatrix &elmat)
 
virtual void AssembleTraceFaceMatrix (int elem, const FiniteElement &trial_face_fe, const FiniteElement &test_fe, FaceElementTransformations &Trans, DenseMatrix &elmat)
 
virtual void AssembleElementVector (const FiniteElement &el, ElementTransformation &Tr, const Vector &elfun, Vector &elvect)
 Perform the local action of the BilinearFormIntegrator. Note that the default implementation in the base class is general but not efficient. More...
 
virtual void AssembleFaceVector (const FiniteElement &el1, const FiniteElement &el2, FaceElementTransformations &Tr, const Vector &elfun, Vector &elvect)
 Perform the local action of the BilinearFormIntegrator resulting from a face integral term. Note that the default implementation in the base class is general but not efficient. 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 void ComputeElementFlux (const FiniteElement &el, ElementTransformation &Trans, Vector &u, const FiniteElement &fluxelem, Vector &flux, bool with_coef=true, const IntegrationRule *ir=NULL)
 Virtual method required for Zienkiewicz-Zhu type error estimators. More...
 
virtual double ComputeFluxEnergy (const FiniteElement &fluxelem, ElementTransformation &Trans, Vector &flux, Vector *d_energy=NULL)
 Virtual method required for Zienkiewicz-Zhu type error estimators. More...
 
virtual ~BilinearFormIntegrator ()
 
virtual void AssemblePA (const FiniteElementSpace &fes)
 Method defining partial assembly. More...
 
virtual void AssemblePA (const FiniteElementSpace &trial_fes, const FiniteElementSpace &test_fes)
 
- Public Member Functions inherited from mfem::NonlinearFormIntegrator
virtual void SetIntRule (const IntegrationRule *ir)
 Prescribe a fixed IntegrationRule to use (when ir != NULL) or let the integrator choose (when ir == NULL). More...
 
void SetIntegrationMode (Mode m)
 
void SetNURBSPatchIntRule (NURBSMeshRules *pr)
 For patchwise integration, SetNURBSPatchIntRule must be called. More...
 
bool HasNURBSPatchIntRule () const
 
bool Patchwise () const
 
void SetIntegrationRule (const IntegrationRule &ir)
 Prescribe a fixed IntegrationRule to use. More...
 
void SetPAMemoryType (MemoryType mt)
 
const IntegrationRuleGetIntegrationRule () const
 Get the integration rule of the integrator (possibly NULL). More...
 
virtual double GetElementEnergy (const FiniteElement &el, ElementTransformation &Tr, const Vector &elfun)
 Compute the local energy. More...
 
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. More...
 
virtual double GetLocalStateEnergyPA (const Vector &x) const
 Compute the local (to the MPI rank) energy with partial assembly. More...
 
virtual void AddMultGradPA (const Vector &x, Vector &y) const
 Method for partially assembled gradient action. More...
 
virtual void AssembleGradDiagonalPA (Vector &diag) const
 Method for computing the diagonal of the gradient with partial assembly. More...
 
virtual bool SupportsCeed () const
 Indicates whether this integrator can use a Ceed backend. More...
 
ceed::OperatorGetCeedOp ()
 
virtual ~NonlinearFormIntegrator ()
 

Static Public Member Functions

static const IntegrationRuleGetRule (Geometry::Type geom, int order, FaceElementTransformations &T)
 

Protected Attributes

Coefficientrho
 
VectorCoefficientu
 
double alpha
 
double beta
 
Vector pa_data
 
const DofToQuadmaps
 Not owned. More...
 
const FaceGeometricFactorsgeom
 Not owned. More...
 
int dim
 
int nf
 
int nq
 
int dofs1D
 
int quad1D
 
- Protected Attributes inherited from mfem::NonlinearFormIntegrator
const IntegrationRuleIntRule
 
Mode integrationMode = Mode::ELEMENTWISE
 
NURBSMeshRulespatchRules = nullptr
 
ceed::OperatorceedOp
 
MemoryType pa_mt = MemoryType::DEFAULT
 

Additional Inherited Members

- Public Types inherited from mfem::NonlinearFormIntegrator
enum  Mode { ELEMENTWISE = 0, PATCHWISE = 1, PATCHWISE_REDUCED = 2 }
 
- Protected Member Functions inherited from mfem::BilinearFormIntegrator
 BilinearFormIntegrator (const IntegrationRule *ir=NULL)
 
- Protected Member Functions inherited from mfem::NonlinearFormIntegrator
 NonlinearFormIntegrator (const IntegrationRule *ir=NULL)
 

Detailed Description

Integrator for the DG form: alpha < rho_u (u.n) {v},[w] > + beta < rho_u |u.n| [v],[w] >, where v and w are the trial and test variables, respectively, and rho/u are given scalar/vector coefficients. {v} represents the average value of v on the face and [v] is the jump such that {v}=(v1+v2)/2 and [v]=(v1-v2) for the face between elements 1 and 2. For boundary elements, v2=0. The vector coefficient, u, is assumed to be continuous across the faces and when given the scalar coefficient, rho, is assumed to be discontinuous. The integrator uses the upwind value of rho, rho_u, which is value from the side into which the vector coefficient, u, points.

One use case for this integrator is to discretize the operator -u.grad(v) with a DG formulation. The resulting formulation uses the ConvectionIntegrator (with coefficient u, and parameter alpha = -1) and the transpose of the DGTraceIntegrator (with coefficient u, and parameters alpha = 1, beta = -1/2 to use the upwind face flux, see also NonconservativeDGTraceIntegrator). This discretization and the handling of the inflow and outflow boundaries is illustrated in Example 9/9p.

Another use case for this integrator is to discretize the operator -div(u v) with a DG formulation. The resulting formulation is conservative and consists of the ConservativeConvectionIntegrator (with coefficient u, and parameter alpha = -1) plus the DGTraceIntegrator (with coefficient u, and parameters alpha = -1, beta = -1/2 to use the upwind face flux).

Definition at line 3065 of file bilininteg.hpp.

Constructor & Destructor Documentation

◆ DGTraceIntegrator() [1/3]

mfem::DGTraceIntegrator::DGTraceIntegrator ( VectorCoefficient u_,
double  a 
)
inline

Construct integrator with rho = 1, b = 0.5*a.

Definition at line 3082 of file bilininteg.hpp.

◆ DGTraceIntegrator() [2/3]

mfem::DGTraceIntegrator::DGTraceIntegrator ( VectorCoefficient u_,
double  a,
double  b 
)
inline

Construct integrator with rho = 1.

Definition at line 3086 of file bilininteg.hpp.

◆ DGTraceIntegrator() [3/3]

mfem::DGTraceIntegrator::DGTraceIntegrator ( Coefficient rho_,
VectorCoefficient u_,
double  a,
double  b 
)
inline

Definition at line 3089 of file bilininteg.hpp.

Member Function Documentation

◆ AddMultPA()

virtual void mfem::DGTraceIntegrator::AddMultPA ( const Vector x,
Vector y 
) const
virtual

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 from mfem::BilinearFormIntegrator.

◆ AddMultTransposePA()

virtual void mfem::DGTraceIntegrator::AddMultTransposePA ( const Vector x,
Vector y 
) const
virtual

Method for partially assembled transposed action.

Perform the transpose 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 from mfem::BilinearFormIntegrator.

◆ AssembleEABoundaryFaces()

virtual void mfem::DGTraceIntegrator::AssembleEABoundaryFaces ( const FiniteElementSpace fes,
Vector ea_data_bdr,
const bool  add 
)
virtual

Reimplemented from mfem::BilinearFormIntegrator.

◆ AssembleEAInteriorFaces()

virtual void mfem::DGTraceIntegrator::AssembleEAInteriorFaces ( const FiniteElementSpace fes,
Vector ea_data_int,
Vector ea_data_ext,
const bool  add 
)
virtual

Reimplemented from mfem::BilinearFormIntegrator.

◆ AssembleFaceMatrix() [1/3]

void mfem::BilinearFormIntegrator::AssembleFaceMatrix

Definition at line 164 of file bilininteg.cpp.

◆ AssembleFaceMatrix() [2/3]

void mfem::BilinearFormIntegrator::AssembleFaceMatrix

Abstract method used for assembling TraceFaceIntegrators in a MixedBilinearForm.

Definition at line 172 of file bilininteg.cpp.

◆ AssembleFaceMatrix() [3/3]

void mfem::DGTraceIntegrator::AssembleFaceMatrix ( const FiniteElement el1,
const FiniteElement el2,
FaceElementTransformations Trans,
DenseMatrix elmat 
)
virtual

Reimplemented from mfem::BilinearFormIntegrator.

Definition at line 3273 of file bilininteg.cpp.

◆ AssemblePA() [1/4]

void mfem::NonlinearFormIntegrator::AssemblePA

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.

Definition at line 31 of file nonlininteg.cpp.

◆ AssemblePA() [2/4]

void mfem::NonlinearFormIntegrator::AssemblePA

Method defining partial assembly.

The result of the partial assembly is stored internally so that it can be used later in the methods AddMultPA().

Definition at line 25 of file nonlininteg.cpp.

◆ AssemblePA() [3/4]

void mfem::BilinearFormIntegrator::AssemblePA

Method defining partial assembly.

The result of the partial assembly is stored internally so that it can be used later in the methods AddMultPA() and AddMultTransposePA().

Definition at line 23 of file bilininteg.cpp.

◆ AssemblePA() [4/4]

void mfem::BilinearFormIntegrator::AssemblePA

Used with BilinearFormIntegrators that have different spaces.

Definition at line 35 of file bilininteg.cpp.

◆ AssemblePABoundaryFaces()

virtual void mfem::DGTraceIntegrator::AssemblePABoundaryFaces ( const FiniteElementSpace fes)
virtual

Reimplemented from mfem::BilinearFormIntegrator.

◆ AssemblePAInteriorFaces()

virtual void mfem::DGTraceIntegrator::AssemblePAInteriorFaces ( const FiniteElementSpace fes)
virtual

Reimplemented from mfem::BilinearFormIntegrator.

◆ GetRule()

const IntegrationRule & mfem::DGTraceIntegrator::GetRule ( Geometry::Type  geom,
int  order,
FaceElementTransformations T 
)
static

Definition at line 3407 of file bilininteg.cpp.

Member Data Documentation

◆ alpha

double mfem::DGTraceIntegrator::alpha
protected

Definition at line 3070 of file bilininteg.hpp.

◆ beta

double mfem::DGTraceIntegrator::beta
protected

Definition at line 3070 of file bilininteg.hpp.

◆ dim

int mfem::DGTraceIntegrator::dim
protected

Definition at line 3075 of file bilininteg.hpp.

◆ dofs1D

int mfem::DGTraceIntegrator::dofs1D
protected

Definition at line 3075 of file bilininteg.hpp.

◆ geom

const FaceGeometricFactors* mfem::DGTraceIntegrator::geom
protected

Not owned.

Definition at line 3074 of file bilininteg.hpp.

◆ maps

const DofToQuad* mfem::DGTraceIntegrator::maps
protected

Not owned.

Definition at line 3073 of file bilininteg.hpp.

◆ nf

int mfem::DGTraceIntegrator::nf
protected

Definition at line 3075 of file bilininteg.hpp.

◆ nq

int mfem::DGTraceIntegrator::nq
protected

Definition at line 3075 of file bilininteg.hpp.

◆ pa_data

Vector mfem::DGTraceIntegrator::pa_data
protected

Definition at line 3072 of file bilininteg.hpp.

◆ quad1D

int mfem::DGTraceIntegrator::quad1D
protected

Definition at line 3075 of file bilininteg.hpp.

◆ rho

Coefficient* mfem::DGTraceIntegrator::rho
protected

Definition at line 3068 of file bilininteg.hpp.

◆ u

VectorCoefficient* mfem::DGTraceIntegrator::u
protected

Definition at line 3069 of file bilininteg.hpp.


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