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

alpha (q . grad u, v) More...

#include <bilininteg.hpp>

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

Public Member Functions

 ConvectionIntegrator (VectorCoefficient &q, double a=1.0)
 
virtual void AssembleElementMatrix (const FiniteElement &, ElementTransformation &, DenseMatrix &)
 Given a particular Finite Element computes the element matrix elmat. More...
 
virtual void AssemblePA (const FiniteElementSpace &)
 Method defining partial assembly. More...
 
virtual void AssembleEA (const FiniteElementSpace &fes, Vector &emat, const bool add)
 Method defining element assembly. More...
 
virtual void AddMultPA (const Vector &, Vector &) const
 Method for partially assembled action. More...
 
- 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 (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 AddMultTransposePA (const Vector &x, Vector &y) const
 Method for partially assembled transposed action. 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 AssembleElementMatrix2 (const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &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)
 
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 &el, ElementTransformation &Trans)
 
static const IntegrationRuleGetRule (const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans)
 

Protected Attributes

VectorCoefficientQ
 
double alpha
 
Vector pa_data
 
const DofToQuadmaps
 Not owned. More...
 
const GeometricFactorsgeom
 Not owned. More...
 
int dim
 
int ne
 
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

alpha (q . grad u, v)

Definition at line 2072 of file bilininteg.hpp.

Constructor & Destructor Documentation

mfem::ConvectionIntegrator::ConvectionIntegrator ( VectorCoefficient q,
double  a = 1.0 
)
inline

Definition at line 2090 of file bilininteg.hpp.

Member Function Documentation

void mfem::ConvectionIntegrator::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 872 of file bilininteg_convection_pa.cpp.

void mfem::ConvectionIntegrator::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 228 of file bilininteg_convection_ea.cpp.

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

void mfem::ConvectionIntegrator::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 769 of file bilininteg_convection_pa.cpp.

const IntegrationRule & mfem::ConvectionIntegrator::GetRule ( const FiniteElement el,
ElementTransformation Trans 
)
static

Definition at line 1131 of file bilininteg.cpp.

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

Definition at line 1121 of file bilininteg.cpp.

Member Data Documentation

double mfem::ConvectionIntegrator::alpha
protected

Definition at line 2076 of file bilininteg.hpp.

int mfem::ConvectionIntegrator::dim
protected

Definition at line 2081 of file bilininteg.hpp.

int mfem::ConvectionIntegrator::dofs1D
protected

Definition at line 2081 of file bilininteg.hpp.

const GeometricFactors* mfem::ConvectionIntegrator::geom
protected

Not owned.

Definition at line 2080 of file bilininteg.hpp.

const DofToQuad* mfem::ConvectionIntegrator::maps
protected

Not owned.

Definition at line 2079 of file bilininteg.hpp.

int mfem::ConvectionIntegrator::ne
protected

Definition at line 2081 of file bilininteg.hpp.

int mfem::ConvectionIntegrator::nq
protected

Definition at line 2081 of file bilininteg.hpp.

Vector mfem::ConvectionIntegrator::pa_data
protected

Definition at line 2078 of file bilininteg.hpp.

VectorCoefficient* mfem::ConvectionIntegrator::Q
protected

Definition at line 2075 of file bilininteg.hpp.

int mfem::ConvectionIntegrator::quad1D
protected

Definition at line 2081 of file bilininteg.hpp.


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