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::MassIntegrator Class Reference

#include <bilininteg.hpp>

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

Public Member Functions

 MassIntegrator (const IntegrationRule *ir=NULL)
 
 MassIntegrator (Coefficient &q, const IntegrationRule *ir=NULL)
 Construct a mass integrator with coefficient q. More...
 
virtual ~MassIntegrator ()
 
virtual void AssembleElementMatrix (const FiniteElement &el, ElementTransformation &Trans, DenseMatrix &elmat)
 
virtual void AssembleElementMatrix2 (const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &elmat)
 
virtual void AssembleMF (const FiniteElementSpace &fes)
 Method defining matrix-free assembly. More...
 
virtual void AssemblePA (const FiniteElementSpace &fes)
 Method defining partial assembly. More...
 
virtual void AssembleEA (const FiniteElementSpace &fes, Vector &emat, const bool add)
 Method defining element assembly. More...
 
virtual void AssembleDiagonalPA (Vector &diag)
 Assemble diagonal and add it to Vector diag. More...
 
virtual void AssembleDiagonalMF (Vector &diag)
 Assemble diagonal and add it to Vector diag. More...
 
virtual void AddMultMF (const Vector &, Vector &) const
 
virtual void AddMultPA (const Vector &, Vector &) const
 Method for partially assembled action. More...
 
void SetupPA (const FiniteElementSpace &fes)
 
- Public Member Functions inherited from mfem::BilinearFormIntegrator
virtual void AssemblePA (const FiniteElementSpace &trial_fes, const FiniteElementSpace &test_fes)
 
virtual void AssemblePAInteriorFaces (const FiniteElementSpace &fes)
 
virtual void AssemblePABoundaryFaces (const FiniteElementSpace &fes)
 
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 AddMultTransposePA (const Vector &x, Vector &y) const
 Method for partially assembled transposed action. More...
 
virtual void AddMultTransposeMF (const Vector &x, Vector &y) const
 
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 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 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 (const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans)
 

Protected Attributes

Vector shape
 
Vector te_shape
 
CoefficientQ
 
const FiniteElementSpacefespace
 
Vector pa_data
 
const DofToQuadmaps
 Not owned. More...
 
const GeometricFactorsgeom
 Not owned. More...
 
int dim
 
int ne
 
int nq
 
int dofs1D
 
int quad1D
 
CeedData * ceedDataPtr
 
- 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

Class for local mass matrix assembling a(u,v) := (Q u, v)

Definition at line 1992 of file bilininteg.hpp.

Constructor & Destructor Documentation

mfem::MassIntegrator::MassIntegrator ( const IntegrationRule ir = NULL)
inline

Definition at line 2010 of file bilininteg.hpp.

mfem::MassIntegrator::MassIntegrator ( Coefficient q,
const IntegrationRule ir = NULL 
)
inline

Construct a mass integrator with coefficient q.

Definition at line 2015 of file bilininteg.hpp.

virtual mfem::MassIntegrator::~MassIntegrator ( )
inlinevirtual

Definition at line 2019 of file bilininteg.hpp.

Member Function Documentation

void mfem::MassIntegrator::AddMultMF ( const Vector x,
Vector y 
) const
virtual

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 AssembleMF() has been called.

Reimplemented from mfem::BilinearFormIntegrator.

Definition at line 43 of file bilininteg_mass_mf.cpp.

void mfem::MassIntegrator::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 1217 of file bilininteg_mass_pa.cpp.

void mfem::MassIntegrator::AssembleDiagonalMF ( Vector diag)
virtual

Assemble diagonal and add it to Vector diag.

Reimplemented from mfem::BilinearFormIntegrator.

Definition at line 57 of file bilininteg_mass_mf.cpp.

void mfem::MassIntegrator::AssembleDiagonalPA ( Vector diag)
virtual

Assemble diagonal and add it to Vector diag.

Reimplemented from mfem::BilinearFormIntegrator.

Definition at line 460 of file bilininteg_mass_pa.cpp.

void mfem::MassIntegrator::AssembleEA ( const FiniteElementSpace fes,
Vector emat,
const bool  add 
)
virtual

Method defining element assembly.

The result of the element assembly is added to the emat Vector if add is true. Otherwise, if add is false, we set emat.

Reimplemented from mfem::BilinearFormIntegrator.

Definition at line 226 of file bilininteg_mass_ea.cpp.

void mfem::MassIntegrator::AssembleElementMatrix ( const FiniteElement el,
ElementTransformation Trans,
DenseMatrix elmat 
)
virtual

Given a particular Finite Element computes the element mass matrix elmat.

Reimplemented from mfem::BilinearFormIntegrator.

Definition at line 891 of file bilininteg.cpp.

void mfem::MassIntegrator::AssembleElementMatrix2 ( const FiniteElement trial_fe,
const FiniteElement test_fe,
ElementTransformation Trans,
DenseMatrix elmat 
)
virtual

Compute the local matrix representation of a bilinear form a(u,v) defined on different trial (given by u) and test (given by v) spaces. The rows in the local matrix correspond to the test dofs and the columns – to the trial dofs.

Reimplemented from mfem::BilinearFormIntegrator.

Definition at line 923 of file bilininteg.cpp.

void mfem::MassIntegrator::AssembleMF ( const FiniteElementSpace fes)
virtual

Method defining matrix-free assembly.

Used with BilinearFormIntegrators that have different spaces.The result of fully matrix-free assembly is stored internally so that it can be used later in the methods AddMultMF() and AddMultTransposeMF().

Reimplemented from mfem::BilinearFormIntegrator.

Definition at line 22 of file bilininteg_mass_mf.cpp.

void mfem::MassIntegrator::AssemblePA ( const FiniteElementSpace fes)
virtual

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

Reimplemented from mfem::BilinearFormIntegrator.

Definition at line 26 of file bilininteg_mass_pa.cpp.

const IntegrationRule & mfem::MassIntegrator::GetRule ( const FiniteElement trial_fe,
const FiniteElement test_fe,
ElementTransformation Trans 
)
static

Definition at line 960 of file bilininteg.cpp.

void mfem::MassIntegrator::SetupPA ( const FiniteElementSpace fes)

Member Data Documentation

CeedData* mfem::MassIntegrator::ceedDataPtr
protected

Definition at line 2007 of file bilininteg.hpp.

int mfem::MassIntegrator::dim
protected

Definition at line 2004 of file bilininteg.hpp.

int mfem::MassIntegrator::dofs1D
protected

Definition at line 2004 of file bilininteg.hpp.

const FiniteElementSpace* mfem::MassIntegrator::fespace
protected

Definition at line 2000 of file bilininteg.hpp.

const GeometricFactors* mfem::MassIntegrator::geom
protected

Not owned.

Definition at line 2003 of file bilininteg.hpp.

const DofToQuad* mfem::MassIntegrator::maps
protected

Not owned.

Definition at line 2002 of file bilininteg.hpp.

int mfem::MassIntegrator::ne
protected

Definition at line 2004 of file bilininteg.hpp.

int mfem::MassIntegrator::nq
protected

Definition at line 2004 of file bilininteg.hpp.

Vector mfem::MassIntegrator::pa_data
protected

Definition at line 2001 of file bilininteg.hpp.

Coefficient* mfem::MassIntegrator::Q
protected

Definition at line 1998 of file bilininteg.hpp.

int mfem::MassIntegrator::quad1D
protected

Definition at line 2004 of file bilininteg.hpp.

Vector mfem::MassIntegrator::shape
protected

Definition at line 1996 of file bilininteg.hpp.

Vector mfem::MassIntegrator::te_shape
protected

Definition at line 1996 of file bilininteg.hpp.


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