MFEM  v4.2.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
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, 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)
 
- 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 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 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 AssembleFaceMatrix (const FiniteElement &trial_face_fe, const FiniteElement &test_fe1, const FiniteElement &test_fe2, 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)
 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 ()
 
- Public Member Functions inherited from mfem::NonlinearFormIntegrator
void SetIntRule (const IntegrationRule *ir)
 Prescribe a fixed IntegrationRule to use (when ir != NULL) or let the integrator choose (when ir == NULL). More...
 
void SetIntegrationRule (const IntegrationRule &irule)
 Prescribe a fixed IntegrationRule to use. More...
 
virtual double GetElementEnergy (const FiniteElement &el, ElementTransformation &Tr, const Vector &elfun)
 Compute the local energy. More...
 
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
 

Additional Inherited Members

- 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. 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.

Definition at line 2662 of file bilininteg.hpp.

Constructor & Destructor Documentation

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

Construct integrator with rho = 1.

Definition at line 2679 of file bilininteg.hpp.

mfem::DGTraceIntegrator::DGTraceIntegrator ( Coefficient _rho,
VectorCoefficient _u,
double  a,
double  b 
)
inline

Definition at line 2682 of file bilininteg.hpp.

Member Function Documentation

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.

Definition at line 1133 of file bilininteg_dgtrace_pa.cpp.

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.

Definition at line 1140 of file bilininteg_dgtrace_pa.cpp.

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

Reimplemented from mfem::BilinearFormIntegrator.

Definition at line 427 of file bilininteg_dgtrace_ea.cpp.

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

Reimplemented from mfem::BilinearFormIntegrator.

Definition at line 348 of file bilininteg_dgtrace_ea.cpp.

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

Reimplemented from mfem::BilinearFormIntegrator.

Definition at line 2597 of file bilininteg.cpp.

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

Reimplemented from mfem::BilinearFormIntegrator.

Definition at line 298 of file bilininteg_dgtrace_pa.cpp.

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

Reimplemented from mfem::BilinearFormIntegrator.

Definition at line 293 of file bilininteg_dgtrace_pa.cpp.

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

Definition at line 2731 of file bilininteg.cpp.

Member Data Documentation

double mfem::DGTraceIntegrator::alpha
protected

Definition at line 2667 of file bilininteg.hpp.

double mfem::DGTraceIntegrator::beta
protected

Definition at line 2667 of file bilininteg.hpp.

int mfem::DGTraceIntegrator::dim
protected

Definition at line 2672 of file bilininteg.hpp.

int mfem::DGTraceIntegrator::dofs1D
protected

Definition at line 2672 of file bilininteg.hpp.

const FaceGeometricFactors* mfem::DGTraceIntegrator::geom
protected

Not owned.

Definition at line 2671 of file bilininteg.hpp.

const DofToQuad* mfem::DGTraceIntegrator::maps
protected

Not owned.

Definition at line 2670 of file bilininteg.hpp.

int mfem::DGTraceIntegrator::nf
protected

Definition at line 2672 of file bilininteg.hpp.

int mfem::DGTraceIntegrator::nq
protected

Definition at line 2672 of file bilininteg.hpp.

Vector mfem::DGTraceIntegrator::pa_data
protected

Definition at line 2669 of file bilininteg.hpp.

int mfem::DGTraceIntegrator::quad1D
protected

Definition at line 2672 of file bilininteg.hpp.

Coefficient* mfem::DGTraceIntegrator::rho
protected

Definition at line 2665 of file bilininteg.hpp.

VectorCoefficient* mfem::DGTraceIntegrator::u
protected

Definition at line 2666 of file bilininteg.hpp.


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