v0.13.0
Public Member Functions | Public Attributes | Friends | List of all members
MoFEM::CubitMeshSets Struct Reference

this struct keeps basic methods for moab meshset about material and boundary conditions More...

#include <src/multi_indices/BCMultiIndices.hpp>

Collaboration diagram for MoFEM::CubitMeshSets:
[legend]

Public Member Functions

 CubitMeshSets (Interface &moab, const EntityHandle meshset)
 
 CubitMeshSets (Interface &moab, const CubitBCType cubit_bc_type, const int msId)
 
int getMeshsetId () const
 get meshset id as it set in preprocessing software More...
 
CubitBCType getBcType () const
 get type of meshset More...
 
EntityHandle getMeshset () const
 get bc meshset More...
 
unsigned long int getBcTypeULong () const
 get bc meshset type More...
 
unsigned long int getMaksedBcTypeULong () const
 get meshset type and mask More...
 
unsigned int getMeshsetEntitiesDimension () const
 Get the meshset entities dimension. More...
 
MoFEMErrorCode getMeshsetIdEntitiesByDimension (Interface &moab, const int dimension, Range &entities, const bool recursive=false) const
 get entities form meshset More...
 
MoFEMErrorCode getMeshsetIdEntitiesByDimension (Interface &moab, Range &entities, const bool recursive=false) const
 get entities form meshset More...
 
MoFEMErrorCode getMeshsetIdEntitiesByType (Interface &moab, const EntityType type, Range &entities, const bool recursive=false) const
 get entities by type More...
 
MoFEMErrorCode getTypeFromBcData (const std::vector< char > &bc_data, CubitBCType &type) const
 Function that returns the CubitBCType type of the contents of bc_data. More...
 
MoFEMErrorCode getTypeFromBcData (CubitBCType &type) const
 Function that returns the CubitBCType type of the contents of bc_data. More...
 
MoFEMErrorCode getBcData (std::vector< char > &bc_data) const
 get bc_data vector from MoFEM database More...
 
MoFEMErrorCode getBlockHeaderData (std::vector< unsigned int > &material_data) const
 get block_headers vector from MoFEM database More...
 
MoFEMErrorCode printBlockHeaderData (std::ostream &os) const
 print material_data int stream given by os More...
 
MoFEMErrorCode printBcData (std::ostream &os) const
 print bc_data int stream given by os More...
 
MoFEMErrorCode getTypeFromName (const std::string &name, CubitBCType &type) const
 Function that returns the CubitBCType type of the block name, sideset name etc. More...
 
MoFEMErrorCode getTypeFromName (CubitBCType &type) const
 Function that returns the CubitBCType type of the block name, sideset name etc. More...
 
MoFEMErrorCode getAttributes (std::vector< double > &attributes) const
 get Cubit block attributes More...
 
MoFEMErrorCode setAttributes (moab::Interface &moab, const std::vector< double > &attributes)
 cet Cubit block attributes More...
 
MoFEMErrorCode printAttributes (std::ostream &os) const
 print the attributes vector More...
 
std::string getName () const
 get name of block, sideset etc. (this is set in Cubit block properties) More...
 
MoFEMErrorCode printName (std::ostream &os) const
 print name of block, sideset etc. (this is set in Cubit setting properties) More...
 
template<class ATTRIBUTE_TYPE >
MoFEMErrorCode getAttributeDataStructure (ATTRIBUTE_TYPE &data) const
 fill data structure with data saved on meshset More...
 
template<class ATTRIBUTE_TYPE >
MoFEMErrorCode setAttributeDataStructure (const ATTRIBUTE_TYPE &data)
 fill meshset data with data on structure More...
 
template<class CUBIT_BC_DATA_TYPE >
MoFEMErrorCode getBcDataStructure (CUBIT_BC_DATA_TYPE &data) const
 
template<class CUBIT_BC_DATA_TYPE >
MoFEMErrorCode setBcDataStructure (CUBIT_BC_DATA_TYPE &data)
 
MoFEMErrorCode getTagsHandlers (Interface &moab)
 

Public Attributes

EntityHandle meshset
 
CubitBCType cubitBcType
 
std::vector< Tag > tag_handles
 vector of tag handles to types of data passed from cubit More...
 
int * msId
 cubit meshset ID More...
 
char * tagBcData
 
int tagBcSize
 
unsigned int * tagBlockHeaderData
 
double * tagBlockAttributes
 
int tagBlockAttributesSize
 
char * tagName
 
const CubitBCType meshsetsMask
 
Tag nsTag
 
Tag ssTag
 
Tag nsTag_data
 
Tag ssTag_data
 
Tag bhTag
 
Tag bhTag_header
 
Tag thBlockAttribs
 
Tag entityNameTag
 

Friends

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

Detailed Description

this struct keeps basic methods for moab meshset about material and boundary conditions

Definition at line 29 of file BCMultiIndices.hpp.

Constructor & Destructor Documentation

◆ CubitMeshSets() [1/2]

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

◆ CubitMeshSets() [2/2]

MoFEM::CubitMeshSets::CubitMeshSets ( Interface moab,
const CubitBCType  cubit_bc_type,
const int  msId 
)

Member Function Documentation

◆ getAttributeDataStructure()

template<class ATTRIBUTE_TYPE >
MoFEMErrorCode MoFEM::CubitMeshSets::getAttributeDataStructure ( ATTRIBUTE_TYPE &  data) const

fill data structure with data saved on meshset

Definition at line 273 of file BCMultiIndices.hpp.

273  {
275  if ((cubitBcType & data.getType()).none()) {
276  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
277  "attributes are not for ATTRIBUTE_TYPE structure");
278  }
279  std::vector<double> attributes;
280  CHKERR getAttributes(attributes);
281  CHKERR data.fill_data(attributes);
283  }
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:359
@ MOFEM_DATA_INCONSISTENCY
Definition: definitions.h:44
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:429
#define CHKERR
Inline error check.
Definition: definitions.h:548
MoFEMErrorCode getAttributes(std::vector< double > &attributes) const
get Cubit block attributes

◆ getAttributes()

MoFEMErrorCode MoFEM::CubitMeshSets::getAttributes ( std::vector< double > &  attributes) const

get Cubit block attributes

Parameters
attributesis the vector where the block attribute data will be stored

Definition at line 287 of file BCMultiIndices.cpp.

287  {
289  attributes.resize(tagBlockAttributesSize);
290  if (tagBlockAttributesSize > 0) {
292  attributes.begin());
293  }
295 }

◆ getBcData()

MoFEMErrorCode MoFEM::CubitMeshSets::getBcData ( std::vector< char > &  bc_data) const

get bc_data vector from MoFEM database

Parameters
bc_datais the in/out vector were bc_data will be stored

Definition at line 176 of file BCMultiIndices.cpp.

176  {
178  bc_data.resize(tagBcSize);
179  copy(&tagBcData[0], &tagBcData[tagBcSize], bc_data.begin());
181 }

◆ getBcDataStructure()

template<class CUBIT_BC_DATA_TYPE >
MoFEMErrorCode MoFEM::CubitMeshSets::getBcDataStructure ( CUBIT_BC_DATA_TYPE &  data) const

Definition at line 301 of file BCMultiIndices.hpp.

301  {
303 
304  if ((cubitBcType & data.tYpe).none()) {
305  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
306  "bc_data are not for CUBIT_BC_DATA_TYPE structure");
307  }
308  std::vector<char> bc_data;
309  getBcData(bc_data);
310  ierr = data.fill_data(bc_data);
311  CHKERRG(ierr);
313  }
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:460
#define CHKERRG(n)
Check error code of MoFEM/MOAB/PETSc function.
Definition: definitions.h:496
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:453
static MoFEMErrorCodeGeneric< PetscErrorCode > ierr
Definition: Exceptions.hpp:87
MoFEMErrorCode getBcData(std::vector< char > &bc_data) const
get bc_data vector from MoFEM database

◆ getBcType()

CubitBCType MoFEM::CubitMeshSets::getBcType ( ) const

get type of meshset

See CubitBC for set of types of meshsets.

Returns
meshset type

Definition at line 66 of file BCMultiIndices.hpp.

66 { return cubitBcType; }

◆ getBcTypeULong()

unsigned long int MoFEM::CubitMeshSets::getBcTypeULong ( ) const

get bc meshset type

Returns
return type as unsigned int

Definition at line 78 of file BCMultiIndices.hpp.

78  {
79  return cubitBcType.to_ulong();
80  }

◆ getBlockHeaderData()

MoFEMErrorCode MoFEM::CubitMeshSets::getBlockHeaderData ( std::vector< unsigned int > &  material_data) const

get block_headers vector from MoFEM database

Parameters
material_datais the in/out vector were the material data will be stored

Definition at line 183 of file BCMultiIndices.cpp.

184  {
186  material_data.resize(3);
187  copy(&tagBlockHeaderData[0], &tagBlockHeaderData[3], material_data.begin());
189 }
unsigned int * tagBlockHeaderData

◆ getMaksedBcTypeULong()

unsigned long int MoFEM::CubitMeshSets::getMaksedBcTypeULong ( ) const

get meshset type and mask

Returns
type is returned as unsigned integer

Definition at line 86 of file BCMultiIndices.hpp.

86  {
87  return (cubitBcType & meshsetsMask).to_ulong();
88  }
const CubitBCType meshsetsMask

◆ getMeshset()

EntityHandle MoFEM::CubitMeshSets::getMeshset ( ) const

get bc meshset

Returns
meshset entity handle

Definition at line 72 of file BCMultiIndices.hpp.

72 { return meshset; }

◆ getMeshsetEntitiesDimension()

unsigned int MoFEM::CubitMeshSets::getMeshsetEntitiesDimension ( ) const

Get the meshset entities dimension.

Note
If dimension is -1, then dimension for meshset ins undetermined.
Returns
unsigned int

Definition at line 97 of file BCMultiIndices.hpp.

97  {
99  return tagBlockHeaderData[2];
100  else
101  return -1;
102  }

◆ getMeshsetId()

int MoFEM::CubitMeshSets::getMeshsetId ( ) const

get meshset id as it set in preprocessing software

Returns
id of meshset

Definition at line 53 of file BCMultiIndices.hpp.

53 { return *msId; }
int * msId
cubit meshset ID

◆ getMeshsetIdEntitiesByDimension() [1/2]

MoFEMErrorCode MoFEM::CubitMeshSets::getMeshsetIdEntitiesByDimension ( Interface moab,
const int  dimension,
Range &  entities,
const bool  recursive = false 
) const

get entities form meshset

Parameters
moabmoab instance
dimensiondimension of entities
entitiesrange of returned entities
recursivetrue if meshset should be searched recursively
Returns
error code

◆ getMeshsetIdEntitiesByDimension() [2/2]

MoFEMErrorCode MoFEM::CubitMeshSets::getMeshsetIdEntitiesByDimension ( Interface moab,
Range &  entities,
const bool  recursive = false 
) const

get entities form meshset

Use if meshset have predefined dimension

Parameters
moabmoab instance
entitiesrange of returned entities
recursivetrue if meshset should be searched recursively
Returns
error code

◆ getMeshsetIdEntitiesByType()

MoFEMErrorCode MoFEM::CubitMeshSets::getMeshsetIdEntitiesByType ( Interface moab,
const EntityType  type,
Range &  entities,
const bool  recursive = false 
) const

get entities by type

Parameters
moabmoab instance
typetype of entity
entitiesreturned entities
recursivetrue if meshset should be searched recursively
Returns
error code

Definition at line 162 of file BCMultiIndices.cpp.

164  {
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 }
static MoFEMErrorCodeGeneric< moab::ErrorCode > rval
Definition: Exceptions.hpp:85

◆ getName()

std::string MoFEM::CubitMeshSets::getName ( ) const

get name of block, sideset etc. (this is set in Cubit block properties)

Block Name Conventions:

Materials are defined with block names starting with MAT_ e.g. MAT_ELASTIC_abcd.

List of materials/solution procedures

Block name / Number of attributes / (1) Attribute 1, (2) Attribute 2 etc.

MAT_ELASTIC / 10 / (1) Young's modulus (2) Poisson's ratio (3) User attribute 8 ... (10) User attribute 8

MAT_ELASTIC_TRANSISO / 5 / (1) Young's modulus in xy plane (Ep) (2) Young's modulus in z-direction (Ez) (3) Poisson's ratio in xy plane (vp) (4) Poisson's ratio in z-direction (vpz) (5) Shear modulus in z-direction (Gzp)

MAT_INTERF / 1 / (1) Elastic modulus multiplier

To be extended as appropriate

Definition at line 213 of file BCMultiIndices.cpp.

213  {
214  if (tagName != nullptr) {
215  return std::string(tagName);
216  } else {
217  return "NoNameSet";
218  }
219 }

◆ getTagsHandlers()

MoFEMErrorCode MoFEM::CubitMeshSets::getTagsHandlers ( Interface moab)

Definition at line 22 of file BCMultiIndices.cpp.

22  {
24  CHKERR moab.tag_get_handle(DIRICHLET_SET_TAG_NAME, nsTag);
25  CHKERR moab.tag_get_handle(NEUMANN_SET_TAG_NAME, ssTag);
26  CHKERR moab.tag_get_handle(
27  (std::string(DIRICHLET_SET_TAG_NAME) + "__BC_DATA").c_str(), nsTag_data);
28  CHKERR moab.tag_get_handle(
29  (std::string(NEUMANN_SET_TAG_NAME) + "__BC_DATA").c_str(), ssTag_data);
30  CHKERR moab.tag_get_handle(MATERIAL_SET_TAG_NAME, bhTag);
31  CHKERR moab.tag_get_handle(BLOCK_HEADER, bhTag_header);
32  CHKERR moab.tag_get_handle(BLOCK_ATTRIBUTES, thBlockAttribs);
33  CHKERR moab.tag_get_handle(NAME_TAG_NAME, entityNameTag);
35 }

◆ getTypeFromBcData() [1/2]

MoFEMErrorCode MoFEM::CubitMeshSets::getTypeFromBcData ( const std::vector< char > &  bc_data,
CubitBCType type 
) const

Function that returns the CubitBCType type of the contents of bc_data.

Definition at line 230 of file BCMultiIndices.cpp.

231  {
233 
234  // See CubitBCType in common.hpp
235  if (bc_data.size() == 0) {
237  }
238 
239  if (strcmp(&bc_data[0], "Displacement") == 0)
241  else if (strcmp(&bc_data[0], "Force") == 0)
242  type |= FORCESET;
243  else if (strcmp(&bc_data[0], "Velocity") == 0)
244  type |= VELOCITYSET;
245  else if (strcmp(&bc_data[0], "Acceleration") == 0)
247  else if (strcmp(&bc_data[0], "Temperature") == 0)
248  type |= TEMPERATURESET;
249  else if (strcmp(&bc_data[0], "Pressure") == 0)
250  type |= PRESSURESET;
251  else if (strcmp(&bc_data[0], "HeatFlux") == 0)
252  type |= HEATFLUXSET;
253  else if (strcmp(&bc_data[0], "cfd_bc") == 0)
254  type |= INTERFACESET;
255  else
256  type |= UNKNOWNNAME;
257 
259 }
@ TEMPERATURESET
Definition: definitions.h:168
@ PRESSURESET
Definition: definitions.h:165
@ ACCELERATIONSET
Definition: definitions.h:167
@ FORCESET
Definition: definitions.h:164
@ HEATFLUXSET
Definition: definitions.h:169
@ UNKNOWNNAME
Definition: definitions.h:171
@ VELOCITYSET
Definition: definitions.h:166
@ DISPLACEMENTSET
Definition: definitions.h:163
@ INTERFACESET
Definition: definitions.h:170

◆ getTypeFromBcData() [2/2]

MoFEMErrorCode MoFEM::CubitMeshSets::getTypeFromBcData ( CubitBCType type) const

Function that returns the CubitBCType type of the contents of bc_data.

Definition at line 260 of file BCMultiIndices.cpp.

260  {
262  std::vector<char> bc_data;
263  CHKERR getBcData(bc_data);
264  CHKERR getTypeFromBcData(bc_data, type);
266 }
MoFEMErrorCode getTypeFromBcData(const std::vector< char > &bc_data, CubitBCType &type) const
Function that returns the CubitBCType type of the contents of bc_data.

◆ getTypeFromName() [1/2]

MoFEMErrorCode MoFEM::CubitMeshSets::getTypeFromName ( const std::string &  name,
CubitBCType type 
) const

Function that returns the CubitBCType type of the block name, sideset name etc.

Definition at line 325 of file BCMultiIndices.cpp.

326  {
328  // See CubitBCType in common.hpp
329  if (name.compare(0, 11, "MAT_ELASTIC") == 0) {
330  type |= MAT_ELASTICSET;
331  } else if (name.compare(0, 11, "MAT_THERMAL") == 0) {
332  type |= MAT_THERMALSET;
333  } else if (name.compare(0, 12, "MAT_MOISTURE") == 0) {
335  } else if (name.compare(0, 10, "MAT_INTERF") == 0) {
336  type |= MAT_INTERFSET;
337  } else if (name.compare(0, 11, "BODY_FORCES") == 0) {
338  type |= BODYFORCESSET;
339  }
340  // To be extended as appropriate
341  else {
342  type |= UNKNOWNNAME;
343  }
345 }
@ BODYFORCESSET
block name is "BODY_FORCES"
Definition: definitions.h:175
@ MAT_ELASTICSET
block name is "MAT_ELASTIC"
Definition: definitions.h:172
@ MAT_INTERFSET
Definition: definitions.h:173
@ MAT_THERMALSET
block name is "MAT_THERMAL"
Definition: definitions.h:174
@ MAT_MOISTURESET
block name is "MAT_MOISTURE"
Definition: definitions.h:176

◆ getTypeFromName() [2/2]

MoFEMErrorCode MoFEM::CubitMeshSets::getTypeFromName ( CubitBCType type) const

Function that returns the CubitBCType type of the block name, sideset name etc.

Definition at line 347 of file BCMultiIndices.cpp.

347  {
349  std::string name = getName();
350  CHKERR getTypeFromName(name, type);
352 }
MoFEMErrorCode getTypeFromName(const std::string &name, CubitBCType &type) const
Function that returns the CubitBCType type of the block name, sideset name etc.
std::string getName() const
get name of block, sideset etc. (this is set in Cubit block properties)

◆ printAttributes()

MoFEMErrorCode MoFEM::CubitMeshSets::printAttributes ( std::ostream &  os) const

print the attributes vector

f.e. it->printAttributes(cout), i.e. printing to standard output f.e. it->printAttributes(std::cerr), i.e. printing to standard error output

Definition at line 311 of file BCMultiIndices.cpp.

311  {
313  std::vector<double> attributes;
314  CHKERR getAttributes(attributes);
315  os << std::endl;
316  os << "Block attributes" << std::endl;
317  os << "----------------" << std::endl;
318  for (unsigned int ii = 0; ii < attributes.size(); ii++) {
319  os << "attr. no: " << ii + 1 << " value: " << attributes[ii] << std::endl;
320  }
321  os << std::endl;
323 }

◆ printBcData()

MoFEMErrorCode MoFEM::CubitMeshSets::printBcData ( std::ostream &  os) const

print bc_data int stream given by os

f.e. it->printBcData(cout), i.e. printing to standard output f.e. it->printBcData(std::cerr), i.e. printing to standard error output

Definition at line 268 of file BCMultiIndices.cpp.

268  {
270  std::vector<char> bc_data;
271  CHKERR getBcData(bc_data);
272  os << "bc_data = ";
273  std::vector<char>::iterator vit = bc_data.begin();
274  for (; vit != bc_data.end(); vit++) {
275  os << std::hex << (int)((unsigned char)*vit) << " ";
276  }
277  os << ": ";
278  vit = bc_data.begin();
279  for (; vit != bc_data.end(); vit++) {
280  os << *vit;
281  }
282  os << std::endl;
284 }

◆ printBlockHeaderData()

MoFEMErrorCode MoFEM::CubitMeshSets::printBlockHeaderData ( std::ostream &  os) const

print material_data int stream given by os

f.e. it->print_Cubit_material_data(cout), i.e. printing to standard output f.e. it->print_Cubit_material_data(std::cerr), i.e. printing to standard error output

Definition at line 191 of file BCMultiIndices.cpp.

191  {
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 }
MoFEMErrorCode getBlockHeaderData(std::vector< unsigned int > &material_data) const
get block_headers vector from MoFEM database

◆ printName()

MoFEMErrorCode MoFEM::CubitMeshSets::printName ( std::ostream &  os) const

print name of block, sideset etc. (this is set in Cubit setting properties)

e.g. it->printName(cout), i.e. printing to standard output e.g it->printName(std::cerr), i.e. printing to standard error output

Definition at line 221 of file BCMultiIndices.cpp.

221  {
223  std::string name = getName();
224  os << std::endl;
225  os << "Block name: " << name << std::endl;
227 }

◆ setAttributeDataStructure()

template<class ATTRIBUTE_TYPE >
MoFEMErrorCode MoFEM::CubitMeshSets::setAttributeDataStructure ( const ATTRIBUTE_TYPE &  data)

fill meshset data with data on structure

Definition at line 289 of file BCMultiIndices.hpp.

289  {
291  if ((cubitBcType & data.getType()).none()) {
292  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
293  "attributes are not for ATTRIBUTE_TYPE structure");
294  }
295  double *ptr = const_cast<double *>(tagBlockAttributes);
296  CHKERR data.set_data(ptr, 8 * tagBlockAttributesSize);
298  }

◆ setAttributes()

MoFEMErrorCode MoFEM::CubitMeshSets::setAttributes ( moab::Interface &  moab,
const std::vector< double > &  attributes 
)

cet Cubit block attributes

Parameters
attributesis the vector where the block attribute data will be stored

Definition at line 298 of file BCMultiIndices.cpp.

299  {
300 
302  int tag_size[] = {(int)attributes.size()};
303  void const *tag_data[] = {&*attributes.begin()};
304  CHKERR moab.tag_set_by_ptr(thBlockAttribs, &meshset, 1, tag_data, tag_size);
305  CHKERR moab.tag_get_by_ptr(thBlockAttribs, &meshset, 1,
306  (const void **)&tagBlockAttributes,
309 }

◆ setBcDataStructure()

template<class CUBIT_BC_DATA_TYPE >
MoFEMErrorCode MoFEM::CubitMeshSets::setBcDataStructure ( CUBIT_BC_DATA_TYPE &  data)

Definition at line 316 of file BCMultiIndices.hpp.

316  {
318 
319  char *ptr = const_cast<char *>(tagBcData);
320  ierr = data.set_data(ptr, tagBcSize);
321  CHKERRG(ierr);
323  }

Friends And Related Function Documentation

◆ operator<<

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

Definition at line 354 of file BCMultiIndices.cpp.

354  {
355  // get name of cubit meshset
356  std::ostringstream ss;
357  unsigned jj = 0;
358  while (1 << jj != LASTSET_BC) {
359  const CubitBCType jj_bc_type = 1 << jj;
360  if ((e.cubitBcType & jj_bc_type).any()) {
361  string bc_type_name;
362  ss << " " << string(CubitBCNames[jj + 1]);
363  }
364  ++jj;
365  }
366 
367  // push data to stream
368  os << "meshset " << e.meshset << " type" << ss.str();
369  if (e.msId != nullptr)
370  os << " msId " << *(e.msId);
371  if (e.tagName != nullptr) {
372  os << " name " << e.getName();
373  }
374  if (e.tagBlockHeaderData != nullptr) {
375  os << " block header: ";
376  os << " blockCol = " << e.tagBlockHeaderData[0];
377  os << " blockMat = " << e.tagBlockHeaderData[1];
378  os << " blockDimension = " << e.tagBlockHeaderData[2];
379  }
380  return os;
381 }
@ LASTSET_BC
Definition: definitions.h:179
static const char *const CubitBCNames[]
Names of types of sets and boundary conditions.
Definition: definitions.h:188
std::bitset< 32 > CubitBCType
Definition: Types.hpp:63

Member Data Documentation

◆ bhTag

Tag MoFEM::CubitMeshSets::bhTag

Definition at line 327 of file BCMultiIndices.hpp.

◆ bhTag_header

Tag MoFEM::CubitMeshSets::bhTag_header

Definition at line 327 of file BCMultiIndices.hpp.

◆ cubitBcType

CubitBCType MoFEM::CubitMeshSets::cubitBcType

type of meshset from cubit NodeSet, BlockSet, SideSet and more

Definition at line 32 of file BCMultiIndices.hpp.

◆ entityNameTag

Tag MoFEM::CubitMeshSets::entityNameTag

Definition at line 328 of file BCMultiIndices.hpp.

◆ meshset

EntityHandle MoFEM::CubitMeshSets::meshset

Definition at line 31 of file BCMultiIndices.hpp.

◆ meshsetsMask

const CubitBCType MoFEM::CubitMeshSets::meshsetsMask

Definition at line 43 of file BCMultiIndices.hpp.

◆ msId

int* MoFEM::CubitMeshSets::msId

cubit meshset ID

Definition at line 36 of file BCMultiIndices.hpp.

◆ nsTag

Tag MoFEM::CubitMeshSets::nsTag

Definition at line 327 of file BCMultiIndices.hpp.

◆ nsTag_data

Tag MoFEM::CubitMeshSets::nsTag_data

Definition at line 327 of file BCMultiIndices.hpp.

◆ ssTag

Tag MoFEM::CubitMeshSets::ssTag

Definition at line 327 of file BCMultiIndices.hpp.

◆ ssTag_data

Tag MoFEM::CubitMeshSets::ssTag_data

Definition at line 327 of file BCMultiIndices.hpp.

◆ tag_handles

std::vector<Tag> MoFEM::CubitMeshSets::tag_handles

vector of tag handles to types of data passed from cubit

Definition at line 35 of file BCMultiIndices.hpp.

◆ tagBcData

char* MoFEM::CubitMeshSets::tagBcData

Definition at line 37 of file BCMultiIndices.hpp.

◆ tagBcSize

int MoFEM::CubitMeshSets::tagBcSize

Definition at line 38 of file BCMultiIndices.hpp.

◆ tagBlockAttributes

double* MoFEM::CubitMeshSets::tagBlockAttributes

Definition at line 40 of file BCMultiIndices.hpp.

◆ tagBlockAttributesSize

int MoFEM::CubitMeshSets::tagBlockAttributesSize

Definition at line 41 of file BCMultiIndices.hpp.

◆ tagBlockHeaderData

unsigned int* MoFEM::CubitMeshSets::tagBlockHeaderData

Definition at line 39 of file BCMultiIndices.hpp.

◆ tagName

char* MoFEM::CubitMeshSets::tagName

Definition at line 42 of file BCMultiIndices.hpp.

◆ thBlockAttribs

Tag MoFEM::CubitMeshSets::thBlockAttribs

Definition at line 327 of file BCMultiIndices.hpp.


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