v0.14.0
BCMultiIndices.cpp
Go to the documentation of this file.
1 /** \file BCMultiIndices.cpp
2  * \brief Structures for managing boundary conditions
3  */
4 
5 namespace MoFEM {
6 
7 // moab base meshsets
10  CHKERR moab.tag_get_handle(DIRICHLET_SET_TAG_NAME, nsTag);
11  CHKERR moab.tag_get_handle(NEUMANN_SET_TAG_NAME, ssTag);
12  CHKERR moab.tag_get_handle(
13  (std::string(DIRICHLET_SET_TAG_NAME) + "__BC_DATA").c_str(), nsTag_data);
14  CHKERR moab.tag_get_handle(
15  (std::string(NEUMANN_SET_TAG_NAME) + "__BC_DATA").c_str(), ssTag_data);
16  CHKERR moab.tag_get_handle(MATERIAL_SET_TAG_NAME, bhTag);
17  CHKERR moab.tag_get_handle(BLOCK_HEADER, bhTag_header);
18  CHKERR moab.tag_get_handle(BLOCK_ATTRIBUTES, thBlockAttribs);
19  CHKERR moab.tag_get_handle(NAME_TAG_NAME, entityNameTag);
21 }
23  : meshset(meshset), cubitBcType(UNKNOWNSET), msId(nullptr),
24  tagBcData(nullptr), tagBcSize(0), tagBlockHeaderData(nullptr),
25  tagBlockAttributes(nullptr), tagBlockAttributesSize(0), tagName(nullptr),
26  meshsetsMask(NODESET | SIDESET | BLOCKSET) {
27 
28  auto construct = [&]() {
30 
31  CHKERR getTagsHandlers(moab);
32  CHKERR moab.tag_get_tags_on_entity(meshset, tag_handles);
33  for (auto tit = tag_handles.begin(); tit != tag_handles.end(); tit++) {
34  if (*tit == nsTag || *tit == ssTag || *tit == bhTag) {
35  CHKERR moab.tag_get_by_ptr(*tit, &meshset, 1, (const void **)&msId);
36  }
37  if (*tit == nsTag) {
38  if (*msId != -1) {
40  }
41  }
42  if (*tit == ssTag) {
43  if (*msId != -1) {
45  }
46  }
47  if (*tit == bhTag) {
48  if (*msId != -1) {
50  }
51  }
52  }
53 
54  for (auto tit = tag_handles.begin(); tit != tag_handles.end(); tit++) {
55  if ((cubitBcType & CubitBCType(NODESET | SIDESET)).any()) {
56  if ((*tit == nsTag_data) || (*tit == ssTag_data)) {
57  CHKERR moab.tag_get_by_ptr(*tit, &meshset, 1,
58  (const void **)&tagBcData, &tagBcSize);
60  }
61  }
62  if ((cubitBcType & CubitBCType(BLOCKSET)).any()) {
63  if (*tit == bhTag_header) {
64  int tag_length;
65  CHKERR moab.tag_get_length(*tit, tag_length);
66  CHKERR moab.tag_get_by_ptr(*tit, &meshset, 1,
67  (const void **)&tagBlockHeaderData);
68  if (tagBlockHeaderData[1] > 0)
70  }
71  if (*tit == thBlockAttribs) {
72  CHKERR moab.tag_get_by_ptr(*tit, &meshset, 1,
73  (const void **)&tagBlockAttributes,
75  }
76  if (*tit == entityNameTag) {
77  CHKERR moab.tag_get_by_ptr(entityNameTag, &meshset, 1,
78  (const void **)&tagName);
80  }
81  }
82  }
83 
84  // If BC set has name, unset UNKNOWNNAME
85  if (cubitBcType.to_ulong() &
88  if ((cubitBcType & CubitBCType(UNKNOWNNAME)).any()) {
90  }
91  }
92 
94  };
95 
96  CHK_THROW_MESSAGE(construct(), "construct CubitMeshSets");
97 }
99  const CubitBCType cubit_bc_type, const int ms_id)
100  : cubitBcType(cubit_bc_type), msId(nullptr), tagBcData(nullptr),
101  tagBcSize(0), tagBlockHeaderData(nullptr), tagBlockAttributes(nullptr),
102  tagBlockAttributesSize(0), tagName(nullptr),
103  meshsetsMask(NODESET | SIDESET | BLOCKSET) {
104 
105  CHKERR getTagsHandlers(moab);
106  CHKERR moab.create_meshset(MESHSET_SET, meshset);
107  switch (cubitBcType.to_ulong()) {
108  case NODESET:
109  CHKERR moab.tag_set_data(nsTag, &meshset, 1, &ms_id);
110  CHKERR moab.tag_get_by_ptr(nsTag, &meshset, 1, (const void **)&msId);
111  break;
112  case SIDESET:
113  CHKERR moab.tag_set_data(ssTag, &meshset, 1, &ms_id);
114  CHKERR moab.tag_get_by_ptr(ssTag, &meshset, 1, (const void **)&msId);
115  break;
116  case BLOCKSET:
117  CHKERR moab.tag_set_data(bhTag, &meshset, 1, &ms_id);
118  CHKERR moab.tag_get_by_ptr(bhTag, &meshset, 1, (const void **)&msId);
119  break;
120  default: {
121  PetscTraceBackErrorHandler(PETSC_COMM_WORLD, __LINE__, PETSC_FUNCTION_NAME,
122  __FILE__, MOFEM_DATA_INCONSISTENCY,
123  PETSC_ERROR_INITIAL, "Unknow meshset type",
124  PETSC_NULL);
125  PetscMPIAbortErrorHandler(PETSC_COMM_WORLD, __LINE__, PETSC_FUNCTION_NAME,
126  __FILE__, MOFEM_DATA_INCONSISTENCY,
127  PETSC_ERROR_INITIAL, "Unknown meshset type",
128  PETSC_NULL);
129  }
130  }
131 }
133  moab::Interface &moab, const int dimension, Range &entities,
134  const bool recursive) const {
136  rval =
137  moab.get_entities_by_dimension(meshset, dimension, entities, recursive);
138  if (rval != MB_SUCCESS) {
139  std::ostringstream ss;
140  ss << "bc set " << *this << std::endl;
141  PetscPrintf(PETSC_COMM_WORLD, ss.str().c_str());
142  }
143  CHKERR rval;
145 }
147  moab::Interface &moab, Range &entities, const bool recursive) const {
149  if ((cubitBcType & CubitBCType(BLOCKSET)).any()) {
150  if (tagBlockHeaderData != nullptr) {
152  entities, recursive);
153  } else {
154  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "dimension unknown");
155  }
156  }
157  if ((cubitBcType & CubitBCType(NODESET)).any()) {
158  return getMeshsetIdEntitiesByDimension(moab, 0, entities, recursive);
159  }
161 }
163  moab::Interface &moab, const EntityType type, Range &entities,
164  const bool recursive) const {
166  rval = moab.get_entities_by_type(meshset, type, entities, recursive);
167  if (rval != MB_SUCCESS) {
168  std::ostringstream ss;
169  ss << "bc set " << *this << std::endl;
170  PetscPrintf(PETSC_COMM_WORLD, ss.str().c_str());
171  }
172  CHKERR rval;
174 }
175 
176 MoFEMErrorCode CubitMeshSets::getBcData(std::vector<char> &bc_data) const {
178  bc_data.resize(tagBcSize);
179  copy(&tagBcData[0], &tagBcData[tagBcSize], bc_data.begin());
181 }
182 
184  std::vector<unsigned int> &material_data) const {
186  material_data.resize(3);
187  copy(&tagBlockHeaderData[0], &tagBlockHeaderData[3], material_data.begin());
189 }
190 
193  if (tagBlockHeaderData != nullptr) {
194  std::vector<unsigned int> material_data;
195  getBlockHeaderData(material_data);
196  os << "block_header_data = ";
197  std::vector<unsigned int>::iterator vit = material_data.begin();
198  for (; vit != material_data.end(); vit++) {
199  os << std::hex << (int)((unsigned int)*vit) << " ";
200  }
201  os << ": ";
202  vit = material_data.begin();
203  for (; vit != material_data.end(); vit++) {
204  os << *vit;
205  }
206  os << std::endl;
207  } else {
208  os << "no block header data" << std::endl;
209  }
211 }
212 
213 std::string CubitMeshSets::getName() const {
214  if (tagName != nullptr) {
215  return std::string(tagName);
216  } else {
217  return "NoNameSet";
218  }
219 }
220 
221 MoFEMErrorCode CubitMeshSets::printName(std::ostream &os) const {
223  std::string name = getName();
224  os << std::endl;
225  os << "Block name: " << name << std::endl;
227 }
228 
230 CubitMeshSets::getTypeFromBcData(const std::vector<char> &bc_data,
231  CubitBCType &type) const {
233 
234  // See CubitBCType in common.hpp
235  if (bc_data.size() == 0) {
237  }
238 
241  UNKNOWNNAME);
242  if (strcmp(&bc_data[0], "Displacement") == 0)
244  else if (strcmp(&bc_data[0], "Force") == 0)
245  type |= FORCESET;
246  else if (strcmp(&bc_data[0], "Velocity") == 0)
247  type |= VELOCITYSET;
248  else if (strcmp(&bc_data[0], "Acceleration") == 0)
250  else if (strcmp(&bc_data[0], "Temperature") == 0)
251  type |= TEMPERATURESET;
252  else if (strcmp(&bc_data[0], "Pressure") == 0)
253  type |= PRESSURESET;
254  else if (strcmp(&bc_data[0], "HeatFlux") == 0)
255  type |= HEATFLUXSET;
256  else if (strcmp(&bc_data[0], "cfd_bc") == 0)
257  type |= INTERFACESET;
258  else
259  type |= UNKNOWNNAME;
260 
262 }
265  std::vector<char> bc_data;
266  CHKERR getBcData(bc_data);
267  CHKERR getTypeFromBcData(bc_data, type);
269 }
270 
273  std::vector<char> bc_data;
274  CHKERR getBcData(bc_data);
275  os << "bc_data = ";
276  std::vector<char>::iterator vit = bc_data.begin();
277  for (; vit != bc_data.end(); vit++) {
278  os << std::hex << (int)((unsigned char)*vit) << " ";
279  }
280  os << ": ";
281  vit = bc_data.begin();
282  for (; vit != bc_data.end(); vit++) {
283  os << *vit;
284  }
285  os << std::endl;
287 }
288 
290 CubitMeshSets::getAttributes(std::vector<double> &attributes) const {
292  attributes.resize(tagBlockAttributesSize);
293  if (tagBlockAttributesSize > 0) {
295  attributes.begin());
296  }
298 }
299 
302  const std::vector<double> &attributes) {
303 
305  int tag_size[] = {(int)attributes.size()};
306  void const *tag_data[] = {&*attributes.begin()};
307  CHKERR moab.tag_set_by_ptr(thBlockAttribs, &meshset, 1, tag_data, tag_size);
308  CHKERR moab.tag_get_by_ptr(thBlockAttribs, &meshset, 1,
309  (const void **)&tagBlockAttributes,
312 }
313 
316  std::vector<double> attributes;
317  CHKERR getAttributes(attributes);
318  os << std::endl;
319  os << "Block attributes" << std::endl;
320  os << "----------------" << std::endl;
321  for (unsigned int ii = 0; ii < attributes.size(); ii++) {
322  os << "attr. no: " << ii + 1 << " value: " << attributes[ii] << std::endl;
323  }
324  os << std::endl;
326 }
327 
329  CubitBCType &type) const {
331  // See CubitBCType in common.hpp
333  BODYFORCESSET);
334 
335  if (name.compare(0, 11, "MAT_ELASTIC") == 0) {
336  type |= MAT_ELASTICSET;
337  } else if (name.compare(0, 11, "MAT_THERMAL") == 0) {
338  type |= MAT_THERMALSET;
339  } else if (name.compare(0, 12, "MAT_MOISTURE") == 0) {
341  } else if (name.compare(0, 10, "MAT_INTERF") == 0) {
342  type |= MAT_INTERFSET;
343  } else if (name.compare(0, 11, "BODY_FORCES") == 0) {
344  type |= BODYFORCESSET;
345  }
346  // To be extended as appropriate
347  else {
348  type |= UNKNOWNNAME;
349  }
351 }
352 
355  std::string name = getName();
356  CHKERR getTypeFromName(name, type);
358 }
359 
360 std::ostream &operator<<(std::ostream &os, const CubitMeshSets &e) {
361  // get name of cubit meshset
362  std::ostringstream ss;
363  unsigned jj = 0;
364  while (1 << jj != LASTSET_BC) {
365  const CubitBCType jj_bc_type = 1 << jj;
366  if ((e.cubitBcType & jj_bc_type).any()) {
367  string bc_type_name;
368  ss << " " << string(CubitBCNames[jj + 1]);
369  }
370  ++jj;
371  }
372 
373  // push data to stream
374  os << "meshset " << e.meshset << " type" << ss.str();
375  if (e.msId != nullptr)
376  os << " msId " << *(e.msId);
377  if (e.tagName != nullptr) {
378  os << " name " << e.getName();
379  }
380  if (e.tagBlockHeaderData != nullptr) {
381  os << " block header: ";
382  os << " blockCol = " << e.tagBlockHeaderData[0];
383  os << " blockMat = " << e.tagBlockHeaderData[1];
384  os << " blockDimension = " << e.tagBlockHeaderData[2];
385  }
386  return os;
387 }
388 
390 
391  switch (e.cubitBcType.to_ulong()) {
392  case BLOCKSET: {
393  nAme.resize(NAME_TAG_SIZE);
394  CHKERR mOab.tag_set_data(e.entityNameTag, &e.meshset, 1, nAme.c_str());
395  CHKERR mOab.tag_get_by_ptr(e.entityNameTag, &e.meshset, 1,
396  (const void **)&e.tagName);
397 
400  e.cubitBcType |= type;
401  }; break;
402  case NODESET:
403  case SIDESET: {
404  nAme.resize(NAME_TAG_SIZE);
405  CHKERR mOab.tag_set_data(e.entityNameTag, &e.meshset, 1, nAme.c_str());
406  CHKERR mOab.tag_get_by_ptr(e.entityNameTag, &e.meshset, 1,
407  (const void **)&e.tagName);
408  }; break;
409  default:
410  THROW_MESSAGE("not implemented for this CubitBC type");
411  }
412 }
413 
415  CubitMeshSets &e) {
416  e.cubitBcType |= bIt;
417 }
418 
421 }
422 
424  CubitMeshSets &e) {
425  // Need to run this to set tag size in number of doubles, don;t know nothing
426  // about structure
427  int tag_size[] = {(int)(aTtr.getSizeOfData() / sizeof(double))};
428  void const *tag_data[] = {aTtr.getDataPtr()};
429  CHKERR mOab.tag_set_by_ptr(e.thBlockAttribs, &e.meshset, 1, tag_data,
430  tag_size);
431  CHKERR mOab.tag_get_by_ptr(e.thBlockAttribs, &e.meshset, 1,
432  (const void **)&e.tagBlockAttributes,
434  // Here I know about structure
436 }
437 
439 
440  // Need to run this to set tag size, don;t know nothing about structure
441  int tag_size[] = {(int)bcData.getSizeOfData()};
442  void const *tag_data[] = {bcData.getDataPtr()};
443  if ((e.cubitBcType & CubitBCType(NODESET)).any()) {
444  CHKERR mOab.tag_set_by_ptr(e.nsTag_data, &e.meshset, 1, tag_data, tag_size);
445  CHKERR mOab.tag_get_by_ptr(e.nsTag_data, &e.meshset, 1,
446  (const void **)&e.tagBcData, &e.tagBcSize);
447  } else if ((e.cubitBcType & CubitBCType(SIDESET)).any()) {
448  CHKERR mOab.tag_set_by_ptr(e.ssTag_data, &e.meshset, 1, tag_data, tag_size);
449  CHKERR mOab.tag_get_by_ptr(e.ssTag_data, &e.meshset, 1,
450  (const void **)&e.tagBcData, &e.tagBcSize);
451  } else {
452  THROW_MESSAGE("You have to have NODESET or SIDESET to apply BC data on it");
453  }
454  // Here I know about structure
456  // Get Type form BC data
458 }
459 
460 } // namespace MoFEM
CubitBCNames
const static char *const CubitBCNames[]
Names of types of sets and boundary conditions.
Definition: definitions.h:188
MoFEMFunctionReturnHot
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:460
MoFEM::operator<<
std::ostream & operator<<(std::ostream &os, const EntitiesFieldData::EntData &e)
Definition: EntitiesFieldData.cpp:240
MoFEM::CubitMeshSets::cubitBcType
CubitBCType cubitBcType
Definition: BCMultiIndices.hpp:22
SIDESET
@ SIDESET
Definition: definitions.h:160
MoFEM::CubitMeshSets::CubitMeshSets
CubitMeshSets(Interface &moab, const EntityHandle meshset)
MoFEM::CubitMeshSets::msId
int * msId
cubit meshset ID
Definition: BCMultiIndices.hpp:26
MoFEM::CubitMeshSets::printAttributes
MoFEMErrorCode printAttributes(std::ostream &os) const
print the attributes vector
Definition: BCMultiIndices.cpp:314
EntityHandle
MoFEM::CubitMeshSets::bhTag
Tag bhTag
Definition: BCMultiIndices.hpp:322
MoFEM::CubitMeshSets_change_attributes_data_structure::mOab
Interface & mOab
Definition: BCMultiIndices.hpp:426
PRESSURESET
@ PRESSURESET
Definition: definitions.h:165
CHK_THROW_MESSAGE
#define CHK_THROW_MESSAGE(err, msg)
Check and throw MoFEM exception.
Definition: definitions.h:609
MoFEM::CubitMeshSets::bhTag_header
Tag bhTag_header
Definition: BCMultiIndices.hpp:322
MoFEM::Exceptions::MoFEMErrorCode
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:56
MoFEM::CubitMeshSets
this struct keeps basic methods for moab meshset about material and boundary conditions
Definition: BCMultiIndices.hpp:19
MoFEM::CubitMeshSets_change_attributes::mOab
Interface & mOab
Definition: BCMultiIndices.hpp:414
MoFEM::CubitMeshSets_change_bc_data_structure::mOab
Interface & mOab
Definition: BCMultiIndices.hpp:438
MoFEM::CubitMeshSets::printBlockHeaderData
MoFEMErrorCode printBlockHeaderData(std::ostream &os) const
print material_data int stream given by os
Definition: BCMultiIndices.cpp:191
MoFEM::CubitMeshSets::tagBcSize
int tagBcSize
Definition: BCMultiIndices.hpp:28
MoFEM::CubitMeshSets::tagBlockAttributes
double * tagBlockAttributes
Definition: BCMultiIndices.hpp:30
MoFEM::GenericCubitBcData::getDataPtr
virtual const void * getDataPtr() const =0
get pointer to data structure
MoFEM::CubitMeshSets::tagName
char * tagName
Definition: BCMultiIndices.hpp:32
MoFEM::CubitMeshSets::getMeshsetIdEntitiesByType
MoFEMErrorCode getMeshsetIdEntitiesByType(Interface &moab, const EntityType type, Range &entities, const bool recursive=false) const
get entities by type
Definition: BCMultiIndices.cpp:162
MoFEM::CubitMeshSets_change_add_bit_to_cubit_bc_type::operator()
void operator()(CubitMeshSets &e)
Definition: BCMultiIndices.cpp:414
MoFEM::CubitMeshSets_change_bc_data_structure::operator()
void operator()(CubitMeshSets &e)
Definition: BCMultiIndices.cpp:438
MoFEM::CubitMeshSets::entityNameTag
Tag entityNameTag
Definition: BCMultiIndices.hpp:323
THROW_MESSAGE
#define THROW_MESSAGE(msg)
Throw MoFEM exception.
Definition: definitions.h:574
MoFEM::CubitMeshSets_change_attributes::aTtr
const std::vector< double > & aTtr
Definition: BCMultiIndices.hpp:415
MoFEM::CubitMeshSets::setAttributes
MoFEMErrorCode setAttributes(moab::Interface &moab, const std::vector< double > &attributes)
cet Cubit block attributes
Definition: BCMultiIndices.cpp:301
MoFEM::CubitMeshSets::getTagsHandlers
MoFEMErrorCode getTagsHandlers(Interface &moab)
Definition: BCMultiIndices.cpp:8
MoFEM::CubitMeshSets::thBlockAttribs
Tag thBlockAttribs
Definition: BCMultiIndices.hpp:322
UNKNOWNNAME
@ UNKNOWNNAME
Definition: definitions.h:171
MoFEM::Interface
DeprecatedCoreInterface Interface
Definition: Interface.hpp:2010
MoFEM::GenericCubitBcData::getSizeOfData
virtual std::size_t getSizeOfData() const =0
get data structure size
MoFEM::CubitMeshSets_change_name::mOab
Interface & mOab
Definition: BCMultiIndices.hpp:403
MoFEM::CubitMeshSets_change_name::operator()
void operator()(CubitMeshSets &e)
Definition: BCMultiIndices.cpp:389
MoFEM::CubitMeshSets::tag_handles
std::vector< Tag > tag_handles
vector of tag handles to types of data passed from cubit
Definition: BCMultiIndices.hpp:25
NODESET
@ NODESET
Definition: definitions.h:159
MATERIALSET
@ MATERIALSET
Definition: definitions.h:162
MoFEM::CubitMeshSets::nsTag
Tag nsTag
Definition: BCMultiIndices.hpp:322
MoFEM::Exceptions::rval
static MoFEMErrorCodeGeneric< moab::ErrorCode > rval
Definition: Exceptions.hpp:74
MoFEM::CubitMeshSets::getName
std::string getName() const
get name of block, sideset etc. (this is set in Cubit block properties)
Definition: BCMultiIndices.cpp:213
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:548
MoFEM
implementation of Data Operators for Forces and Sources
Definition: Common.hpp:10
MoFEM::CubitMeshSets::setBcDataStructure
MoFEMErrorCode setBcDataStructure(CUBIT_BC_DATA_TYPE &data)
Definition: BCMultiIndices.hpp:311
MoFEM::CubitMeshSets::printName
MoFEMErrorCode printName(std::ostream &os) const
print name of block, sideset etc. (this is set in Cubit setting properties)
Definition: BCMultiIndices.cpp:221
FORCESET
@ FORCESET
Definition: definitions.h:164
MAT_THERMALSET
@ MAT_THERMALSET
block name is "MAT_THERMAL"
Definition: definitions.h:174
BODYFORCESSET
@ BODYFORCESSET
block name is "BODY_FORCES"
Definition: definitions.h:175
double
convert.type
type
Definition: convert.py:64
ACCELERATIONSET
@ ACCELERATIONSET
Definition: definitions.h:167
MAT_ELASTICSET
@ MAT_ELASTICSET
block name is "MAT_ELASTIC"
Definition: definitions.h:172
DISPLACEMENTSET
@ DISPLACEMENTSET
Definition: definitions.h:163
INTERFACESET
@ INTERFACESET
Definition: definitions.h:170
MoFEM::CubitMeshSets::meshset
EntityHandle meshset
Definition: BCMultiIndices.hpp:21
MoFEM::CubitMeshSets::setAttributeDataStructure
MoFEMErrorCode setAttributeDataStructure(const ATTRIBUTE_TYPE &data)
fill meshset data with data on structure
Definition: BCMultiIndices.hpp:284
MoFEM::CubitMeshSets_change_name::nAme
std::string nAme
Definition: BCMultiIndices.hpp:404
MoFEM::CubitMeshSets_change_add_bit_to_cubit_bc_type::bIt
CubitBCType bIt
Definition: BCMultiIndices.hpp:393
MoFEM::CubitMeshSets_change_attributes_data_structure::operator()
void operator()(CubitMeshSets &e)
Definition: BCMultiIndices.cpp:423
MoFEM::CubitMeshSets::ssTag
Tag ssTag
Definition: BCMultiIndices.hpp:322
UNKNOWNSET
@ UNKNOWNSET
Definition: definitions.h:158
MoFEM::GenericAttributeData::getSizeOfData
virtual std::size_t getSizeOfData() const =0
get data structure size
MoFEM::CubitMeshSets::tagBlockAttributesSize
int tagBlockAttributesSize
Definition: BCMultiIndices.hpp:31
VELOCITYSET
@ VELOCITYSET
Definition: definitions.h:166
MoFEM::CubitMeshSets_change_bc_data_structure::bcData
const GenericCubitBcData & bcData
Definition: BCMultiIndices.hpp:439
LASTSET_BC
@ LASTSET_BC
Definition: definitions.h:179
MoFEM::CubitMeshSets::tagBlockHeaderData
unsigned int * tagBlockHeaderData
Definition: BCMultiIndices.hpp:29
Range
MoFEM::Types::CubitBCType
std::bitset< 32 > CubitBCType
Definition: Types.hpp:52
MoFEM::CubitMeshSets::getAttributes
MoFEMErrorCode getAttributes(std::vector< double > &attributes) const
get Cubit block attributes
Definition: BCMultiIndices.cpp:290
BLOCKSET
@ BLOCKSET
Definition: definitions.h:161
MoFEM::CubitMeshSets::getBcData
MoFEMErrorCode getBcData(std::vector< char > &bc_data) const
get bc_data vector from MoFEM database
Definition: BCMultiIndices.cpp:176
MOFEM_DATA_INCONSISTENCY
@ MOFEM_DATA_INCONSISTENCY
Definition: definitions.h:31
TEMPERATURESET
@ TEMPERATURESET
Definition: definitions.h:168
MoFEM::CubitMeshSets::getTypeFromName
MoFEMErrorCode getTypeFromName(const std::string &name, CubitBCType &type) const
Function that returns the CubitBCType type of the block name, sideset name etc.
Definition: BCMultiIndices.cpp:328
MoFEM::CubitMeshSets_change_attributes_data_structure::aTtr
const GenericAttributeData & aTtr
Definition: BCMultiIndices.hpp:427
MoFEM::CubitMeshSets::nsTag_data
Tag nsTag_data
Definition: BCMultiIndices.hpp:322
MoFEM::CubitMeshSets_change_attributes::operator()
void operator()(CubitMeshSets &e)
Definition: BCMultiIndices.cpp:419
MAT_INTERFSET
@ MAT_INTERFSET
Definition: definitions.h:173
MoFEM::CubitMeshSets::getTypeFromBcData
MoFEMErrorCode getTypeFromBcData(const std::vector< char > &bc_data, CubitBCType &type) const
Function that returns the CubitBCType type of the contents of bc_data.
Definition: BCMultiIndices.cpp:230
MoFEM::GenericAttributeData::getDataPtr
virtual const void * getDataPtr() const =0
get pointer to data structure
HEATFLUXSET
@ HEATFLUXSET
Definition: definitions.h:169
MoFEM::CubitMeshSets::ssTag_data
Tag ssTag_data
Definition: BCMultiIndices.hpp:322
MoFEM::CubitMeshSets::getMeshsetIdEntitiesByDimension
MoFEMErrorCode getMeshsetIdEntitiesByDimension(Interface &moab, const int dimension, Range &entities, const bool recursive=false) const
get entities form meshset
convert.int
int
Definition: convert.py:64
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:429
MoFEM::CubitMeshSets::getBlockHeaderData
MoFEMErrorCode getBlockHeaderData(std::vector< unsigned int > &material_data) const
get block_headers vector from MoFEM database
Definition: BCMultiIndices.cpp:183
MAT_MOISTURESET
@ MAT_MOISTURESET
block name is "MAT_MOISTURE"
Definition: definitions.h:176
MoFEM::CubitMeshSets::tagBcData
char * tagBcData
Definition: BCMultiIndices.hpp:27
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:359
MoFEM::CubitMeshSets::printBcData
MoFEMErrorCode printBcData(std::ostream &os) const
print bc_data int stream given by os
Definition: BCMultiIndices.cpp:271