MFEM v4.7.0
Finite element discretization library
|
Integrator for \((\mathrm{curl}(u), \mathrm{curl}(v))\) for Nedelec elements. More...
#include <bilininteg.hpp>
Public Member Functions | |
CurlCurlIntegrator () | |
CurlCurlIntegrator (Coefficient &q, const IntegrationRule *ir=NULL) | |
Construct a bilinear form integrator for Nedelec elements. | |
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. | |
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. | |
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 void | AssemblePA (const FiniteElementSpace &fes) |
Method defining partial assembly. | |
virtual void | AddMultPA (const Vector &x, Vector &y) const |
Method for partially assembled action. | |
virtual void | AssembleDiagonalPA (Vector &diag) |
Assemble diagonal and add it to Vector diag. | |
const Coefficient * | GetCoefficient () const |
virtual void | AssemblePA (const FiniteElementSpace &trial_fes, const FiniteElementSpace &test_fes) |
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 | 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 | 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 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 () |
Protected Attributes | |
Coefficient * | Q |
DiagonalMatrixCoefficient * | DQ |
MatrixCoefficient * | MQ |
Vector | pa_data |
const DofToQuad * | mapsO |
Not owned. DOF-to-quad map, open. | |
const DofToQuad * | mapsC |
Not owned. DOF-to-quad map, closed. | |
const GeometricFactors * | geom |
Not owned. | |
int | dim |
int | ne |
int | nq |
int | dofs1D |
int | quad1D |
bool | symmetric = true |
False if using a nonsymmetric matrix coefficient. | |
Protected Attributes inherited from mfem::NonlinearFormIntegrator | |
const IntegrationRule * | IntRule |
Mode | integrationMode = Mode::ELEMENTWISE |
NURBSMeshRules * | patchRules = nullptr |
ceed::Operator * | ceedOp |
MemoryType | pa_mt = MemoryType::DEFAULT |
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) | |
Integrator for \((\mathrm{curl}(u), \mathrm{curl}(v))\) for Nedelec elements.
Definition at line 2639 of file bilininteg.hpp.
|
inline |
Definition at line 2664 of file bilininteg.hpp.
|
inline |
Construct a bilinear form integrator for Nedelec elements.
Definition at line 2666 of file bilininteg.hpp.
|
inline |
Definition at line 2668 of file bilininteg.hpp.
|
inline |
Definition at line 2671 of file bilininteg.hpp.
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.
|
virtual |
Assemble diagonal and add it to Vector diag.
Reimplemented from mfem::BilinearFormIntegrator.
|
virtual |
Given a particular Finite Element computes the element matrix elmat.
Reimplemented from mfem::BilinearFormIntegrator.
Definition at line 1996 of file bilininteg.cpp.
|
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 2067 of file bilininteg.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.
|
virtual |
Used with BilinearFormIntegrators that have different spaces.
Reimplemented from mfem::BilinearFormIntegrator.
Definition at line 50 of file bilininteg.cpp.
|
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.
[in] | el | FiniteElement of the solution. |
[in] | Trans | The ElementTransformation describing the physical position of the mesh element. |
[in] | u | Solution coefficients representing the expansion of the solution function in the basis of el. |
[in] | fluxelem | FiniteElement 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_coef | If 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] | ir | If 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 2144 of file bilininteg.cpp.
|
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.
[in] | fluxelem | FiniteElement of the "flux". |
[in] | Trans | The 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_energy | If not NULL, the given Vector should be set to represent directional energy split that can be used for anisotropic error estimation. |
Reimplemented from mfem::BilinearFormIntegrator.
Definition at line 2163 of file bilininteg.cpp.
|
inline |
Definition at line 2700 of file bilininteg.hpp.
|
protected |
Definition at line 2660 of file bilininteg.hpp.
|
protected |
Definition at line 2660 of file bilininteg.hpp.
|
protected |
Definition at line 2652 of file bilininteg.hpp.
|
protected |
Not owned.
Definition at line 2659 of file bilininteg.hpp.
|
protected |
Not owned. DOF-to-quad map, closed.
Definition at line 2658 of file bilininteg.hpp.
|
protected |
Not owned. DOF-to-quad map, open.
Definition at line 2657 of file bilininteg.hpp.
|
protected |
Definition at line 2653 of file bilininteg.hpp.
|
protected |
Definition at line 2660 of file bilininteg.hpp.
|
protected |
Definition at line 2660 of file bilininteg.hpp.
|
protected |
Definition at line 2656 of file bilininteg.hpp.
|
protected |
Definition at line 2651 of file bilininteg.hpp.
|
protected |
Definition at line 2660 of file bilininteg.hpp.
|
protected |
False if using a nonsymmetric matrix coefficient.
Definition at line 2661 of file bilininteg.hpp.