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

Integrator for (curl u, curl v) for Nedelec elements. More...

#include <bilininteg.hpp>

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

Public Member Functions

 CurlCurlIntegrator ()
 
 CurlCurlIntegrator (Coefficient &q, const IntegrationRule *ir=NULL)
 Construct a bilinear form integrator for Nedelec elements. More...
 
 CurlCurlIntegrator (DiagonalMatrixCoefficient &dq, const IntegrationRule *ir=NULL)
 
 CurlCurlIntegrator (MatrixCoefficient &mq, const IntegrationRule *ir=NULL)
 
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 ComputeElementFlux (const FiniteElement &el, ElementTransformation &Trans, Vector &u, const FiniteElement &fluxelem, Vector &flux, bool with_coef, 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 void AssemblePA (const FiniteElementSpace &fes)
 Method defining partial assembly. More...
 
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...
 
const CoefficientGetCoefficient () const
 
virtual void AssemblePA (const FiniteElementSpace &fes)
 Method defining partial assembly. More...
 
virtual void AssemblePA (const FiniteElementSpace &trial_fes, const FiniteElementSpace &test_fes)
 
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::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 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 ~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 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
 
DiagonalMatrixCoefficientDQ
 
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 GeometricFactorsgeom
 Not owned. More...
 
int dim
 
int ne
 
int nq
 
int dofs1D
 
int quad1D
 
bool symmetric = true
 False if using a nonsymmetric matrix coefficient. More...
 
- Protected Attributes inherited from mfem::NonlinearFormIntegrator
const IntegrationRuleIntRule
 
ceed::OperatorceedOp
 
MemoryType pa_mt = MemoryType::DEFAULT
 

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 (curl u, curl v) for Nedelec elements.

Definition at line 2521 of file bilininteg.hpp.

Constructor & Destructor Documentation

◆ CurlCurlIntegrator() [1/4]

mfem::CurlCurlIntegrator::CurlCurlIntegrator ( )
inline

Definition at line 2546 of file bilininteg.hpp.

◆ CurlCurlIntegrator() [2/4]

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

Construct a bilinear form integrator for Nedelec elements.

Definition at line 2548 of file bilininteg.hpp.

◆ CurlCurlIntegrator() [3/4]

mfem::CurlCurlIntegrator::CurlCurlIntegrator ( DiagonalMatrixCoefficient dq,
const IntegrationRule ir = NULL 
)
inline

Definition at line 2550 of file bilininteg.hpp.

◆ CurlCurlIntegrator() [4/4]

mfem::CurlCurlIntegrator::CurlCurlIntegrator ( MatrixCoefficient mq,
const IntegrationRule ir = NULL 
)
inline

Definition at line 2553 of file bilininteg.hpp.

Member Function Documentation

◆ AddMultPA()

void mfem::CurlCurlIntegrator::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 2188 of file bilininteg_hcurl.cpp.

◆ AssembleDiagonalPA()

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

Assemble diagonal and add it to Vector diag.

Reimplemented from mfem::BilinearFormIntegrator.

Definition at line 2661 of file bilininteg_hcurl.cpp.

◆ AssembleElementMatrix()

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

◆ AssembleElementMatrix2()

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

◆ AssemblePA() [1/5]

void mfem::BilinearFormIntegrator::AssemblePA

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

Definition at line 23 of file bilininteg.cpp.

◆ AssemblePA() [2/5]

void mfem::BilinearFormIntegrator::AssemblePA

Used with BilinearFormIntegrators that have different spaces.

Definition at line 29 of file bilininteg.cpp.

◆ AssemblePA() [3/5]

void mfem::NonlinearFormIntegrator::AssemblePA

Method defining partial assembly.

The result of the partial assembly is stored internally so that it can be used later in the methods AddMultPA().

Definition at line 25 of file nonlininteg.cpp.

◆ AssemblePA() [4/5]

void mfem::NonlinearFormIntegrator::AssemblePA

The result of the partial assembly is stored internally so that it can be used later in the methods AddMultPA(). Used with BilinearFormIntegrators that have different spaces.

Definition at line 31 of file nonlininteg.cpp.

◆ AssemblePA() [5/5]

void mfem::CurlCurlIntegrator::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 950 of file bilininteg_hcurl.cpp.

◆ ComputeElementFlux()

void mfem::CurlCurlIntegrator::ComputeElementFlux ( const FiniteElement el,
ElementTransformation Trans,
Vector u,
const FiniteElement fluxelem,
Vector flux,
bool  with_coef,
const IntegrationRule ir = NULL 
)
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.
[in]irIf passed (the default value is NULL), the implementation of the method will ignore the integration rule provided by the fluxelem parameter and, instead, compute the discrete flux at the points specified by the integration rule ir.

Reimplemented from mfem::BilinearFormIntegrator.

Definition at line 2084 of file bilininteg.cpp.

◆ ComputeFluxEnergy()

double mfem::CurlCurlIntegrator::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 2102 of file bilininteg.cpp.

◆ GetCoefficient()

const Coefficient* mfem::CurlCurlIntegrator::GetCoefficient ( ) const
inline

Definition at line 2582 of file bilininteg.hpp.

Member Data Documentation

◆ dim

int mfem::CurlCurlIntegrator::dim
protected

Definition at line 2542 of file bilininteg.hpp.

◆ dofs1D

int mfem::CurlCurlIntegrator::dofs1D
protected

Definition at line 2542 of file bilininteg.hpp.

◆ DQ

DiagonalMatrixCoefficient* mfem::CurlCurlIntegrator::DQ
protected

Definition at line 2534 of file bilininteg.hpp.

◆ geom

const GeometricFactors* mfem::CurlCurlIntegrator::geom
protected

Not owned.

Definition at line 2541 of file bilininteg.hpp.

◆ mapsC

const DofToQuad* mfem::CurlCurlIntegrator::mapsC
protected

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

Definition at line 2540 of file bilininteg.hpp.

◆ mapsO

const DofToQuad* mfem::CurlCurlIntegrator::mapsO
protected

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

Definition at line 2539 of file bilininteg.hpp.

◆ MQ

MatrixCoefficient* mfem::CurlCurlIntegrator::MQ
protected

Definition at line 2535 of file bilininteg.hpp.

◆ ne

int mfem::CurlCurlIntegrator::ne
protected

Definition at line 2542 of file bilininteg.hpp.

◆ nq

int mfem::CurlCurlIntegrator::nq
protected

Definition at line 2542 of file bilininteg.hpp.

◆ pa_data

Vector mfem::CurlCurlIntegrator::pa_data
protected

Definition at line 2538 of file bilininteg.hpp.

◆ Q

Coefficient* mfem::CurlCurlIntegrator::Q
protected

Definition at line 2533 of file bilininteg.hpp.

◆ quad1D

int mfem::CurlCurlIntegrator::quad1D
protected

Definition at line 2542 of file bilininteg.hpp.

◆ symmetric

bool mfem::CurlCurlIntegrator::symmetric = true
protected

False if using a nonsymmetric matrix coefficient.

Definition at line 2543 of file bilininteg.hpp.


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