v0.15.0
Loading...
Searching...
No Matches
Public Member Functions | Public Attributes | List of all members
ApplyDirichletBc Struct Reference
Inheritance diagram for ApplyDirichletBc:
[legend]
Collaboration diagram for ApplyDirichletBc:
[legend]

Public Member Functions

 ApplyDirichletBc (const Range &fix_faces, const Range &fix_nodes, const Range &fix_second_node)
 
MoFEMErrorCode postProcess ()
 Post-processing function executed at loop completion.
 
- Public Member Functions inherited from MoFEM::FEMethod
MoFEMErrorCode query_interface (boost::typeindex::type_index type_index, UnknownInterface **iface) const
 Query interface for type casting.
 
 FEMethod ()=default
 Default constructor.
 
auto getFEName () const
 Get the name of the current finite element.
 
auto getDataDofsPtr () const
 Get pointer to DOF data for the current finite element.
 
auto getDataVectorDofsPtr () const
 Get pointer to vector DOF data for the current finite element.
 
const FieldEntity_vector_viewgetDataFieldEnts () const
 Get reference to data field entities for the current finite element.
 
boost::shared_ptr< FieldEntity_vector_view > & getDataFieldEntsPtr () const
 Get shared pointer to data field entities for the current finite element.
 
auto & getRowFieldEnts () const
 Get reference to row field entities for the current finite element.
 
auto & getRowFieldEntsPtr () const
 Get shared pointer to row field entities for the current finite element.
 
auto & getColFieldEnts () const
 Get reference to column field entities for the current finite element.
 
auto & getColFieldEntsPtr () const
 Get shared pointer to column field entities for the current finite element.
 
auto getRowDofsPtr () const
 Get pointer to row DOFs for the current finite element.
 
auto getColDofsPtr () const
 Get pointer to column DOFs for the current finite element.
 
auto getNumberOfNodes () const
 Get the number of nodes in the current finite element.
 
EntityHandle getFEEntityHandle () const
 Get the entity handle of the current finite element.
 
MoFEMErrorCode getNodeData (const std::string field_name, VectorDouble &data, const bool reset_dofs=true)
 Get nodal data for a specific field.
 
- Public Member Functions inherited from MoFEM::BasicMethod
 BasicMethod ()
 Default constructor.
 
virtual ~BasicMethod ()=default
 Virtual destructor.
 
int getNinTheLoop () const
 Get current loop iteration index.
 
int getLoopSize () const
 Get total loop size.
 
auto getLoHiFERank () const
 Get processor rank range for finite element iteration.
 
auto getLoFERank () const
 Get lower processor rank for finite element iteration.
 
auto getHiFERank () const
 Get upper processor rank for finite element iteration.
 
unsigned int getFieldBitNumber (std::string field_name) const
 Get bit number for a specific field by name.
 
MoFEMErrorCode copyBasicMethod (const BasicMethod &basic)
 Copy data from another BasicMethod instance.
 
virtual MoFEMErrorCode preProcess ()
 Pre-processing function executed at loop initialization.
 
virtual MoFEMErrorCode operator() ()
 Main operator function executed for each loop iteration.
 
boost::weak_ptr< CacheTuplegetCacheWeakPtr () const
 Get the cache weak pointer object.
 
- Public Member Functions inherited from MoFEM::KspMethod
MoFEMErrorCode query_interface (boost::typeindex::type_index type_index, UnknownInterface **iface) const
 Query interface for type casting.
 
 KspMethod ()
 Default constructor.
 
virtual ~KspMethod ()=default
 Virtual destructor.
 
MoFEMErrorCode copyKsp (const KspMethod &ksp)
 Copy data from another KSP method.
 
- Public Member Functions inherited from MoFEM::PetscData
 PetscData ()
 Default constructor.
 
virtual ~PetscData ()=default
 Virtual destructor.
 
MoFEMErrorCode copyPetscData (const PetscData &petsc_data)
 Copy PETSc data from another instance.
 
- Public Member Functions inherited from MoFEM::UnknownInterface
template<class IFACE >
MoFEMErrorCode registerInterface (bool error_if_registration_failed=true)
 Register interface.
 
template<class IFACE >
MoFEMErrorCode getInterface (IFACE *&iface) const
 Get interface reference to pointer of interface.
 
template<class IFACE >
MoFEMErrorCode getInterface (IFACE **const iface) const
 Get interface pointer to pointer of interface.
 
template<class IFACE , typename boost::enable_if< boost::is_pointer< IFACE >, int >::type = 0>
IFACE getInterface () const
 Get interface pointer to pointer of interface.
 
template<class IFACE , typename boost::enable_if< boost::is_reference< IFACE >, int >::type = 0>
IFACE getInterface () const
 Get reference to interface.
 
template<class IFACE >
IFACE * getInterface () const
 Function returning pointer to interface.
 
virtual ~UnknownInterface ()=default
 
- Public Member Functions inherited from MoFEM::SnesMethod
MoFEMErrorCode query_interface (boost::typeindex::type_index type_index, UnknownInterface **iface) const
 Query interface for type casting.
 
 SnesMethod ()
 Default constructor.
 
virtual ~SnesMethod ()=default
 Virtual destructor.
 
MoFEMErrorCode copySnes (const SnesMethod &snes)
 Copy SNES data from another instance.
 
- Public Member Functions inherited from MoFEM::TSMethod
MoFEMErrorCode query_interface (boost::typeindex::type_index type_index, UnknownInterface **iface) const
 Query interface for type casting.
 
 TSMethod ()
 Default constructor.
 
virtual ~TSMethod ()=default
 Virtual destructor.
 
MoFEMErrorCode copyTs (const TSMethod &ts)
 Copy TS solver data from another instance.
 
- Public Member Functions inherited from MoFEM::TaoMethod
MoFEMErrorCode query_interface (boost::typeindex::type_index type_index, UnknownInterface **iface) const
 Query interface for type casting.
 
 TaoMethod ()
 Default constructor.
 
virtual ~TaoMethod ()=default
 Virtual destructor.
 
MoFEMErrorCode copyTao (const TaoMethod &tao)
 Copy TAO data from another instance.
 

Public Attributes

Range fixFaces
 
Range fixNodes
 
Range fixSecondNode
 
- Public Attributes inherited from MoFEM::FEMethod
std::string feName
 Name of the finite element being processed.
 
boost::shared_ptr< const NumeredEntFiniteElementnumeredEntFiniteElementPtr
 Shared pointer to finite element database structure.
 
boost::function< bool(FEMethod *fe_method_ptr)> exeTestHook
 Test function to determine if element should be skipped.
 
- Public Attributes inherited from MoFEM::BasicMethod
int nInTheLoop
 Current index of processed method in the loop.
 
int loopSize
 Total number of methods to process in the loop.
 
std::pair< int, int > loHiFERank
 Processor rank range for distributed finite element iteration.
 
int rAnk
 Current processor rank in MPI communicator.
 
int sIze
 Total number of processors in MPI communicator.
 
const RefEntity_multiIndexrefinedEntitiesPtr
 Pointer to container of refined MoFEM DOF entities.
 
const RefElement_multiIndexrefinedFiniteElementsPtr
 Pointer to container of refined finite element entities.
 
const ProblemproblemPtr
 Raw pointer to current MoFEM problem instance.
 
const Field_multiIndexfieldsPtr
 Raw pointer to fields multi-index container.
 
const FieldEntity_multiIndexentitiesPtr
 Raw pointer to container of field entities.
 
const DofEntity_multiIndexdofsPtr
 Raw pointer to container of degree of freedom entities.
 
const FiniteElement_multiIndexfiniteElementsPtr
 Raw pointer to container of finite elements.
 
const EntFiniteElement_multiIndexfiniteElementsEntitiesPtr
 Raw pointer to container of finite element entities.
 
const FieldEntityEntFiniteElementAdjacencyMap_multiIndexadjacenciesPtr
 Raw pointer to container of adjacencies between DOFs and finite elements.
 
boost::function< MoFEMErrorCode()> preProcessHook
 Hook function for pre-processing operations.
 
boost::function< MoFEMErrorCode()> postProcessHook
 Hook function for post-processing operations.
 
boost::function< MoFEMErrorCode()> operatorHook
 Hook function for main operator execution.
 
boost::movelib::unique_ptr< boolvecAssembleSwitch
 Switch for vector assembly operations.
 
boost::movelib::unique_ptr< boolmatAssembleSwitch
 Switch for matrix assembly operations.
 
boost::weak_ptr< CacheTuplecacheWeakPtr
 Weak pointer to cached entity data.
 
- Public Attributes inherited from MoFEM::KspMethod
KSPContext ksp_ctx
 Current KSP computation context.
 
KSP ksp
 PETSc KSP linear solver object.
 
Vec & ksp_f
 Reference to residual vector in KSP context.
 
Mat & ksp_A
 Reference to system matrix in KSP context.
 
Mat & ksp_B
 Reference to preconditioner matrix in KSP context.
 
- Public Attributes inherited from MoFEM::PetscData
Switches data_ctx
 Current data context switches.
 
Vec f
 PETSc residual vector.
 
Mat A
 PETSc Jacobian matrix.
 
Mat B
 PETSc preconditioner matrix.
 
Vec x
 PETSc solution vector.
 
Vec dx
 PETSc solution increment vector.
 
Vec x_t
 PETSc first time derivative vector.
 
Vec x_tt
 PETSc second time derivative vector.
 
- Public Attributes inherited from MoFEM::SnesMethod
SNESContext snes_ctx
 Current SNES computation context.
 
SNES snes
 PETSc SNES nonlinear solver object.
 
Vec & snes_x
 Reference to current solution state vector.
 
Vec & snes_dx
 Reference to solution update/increment vector.
 
Vec & snes_f
 Reference to residual vector.
 
Mat & snes_A
 Reference to Jacobian matrix.
 
Mat & snes_B
 Reference to preconditioner of Jacobian matrix.
 
- Public Attributes inherited from MoFEM::TSMethod
TS ts
 PETSc time stepping solver object.
 
TSContext ts_ctx
 Current TS computation context.
 
PetscInt ts_step
 Current time step number.
 
PetscReal ts_a
 Shift parameter for U_t (see PETSc Time Solver documentation)
 
PetscReal ts_aa
 Shift parameter for U_tt (second time derivative)
 
PetscReal ts_t
 Current time value.
 
PetscReal ts_dt
 Current time step size.
 
Vec & ts_u
 Reference to current state vector U(t)
 
Vec & ts_u_t
 Reference to first time derivative of state vector dU/dt.
 
Vec & ts_u_tt
 Reference to second time derivative of state vector d²U/dt²
 
Vec & ts_F
 Reference to residual vector F(t,U,U_t,U_tt)
 
Mat & ts_A
 Reference to Jacobian matrix: dF/dU + a*dF/dU_t + aa*dF/dU_tt.
 
Mat & ts_B
 Reference to preconditioner matrix for ts_A.
 
- Public Attributes inherited from MoFEM::TaoMethod
TAOContext tao_ctx
 Current TAO computation context.
 
Tao tao
 PETSc TAO optimization solver object.
 
Vec & tao_x
 Reference to optimization variables vector.
 
Vec & tao_f
 Reference to gradient vector.
 
Mat & tao_A
 Reference to Hessian matrix.
 
Mat & tao_B
 Reference to preconditioner matrix for Hessian.
 

Additional Inherited Members

- Public Types inherited from MoFEM::KspMethod
enum  KSPContext { CTX_SETFUNCTION , CTX_OPERATORS , CTX_KSPNONE }
 Context enumeration for KSP solver phases. More...
 
- Public Types inherited from MoFEM::PetscData
enum  DataContext {
  CTX_SET_NONE = 0 , CTX_SET_F = 1 << 0 , CTX_SET_A = 1 << 1 , CTX_SET_B = 1 << 2 ,
  CTX_SET_X = 1 << 3 , CTX_SET_DX = 1 << 4 , CTX_SET_X_T = 1 << 5 , CTX_SET_X_TT = 1 << 6 ,
  CTX_SET_TIME = 1 << 7
}
 Enumeration for data context flags. More...
 
using Switches = std::bitset< 8 >
 Bitset type for context switches.
 
- Public Types inherited from MoFEM::SnesMethod
enum  SNESContext { CTX_SNESSETFUNCTION , CTX_SNESSETJACOBIAN , CTX_SNESNONE }
 Context enumeration for SNES solver phases. More...
 
- Public Types inherited from MoFEM::TSMethod
enum  TSContext {
  CTX_TSSETRHSFUNCTION , CTX_TSSETRHSJACOBIAN , CTX_TSSETIFUNCTION , CTX_TSSETIJACOBIAN ,
  CTX_TSTSMONITORSET , CTX_TSNONE
}
 Context enumeration for TS solver phases. More...
 
- Public Types inherited from MoFEM::TaoMethod
enum  TAOContext { CTX_TAO_OBJECTIVE , CTX_TAO_GRADIENT , CTX_TAO_HESSIAN , CTX_TAO_NONE }
 Context enumeration for TAO solver phases. More...
 
- Static Public Member Functions inherited from MoFEM::UnknownInterface
static MoFEMErrorCode getLibVersion (Version &version)
 Get library version.
 
static MoFEMErrorCode getFileVersion (moab::Interface &moab, Version &version)
 Get database major version.
 
static MoFEMErrorCode setFileVersion (moab::Interface &moab, Version version=Version(MoFEM_VERSION_MAJOR, MoFEM_VERSION_MINOR, MoFEM_VERSION_BUILD))
 Get database major version.
 
static MoFEMErrorCode getInterfaceVersion (Version &version)
 Get database major version.
 
- Static Public Attributes inherited from MoFEM::PetscData
static constexpr Switches CtxSetNone = PetscData::Switches(CTX_SET_NONE)
 No data switch.
 
static constexpr Switches CtxSetF = PetscData::Switches(CTX_SET_F)
 Residual vector switch.
 
static constexpr Switches CtxSetA = PetscData::Switches(CTX_SET_A)
 Jacobian matrix switch.
 
static constexpr Switches CtxSetB = PetscData::Switches(CTX_SET_B)
 Preconditioner matrix switch.
 
static constexpr Switches CtxSetX = PetscData::Switches(CTX_SET_X)
 Solution vector switch.
 
static constexpr Switches CtxSetDX = PetscData::Switches(CTX_SET_DX)
 Solution increment switch.
 
static constexpr Switches CtxSetX_T = PetscData::Switches(CTX_SET_X_T)
 First time derivative switch.
 
static constexpr Switches CtxSetX_TT = PetscData::Switches(CTX_SET_X_TT)
 Second time derivative switch.
 
static constexpr Switches CtxSetTime = PetscData::Switches(CTX_SET_TIME)
 Time value switch.
 

Detailed Description

Examples
mofem/tutorials/cor-6/simple_elasticity.cpp.

Definition at line 301 of file simple_elasticity.cpp.

Constructor & Destructor Documentation

◆ ApplyDirichletBc()

ApplyDirichletBc::ApplyDirichletBc ( const Range fix_faces,
const Range fix_nodes,
const Range fix_second_node 
)
inline

Definition at line 305 of file simple_elasticity.cpp.

307 : MoFEM::FEMethod(), fixFaces(fix_faces), fixNodes(fix_nodes),
308 fixSecondNode(fix_second_node) {
309 // constructor
310 }
Structure for user loop methods on finite elements.

Member Function Documentation

◆ postProcess()

MoFEMErrorCode ApplyDirichletBc::postProcess ( )
inlinevirtual

Post-processing function executed at loop completion.

This virtual function is called once at the end of the loop. It is typically used for:

  • Assembling matrices and vectors to global system
  • Computing global variables (e.g., total internal energy)
  • Finalizing computation results
  • Cleanup operations

Example of iterating over DOFs:

for(_IT_GET_FEROW_BY_NAME_DOFS_FOR_LOOP_(this,"DISPLACEMENT",it)) {
// Process each DOF
}
Returns
Error code indicating success or failure

Reimplemented from MoFEM::BasicMethod.

Examples
mofem/tutorials/cor-6/simple_elasticity.cpp.

Definition at line 312 of file simple_elasticity.cpp.

312 {
313
315 std::set<int> set_fix_dofs;
316
318 if (dit->get()->getDofCoeffIdx() == 2) {
319 if (fixFaces.find(dit->get()->getEnt()) != fixFaces.end()) {
320 set_fix_dofs.insert(dit->get()->getPetscGlobalDofIdx());
321 }
322 }
323
324 if (fixSecondNode.find(dit->get()->getEnt()) != fixSecondNode.end()) {
325 if (dit->get()->getDofCoeffIdx() == 1) {
326 set_fix_dofs.insert(dit->get()->getPetscGlobalDofIdx());
327 }
328 }
329
330 if (fixNodes.find(dit->get()->getEnt()) != fixNodes.end()) {
331 set_fix_dofs.insert(dit->get()->getPetscGlobalDofIdx());
332 }
333 }
334
335 std::vector<int> fix_dofs(set_fix_dofs.size());
336
337 std::copy(set_fix_dofs.begin(), set_fix_dofs.end(), fix_dofs.begin());
338
339 CHKERR MatAssemblyBegin(ksp_B, MAT_FINAL_ASSEMBLY);
340 CHKERR MatAssemblyEnd(ksp_B, MAT_FINAL_ASSEMBLY);
341 CHKERR VecAssemblyBegin(ksp_f);
342 CHKERR VecAssemblyEnd(ksp_f);
343
344 Vec x;
345
346 CHKERR VecDuplicate(ksp_f, &x);
347 CHKERR VecZeroEntries(x);
348 CHKERR MatZeroRowsColumns(ksp_B, fix_dofs.size(), &fix_dofs[0], 1, x,
349 ksp_f);
350
351 CHKERR VecDestroy(&x);
352
354 }
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define CHKERR
Inline error check.
#define _IT_NUMEREDDOF_ROW_FOR_LOOP_(PROBLEMPTR, IT)
use with loops to iterate row DOFs
const FTensor::Tensor2< T, Dim, Dim > Vec
const Problem * problemPtr
Raw pointer to current MoFEM problem instance.
Mat & ksp_B
Reference to preconditioner matrix in KSP context.
Vec & ksp_f
Reference to residual vector in KSP context.
Vec x
PETSc solution vector.

Member Data Documentation

◆ fixFaces

Range ApplyDirichletBc::fixFaces

◆ fixNodes

Range ApplyDirichletBc::fixNodes

◆ fixSecondNode

Range ApplyDirichletBc::fixSecondNode

The documentation for this struct was generated from the following file: