MFEM v4.8.0
Finite element discretization library
Loading...
Searching...
No Matches
mfem::SBM2NeumannIntegrator Class Reference

#include <sbm_solver.hpp>

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

Public Member Functions

 SBM2NeumannIntegrator (const ParMesh *pmesh, VectorCoefficient &vD_, ShiftedVectorFunctionCoefficient &vN_, Array< int > &elem_marker_, Array< int > &cut_marker_, bool include_cut_cell_=false, int nterms_=1)
 
void AssembleFaceMatrix (const FiniteElement &el1, const FiniteElement &el2, FaceElementTransformations &Trans, DenseMatrix &elmat) override
 
bool GetTrimFlag () const
 
virtual ~SBM2NeumannIntegrator ()
 
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)
 
- Public Member Functions inherited from mfem::BilinearFormIntegrator
void AssemblePA (const FiniteElementSpace &fes) override
 Method defining partial assembly.
 
void AssemblePA (const FiniteElementSpace &trial_fes, const FiniteElementSpace &test_fes) override
 
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.
 
void AddMultPA (const Vector &x, Vector &y) const override
 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.
 
void AssembleMF (const FiniteElementSpace &fes) override
 Method defining matrix-free assembly.
 
void AddMultMF (const Vector &x, Vector &y) const override
 
virtual void AddMultTransposeMF (const Vector &x, Vector &y) const
 
virtual void AssembleDiagonalMF (Vector &diag)
 Assemble diagonal and add it to Vector diag.
 
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 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 AssemblePatchMatrix (const int patch, const FiniteElementSpace &fes, SparseMatrix *&smat)
 
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.
 
virtual bool SupportsCeed () const
 Indicates whether this integrator can use a Ceed backend.
 
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.
 

Protected Attributes

ShiftedVectorFunctionCoefficientvN
 
VectorCoefficientvD
 
Array< int > * elem_marker
 
bool include_cut_cell
 
int nterms
 
int NEproc
 
int par_shared_face_count
 
Array< int > cut_marker
 
Vector shape
 
Vector dshapedn
 
Vector nor
 
Vector nh
 
Vector ni
 
DenseMatrix dshape
 
DenseMatrix adjJ
 
- 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
 

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 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)
 
virtual const IntegrationRuleGetDefaultIntegrationRule (const FiniteElement &trial_fe, const FiniteElement &test_fe, const ElementTransformation &trans) const
 Subclasses should override to choose a default integration rule.
 

Detailed Description

BilinearFormIntegrator for Neumann boundaries using the shifted boundary method.

\[ A(u,w) = \langle [\nabla u + \nabla(\nabla u) \cdot d + h.o.t.] \cdot \hat{n} \, (n \cdot \hat{n}),w ⟩ - \langle \nabla u \cdot n,w \rangle \]

where h.o.t are the high-order terms due to Taylor expansion for \(\nabla u\), \(\hat{n}\) is the normal vector at the true boundary, \(n\) is the normal vector at the surrogate boundary. Since this interior face integrator is applied to the surrogate boundary (see marking.hpp for notes on how the surrogate faces are determined and elements are marked), this integrator adds contribution to only the element that is adjacent to that face (Trans.Elem1 or Trans.Elem2) and is part of the surrogate domain.

Definition at line 214 of file sbm_solver.hpp.

Constructor & Destructor Documentation

◆ SBM2NeumannIntegrator()

mfem::SBM2NeumannIntegrator::SBM2NeumannIntegrator ( const ParMesh * pmesh,
VectorCoefficient & vD_,
ShiftedVectorFunctionCoefficient & vN_,
Array< int > & elem_marker_,
Array< int > & cut_marker_,
bool include_cut_cell_ = false,
int nterms_ = 1 )
inline

Definition at line 235 of file sbm_solver.hpp.

◆ ~SBM2NeumannIntegrator()

virtual mfem::SBM2NeumannIntegrator::~SBM2NeumannIntegrator ( )
inlinevirtual

Definition at line 258 of file sbm_solver.hpp.

Member Function Documentation

◆ AssembleFaceMatrix() [1/3]

void mfem::SBM2NeumannIntegrator::AssembleFaceMatrix ( const FiniteElement & el1,
const FiniteElement & el2,
FaceElementTransformations & Trans,
DenseMatrix & elmat )
overridevirtual

Reimplemented from mfem::BilinearFormIntegrator.

Definition at line 646 of file sbm_solver.cpp.

◆ AssembleFaceMatrix() [2/3]

void mfem::BilinearFormIntegrator::AssembleFaceMatrix ( const FiniteElement & trial_face_fe,
const FiniteElement & test_fe1,
const FiniteElement & test_fe2,
FaceElementTransformations & Trans,
DenseMatrix & elmat )
virtual

Abstract method used for assembling TraceFaceIntegrators in a MixedBilinearForm.

Reimplemented from mfem::BilinearFormIntegrator.

Definition at line 185 of file bilininteg.cpp.

◆ AssembleFaceMatrix() [3/3]

void mfem::BilinearFormIntegrator::AssembleFaceMatrix ( const FiniteElement & trial_fe1,
const FiniteElement & test_fe1,
const FiniteElement & trial_fe2,
const FiniteElement & test_fe2,
FaceElementTransformations & Trans,
DenseMatrix & elmat )
virtual

Reimplemented from mfem::BilinearFormIntegrator.

Definition at line 176 of file bilininteg.cpp.

◆ GetTrimFlag()

bool mfem::SBM2NeumannIntegrator::GetTrimFlag ( ) const
inline

Definition at line 256 of file sbm_solver.hpp.

Member Data Documentation

◆ adjJ

DenseMatrix mfem::SBM2NeumannIntegrator::adjJ
protected

Definition at line 231 of file sbm_solver.hpp.

◆ cut_marker

Array<int> mfem::SBM2NeumannIntegrator::cut_marker
protected

Definition at line 226 of file sbm_solver.hpp.

◆ dshape

DenseMatrix mfem::SBM2NeumannIntegrator::dshape
protected

Definition at line 231 of file sbm_solver.hpp.

◆ dshapedn

Vector mfem::SBM2NeumannIntegrator::dshapedn
protected

Definition at line 230 of file sbm_solver.hpp.

◆ elem_marker

Array<int>* mfem::SBM2NeumannIntegrator::elem_marker
protected

Definition at line 219 of file sbm_solver.hpp.

◆ include_cut_cell

bool mfem::SBM2NeumannIntegrator::include_cut_cell
protected

Definition at line 221 of file sbm_solver.hpp.

◆ NEproc

int mfem::SBM2NeumannIntegrator::NEproc
protected

Definition at line 224 of file sbm_solver.hpp.

◆ nh

Vector mfem::SBM2NeumannIntegrator::nh
protected

Definition at line 230 of file sbm_solver.hpp.

◆ ni

Vector mfem::SBM2NeumannIntegrator::ni
protected

Definition at line 230 of file sbm_solver.hpp.

◆ nor

Vector mfem::SBM2NeumannIntegrator::nor
protected

Definition at line 230 of file sbm_solver.hpp.

◆ nterms

int mfem::SBM2NeumannIntegrator::nterms
protected

Definition at line 222 of file sbm_solver.hpp.

◆ par_shared_face_count

int mfem::SBM2NeumannIntegrator::par_shared_face_count
protected

Definition at line 225 of file sbm_solver.hpp.

◆ shape

Vector mfem::SBM2NeumannIntegrator::shape
protected

Definition at line 230 of file sbm_solver.hpp.

◆ vD

VectorCoefficient* mfem::SBM2NeumannIntegrator::vD
protected

Definition at line 218 of file sbm_solver.hpp.

◆ vN

ShiftedVectorFunctionCoefficient* mfem::SBM2NeumannIntegrator::vN
protected

Definition at line 217 of file sbm_solver.hpp.


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