|
| v0.14.0
|
keeps basic data about problem
More...
#include <src/multi_indices/ProblemsMultiIndices.hpp>
|
auto & | getNumeredRowDofsPtr () const |
| get access to numeredRowDofsPtr storing DOFs on rows More...
|
|
auto & | getNumeredColDofsPtr () const |
| get access to numeredColDofsPtr storing DOFs on cols More...
|
|
const auto & | getNumeredFiniteElementsPtr () const |
| get access to reference for multi-index storing finite elements More...
|
|
MoFEMErrorCode | eraseElements (Range entities) const |
| Erase elements by entities. More...
|
|
boost::shared_ptr< SubProblemData > & | getSubData () const |
| Get main problem of sub-problem is. More...
|
|
auto & | getComposedProblemsData () const |
| Het composed problems data structure. More...
|
|
MoFEMErrorCode | getDofByNameEntAndEntDofIdx (const int field_bit_number, const EntityHandle ent, const int ent_dof_idx, const RowColData row_or_col, boost::shared_ptr< NumeredDofEntity > &dof_ptr) const |
| get DOFs from problem More...
|
|
NumeredDofEntity_multiIndex::iterator | getNumeredRowDofsBegin () const |
|
NumeredDofEntity_multiIndex::iterator | getNumeredRowDofsEnd () const |
|
NumeredDofEntity_multiIndex::iterator | getNumeredColDofsBegin () const |
|
NumeredDofEntity_multiIndex::iterator | getNumeredColDofsEnd () const |
|
auto | getNumeredRowDofsByLocIdxBegin (const DofIdx locidx) const |
|
auto | getNumeredRowDofsByLocIdxEnd (const DofIdx locidx) const |
|
auto | getNumeredColDofsByLocIdxBegin (const DofIdx locidx) const |
|
auto | getNumeredColDofsByLocIdxEnd (const DofIdx locidx) const |
|
auto | getNumeredRowDofsByEntBegin (const EntityHandle ent) const |
|
auto | getNumeredRowDofsByEntEnd (const EntityHandle ent) const |
|
auto | getNumeredColDofsByEntBegin (const EntityHandle ent) const |
|
auto | getNumeredColDofsByEntEnd (const EntityHandle ent) const |
|
| Problem (Interface &moab, const EntityHandle meshset) |
|
virtual | ~Problem ()=default |
|
BitProblemId | getId () const |
|
auto | getName () const |
|
DofIdx | getNbDofsRow () const |
|
DofIdx | getNbDofsCol () const |
|
DofIdx | getNbLocalDofsRow () const |
|
DofIdx | getNbLocalDofsCol () const |
|
DofIdx | getNbGhostDofsRow () const |
|
DofIdx | getNbGhostDofsCol () const |
|
BitRefLevel | getBitRefLevel () const |
|
BitRefLevel | getBitRefLevelMask () const |
|
MoFEMErrorCode | getRowDofsByPetscGlobalDofIdx (DofIdx idx, const NumeredDofEntity **dof_ptr, MoFEMTypes bh=MF_EXIST) const |
| Get the Row Dofs By Petsc Global Dof Idx object. More...
|
|
MoFEMErrorCode | getColDofsByPetscGlobalDofIdx (DofIdx idx, const NumeredDofEntity **dof_ptr, MoFEMTypes bh=MF_EXIST) const |
| Get the Col Dofs By Petsc Global Dof Idx object. More...
|
|
boost::weak_ptr< NumeredDofEntity > | getRowDofsByPetscGlobalDofIdx (DofIdx idx) const |
| Get the Row Dofs By Petsc Global Dof Idx object. More...
|
|
boost::weak_ptr< NumeredDofEntity > | getColDofsByPetscGlobalDofIdx (DofIdx idx) const |
| Get the Col Dofs By Petsc Global Dof Idx object. More...
|
|
BitFEId | getBitFEId () const |
| Get the BitFEIDs in problem
More...
|
|
MoFEMErrorCode | getNumberOfElementsByNameAndPart (MPI_Comm comm, const std::string name, PetscLayout *layout) const |
| Get number of finite elements by name on processors. More...
|
|
MoFEMErrorCode | getNumberOfElementsByPart (MPI_Comm comm, PetscLayout *layout) const |
| Get number of finite elements on processors. More...
|
|
auto & | getRowDofsSequence () const |
| Get reference to sequence data numbered dof container. More...
|
|
auto & | getColDofsSequence () const |
| Get reference to sequence data numbered dof container. More...
|
|
EmptyFieldBlocks & | getEmptyFieldBlocks () const |
| Get the empty field blocks. More...
|
|
BlockFieldPair & | addFieldToEmptyFieldBlocks (const BlockFieldPair add_fields) const |
| Add fields to the empty field blocks. More...
|
|
keeps basic data about problem
This is low level structure with information about problem, what elements compose problem and what DOFs are on rows and columns.
- Todo:
- fix names following name convention
- Examples
- cell_forces.cpp, field_evaluator.cpp, loop_entities.cpp, lorentz_force.cpp, MagneticElement.hpp, nonlinear_dynamics.cpp, partition_mesh.cpp, poisson_2d_dis_galerkin.cpp, Remodeling.cpp, remove_entities_from_problem.cpp, remove_entities_from_problem_not_partitioned.cpp, and test_cache_on_entities.cpp.
Definition at line 54 of file ProblemsMultiIndices.hpp.
◆ BlockFieldPair
◆ EmptyFieldBlocks
◆ SequenceDofContainer
◆ Problem()
Definition at line 9 of file ProblemsMultiIndices.cpp.
16 MOAB_THROW(moab.tag_get_handle(
"_ProblemId", th_ProblemId));
18 moab.tag_get_by_ptr(th_ProblemId, &
meshset, 1, (
const void **)&
tagId));
20 MOAB_THROW(moab.tag_get_handle(
"_ProblemName", th_ProblemName));
24 MOAB_THROW(moab.tag_get_handle(
"_ProblemFEId", th_ProblemFEId));
28 MOAB_THROW(moab.tag_get_handle(
"_RefBitLevel", th_RefBitLevel));
31 Tag th_RefBitLevel_Mask;
32 MOAB_THROW(moab.tag_get_handle(
"_RefBitLevelMask", th_RefBitLevel_Mask));
◆ ~Problem()
virtual MoFEM::Problem::~Problem |
( |
| ) |
|
|
virtualdefault |
◆ addFieldToEmptyFieldBlocks()
Add fields to the empty field blocks.
Emtpy field blocks is a pair contains IDs of the fields for which matrix has zero entries.
- Parameters
-
- Returns
- EmptyFieldBlocks&
Definition at line 553 of file ProblemsMultiIndices.hpp.
◆ eraseElements()
Erase elements by entities.
- Parameters
-
- Returns
- MoFEMErrorCode
Definition at line 170 of file ProblemsMultiIndices.cpp.
172 for (
auto p = entities.pair_begin(); p != entities.pair_end(); ++p) {
◆ getBitFEId()
BitFEId MoFEM::Problem::getBitFEId |
( |
| ) |
const |
◆ getBitRefLevel()
◆ getBitRefLevelMask()
BitRefLevel MoFEM::Problem::getBitRefLevelMask |
( |
| ) |
const |
|
inline |
◆ getColDofsByPetscGlobalDofIdx() [1/2]
Get the Col Dofs By Petsc Global Dof Idx object.
- Parameters
-
- Returns
- boost::weak_ptr<NumeredDofEntity>
Definition at line 59 of file ProblemsMultiIndices.cpp.
61 boost::weak_ptr<NumeredDofEntity> dof_weak_ptr;
62 NumeredDofEntity_multiIndex::index<PetscGlobalIdx_mi_tag>::type::iterator dit;
◆ getColDofsByPetscGlobalDofIdx() [2/2]
Get the Col Dofs By Petsc Global Dof Idx object.
- Parameters
-
- Returns
- MoFEMErrorCode
Definition at line 83 of file ProblemsMultiIndices.cpp.
87 if (
auto shared_dof_ptr = weak_dof_ptr.lock()) {
88 *dof_ptr = shared_dof_ptr.get();
91 SETERRQ1(PETSC_COMM_SELF,
MOFEM_NOT_FOUND,
"row dof <%d> not found", idx);
◆ getColDofsSequence()
auto& MoFEM::Problem::getColDofsSequence |
( |
| ) |
const |
|
inline |
Get reference to sequence data numbered dof container.
In sequence data container data are physically stored. The purpose of this is to allocate NumeredDofEntity data in bulk, having only one allocation instead each time entity is inserted. That makes code efficient.
The vector in sequence is destroyed if last entity inside that vector is destroyed. All MoFEM::NumeredDofEntity have aliased shared_ptr which points to the vector.
- Returns
- MoFEM::Problem::SequenceDofContainer
Definition at line 526 of file ProblemsMultiIndices.hpp.
◆ getComposedProblemsData()
auto& MoFEM::Problem::getComposedProblemsData |
( |
| ) |
const |
|
inline |
◆ getDofByNameEntAndEntDofIdx()
get DOFs from problem
Note that ent_dof_idx is not coefficient number, is local number of DOFs on the entity. The coefficient number and local index of DOFs or entity are the same on vertices and H1 approximation.
- Parameters
-
field_bit_number | field name field_bit_number = (use m_field.get_field_bit_number(field_name); |
ent | entity handle |
ent_dof_idx | index of DOFs on entity |
row_or_col | ROW or COL |
dof_ptr | shared pointer to DOFs if found |
- Returns
- error code
- Examples
- cell_forces.cpp.
Definition at line 132 of file ProblemsMultiIndices.cpp.
138 switch (row_or_col) {
142 "Row numbered index in problem not allocated");
149 "Col numbered index in problem not allocated");
155 "Only ROW and COL is possible for 3rd argument");
162 if (it != numered_dofs->get<Unique_mi_tag>().end()) {
165 dof_ptr = boost::shared_ptr<NumeredDofEntity>();
◆ getEmptyFieldBlocks()
Get the empty field blocks.
Emtpy field blocks is a pair contains IDs of the fields for which matrix has zero entries.
- Returns
- EmptyFieldBlocks&
Definition at line 539 of file ProblemsMultiIndices.hpp.
◆ getId()
◆ getName()
auto MoFEM::Problem::getName |
( |
| ) |
const |
|
inline |
◆ getNbDofsCol()
DofIdx MoFEM::Problem::getNbDofsCol |
( |
| ) |
const |
|
inline |
◆ getNbDofsRow()
DofIdx MoFEM::Problem::getNbDofsRow |
( |
| ) |
const |
|
inline |
◆ getNbGhostDofsCol()
DofIdx MoFEM::Problem::getNbGhostDofsCol |
( |
| ) |
const |
|
inline |
◆ getNbGhostDofsRow()
DofIdx MoFEM::Problem::getNbGhostDofsRow |
( |
| ) |
const |
|
inline |
◆ getNbLocalDofsCol()
DofIdx MoFEM::Problem::getNbLocalDofsCol |
( |
| ) |
const |
|
inline |
◆ getNbLocalDofsRow()
DofIdx MoFEM::Problem::getNbLocalDofsRow |
( |
| ) |
const |
|
inline |
◆ getNumberOfElementsByNameAndPart()
MoFEMErrorCode MoFEM::Problem::getNumberOfElementsByNameAndPart |
( |
MPI_Comm |
comm, |
|
|
const std::string |
name, |
|
|
PetscLayout * |
layout |
|
) |
| const |
Get number of finite elements by name on processors.
Size retuned IS is equal to size of processors.
What PetscLayout see, http://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/IS/index.html
Example of usage for layout
PetscInt rstart, rend;
int global_size;
- Parameters
-
comm | Communicator |
name | Finite element name |
layout | Get number of elements on each processor |
- Returns
- error code
Definition at line 98 of file ProblemsMultiIndices.cpp.
102 MPI_Comm_size(comm, &size);
103 MPI_Comm_rank(comm, &rank);
104 CHKERR PetscLayoutCreate(comm, layout);
105 CHKERR PetscLayoutSetBlockSize(*layout, 1);
106 const auto &fe_by_name_and_part =
109 nb_elems = fe_by_name_and_part.count(boost::make_tuple(name, rank));
110 CHKERR PetscLayoutSetLocalSize(*layout, nb_elems);
111 CHKERR PetscLayoutSetUp(*layout);
◆ getNumberOfElementsByPart()
MoFEMErrorCode MoFEM::Problem::getNumberOfElementsByPart |
( |
MPI_Comm |
comm, |
|
|
PetscLayout * |
layout |
|
) |
| const |
Get number of finite elements on processors.
Size retuned IS is equal to size of processors.
What PetscLayout see, http://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/IS/index.html
Example of usage for layout
PetscInt rstart, rend;
int global_size;
- Parameters
-
comm | Communicator |
layout | Get number of elements on each processor |
- Returns
- error code
Definition at line 115 of file ProblemsMultiIndices.cpp.
119 MPI_Comm_size(comm, &size);
120 MPI_Comm_rank(comm, &rank);
121 CHKERR PetscLayoutCreate(comm, layout);
122 CHKERR PetscLayoutSetBlockSize(*layout, 1);
126 nb_elems = fe_by_part.count(rank);
127 CHKERR PetscLayoutSetLocalSize(*layout, nb_elems);
128 CHKERR PetscLayoutSetUp(*layout);
◆ getNumeredColDofsBegin()
NumeredDofEntity_multiIndex::iterator MoFEM::Problem::getNumeredColDofsBegin |
( |
| ) |
const |
|
inline |
get begin iterator for numeredColDofsPtr (insted you can use #_IT_NUMEREDDOFMOFEMENTITY_COL_FOR_LOOP_ for loops)
Definition at line 201 of file ProblemsMultiIndices.hpp.
◆ getNumeredColDofsByEntBegin()
auto MoFEM::Problem::getNumeredColDofsByEntBegin |
( |
const EntityHandle |
ent | ) |
const |
|
inline |
◆ getNumeredColDofsByEntEnd()
auto MoFEM::Problem::getNumeredColDofsByEntEnd |
( |
const EntityHandle |
ent | ) |
const |
|
inline |
◆ getNumeredColDofsByLocIdxBegin()
auto MoFEM::Problem::getNumeredColDofsByLocIdxBegin |
( |
const DofIdx |
locidx | ) |
const |
|
inline |
◆ getNumeredColDofsByLocIdxEnd()
auto MoFEM::Problem::getNumeredColDofsByLocIdxEnd |
( |
const DofIdx |
locidx | ) |
const |
|
inline |
◆ getNumeredColDofsEnd()
NumeredDofEntity_multiIndex::iterator MoFEM::Problem::getNumeredColDofsEnd |
( |
| ) |
const |
|
inline |
get end iterator for numeredColDofsPtr (insted you can use #_IT_NUMEREDDOFMOFEMENTITY_COL_FOR_LOOP_ for loops)
Definition at line 207 of file ProblemsMultiIndices.hpp.
◆ getNumeredColDofsPtr()
auto& MoFEM::Problem::getNumeredColDofsPtr |
( |
| ) |
const |
|
inline |
◆ getNumeredFiniteElementsPtr()
const auto& MoFEM::Problem::getNumeredFiniteElementsPtr |
( |
| ) |
const |
|
inline |
◆ getNumeredRowDofsBegin()
NumeredDofEntity_multiIndex::iterator MoFEM::Problem::getNumeredRowDofsBegin |
( |
| ) |
const |
|
inline |
get begin iterator for numeredRowDofsPtr (insted you can use #_IT_NUMEREDDOFMOFEMENTITY_ROW_FOR_LOOP_ for loops)
Definition at line 189 of file ProblemsMultiIndices.hpp.
◆ getNumeredRowDofsByEntBegin()
auto MoFEM::Problem::getNumeredRowDofsByEntBegin |
( |
const EntityHandle |
ent | ) |
const |
|
inline |
◆ getNumeredRowDofsByEntEnd()
auto MoFEM::Problem::getNumeredRowDofsByEntEnd |
( |
const EntityHandle |
ent | ) |
const |
|
inline |
◆ getNumeredRowDofsByLocIdxBegin()
auto MoFEM::Problem::getNumeredRowDofsByLocIdxBegin |
( |
const DofIdx |
locidx | ) |
const |
|
inline |
◆ getNumeredRowDofsByLocIdxEnd()
auto MoFEM::Problem::getNumeredRowDofsByLocIdxEnd |
( |
const DofIdx |
locidx | ) |
const |
|
inline |
◆ getNumeredRowDofsEnd()
NumeredDofEntity_multiIndex::iterator MoFEM::Problem::getNumeredRowDofsEnd |
( |
| ) |
const |
|
inline |
get end iterator for numeredRowDofsPtr (insted you can use #_IT_NUMEREDDOFMOFEMENTITY_ROW_FOR_LOOP_ for loops)
Definition at line 195 of file ProblemsMultiIndices.hpp.
◆ getNumeredRowDofsPtr()
auto& MoFEM::Problem::getNumeredRowDofsPtr |
( |
| ) |
const |
|
inline |
◆ getRowDofsByPetscGlobalDofIdx() [1/2]
Get the Row Dofs By Petsc Global Dof Idx object.
- Parameters
-
- Returns
- boost::weak_ptr<NumeredDofEntity>
Definition at line 48 of file ProblemsMultiIndices.cpp.
50 boost::weak_ptr<NumeredDofEntity> dof_weak_ptr;
51 NumeredDofEntity_multiIndex::index<PetscGlobalIdx_mi_tag>::type::iterator dit;
◆ getRowDofsByPetscGlobalDofIdx() [2/2]
Get the Row Dofs By Petsc Global Dof Idx object.
- Parameters
-
- Returns
- MoFEMErrorCode
Definition at line 69 of file ProblemsMultiIndices.cpp.
73 if (
auto shared_dof_ptr = weak_dof_ptr.lock()) {
74 *dof_ptr = shared_dof_ptr.get();
77 SETERRQ1(PETSC_COMM_SELF,
MOFEM_NOT_FOUND,
"row dof <%d> not found", idx);
◆ getRowDofsSequence()
auto& MoFEM::Problem::getRowDofsSequence |
( |
| ) |
const |
|
inline |
Get reference to sequence data numbered dof container.
In sequence data container data are physically stored. The purpose of this is to allocate NumeredDofEntity data in bulk, having only one allocation instead each time entity is inserted. That makes code efficient.
The vector in sequence is destroyed if last entity inside that vector is destroyed. All MoFEM::NumeredDofEntity have aliased shared_ptr which points to the vector.
- Returns
- MoFEM::Problem::SequenceDofContainer
Definition at line 510 of file ProblemsMultiIndices.hpp.
◆ getSubData()
boost::shared_ptr<SubProblemData>& MoFEM::Problem::getSubData |
( |
| ) |
const |
|
inline |
◆ operator<<
std::ostream& operator<< |
( |
std::ostream & |
os, |
|
|
const Problem & |
e |
|
) |
| |
|
friend |
Definition at line 37 of file ProblemsMultiIndices.cpp.
38 os <<
"Problem id " << e.getId() <<
" name " << e.getName() << endl;
39 os <<
"FiniteElement id " << e.getBitFEId() << endl;
40 os <<
"BitRefLevel " << e.getBitRefLevel() << endl;
41 os <<
"BitRefLevelMask " << e.getBitRefLevelMask();
◆ composedProblemsData
◆ emptyFieldBlocks
◆ meshset
◆ nbDofsCol
DofIdx MoFEM::Problem::nbDofsCol |
|
mutable |
◆ nbDofsRow
DofIdx MoFEM::Problem::nbDofsRow |
|
mutable |
◆ nbGhostDofsCol
DofIdx MoFEM::Problem::nbGhostDofsCol |
|
mutable |
◆ nbGhostDofsRow
DofIdx MoFEM::Problem::nbGhostDofsRow |
|
mutable |
◆ nbLocDofsCol
DofIdx MoFEM::Problem::nbLocDofsCol |
|
mutable |
◆ nbLocDofsRow
DofIdx MoFEM::Problem::nbLocDofsRow |
|
mutable |
◆ numeredColDofsPtr
◆ numeredFiniteElementsPtr
◆ numeredRowDofsPtr
◆ sequenceColDofContainer
◆ sequenceRowDofContainer
◆ subProblemData
Pointer to data structure. This pointer has allocated data only for sub problems.
Definition at line 107 of file ProblemsMultiIndices.hpp.
◆ tagBitFEId
BitFEId* MoFEM::Problem::tagBitFEId |
◆ tagBitRefLevel
◆ tagBitRefLevelMask
◆ tagId
◆ tagName
const char* MoFEM::Problem::tagName |
◆ tagNameSize
int MoFEM::Problem::tagNameSize |
The documentation for this struct was generated from the following files:
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
BitRefLevel * tagBitRefLevel
BitRef level of finite elements in problem.
DofIdx nbGhostDofsRow
Number of ghost DOFs in row.
static UId getUniqueIdCalculate(const DofIdx dof, UId ent_uid)
std::bitset< BITPROBLEMID_SIZE > BitProblemId
Problem Id.
BitFEId * tagBitFEId
IDs of finite elements in problem.
#define MOAB_THROW(err)
Check error code of MoAB function and throw MoFEM exception.
int tagNameSize
Size of problem name.
multi_index_container< boost::shared_ptr< NumeredDofEntity >, indexed_by< ordered_unique< tag< Unique_mi_tag >, const_mem_fun< NumeredDofEntity::interface_type_DofEntity, UId, &NumeredDofEntity::getLocalUniqueId > >, ordered_non_unique< tag< Part_mi_tag >, member< NumeredDofEntity, unsigned int, &NumeredDofEntity::pArt > >, ordered_non_unique< tag< Idx_mi_tag >, member< NumeredDofEntity, DofIdx, &NumeredDofEntity::dofIdx > >, ordered_non_unique< tag< PetscGlobalIdx_mi_tag >, member< NumeredDofEntity, DofIdx, &NumeredDofEntity::petscGloablDofIdx > >, ordered_non_unique< tag< PetscLocalIdx_mi_tag >, member< NumeredDofEntity, DofIdx, &NumeredDofEntity::petscLocalDofIdx > >, ordered_non_unique< tag< Ent_mi_tag >, const_mem_fun< NumeredDofEntity::interface_type_DofEntity, EntityHandle, &NumeredDofEntity::getEnt > > > > NumeredDofEntity_multiIndex
MultiIndex container keeps NumeredDofEntity.
#define CHKERR
Inline error check.
boost::shared_ptr< ComposedProblemsData > composedProblemsData
multi_index_container< boost::weak_ptr< std::vector< NumeredDofEntity > >, indexed_by< sequenced<> > > SequenceDofContainer
boost::shared_ptr< NumeredDofEntity_multiIndex > numeredColDofsPtr
store DOFs on columns for this problem
UId getLocalUniqueIdCalculate()
Get the Local Unique Id Calculate object.
DofIdx nbDofsCol
Global number of DOFs in col.
MoFEMErrorCode getColDofsByPetscGlobalDofIdx(DofIdx idx, const NumeredDofEntity **dof_ptr, MoFEMTypes bh=MF_EXIST) const
Get the Col Dofs By Petsc Global Dof Idx object.
DofIdx nbLocDofsRow
Local number of DOFs in row.
const char * tagName
Problem name.
BitProblemId * tagId
Unique problem ID.
boost::shared_ptr< NumeredEntFiniteElement_multiIndex > numeredFiniteElementsPtr
store finite elements
multi_index_container< boost::shared_ptr< NumeredEntFiniteElement >, indexed_by< ordered_unique< tag< Unique_mi_tag >, const_mem_fun< NumeredEntFiniteElement::interface_type_EntFiniteElement, UId, &NumeredEntFiniteElement::getLocalUniqueId > >, ordered_non_unique< tag< Part_mi_tag >, member< NumeredEntFiniteElement, unsigned int, &NumeredEntFiniteElement::part > >, ordered_non_unique< tag< Ent_mi_tag >, const_mem_fun< NumeredEntFiniteElement::interface_type_RefEntity, EntityHandle, &NumeredEntFiniteElement::getEnt > >, ordered_non_unique< tag< Composite_Name_And_Part_mi_tag >, composite_key< NumeredEntFiniteElement, const_mem_fun< NumeredEntFiniteElement::interface_type_FiniteElement, boost::string_ref, &NumeredEntFiniteElement::getNameRef >, member< NumeredEntFiniteElement, unsigned int, &NumeredEntFiniteElement::part > >> >> NumeredEntFiniteElement_multiIndex
MultiIndex for entities for NumeredEntFiniteElement.
DofIdx nbGhostDofsCol
Number of ghost DOFs in col.
MoFEMErrorCode getRowDofsByPetscGlobalDofIdx(DofIdx idx, const NumeredDofEntity **dof_ptr, MoFEMTypes bh=MF_EXIST) const
Get the Row Dofs By Petsc Global Dof Idx object.
boost::shared_ptr< NumeredDofEntity_multiIndex > numeredRowDofsPtr
store DOFs on rows for this problem
DofIdx nbLocDofsCol
Local number of DOFs in colIs.
boost::shared_ptr< SequenceDofContainer > sequenceRowDofContainer
boost::shared_ptr< SequenceDofContainer > sequenceColDofContainer
static MoFEMErrorCodeGeneric< PetscErrorCode > ierr
DofIdx nbDofsRow
Global number of DOFs in row.
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
EmptyFieldBlocks emptyFieldBlocks
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define CHKERRG(n)
Check error code of MoFEM/MOAB/PETSc function.
BitRefLevel * tagBitRefLevelMask
BItRefMask of elements in problem.
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
boost::shared_ptr< SubProblemData > subProblemData