![]() |
MFEM v4.8.0
Finite element discretization library
|
#include <sbm_solver.hpp>
Public Member Functions | |
| SBM2NeumannLFIntegrator (const ParMesh *pmesh, ShiftedFunctionCoefficient &u, VectorCoefficient &vD_, ShiftedVectorFunctionCoefficient &vN_, Array< int > &elem_marker_, int nterms_=0, bool include_cut_cell_=false, int ls_cut_marker_=ShiftedFaceMarker::SBElementType::CUT) | |
| void | AssembleRHSElementVect (const FiniteElement &el, ElementTransformation &Tr, Vector &elvect) override |
| void | AssembleRHSElementVect (const FiniteElement &el, FaceElementTransformations &Tr, Vector &elvect) override |
| void | AssembleRHSElementVect (const FiniteElement &el1, const FiniteElement &el2, FaceElementTransformations &Tr, Vector &elvect) override |
| bool | GetTrimFlag () const |
Public Member Functions inherited from mfem::LinearFormIntegrator | |
| virtual bool | SupportsDevice () const |
| Method probing for assembly on device. | |
| virtual void | AssembleDevice (const FiniteElementSpace &fes, const Array< int > &markers, Vector &b) |
| Method defining assembly on device. | |
| virtual | ~LinearFormIntegrator () |
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 IntegrationRule * | GetIntRule () const |
| Directly return the IntRule pointer (possibly null) without checking for NURBS patch rules or falling back on a default. | |
| const IntegrationRule * | GetIntegrationRule () const |
| Equivalent to GetIntRule, but retained for backward compatibility with applications. | |
Protected Attributes | |
| ShiftedVectorFunctionCoefficient * | vN |
| ShiftedFunctionCoefficient * | uN |
| VectorCoefficient * | vD |
| Array< int > * | elem_marker |
| int | nterms |
| bool | include_cut_cell |
| int | NEproc |
| int | par_shared_face_count |
| int | ls_cut_marker |
| Vector | shape |
| Vector | nor |
Protected Attributes inherited from mfem::Integrator | |
| const IntegrationRule * | IntRule |
| NURBSMeshRules * | patchRules = nullptr |
Additional Inherited Members | |
Protected Member Functions inherited from mfem::LinearFormIntegrator | |
| LinearFormIntegrator (const IntegrationRule *ir=NULL) | |
Protected Member Functions inherited from mfem::Integrator | |
| const IntegrationRule * | GetIntegrationRule (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 IntegrationRule * | GetIntegrationRule (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 IntegrationRule * | GetDefaultIntegrationRule (const FiniteElement &trial_fe, const FiniteElement &test_fe, const ElementTransformation &trans) const |
| Subclasses should override to choose a default integration rule. | |
LinearFormIntegrator for Neumann boundaries using the shifted boundary method.
\[ (u, w) = \langle \hat{n} \cdot n \, t_n, w \rangle \]
where \(\hat{n}\) is the normal vector at the true boundary, \(n\) is the normal vector at the surrogate boundary, and \(t_n\) is the traction boundary condition. 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. Note that \(t_n\) is evaluated at the true boundary using the distance function and ShiftedFunctionCoefficient, i.e. \(t_n(x_{true}) = t_N(x_{surrogate} + D)\), where \(x_{surrogate}\) is the location of the integration point on the surrogate boundary and \(D\) is the distance vector from the surrogate boundary to the true boundary.
Definition at line 278 of file sbm_solver.hpp.
|
inline |
Definition at line 296 of file sbm_solver.hpp.
|
overridevirtual |
Given a particular Finite Element and a transformation (Tr) computes the element vector, elvect.
Implements mfem::LinearFormIntegrator.
Definition at line 932 of file sbm_solver.cpp.
|
overridevirtual |
Reimplemented from mfem::LinearFormIntegrator.
Definition at line 938 of file sbm_solver.cpp.
|
overridevirtual |
Reimplemented from mfem::LinearFormIntegrator.
Definition at line 944 of file sbm_solver.cpp.
|
inline |
Definition at line 322 of file sbm_solver.hpp.
|
protected |
Definition at line 284 of file sbm_solver.hpp.
|
protected |
Definition at line 287 of file sbm_solver.hpp.
|
protected |
Definition at line 290 of file sbm_solver.hpp.
|
protected |
Definition at line 288 of file sbm_solver.hpp.
|
protected |
Definition at line 293 of file sbm_solver.hpp.
|
protected |
Definition at line 285 of file sbm_solver.hpp.
|
protected |
Definition at line 289 of file sbm_solver.hpp.
|
protected |
Definition at line 293 of file sbm_solver.hpp.
|
protected |
Definition at line 282 of file sbm_solver.hpp.
|
protected |
Definition at line 283 of file sbm_solver.hpp.
|
protected |
Definition at line 281 of file sbm_solver.hpp.