v0.14.0
BitRefManager.cpp
Go to the documentation of this file.
1 /** \file BitRefManager.cpp
2  * \brief Managing BitRefLevels
3  * \mofem_bit_ref
4  */
5 
6 namespace MoFEM {
7 
9 BitRefManager::query_interface(boost::typeindex::type_index type_index,
10  UnknownInterface **iface) const {
12  *iface = const_cast<BitRefManager *>(this);
14 }
15 
17  : cOre(const_cast<MoFEM::Core &>(core)), dEbug(false) {
18 
19  if (!LogManager::checkIfChannelExist("BitRefSelf")) {
20  auto core_log = logging::core::get();
21  core_log->add_sink(
23  LogManager::setLog("BitRefSelf");
24  core_log->add_sink(
26  LogManager::setLog("BitRefWorld");
27  core_log->add_sink(
29  LogManager::setLog("BitRefSync");
30  MOFEM_LOG_TAG("BitRefSelf", "BitRefManager");
31  MOFEM_LOG_TAG("BitRefWorld", "BitRefManager");
32  MOFEM_LOG_TAG("BitRefSync", "BitRefManager");
33  }
34 
36  MOFEM_LOG("BitRefWorld", Sev::noisy) << "BitRefManager interface created";
37 }
38 
39 /// tool class with methods used more than twp times
41 
43 
44  const BitRefLevel &bIt; ///< bit to set
45  const RefEntity_multiIndex *refEntsPtr; ///< access to ents database
46  const RefElement_multiIndex *refElementPtr; ///< access to fe database
47 
48  boost::shared_ptr<BasicEntityData> &baseEntData; ///< base entity data
49 
50  /// constrictor
52  const RefEntity_multiIndex *ref_ents_ptr,
53  const RefElement_multiIndex *ref_element_ptr)
54  : mField(m_field), bIt(bit), refEntsPtr(ref_ents_ptr),
55  refElementPtr(ref_element_ptr),
56  baseEntData(m_field.get_basic_entity_data_ptr()) {}
57 
58  /// find entities and change entity bit if in database
60  Range &seed_ents_range) const {
62 
63  seed_ents_range.insert(f, s);
64  // get lower bound of multi-index
65  auto rit = refEntsPtr->lower_bound(f);
66  if (rit == refEntsPtr->end()) {
67  // all enties in range are added to database
69  } else {
70 
71  // some entities from range are in database
72  auto hi_rit = refEntsPtr->upper_bound(s);
73 
74  Range to_erase;
75  insertOrdered(to_erase, RefEntExtractor(), rit, hi_rit);
76  if (bIt.any())
77  for (; rit != hi_rit; ++rit)
78  const_cast<BitRefLevel &>((*rit)->getBitRefLevel()) |= bIt;
79 
80  Range::iterator lo, hi = seed_ents_range.begin();
81  for (auto pt = to_erase.pair_begin(); pt != to_erase.pair_end(); ++pt) {
82  lo = seed_ents_range.lower_bound(hi, seed_ents_range.end(), pt->first);
83  if (lo != seed_ents_range.end()) {
84  hi = seed_ents_range.upper_bound(lo, seed_ents_range.end(),
85  pt->second);
86  seed_ents_range.erase(lo, hi);
87  } else
88  break;
89  }
90  }
91 
93  }
94 
95  /// add entities to database
96  template <int N>
97  MoFEMErrorCode addEntsToDatabaseImpl(const Range &seed_ents_range) const {
99  std::vector<boost::shared_ptr<RefEntity>> shared_ref_ents_vec;
100  shared_ref_ents_vec.reserve(seed_ents_range.size());
101  std::vector<const void *> tag_by_ptr;
102  for (Range::const_pair_iterator pit = seed_ents_range.pair_begin();
103  pit != seed_ents_range.pair_end(); pit++) {
104  // add entities to database
105  EntityHandle f = pit->first;
106  EntityHandle s = pit->second;
107 
108  boost::shared_ptr<std::vector<RefEntityTmp<N>>> ref_ents_vec(
109  new std::vector<RefEntityTmp<N>>());
110  ref_ents_vec->reserve(s - f + 1);
111 
112  tag_by_ptr.resize(s - f + 1);
113  CHKERR baseEntData->moab.tag_get_by_ptr(
114  baseEntData->th_RefParentHandle, Range(f, s), &*tag_by_ptr.begin());
115  auto tag_parent_it = tag_by_ptr.begin();
116  for (auto f : Range(f, s)) {
117  ref_ents_vec->emplace_back(
118  baseEntData, f,
119  const_cast<EntityHandle *>(
120  static_cast<const EntityHandle *>(*tag_parent_it)));
121  ++tag_parent_it;
122  }
123 
124  // Set bits to range
125  if (bIt.any()) {
126  tag_by_ptr.resize(s - f + 1);
127  CHKERR baseEntData->moab.tag_get_by_ptr(
128  baseEntData->th_RefBitLevel, Range(f, s), &*tag_by_ptr.begin());
129  for (auto &v_bit_ptr : tag_by_ptr)
130  const_cast<BitRefLevel &>(
131  *(static_cast<const BitRefLevel *>(v_bit_ptr))) |= bIt;
132  }
133 
134  for (auto &re : *ref_ents_vec)
135  shared_ref_ents_vec.emplace_back(ref_ents_vec,
136  static_cast<RefEntity *>(&re));
137  }
138  if (!shared_ref_ents_vec.empty()) {
139  int s0 = refEntsPtr->size();
140  const_cast<RefEntity_multiIndex *>(refEntsPtr)
141  ->insert(shared_ref_ents_vec.begin(), shared_ref_ents_vec.end());
142  if ((refEntsPtr->size() - s0) != shared_ref_ents_vec.size()) {
143  SETERRQ2(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
144  "Data inconsistency %d != %d", refEntsPtr->size() - s0,
145  shared_ref_ents_vec.size());
146  }
147  }
149  }
150 
151  MoFEMErrorCode addEntsToDatabase(const Range &seed_ents_range) const {
153 
154  switch (mField.getValue()) {
155  case -1:
156  return addEntsToDatabaseImpl<-1>(seed_ents_range);
157  case 0:
158  return addEntsToDatabaseImpl<0>(seed_ents_range);
159  case 1:
160  return addEntsToDatabaseImpl<1>(seed_ents_range);
161  default:
162  SETERRQ1(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
163  "Core index can vary from -1 to %d", MAX_CORE_TMP);
164  }
165 
167  }
168 
170  Range &seed_fe_range) const {
172  seed_fe_range.insert(f, s);
173  RefElement_multiIndex::iterator rit, hi_rit;
174  // get lower bound of multi-index
175  rit = refElementPtr->lower_bound(f);
176  if (rit == refElementPtr->end()) {
177  // all enties in range are added to database
179  } else {
180  // some entities from range are in database
181  hi_rit = refElementPtr->upper_bound(s);
182  for (; rit != hi_rit; ++rit) {
183  seed_fe_range.erase(rit->get()->getEnt());
184  }
185  }
187  }
188 
191  std::vector<boost::shared_ptr<RefElement>> shared_ref_fe_vec;
192  shared_ref_fe_vec.reserve(seed_fe_range.size());
193  // create ref entity instances
194  for (Range::const_pair_iterator pit = seed_fe_range.const_pair_begin();
195  pit != seed_fe_range.const_pair_end(); ++pit) {
196  RefEntity_multiIndex::iterator rit, hi_rit;
197  rit = refEntsPtr->lower_bound(pit->first);
198  hi_rit = refEntsPtr->upper_bound(pit->second);
199  if (static_cast<int>(std::distance(rit, hi_rit)) !=
200  static_cast<int>(pit->second - pit->first + 1)) {
201  SETERRQ2(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
202  "data inconsistency %d != %d", std::distance(rit, hi_rit),
203  pit->second - pit->first + 1);
204  }
205  switch ((*rit)->getEntType()) {
206  case MBVERTEX: {
207  boost::shared_ptr<std::vector<RefElement_VERTEX>> ref_fe_vec =
208  boost::make_shared<std::vector<RefElement_VERTEX>>();
209  ref_fe_vec->reserve(pit->second - pit->first + 1);
210  for (; rit != hi_rit; ++rit) {
211  ref_fe_vec->push_back(RefElement_VERTEX(*rit));
212  shared_ref_fe_vec.push_back(
213  boost::shared_ptr<RefElement>(ref_fe_vec, &ref_fe_vec->back()));
214  }
215  } break;
216  case MBEDGE: {
217  boost::shared_ptr<std::vector<RefElement_EDGE>> ref_fe_vec =
218  boost::make_shared<std::vector<RefElement_EDGE>>();
219  ref_fe_vec->reserve(pit->second - pit->first + 1);
220  for (; rit != hi_rit; ++rit) {
221  ref_fe_vec->push_back(RefElement_EDGE(*rit));
222  shared_ref_fe_vec.push_back(
223  boost::shared_ptr<RefElement>(ref_fe_vec, &ref_fe_vec->back()));
224  }
225  } break;
226  case MBTRI:
227  case MBQUAD: {
228  boost::shared_ptr<std::vector<RefElementFace>> ref_fe_vec =
229  boost::make_shared<std::vector<RefElementFace>>();
230  ref_fe_vec->reserve(pit->second - pit->first + 1);
231  for (; rit != hi_rit; ++rit) {
232  ref_fe_vec->push_back(RefElementFace(*rit));
233  shared_ref_fe_vec.push_back(
234  boost::shared_ptr<RefElement>(ref_fe_vec, &ref_fe_vec->back()));
235  }
236  } break;
237  case MBTET:
238  case MBHEX: {
239  boost::shared_ptr<std::vector<RefElementVolume>> ref_fe_vec =
240  boost::make_shared<std::vector<RefElementVolume>>();
241  ref_fe_vec->reserve(pit->second - pit->first + 1);
242  for (; rit != hi_rit; ++rit) {
243  ref_fe_vec->push_back(RefElementVolume(*rit));
244  shared_ref_fe_vec.push_back(
245  boost::shared_ptr<RefElement>(ref_fe_vec, &ref_fe_vec->back()));
246  }
247  } break;
248  case MBPRISM: {
249  boost::shared_ptr<std::vector<RefElement_PRISM>> ref_fe_vec =
250  boost::make_shared<std::vector<RefElement_PRISM>>();
251  ref_fe_vec->reserve(pit->second - pit->first + 1);
252  for (; rit != hi_rit; ++rit) {
253  ref_fe_vec->push_back(RefElement_PRISM(*rit));
254  shared_ref_fe_vec.push_back(
255  boost::shared_ptr<RefElement>(ref_fe_vec, &ref_fe_vec->back()));
256  }
257  } break;
258  case MBENTITYSET: {
259  boost::shared_ptr<std::vector<RefElement_MESHSET>> ref_fe_vec =
260  boost::make_shared<std::vector<RefElement_MESHSET>>();
261  ref_fe_vec->reserve(pit->second - pit->first + 1);
262  for (; rit != hi_rit; ++rit) {
263  ref_fe_vec->push_back(RefElement_MESHSET(*rit));
264  shared_ref_fe_vec.push_back(
265  boost::shared_ptr<RefElement>(ref_fe_vec, &ref_fe_vec->back()));
266  }
267  } break;
268  default:
269  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "not implemented");
270  }
271  }
272  // add shared pointers to database
273  const_cast<RefElement_multiIndex *>(refElementPtr)
274  ->insert(shared_ref_fe_vec.begin(), shared_ref_fe_vec.end());
276  }
277 };
278 
280  const BitRefLevel bit,
281  const bool only_tets, int verb,
282  Range *adj_ents_ptr) const {
283  MoFEM::Interface &m_field = cOre;
284  auto ref_ents_ptr = m_field.get_ref_ents();
285  auto ref_fe_ptr = m_field.get_ref_finite_elements();
287 
289  MOFEM_LOG_C("BitRefSelf", Sev::noisy, "Number of entities to add %d",
290  ents.size());
291 
292  CHKERR setElementsBitRefLevel(ents, bit, verb);
293 
294  if (!ents.empty()) {
295  for (int d = 3; d >= 1; --d) {
296  Range dim_ents;
297  if (only_tets && d == 3) {
298  dim_ents = ents.subset_by_type(MBTET);
299  } else {
300  dim_ents = ents.subset_by_dimension(d);
301  }
302 
304  MOFEM_LOG_C("BitRefSelf", Sev::noisy,
305  " Number of dim %d entities to add %d", d, dim_ents.size());
306 
307  if (!dim_ents.empty()) {
308  for (int dd = 0; dd < d; ++dd) {
309  Range adj_ents;
310 
311  if (dd == 0) {
312  rval = m_field.get_moab().get_connectivity(ents, adj_ents, true);
313  // rval = m_field.get_moab().get_adjacencies(
314  // dim_ents, dd, true, adj_ents, moab::Interface::UNION);
315 
316  } else {
317  if (adj_ents_ptr) {
318  if (dd == 1) {
319  adj_ents = adj_ents_ptr->subset_by_dimension(MBEDGE);
320  } else if (dd == 2) {
321  adj_ents = adj_ents_ptr->subset_by_dimension(MBTRI);
322  }
323  } else {
324  rval = m_field.get_moab().get_adjacencies(
325  dim_ents, dd, true, adj_ents, moab::Interface::UNION);
326  }
327  }
328 
329  // rval = m_field.get_moab().get_adjacencies(
330  // dim_ents, dd, true, adj_ents, moab::Interface::UNION);
331 
333  MOFEM_LOG_C("BitRefSelf", Sev::noisy,
334  " Number of dim %d adj entities for dim %d to add %d", d,
335  dd, adj_ents.size());
336 
337  if (rval == MB_MULTIPLE_ENTITIES_FOUND) {
338  auto log_message = [&](const auto sev) {
341  MOFEM_LOG("BitRefSelf", sev)
342  << "When get adjacencies moab return MB_MULTIPLE_ENTITIES_ "
343  "FOUND for dim = "
344  << dd << " and dim of entities " << d;
345  MOFEM_LOG_CHANNEL("BitRefSelf"); // reset channel
346  };
347 
348  if (verb <= QUIET)
349  log_message(Sev::noisy);
350  else
351  log_message(Sev::warning);
352 
353  rval = MB_SUCCESS;
354  }
355  MOAB_THROW(rval);
356  for (Range::pair_iterator pit = adj_ents.pair_begin();
357  pit != adj_ents.pair_end(); ++pit) {
358  Range seed_ents_range;
359  // get first and last element of range
360  EntityHandle f = pit->first;
361  EntityHandle s = pit->second;
362  CHKERR SetBitRefLevelTool(m_field, bit, ref_ents_ptr, ref_fe_ptr)
363  .findEntsToAdd(f, s, seed_ents_range);
364  if (!seed_ents_range.empty())
365  CHKERR SetBitRefLevelTool(m_field, bit, ref_ents_ptr, ref_fe_ptr)
366  .addEntsToDatabase(seed_ents_range);
367  }
368  }
369  }
370  }
371  }
372 
374 }
375 
377  const BitRefLevel bit,
378  int verb) const {
379  MoFEM::Interface &m_field = cOre;
380  auto ref_ents_ptr = m_field.get_ref_ents();
381  auto ref_fe_ptr = m_field.get_ref_finite_elements();
383 
384  for (Range::const_pair_iterator pit = ents.pair_begin();
385  pit != ents.pair_end(); pit++) {
386  // get first and last element of range
387  EntityHandle f = pit->first;
388  EntityHandle s = pit->second;
389  Range seed_ents_range; // entities seeded not in database
390  // find ents to add
391  CHKERR SetBitRefLevelTool(m_field, bit, ref_ents_ptr, ref_fe_ptr)
392  .findEntsToAdd(f, s, seed_ents_range);
393  // add elements
394  if (!seed_ents_range.empty())
395  CHKERR SetBitRefLevelTool(m_field, bit, ref_ents_ptr, ref_fe_ptr)
396  .addEntsToDatabase(seed_ents_range);
397 
398  Range seed_fe_range;
399  CHKERR SetBitRefLevelTool(m_field, bit, ref_ents_ptr, ref_fe_ptr)
400  .findElementsToAdd(f, s, seed_fe_range);
401  if (!seed_fe_range.empty()) {
402  CHKERR SetBitRefLevelTool(m_field, bit, ref_ents_ptr, ref_fe_ptr)
403  .addElementsToDatabase(seed_fe_range);
404  }
405  }
406 
408  MOFEM_LOG("BitRefSelf", Sev::noisy)
409  << "Number of entities in databse " << ref_ents_ptr->size();
410  MOFEM_LOG("BitRefSelf", Sev::noisy)
411  << "Number of finite element entities in databse " << ref_fe_ptr->size();
412 
414 }
415 
417  const BitRefLevel bit,
418  int verb) const {
419  MoFEM::Interface &m_field = cOre;
420  auto ref_ents_ptr = m_field.get_ref_ents();
421  auto ref_fe_ptr = m_field.get_ref_finite_elements();
423 
424  for (Range::const_pair_iterator pit = ents.pair_begin();
425  pit != ents.pair_end(); pit++) {
426  // get first and last element of range
427  EntityHandle f = pit->first;
428  EntityHandle s = pit->second;
429  Range seed_ents_range; // entities seeded not in database
430  // find ents to add
431  CHKERR SetBitRefLevelTool(m_field, bit, ref_ents_ptr, ref_fe_ptr)
432  .findEntsToAdd(f, s, seed_ents_range);
433  // add elements
434  if (!seed_ents_range.empty())
435  CHKERR SetBitRefLevelTool(m_field, bit, ref_ents_ptr, ref_fe_ptr)
436  .addEntsToDatabase(seed_ents_range);
437  }
439 }
440 
442  const std::string field_name, const BitRefLevel bit, int verb) const {
443  MoFEM::Interface &m_field = cOre;
445  EntityHandle field_meshset = m_field.get_field_meshset(field_name);
446  Range field_ents;
447  CHKERR m_field.get_moab().get_entities_by_handle(field_meshset, field_ents,
448  true);
449  CHKERR setEntitiesBitRefLevel(field_ents, bit, verb);
451 }
452 
454  const EntityType type, const BitRefLevel bit, const BitRefLevel mask,
455  int verb) const {
457  Range ents;
459  CHKERR setBitRefLevel(ents, BitRefLevel(), false, verb);
461 }
462 
464  const int dim, const BitRefLevel bit, const BitRefLevel mask,
465  int verb) const {
467  Range ents;
468  CHKERR getEntitiesByDimAndRefLevel(bit, mask, dim, ents);
469  CHKERR setBitRefLevel(ents, BitRefLevel(), false, verb);
471 }
472 
474  const BitRefLevel bit,
475  int verb) const {
476  MoFEM::Interface &m_field = cOre;
477  auto ref_ents_ptr = m_field.get_ref_ents();
478  auto ref_fe_ptr = m_field.get_ref_finite_elements();
480  // Add ref entity
481  std::pair<RefEntity_multiIndex::iterator, bool> p_ent =
482  const_cast<RefEntity_multiIndex *>(ref_ents_ptr)
483  ->insert(boost::shared_ptr<RefEntity>(
484  new RefEntity(m_field.get_basic_entity_data_ptr(), meshset)));
485  *(const_cast<RefEntity *>(p_ent.first->get())->getBitRefLevelPtr()) |= bit;
486  // Add ref element
487  boost::shared_ptr<RefElement> fe_ptr =
488  boost::shared_ptr<RefElement>(new RefElement_MESHSET(*p_ent.first));
489  std::pair<RefElement_multiIndex::iterator, bool> p_fe =
490  const_cast<RefElement_multiIndex *>(ref_fe_ptr)->insert(fe_ptr);
491 
493  MOFEM_LOG("BitRefSelf", Sev::noisy)
494  << "Add meshset as ref_ent " << **p_fe.first;
495 
497 }
498 
500  const int dim,
501  const BitRefLevel bit,
502  int verb) const {
503  MoFEM::Interface &m_field = cOre;
505  Range ents;
506  CHKERR m_field.get_moab().get_entities_by_dimension(meshset, dim, ents,
507  false);
508  CHKERR setBitRefLevel(ents, bit, false, verb);
510 }
511 
513  const EntityType type,
514  const BitRefLevel bit,
515  int verb) const {
516  MoFEM::Interface &m_field = cOre;
518  Range ents;
519  CHKERR m_field.get_moab().get_entities_by_type(meshset, type, ents, false);
520  CHKERR setBitRefLevel(ents, bit, false, verb);
522 }
523 
525  boost::function<void(EntityHandle ent, BitRefLevel &bit)> fun) const {
526  MoFEM::Interface &m_field = cOre;
528  auto get_ents = [&]() {
529  Range ents;
530  CHKERR m_field.get_moab().get_entities_by_handle(
531  m_field.get_moab().get_root_set(), ents, true);
532  ents = subtract(ents, ents.subset_by_type(MBENTITYSET));
533  return ents;
534  };
535  CHKERR lambdaBitRefLevel(get_ents(), fun);
537 }
538 
540  const Range &ents,
541  boost::function<void(EntityHandle ent, BitRefLevel &bit)> fun) const {
542  MoFEM::Interface &m_field = cOre;
544  std::vector<const BitRefLevel *> ents_bits_vec;
545  CHKERR RefEntity::getBitRefLevel(m_field.get_moab(), ents, ents_bits_vec);
546  auto eit = ents.begin();
547  for (auto &it : ents_bits_vec) {
548  fun(*eit, const_cast<BitRefLevel &>(*it));
549  ++eit;
550  }
552 };
553 
555  const BitRefLevel &bit,
556  int verb) const {
557  return lambdaBitRefLevel(
558  ents, [&](EntityHandle ent, BitRefLevel &ent_bit) { ent_bit |= bit; });
559 }
560 
562  const int dim,
563  const BitRefLevel bit,
564  int verb) const {
565  MoFEM::Interface &m_field = cOre;
566  moab::Interface &moab = m_field.get_moab();
567  Range ents, adj;
569  CHKERR moab.get_entities_by_dimension(meshset, dim, ents, true);
570  for (int dd = dim - 1; dd >= 0; dd--)
571  CHKERR moab.get_adjacencies(ents, dd, false, adj, moab::Interface::UNION);
572  ents.merge(adj);
573  if (verb == VERY_NOISY)
574  MOFEM_LOG("BitRefSelf", Sev::noisy) << "Add add bit ref level by dim";
575  CHKERR addBitRefLevel(ents, bit, verb);
577 }
578 
580  const bool b, int verb) const {
581  if (verb == VERY_NOISY)
582  MOFEM_LOG("BitRefSelf", Sev::noisy) << "Set bit to " << ents;
583  return lambdaBitRefLevel(
584  ents, [&](EntityHandle ent, BitRefLevel &ent_bit) { ent_bit[n] = b; });
585 }
586 
588  int verb) const {
589  if (verb == VERY_NOISY)
590  MOFEM_LOG("BitRefSelf", Sev::noisy) << "Set bit to all entities";
591  return lambdaBitRefLevel(
592  [&](EntityHandle ent, BitRefLevel &ent_bit) { ent_bit[n] = b; });
593 }
594 
596  const BitRefLevel mask,
597  int verb) const {
599  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "not implemented");
601 }
602 
604  const BitRefLevel mask, int verb,
605  MoFEMTypes mf) const {
606  MoFEM::Interface &m_field = cOre;
607  auto ref_ents_ptr = m_field.get_ref_ents();
609  RefEntity_change_right_shift right_shift(1, mask);
610  for (int ii = 0; ii < shift; ii++) {
611  // delete bits on the right which are shifted to zero
612  BitRefLevel delete_bits = BitRefLevel().set(0) & mask;
613  if (delete_bits.any()) {
614  CHKERR m_field.delete_ents_by_bit_ref(delete_bits, delete_bits, true,
615  verb, mf);
616  }
617  for (RefEntity_multiIndex::iterator ent_it = ref_ents_ptr->begin();
618  ent_it != ref_ents_ptr->end(); ent_it++) {
619  if (verb >= NOISY) {
621  MOFEM_LOG("BitRefSelf", Sev::noisy)
622  << (*ent_it)->getBitRefLevel() << " : ";
623  }
624  right_shift(const_cast<boost::shared_ptr<RefEntity> &>(*ent_it));
625  if (verb == VERY_NOISY) {
627  MOFEM_LOG("BitRefSelf", Sev::noisy) << (*ent_it)->getBitRefLevel();
628  }
629  }
630  }
632 }
633 
635  const BitRefLevel mask,
636  const char *file_name,
637  const char *file_type,
638  const char *options,
639  const bool check_for_empty) const {
640  MoFEM::Interface &m_field = cOre;
641  moab::Interface &moab(m_field.get_moab());
643  EntityHandle meshset;
644  CHKERR moab.create_meshset(MESHSET_SET, meshset);
645  CHKERR getEntitiesByRefLevel(bit, mask, meshset);
646  int nb_ents;
647  CHKERR moab.get_number_entities_by_handle(meshset, nb_ents, true);
648  if (check_for_empty && !nb_ents) {
649  MOFEM_LOG("SELF", Sev::warning)
650  << "No entities to save < " << file_name << " > in writeBitLevel";
652  }
653 
654  CHKERR moab.write_file(file_name, file_type, options, &meshset, 1);
655  CHKERR moab.delete_entities(&meshset, 1);
657 }
658 
661  const int dim, const char *file_name,
662  const char *file_type, const char *options,
663  const bool check_for_empty) const {
664  MoFEM::Interface &m_field = cOre;
665  moab::Interface &moab(m_field.get_moab());
667  Range ents;
668  CHKERR getEntitiesByDimAndRefLevel(bit, mask, dim, ents);
669  if (check_for_empty && ents.empty()) {
670  MOFEM_LOG("SELF", Sev::warning)
671  << "No entities to save < " << file_name << " > in writeBitLevelByDim";
673  }
674  EntityHandle meshset;
675  CHKERR moab.create_meshset(MESHSET_SET, meshset);
676  CHKERR moab.add_entities(meshset, ents);
677  CHKERR moab.write_file(file_name, file_type, options, &meshset, 1);
678  CHKERR moab.delete_entities(&meshset, 1);
680 }
681 
683  const BitRefLevel bit, const BitRefLevel mask, const EntityType type,
684  const char *file_name, const char *file_type, const char *options,
685  const bool check_for_empty) const {
686  MoFEM::Interface &m_field = cOre;
687  moab::Interface &moab(m_field.get_moab());
689  Range ents;
691  if (check_for_empty && ents.empty()) {
692  MOFEM_LOG("SELF", Sev::warning)
693  << "No entities to save < " << file_name << " > in writeBitLevelByType";
695  }
696  EntityHandle meshset;
697  CHKERR moab.create_meshset(MESHSET_SET, meshset);
698  CHKERR moab.add_entities(meshset, ents);
699  CHKERR moab.write_file(file_name, file_type, options, &meshset, 1);
700  CHKERR moab.delete_entities(&meshset, 1);
702 }
703 
705  const char *file_name, const char *file_type, const char *options,
706  const bool check_for_empty) const {
707  MoFEM::Interface &m_field = cOre;
708  moab::Interface &moab(m_field.get_moab());
710  EntityHandle meshset;
711  Range ents;
713  if (check_for_empty && ents.empty())
715  CHKERR moab.create_meshset(MESHSET_SET, meshset);
716  CHKERR moab.add_entities(meshset, ents);
717  CHKERR moab.write_file(file_name, file_type, options, &meshset, 1);
718  CHKERR moab.delete_entities(&meshset, 1);
720 }
721 
723  const BitRefLevel mask, const EntityType type, const char *file_name,
724  const char *file_type, const char *options) {
726  for (int ll = 0; ll != BITREFLEVEL_SIZE; ++ll) {
727  std::string name = boost::lexical_cast<std::string>(ll) + "_" + file_name;
728  CHKERR writeBitLevelByType(BitRefLevel().set(ll), mask, type, name.c_str(),
729  file_type, options, true);
730  }
732 }
733 
735  const BitRefLevel bit, const BitRefLevel mask, const EntityType type,
736  const EntityHandle meshset, int verb) const {
737  MoFEM::Interface &m_field = cOre;
738  moab::Interface &moab(m_field.get_moab());
740  Range ents;
741  CHKERR getEntitiesByTypeAndRefLevel(bit, mask, type, ents, verb);
742  CHKERR moab.add_entities(meshset, ents);
744 }
745 
747  const BitRefLevel mask,
748  Range &ents,
749  int verb) const {
750  MoFEM::Interface &m_field = cOre;
751  moab::Interface &moab(m_field.get_moab());
753 
754  std::vector<EntityHandle> ents_vec;
755  ents_vec.reserve(ents.size());
756 
757  std::vector<BitRefLevel *> tags_bits_ptr_vec(ents.size());
758 
759  Range swap_ents;
760  auto hint = swap_ents.begin();
761 
762  for (Range::pair_iterator p_eit = ents.pair_begin(); p_eit != ents.pair_end();
763  ++p_eit) {
764 
765  EntityHandle f = p_eit->first;
766  const EntityHandle s = p_eit->second;
767 
768  // get bits on entities
769  rval = moab.tag_get_by_ptr(cOre.get_th_RefBitLevel(), Range(f, s),
770  (const void **)(&*tags_bits_ptr_vec.begin()));
771 
772  if (rval == MB_SUCCESS) {
773 
774  auto bit_it = tags_bits_ptr_vec.begin();
775 
776  auto check = [&bit, &mask](const auto &entity_bit) -> bool {
777  return
778 
779  (entity_bit & bit).any() &&
780 
781  ((entity_bit & mask) == entity_bit);
782  };
783 
784  while (f != s + 1) {
785 
786  while (f != s + 1 && !check(**bit_it)) {
787  ++bit_it;
788  ++f;
789  }
790 
791  if (f != s + 1) {
792 
793  const EntityHandle start = f;
794 
795  while (f != (s + 1) && check(**bit_it)) {
796  ++bit_it;
797  ++f;
798  };
799 
800  hint = swap_ents.insert(hint, start, f - 1);
801  }
802  }
803  }
804  }
805 
806  ents.swap(swap_ents);
807 
809 }
810 
812  const BitRefLevel bit, const BitRefLevel mask, const EntityType type,
813  Range &ents, int verb) const {
814  MoFEM::Interface &m_field = cOre;
815  moab::Interface &moab(m_field.get_moab());
817  CHKERR moab.get_entities_by_type(0, type, ents, false);
818  CHKERR filterEntitiesByRefLevel(bit, mask, ents, verb);
820 }
821 
823  const BitRefLevel bit, const BitRefLevel mask, const int dim,
824  const EntityHandle meshset, int verb) const {
825  MoFEM::Interface &m_field = cOre;
826  moab::Interface &moab(m_field.get_moab());
828  Range ents;
829  CHKERR getEntitiesByDimAndRefLevel(bit, mask, dim, ents, verb);
830  CHKERR moab.add_entities(meshset, ents);
832 }
833 
835  const BitRefLevel bit, const BitRefLevel mask, const int dim, Range &ents,
836  int verb) const {
837  MoFEM::Interface &m_field = cOre;
838  moab::Interface &moab(m_field.get_moab());
840  CHKERR moab.get_entities_by_dimension(0, dim, ents, false);
841  CHKERR filterEntitiesByRefLevel(bit, mask, ents, verb);
843 }
844 
846  const BitRefLevel mask,
847  const EntityHandle meshset,
848  const int verb) const {
849  MoFEM::Interface &m_field = cOre;
850  moab::Interface &moab(m_field.get_moab());
852  Range ents;
853  CHKERR getEntitiesByRefLevel(bit, mask, ents, verb);
854  CHKERR moab.add_entities(meshset, ents);
856 }
857 
859  const BitRefLevel mask,
860  Range &ents,
861  const int verb) const {
862  MoFEM::Interface &m_field = cOre;
863  moab::Interface &moab(m_field.get_moab());
865  Range meshset_ents;
866  CHKERR moab.get_entities_by_type(0, MBENTITYSET, meshset_ents, false);
867  CHKERR moab.get_entities_by_handle(0, ents, false);
868  ents.merge(meshset_ents);
869  CHKERR filterEntitiesByRefLevel(bit, mask, ents, verb);
871 }
872 
874  const BitRefLevel mask,
875  const EntityType type,
876  Range &ents,
877  const int verb) const {
878  MoFEM::Interface &m_field = cOre;
879  auto ref_ents_ptr = m_field.get_ref_ents();
881  auto &ref_ents = ref_ents_ptr->get<Ent_Ent_mi_tag>();
882  auto it = ref_ents.lower_bound(get_id_for_min_type(type));
883  auto hi_it = ref_ents.upper_bound(get_id_for_max_type(type));
884  std::vector<EntityHandle> ents_vec;
885  ents_vec.reserve(std::distance(it, hi_it));
886  for (; it != hi_it; it++) {
887  const BitRefLevel &ent_bit = it->get()->getBitRefLevel();
888  if ((ent_bit & mask) == ent_bit && (ent_bit & bit).any())
889  ents_vec.emplace_back(it->get()->getEnt());
890  }
891  ents.insert_list(ents_vec.begin(), ents_vec.end());
892  if (verb > NOISY)
893  MOFEM_LOG("BitRefSelf", Sev::noisy)
894  << "getEntitiesByParentType: " << ents << endl;
896 }
897 
899  MoFEM::Interface &m_field = cOre;
900  moab::Interface &moab = m_field.get_moab();
902  CHKERR moab.get_entities_by_handle(0, ents, false);
903  ents = subtract(ents, ents.subset_by_type(MBENTITYSET));
906 }
907 
909  MoFEM::Interface &m_field = cOre;
910  auto ref_ents_ptr = m_field.get_ref_ents();
912  auto eit = ents.begin();
913  for (; eit != ents.end();) {
914  auto rit = ref_ents_ptr->get<Ent_mi_tag>().find(*eit);
915  if (rit != ref_ents_ptr->get<Ent_mi_tag>().end()) {
916  eit = ents.erase(eit);
917  } else {
918  eit++;
919  }
920  }
922 }
923 
926  const int to_dimension,
927  Range &adj_entities) const {
928  MoFEM::Interface &m_field = cOre;
929  moab::Interface &moab(m_field.get_moab());
931  BitRefLevel bit_from_entity;
932  CHKERR moab.tag_get_data(cOre.get_th_RefBitLevel(), &from_entity, 1,
933  &bit_from_entity);
934  CHKERR moab.get_adjacencies(&from_entity, 1, to_dimension, false,
935  adj_entities);
936  std::vector<BitRefLevel> bit_levels(adj_entities.size());
937  CHKERR moab.tag_get_data(cOre.get_th_RefBitLevel(), adj_entities,
938  &*bit_levels.begin());
939  auto b_it = bit_levels.begin();
940  auto eit = adj_entities.begin();
941  for (; eit != adj_entities.end(); b_it++) {
942  if (bit_from_entity != *b_it) {
943  eit = adj_entities.erase(eit);
944  } else {
945  eit++;
946  }
947  }
949 }
950 
952  const int to_dimension,
953  Range &adj_entities) const {
954  MoFEM::Interface &m_field = cOre;
955  moab::Interface &moab(m_field.get_moab());
957  BitRefLevel bit_from_entity;
958  CHKERR moab.tag_get_data(cOre.get_th_RefBitLevel(), &from_entity, 1,
959  &bit_from_entity);
960  CHKERR moab.get_adjacencies(&from_entity, 1, to_dimension, false,
961  adj_entities);
962  std::vector<BitRefLevel> bit_levels(adj_entities.size());
963  CHKERR moab.tag_get_data(cOre.get_th_RefBitLevel(), adj_entities,
964  &*bit_levels.begin());
965  std::vector<BitRefLevel>::iterator b_it = bit_levels.begin();
966  Range::iterator eit = adj_entities.begin();
967  for (; eit != adj_entities.end(); b_it++) {
968  if (!(bit_from_entity & (*b_it)).any()) {
969  eit = adj_entities.erase(eit);
970  } else {
971  eit++;
972  }
973  }
975 }
976 
978  const Problem *problem_ptr, const EntityHandle *from_entities,
979  const int num_entities, const int to_dimension, Range &adj_entities,
980  const int operation_type, const int verb) const {
982  BitRefLevel bit = problem_ptr->getBitRefLevel();
983  CHKERR getAdjacencies(bit, from_entities, num_entities, to_dimension,
984  adj_entities, operation_type);
986 }
987 
989  const BitRefLevel bit, const EntityHandle *from_entities,
990  const int num_entities, const int to_dimension, Range &adj_entities,
991  const int operation_type, const int verb) const {
992  MoFEM::Interface &m_field = cOre;
993  moab::Interface &moab(m_field.get_moab());
995 
996  if (verb > VERBOSE) {
998  MOFEM_LOG("BitRef", Sev::noisy) << "from: " << bit;
999  }
1000 
1001  CHKERR moab.get_adjacencies(from_entities, num_entities, to_dimension, false,
1002  adj_entities, operation_type);
1003  std::vector<BitRefLevel> bit_levels(adj_entities.size());
1004  CHKERR moab.tag_get_data(cOre.get_th_RefBitLevel(), adj_entities,
1005  &*bit_levels.begin());
1006  std::vector<BitRefLevel>::iterator b_it = bit_levels.begin();
1007  for (Range::iterator eit = adj_entities.begin(); eit != adj_entities.end();
1008  b_it++) {
1009 
1010  if (verb > VERBOSE) {
1011  RefEntity adj_entity(m_field.get_basic_entity_data_ptr(), *eit);
1013  MOFEM_LOG("BitRef", Sev::noisy)
1014  << "to: " << adj_entity.getBitRefLevel() << " : " << adj_entity;
1015  }
1016 
1017  if (!((*b_it) & bit).any()) {
1018  eit = adj_entities.erase(eit);
1019  } else {
1020  eit++;
1021  }
1022  }
1023  if (b_it != bit_levels.end()) {
1024  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "Data inconsistency");
1025  }
1027 }
1028 
1030  const EntityHandle parent, const BitRefLevel &parent_bit,
1031  const BitRefLevel &parent_mask, const BitRefLevel &child_bit,
1032  const BitRefLevel &child_mask, const EntityHandle child,
1033  EntityType child_type, const bool recursive, int verb) {
1034  MoFEM::Interface &m_field = cOre;
1035  moab::Interface &moab = m_field.get_moab();
1037  Range parent_ents;
1038  CHKERR moab.get_entities_by_handle(parent, parent_ents, recursive);
1039  CHKERR filterEntitiesByRefLevel(parent_bit, parent_mask, parent_ents, verb);
1040  if (verb >= VERY_VERBOSE) {
1042  MOFEM_LOG("BitRefSelf", Sev::noisy) << "Parents: " << parent;
1043  }
1044  Range children_ents;
1045  CHKERR updateRangeByChildren(parent_ents, children_ents);
1046  if (child_type < MBMAXTYPE)
1047  children_ents = children_ents.subset_by_type(child_type);
1048  CHKERR filterEntitiesByRefLevel(child_bit, BitRefLevel().set(), children_ents,
1049  verb);
1050  CHKERR moab.add_entities(child, children_ents);
1052 }
1053 
1055  const EntityHandle parent, const BitRefLevel &child_bit,
1056  const EntityHandle child, EntityType child_type, const bool recursive,
1057  int verb) {
1060  parent, BitRefLevel().set(), BitRefLevel().set(), child_bit, child_bit,
1061  child, child_type, recursive, verb);
1063 }
1064 
1066  const BitRefLevel &child_bit, int verb) {
1067  MoFEM::Interface &m_field = cOre;
1068  moab::Interface &moab = m_field.get_moab();
1069  auto fields_ptr = m_field.get_fields();
1071 
1072  for (auto &fit : (*fields_ptr)) {
1073 
1074  const EntityHandle meshset = fit->getMeshset();
1075  Range parent_ents;
1076  CHKERR moab.get_entities_by_handle(meshset, parent_ents, true);
1077 
1078  if (verb >= VERY_VERBOSE) {
1080  MOFEM_LOG("BitRefSelf", Sev::noisy) << "Parnets:" << std::endl
1081  << parent_ents << std::endl;
1082  }
1083 
1084  Range children_ents;
1085  CHKERR updateRangeByChildren(parent_ents, children_ents);
1086  CHKERR filterEntitiesByRefLevel(child_bit, BitRefLevel().set(),
1087  children_ents, verb);
1088 
1089  CHKERR moab.add_entities(meshset, children_ents);
1090 
1091  if (verb >= VERY_VERBOSE) {
1093  MOFEM_LOG("BitRefSelf", Sev::noisy) << "Children: " << std::endl
1094  << children_ents << std::endl;
1095  }
1096  }
1098 }
1099 
1101  const std::string name, const BitRefLevel &child_bit, int verb) {
1102  MoFEM::Interface &m_field = cOre;
1103  moab::Interface &moab = m_field.get_moab();
1105 
1106  EntityHandle field_meshset = m_field.get_field_structure(name)->getMeshset();
1107 
1108  Range parent_ents;
1109  CHKERR moab.get_entities_by_handle(field_meshset, parent_ents, true);
1110 
1111  if (verb >= VERBOSE) {
1113  MOFEM_LOG("BitRefSelf", Sev::noisy) << "Parnets:" << endl
1114  << parent_ents << std::endl;
1115  }
1116 
1117  Range children_ents;
1118  CHKERR updateRangeByChildren(parent_ents, children_ents);
1119  CHKERR filterEntitiesByRefLevel(child_bit, BitRefLevel().set(), children_ents,
1120  verb);
1121 
1122  CHKERR moab.add_entities(field_meshset, children_ents);
1123 
1124  if (verb >= VERBOSE) {
1126  MOFEM_LOG("BitRefSelf", Sev::noisy) << "Children: " << endl
1127  << children_ents << std::endl;
1128  }
1129 
1131 }
1132 
1134  const std::string name, const BitRefLevel &child_bit,
1135  const EntityType fe_ent_type, int verb) {
1136  MoFEM::Interface &m_field = cOre;
1138  EntityHandle meshset = m_field.get_finite_element_meshset(name);
1139  CHKERR updateMeshsetByEntitiesChildren(meshset, child_bit, meshset,
1140  fe_ent_type, false, verb);
1142 }
1143 
1145  Range &child_ents,
1146  MoFEMTypes bh) {
1147  MoFEM::Interface &m_field = cOre;
1148  auto ref_ents_ptr = m_field.get_ref_ents();
1150  auto &ref_ents =
1151  const_cast<RefEntity_multiIndex *>(ref_ents_ptr)->get<Ent_Ent_mi_tag>();
1152  std::vector<EntityHandle> child_ents_vec;
1153  child_ents_vec.reserve(ref_ents.size());
1154  for (Range::const_pair_iterator pit = parent_ents.const_pair_begin();
1155  pit != parent_ents.const_pair_end(); pit++) {
1156  const auto f = pit->first;
1157  auto it = ref_ents.lower_bound(f);
1158  if (it != ref_ents.end()) {
1159  const auto s = pit->second;
1160  auto hi_it = ref_ents.upper_bound(s);
1161  if (bh == MF_EXIST) {
1162  if (std::distance(it, hi_it) != (s - f) + 1) {
1163  SETERRQ2(
1164  PETSC_COMM_SELF, MOFEM_NOT_FOUND,
1165  "Number of entities and entities parents is different %d != %d ",
1166  std::distance(it, hi_it), (s - f) + 1);
1167  }
1168  }
1169  for (; it != hi_it; it++) {
1170 #ifndef NDEBUG
1171  if (it->get()->getEntType() == MBENTITYSET) {
1172  SETERRQ(PETSC_COMM_SELF, MOFEM_IMPOSSIBLE_CASE,
1173  "This should not happen; Entity should not have part of the "
1174  "meshset. It has no children.");
1175  }
1176 #endif
1177  child_ents_vec.emplace_back((*it)->getEnt());
1178  }
1179  } else if (bh == MF_EXIST) {
1180  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_FOUND, "Entities not found");
1181  }
1182  }
1183  child_ents.insert_list(child_ents_vec.begin(), child_ents_vec.end());
1185 }
1186 
1188  Range &parent_ents,
1189  MoFEMTypes bh) {
1190  MoFEM::Interface &m_field = cOre;
1191  auto ref_ents_ptr = m_field.get_ref_ents();
1193  auto &ref_ents =
1194  const_cast<RefEntity_multiIndex *>(ref_ents_ptr)->get<Ent_mi_tag>();
1195  std::vector<EntityHandle> parent_ents_vec;
1196  parent_ents_vec.reserve(ref_ents.size());
1197  for (Range::const_pair_iterator pit = child_ents.const_pair_begin();
1198  pit != child_ents.const_pair_end(); pit++) {
1199  const auto f = pit->first;
1200  auto it = ref_ents.lower_bound(f);
1201  if (it != ref_ents.end()) {
1202  const auto s = pit->second;
1203  auto hi_it = ref_ents.upper_bound(s);
1204  if (bh == MF_EXIST) {
1205  if (std::distance(it, hi_it) != (s - f) + 1) {
1206  SETERRQ2(
1207  PETSC_COMM_SELF, MOFEM_NOT_FOUND,
1208  "Number of entities and entities parents is different %d != %d ",
1209  std::distance(it, hi_it), (s - f) + 1);
1210  }
1211  }
1212  for (; it != hi_it; it++) {
1213 #ifndef NDEBUG
1214  if (it->get()->getEntType() == MBENTITYSET) {
1215  SETERRQ(PETSC_COMM_SELF, MOFEM_IMPOSSIBLE_CASE,
1216  "This should not happen; Entity should not have part of the "
1217  "meshset. It has no children.");
1218  }
1219 #endif
1220  auto parent = (*it)->getParentEnt();
1221  if (parent)
1222  parent_ents_vec.emplace_back(parent);
1223  }
1224  } else if (bh == MF_EXIST) {
1225  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_FOUND, "Entities not found");
1226  }
1227  }
1228  parent_ents.insert_list(parent_ents_vec.begin(), parent_ents_vec.end());
1230 }
1231 
1234  MOFEM_LOG_CHANNEL("WORLD");
1235 
1236  if (changes)
1237  *changes = false;
1238 
1239  if (Tag th = 0; moab.tag_get_handle("_RefBitLevel", th) == MB_SUCCESS) {
1240 
1241  MOFEM_TAG_AND_LOG("WORLD", Sev::verbose, "BitRefManager") << "Tag found";
1242 
1243  auto get_old_tag = [&](auto &&name) {
1244  Tag th;
1245  CHK_MOAB_THROW(moab.tag_get_handle(name, th),
1246  "bit ref level handle does not exist");
1247  return th;
1248  };
1249 
1250  auto get_new_tag = [&](auto &&name, auto &&def_val) {
1251  Tag th;
1252  CHK_MOAB_THROW(moab.tag_get_handle(
1253  name, sizeof(BitRefLevel), MB_TYPE_OPAQUE, th,
1254  MB_TAG_CREAT | MB_TAG_BYTES | MB_TAG_SPARSE, &def_val),
1255  "can not create tag");
1256  return th;
1257  };
1258 
1259  int length;
1260  CHKERR moab.tag_get_length(get_old_tag("_RefBitLevel"), length);
1261 
1262  if (sizeof(BitRefLevel) != length) {
1263 
1264  if(changes)
1265  *changes = true;
1266 
1267  MOFEM_TAG_AND_LOG("WORLD", Sev::verbose, "BitRefManager")
1268  << "Fixing tag length";
1269 
1270  Range all_ents;
1271  CHKERR moab.get_entities_by_type(0, MBENTITYSET, all_ents, true);
1272  CHKERR moab.get_entities_by_handle(0, all_ents, true);
1273 
1274  auto process_tag = [&](auto &&name, auto &&def_val) {
1276  auto tag_old = get_old_tag(name);
1277  auto get_bits = [&]() {
1278  std::vector<BitRefLevel> bits;
1279  bits.reserve(all_ents.size());
1280  auto it_bit = bits.begin();
1281  for (auto e : all_ents) {
1282  const void *data;
1283  int data_size;
1284  CHKERR moab.tag_get_by_ptr(tag_old, &e, 1, (const void **)&data,
1285  &data_size);
1286  bcopy(
1287  data, &*it_bit,
1288  std::min(sizeof(BitRefLevel), static_cast<size_t>(data_size)));
1289  ++it_bit;
1290  }
1291  return bits;
1292  };
1293  auto bits = get_bits();
1294  CHKERR moab.tag_delete(tag_old);
1295  auto tag_new = get_new_tag(name, def_val);
1296  auto it_bit = bits.begin();
1297  for (auto e : all_ents) {
1298  if (it_bit->any()) {
1299  CHKERR moab.tag_set_data(tag_new, &e, 1, &*it_bit);
1300  }
1301  ++it_bit;
1302  }
1304  };
1305 
1306  CHKERR process_tag("_RefBitLevel", BitRefLevel() = 0);
1307  CHKERR process_tag("_RefBitLevelMask", BitRefLevel().set());
1308  }
1309  }
1310 
1311  MOFEM_LOG_CHANNEL("WORLD");
1313 }
1314 
1315 } // namespace MoFEM
MoFEMFunctionReturnHot
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:447
CHK_MOAB_THROW
#define CHK_MOAB_THROW(err, msg)
Check error code of MoAB function and throw MoFEM exception.
Definition: definitions.h:576
MoFEM::BitRefManager::addBitRefLevelByDim
MoFEMErrorCode addBitRefLevelByDim(const EntityHandle meshset, const int dim, const BitRefLevel bit, int verb=QUIET) const
add bit ref level by dimension
Definition: BitRefManager.cpp:561
MoFEM::Ent_mi_tag
Definition: TagMultiIndices.hpp:21
MoFEM::LogManager::checkIfChannelExist
static bool checkIfChannelExist(const std::string channel)
Check if channel exist.
Definition: LogManager.cpp:404
MoFEM::CoreInterface::get_finite_element_meshset
virtual EntityHandle get_finite_element_meshset(const std::string name) const =0
MoFEM::BitRefManager::getEntitiesByRefLevel
MoFEMErrorCode getEntitiesByRefLevel(const BitRefLevel bit, const BitRefLevel mask, const EntityHandle meshset, const int verb=QUIET) const
add all ents from ref level given by bit to meshset
Definition: BitRefManager.cpp:845
MoFEM::CoreTmp< 0 >
Core (interface) class.
Definition: Core.hpp:82
MoFEM::BitRefManager::shiftLeftBitRef
MoFEMErrorCode shiftLeftBitRef(const int shift, const BitRefLevel mask=BitRefLevel().set(), int verb=DEFAULT_VERBOSITY) const
left shift bit ref level
Definition: BitRefManager.cpp:595
sdf_hertz.d
float d
Definition: sdf_hertz.py:5
MOFEM_LOG_CHANNEL
#define MOFEM_LOG_CHANNEL(channel)
Set and reset channel.
Definition: LogManager.hpp:284
EntityHandle
MoFEM::BitRefManager::setElementsBitRefLevel
MoFEMErrorCode setElementsBitRefLevel(const Range &ents, const BitRefLevel bit=BitRefLevel(), int verb=QUIET) const
add entities to database and set bit ref level
Definition: BitRefManager.cpp:376
NOISY
@ NOISY
Definition: definitions.h:211
MoFEM::BitRefManager::setBitRefLevelByType
MoFEMErrorCode setBitRefLevelByType(const EntityHandle meshset, const EntityType type, const BitRefLevel bit, int verb=QUIET) const
Set the Bit Ref Level By Type object.
Definition: BitRefManager.cpp:512
MoFEM::BitRefManager::getEntitiesByParentType
MoFEMErrorCode getEntitiesByParentType(const BitRefLevel bit, const BitRefLevel mask, const EntityType type, Range &ents, const int verb=QUIET) const
get entities by bit ref level and type of parent
Definition: BitRefManager.cpp:873
MoFEM::BitRefManager::filterEntitiesNotInDatabase
MoFEMErrorCode filterEntitiesNotInDatabase(Range &ents) const
Get entities not in database.
Definition: BitRefManager.cpp:908
MoFEM::Field::getMeshset
EntityHandle getMeshset() const
Get field meshset.
Definition: FieldMultiIndices.hpp:122
MoFEM::RefElement_multiIndex
multi_index_container< boost::shared_ptr< RefElement >, indexed_by< ordered_unique< tag< Ent_mi_tag >, const_mem_fun< RefElement::interface_type_RefEntity, EntityHandle, &RefElement::getEnt > > > > RefElement_multiIndex
Definition: RefElementMultiIndices.hpp:180
MoFEM::CoreInterface::get_field_structure
virtual const Field * get_field_structure(const std::string &name, enum MoFEMTypes bh=MF_EXIST) const =0
get field structure
MoFEM::Exceptions::MoFEMErrorCode
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:56
MoFEM::RefEntityTmp< 0 >
Struct keeps handle to refined handle.
Definition: RefEntsMultiIndices.hpp:141
MoFEM::BitRefManager::setFieldEntitiesBitRefLevel
MoFEMErrorCode setFieldEntitiesBitRefLevel(const std::string field_name, const BitRefLevel bit=BitRefLevel(), int verb=QUIET) const
Set the bit ref level to entities in the field meshset.
Definition: BitRefManager.cpp:441
VERY_NOISY
@ VERY_NOISY
Definition: definitions.h:212
MoFEM::th
Tag th
Definition: Projection10NodeCoordsOnField.cpp:122
MOFEM_LOG_ATTRIBUTES
#define MOFEM_LOG_ATTRIBUTES(channel, bit)
Add attributes to channel.
Definition: LogManager.hpp:296
MoFEM::SetBitRefLevelTool::refEntsPtr
const RefEntity_multiIndex * refEntsPtr
access to ents database
Definition: BitRefManager.cpp:45
MoFEM::BitRefManager::cOre
MoFEM::Core & cOre
Definition: BitRefManager.hpp:26
MOAB_THROW
#define MOAB_THROW(err)
Check error code of MoAB function and throw MoFEM exception.
Definition: definitions.h:541
MoFEM::Ent_Ent_mi_tag
Definition: TagMultiIndices.hpp:55
MoFEM::BitRefManager::lambdaBitRefLevel
MoFEMErrorCode lambdaBitRefLevel(boost::function< void(EntityHandle ent, BitRefLevel &bit)> fun) const
Process bit ref level by lambda function.
Definition: BitRefManager.cpp:524
MoFEM::CoreInterface::get_field_meshset
virtual EntityHandle get_field_meshset(const std::string name) const =0
get field meshset
MoFEM::LogManager::createSink
static boost::shared_ptr< SinkType > createSink(boost::shared_ptr< std::ostream > stream_ptr, std::string comm_filter)
Create a sink object.
Definition: LogManager.cpp:298
MoFEM::RefEntity_multiIndex
multi_index_container< boost::shared_ptr< RefEntity >, indexed_by< ordered_unique< tag< Ent_mi_tag >, const_mem_fun< RefEntity, EntityHandle, &RefEntity::getEnt > >, ordered_non_unique< tag< Ent_Ent_mi_tag >, const_mem_fun< RefEntity, EntityHandle, &RefEntity::getParentEnt > >, ordered_non_unique< tag< Composite_EntType_and_ParentEntType_mi_tag >, composite_key< RefEntity, const_mem_fun< RefEntity, EntityType, &RefEntity::getEntType >, const_mem_fun< RefEntity, EntityType, &RefEntity::getParentEntType > > >, ordered_non_unique< tag< Composite_ParentEnt_And_EntType_mi_tag >, composite_key< RefEntity, const_mem_fun< RefEntity, EntityType, &RefEntity::getEntType >, const_mem_fun< RefEntity, EntityHandle, &RefEntity::getParentEnt > > > > > RefEntity_multiIndex
Definition: RefEntsMultiIndices.hpp:760
MOFEM_IMPOSSIBLE_CASE
@ MOFEM_IMPOSSIBLE_CASE
Definition: definitions.h:35
MoFEM::BitRefManager::writeBitLevel
MoFEMErrorCode writeBitLevel(const BitRefLevel bit, const BitRefLevel mask, const char *file_name, const char *file_type, const char *options, const bool check_for_empty=true) const
Write bit ref level to file.
Definition: BitRefManager.cpp:634
MoFEM::RefEntityTmp< 0 >::getBitRefLevel
const BitRefLevel & getBitRefLevel() const
Get entity ref bit refinement signature.
Definition: RefEntsMultiIndices.hpp:509
MoFEM::SetBitRefLevelTool::addElementsToDatabase
MoFEMErrorCode addElementsToDatabase(Range &seed_fe_range) const
Definition: BitRefManager.cpp:189
MoFEM::BitRefManager::filterEntitiesByRefLevel
MoFEMErrorCode filterEntitiesByRefLevel(const BitRefLevel bit, const BitRefLevel mask, Range &ents, int verb=QUIET) const
filter entities by bit ref level
Definition: BitRefManager.cpp:746
MoFEM::DeprecatedCoreInterface
Deprecated interface functions.
Definition: DeprecatedCoreInterface.hpp:16
MoFEM::BitRefManager::fixTagSize
static MoFEMErrorCode fixTagSize(moab::Interface &moab, bool *changed=nullptr)
Fix tag size when BITREFLEVEL_SIZE of core library is different than file BITREFLEVEL_SIZE.
Definition: BitRefManager.cpp:1232
MoFEM::Interface
DeprecatedCoreInterface Interface
Definition: Interface.hpp:1975
MOFEM_LOG_FUNCTION
#define MOFEM_LOG_FUNCTION()
Set scope.
Definition: LogManager.hpp:325
MoFEM::BitRefManager::addToDatabaseBitRefLevelByType
MoFEMErrorCode addToDatabaseBitRefLevelByType(const EntityType type, const BitRefLevel bit, const BitRefLevel mask=BitRefLevel().set(), int verb=QUIET) const
Add entities which exist in MoAB database, and have set appropiate BitRef level tag,...
Definition: BitRefManager.cpp:453
MoFEM::LogManager::BitScope
@ BitScope
Definition: LogManager.hpp:49
MoFEM::Exceptions::rval
static MoFEMErrorCodeGeneric< moab::ErrorCode > rval
Definition: Exceptions.hpp:74
MoFEM::BitRefManager::setEntitiesBitRefLevel
MoFEMErrorCode setEntitiesBitRefLevel(const Range &ents, const BitRefLevel bit=BitRefLevel(), int verb=QUIET) const
add entities to database and set bit ref level
Definition: BitRefManager.cpp:416
VERBOSE
@ VERBOSE
Definition: definitions.h:209
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:535
MoFEM::BitRefManager::writeEntitiesNotInDatabase
MoFEMErrorCode writeEntitiesNotInDatabase(const char *file_name, const char *file_type, const char *options, const bool check_for_empty=true) const
Write ents not in database.
Definition: BitRefManager.cpp:704
MoFEM::BitRefManager::setBitRefLevel
MoFEMErrorCode setBitRefLevel(const Range &ents, const BitRefLevel bit, const bool only_tets=true, int verb=0, Range *adj_ents_ptr=nullptr) const
add entities to database and set bit ref level
Definition: BitRefManager.cpp:279
MoFEM::CoreInterface::get_moab
virtual moab::Interface & get_moab()=0
MoFEM
implementation of Data Operators for Forces and Sources
Definition: Common.hpp:10
MoFEM::BitRefManager::setBitRefLevelByDim
MoFEMErrorCode setBitRefLevelByDim(const EntityHandle meshset, const int dim, const BitRefLevel bit, int verb=QUIET) const
Set the Bit Ref Level By Dim object.
Definition: BitRefManager.cpp:499
MoFEM::insertOrdered
moab::Range::iterator insertOrdered(Range &r, Extractor, Iterator begin_iter, Iterator end_iter)
Insert ordered mofem multi-index into range.
Definition: Templates.hpp:1788
MoFEM::BitRefManager::getEntitiesByDimAndRefLevel
MoFEMErrorCode getEntitiesByDimAndRefLevel(const BitRefLevel bit, const BitRefLevel mask, const int dim, const EntityHandle meshset, int verb=0) const
add all ents from ref level given by bit to meshset
Definition: BitRefManager.cpp:822
MoFEM::BitRefManager::updateFiniteElementMeshsetByEntitiesChildren
MoFEMErrorCode updateFiniteElementMeshsetByEntitiesChildren(const std::string name, const BitRefLevel &child_bit, const EntityType fe_ent_type, int verb=0)
update finite element meshset by child entities
Definition: BitRefManager.cpp:1133
MOFEM_LOG_C
#define MOFEM_LOG_C(channel, severity, format,...)
Definition: LogManager.hpp:311
MoFEM::RefElement_EDGE
keeps data about abstract EDGE finite element
Definition: RefElementMultiIndices.hpp:111
bit
auto bit
set bit
Definition: hanging_node_approx.cpp:75
convert.type
type
Definition: convert.py:64
MoFEM::CoreInterface::get_basic_entity_data_ptr
virtual boost::shared_ptr< BasicEntityData > & get_basic_entity_data_ptr()=0
Get pointer to basic entity data.
MoFEM::SetBitRefLevelTool::mField
MoFEM::Interface & mField
Definition: BitRefManager.cpp:42
MoFEM::BitRefManager::updateRangeByParent
MoFEMErrorCode updateRangeByParent(const Range &child_ents, Range &parent_ents, MoFEMTypes bh=MF_ZERO)
Update range by parents.
Definition: BitRefManager.cpp:1187
MoFEM::LogManager::getStrmSync
static boost::shared_ptr< std::ostream > getStrmSync()
Get the strm sync object.
Definition: LogManager.cpp:348
MoFEM::CoreTmp< 0 >::get_th_RefBitLevel
Tag get_th_RefBitLevel() const
Definition: Core.hpp:198
MoFEM::LogManager::getStrmWorld
static boost::shared_ptr< std::ostream > getStrmWorld()
Get the strm world object.
Definition: LogManager.cpp:344
MoFEM::BitRefManager::updateMeshsetByEntitiesChildren
MoFEMErrorCode updateMeshsetByEntitiesChildren(const EntityHandle parent, const BitRefLevel &parent_bit, const BitRefLevel &parent_mask, const BitRefLevel &child_bit, const BitRefLevel &child_mask, const EntityHandle child, EntityType child_type, const bool recursive=false, int verb=0)
Get child entities form meshset containing parent entities.
Definition: BitRefManager.cpp:1029
MoFEM::RefElement_MESHSET
keeps data about abstract MESHSET finite element
Definition: RefElementMultiIndices.hpp:56
MoFEM::BitRefManager::addBitRefLevel
MoFEMErrorCode addBitRefLevel(const Range &ents, const BitRefLevel &bit, int verb=QUIET) const
add bit ref level to ref entity
Definition: BitRefManager.cpp:554
MoFEM::get_id_for_max_type
EntityHandle get_id_for_max_type()
Definition: RefEntsMultiIndices.hpp:13
MoFEM::SetBitRefLevelTool
tool class with methods used more than twp times
Definition: BitRefManager.cpp:40
MoFEM::BitRefManager::setBitLevelToMeshset
MoFEMErrorCode setBitLevelToMeshset(const EntityHandle meshset, const BitRefLevel bit, int verb=0) const
Definition: BitRefManager.cpp:473
MoFEM::RefEntExtractor
Extract entity handle form multi-index container.
Definition: Templates.hpp:1762
MoFEM::LogManager::getStrmSelf
static boost::shared_ptr< std::ostream > getStrmSelf()
Get the strm self object.
Definition: LogManager.cpp:340
MoFEM::BitRefManager::writeBitLevelByDim
MoFEMErrorCode writeBitLevelByDim(const BitRefLevel bit, const BitRefLevel mask, const int dim, const char *file_name, const char *file_type, const char *options, const bool check_for_empty=true) const
Write bit ref level to file.
Definition: BitRefManager.cpp:660
MoFEM::CoreInterface::delete_ents_by_bit_ref
virtual MoFEMErrorCode delete_ents_by_bit_ref(const BitRefLevel bit, const BitRefLevel mask, const bool remove_parent=false, int verb=DEFAULT_VERBOSITY, MoFEMTypes mf=MF_ZERO)=0
delete entities form mofem and moab database
MoFEM::BitRefManager::BitRefManager
BitRefManager(const MoFEM::Core &core)
Definition: BitRefManager.cpp:16
MOFEM_LOG_TAG
#define MOFEM_LOG_TAG(channel, tag)
Tag channel.
Definition: LogManager.hpp:339
MoFEM::RefEntity
RefEntityTmp< 0 > RefEntity
Definition: RefEntsMultiIndices.hpp:566
MoFEM::BitRefManager::addToDatabaseBitRefLevelByDim
MoFEMErrorCode addToDatabaseBitRefLevelByDim(const int dim, const BitRefLevel bit, const BitRefLevel mask=BitRefLevel().set(), int verb=QUIET) const
Add entities which exist in MoAB database, and have set appropiate BitRef level tag,...
Definition: BitRefManager.cpp:463
MOFEM_NOT_FOUND
@ MOFEM_NOT_FOUND
Definition: definitions.h:33
field_name
constexpr auto field_name
Definition: poisson_2d_homogeneous.cpp:13
MoFEM::SetBitRefLevelTool::SetBitRefLevelTool
SetBitRefLevelTool(MoFEM::Interface &m_field, const BitRefLevel &bit, const RefEntity_multiIndex *ref_ents_ptr, const RefElement_multiIndex *ref_element_ptr)
constrictor
Definition: BitRefManager.cpp:51
convert.n
n
Definition: convert.py:82
MoFEM::CoreInterface::get_fields
virtual const Field_multiIndex * get_fields() const =0
Get the fields object.
MoFEM::BitRefManager::updateRangeByChildren
MoFEMErrorCode updateRangeByChildren(const Range &parent, Range &child, MoFEMTypes bh=MF_ZERO)
Update range by childrens.
Definition: BitRefManager.cpp:1144
MoFEM::UnknownInterface
base class for all interface classes
Definition: UnknownInterface.hpp:34
MoFEM::RefEntity_change_right_shift
ref mofem entity, right shift
Definition: RefEntsMultiIndices.hpp:840
MoFEM::SetBitRefLevelTool::findEntsToAdd
MoFEMErrorCode findEntsToAdd(EntityHandle f, EntityHandle s, Range &seed_ents_range) const
find entities and change entity bit if in database
Definition: BitRefManager.cpp:59
Range
FTensor::dd
const Tensor2_symmetric_Expr< const ddTensor0< T, Dim, i, j >, typename promote< T, double >::V, Dim, i, j > dd(const Tensor0< T * > &a, const Index< i, Dim > index1, const Index< j, Dim > index2, const Tensor1< int, Dim > &d_ijk, const Tensor1< double, Dim > &d_xyz)
Definition: ddTensor0.hpp:33
MOFEM_TAG_AND_LOG
#define MOFEM_TAG_AND_LOG(channel, severity, tag)
Tag and log in channel.
Definition: LogManager.hpp:362
MOFEM_LOG
#define MOFEM_LOG(channel, severity)
Log.
Definition: LogManager.hpp:308
MoFEM::RefElement_VERTEX
keeps data about abstract VERTEX finite element
Definition: RefElementMultiIndices.hpp:123
MoFEM::RefElement_PRISM
keeps data about abstract PRISM finite element
Definition: RefElementMultiIndices.hpp:66
MoFEM::BitRefManager::getAdjacencies
virtual MoFEMErrorCode getAdjacencies(const Problem *problem_ptr, const EntityHandle *from_entities, const int num_entities, const int to_dimension, Range &adj_entities, const int operation_type=moab::Interface::INTERSECT, const int verb=0) const
Get the adjacencies associated with a entity to entities of a specified dimension.
Definition: BitRefManager.cpp:977
MAX_CORE_TMP
#define MAX_CORE_TMP
maximal number of cores
Definition: definitions.h:217
MoFEM::BitRefManager::getAdjacenciesAny
virtual MoFEMErrorCode getAdjacenciesAny(const EntityHandle from_entity, const int to_dimension, Range &adj_entities) const
Get the adjacencies associated with a entity to entities of a specified dimension.
Definition: BitRefManager.cpp:951
MoFEM::BitRefManager::query_interface
MoFEMErrorCode query_interface(boost::typeindex::type_index type_index, UnknownInterface **iface) const
Definition: BitRefManager.cpp:9
HenckyOps::f
auto f
Definition: HenckyOps.hpp:15
MoFEM::BitRefManager::updateFieldMeshsetByEntitiesChildren
MoFEMErrorCode updateFieldMeshsetByEntitiesChildren(const BitRefLevel &child_bit, int verb=0)
update fields meshesets by child entities
Definition: BitRefManager.cpp:1065
MoFEM::SetBitRefLevelTool::bIt
const BitRefLevel & bIt
bit to set
Definition: BitRefManager.cpp:44
BITREFLEVEL_SIZE
#define BITREFLEVEL_SIZE
max number of refinements
Definition: definitions.h:219
MoFEM::BitRefManager::writeBitLevelByType
MoFEMErrorCode writeBitLevelByType(const BitRefLevel bit, const BitRefLevel mask, const EntityType type, const char *file_name, const char *file_type, const char *options, const bool check_for_empty=true) const
Write bit ref level to file.
Definition: BitRefManager.cpp:682
MoFEM::CoreInterface::get_ref_ents
virtual const RefEntity_multiIndex * get_ref_ents() const =0
Get the ref ents object.
MoFEM::BitRefManager::setNthBitRefLevel
MoFEMErrorCode setNthBitRefLevel(const Range &ents, const int n, const bool b, int verb=QUIET) const
Set nth bit ref level.
Definition: BitRefManager.cpp:579
MOFEM_DATA_INCONSISTENCY
@ MOFEM_DATA_INCONSISTENCY
Definition: definitions.h:31
MoFEM::CoreInterface::get_ref_finite_elements
virtual const RefElement_multiIndex * get_ref_finite_elements() const =0
Get the ref finite elements object.
MoFEM::Problem::getBitRefLevel
BitRefLevel getBitRefLevel() const
Definition: ProblemsMultiIndices.hpp:383
MoFEM::BitRefManager::getAdjacenciesEquality
virtual MoFEMErrorCode getAdjacenciesEquality(const EntityHandle from_entity, const int to_dimension, Range &adj_entities) const
Get the adjacencies associated with a entity to entities of a specified dimension.
Definition: BitRefManager.cpp:925
MoFEM::RefElementVolume
keeps data about abstract TET finite element
Definition: RefElementMultiIndices.hpp:81
MoFEMTypes
MoFEMTypes
Those types control how functions respond on arguments, f.e. error handling.
Definition: definitions.h:97
MoFEM::RefEntityTmp< 0 >::getBitRefLevel
static MoFEMErrorCode getBitRefLevel(Interface &moab, Range ents, std::vector< BitRefLevel > &vec_bit_ref_level)
Definition: RefEntsMultiIndices.hpp:433
MoFEM::BitRefManager
Managing BitRefLevels.
Definition: BitRefManager.hpp:21
MoFEM::Types::BitRefLevel
std::bitset< BITREFLEVEL_SIZE > BitRefLevel
Bit structure attached to each entity identifying to what mesh entity is attached.
Definition: Types.hpp:40
MoFEMFunctionBeginHot
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:440
QUIET
@ QUIET
Definition: definitions.h:208
MoFEM::Problem
keeps basic data about problem
Definition: ProblemsMultiIndices.hpp:54
MoFEM::SetBitRefLevelTool::baseEntData
boost::shared_ptr< BasicEntityData > & baseEntData
base entity data
Definition: BitRefManager.cpp:48
MoFEM::SetBitRefLevelTool::refElementPtr
const RefElement_multiIndex * refElementPtr
access to fe database
Definition: BitRefManager.cpp:46
MoFEM::LogManager::setLog
static LoggerType & setLog(const std::string channel)
Set ans resset chanel logger.
Definition: LogManager.cpp:389
MoFEM::BitRefManager::writeEntitiesAllBitLevelsByType
MoFEMErrorCode writeEntitiesAllBitLevelsByType(const BitRefLevel mask, const EntityType type, const char *file_name, const char *file_type, const char *options)
Write all entities by bit levels and type.
Definition: BitRefManager.cpp:722
MoFEM::BitRefManager::getEntitiesByTypeAndRefLevel
MoFEMErrorCode getEntitiesByTypeAndRefLevel(const BitRefLevel bit, const BitRefLevel mask, const EntityType type, const EntityHandle meshset, int verb=0) const
add all ents from ref level given by bit to meshset
Definition: BitRefManager.cpp:734
fun
auto fun
Function to approximate.
Definition: dg_projection.cpp:36
MoFEM::get_id_for_min_type
EntityHandle get_id_for_min_type()
Definition: RefEntsMultiIndices.hpp:18
MF_EXIST
@ MF_EXIST
Definition: definitions.h:100
MoFEM::BitRefManager::getAllEntitiesNotInDatabase
MoFEMErrorCode getAllEntitiesNotInDatabase(Range &ents) const
Get all entities not in database.
Definition: BitRefManager.cpp:898
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:416
MoFEM::SetBitRefLevelTool::addEntsToDatabase
MoFEMErrorCode addEntsToDatabase(const Range &seed_ents_range) const
Definition: BitRefManager.cpp:151
MOFEM_NOT_IMPLEMENTED
@ MOFEM_NOT_IMPLEMENTED
Definition: definitions.h:32
MoFEM::RefEntityTmp
Definition: RefEntsMultiIndices.hpp:118
VERY_VERBOSE
@ VERY_VERBOSE
Definition: definitions.h:210
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:346
MoFEM::CoreInterface::getValue
virtual const int getValue() const =0
Get the core.
MoFEM::BitRefManager::shiftRightBitRef
MoFEMErrorCode shiftRightBitRef(const int shift, const BitRefLevel mask=BitRefLevel().set(), int verb=DEFAULT_VERBOSITY, MoFEMTypes mf=MF_ZERO) const
right shift bit ref level
Definition: BitRefManager.cpp:603
MoFEM::RefElementFace
keeps data about abstract TRI finite element
Definition: RefElementMultiIndices.hpp:99
MoFEM::SetBitRefLevelTool::addEntsToDatabaseImpl
MoFEMErrorCode addEntsToDatabaseImpl(const Range &seed_ents_range) const
add entities to database
Definition: BitRefManager.cpp:97
MoFEM::SetBitRefLevelTool::findElementsToAdd
MoFEMErrorCode findElementsToAdd(EntityHandle f, EntityHandle s, Range &seed_fe_range) const
Definition: BitRefManager.cpp:169