MFEM  v4.2.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
Public Member Functions | Protected Attributes | List of all members
mfem::VectorFEMassIntegrator Class Reference

#include <bilininteg.hpp>

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

Public Member Functions

 VectorFEMassIntegrator ()
 
 VectorFEMassIntegrator (Coefficient *_q)
 
 VectorFEMassIntegrator (Coefficient &q)
 
 VectorFEMassIntegrator (VectorCoefficient *_vq)
 
 VectorFEMassIntegrator (VectorCoefficient &vq)
 
 VectorFEMassIntegrator (MatrixCoefficient *_mq)
 
 VectorFEMassIntegrator (MatrixCoefficient &mq)
 
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 AssemblePA (const FiniteElementSpace &fes)
 Method defining partial assembly. More...
 
virtual void AssemblePA (const FiniteElementSpace &trial_fes, const FiniteElementSpace &test_fes)
 
virtual void AddMultPA (const Vector &x, Vector &y) const
 Method for partially assembled action. More...
 
virtual void AssembleDiagonalPA (Vector &diag)
 Assemble diagonal and add it to Vector diag. More...
 
- Public Member Functions inherited from mfem::BilinearFormIntegrator
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 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 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 ()
 

Protected Attributes

CoefficientQ
 
VectorCoefficientVQ
 
MatrixCoefficientMQ
 
Vector pa_data
 
const DofToQuadmapsO
 Not owned. DOF-to-quad map, open. More...
 
const DofToQuadmapsC
 Not owned. DOF-to-quad map, closed. More...
 
const DofToQuadmapsOtest
 Not owned. DOF-to-quad map, open. More...
 
const DofToQuadmapsCtest
 Not owned. DOF-to-quad map, closed. More...
 
const GeometricFactorsgeom
 Not owned. More...
 
int dim
 
int ne
 
int nq
 
int dofs1D
 
int dofs1Dtest
 
int quad1D
 
int trial_fetype
 
int test_fetype
 
bool symmetric = true
 False if using a nonsymmetric matrix coefficient. More...
 
- 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 (Q u, v), where Q is an optional coefficient (of type scalar, vector (diagonal matrix), or matrix), trial function u is in H(Curl) or H(Div), and test function v is in H(Curl), H(Div), or v=(v1,...,vn), where vi are in H1.

Definition at line 2411 of file bilininteg.hpp.

Constructor & Destructor Documentation

mfem::VectorFEMassIntegrator::VectorFEMassIntegrator ( )
inline

Definition at line 2442 of file bilininteg.hpp.

mfem::VectorFEMassIntegrator::VectorFEMassIntegrator ( Coefficient _q)
inline

Definition at line 2443 of file bilininteg.hpp.

mfem::VectorFEMassIntegrator::VectorFEMassIntegrator ( Coefficient q)
inline

Definition at line 2444 of file bilininteg.hpp.

mfem::VectorFEMassIntegrator::VectorFEMassIntegrator ( VectorCoefficient _vq)
inline

Definition at line 2445 of file bilininteg.hpp.

mfem::VectorFEMassIntegrator::VectorFEMassIntegrator ( VectorCoefficient vq)
inline

Definition at line 2446 of file bilininteg.hpp.

mfem::VectorFEMassIntegrator::VectorFEMassIntegrator ( MatrixCoefficient _mq)
inline

Definition at line 2447 of file bilininteg.hpp.

mfem::VectorFEMassIntegrator::VectorFEMassIntegrator ( MatrixCoefficient mq)
inline

Definition at line 2448 of file bilininteg.hpp.

Member Function Documentation

void mfem::VectorFEMassIntegrator::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 962 of file bilininteg_vectorfe.cpp.

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

Assemble diagonal and add it to Vector diag.

Reimplemented from mfem::BilinearFormIntegrator.

Definition at line 900 of file bilininteg_vectorfe.cpp.

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

Given a particular Finite Element computes the element matrix elmat.

Reimplemented from mfem::BilinearFormIntegrator.

Definition at line 1905 of file bilininteg.cpp.

void mfem::VectorFEMassIntegrator::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 1970 of file bilininteg.cpp.

void mfem::VectorFEMassIntegrator::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 714 of file bilininteg_vectorfe.cpp.

void mfem::VectorFEMassIntegrator::AssemblePA ( const FiniteElementSpace trial_fes,
const FiniteElementSpace test_fes 
)
virtual

Used with BilinearFormIntegrators that have different spaces.

Reimplemented from mfem::BilinearFormIntegrator.

Definition at line 719 of file bilininteg_vectorfe.cpp.

Member Data Documentation

int mfem::VectorFEMassIntegrator::dim
protected

Definition at line 2438 of file bilininteg.hpp.

int mfem::VectorFEMassIntegrator::dofs1D
protected

Definition at line 2438 of file bilininteg.hpp.

int mfem::VectorFEMassIntegrator::dofs1Dtest
protected

Definition at line 2438 of file bilininteg.hpp.

const GeometricFactors* mfem::VectorFEMassIntegrator::geom
protected

Not owned.

Definition at line 2437 of file bilininteg.hpp.

const DofToQuad* mfem::VectorFEMassIntegrator::mapsC
protected

Not owned. DOF-to-quad map, closed.

Definition at line 2434 of file bilininteg.hpp.

const DofToQuad* mfem::VectorFEMassIntegrator::mapsCtest
protected

Not owned. DOF-to-quad map, closed.

Definition at line 2436 of file bilininteg.hpp.

const DofToQuad* mfem::VectorFEMassIntegrator::mapsO
protected

Not owned. DOF-to-quad map, open.

Definition at line 2433 of file bilininteg.hpp.

const DofToQuad* mfem::VectorFEMassIntegrator::mapsOtest
protected

Not owned. DOF-to-quad map, open.

Definition at line 2435 of file bilininteg.hpp.

MatrixCoefficient* mfem::VectorFEMassIntegrator::MQ
protected

Definition at line 2429 of file bilininteg.hpp.

int mfem::VectorFEMassIntegrator::ne
protected

Definition at line 2438 of file bilininteg.hpp.

int mfem::VectorFEMassIntegrator::nq
protected

Definition at line 2438 of file bilininteg.hpp.

Vector mfem::VectorFEMassIntegrator::pa_data
protected

Definition at line 2432 of file bilininteg.hpp.

Coefficient* mfem::VectorFEMassIntegrator::Q
protected

Definition at line 2427 of file bilininteg.hpp.

int mfem::VectorFEMassIntegrator::quad1D
protected

Definition at line 2438 of file bilininteg.hpp.

bool mfem::VectorFEMassIntegrator::symmetric = true
protected

False if using a nonsymmetric matrix coefficient.

Definition at line 2439 of file bilininteg.hpp.

int mfem::VectorFEMassIntegrator::test_fetype
protected

Definition at line 2438 of file bilininteg.hpp.

int mfem::VectorFEMassIntegrator::trial_fetype
protected

Definition at line 2438 of file bilininteg.hpp.

VectorCoefficient* mfem::VectorFEMassIntegrator::VQ
protected

Definition at line 2428 of file bilininteg.hpp.


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