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

#include <bilininteg.hpp>

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

Public Member Functions

 DiffusionIntegrator ()
 Construct a diffusion integrator with coefficient Q = 1. More...
 
 DiffusionIntegrator (Coefficient &q)
 Construct a diffusion integrator with a scalar coefficient q. More...
 
 DiffusionIntegrator (VectorCoefficient &q)
 Construct a diffusion integrator with a vector coefficient q. More...
 
 DiffusionIntegrator (MatrixCoefficient &q)
 Construct a diffusion integrator with a matrix coefficient q. More...
 
virtual ~DiffusionIntegrator ()
 
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 AssembleElementVector (const FiniteElement &el, ElementTransformation &Tr, const Vector &elfun, Vector &elvect)
 Perform the local action of the BilinearFormIntegrator. 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 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...
 
- 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 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 ~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)
 

Protected Attributes

CoefficientQ
 
VectorCoefficientVQ
 
MatrixCoefficientMQ
 
- 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 integrating the bilinear form a(u,v) := (Q grad u, grad v) where Q can be a scalar or a matrix coefficient.

Definition at line 1897 of file bilininteg.hpp.

Constructor & Destructor Documentation

mfem::DiffusionIntegrator::DiffusionIntegrator ( )
inline

Construct a diffusion integrator with coefficient Q = 1.

Definition at line 1924 of file bilininteg.hpp.

mfem::DiffusionIntegrator::DiffusionIntegrator ( Coefficient q)
inline

Construct a diffusion integrator with a scalar coefficient q.

Definition at line 1928 of file bilininteg.hpp.

mfem::DiffusionIntegrator::DiffusionIntegrator ( VectorCoefficient q)
inline

Construct a diffusion integrator with a vector coefficient q.

Definition at line 1932 of file bilininteg.hpp.

mfem::DiffusionIntegrator::DiffusionIntegrator ( MatrixCoefficient q)
inline

Construct a diffusion integrator with a matrix coefficient q.

Definition at line 1936 of file bilininteg.hpp.

virtual mfem::DiffusionIntegrator::~DiffusionIntegrator ( )
inlinevirtual

Definition at line 1939 of file bilininteg.hpp.

Member Function Documentation

void mfem::DiffusionIntegrator::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 56 of file bilininteg_diffusion_mf.cpp.

void mfem::DiffusionIntegrator::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 1890 of file bilininteg_diffusion_pa.cpp.

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

Assemble diagonal and add it to Vector diag.

Reimplemented from mfem::BilinearFormIntegrator.

Definition at line 42 of file bilininteg_diffusion_mf.cpp.

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

Assemble diagonal and add it to Vector diag.

Reimplemented from mfem::BilinearFormIntegrator.

Definition at line 914 of file bilininteg_diffusion_pa.cpp.

void mfem::DiffusionIntegrator::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 245 of file bilininteg_diffusion_ea.cpp.

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

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

Reimplemented from mfem::BilinearFormIntegrator.

Definition at line 546 of file bilininteg.cpp.

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

Given a trial and test Finite Element computes the element stiffness matrix elmat.

Reimplemented from mfem::BilinearFormIntegrator.

Definition at line 604 of file bilininteg.cpp.

void mfem::DiffusionIntegrator::AssembleElementVector ( const FiniteElement el,
ElementTransformation Tr,
const Vector elfun,
Vector elvect 
)
virtual

Perform the local action of the BilinearFormIntegrator.

Reimplemented from mfem::BilinearFormIntegrator.

Definition at line 671 of file bilininteg.cpp.

void mfem::DiffusionIntegrator::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_diffusion_mf.cpp.

void mfem::DiffusionIntegrator::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 352 of file bilininteg_diffusion_pa.cpp.

void mfem::DiffusionIntegrator::ComputeElementFlux ( const FiniteElement el,
ElementTransformation Trans,
Vector u,
const FiniteElement fluxelem,
Vector flux,
bool  with_coef = true 
)
virtual

Virtual method required for Zienkiewicz-Zhu type error estimators.

The purpose of the method is to compute a local "flux" finite element function given a local finite element solution. The "flux" function has to be computed in terms of its coefficients (represented by the Vector flux) which multiply the basis functions defined by the FiniteElement fluxelem. Typically, the "flux" function will have more than one component and consequently flux should be store the coefficients of all components: first all coefficient for component 0, then all coefficients for component 1, etc. What the "flux" function represents depends on the specific integrator. For example, in the case of DiffusionIntegrator, the flux is the gradient of the solution multiplied by the diffusion coefficient.

Parameters
[in]elFiniteElement of the solution.
[in]TransThe ElementTransformation describing the physical position of the mesh element.
[in]uSolution coefficients representing the expansion of the solution function in the basis of el.
[in]fluxelemFiniteElement of the "flux".
[out]flux"Flux" coefficients representing the expansion of the "flux" function in the basis of fluxelem. The size of flux as a Vector has to be set by this method, e.g. using Vector::SetSize().
[in]with_coefIf zero (the default value is 1) the implementation of the method may choose not to scale the "flux" function by any coefficients describing the integrator.

Reimplemented from mfem::BilinearFormIntegrator.

Definition at line 744 of file bilininteg.cpp.

double mfem::DiffusionIntegrator::ComputeFluxEnergy ( const FiniteElement fluxelem,
ElementTransformation Trans,
Vector flux,
Vector d_energy = NULL 
)
virtual

Virtual method required for Zienkiewicz-Zhu type error estimators.

The purpose of this method is to compute a local number that measures the energy of a given "flux" function (see ComputeElementFlux() for a description of the "flux" function). Typically, the energy of a "flux" function should be equal to a_local(u,u), if the "flux" is defined from a solution u; here a_local(.,.) denotes the element-local bilinear form represented by the integrator.

Parameters
[in]fluxelemFiniteElement of the "flux".
[in]TransThe ElementTransformation describing the physical position of the mesh element.
[in]flux"Flux" coefficients representing the expansion of the "flux" function in the basis of fluxelem.
[out]d_energyIf not NULL, the given Vector should be set to represent directional energy split that can be used for anisotropic error estimation.
Returns
The computed energy.

Reimplemented from mfem::BilinearFormIntegrator.

Definition at line 802 of file bilininteg.cpp.

const IntegrationRule & mfem::DiffusionIntegrator::GetRule ( const FiniteElement trial_fe,
const FiniteElement test_fe 
)
static

Definition at line 868 of file bilininteg.cpp.

Member Data Documentation

MatrixCoefficient* mfem::DiffusionIntegrator::MQ
protected

Definition at line 1902 of file bilininteg.hpp.

Coefficient* mfem::DiffusionIntegrator::Q
protected

Definition at line 1900 of file bilininteg.hpp.

VectorCoefficient* mfem::DiffusionIntegrator::VQ
protected

Definition at line 1901 of file bilininteg.hpp.


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