MFEM v4.7.0
Finite element discretization library
|
#include <bilininteg.hpp>
Public Member Functions | |
virtual void | AssembleElementMatrix2 (const FiniteElement &dom_fe, const FiniteElement &ran_fe, ElementTransformation &Trans, DenseMatrix &elmat) |
Public Member Functions inherited from mfem::BilinearFormIntegrator | |
virtual void | AssemblePA (const FiniteElementSpace &fes) |
Method defining partial assembly. | |
virtual void | AssemblePA (const FiniteElementSpace &trial_fes, const FiniteElementSpace &test_fes) |
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 (Vector &diag) |
Assemble diagonal and add it to Vector diag. | |
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 | AddMultPA (const Vector &x, Vector &y) const |
Method for partially assembled action. | |
virtual void | AddMultNURBSPA (const Vector &x, Vector &y) const |
Method for partially assembled action on NURBS patches. | |
virtual void | AddMultTransposePA (const Vector &x, Vector &y) const |
Method for partially assembled transposed action. | |
virtual void | AssembleEA (const FiniteElementSpace &fes, Vector &emat, const bool add=true) |
Method defining element assembly. | |
virtual void | AssembleMF (const FiniteElementSpace &fes) |
Method defining matrix-free assembly. | |
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. | |
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 | AssembleElementMatrix (const FiniteElement &el, ElementTransformation &Trans, DenseMatrix &elmat) |
Given a particular Finite Element computes the element matrix 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_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) |
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. | |
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. | |
virtual void | AssembleElementGrad (const FiniteElement &el, ElementTransformation &Tr, const Vector &elfun, DenseMatrix &elmat) |
Assemble the local gradient matrix. | |
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. | |
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 | |
virtual void | SetIntRule (const IntegrationRule *ir) |
Prescribe a fixed IntegrationRule to use (when ir != NULL) or let the integrator choose (when ir == NULL). | |
void | SetIntegrationMode (Mode m) |
void | SetNURBSPatchIntRule (NURBSMeshRules *pr) |
For patchwise integration, SetNURBSPatchIntRule must be called. | |
bool | HasNURBSPatchIntRule () const |
bool | Patchwise () const |
void | SetIntegrationRule (const IntegrationRule &ir) |
Prescribe a fixed IntegrationRule to use. | |
void | SetPAMemoryType (MemoryType mt) |
const IntegrationRule * | GetIntegrationRule () const |
Get the integration rule of the integrator (possibly NULL). | |
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. | |
virtual bool | SupportsCeed () const |
Indicates whether this integrator can use a Ceed backend. | |
ceed::Operator & | GetCeedOp () |
virtual | ~NonlinearFormIntegrator () |
Additional Inherited Members | |
Public Types inherited from mfem::NonlinearFormIntegrator | |
enum | Mode { ELEMENTWISE = 0 , PATCHWISE = 1 , PATCHWISE_REDUCED = 2 } |
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 Attributes inherited from mfem::NonlinearFormIntegrator | |
const IntegrationRule * | IntRule |
Mode | integrationMode = Mode::ELEMENTWISE |
NURBSMeshRules * | patchRules = nullptr |
ceed::Operator * | ceedOp |
MemoryType | pa_mt = MemoryType::DEFAULT |
A trace face interpolator class for interpolating the normal component of the domain space, e.g. vector \(H^1\), into the range space, e.g. the trace of \(H(div)\) which uses FiniteElement::INTEGRAL map type.
Definition at line 3705 of file bilininteg.hpp.
|
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 4277 of file bilininteg.cpp.