MFEM
v4.5.1
Finite element discretization library
|
alpha (q . grad u, v) More...
#include <bilininteg.hpp>
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 | AssembleMF (const FiniteElementSpace &fes) |
Method defining matrix-free assembly. 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 | 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... | |
virtual void | AddMultTransposePA (const Vector &x, Vector &y) const |
Method for partially assembled transposed action. More... | |
bool | SupportsCeed () const |
Indicates whether this integrator can use a Ceed backend. 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 | 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 | 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, 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 | ~BilinearFormIntegrator () |
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 IntegrationRule * | GetIntegrationRule () 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... | |
ceed::Operator & | GetCeedOp () |
virtual | ~NonlinearFormIntegrator () |
Static Public Member Functions | |
static const IntegrationRule & | GetRule (const FiniteElement &el, ElementTransformation &Trans) |
static const IntegrationRule & | GetRule (const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans) |
Protected Attributes | |
VectorCoefficient * | Q |
double | alpha |
Vector | pa_data |
const DofToQuad * | maps |
Not owned. More... | |
const GeometricFactors * | geom |
Not owned. More... | |
int | dim |
int | ne |
int | nq |
int | dofs1D |
int | quad1D |
Protected Attributes inherited from mfem::NonlinearFormIntegrator | |
const IntegrationRule * | IntRule |
ceed::Operator * | ceedOp |
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) | |
alpha (q . grad u, v)
Definition at line 2251 of file bilininteg.hpp.
|
inline |
Definition at line 2269 of file bilininteg.hpp.
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 62 of file bilininteg_convection_mf.cpp.
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 1525 of file bilininteg_convection_pa.cpp.
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.
Definition at line 1540 of file bilininteg_convection_pa.cpp.
|
virtual |
Assemble diagonal and add it to Vector diag.
Reimplemented from mfem::BilinearFormIntegrator.
Definition at line 49 of file bilininteg_convection_mf.cpp.
|
virtual |
Assemble diagonal and add it to Vector diag.
Reimplemented from mfem::BilinearFormIntegrator.
Definition at line 1555 of file bilininteg_convection_pa.cpp.
|
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.
|
virtual |
Given a particular Finite Element computes the element matrix elmat.
Reimplemented from mfem::BilinearFormIntegrator.
Definition at line 1366 of file bilininteg.cpp.
|
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_convection_mf.cpp.
|
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 1378 of file bilininteg_convection_pa.cpp.
|
static |
Definition at line 1475 of file bilininteg.cpp.
|
static |
Definition at line 1466 of file bilininteg.cpp.
|
inlinevirtual |
Indicates whether this integrator can use a Ceed backend.
Reimplemented from mfem::NonlinearFormIntegrator.
Definition at line 2302 of file bilininteg.hpp.
|
protected |
Definition at line 2255 of file bilininteg.hpp.
|
protected |
Definition at line 2260 of file bilininteg.hpp.
|
protected |
Definition at line 2260 of file bilininteg.hpp.
|
protected |
Not owned.
Definition at line 2259 of file bilininteg.hpp.
|
protected |
Not owned.
Definition at line 2258 of file bilininteg.hpp.
|
protected |
Definition at line 2260 of file bilininteg.hpp.
|
protected |
Definition at line 2260 of file bilininteg.hpp.
|
protected |
Definition at line 2257 of file bilininteg.hpp.
|
protected |
Definition at line 2254 of file bilininteg.hpp.
|
protected |
Definition at line 2260 of file bilininteg.hpp.