v0.15.0
Loading...
Searching...
No Matches
Files | Classes | Functions
Simple interface

Implementation of simple interface for fast problem set-up. More...

Collaboration diagram for Simple interface:

Files

file  continuity_check_on_skeleton_with_simple_2d_for_h1.cpp
 
file  continuity_check_on_skeleton_with_simple_2d_for_hcurl.cpp
 
file  continuity_check_on_skeleton_with_simple_2d_for_hdiv.cpp
 
file  simple_interface.cpp
 
file  simple_l2_only.cpp
 
file  Simple.cpp
 Implementation of simple interface.
 
file  Simple.hpp
 Header file for simple interface.
 
file  analytical_nonlinear_poisson.cpp
 
file  analytical_poisson.cpp
 
file  analytical_poisson_field_split.cpp
 
file  hello_world.cpp
 

Classes

struct  MoFEM::BcManager
 Boundary condition manager for finite element problem setup. More...
 
struct  MoFEM::BcManager::BCs
 Data structure storing boundary condition markers and attributes. More...
 
struct  MoFEM::Simple
 Simple interface for fast problem set-up. More...
 

Functions

MoFEMErrorCode MoFEM::BcManager::getOptions ()
 Get command line options for BcManager.
 
MoFEMErrorCode MoFEM::BcManager::pushMarkSideDofs (const std::string problem_name, const std::string block_name, const std::string field_name, int bridge_dim, int lo, int hi)
 Mark DOFs on side entities for boundary conditions.
 
MoFEMErrorCode MoFEM::BcManager::removeSideDOFs (const std::string problem_name, const std::string block_name, const std::string field_name, int bridge_dim, int lo, int hi, bool is_distributed_mesh=true)
 Remove DOFs on side entities from problem.
 
MoFEMErrorCode MoFEM::BcManager::removeBlockDOFsOnEntities (const std::string problem_name, const std::string block_name, const std::string field_name, int lo, int hi, bool get_low_dim_ents=true, bool is_distributed_mesh=true)
 Remove DOFs from problem based on block entities.
 
MoFEMErrorCode MoFEM::BcManager::pushMarkDOFsOnEntities (const std::string problem_name, const std::string block_name, const std::string field_name, int lo, int hi, bool get_low_dim_ents=true)
 Mark DOFs on block entities for boundary conditions.
 
template<typename T >
MoFEMErrorCode MoFEM::BcManager::removeBlockDOFsOnEntities (const std::string problem_name, const std::string field_name, bool get_low_dim_ents=true, bool block_name_field_prefix=false, bool is_distributed_mesh=true)
 Remove DOFs based on boundary condition type template.
 
template<typename T >
MoFEMErrorCode MoFEM::BcManager::pushMarkDOFsOnEntities (const std::string problem_name, const std::string field_name, bool get_low_dim_ents=true, bool block_name_field_prefix=false)
 Mark DOFs based on boundary condition type template.
 
template<typename T >
MoFEMErrorCode MoFEM::BcManager::removeBlockDOFsOnEntities (const std::string problem_name, const std::string block_name, const std::string field_name, bool get_low_dim_ents=true, bool is_distributed_mesh=true)
 Remove DOFs with explicit block specification.
 
template<typename T >
MoFEMErrorCode MoFEM::BcManager::pushMarkDOFsOnEntities (const std::string problem_name, const std::string block_name, const std::string field_name, bool get_low_dim_ents=true)
 Mark DOFs with explicit block specification.
 
boost::shared_ptr< BCsMoFEM::BcManager::popMarkDOFsOnEntities (const std::string block_name)
 Retrieve and remove boundary condition data.
 
auto MoFEM::BcManager::getBcStructure (const std::string bc_id)
 Get boundary condition structure by identifier.
 
BcMapByBlockNameMoFEM::BcManager::getBcMapByBlockName ()
 Get the boundary condition map.
 
Range MoFEM::BcManager::getMergedBlocksRange (std::vector< std::regex > bc_regex_vec)
 Merge entity ranges from multiple boundary condition blocks.
 
auto MoFEM::BcManager::getMergedBlocksRange (std::vector< string > bc_names)
 Merge entity ranges by boundary condition names.
 
BcMarkerPtr MoFEM::BcManager::getMergedBlocksMarker (std::vector< std::regex > bc_regex_vec)
 Merge DOF markers from multiple boundary condition blocks.
 
auto MoFEM::BcManager::getMergedBlocksMarker (std::vector< string > bc_names)
 Merge DOF markers by boundary condition names.
 
BcMarkerPtr MoFEM::BcManager::getMergedBlocksMarker (const std::vector< BcMarkerPtr > &boundary_markers_ptr_vec)
 Merge pre-existing marker vectors.
 
auto MoFEM::BcManager::checkBlock (const std::pair< string, boost::shared_ptr< BCs > > &bc, std::regex reg)
 Check if boundary condition matches regular expression.
 
auto MoFEM::BcManager::checkBlock (const std::pair< std::string, boost::shared_ptr< BCs > > &bc, std::string name)
 Check if boundary condition matches name pattern.
 
SmartPetscObj< IS > MoFEM::BcManager::getBlockIS (const std::string block_prefix, const std::string block_name, const std::string field_name, const std::string problem_name, int lo, int hi, SmartPetscObj< IS > is_expand=SmartPetscObj< IS >())
 Create PETSc Index Set for boundary condition block.
 
SmartPetscObj< IS > MoFEM::BcManager::getBlockIS (const std::string problem_name, const std::string block_name, const std::string field_name, int lo, int hi, SmartPetscObj< IS > is_expand=SmartPetscObj< IS >())
 Create PETSc Index Set for boundary condition block (simplified)
 

Detailed Description

Implementation of simple interface for fast problem set-up.

Function Documentation

◆ checkBlock() [1/2]

auto MoFEM::BcManager::checkBlock ( const std::pair< std::string, boost::shared_ptr< BCs > > &  bc,
std::string  name 
)
inline

#include <src/interfaces/BcManager.hpp>

Check if boundary condition matches name pattern.

Tests whether a boundary condition block name contains the specified substring pattern. Automatically constructs a flexible regular expression for convenient name-based filtering.

Parameters
bcmap element containing block name and boundary condition data
namesubstring pattern to search for in the block name
Returns
bool true if the block name contains the pattern

Definition at line 813 of file BcManager.hpp.

814 {
815 auto full_name = std::string("(.*)_") + name + std::string("(.*)");
816 return checkBlock(bc, std::regex(full_name));
817 }
auto checkBlock(const std::pair< string, boost::shared_ptr< BCs > > &bc, std::regex reg)
Check if boundary condition matches regular expression.

◆ checkBlock() [2/2]

auto MoFEM::BcManager::checkBlock ( const std::pair< string, boost::shared_ptr< BCs > > &  bc,
std::regex  reg 
)
inline

#include <src/interfaces/BcManager.hpp>

Check if boundary condition matches regular expression.

Tests whether a boundary condition block name matches the provided regular expression pattern. This is used internally for filtering and selecting boundary condition blocks based on naming patterns.

Parameters
bcmap element containing block name and boundary condition data
regregular expression pattern for matching
Returns
bool true if the block name matches the pattern

Definition at line 795 of file BcManager.hpp.

796 {
797 return std::regex_match(bc.first, reg);
798 }

◆ getBcMapByBlockName()

BcMapByBlockName & MoFEM::BcManager::getBcMapByBlockName ( )
inline

#include <src/interfaces/BcManager.hpp>

Get the boundary condition map.

Provides access to the complete map of boundary condition blocks managed by this BcManager instance. Useful for iteration and bulk operations on multiple boundary conditions.

Returns
BcMapByBlockName& Reference to the boundary condition map
Examples
mofem/tutorials/scl-10/initial_diffusion.cpp, and mofem/tutorials/scl-10/photon_diffusion.cpp.

Definition at line 700 of file BcManager.hpp.

700{ return bcMapByBlockName; }
BcMapByBlockName bcMapByBlockName

◆ getBcStructure()

auto MoFEM::BcManager::getBcStructure ( const std::string  bc_id)
inline

#include <src/interfaces/BcManager.hpp>

Get boundary condition structure by identifier.

Retrieves the boundary condition data structure for a specific block identifier. Provides direct access to the internal boundary condition management system for advanced operations.

Parameters
bc_ididentifier of the boundary condition block
Returns
auto Reference to the boundary condition structure

Definition at line 686 of file BcManager.hpp.

686 {
687 return bcMapByBlockName.at(bc_id);
688 }

◆ getBlockIS() [1/2]

SmartPetscObj< IS > MoFEM::BcManager::getBlockIS ( const std::string  block_prefix,
const std::string  block_name,
const std::string  field_name,
const std::string  problem_name,
int  lo,
int  hi,
SmartPetscObj< IS >  is_expand = SmartPetscObj<IS>() 
)

#include <src/interfaces/BcManager.hpp>

Create PETSc Index Set for boundary condition block.

Generates a PETSc Index Set (IS) containing DOF indices for a specified boundary condition block. The IS can be used for matrix/vector operations, subproblem creation, and solver configuration in PETSc-based computations.

Parameters
block_prefixprefix for block name construction
block_namename of the boundary condition block
field_namename of the field containing the DOFs
problem_namename of the finite element problem
lolowest coefficient (rank) to include
hihighest coefficient (rank) to include
is_expandoptional existing IS to expand with new indices
Returns
SmartPetscObj<IS> Smart pointer to PETSc Index Set
Examples
mofem/tutorials/scl-6/heat_equation.cpp, and mofem/tutorials/scl-7/wave_equation.cpp.

Definition at line 241 of file BcManager.cpp.

245 {
246 Interface &m_field = cOre;
247
248 const std::string bc_id =
249 block_prefix + "_" + field_name + "_" + block_name + "(.*)";
250
251 Range bc_ents;
252 for (auto bc : getBcMapByBlockName()) {
253 if (std::regex_match(bc.first, std::regex(bc_id))) {
254 bc_ents.merge(*(bc.second->getBcEntsPtr()));
255 MOFEM_LOG("BcMngWorld", Sev::verbose)
256 << "Get entities from block and add to IS. Block name " << bc.first;
257 }
258 }
259
260 SmartPetscObj<IS> is_bc;
261 auto get_is = [&]() {
263 CHKERR m_field.getInterface<CommInterface>()->synchroniseEntities(bc_ents);
264 CHKERR m_field.getInterface<ISManager>()->isCreateProblemFieldAndRank(
265 problem_name, ROW, field_name, lo, hi, is_bc, &bc_ents);
266 if (is_expand) {
267 IS is_tmp;
268 CHKERR ISExpand(is_bc, is_expand, &is_tmp);
269 is_bc = SmartPetscObj<IS>(is_tmp);
270 }
271 CHKERR ISSort(is_bc);
273 };
274
275 if (get_is())
276 CHK_THROW_MESSAGE(MOFEM_DATA_INCONSISTENCY, "IS is not created");
277
278 return is_bc;
279}
@ ROW
#define CHK_THROW_MESSAGE(err, msg)
Check and throw MoFEM exception.
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
@ MOFEM_DATA_INCONSISTENCY
Definition definitions.h:31
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define CHKERR
Inline error check.
#define MOFEM_LOG(channel, severity)
Log.
BcMapByBlockName & getBcMapByBlockName()
Get the boundary condition map.
DeprecatedCoreInterface Interface
constexpr auto field_name
MoFEM::Core & cOre

◆ getBlockIS() [2/2]

SmartPetscObj< IS > MoFEM::BcManager::getBlockIS ( const std::string  problem_name,
const std::string  block_name,
const std::string  field_name,
int  lo,
int  hi,
SmartPetscObj< IS >  is_expand = SmartPetscObj<IS>() 
)

#include <src/interfaces/BcManager.hpp>

Create PETSc Index Set for boundary condition block (simplified)

Simplified version that generates a PETSc Index Set for a boundary condition block without requiring a separate block prefix. This is the most commonly used interface for creating IS objects for boundary condition processing.

Parameters
problem_namename of the finite element problem
block_namename of the boundary condition block
field_namename of the field containing the DOFs
lolowest coefficient (rank) to include
hihighest coefficient (rank) to include
is_expandoptional existing IS to expand with new indices
Returns
SmartPetscObj<IS> Smart pointer to PETSc Index Set

Definition at line 281 of file BcManager.cpp.

284 {
285 return getBlockIS(problem_name, block_name, field_name, problem_name, lo, hi,
286 is_expand);
287}
SmartPetscObj< IS > getBlockIS(const std::string block_prefix, const std::string block_name, const std::string field_name, const std::string problem_name, int lo, int hi, SmartPetscObj< IS > is_expand=SmartPetscObj< IS >())
Create PETSc Index Set for boundary condition block.

◆ getMergedBlocksMarker() [1/3]

BcMarkerPtr MoFEM::BcManager::getMergedBlocksMarker ( const std::vector< BcMarkerPtr > &  boundary_markers_ptr_vec)

#include <src/interfaces/BcManager.hpp>

Merge pre-existing marker vectors.

Combines multiple existing DOF marker vectors into a single unified marker vector. This allows for flexible composition of boundary condition constraints from different sources or processing stages.

Parameters
boundary_markers_ptr_vecvector of marker pointers to merge
Returns
BcMarkerPtr Shared pointer to merged marker vector

◆ getMergedBlocksMarker() [2/3]

BcManager::BcMarkerPtr MoFEM::BcManager::getMergedBlocksMarker ( std::vector< std::regex >  bc_regex_vec)

#include <src/interfaces/BcManager.hpp>

Merge DOF markers from multiple boundary condition blocks.

Combines DOF marker vectors from multiple boundary condition blocks that match the provided regular expression patterns. The resulting marker vector represents the union of all DOF constraints from the matching blocks.

Parameters
bc_regex_vecvector of regular expressions for block name matching
Returns
BcMarkerPtr Shared pointer to merged marker vector

Definition at line 210 of file BcManager.cpp.

210 {
211 BcManager::BcMarkerPtr boundary_marker_ptr;
212 if (bcMapByBlockName.size()) {
213 for (auto b : bcMapByBlockName) {
214 for (auto &reg_name : bc_regex_vec) {
215 if (std::regex_match(b.first, reg_name)) {
216 if (!boundary_marker_ptr)
217 boundary_marker_ptr =
218 boost::make_shared<std::vector<char unsigned>>();
219 boundary_marker_ptr->resize(b.second->bcMarkers.size(), 0);
220 for (int i = 0; i != b.second->bcMarkers.size(); ++i) {
221 (*boundary_marker_ptr)[i] |= b.second->bcMarkers[i];
222 }
223 }
224 }
225 }
226 }
227 return boundary_marker_ptr;
228}
FTensor::Index< 'i', SPACE_DIM > i
boost::shared_ptr< std::vector< char unsigned > > BcMarkerPtr
Shared pointer to marker vector.

◆ getMergedBlocksMarker() [3/3]

auto MoFEM::BcManager::getMergedBlocksMarker ( std::vector< string >  bc_names)
inline

#include <src/interfaces/BcManager.hpp>

Merge DOF markers by boundary condition names.

Convenience function that merges DOF marker vectors from boundary condition blocks with names containing the specified substrings. Useful for combining markers from related boundary conditions into a single constraint system.

Parameters
bc_namesvector of boundary condition name patterns to match
Returns
BcMarkerPtr Shared pointer to merged marker vector

Definition at line 760 of file BcManager.hpp.

760 {
761 std::vector<std::regex> reg_vec(bc_names.size());
762 for (int i = 0; i != bc_names.size(); ++i) {
763 auto full_name = std::string("(.*)_") + bc_names[i] + std::string("(.*)");
764 reg_vec[i] = std::regex(full_name);
765 }
766 return getMergedBlocksMarker(reg_vec);
767 }
BcMarkerPtr getMergedBlocksMarker(std::vector< std::regex > bc_regex_vec)
Merge DOF markers from multiple boundary condition blocks.

◆ getMergedBlocksRange() [1/2]

Range MoFEM::BcManager::getMergedBlocksRange ( std::vector< std::regex >  bc_regex_vec)

#include <src/interfaces/BcManager.hpp>

Merge entity ranges from multiple boundary condition blocks.

Combines entity ranges from multiple boundary condition blocks that match the provided regular expression patterns. This is useful for creating composite boundary conditions or applying operations to multiple related boundary condition blocks simultaneously.

Parameters
bc_regex_vecvector of regular expressions for block name matching
Returns
Range Merged range containing entities from all matching blocks

Definition at line 195 of file BcManager.cpp.

195 {
196 Range ents;
197 if (bcMapByBlockName.size()) {
198 for (auto b : bcMapByBlockName) {
199 for (auto &reg_name : bc_regex_vec) {
200 if (std::regex_match(b.first, reg_name)) {
201 ents.merge(b.second->bcEnts);
202 }
203 }
204 }
205 }
206 return ents;
207}

◆ getMergedBlocksRange() [2/2]

auto MoFEM::BcManager::getMergedBlocksRange ( std::vector< string >  bc_names)
inline

#include <src/interfaces/BcManager.hpp>

Merge entity ranges by boundary condition names.

Convenience function that merges entity ranges from boundary condition blocks with names containing the specified substrings. Automatically constructs regular expressions for flexible name matching.

Parameters
bc_namesvector of boundary condition name patterns to match
Returns
Range Merged range containing entities from all matching blocks

Definition at line 727 of file BcManager.hpp.

727 {
728 std::vector<std::regex> reg_vec(bc_names.size());
729 for (int i = 0; i != bc_names.size(); ++i) {
730 auto full_name = std::string("(.*)_") + bc_names[i] + std::string("(.*)");
731 reg_vec[i] = std::regex(full_name);
732 }
733 return getMergedBlocksRange(reg_vec);
734 }
Range getMergedBlocksRange(std::vector< std::regex > bc_regex_vec)
Merge entity ranges from multiple boundary condition blocks.

◆ getOptions()

MoFEMErrorCode MoFEM::BcManager::getOptions ( )

#include <src/interfaces/BcManager.hpp>

Get command line options for BcManager.

Retrieves configuration options from command line arguments that control boundary condition processing behavior and integration settings.

Returns
MoFEMErrorCode Error code indicating success or failure

Definition at line 65 of file BcManager.cpp.

65 {
67 PetscOptionsBegin(PETSC_COMM_WORLD, "", "BcManager options", "none");
68 PetscOptionsEnd();
70}
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...

◆ popMarkDOFsOnEntities()

boost::shared_ptr< BcManager::BCs > MoFEM::BcManager::popMarkDOFsOnEntities ( const std::string  block_name)

#include <src/interfaces/BcManager.hpp>

Retrieve and remove boundary condition data.

Retrieves the boundary condition data structure for a specified block and removes it from the internal management system. This is typically used when transferring boundary condition data to other processing stages or when cleaning up after boundary condition application.

Parameters
block_namename of the boundary condition block to retrieve
Returns
boost::shared_ptr<BCs> Shared pointer to boundary condition data structure

Definition at line 185 of file BcManager.cpp.

185 {
186 auto bc_it = bcMapByBlockName.find(block_name);
187 if (bc_it != bcMapByBlockName.end()) {
188 auto bc = bc_it->second;
189 bcMapByBlockName.erase(bc_it);
190 return bc;
191 }
192 return boost::shared_ptr<BCs>();
193}

◆ pushMarkDOFsOnEntities() [1/3]

template<typename T >
MoFEMErrorCode MoFEM::BcManager::pushMarkDOFsOnEntities ( const std::string  problem_name,
const std::string  block_name,
const std::string  field_name,
bool  get_low_dim_ents = true 
)

#include <src/interfaces/BcManager.hpp>

Mark DOFs with explicit block specification.

Template-based DOF marking with explicit block name specification. Creates internal data structures for the specified boundary condition block while maintaining type-safe boundary condition handling.

Template Parameters
Tboundary condition type template parameter
Parameters
problem_namename of the finite element problem
block_namespecific name of the boundary condition block
field_namename of the field containing the DOFs
get_low_dim_entsinclude lower dimensional entities
Returns
MoFEMErrorCode Error code indicating success or failure

◆ pushMarkDOFsOnEntities() [2/3]

MoFEMErrorCode MoFEM::BcManager::pushMarkDOFsOnEntities ( const std::string  problem_name,
const std::string  block_name,
const std::string  field_name,
int  lo,
int  hi,
bool  get_low_dim_ents = true 
)

#include <src/interfaces/BcManager.hpp>

Mark DOFs on block entities for boundary conditions.

Marks DOFs associated with entities in a specified block for boundary condition processing. This creates internal data structures needed for subsequent boundary condition application or DOF removal operations.

Parameters
problem_namename of the finite element problem
block_namename of the boundary condition block/meshset
field_namename of the field containing the DOFs to mark
lolowest coefficient (rank) to mark
hihighest coefficient (rank) to mark
get_low_dim_entsinclude lower dimensional entities in the operation
Returns
MoFEMErrorCode Error code indicating success or failure
Examples
mofem/users_modules/eshelbian_plasticity/src/impl/EshelbianPlasticity.cpp.

Definition at line 110 of file BcManager.cpp.

114 {
115 Interface &m_field = cOre;
116 auto prb_mng = m_field.getInterface<ProblemsManager>();
118
119 // if(problem_name.size())
120 // MOFEM_LOG("BcMngWorld", Sev::warning)
121 // << "Argument problem_name has no effect";
122
123 auto fix_disp = [&]() {
125
126 auto mark_fix_dofs = [&](std::vector<unsigned char> &marked_field_dofs,
127 const auto lo, const auto hi) {
128 return prb_mng->modifyMarkDofs(problem_name, ROW, field_name, lo, hi,
130 marked_field_dofs);
131 };
132
133 auto iterate_meshsets = [&](auto &&meshset_vec_ptr) {
135 for (auto m : meshset_vec_ptr) {
136 auto bc = boost::make_shared<BCs>();
137 CHKERR m_field.get_moab().get_entities_by_handle(m->getMeshset(),
138 bc->bcEnts, true);
139 CHKERR m->getAttributes(bc->bcAttributes);
140
141 if (problem_name.size())
142 CHKERR mark_fix_dofs(bc->bcMarkers, lo, hi);
143 MOFEM_LOG("BcMngWorld", Sev::verbose)
144 << "Found block " << m->getName() << " number of attributes "
145 << bc->bcAttributes.size();
146
147 if (get_low_dim_ents) {
148 auto low_dim_ents =
149 BcManagerImplTools::get_adj_ents(m_field.get_moab(), bc->bcEnts);
150 bc->bcEnts.swap(low_dim_ents);
151 }
152
153 CHKERR m_field.getInterface<CommInterface>()->synchroniseEntities(
154 bc->bcEnts);
155 if (problem_name.size())
156 CHKERR prb_mng->markDofs(problem_name, ROW, ProblemsManager::AND,
157 bc->bcEnts, bc->bcMarkers);
158
159 const std::string bc_id =
160 problem_name + "_" + field_name + "_" + m->getName();
161 bcMapByBlockName[bc_id] = bc;
162 }
164 };
165
166 CHKERR iterate_meshsets(
167
168 m_field.getInterface<MeshsetsManager>()->getCubitMeshsetPtr(std::regex(
169
170 (boost::format("%s(.*)") % block_name).str()
171
172 ))
173
174 );
175
177 };
178
179 CHKERR fix_disp();
180
182}
auto get_adj_ents(moab::Interface &moab, const Range &ents)
Definition BcManager.cpp:17
FTensor::Index< 'm', 3 > m
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.

◆ pushMarkDOFsOnEntities() [3/3]

template<typename T >
MoFEMErrorCode MoFEM::BcManager::pushMarkDOFsOnEntities ( const std::string  problem_name,
const std::string  field_name,
bool  get_low_dim_ents = true,
bool  block_name_field_prefix = false 
)

#include <src/interfaces/BcManager.hpp>

Mark DOFs based on boundary condition type template.

Template-based DOF marking that automatically determines the appropriate meshset type and boundary condition handling. Creates internal data structures for subsequent boundary condition processing operations.

Template Parameters
Tboundary condition type template parameter
Parameters
problem_namename of the finite element problem
field_namename of the field containing the DOFs
get_low_dim_entsinclude lower dimensional entities
block_name_field_prefixuse field name as block name prefix
Returns
MoFEMErrorCode Error code indicating success or failure

◆ pushMarkSideDofs()

MoFEMErrorCode MoFEM::BcManager::pushMarkSideDofs ( const std::string  problem_name,
const std::string  block_name,
const std::string  field_name,
int  bridge_dim,
int  lo,
int  hi 
)

#include <src/interfaces/BcManager.hpp>

Mark DOFs on side entities for boundary conditions.

Marks DOFs associated with side entities (faces, edges) of a specified block for boundary condition application. This is commonly used for surface-based boundary conditions like pressure, traction, or flux conditions.

Parameters
problem_namename of the finite element problem
block_namename of the boundary condition block/meshset
field_namename of the field containing the DOFs
bridge_dimdimension of bridge entities connecting DOFs to sides
lolowest coefficient (rank) to mark
hihighest coefficient (rank) to mark
Returns
MoFEMErrorCode Error code indicating success or failure

Definition at line 1630 of file BcManager.cpp.

1635 {
1636 Interface &m_field = cOre;
1638
1639 if (problem_name.empty())
1641
1642 auto iterate_meshsets = [&](auto &&meshset_vec_ptr) {
1643 auto prb_mng = m_field.getInterface<ProblemsManager>();
1645 for (auto m : meshset_vec_ptr) {
1646 auto bc = boost::make_shared<BCs>();
1647 CHKERR m_field.get_moab().get_entities_by_handle(m->getMeshset(),
1648 bc->bcEnts, true);
1649 CHKERR m->getAttributes(bc->bcAttributes);
1650
1651 bc->dofsViewPtr = boost::make_shared<BCs::DofsView>();
1652
1653 CHKERR m_field.getInterface<CommInterface>()->synchroniseEntities(
1654 bc->bcEnts);
1655 CHKERR prb_mng->getSideDofsOnBrokenSpaceEntities(
1656 *(bc->dofsViewPtr), problem_name, ROW, field_name, bc->bcEnts,
1657 bridge_dim, lo, hi);
1658 CHKERR prb_mng->markDofs(problem_name, ROW, *(bc->dofsViewPtr),
1659 ProblemsManager::OR, bc->bcMarkers);
1660
1661 MOFEM_LOG("BcMngWorld", Sev::inform)
1662 << "Found block " << m->getName() << " number of attributes "
1663 << bc->bcAttributes.size() << " number of entities "
1664 << bc->bcEnts.size();
1665
1666 const std::string bc_id =
1667 problem_name + "_" + field_name + "_" + m->getName();
1668 bcMapByBlockName[bc_id] = bc;
1669 }
1671 };
1672
1673 CHKERR iterate_meshsets(
1674
1675 m_field.getInterface<MeshsetsManager>()->getCubitMeshsetPtr(std::regex(
1676
1677 (boost::format("%s(.*)") % block_name).str()
1678
1679 ))
1680
1681 );
1682
1684}

◆ removeBlockDOFsOnEntities() [1/3]

template<typename T >
MoFEMErrorCode MoFEM::BcManager::removeBlockDOFsOnEntities ( const std::string  problem_name,
const std::string  block_name,
const std::string  field_name,
bool  get_low_dim_ents = true,
bool  is_distributed_mesh = true 
)

#include <src/interfaces/BcManager.hpp>

Remove DOFs with explicit block specification.

Template-based DOF removal with explicit block name specification. Provides fine-grained control over which specific block is processed while maintaining type-safe boundary condition handling.

Template Parameters
Tboundary condition type template parameter
Parameters
problem_namename of the finite element problem
block_namespecific name of the boundary condition block
field_namename of the field containing the DOFs
get_low_dim_entsinclude lower dimensional entities
is_distributed_meshtrue if mesh is distributed across processes
Returns
MoFEMErrorCode Error code indicating success or failure

◆ removeBlockDOFsOnEntities() [2/3]

MoFEMErrorCode MoFEM::BcManager::removeBlockDOFsOnEntities ( const std::string  problem_name,
const std::string  block_name,
const std::string  field_name,
int  lo,
int  hi,
bool  get_low_dim_ents = true,
bool  is_distributed_mesh = true 
)

#include <src/interfaces/BcManager.hpp>

Remove DOFs from problem based on block entities.

Removes DOFs associated with entities in a specified block from the finite element problem. This is the primary method for enforcing Dirichlet boundary conditions by eliminating constrained DOFs from the system.

Parameters
problem_namename of the finite element problem
block_namename of the boundary condition block/meshset
field_namename of the field containing the DOFs to remove
lolowest coefficient (rank) to remove
hihighest coefficient (rank) to remove
get_low_dim_entsinclude lower dimensional entities in the operation
is_distributed_meshtrue if mesh is distributed across processes
Returns
MoFEMErrorCode Error code indicating success or failure
Examples
mofem/atom_tests/schur_test_diag_mat.cpp, mofem/atom_tests/tensor_divergence_operator.cpp, and mofem/tutorials/scl-1/poisson_2d_homogeneous.cpp.

Definition at line 72 of file BcManager.cpp.

75 {
76 Interface &m_field = cOre;
77 auto prb_mng = m_field.getInterface<ProblemsManager>();
79
80 CHKERR pushMarkDOFsOnEntities(problem_name, block_name, field_name, lo, hi,
81 get_low_dim_ents);
82
83 auto remove_dofs_on_ents = [&](const Range &ents, const int lo,
84 const int hi) {
85 if (is_distributed_mesh)
86 return prb_mng->removeDofsOnEntities(problem_name, field_name, ents, lo,
87 hi);
88 else
89 return prb_mng->removeDofsOnEntitiesNotDistributed(
90 problem_name, field_name, ents, lo, hi);
91 };
92
93 for (auto m :
94 m_field.getInterface<MeshsetsManager>()->getCubitMeshsetPtr(std::regex(
95
96 (boost::format("%s(.*)") % block_name).str()
97
98 ))
99
100 ) {
101 const std::string bc_id =
102 problem_name + "_" + field_name + "_" + m->getName();
103 CHKERR remove_dofs_on_ents(bcMapByBlockName.at(bc_id)->bcEnts, lo, hi);
104 bcMapByBlockName.at(bc_id)->bcMarkers = std::vector<unsigned char>();
105 }
106
108}
MoFEMErrorCode pushMarkDOFsOnEntities(const std::string problem_name, const std::string block_name, const std::string field_name, int lo, int hi, bool get_low_dim_ents=true)
Mark DOFs on block entities for boundary conditions.
IFACE getInterface() const
Get interface pointer to pointer of interface.

◆ removeBlockDOFsOnEntities() [3/3]

template<typename T >
MoFEMErrorCode MoFEM::BcManager::removeBlockDOFsOnEntities ( const std::string  problem_name,
const std::string  field_name,
bool  get_low_dim_ents = true,
bool  block_name_field_prefix = false,
bool  is_distributed_mesh = true 
)

#include <src/interfaces/BcManager.hpp>

Remove DOFs based on boundary condition type template.

Template-based DOF removal that automatically determines the appropriate meshset type and boundary condition handling based on the template parameter. Supports various boundary condition types including displacement, temperature, heat flux, force, and pressure conditions.

Template Parameters
Tboundary condition type template parameter
Parameters
problem_namename of the finite element problem
field_namename of the field containing the DOFs
get_low_dim_entsinclude lower dimensional entities
block_name_field_prefixuse field name as block name prefix
is_distributed_meshtrue if mesh is distributed across processes
Returns
MoFEMErrorCode Error code indicating success or failure

◆ removeSideDOFs()

MoFEMErrorCode MoFEM::BcManager::removeSideDOFs ( const std::string  problem_name,
const std::string  block_name,
const std::string  field_name,
int  bridge_dim,
int  lo,
int  hi,
bool  is_distributed_mesh = true 
)

#include <src/interfaces/BcManager.hpp>

Remove DOFs on side entities from problem.

Removes DOFs associated with side entities from the finite element problem. This is typically used to enforce boundary conditions by eliminating constrained DOFs from the solution system.

Parameters
problem_namename of the finite element problem
block_namename of the boundary condition block/meshset
field_namename of the field containing the DOFs
bridge_dimdimension of bridge entities connecting DOFs to sides
lolowest coefficient (rank) to remove
hihighest coefficient (rank) to remove
is_distributed_meshtrue if mesh is distributed across processes
Returns
MoFEMErrorCode Error code indicating success or failure
Examples
mofem/atom_tests/test_broken_space.cpp, and mofem/tutorials/fun-2/plot_base.cpp.

Definition at line 1686 of file BcManager.cpp.

1690 {
1691 Interface &m_field = cOre;
1693
1694 CHKERR pushMarkSideDofs(problem_name, block_name, field_name, bridge_dim, lo,
1695 hi);
1696
1697 auto iterate_meshsets = [&](auto &&meshset_vec_ptr) {
1699 auto prb_mng = m_field.getInterface<ProblemsManager>();
1700 for (auto m : meshset_vec_ptr) {
1701 const std::string bc_id =
1702 problem_name + "_" + field_name + "_" + m->getName();
1703 auto &bc = bcMapByBlockName.at(bc_id);
1704 CHKERR prb_mng->removeDofs(problem_name, ROW, *(bc->dofsViewPtr), lo, hi);
1705 CHKERR prb_mng->removeDofs(problem_name, COL, *(bc->dofsViewPtr), lo, hi);
1706 }
1708 };
1709
1710 if (is_distributed_mesh) {
1711
1712 CHKERR iterate_meshsets(
1713
1714 m_field.getInterface<MeshsetsManager>()->getCubitMeshsetPtr(std::regex(
1715
1716 (boost::format("%s(.*)") % block_name).str()
1717
1718 ))
1719
1720 );
1721 } else {
1722 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "Not implemented");
1723 }
1724
1726}
@ COL
MoFEMErrorCode pushMarkSideDofs(const std::string problem_name, const std::string block_name, const std::string field_name, int bridge_dim, int lo, int hi)
Mark DOFs on side entities for boundary conditions.