v0.14.0
Classes | Public Types | Public Member Functions | Public Attributes | Private Attributes | Friends | List of all members
MoFEM::Problem Struct Reference

keeps basic data about problem More...

#include <src/multi_indices/ProblemsMultiIndices.hpp>

Collaboration diagram for MoFEM::Problem:
[legend]

Classes

struct  SubProblemData
 Subproblem problem data. More...
 

Public Types

typedef multi_index_container< boost::weak_ptr< std::vector< NumeredDofEntity > >, indexed_by< sequenced<> > > SequenceDofContainer
 
using BlockFieldPair = std::pair< BitFieldId, BitFieldId >
 
using EmptyFieldBlocks = std::vector< BlockFieldPair >
 

Public Member Functions

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< NumeredDofEntitygetRowDofsByPetscGlobalDofIdx (DofIdx idx) const
 Get the Row Dofs By Petsc Global Dof Idx object. More...
 
boost::weak_ptr< NumeredDofEntitygetColDofsByPetscGlobalDofIdx (DofIdx idx) const
 Get the Col Dofs By Petsc Global Dof Idx object. More...
 
BitFEId getBitFEId () const
 
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...
 
EmptyFieldBlocksgetEmptyFieldBlocks () const
 Get the empty field blocks. More...
 
BlockFieldPairaddFieldToEmptyFieldBlocks (const BlockFieldPair add_fields) const
 Add fields to the empty field blocks. More...
 

Public Attributes

EntityHandle meshset
 
BitProblemIdtagId
 Unique problem ID. More...
 
const char * tagName
 Problem name. More...
 
int tagNameSize
 Size of problem name. More...
 
BitFEIdtagBitFEId
 IDs of finite elements in problem. More...
 
BitRefLeveltagBitRefLevel
 BitRef level of finite elements in problem. More...
 
BitRefLeveltagBitRefLevelMask
 BItRefMask of elements in problem. More...
 
DofIdx nbDofsRow
 Global number of DOFs in row. More...
 
DofIdx nbDofsCol
 Global number of DOFs in col. More...
 
DofIdx nbLocDofsRow
 Local number of DOFs in row. More...
 
DofIdx nbLocDofsCol
 Local number of DOFs in colIs. More...
 
DofIdx nbGhostDofsRow
 Number of ghost DOFs in row. More...
 
DofIdx nbGhostDofsCol
 Number of ghost DOFs in col. More...
 
boost::shared_ptr< NumeredDofEntity_multiIndexnumeredRowDofsPtr
 store DOFs on rows for this problem More...
 
boost::shared_ptr< NumeredDofEntity_multiIndexnumeredColDofsPtr
 store DOFs on columns for this problem More...
 
boost::shared_ptr< NumeredEntFiniteElement_multiIndexnumeredFiniteElementsPtr
 store finite elements More...
 
boost::shared_ptr< SubProblemDatasubProblemData
 
boost::shared_ptr< ComposedProblemsDatacomposedProblemsData
 

Private Attributes

boost::shared_ptr< SequenceDofContainersequenceRowDofContainer
 
boost::shared_ptr< SequenceDofContainersequenceColDofContainer
 
EmptyFieldBlocks emptyFieldBlocks
 

Friends

std::ostream & operator<< (std::ostream &os, const Problem &e)
 

Detailed Description

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.

Member Typedef Documentation

◆ BlockFieldPair

Definition at line 524 of file ProblemsMultiIndices.hpp.

◆ EmptyFieldBlocks

Definition at line 525 of file ProblemsMultiIndices.hpp.

◆ SequenceDofContainer

typedef multi_index_container<boost::weak_ptr<std::vector<NumeredDofEntity> >, indexed_by<sequenced<> > > MoFEM::Problem::SequenceDofContainer

Definition at line 490 of file ProblemsMultiIndices.hpp.

Constructor & Destructor Documentation

◆ Problem()

MoFEM::Problem::Problem ( Interface moab,
const EntityHandle  meshset 
)

Definition at line 9 of file ProblemsMultiIndices.cpp.

15  Tag th_ProblemId;
16  MOAB_THROW(moab.tag_get_handle("_ProblemId", th_ProblemId));
17  MOAB_THROW(
18  moab.tag_get_by_ptr(th_ProblemId, &meshset, 1, (const void **)&tagId));
19  Tag th_ProblemName;
20  MOAB_THROW(moab.tag_get_handle("_ProblemName", th_ProblemName));
21  MOAB_THROW(moab.tag_get_by_ptr(th_ProblemName, &meshset, 1,
22  (const void **)&tagName, &tagNameSize));
23  Tag th_ProblemFEId;
24  MOAB_THROW(moab.tag_get_handle("_ProblemFEId", th_ProblemFEId));
25  MOAB_THROW(moab.tag_get_by_ptr(th_ProblemFEId, &meshset, 1,
26  (const void **)&tagBitFEId));
27  Tag th_RefBitLevel;
28  MOAB_THROW(moab.tag_get_handle("_RefBitLevel", th_RefBitLevel));
29  MOAB_THROW(moab.tag_get_by_ptr(th_RefBitLevel, &meshset, 1,
30  (const void **)&tagBitRefLevel));
31  Tag th_RefBitLevel_Mask;
32  MOAB_THROW(moab.tag_get_handle("_RefBitLevelMask", th_RefBitLevel_Mask));
33  MOAB_THROW(moab.tag_get_by_ptr(th_RefBitLevel_Mask, &meshset, 1,
34  (const void **)&tagBitRefLevelMask));
35 }

◆ ~Problem()

virtual MoFEM::Problem::~Problem ( )
virtualdefault

Member Function Documentation

◆ addFieldToEmptyFieldBlocks()

BlockFieldPair& MoFEM::Problem::addFieldToEmptyFieldBlocks ( const BlockFieldPair  add_fields) const
inline

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
add_fields
Returns
EmptyFieldBlocks&

Definition at line 549 of file ProblemsMultiIndices.hpp.

549  {
550  emptyFieldBlocks.push_back(add_fields);
551  return emptyFieldBlocks.back();
552  }

◆ eraseElements()

MoFEMErrorCode MoFEM::Problem::eraseElements ( Range  entities) const

Erase elements by entities.

Parameters
entities
Returns
MoFEMErrorCode

Definition at line 170 of file ProblemsMultiIndices.cpp.

170  {
172  for (auto p = entities.pair_begin(); p != entities.pair_end(); ++p) {
173  auto lo = numeredFiniteElementsPtr->get<Ent_mi_tag>().lower_bound(p->first);
174  auto hi =
175  numeredFiniteElementsPtr->get<Ent_mi_tag>().upper_bound(p->second);
176  numeredFiniteElementsPtr->get<Ent_mi_tag>().erase(lo, hi);
177  }
179 }

◆ getBitFEId()

BitFEId MoFEM::Problem::getBitFEId ( ) const

Definition at line 45 of file ProblemsMultiIndices.cpp.

45 { return *tagBitFEId; }

◆ getBitRefLevel()

BitRefLevel MoFEM::Problem::getBitRefLevel ( ) const
inline

Definition at line 383 of file ProblemsMultiIndices.hpp.

383 { return *tagBitRefLevel; }

◆ getBitRefLevelMask()

BitRefLevel MoFEM::Problem::getBitRefLevelMask ( ) const
inline

Definition at line 384 of file ProblemsMultiIndices.hpp.

384 { return *tagBitRefLevelMask; }

◆ getColDofsByPetscGlobalDofIdx() [1/2]

boost::weak_ptr< NumeredDofEntity > MoFEM::Problem::getColDofsByPetscGlobalDofIdx ( DofIdx  idx) const

Get the Col Dofs By Petsc Global Dof Idx object.

Parameters
idx
Returns
boost::weak_ptr<NumeredDofEntity>

Definition at line 59 of file ProblemsMultiIndices.cpp.

59  {
61  boost::weak_ptr<NumeredDofEntity> dof_weak_ptr;
62  NumeredDofEntity_multiIndex::index<PetscGlobalIdx_mi_tag>::type::iterator dit;
63  dit = numeredColDofsPtr->get<PetscGlobalIdx_mi_tag>().find(idx);
64  if (dit != numeredColDofsPtr->get<PetscGlobalIdx_mi_tag>().end())
65  dof_weak_ptr = *dit;
66  return dof_weak_ptr;
67 }

◆ getColDofsByPetscGlobalDofIdx() [2/2]

MoFEMErrorCode MoFEM::Problem::getColDofsByPetscGlobalDofIdx ( DofIdx  idx,
const NumeredDofEntity **  dof_ptr,
MoFEMTypes  bh = MF_EXIST 
) const

Get the Col Dofs By Petsc Global Dof Idx object.

Parameters
idx
dof_ptr
bh
Returns
MoFEMErrorCode

Definition at line 83 of file ProblemsMultiIndices.cpp.

84  {
86  auto weak_dof_ptr = getColDofsByPetscGlobalDofIdx(idx);
87  if (auto shared_dof_ptr = weak_dof_ptr.lock()) {
88  *dof_ptr = shared_dof_ptr.get();
89  } else {
90  if (bh == MF_EXIST)
91  SETERRQ1(PETSC_COMM_SELF, MOFEM_NOT_FOUND, "row dof <%d> not found", idx);
92  *dof_ptr = nullptr;
93  }
95 }

◆ 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 522 of file ProblemsMultiIndices.hpp.

522 { return sequenceColDofContainer; }

◆ getComposedProblemsData()

auto& MoFEM::Problem::getComposedProblemsData ( ) const
inline

Het composed problems data structure.

Definition at line 131 of file ProblemsMultiIndices.hpp.

131 { return composedProblemsData; }

◆ getDofByNameEntAndEntDofIdx()

MoFEMErrorCode MoFEM::Problem::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

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_numberfield name field_bit_number = (use m_field.get_field_bit_number(field_name);
ententity handle
ent_dof_idxindex of DOFs on entity
row_or_colROW or COL
dof_ptrshared pointer to DOFs if found
Returns
error code
Examples
cell_forces.cpp.

Definition at line 132 of file ProblemsMultiIndices.cpp.

135  {
137  decltype(numeredRowDofsPtr) numered_dofs;
138  switch (row_or_col) {
139  case ROW:
140  if (!numeredRowDofsPtr) {
141  SETERRQ(PETSC_COMM_SELF, MOFEM_INVALID_DATA,
142  "Row numbered index in problem not allocated");
143  }
144  numered_dofs = numeredRowDofsPtr;
145  break;
146  case COL:
147  if (!numeredColDofsPtr) {
148  SETERRQ(PETSC_COMM_SELF, MOFEM_INVALID_DATA,
149  "Col numbered index in problem not allocated");
150  }
151  numered_dofs = numeredColDofsPtr;
152  break;
153  default:
154  SETERRQ(PETSC_COMM_SELF, MOFEM_INVALID_DATA,
155  "Only ROW and COL is possible for 3rd argument");
156  }
157 
158  auto it =
159  numered_dofs->get<Unique_mi_tag>().find(DofEntity::getUniqueIdCalculate(
160  ent_dof_idx,
161  FieldEntity::getLocalUniqueIdCalculate(field_bit_number, ent)));
162  if (it != numered_dofs->get<Unique_mi_tag>().end()) {
163  dof_ptr = *it;
164  } else {
165  dof_ptr = boost::shared_ptr<NumeredDofEntity>();
166  }
168 }

◆ getEmptyFieldBlocks()

EmptyFieldBlocks& MoFEM::Problem::getEmptyFieldBlocks ( ) const
inline

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 535 of file ProblemsMultiIndices.hpp.

535  {
536  return emptyFieldBlocks;
537  }

◆ getId()

BitProblemId MoFEM::Problem::getId ( ) const
inline

Definition at line 370 of file ProblemsMultiIndices.hpp.

370 { return *((BitProblemId *)tagId); }

◆ getName()

auto MoFEM::Problem::getName ( ) const
inline
Examples
field_evaluator.cpp, lorentz_force.cpp, MagneticElement.hpp, and Remodeling.cpp.

Definition at line 372 of file ProblemsMultiIndices.hpp.

372  {
373  return std::string((char *)tagName, tagNameSize);
374  }

◆ getNbDofsCol()

DofIdx MoFEM::Problem::getNbDofsCol ( ) const
inline

Definition at line 377 of file ProblemsMultiIndices.hpp.

377 { return nbDofsCol; }

◆ getNbDofsRow()

DofIdx MoFEM::Problem::getNbDofsRow ( ) const
inline
Examples
nonlinear_dynamics.cpp, and poisson_2d_dis_galerkin.cpp.

Definition at line 376 of file ProblemsMultiIndices.hpp.

376 { return nbDofsRow; }

◆ getNbGhostDofsCol()

DofIdx MoFEM::Problem::getNbGhostDofsCol ( ) const
inline

Definition at line 381 of file ProblemsMultiIndices.hpp.

381 { return nbGhostDofsCol; }

◆ getNbGhostDofsRow()

DofIdx MoFEM::Problem::getNbGhostDofsRow ( ) const
inline

Definition at line 380 of file ProblemsMultiIndices.hpp.

380 { return nbGhostDofsRow; }

◆ getNbLocalDofsCol()

DofIdx MoFEM::Problem::getNbLocalDofsCol ( ) const
inline
Examples
nonlinear_dynamics.cpp.

Definition at line 379 of file ProblemsMultiIndices.hpp.

379 { return nbLocDofsCol; }

◆ getNbLocalDofsRow()

DofIdx MoFEM::Problem::getNbLocalDofsRow ( ) const
inline
Examples
nonlinear_dynamics.cpp.

Definition at line 378 of file ProblemsMultiIndices.hpp.

378 { return nbLocDofsRow; }

◆ 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;
ierr = PetscLayoutGetRange(layout, &rstart, &rend); CHKERRG(ierr);
int global_size;
ierr = PetscLayoutGetSize(layout,&global_size); CHKERRG(ierr);
Parameters
commCommunicator
nameFinite element name
layoutGet number of elements on each processor
Returns
error code

Definition at line 98 of file ProblemsMultiIndices.cpp.

99  {
101  int size, rank;
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 =
107  numeredFiniteElementsPtr->get<Composite_Name_And_Part_mi_tag>();
108  int nb_elems;
109  nb_elems = fe_by_name_and_part.count(boost::make_tuple(name, rank));
110  CHKERR PetscLayoutSetLocalSize(*layout, nb_elems);
111  CHKERR PetscLayoutSetUp(*layout);
113 }

◆ 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;
ierr = PetscLayoutGetRange(layout, &rstart, &rend); CHKERRG(ierr);
int global_size;
ierr = PetscLayoutGetSize(layout,&global_size); CHKERRG(ierr);
Parameters
commCommunicator
layoutGet number of elements on each processor
Returns
error code

Definition at line 115 of file ProblemsMultiIndices.cpp.

116  {
118  int size, rank;
119  MPI_Comm_size(comm, &size);
120  MPI_Comm_rank(comm, &rank);
121  CHKERR PetscLayoutCreate(comm, layout);
122  CHKERR PetscLayoutSetBlockSize(*layout, 1);
124  const FeByPart &fe_by_part = numeredFiniteElementsPtr->get<Part_mi_tag>();
125  int nb_elems;
126  nb_elems = fe_by_part.count(rank);
127  CHKERR PetscLayoutSetLocalSize(*layout, nb_elems);
128  CHKERR PetscLayoutSetUp(*layout);
130 }

◆ 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.

201  {
202  return numeredColDofsPtr->begin();
203  }

◆ getNumeredColDofsByEntBegin()

auto MoFEM::Problem::getNumeredColDofsByEntBegin ( const EntityHandle  ent) const
inline

get begin iterator for numeredColDofsPtr (insted you can use _IT_NUMEREDDOF_COL_BY_ENT_FOR_LOOP_ for loops)

Definition at line 316 of file ProblemsMultiIndices.hpp.

316  {
317  return numeredColDofsPtr->get<Ent_mi_tag>().lower_bound(ent);
318  }

◆ getNumeredColDofsByEntEnd()

auto MoFEM::Problem::getNumeredColDofsByEntEnd ( const EntityHandle  ent) const
inline

get end iterator for numeredColDofsPtr (insted you can use _IT_NUMEREDDOF_COL_BY_ENT_FOR_LOOP_ for loops)

Definition at line 322 of file ProblemsMultiIndices.hpp.

322  {
323  return numeredColDofsPtr->get<Ent_mi_tag>().upper_bound(ent);
324  }

◆ getNumeredColDofsByLocIdxBegin()

auto MoFEM::Problem::getNumeredColDofsByLocIdxBegin ( const DofIdx  locidx) const
inline

get begin iterator for numeredColDofsPtr (insted you can use _IT_NUMEREDDOF_COL_FOR_LOOP_ for loops)

Definition at line 260 of file ProblemsMultiIndices.hpp.

260  {
261  return numeredColDofsPtr->get<PetscLocalIdx_mi_tag>().lower_bound(locidx);
262  }

◆ getNumeredColDofsByLocIdxEnd()

auto MoFEM::Problem::getNumeredColDofsByLocIdxEnd ( const DofIdx  locidx) const
inline

get end iterator for numeredColDofsPtr (insted you can use _IT_NUMEREDDOF_COL_FOR_LOOP_ for loops)

Definition at line 266 of file ProblemsMultiIndices.hpp.

266  {
267  return numeredColDofsPtr->get<PetscLocalIdx_mi_tag>().upper_bound(locidx);
268  }

◆ 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.

207  {
208  return numeredColDofsPtr->end();
209  }

◆ getNumeredColDofsPtr()

auto& MoFEM::Problem::getNumeredColDofsPtr ( ) const
inline

get access to numeredColDofsPtr storing DOFs on cols

Definition at line 87 of file ProblemsMultiIndices.hpp.

87 { return numeredColDofsPtr; }

◆ getNumeredFiniteElementsPtr()

const auto& MoFEM::Problem::getNumeredFiniteElementsPtr ( ) const
inline

get access to reference for multi-index storing finite elements

Definition at line 92 of file ProblemsMultiIndices.hpp.

92  {
94  }

◆ 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.

189  {
190  return numeredRowDofsPtr->begin();
191  }

◆ getNumeredRowDofsByEntBegin()

auto MoFEM::Problem::getNumeredRowDofsByEntBegin ( const EntityHandle  ent) const
inline

get begin iterator for numeredRowDofsPtr (insted you can use _IT_NUMEREDDOF_ROW_BY_ENT_FOR_LOOP_ for loops)

Definition at line 304 of file ProblemsMultiIndices.hpp.

304  {
305  return numeredRowDofsPtr->get<Ent_mi_tag>().lower_bound(ent);
306  }

◆ getNumeredRowDofsByEntEnd()

auto MoFEM::Problem::getNumeredRowDofsByEntEnd ( const EntityHandle  ent) const
inline

get end iterator for numeredRowDofsPtr (insted you can use _IT_NUMEREDDOF_ROW_BY_ENT_FOR_LOOP_ for loops)

Definition at line 310 of file ProblemsMultiIndices.hpp.

310  {
311  return numeredRowDofsPtr->get<Ent_mi_tag>().upper_bound(ent);
312  }

◆ getNumeredRowDofsByLocIdxBegin()

auto MoFEM::Problem::getNumeredRowDofsByLocIdxBegin ( const DofIdx  locidx) const
inline

get begin iterator for numeredRowDofsPtr (insted you can use _IT_NUMEREDDOF_ROW_FOR_LOOP_ for loops)

Definition at line 248 of file ProblemsMultiIndices.hpp.

248  {
249  return numeredRowDofsPtr->get<PetscLocalIdx_mi_tag>().lower_bound(locidx);
250  }

◆ getNumeredRowDofsByLocIdxEnd()

auto MoFEM::Problem::getNumeredRowDofsByLocIdxEnd ( const DofIdx  locidx) const
inline

get end iterator for numeredRowDofsPtr (insted you can use _IT_NUMEREDDOF_ROW_FOR_LOOP_ for loops)

Definition at line 254 of file ProblemsMultiIndices.hpp.

254  {
255  return numeredRowDofsPtr->get<PetscLocalIdx_mi_tag>().upper_bound(locidx);
256  }

◆ 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.

195  {
196  return numeredRowDofsPtr->end();
197  }

◆ getNumeredRowDofsPtr()

auto& MoFEM::Problem::getNumeredRowDofsPtr ( ) const
inline

get access to numeredRowDofsPtr storing DOFs on rows

Examples
test_cache_on_entities.cpp.

Definition at line 82 of file ProblemsMultiIndices.hpp.

82 { return numeredRowDofsPtr; }

◆ getRowDofsByPetscGlobalDofIdx() [1/2]

boost::weak_ptr< NumeredDofEntity > MoFEM::Problem::getRowDofsByPetscGlobalDofIdx ( DofIdx  idx) const

Get the Row Dofs By Petsc Global Dof Idx object.

Parameters
idx
Returns
boost::weak_ptr<NumeredDofEntity>

Definition at line 48 of file ProblemsMultiIndices.cpp.

48  {
50  boost::weak_ptr<NumeredDofEntity> dof_weak_ptr;
51  NumeredDofEntity_multiIndex::index<PetscGlobalIdx_mi_tag>::type::iterator dit;
52  dit = numeredRowDofsPtr->get<PetscGlobalIdx_mi_tag>().find(idx);
53  if (dit != numeredRowDofsPtr->get<PetscGlobalIdx_mi_tag>().end())
54  dof_weak_ptr = *dit;
55  return dof_weak_ptr;
56 }

◆ getRowDofsByPetscGlobalDofIdx() [2/2]

MoFEMErrorCode MoFEM::Problem::getRowDofsByPetscGlobalDofIdx ( DofIdx  idx,
const NumeredDofEntity **  dof_ptr,
MoFEMTypes  bh = MF_EXIST 
) const

Get the Row Dofs By Petsc Global Dof Idx object.

Parameters
idx
dof_ptr
bh
Returns
MoFEMErrorCode

Definition at line 69 of file ProblemsMultiIndices.cpp.

70  {
72  auto weak_dof_ptr = getRowDofsByPetscGlobalDofIdx(idx);
73  if (auto shared_dof_ptr = weak_dof_ptr.lock()) {
74  *dof_ptr = shared_dof_ptr.get();
75  } else {
76  if (bh == MF_EXIST)
77  SETERRQ1(PETSC_COMM_SELF, MOFEM_NOT_FOUND, "row dof <%d> not found", idx);
78  *dof_ptr = nullptr;
79  }
81 }

◆ 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 506 of file ProblemsMultiIndices.hpp.

506 { return sequenceRowDofContainer; }

◆ getSubData()

boost::shared_ptr<SubProblemData>& MoFEM::Problem::getSubData ( ) const
inline

Get main problem of sub-problem is.

Returns
sub problem data structure
Examples
cell_forces.cpp, remove_entities_from_problem.cpp, and remove_entities_from_problem_not_partitioned.cpp.

Definition at line 119 of file ProblemsMultiIndices.hpp.

119  {
120  return subProblemData;
121  }

Friends And Related Function Documentation

◆ operator<<

std::ostream& operator<< ( std::ostream &  os,
const Problem e 
)
friend

Definition at line 37 of file ProblemsMultiIndices.cpp.

37  {
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();
42  return os;
43 }

Member Data Documentation

◆ composedProblemsData

boost::shared_ptr<ComposedProblemsData> MoFEM::Problem::composedProblemsData
mutable

Pointer to data structure from which this problem is composed

Definition at line 126 of file ProblemsMultiIndices.hpp.

◆ emptyFieldBlocks

EmptyFieldBlocks MoFEM::Problem::emptyFieldBlocks
mutableprivate

Definition at line 559 of file ProblemsMultiIndices.hpp.

◆ meshset

EntityHandle MoFEM::Problem::meshset

Problem meshset (on tag of this meshset all data related to problem are stored)

Definition at line 56 of file ProblemsMultiIndices.hpp.

◆ nbDofsCol

DofIdx MoFEM::Problem::nbDofsCol
mutable

Global number of DOFs in col.

Definition at line 66 of file ProblemsMultiIndices.hpp.

◆ nbDofsRow

DofIdx MoFEM::Problem::nbDofsRow
mutable

Global number of DOFs in row.

Examples
partition_mesh.cpp.

Definition at line 65 of file ProblemsMultiIndices.hpp.

◆ nbGhostDofsCol

DofIdx MoFEM::Problem::nbGhostDofsCol
mutable

Number of ghost DOFs in col.

Definition at line 70 of file ProblemsMultiIndices.hpp.

◆ nbGhostDofsRow

DofIdx MoFEM::Problem::nbGhostDofsRow
mutable

Number of ghost DOFs in row.

Definition at line 69 of file ProblemsMultiIndices.hpp.

◆ nbLocDofsCol

DofIdx MoFEM::Problem::nbLocDofsCol
mutable

Local number of DOFs in colIs.

Definition at line 68 of file ProblemsMultiIndices.hpp.

◆ nbLocDofsRow

DofIdx MoFEM::Problem::nbLocDofsRow
mutable

Local number of DOFs in row.

Definition at line 67 of file ProblemsMultiIndices.hpp.

◆ numeredColDofsPtr

boost::shared_ptr<NumeredDofEntity_multiIndex> MoFEM::Problem::numeredColDofsPtr
mutable

store DOFs on columns for this problem

Definition at line 75 of file ProblemsMultiIndices.hpp.

◆ numeredFiniteElementsPtr

boost::shared_ptr<NumeredEntFiniteElement_multiIndex> MoFEM::Problem::numeredFiniteElementsPtr
mutable

store finite elements

Definition at line 77 of file ProblemsMultiIndices.hpp.

◆ numeredRowDofsPtr

boost::shared_ptr<NumeredDofEntity_multiIndex> MoFEM::Problem::numeredRowDofsPtr
mutable

store DOFs on rows for this problem

Examples
cell_forces.cpp.

Definition at line 73 of file ProblemsMultiIndices.hpp.

◆ sequenceColDofContainer

boost::shared_ptr<SequenceDofContainer> MoFEM::Problem::sequenceColDofContainer
mutableprivate

Definition at line 557 of file ProblemsMultiIndices.hpp.

◆ sequenceRowDofContainer

boost::shared_ptr<SequenceDofContainer> MoFEM::Problem::sequenceRowDofContainer
mutableprivate

Definition at line 556 of file ProblemsMultiIndices.hpp.

◆ subProblemData

boost::shared_ptr<SubProblemData> MoFEM::Problem::subProblemData
mutable

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

IDs of finite elements in problem.

Definition at line 61 of file ProblemsMultiIndices.hpp.

◆ tagBitRefLevel

BitRefLevel* MoFEM::Problem::tagBitRefLevel

BitRef level of finite elements in problem.

Definition at line 62 of file ProblemsMultiIndices.hpp.

◆ tagBitRefLevelMask

BitRefLevel* MoFEM::Problem::tagBitRefLevelMask

BItRefMask of elements in problem.

Definition at line 63 of file ProblemsMultiIndices.hpp.

◆ tagId

BitProblemId* MoFEM::Problem::tagId

Unique problem ID.

Definition at line 58 of file ProblemsMultiIndices.hpp.

◆ tagName

const char* MoFEM::Problem::tagName

Problem name.

Definition at line 59 of file ProblemsMultiIndices.hpp.

◆ tagNameSize

int MoFEM::Problem::tagNameSize

Size of problem name.

Definition at line 60 of file ProblemsMultiIndices.hpp.


The documentation for this struct was generated from the following files:
MoFEMFunctionReturnHot
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:447
MoFEM::Problem::tagBitRefLevel
BitRefLevel * tagBitRefLevel
BitRef level of finite elements in problem.
Definition: ProblemsMultiIndices.hpp:62
MoFEM::Problem::nbGhostDofsRow
DofIdx nbGhostDofsRow
Number of ghost DOFs in row.
Definition: ProblemsMultiIndices.hpp:69
MoFEM::DofEntity::getUniqueIdCalculate
static UId getUniqueIdCalculate(const DofIdx dof, UId ent_uid)
Definition: DofsMultiIndices.hpp:54
MoFEM::Types::BitProblemId
std::bitset< BITPROBLEMID_SIZE > BitProblemId
Problem Id.
Definition: Types.hpp:44
MoFEM::Problem::tagBitFEId
BitFEId * tagBitFEId
IDs of finite elements in problem.
Definition: ProblemsMultiIndices.hpp:61
MOAB_THROW
#define MOAB_THROW(err)
Check error code of MoAB function and throw MoFEM exception.
Definition: definitions.h:541
MoFEM::Problem::tagNameSize
int tagNameSize
Size of problem name.
Definition: ProblemsMultiIndices.hpp:60
ROW
@ ROW
Definition: definitions.h:123
NumeredDofEntity_multiIndex
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.
Definition: DofsMultiIndices.hpp:469
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:535
MoFEM::Problem::composedProblemsData
boost::shared_ptr< ComposedProblemsData > composedProblemsData
Definition: ProblemsMultiIndices.hpp:126
MoFEM::Problem::SequenceDofContainer
multi_index_container< boost::weak_ptr< std::vector< NumeredDofEntity > >, indexed_by< sequenced<> > > SequenceDofContainer
Definition: ProblemsMultiIndices.hpp:490
MoFEM::Problem::numeredColDofsPtr
boost::shared_ptr< NumeredDofEntity_multiIndex > numeredColDofsPtr
store DOFs on columns for this problem
Definition: ProblemsMultiIndices.hpp:75
convert.type
type
Definition: convert.py:64
MoFEM::FieldEntity::getLocalUniqueIdCalculate
UId getLocalUniqueIdCalculate()
Get the Local Unique Id Calculate object.
Definition: FieldEntsMultiIndices.hpp:136
COL
@ COL
Definition: definitions.h:123
MoFEM::Problem::nbDofsCol
DofIdx nbDofsCol
Global number of DOFs in col.
Definition: ProblemsMultiIndices.hpp:66
MoFEM::Problem::getColDofsByPetscGlobalDofIdx
MoFEMErrorCode getColDofsByPetscGlobalDofIdx(DofIdx idx, const NumeredDofEntity **dof_ptr, MoFEMTypes bh=MF_EXIST) const
Get the Col Dofs By Petsc Global Dof Idx object.
Definition: ProblemsMultiIndices.cpp:83
MoFEM::Problem::nbLocDofsRow
DofIdx nbLocDofsRow
Local number of DOFs in row.
Definition: ProblemsMultiIndices.hpp:67
MoFEM::Problem::tagName
const char * tagName
Problem name.
Definition: ProblemsMultiIndices.hpp:59
MoFEM::Problem::tagId
BitProblemId * tagId
Unique problem ID.
Definition: ProblemsMultiIndices.hpp:58
MoFEM::Problem::numeredFiniteElementsPtr
boost::shared_ptr< NumeredEntFiniteElement_multiIndex > numeredFiniteElementsPtr
store finite elements
Definition: ProblemsMultiIndices.hpp:77
MOFEM_NOT_FOUND
@ MOFEM_NOT_FOUND
Definition: definitions.h:33
NumeredEntFiniteElement_multiIndex
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.
Definition: FEMultiIndices.hpp:830
MoFEM::Problem::nbGhostDofsCol
DofIdx nbGhostDofsCol
Number of ghost DOFs in col.
Definition: ProblemsMultiIndices.hpp:70
MoFEM::Problem::getRowDofsByPetscGlobalDofIdx
MoFEMErrorCode getRowDofsByPetscGlobalDofIdx(DofIdx idx, const NumeredDofEntity **dof_ptr, MoFEMTypes bh=MF_EXIST) const
Get the Row Dofs By Petsc Global Dof Idx object.
Definition: ProblemsMultiIndices.cpp:69
MoFEM::Problem::numeredRowDofsPtr
boost::shared_ptr< NumeredDofEntity_multiIndex > numeredRowDofsPtr
store DOFs on rows for this problem
Definition: ProblemsMultiIndices.hpp:73
MoFEM::Problem::nbLocDofsCol
DofIdx nbLocDofsCol
Local number of DOFs in colIs.
Definition: ProblemsMultiIndices.hpp:68
MoFEM::Problem::sequenceRowDofContainer
boost::shared_ptr< SequenceDofContainer > sequenceRowDofContainer
Definition: ProblemsMultiIndices.hpp:556
MoFEM::Problem::sequenceColDofContainer
boost::shared_ptr< SequenceDofContainer > sequenceColDofContainer
Definition: ProblemsMultiIndices.hpp:557
MoFEM::Exceptions::ierr
static MoFEMErrorCodeGeneric< PetscErrorCode > ierr
Definition: Exceptions.hpp:76
MoFEM::Problem::nbDofsRow
DofIdx nbDofsRow
Global number of DOFs in row.
Definition: ProblemsMultiIndices.hpp:65
MoFEM::Problem::meshset
EntityHandle meshset
Definition: ProblemsMultiIndices.hpp:56
MoFEMFunctionBeginHot
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:440
MoFEM::Problem::emptyFieldBlocks
EmptyFieldBlocks emptyFieldBlocks
Definition: ProblemsMultiIndices.hpp:559
MF_EXIST
@ MF_EXIST
Definition: definitions.h:100
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:416
CHKERRG
#define CHKERRG(n)
Check error code of MoFEM/MOAB/PETSc function.
Definition: definitions.h:483
MoFEM::Problem::tagBitRefLevelMask
BitRefLevel * tagBitRefLevelMask
BItRefMask of elements in problem.
Definition: ProblemsMultiIndices.hpp:63
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:346
MOFEM_INVALID_DATA
@ MOFEM_INVALID_DATA
Definition: definitions.h:36
MoFEM::Problem::subProblemData
boost::shared_ptr< SubProblemData > subProblemData
Definition: ProblemsMultiIndices.hpp:107