MFEM v4.9.0
Finite element discretization library
Loading...
Searching...
No Matches
mfem::ConvectionIntegrator Class Reference

\(\alpha (Q \cdot \nabla u, v)\) More...

#include <bilininteg.hpp>

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

Classes

struct  Kernels
 

Public Types

using ApplyKernelType
 arguments: NE, B, G, Bt, Gt, pa_data, x, y, D1D, Q1D
 
- Public Types inherited from mfem::NonlinearFormIntegrator
enum  Mode { ELEMENTWISE = 0 , PATCHWISE = 1 , PATCHWISE_REDUCED = 2 }
 

Public Member Functions

 ConvectionIntegrator (VectorCoefficient &q, real_t a=1.0)
 
void AssembleElementMatrix (const FiniteElement &, ElementTransformation &, DenseMatrix &) override
 Given a particular Finite Element computes the element matrix elmat.
 
void AssembleMF (const FiniteElementSpace &fes) override
 Method defining matrix-free assembly.
 
void AssemblePA (const FiniteElementSpace &) override
 Method defining partial assembly.
 
void AssembleEA (const FiniteElementSpace &fes, Vector &emat, const bool add) override
 Method defining element assembly.
 
void AssembleDiagonalPA (Vector &diag) override
 Assemble diagonal and add it to Vector diag.
 
void AssembleDiagonalMF (Vector &diag) override
 Assemble diagonal and add it to Vector diag.
 
void AddMultMF (const Vector &, Vector &) const override
 
void AddMultPA (const Vector &, Vector &) const override
 Method for partially assembled action.
 
void AddMultTransposePA (const Vector &x, Vector &y) const override
 Method for partially assembled transposed action.
 
bool SupportsCeed () const override
 Indicates whether this integrator can use a Ceed backend.
 
 MFEM_REGISTER_KERNELS (ApplyPAKernels, ApplyKernelType,(int, int, int))
 arguments: DIMS, D1D, Q1D
 
 MFEM_REGISTER_KERNELS (ApplyPATKernels, ApplyKernelType,(int, int, int))
 arguments: DIMS, D1D, Q1D
 
void AssemblePA (const FiniteElementSpace &trial_fes, const FiniteElementSpace &test_fes) override
 
- Public Member Functions inherited from mfem::BilinearFormIntegrator
virtual void AssembleNURBSPA (const FiniteElementSpace &fes)
 Method defining partial assembly on NURBS patches.
 
virtual void AssemblePABoundary (const FiniteElementSpace &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 \(A D A^T\) ( \(A\) is this integrator) and add it to diag.
 
virtual void AddAbsMultPA (const Vector &x, Vector &y) const
 
virtual void AddMultNURBSPA (const Vector &x, Vector &y) const
 Method for partially assembled action on NURBS patches.
 
virtual void AddAbsMultTransposePA (const Vector &x, Vector &y) const
 
virtual void AddMultTransposeMF (const Vector &x, Vector &y) const
 
virtual void AssembleEABoundary (const FiniteElementSpace &fes, Vector &ea_data_bdr, const bool add=true)
 
virtual void AssembleEAInteriorFaces (const FiniteElementSpace &fes, Vector &ea_data_int, Vector &ea_data_ext, const bool add=true)
 
virtual void AssembleEAInteriorFaces (const FiniteElementSpace &trial_fes, const FiniteElementSpace &test_fes, Vector &emat, const bool add=true)
 Method defining element assembly for mixed trace integrators.
 
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 AssemblePatchMatrix (const int patch, const FiniteElementSpace &fes, SparseMatrix *&smat)
 
virtual void AssembleFaceMatrix (const FiniteElement &el1, const FiniteElement &el2, FaceElementTransformations &Trans, DenseMatrix &elmat)
 
virtual void AssembleFaceMatrix (const FiniteElement &trial_fe1, const FiniteElement &test_fe1, const FiniteElement &trial_fe2, const FiniteElement &test_fe2, 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 AssembleTraceFaceMatrix (int elem, const FiniteElement &trial_face_fe, const FiniteElement &test_fe, FaceElementTransformations &Trans, DenseMatrix &elmat)
 
void AssembleElementVector (const FiniteElement &el, ElementTransformation &Tr, const Vector &elfun, Vector &elvect) override
 Perform the local action of the BilinearFormIntegrator. Note that the default implementation in the base class is general but not efficient.
 
void AssembleFaceVector (const FiniteElement &el1, const FiniteElement &el2, FaceElementTransformations &Tr, const Vector &elfun, Vector &elvect) override
 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.
 
void AssembleElementGrad (const FiniteElement &el, ElementTransformation &Tr, const Vector &elfun, DenseMatrix &elmat) override
 Assemble the local gradient matrix.
 
void AssembleFaceGrad (const FiniteElement &el1, const FiniteElement &el2, FaceElementTransformations &Tr, const Vector &elfun, DenseMatrix &elmat) override
 Assemble the local action of the gradient of the NonlinearFormIntegrator resulting from a face integral term.
 
virtual void ComputeElementFlux (const FiniteElement &el, ElementTransformation &Trans, Vector &u, const FiniteElement &fluxelem, Vector &flux, bool with_coef=true, const IntegrationRule *ir=NULL)
 Virtual method required for Zienkiewicz-Zhu type error estimators.
 
virtual real_t ComputeFluxEnergy (const FiniteElement &fluxelem, ElementTransformation &Trans, Vector &flux, Vector *d_energy=NULL)
 Virtual method required for Zienkiewicz-Zhu type error estimators.
 
virtual bool RequiresFaceNormalDerivatives () const
 For bilinear forms on element faces, specifies if the normal derivatives are needed on the faces or just the face restriction.
 
virtual void AddMultPAFaceNormalDerivatives (const Vector &x, const Vector &dxdn, Vector &y, Vector &dydn) const
 Method for partially assembled action.
 
virtual ~BilinearFormIntegrator ()
 
- Public Member Functions inherited from mfem::NonlinearFormIntegrator
void SetIntegrationMode (Mode m)
 
bool Patchwise () const
 
void SetPAMemoryType (MemoryType mt)
 
virtual real_t GetElementEnergy (const FiniteElement &el, ElementTransformation &Tr, const Vector &elfun)
 Compute the local energy.
 
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.
 
virtual real_t GetLocalStateEnergyPA (const Vector &x) const
 Compute the local (to the MPI rank) energy with partial assembly.
 
virtual void AddMultGradPA (const Vector &x, Vector &y) const
 Method for partially assembled gradient action.
 
virtual void AssembleGradDiagonalPA (Vector &diag) const
 Method for computing the diagonal of the gradient with partial assembly.
 
ceed::OperatorGetCeedOp ()
 
virtual ~NonlinearFormIntegrator ()
 
- Public Member Functions inherited from mfem::Integrator
 Integrator (const IntegrationRule *ir=NULL)
 Create a new Integrator, optionally providing a prescribed quadrature rule to use in assembly.
 
virtual void SetIntRule (const IntegrationRule *ir)
 Prescribe a fixed IntegrationRule to use, or set to null to let the integrator choose an appropriate rule.
 
void SetIntegrationRule (const IntegrationRule &ir)
 Prescribe a fixed IntegrationRule to use. Sets the NURBS patch integration rule to null.
 
void SetNURBSPatchIntRule (NURBSMeshRules *pr)
 Sets an integration rule for use on NURBS patches.
 
bool HasNURBSPatchIntRule () const
 Check if a NURBS patch integration rule has been set.
 
const IntegrationRuleGetIntRule () const
 Directly return the IntRule pointer (possibly null) without checking for NURBS patch rules or falling back on a default.
 
const IntegrationRuleGetIntegrationRule () const
 Equivalent to GetIntRule, but retained for backward compatibility with applications.
 

Static Public Member Functions

static const IntegrationRuleGetRule (const FiniteElement &el, const ElementTransformation &Trans)
 
static const IntegrationRuleGetRule (const FiniteElement &trial_fe, const FiniteElement &test_fe, const ElementTransformation &Trans)
 
template<int DIM, int D1D, int Q1D>
static void AddSpecialization ()
 

Protected Member Functions

const IntegrationRuleGetDefaultIntegrationRule (const FiniteElement &trial_fe, const FiniteElement &test_fe, const ElementTransformation &trans) const override
 Subclasses should override to choose a default integration rule.
 
- Protected Member Functions inherited from mfem::BilinearFormIntegrator
 BilinearFormIntegrator (const IntegrationRule *ir=NULL)
 
- Protected Member Functions inherited from mfem::NonlinearFormIntegrator
 NonlinearFormIntegrator (const IntegrationRule *ir=NULL)
 
- Protected Member Functions inherited from mfem::Integrator
const IntegrationRuleGetIntegrationRule (const FiniteElement &trial_fe, const FiniteElement &test_fe, const ElementTransformation &trans) const
 Returns an integration rule based on the arguments and internal state of the Integrator object.
 
const IntegrationRuleGetIntegrationRule (const FiniteElement &el, const ElementTransformation &trans) const
 Returns an integration rule based on the arguments and internal state. (Version for identical trial_fe and test_fe)
 

Protected Attributes

VectorCoefficientQ
 
real_t alpha
 
Vector pa_data
 
const DofToQuadmaps
 Not owned.
 
const GeometricFactorsgeom
 Not owned.
 
int dim
 
int ne
 
int nq
 
int dofs1D
 
int quad1D
 
- Protected Attributes inherited from mfem::NonlinearFormIntegrator
Mode integrationMode = Mode::ELEMENTWISE
 
ceed::OperatorceedOp
 
MemoryType pa_mt = MemoryType::DEFAULT
 
- Protected Attributes inherited from mfem::Integrator
const IntegrationRuleIntRule
 
NURBSMeshRulespatchRules = nullptr
 

Detailed Description

\(\alpha (Q \cdot \nabla u, v)\)

Definition at line 2481 of file bilininteg.hpp.

Member Typedef Documentation

◆ ApplyKernelType

Initial value:
void (*)(const int, const Array<real_t> &,
const Array<real_t> &,
const Array<real_t> &,
const Array<real_t> &, const Vector &,
const Vector &, Vector &, const int,
const int)

arguments: NE, B, G, Bt, Gt, pa_data, x, y, D1D, Q1D

Definition at line 2533 of file bilininteg.hpp.

Constructor & Destructor Documentation

◆ ConvectionIntegrator()

mfem::ConvectionIntegrator::ConvectionIntegrator ( VectorCoefficient & q,
real_t a = 1.0 )

Member Function Documentation

◆ AddMultMF()

void mfem::ConvectionIntegrator::AddMultMF ( const Vector & x,
Vector & y ) const
overridevirtual

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 59 of file bilininteg_convection_mf.cpp.

◆ AddMultPA()

void mfem::ConvectionIntegrator::AddMultPA ( const Vector & x,
Vector & y ) const
overridevirtual

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.

◆ AddMultTransposePA()

void mfem::ConvectionIntegrator::AddMultTransposePA ( const Vector & x,
Vector & y ) const
overridevirtual

Method for partially assembled transposed action.

Perform the transpose 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.

◆ AddSpecialization()

template<int DIM, int D1D, int Q1D>
static void mfem::ConvectionIntegrator::AddSpecialization ( )
inlinestatic

Definition at line 2546 of file bilininteg.hpp.

◆ AssembleDiagonalMF()

void mfem::ConvectionIntegrator::AssembleDiagonalMF ( Vector & diag)
overridevirtual

Assemble diagonal and add it to Vector diag.

Reimplemented from mfem::BilinearFormIntegrator.

Definition at line 46 of file bilininteg_convection_mf.cpp.

◆ AssembleDiagonalPA()

void mfem::ConvectionIntegrator::AssembleDiagonalPA ( Vector & diag)
overridevirtual

Assemble diagonal and add it to Vector diag.

Reimplemented from mfem::BilinearFormIntegrator.

◆ AssembleEA()

void mfem::ConvectionIntegrator::AssembleEA ( const FiniteElementSpace & fes,
Vector & emat,
const bool add )
overridevirtual

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.

◆ AssembleElementMatrix()

void mfem::ConvectionIntegrator::AssembleElementMatrix ( const FiniteElement & el,
ElementTransformation & Trans,
DenseMatrix & elmat )
overridevirtual

Given a particular Finite Element computes the element matrix elmat.

Reimplemented from mfem::BilinearFormIntegrator.

Definition at line 1509 of file bilininteg.cpp.

◆ AssembleMF()

void mfem::ConvectionIntegrator::AssembleMF ( const FiniteElementSpace & fes)
overridevirtual

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 19 of file bilininteg_convection_mf.cpp.

◆ AssemblePA() [1/2]

void mfem::ConvectionIntegrator::AssemblePA ( const FiniteElementSpace & fes)
overridevirtual

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.

◆ AssemblePA() [2/2]

void mfem::BilinearFormIntegrator::AssemblePA ( const FiniteElementSpace & trial_fes,
const FiniteElementSpace & test_fes )
overridevirtual

Used with BilinearFormIntegrators that have different spaces.

Reimplemented from mfem::BilinearFormIntegrator.

Definition at line 54 of file bilininteg.cpp.

◆ GetDefaultIntegrationRule()

const IntegrationRule * mfem::ConvectionIntegrator::GetDefaultIntegrationRule ( const FiniteElement & trial_fe,
const FiniteElement & test_fe,
const ElementTransformation & trans ) const
inlineoverrideprotectedvirtual

Subclasses should override to choose a default integration rule.

This method is intended to be overridden by subclasses to choose an appropriate integration rule based on the finite element spaces and/or element transformation. The trial_fe and test_fe should be equal for linear forms. The default base-class implementation returns null, which assumes that an appropriate rule is provided by another means, or that null integration rules are handled appropriately by the caller.

Reimplemented from mfem::Integrator.

Definition at line 2555 of file bilininteg.hpp.

◆ GetRule() [1/2]

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

Definition at line 1619 of file bilininteg.cpp.

◆ GetRule() [2/2]

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

Definition at line 1610 of file bilininteg.cpp.

◆ MFEM_REGISTER_KERNELS() [1/2]

mfem::ConvectionIntegrator::MFEM_REGISTER_KERNELS ( ApplyPAKernels ,
ApplyKernelType ,
(int, int, int)  )

arguments: DIMS, D1D, Q1D

◆ MFEM_REGISTER_KERNELS() [2/2]

mfem::ConvectionIntegrator::MFEM_REGISTER_KERNELS ( ApplyPATKernels ,
ApplyKernelType ,
(int, int, int)  )

arguments: DIMS, D1D, Q1D

◆ SupportsCeed()

bool mfem::ConvectionIntegrator::SupportsCeed ( ) const
inlineoverridevirtual

Indicates whether this integrator can use a Ceed backend.

Reimplemented from mfem::NonlinearFormIntegrator.

Definition at line 2530 of file bilininteg.hpp.

Member Data Documentation

◆ alpha

real_t mfem::ConvectionIntegrator::alpha
protected

Definition at line 2485 of file bilininteg.hpp.

◆ dim

int mfem::ConvectionIntegrator::dim
protected

Definition at line 2490 of file bilininteg.hpp.

◆ dofs1D

int mfem::ConvectionIntegrator::dofs1D
protected

Definition at line 2490 of file bilininteg.hpp.

◆ geom

const GeometricFactors* mfem::ConvectionIntegrator::geom
protected

Not owned.

Definition at line 2489 of file bilininteg.hpp.

◆ maps

const DofToQuad* mfem::ConvectionIntegrator::maps
protected

Not owned.

Definition at line 2488 of file bilininteg.hpp.

◆ ne

int mfem::ConvectionIntegrator::ne
protected

Definition at line 2490 of file bilininteg.hpp.

◆ nq

int mfem::ConvectionIntegrator::nq
protected

Definition at line 2490 of file bilininteg.hpp.

◆ pa_data

Vector mfem::ConvectionIntegrator::pa_data
protected

Definition at line 2487 of file bilininteg.hpp.

◆ Q

VectorCoefficient* mfem::ConvectionIntegrator::Q
protected

Definition at line 2484 of file bilininteg.hpp.

◆ quad1D

int mfem::ConvectionIntegrator::quad1D
protected

Definition at line 2490 of file bilininteg.hpp.


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