v0.13.2
Loading...
Searching...
No Matches
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 EmptyFieldBlocks = std::pair< BitFieldId, BitFieldId >
 

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...
 
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...
 
EmptyFieldBlocksaddFieldToEmptyFieldBlocks (const EmptyFieldBlocks 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
EshelbianPlasticity.cpp, MagneticElement.hpp, Remodeling.cpp, cell_forces.cpp, field_evaluator.cpp, loop_entities.cpp, lorentz_force.cpp, nonlinear_dynamics.cpp, partition_mesh.cpp, poisson_2d_dis_galerkin.cpp, remove_entities_from_problem.cpp, remove_entities_from_problem_not_partitioned.cpp, and test_cache_on_entities.cpp.

Definition at line 40 of file ProblemsMultiIndices.hpp.

Member Typedef Documentation

◆ EmptyFieldBlocks

Definition at line 502 of file ProblemsMultiIndices.hpp.

◆ SequenceDofContainer

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

Definition at line 468 of file ProblemsMultiIndices.hpp.

Constructor & Destructor Documentation

◆ Problem()

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

Definition at line 18 of file ProblemsMultiIndices.cpp.

24 ErrorCode rval;
25 Tag th_ProblemId;
26 rval = moab.tag_get_handle("_ProblemId", th_ProblemId);
28 rval = moab.tag_get_by_ptr(th_ProblemId, &meshset, 1, (const void **)&tagId);
30 Tag th_ProblemName;
31 rval = moab.tag_get_handle("_ProblemName", th_ProblemName);
33 rval = moab.tag_get_by_ptr(th_ProblemName, &meshset, 1,
34 (const void **)&tagName, &tagNameSize);
36 Tag th_ProblemFEId;
37 rval = moab.tag_get_handle("_ProblemFEId", th_ProblemFEId);
39 rval = moab.tag_get_by_ptr(th_ProblemFEId, &meshset, 1,
40 (const void **)&tagBitFEId);
42 Tag th_RefBitLevel;
43 rval = moab.tag_get_handle("_RefBitLevel", th_RefBitLevel);
45 rval = moab.tag_get_by_ptr(th_RefBitLevel, &meshset, 1,
46 (const void **)&tagBitRefLevel);
48 Tag th_RefBitLevel_Mask;
49 rval = moab.tag_get_handle("_RefBitLevelMask", th_RefBitLevel_Mask);
51 rval = moab.tag_get_by_ptr(th_RefBitLevel_Mask, &meshset, 1,
52 (const void **)&tagBitRefLevelMask);
54}
#define MOAB_THROW(err)
Check error code of MoAB function and throw MoFEM exception.
Definition: definitions.h:541
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.
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.
static MoFEMErrorCodeGeneric< moab::ErrorCode > rval
Definition: Exceptions.hpp:74
const char * tagName
Problem name.
boost::shared_ptr< SequenceDofContainer > sequenceColDofContainer
BitProblemId * tagId
Unique problem ID.
boost::shared_ptr< NumeredDofEntity_multiIndex > numeredRowDofsPtr
store DOFs on rows for this problem
boost::shared_ptr< NumeredEntFiniteElement_multiIndex > numeredFiniteElementsPtr
store finite elements
boost::shared_ptr< SequenceDofContainer > sequenceRowDofContainer
int tagNameSize
Size of problem name.
BitFEId * tagBitFEId
IDs of finite elements in problem.
BitRefLevel * tagBitRefLevelMask
BItRefMask of elements in problem.
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
BitRefLevel * tagBitRefLevel
BitRef level of finite elements in problem.

◆ ~Problem()

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

Member Function Documentation

◆ addFieldToEmptyFieldBlocks()

EmptyFieldBlocks & MoFEM::Problem::addFieldToEmptyFieldBlocks ( const EmptyFieldBlocks  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 526 of file ProblemsMultiIndices.hpp.

526 {
527 emptyFieldBlocks.first |= add_fields.first;
528 emptyFieldBlocks.second |= add_fields.second;
529 return emptyFieldBlocks;
530 }
EmptyFieldBlocks emptyFieldBlocks

◆ getBitFEId()

BitFEId MoFEM::Problem::getBitFEId ( ) const

Definition at line 64 of file ProblemsMultiIndices.cpp.

64{ return *tagBitFEId; }

◆ getBitRefLevel()

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

Definition at line 361 of file ProblemsMultiIndices.hpp.

361{ return *tagBitRefLevel; }

◆ getBitRefLevelMask()

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

Definition at line 362 of file ProblemsMultiIndices.hpp.

362{ 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 78 of file ProblemsMultiIndices.cpp.

78 {
80 boost::weak_ptr<NumeredDofEntity> dof_weak_ptr;
81 NumeredDofEntity_multiIndex::index<PetscGlobalIdx_mi_tag>::type::iterator dit;
82 dit = numeredColDofsPtr->get<PetscGlobalIdx_mi_tag>().find(idx);
83 if (dit != numeredColDofsPtr->get<PetscGlobalIdx_mi_tag>().end())
84 dof_weak_ptr = *dit;
85 return dof_weak_ptr;
86}
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:440

◆ 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 102 of file ProblemsMultiIndices.cpp.

103 {
105 auto weak_dof_ptr = getColDofsByPetscGlobalDofIdx(idx);
106 if (auto shared_dof_ptr = weak_dof_ptr.lock()) {
107 *dof_ptr = shared_dof_ptr.get();
108 } else {
109 if (bh == MF_EXIST)
110 SETERRQ1(PETSC_COMM_SELF, MOFEM_NOT_FOUND, "row dof <%d> not found", idx);
111 *dof_ptr = nullptr;
112 }
114}
@ MF_EXIST
Definition: definitions.h:100
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:447
@ MOFEM_NOT_FOUND
Definition: definitions.h:33
MoFEMErrorCode getColDofsByPetscGlobalDofIdx(DofIdx idx, const NumeredDofEntity **dof_ptr, MoFEMTypes bh=MF_EXIST) const
Get the Col Dofs By Petsc Global Dof Idx object.

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

500{ return sequenceColDofContainer; }

◆ getComposedProblemsData()

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

Het composed problems data structure.

Definition at line 109 of file ProblemsMultiIndices.hpp.

109{ return composedProblemsData; }
boost::shared_ptr< ComposedProblemsData > 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 151 of file ProblemsMultiIndices.cpp.

154 {
156 decltype(numeredRowDofsPtr) numered_dofs;
157 switch (row_or_col) {
158 case ROW:
159 if (!numeredRowDofsPtr) {
160 SETERRQ(PETSC_COMM_SELF, MOFEM_INVALID_DATA,
161 "Row numbered index in problem not allocated");
162 }
163 numered_dofs = numeredRowDofsPtr;
164 break;
165 case COL:
166 if (!numeredColDofsPtr) {
167 SETERRQ(PETSC_COMM_SELF, MOFEM_INVALID_DATA,
168 "Col numbered index in problem not allocated");
169 }
170 numered_dofs = numeredColDofsPtr;
171 break;
172 default:
173 SETERRQ(PETSC_COMM_SELF, MOFEM_INVALID_DATA,
174 "Only ROW and COL is possible for 3rd argument");
175 }
176
177 auto it =
178 numered_dofs->get<Unique_mi_tag>().find(DofEntity::getUniqueIdCalculate(
179 ent_dof_idx,
180 FieldEntity::getLocalUniqueIdCalculate(field_bit_number, ent)));
181 if (it != numered_dofs->get<Unique_mi_tag>().end()) {
182 dof_ptr = *it;
183 } else {
184 dof_ptr = boost::shared_ptr<NumeredDofEntity>();
185 }
187}
@ COL
Definition: definitions.h:123
@ ROW
Definition: definitions.h:123
#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
Definition: definitions.h:36
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:416
static UId getUniqueIdCalculate(const DofIdx dof, UId ent_uid)
UId getLocalUniqueIdCalculate()
Get the Local Unique Id Calculate object.

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

512 {
513 return emptyFieldBlocks;
514 }

◆ getId()

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

Definition at line 348 of file ProblemsMultiIndices.hpp.

348{ return *((BitProblemId *)tagId); }
std::bitset< BITPROBLEMID_SIZE > BitProblemId
Problem Id.
Definition: Types.hpp:44

◆ getName()

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

Definition at line 350 of file ProblemsMultiIndices.hpp.

350 {
351 return std::string((char *)tagName, tagNameSize);
352 }

◆ getNbDofsCol()

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

Definition at line 355 of file ProblemsMultiIndices.hpp.

355{ return nbDofsCol; }
DofIdx nbDofsCol
Global number of DOFs in col.

◆ getNbDofsRow()

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

Definition at line 354 of file ProblemsMultiIndices.hpp.

354{ return nbDofsRow; }
DofIdx nbDofsRow
Global number of DOFs in row.

◆ getNbGhostDofsCol()

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

Definition at line 359 of file ProblemsMultiIndices.hpp.

359{ return nbGhostDofsCol; }
DofIdx nbGhostDofsCol
Number of ghost DOFs in col.

◆ getNbGhostDofsRow()

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

Definition at line 358 of file ProblemsMultiIndices.hpp.

358{ return nbGhostDofsRow; }
DofIdx nbGhostDofsRow
Number of ghost DOFs in row.

◆ getNbLocalDofsCol()

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

Definition at line 357 of file ProblemsMultiIndices.hpp.

357{ return nbLocDofsCol; }
DofIdx nbLocDofsCol
Local number of DOFs in colIs.

◆ getNbLocalDofsRow()

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

Definition at line 356 of file ProblemsMultiIndices.hpp.

356{ return nbLocDofsRow; }
DofIdx nbLocDofsRow
Local number of DOFs in row.

◆ 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);
#define CHKERRG(n)
Check error code of MoFEM/MOAB/PETSc function.
Definition: definitions.h:483
static MoFEMErrorCodeGeneric< PetscErrorCode > ierr
Definition: Exceptions.hpp:76
Parameters
commCommunicator
nameFinite element name
layoutGet number of elements on each processor
Returns
error code

Definition at line 117 of file ProblemsMultiIndices.cpp.

118 {
120 int size, rank;
121 MPI_Comm_size(comm, &size);
122 MPI_Comm_rank(comm, &rank);
123 CHKERR PetscLayoutCreate(comm, layout);
124 CHKERR PetscLayoutSetBlockSize(*layout, 1);
125 const auto &fe_by_name_and_part =
126 numeredFiniteElementsPtr->get<Composite_Name_And_Part_mi_tag>();
127 int nb_elems;
128 nb_elems = fe_by_name_and_part.count(boost::make_tuple(name, rank));
129 CHKERR PetscLayoutSetLocalSize(*layout, nb_elems);
130 CHKERR PetscLayoutSetUp(*layout);
132}
#define CHKERR
Inline error check.
Definition: definitions.h:535

◆ 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 134 of file ProblemsMultiIndices.cpp.

135 {
137 int size, rank;
138 MPI_Comm_size(comm, &size);
139 MPI_Comm_rank(comm, &rank);
140 CHKERR PetscLayoutCreate(comm, layout);
141 CHKERR PetscLayoutSetBlockSize(*layout, 1);
142 typedef NumeredEntFiniteElement_multiIndex::index<Part_mi_tag>::type FeByPart;
143 const FeByPart &fe_by_part = numeredFiniteElementsPtr->get<Part_mi_tag>();
144 int nb_elems;
145 nb_elems = fe_by_part.count(rank);
146 CHKERR PetscLayoutSetLocalSize(*layout, nb_elems);
147 CHKERR PetscLayoutSetUp(*layout);
149}

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

179 {
180 return numeredColDofsPtr->begin();
181 }

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

294 {
295 return numeredColDofsPtr->get<Ent_mi_tag>().lower_bound(ent);
296 }

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

300 {
301 return numeredColDofsPtr->get<Ent_mi_tag>().upper_bound(ent);
302 }

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

238 {
239 return numeredColDofsPtr->get<PetscLocalIdx_mi_tag>().lower_bound(locidx);
240 }

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

244 {
245 return numeredColDofsPtr->get<PetscLocalIdx_mi_tag>().upper_bound(locidx);
246 }

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

185 {
186 return numeredColDofsPtr->end();
187 }

◆ getNumeredColDofsPtr()

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

get access to numeredColDofsPtr storing DOFs on cols

Definition at line 73 of file ProblemsMultiIndices.hpp.

73{ return numeredColDofsPtr; }

◆ getNumeredFiniteElementsPtr()

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

get access to reference for multi-index storing finite elements

Definition at line 78 of file ProblemsMultiIndices.hpp.

78 {
80 }

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

167 {
168 return numeredRowDofsPtr->begin();
169 }

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

282 {
283 return numeredRowDofsPtr->get<Ent_mi_tag>().lower_bound(ent);
284 }

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

288 {
289 return numeredRowDofsPtr->get<Ent_mi_tag>().upper_bound(ent);
290 }

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

226 {
227 return numeredRowDofsPtr->get<PetscLocalIdx_mi_tag>().lower_bound(locidx);
228 }

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

232 {
233 return numeredRowDofsPtr->get<PetscLocalIdx_mi_tag>().upper_bound(locidx);
234 }

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

173 {
174 return numeredRowDofsPtr->end();
175 }

◆ getNumeredRowDofsPtr()

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

get access to numeredRowDofsPtr storing DOFs on rows

Examples
eigen_elastic.cpp, and test_cache_on_entities.cpp.

Definition at line 68 of file ProblemsMultiIndices.hpp.

68{ 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 67 of file ProblemsMultiIndices.cpp.

67 {
69 boost::weak_ptr<NumeredDofEntity> dof_weak_ptr;
70 NumeredDofEntity_multiIndex::index<PetscGlobalIdx_mi_tag>::type::iterator dit;
71 dit = numeredRowDofsPtr->get<PetscGlobalIdx_mi_tag>().find(idx);
72 if (dit != numeredRowDofsPtr->get<PetscGlobalIdx_mi_tag>().end())
73 dof_weak_ptr = *dit;
74 return dof_weak_ptr;
75}

◆ 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 88 of file ProblemsMultiIndices.cpp.

89 {
91 auto weak_dof_ptr = getRowDofsByPetscGlobalDofIdx(idx);
92 if (auto shared_dof_ptr = weak_dof_ptr.lock()) {
93 *dof_ptr = shared_dof_ptr.get();
94 } else {
95 if (bh == MF_EXIST)
96 SETERRQ1(PETSC_COMM_SELF, MOFEM_NOT_FOUND, "row dof <%d> not found", idx);
97 *dof_ptr = nullptr;
98 }
100}
MoFEMErrorCode getRowDofsByPetscGlobalDofIdx(DofIdx idx, const NumeredDofEntity **dof_ptr, MoFEMTypes bh=MF_EXIST) const
Get the Row Dofs By Petsc Global Dof Idx object.

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

484{ 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 97 of file ProblemsMultiIndices.hpp.

97 {
98 return subProblemData;
99 }
boost::shared_ptr< SubProblemData > subProblemData

Friends And Related Function Documentation

◆ operator<<

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

Definition at line 56 of file ProblemsMultiIndices.cpp.

56 {
57 os << "Problem id " << e.getId() << " name " << e.getName() << endl;
58 os << "FiniteElement id " << e.getBitFEId() << endl;
59 os << "BitRefLevel " << e.getBitRefLevel() << endl;
60 os << "BitRefLevelMask " << e.getBitRefLevelMask();
61 return os;
62}

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

◆ emptyFieldBlocks

EmptyFieldBlocks MoFEM::Problem::emptyFieldBlocks
mutableprivate

Definition at line 537 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 42 of file ProblemsMultiIndices.hpp.

◆ nbDofsCol

DofIdx MoFEM::Problem::nbDofsCol
mutable

Global number of DOFs in col.

Definition at line 52 of file ProblemsMultiIndices.hpp.

◆ nbDofsRow

DofIdx MoFEM::Problem::nbDofsRow
mutable

Global number of DOFs in row.

Examples
partition_mesh.cpp.

Definition at line 51 of file ProblemsMultiIndices.hpp.

◆ nbGhostDofsCol

DofIdx MoFEM::Problem::nbGhostDofsCol
mutable

Number of ghost DOFs in col.

Definition at line 56 of file ProblemsMultiIndices.hpp.

◆ nbGhostDofsRow

DofIdx MoFEM::Problem::nbGhostDofsRow
mutable

Number of ghost DOFs in row.

Definition at line 55 of file ProblemsMultiIndices.hpp.

◆ nbLocDofsCol

DofIdx MoFEM::Problem::nbLocDofsCol
mutable

Local number of DOFs in colIs.

Definition at line 54 of file ProblemsMultiIndices.hpp.

◆ nbLocDofsRow

DofIdx MoFEM::Problem::nbLocDofsRow
mutable

Local number of DOFs in row.

Definition at line 53 of file ProblemsMultiIndices.hpp.

◆ numeredColDofsPtr

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

store DOFs on columns for this problem

Definition at line 61 of file ProblemsMultiIndices.hpp.

◆ numeredFiniteElementsPtr

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

store finite elements

Definition at line 63 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 59 of file ProblemsMultiIndices.hpp.

◆ sequenceColDofContainer

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

Definition at line 535 of file ProblemsMultiIndices.hpp.

◆ sequenceRowDofContainer

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

Definition at line 534 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.

Examples
EshelbianPlasticity.cpp.

Definition at line 91 of file ProblemsMultiIndices.hpp.

◆ tagBitFEId

BitFEId* MoFEM::Problem::tagBitFEId

IDs of finite elements in problem.

Definition at line 47 of file ProblemsMultiIndices.hpp.

◆ tagBitRefLevel

BitRefLevel* MoFEM::Problem::tagBitRefLevel

BitRef level of finite elements in problem.

Definition at line 48 of file ProblemsMultiIndices.hpp.

◆ tagBitRefLevelMask

BitRefLevel* MoFEM::Problem::tagBitRefLevelMask

BItRefMask of elements in problem.

Definition at line 49 of file ProblemsMultiIndices.hpp.

◆ tagId

BitProblemId* MoFEM::Problem::tagId

Unique problem ID.

Definition at line 44 of file ProblemsMultiIndices.hpp.

◆ tagName

const char* MoFEM::Problem::tagName

Problem name.

Definition at line 45 of file ProblemsMultiIndices.hpp.

◆ tagNameSize

int MoFEM::Problem::tagNameSize

Size of problem name.

Definition at line 46 of file ProblemsMultiIndices.hpp.


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