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

#include <bilininteg.hpp>

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

Public Member Functions

 DGDiffusionIntegrator (const double s, const double k)
 
 DGDiffusionIntegrator (Coefficient &q, const double s, const double k)
 
 DGDiffusionIntegrator (MatrixCoefficient &q, const double s, const double k)
 
virtual void AssembleFaceMatrix (const FiniteElement &el1, const FiniteElement &el2, FaceElementTransformations &Trans, DenseMatrix &elmat)
 
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)
 
- 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 AssemblePAInteriorFaces (const FiniteElementSpace &fes)
 
virtual void AssemblePABoundaryFaces (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 AddMultPA (const Vector &x, Vector &y) const
 Method for partially assembled action. More...
 
virtual void AddMultNURBSPA (const Vector &x, Vector &y) const
 Method for partially assembled action on NURBS patches. More...
 
virtual void AddMultTransposePA (const Vector &x, Vector &y) const
 Method for partially assembled transposed action. 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 AssembleEAInteriorFaces (const FiniteElementSpace &fes, Vector &ea_data_int, Vector &ea_data_ext, const bool add=true)
 
virtual void AssembleEABoundaryFaces (const FiniteElementSpace &fes, Vector &ea_data_bdr, const bool add=true)
 
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 ()
 

Protected Attributes

CoefficientQ
 
MatrixCoefficientMQ
 
double sigma
 
double kappa
 
Vector shape1
 
Vector shape2
 
Vector dshape1dn
 
Vector dshape2dn
 
Vector nor
 
Vector nh
 
Vector ni
 
DenseMatrix jmat
 
DenseMatrix dshape1
 
DenseMatrix dshape2
 
DenseMatrix mq
 
DenseMatrix adjJ
 
- 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:

- < {(Q grad(u)).n}, [v] > + sigma < [u], {(Q grad(v)).n} >
+ kappa < {h^{-1} Q} [u], [v] >

where Q is a scalar or matrix diffusion coefficient and u, v are the trial and test spaces, respectively. The parameters sigma and kappa determine the DG method to be used (when this integrator is added to the "broken" DiffusionIntegrator): sigma = -1, kappa >= kappa0: symm. interior penalty (IP or SIPG) method, sigma = +1, kappa > 0: non-symmetric interior penalty (NIPG) method, sigma = +1, kappa = 0: the method of Baumann and Oden.

Definition at line 3161 of file bilininteg.hpp.

Constructor & Destructor Documentation

◆ DGDiffusionIntegrator() [1/3]

mfem::DGDiffusionIntegrator::DGDiffusionIntegrator ( const double  s,
const double  k 
)
inline

Definition at line 3173 of file bilininteg.hpp.

◆ DGDiffusionIntegrator() [2/3]

mfem::DGDiffusionIntegrator::DGDiffusionIntegrator ( Coefficient q,
const double  s,
const double  k 
)
inline

Definition at line 3175 of file bilininteg.hpp.

◆ DGDiffusionIntegrator() [3/3]

mfem::DGDiffusionIntegrator::DGDiffusionIntegrator ( MatrixCoefficient q,
const double  s,
const double  k 
)
inline

Definition at line 3177 of file bilininteg.hpp.

Member Function Documentation

◆ 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::DGDiffusionIntegrator::AssembleFaceMatrix ( const FiniteElement el1,
const FiniteElement el2,
FaceElementTransformations Trans,
DenseMatrix elmat 
)
virtual

Reimplemented from mfem::BilinearFormIntegrator.

Definition at line 3414 of file bilininteg.cpp.

Member Data Documentation

◆ adjJ

DenseMatrix mfem::DGDiffusionIntegrator::adjJ
protected

Definition at line 3170 of file bilininteg.hpp.

◆ dshape1

DenseMatrix mfem::DGDiffusionIntegrator::dshape1
protected

Definition at line 3170 of file bilininteg.hpp.

◆ dshape1dn

Vector mfem::DGDiffusionIntegrator::dshape1dn
protected

Definition at line 3169 of file bilininteg.hpp.

◆ dshape2

DenseMatrix mfem::DGDiffusionIntegrator::dshape2
protected

Definition at line 3170 of file bilininteg.hpp.

◆ dshape2dn

Vector mfem::DGDiffusionIntegrator::dshape2dn
protected

Definition at line 3169 of file bilininteg.hpp.

◆ jmat

DenseMatrix mfem::DGDiffusionIntegrator::jmat
protected

Definition at line 3170 of file bilininteg.hpp.

◆ kappa

double mfem::DGDiffusionIntegrator::kappa
protected

Definition at line 3166 of file bilininteg.hpp.

◆ MQ

MatrixCoefficient* mfem::DGDiffusionIntegrator::MQ
protected

Definition at line 3165 of file bilininteg.hpp.

◆ mq

DenseMatrix mfem::DGDiffusionIntegrator::mq
protected

Definition at line 3170 of file bilininteg.hpp.

◆ nh

Vector mfem::DGDiffusionIntegrator::nh
protected

Definition at line 3169 of file bilininteg.hpp.

◆ ni

Vector mfem::DGDiffusionIntegrator::ni
protected

Definition at line 3169 of file bilininteg.hpp.

◆ nor

Vector mfem::DGDiffusionIntegrator::nor
protected

Definition at line 3169 of file bilininteg.hpp.

◆ Q

Coefficient* mfem::DGDiffusionIntegrator::Q
protected

Definition at line 3164 of file bilininteg.hpp.

◆ shape1

Vector mfem::DGDiffusionIntegrator::shape1
protected

Definition at line 3169 of file bilininteg.hpp.

◆ shape2

Vector mfem::DGDiffusionIntegrator::shape2
protected

Definition at line 3169 of file bilininteg.hpp.

◆ sigma

double mfem::DGDiffusionIntegrator::sigma
protected

Definition at line 3166 of file bilininteg.hpp.


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