v0.14.0
BcManager.cpp
Go to the documentation of this file.
1 /** \file BcManager.cpp
2  * \brief Manages boundary conditions
3  * \ingroup bc_manager
4  */
5 
6 namespace MoFEM {
7 
8 namespace BcManagerImplTools {
9 
10 auto get_dim(const Range &ents) {
11  for (auto d : {3, 2, 1})
12  if (ents.num_of_dimension(d))
13  return d;
14  return 0;
15 };
16 
17 auto get_adj_ents(moab::Interface &moab, const Range &ents) {
18  Range verts;
19  CHK_MOAB_THROW(moab.get_connectivity(ents, verts, true), "get verts");
20  const auto dim = get_dim(ents);
21  for (size_t d = 1; d < dim; ++d) {
22  for (auto dd = d + 1; dd <= dim; ++dd) {
23  CHK_MOAB_THROW(moab.get_adjacencies(ents.subset_by_dimension(dd), d,
24  false, verts, moab::Interface::UNION),
25  "get adj");
26  }
27  }
28  verts.merge(ents);
29  return verts;
30 }
31 } // namespace BcManagerImplTools
32 
34 BcManager::query_interface(boost::typeindex::type_index type_index,
35  UnknownInterface **iface) const {
37  *iface = const_cast<BcManager *>(this);
39 }
40 
41 BcManager::BcManager(const Core &core) : cOre(const_cast<Core &>(core)) {
42 
43  if (!LogManager::checkIfChannelExist("BcMngWorld")) {
44  auto core_log = logging::core::get();
45 
46  core_log->add_sink(
48  core_log->add_sink(
50  core_log->add_sink(
52 
53  LogManager::setLog("BcMngWorld");
54  LogManager::setLog("BcMngSync");
55  LogManager::setLog("BcMngSelf");
56 
57  MOFEM_LOG_TAG("BcMngWorld", "BcMng");
58  MOFEM_LOG_TAG("BcMngSync", "BcMng");
59  MOFEM_LOG_TAG("BcMngSelf", "BcMng");
60  }
61 
62  MOFEM_LOG("BcMngWorld", Sev::noisy) << "BC manager created";
63 }
64 
67  ierr = PetscOptionsBegin(PETSC_COMM_WORLD, "", "BcManager options", "none");
68  ierr = PetscOptionsEnd();
69  CHKERRG(ierr);
71 }
72 
74  const std::string problem_name, const std::string block_name,
75  const std::string field_name, int lo, int hi, bool get_low_dim_ents,
76  bool is_distributed_mesh) {
77  Interface &m_field = cOre;
78  auto prb_mng = m_field.getInterface<ProblemsManager>();
80 
81  CHKERR pushMarkDOFsOnEntities(problem_name, block_name, field_name, lo, hi,
82  get_low_dim_ents);
83 
84  auto remove_dofs_on_ents = [&](const Range &ents, const int lo,
85  const int hi) {
86  if (is_distributed_mesh)
87  return prb_mng->removeDofsOnEntities(problem_name, field_name, ents, lo,
88  hi);
89  else
90  return prb_mng->removeDofsOnEntitiesNotDistributed(
91  problem_name, field_name, ents, lo, hi);
92  };
93 
94  for (auto m :
95  m_field.getInterface<MeshsetsManager>()->getCubitMeshsetPtr(std::regex(
96 
97  (boost::format("%s(.*)") % block_name).str()
98 
99  ))
100 
101  ) {
102  const std::string bc_id =
103  problem_name + "_" + field_name + "_" + m->getName();
104  CHKERR remove_dofs_on_ents(bcMapByBlockName.at(bc_id)->bcEnts, lo, hi);
105  bcMapByBlockName.at(bc_id)->bcMarkers = std::vector<unsigned char>();
106  }
107 
109 }
110 
111 MoFEMErrorCode BcManager::pushMarkDOFsOnEntities(const std::string problem_name,
112  const std::string block_name,
113  const std::string field_name,
114  int lo, int hi,
115  bool get_low_dim_ents) {
116  Interface &m_field = cOre;
117  auto prb_mng = m_field.getInterface<ProblemsManager>();
119 
120  // if(problem_name.size())
121  // MOFEM_LOG("BcMngWorld", Sev::warning)
122  // << "Argument problem_name has no effect";
123 
124  auto fix_disp = [&]() {
126 
127  auto mark_fix_dofs = [&](std::vector<unsigned char> &marked_field_dofs,
128  const auto lo, const auto hi) {
129  return prb_mng->modifyMarkDofs(problem_name, ROW, field_name, lo, hi,
130  ProblemsManager::MarkOP::OR, 1,
131  marked_field_dofs);
132  };
133 
134  auto iterate_meshsets = [&](auto &&meshset_vec_ptr) {
136  for (auto m : meshset_vec_ptr) {
137  auto bc = boost::make_shared<BCs>();
138  CHKERR m_field.get_moab().get_entities_by_handle(m->getMeshset(),
139  bc->bcEnts, true);
140  CHKERR m->getAttributes(bc->bcAttributes);
141 
142  if (problem_name.size())
143  CHKERR mark_fix_dofs(bc->bcMarkers, lo, hi);
144  MOFEM_LOG("BcMngWorld", Sev::verbose)
145  << "Found block " << m->getName() << " number of attributes "
146  << bc->bcAttributes.size();
147 
148  if (get_low_dim_ents) {
149  auto low_dim_ents =
150  BcManagerImplTools::get_adj_ents(m_field.get_moab(), bc->bcEnts);
151  bc->bcEnts.swap(low_dim_ents);
152  }
153 
154  CHKERR m_field.getInterface<CommInterface>()->synchroniseEntities(
155  bc->bcEnts);
156  if (problem_name.size())
157  CHKERR prb_mng->markDofs(problem_name, ROW, ProblemsManager::AND,
158  bc->bcEnts, bc->bcMarkers);
159 
160  const std::string bc_id =
161  problem_name + "_" + field_name + "_" + m->getName();
162  bcMapByBlockName[bc_id] = bc;
163  }
165  };
166 
167  CHKERR iterate_meshsets(
168 
169  m_field.getInterface<MeshsetsManager>()->getCubitMeshsetPtr(std::regex(
170 
171  (boost::format("%s(.*)") % block_name).str()
172 
173  ))
174 
175  );
176 
178  };
179 
180  CHKERR fix_disp();
181 
183 }
184 
185 MoFEMErrorCode BcManager::addBlockDOFsToMPCs(const std::string problem_name,
186  const std::string field_name,
187  bool get_low_dim_ents,
188  bool block_name_field_prefix,
189  bool is_distributed_mesh) {
190  Interface &m_field = cOre;
192 
193  if (block_name_field_prefix)
194  MOFEM_LOG("BcMngWorld", Sev::warning)
195  << "Argument block_name_field_prefix=true has no effect";
196  if (is_distributed_mesh)
197  MOFEM_LOG("BcMngWorld", Sev::warning)
198  << "Argument is_distributed_mesh=true has no effect";
199  if (get_low_dim_ents)
200  MOFEM_LOG("BcMngWorld", Sev::warning)
201  << "Argument get_low_dim_ents=true has no effect";
202 
203  auto get_dim = [&](const Range &ents) {
204  for (auto d : {3, 2, 1})
205  if (ents.num_of_dimension(d))
206  return d;
207  return 0;
208  };
209 
210  auto iterate_mpc_meshsets = [&]() {
212 
213  auto mpc_meshset_ptr =
214  m_field.getInterface<MeshsetsManager>()->getCubitMeshsetPtr(
215  std::regex((boost::format("%s(.*)") % "MPC_(.*)").str()));
216 
217  for (auto m : mpc_meshset_ptr) {
218 
219  if (std::regex_match(m->getName(),
220  std::regex("(.*)COUPLING_LINKS(.*)"))) {
221 
222  auto bc = boost::make_shared<BCs>();
223  bc->mpcPtr = boost::make_shared<MPCsType>();
224  bc->mpcPtr->mpcType = MPC::COUPLING;
225 
226  std::string const corresponding_master_ms =
227  std::regex_replace(m->getName(), std::regex("LINKS"), "MASTER");
228 
229  Range links_ents;
230  CHKERR m_field.get_moab().get_entities_by_handle(m->getMeshset(),
231  links_ents, true);
232 
233  Range master_nodes;
234  if (m_field.getInterface<MeshsetsManager>()->checkMeshset(
235  corresponding_master_ms)) {
236  const CubitMeshSets *l;
237  CHKERR m_field.getInterface<MeshsetsManager>()->getCubitMeshsetPtr(
238  corresponding_master_ms, &l);
239  bc->mpcPtr->isReprocitical = false;
240 
241  CHKERR m_field.get_moab().get_entities_by_handle(l->getMeshset(),
242  master_nodes, true);
243  // if (master_nodes.subset_by_dimension(0).size() <
244  // links_ents.subset_by_dimension(1))
245  {
246  auto low_dim_ents = BcManagerImplTools::get_adj_ents(
247  m_field.get_moab(), master_nodes);
248  low_dim_ents = low_dim_ents.subset_by_dimension(0);
249  master_nodes.swap(low_dim_ents);
250  }
251 
252  MOFEM_LOG("BcMngWorld", Sev::verbose)
253  << "Found block MASTER LINKS block: " << l->getName()
254  << " Entities size: " << master_nodes.size();
255 
256  } else {
257  MOFEM_LOG("BcMngWorld", Sev::warning)
258  << "MASTER LINKS block not found: " << corresponding_master_ms
259  << " setting reprocitical constraint. ("
260  << bc->mpcPtr->isReprocitical << ").";
261  }
262 
263  // if (std::regex_match(bc_id, std::regex("(.*)TIE(.*)")))
264  // bc->mpcPtr->mpcType = MPC::TIE;
265  // if (std::regex_match(bc_id, std::regex("(.*)RIGID_BODY(.*)")))
266  // bc->mpcPtr->mpcType = MPC::RIGID_BODY;
267 
268  // std::cout << "links_ents range is: " << links_ents << "\n";
269  for (auto &link : links_ents.subset_by_dimension(1)) {
270  Range verts;
271  CHKERR m_field.get_moab().get_connectivity(&link, 1, verts, true);
272  // std::cout << "verts range is: " << verts << "\n";
273  if (bc->mpcPtr->isReprocitical) {
274  bc->mpcPtr->mpcMasterEnts.insert(verts[0]);
275  bc->mpcPtr->mpcSlaveEnts.insert(verts[1]);
276  } else {
277  for (auto &m_node : verts)
278  if (master_nodes.find(m_node) != master_nodes.end()) {
279  // std::cout << "found ent: " << m_node << "\n";
280  bc->mpcPtr->mpcMasterEnts.insert(m_node);
281  bc->mpcPtr->mpcSlaveEnts.merge(
282  subtract(verts, Range(m_node, m_node)));
283  break;
284  }
285  }
286  }
287 
288  MOFEM_LOG("BcMngWorld", Sev::verbose)
289  << "Found block MPC LINKS block: " << m->getName()
290  << " Entities size (edges): " << links_ents.size()
291  << " Entities size (nodes): "
292  << bc->mpcPtr->mpcMasterEnts.size() +
293  bc->mpcPtr->mpcSlaveEnts.size()
294  << " (" << bc->mpcPtr->mpcMasterEnts.size() << " "
295  << bc->mpcPtr->mpcSlaveEnts.size() << ")";
296 
297  MOFEM_LOG("BcMngSync", Sev::noisy) << *bc->mpcPtr;
298  MOFEM_LOG_SEVERITY_SYNC(m_field.get_comm(), Sev::noisy);
299 
300  // CHKERR m_field.getInterface<CommInterface>()->synchroniseEntities(
301  // bc->mpcPtr->mpcMasterEnts);
302  // CHKERR m_field.getInterface<CommInterface>()->synchroniseEntities(
303  // bc->mpcPtr->mpcSlaveEnts);
304  vector<double> mAttributes;
305  CHKERR m->getAttributes(mAttributes);
306 
307  auto setFlags = [&](const auto &flags) {
308  auto &d = bc->mpcPtr->data;
309  if (flags.empty()) {
310  d.flag1 = d.flag2 = d.flag3 = d.flag4 = d.flag5 = d.flag6 = true;
311  return;
312  }
313  for (size_t i = 0; i < std::min(flags.size(), size_t(6)); ++i)
314  (&d.flag1)[i] = flags[i] > 0.0;
315  };
316 
317  setFlags(mAttributes);
318 
319  MOFEM_LOG("BcMngWorld", Sev::verbose)
320  << "Found block " << m->getName() << " number of entities "
321  << bc->mpcPtr->mpcMasterEnts.size() << " number of attributes "
322  << mAttributes.size() << " highest dim of entities "
323  << get_dim(bc->mpcPtr->mpcMasterEnts);
324 
325  // NOTE: we are not using markers at the moment
326  // for (int i = 0; i != mAttributes.size(); ++i) {
327  // if (mAttributes[i] > 0.0)
328  // CHKERR mark_fix_dofs(bc->mpcPtr->bcMasterMarkers, i, i);
329  // }
330 
331  // if (mAttributes.empty())
332  // CHKERR mark_fix_dofs(bc->mpcPtr->bcMasterMarkers, 0,
333  // MAX_DOFS_ON_ENTITY);
334 
335  const std::string bc_id =
336  problem_name + "_" + field_name + "_" + m->getName();
337  bcMapByBlockName[bc_id] = bc;
338  }
339  }
341  };
342 
343  CHKERR iterate_mpc_meshsets();
344 
346 }
347 
348 boost::shared_ptr<BcManager::BCs>
349 BcManager::popMarkDOFsOnEntities(const std::string block_name) {
350  auto bc_it = bcMapByBlockName.find(block_name);
351  if (bc_it != bcMapByBlockName.end()) {
352  auto bc = bc_it->second;
353  bcMapByBlockName.erase(bc_it);
354  return bc;
355  }
356  return boost::shared_ptr<BCs>();
357 }
358 
359 Range BcManager::getMergedBlocksRange(std::vector<std::regex> bc_regex_vec) {
360  Range ents;
361  if (bcMapByBlockName.size()) {
362  for (auto b : bcMapByBlockName) {
363  for (auto &reg_name : bc_regex_vec) {
364  if (std::regex_match(b.first, reg_name)) {
365  ents.merge(b.second->bcEnts);
366  }
367  }
368  }
369  }
370  return ents;
371 }
372 
374 BcManager::getMergedBlocksMarker(std::vector<std::regex> bc_regex_vec) {
375  BcManager::BcMarkerPtr boundary_marker_ptr;
376  if (bcMapByBlockName.size()) {
377  for (auto b : bcMapByBlockName) {
378  for (auto &reg_name : bc_regex_vec) {
379  if (std::regex_match(b.first, reg_name)) {
380  if (!boundary_marker_ptr)
381  boundary_marker_ptr =
382  boost::make_shared<std::vector<char unsigned>>();
383  boundary_marker_ptr->resize(b.second->bcMarkers.size(), 0);
384  for (int i = 0; i != b.second->bcMarkers.size(); ++i) {
385  (*boundary_marker_ptr)[i] |= b.second->bcMarkers[i];
386  }
387  }
388  }
389  }
390  }
391  return boundary_marker_ptr;
392 }
393 
395  const std::vector<BcManager::BcMarkerPtr> &boundary_markers_ptr_vec) {
396  auto boundary_marker_ptr = boost::make_shared<std::vector<char unsigned>>();
397  for (auto &bcm : boundary_markers_ptr_vec) {
398  boundary_marker_ptr->resize(bcm->size(), 0);
399  for (int i = 0; i != bcm->size(); ++i)
400  (*boundary_marker_ptr)[i] |= (*bcm)[i];
401  }
402  return boundary_marker_ptr;
403 }
404 
405 SmartPetscObj<IS> BcManager::getBlockIS(const std::string block_prefix,
406  const std::string block_name,
407  const std::string field_name,
408  const std::string problem_name, int lo,
409  int hi, SmartPetscObj<IS> is_expand) {
410  Interface &m_field = cOre;
411 
412  const std::string bc_id =
413  block_prefix + "_" + field_name + "_" + block_name + "(.*)";
414 
415  Range bc_ents;
416  for (auto bc : getBcMapByBlockName()) {
417  if (std::regex_match(bc.first, std::regex(bc_id))) {
418  bc_ents.merge(*(bc.second->getBcEntsPtr()));
419  MOFEM_LOG("BcMngWorld", Sev::verbose)
420  << "Get entities from block and add to IS. Block name " << bc.first;
421  }
422  }
423 
424  SmartPetscObj<IS> is_bc;
425  auto get_is = [&]() {
427  CHKERR m_field.getInterface<CommInterface>()->synchroniseEntities(bc_ents);
428  CHKERR m_field.getInterface<ISManager>()->isCreateProblemFieldAndRank(
429  problem_name, ROW, field_name, lo, hi, is_bc, &bc_ents);
430  if (is_expand) {
431  IS is_tmp;
432  CHKERR ISExpand(is_bc, is_expand, &is_tmp);
433  is_bc = SmartPetscObj<IS>(is_tmp);
434  }
435  CHKERR ISSort(is_bc);
437  };
438 
439  if (get_is())
440  CHK_THROW_MESSAGE(MOFEM_DATA_INCONSISTENCY, "IS is not created");
441 
442  return is_bc;
443 }
444 
445 SmartPetscObj<IS> BcManager::getBlockIS(const std::string problem_name,
446  const std::string block_name,
447  const std::string field_name, int lo,
448  int hi, SmartPetscObj<IS> is_expand) {
449  return getBlockIS(problem_name, block_name, field_name, problem_name, lo, hi,
450  is_expand);
451 }
452 
453 template <>
455 BcManager::removeBlockDOFsOnEntities<BcMeshsetType<DISPLACEMENTSET>>(
456  const std::string problem_name, const std::string field_name,
457  bool get_low_dim_ents, bool block_name_field_prefix,
458  bool is_distributed_mesh) {
459  Interface &m_field = cOre;
460  auto prb_mng = m_field.getInterface<ProblemsManager>();
462 
463  CHKERR pushMarkDOFsOnEntities<BcMeshsetType<DISPLACEMENTSET>>(
464  problem_name, field_name, get_low_dim_ents, block_name_field_prefix);
465 
466  std::array<Range, 3> ents_to_remove;
467 
468  for (auto m :
469 
472 
473  const std::string bc_id =
474  problem_name + "_" + field_name + "_DISPLACEMENTSET" +
475  boost::lexical_cast<std::string>(m->getMeshsetId());
476 
477  auto bc = bcMapByBlockName.at(bc_id);
478 
479  if (bc->dispBcPtr) {
480  if (bc->dispBcPtr->data.flag1) {
481  ents_to_remove[0].merge(bc->bcEnts);
482  }
483  if (bc->dispBcPtr->data.flag2) {
484  ents_to_remove[1].merge(bc->bcEnts);
485  }
486  if (bc->dispBcPtr->data.flag3) {
487  ents_to_remove[2].merge(bc->bcEnts);
488  }
489  if (bc->dispBcPtr->data.flag4) {
490  ents_to_remove[1].merge(bc->bcEnts);
491  ents_to_remove[2].merge(bc->bcEnts);
492  }
493  if (bc->dispBcPtr->data.flag5) {
494  ents_to_remove[0].merge(bc->bcEnts);
495  ents_to_remove[2].merge(bc->bcEnts);
496  }
497  if (bc->dispBcPtr->data.flag6) {
498  ents_to_remove[0].merge(bc->bcEnts);
499  ents_to_remove[1].merge(bc->bcEnts);
500  }
501  }
502  bc->bcMarkers = std::vector<unsigned char>();
503  }
504 
505  auto remove_dofs_on_ents = [&](const Range &ents, const int lo,
506  const int hi) {
507  if (is_distributed_mesh)
508  return prb_mng->removeDofsOnEntities(problem_name, field_name, ents, lo,
509  hi);
510  else
511  return prb_mng->removeDofsOnEntitiesNotDistributed(
512  problem_name, field_name, ents, lo, hi);
513  };
514 
515  CHKERR remove_dofs_on_ents(ents_to_remove[0], 0, 0);
516  CHKERR remove_dofs_on_ents(ents_to_remove[1], 1, 1);
517  CHKERR remove_dofs_on_ents(ents_to_remove[2], 2, 2);
518 
520 }
521 
522 template <>
524 BcManager::removeBlockDOFsOnEntities<BcMeshsetType<TEMPERATURESET>>(
525  const std::string problem_name, const std::string field_name,
526  bool get_low_dim_ents, bool block_name_field_prefix,
527  bool is_distributed_mesh) {
528  Interface &m_field = cOre;
529  auto prb_mng = m_field.getInterface<ProblemsManager>();
531 
532  CHKERR pushMarkDOFsOnEntities<BcMeshsetType<TEMPERATURESET>>(
533  problem_name, field_name, get_low_dim_ents, block_name_field_prefix);
534 
535  Range ents_to_remove;
536 
537  for (auto m :
538 
541 
542  ) {
543  const std::string bc_id =
544  problem_name + "_" + field_name + "_TEMPERATURESET" +
545  boost::lexical_cast<std::string>(m->getMeshsetId());
546  auto bc = bcMapByBlockName.at(bc_id);
547  ents_to_remove.merge(bc->bcEnts);
548  bc->bcMarkers = std::vector<unsigned char>();
549  }
550 
551  auto remove_dofs_on_ents = [&](const Range &ents, const int lo,
552  const int hi) {
553  if (is_distributed_mesh)
554  return prb_mng->removeDofsOnEntities(problem_name, field_name, ents, lo,
555  hi);
556  else
557  return prb_mng->removeDofsOnEntitiesNotDistributed(
558  problem_name, field_name, ents, lo, hi);
559  };
560 
561  CHKERR remove_dofs_on_ents(ents_to_remove, 0, MAX_DOFS_ON_ENTITY);
562 
564 }
565 
566 template <>
567 MoFEMErrorCode BcManager::removeBlockDOFsOnEntities<BcMeshsetType<HEATFLUXSET>>(
568  const std::string problem_name, const std::string field_name,
569  bool get_low_dim_ents, bool block_name_field_prefix,
570  bool is_distributed_mesh) {
571  Interface &m_field = cOre;
572  auto prb_mng = m_field.getInterface<ProblemsManager>();
574 
575  CHKERR pushMarkDOFsOnEntities<BcMeshsetType<HEATFLUXSET>>(
576  problem_name, field_name, get_low_dim_ents, block_name_field_prefix);
577 
578  Range ents_to_remove;
579 
580  for (auto m :
581 
583  HEATFLUXSET)
584 
585  ) {
586  const std::string bc_id =
587  problem_name + "_" + field_name + "_HEATFLUXSET" +
588  boost::lexical_cast<std::string>(m->getMeshsetId());
589  auto bc = bcMapByBlockName.at(bc_id);
590  ents_to_remove.merge(bc->bcEnts);
591  bc->bcMarkers = std::vector<unsigned char>();
592  }
593 
594  auto remove_dofs_on_ents = [&](const Range &ents, const int lo,
595  const int hi) {
596  if (is_distributed_mesh)
597  return prb_mng->removeDofsOnEntities(problem_name, field_name, ents, lo,
598  hi);
599  else
600  return prb_mng->removeDofsOnEntitiesNotDistributed(
601  problem_name, field_name, ents, lo, hi);
602  };
603 
604  CHKERR remove_dofs_on_ents(ents_to_remove, 0, MAX_DOFS_ON_ENTITY);
605 
607 }
608 
609 template <>
611 BcManager::removeBlockDOFsOnEntities<BcDisplacementMeshsetType<BLOCKSET>>(
612  const std::string problem_name, const std::string field_name,
613  bool get_low_dim_ents, bool block_name_field_prefix,
614  bool is_distributed_mesh) {
615  Interface &m_field = cOre;
616  auto prb_mng = m_field.getInterface<ProblemsManager>();
618 
619  CHKERR pushMarkDOFsOnEntities<BcDisplacementMeshsetType<BLOCKSET>>(
620  problem_name, field_name, get_low_dim_ents, block_name_field_prefix);
621 
622  std::array<Range, 3> ents_to_remove;
623 
624  for (auto m :
625 
627  BLOCKSET | UNKNOWNNAME)) {
628 
629  const auto block_name = m->getName();
630 
631  std::string bc_id = problem_name + "_" + field_name + "_" + block_name;
632  std::string regex_str;
633  if (block_name_field_prefix) {
634  regex_str = (boost::format("%s_%s_%s_((FIX_(ALL|X|Y|Z))|("
635  "DISPLACEMENT|ROTATE))(.*)") %
636  problem_name % field_name % field_name)
637  .str();
638  } else {
639  regex_str = (boost::format("%s_%s_((FIX_(ALL|X|Y|Z))|("
640  "DISPLACEMENT|ROTATE))(.*)") %
641  problem_name % field_name)
642  .str();
643  }
644 
645  if (std::regex_match(bc_id, std::regex(regex_str))) {
646 
647  auto bc = bcMapByBlockName.at(bc_id);
648 
649  if (auto disp_bc = bc->dispBcPtr) {
650  if (disp_bc->data.flag1) {
651  ents_to_remove[0].merge(bc->bcEnts);
652  }
653  if (disp_bc->data.flag2) {
654  ents_to_remove[1].merge(bc->bcEnts);
655  }
656  if (disp_bc->data.flag3) {
657  ents_to_remove[2].merge(bc->bcEnts);
658  }
659  if (disp_bc->data.flag4) {
660  ents_to_remove[1].merge(bc->bcEnts);
661  ents_to_remove[2].merge(bc->bcEnts);
662  }
663  if (disp_bc->data.flag5) {
664  ents_to_remove[0].merge(bc->bcEnts);
665  ents_to_remove[2].merge(bc->bcEnts);
666  }
667  if (disp_bc->data.flag6) {
668  ents_to_remove[0].merge(bc->bcEnts);
669  ents_to_remove[1].merge(bc->bcEnts);
670  }
671  } else {
672  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
673  "BC type not implemented");
674  }
675  }
676  }
677 
678  auto remove_dofs_on_ents = [&](const Range &ents, const int lo,
679  const int hi) {
680  if (is_distributed_mesh)
681  return prb_mng->removeDofsOnEntities(problem_name, field_name, ents, lo,
682  hi);
683  else
684  return prb_mng->removeDofsOnEntitiesNotDistributed(
685  problem_name, field_name, ents, lo, hi);
686  };
687 
688  CHKERR remove_dofs_on_ents(ents_to_remove[0], 0, 0);
689  CHKERR remove_dofs_on_ents(ents_to_remove[1], 1, 1);
690  CHKERR remove_dofs_on_ents(ents_to_remove[2], 2, 2);
691 
693 }
694 
695 template <>
697 BcManager::removeBlockDOFsOnEntities<BcScalarMeshsetType<BLOCKSET>>(
698  const std::string problem_name, const std::string block_name,
699  const std::string field_name, bool get_low_dim_ents,
700  bool is_distributed_mesh) {
701  Interface &m_field = cOre;
702  auto prb_mng = m_field.getInterface<ProblemsManager>();
704 
705  CHKERR pushMarkDOFsOnEntities<BcScalarMeshsetType<BLOCKSET>>(
706  problem_name, block_name, field_name, get_low_dim_ents);
707 
708  Range ents_to_remove;
709 
710  for (auto m :
711 
713  BLOCKSET | UNKNOWNNAME)) {
714 
715  std::string bc_id = problem_name + "_" + field_name + "_" + m->getName();
716 
717  auto str = boost::format("%s_%s_%s(.*)")
718 
719  % problem_name % field_name % block_name;
720 
721  if (std::regex_match(bc_id, std::regex(str.str()))) {
722 
723  auto bc = bcMapByBlockName.at(bc_id);
724 
725  if (auto temp_bc = bc->tempBcPtr) {
726  if (temp_bc->data.flag1) {
727  ents_to_remove.merge(bc->bcEnts);
728  }
729  } else {
730  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
731  "BC type not implemented");
732  }
733  }
734  }
735 
736  auto remove_dofs_on_ents = [&](const Range &ents, const int lo,
737  const int hi) {
738  if (is_distributed_mesh)
739  return prb_mng->removeDofsOnEntities(problem_name, field_name, ents, lo,
740  hi);
741  else
742  return prb_mng->removeDofsOnEntitiesNotDistributed(
743  problem_name, field_name, ents, lo, hi);
744  };
745 
746  CHKERR remove_dofs_on_ents(ents_to_remove, 0, MAX_DOFS_ON_ENTITY);
747 
749 }
750 
751 template <>
753 BcManager::pushMarkDOFsOnEntities<BcMeshsetType<DISPLACEMENTSET>>(
754  const std::string problem_name, const std::string field_name,
755  bool get_low_dim_ents, bool block_name_field_prefix) {
756  Interface &m_field = cOre;
757  auto prb_mng = m_field.getInterface<ProblemsManager>();
759 
760  // if(problem_name.size()==0)
761  // MOFEM_LOG("BcMngWorld", Sev::warning)
762  // << "Argument problem_name has no effect";
763 
764  if (block_name_field_prefix)
765  MOFEM_LOG("BcMngWorld", Sev::warning)
766  << "Argument block_name_field_prefix=true has no effect";
767 
768  auto fix_disp = [&]() {
770 
771  auto iterate_meshsets = [&](auto &&meshset_vec_ptr) {
773  for (auto m : meshset_vec_ptr) {
774  auto bc = boost::make_shared<BCs>();
775  CHKERR m_field.get_moab().get_entities_by_handle(m->getMeshset(),
776  bc->bcEnts, true);
777  bc->dispBcPtr = boost::make_shared<DisplacementCubitBcData>();
778  CHKERR m->getBcDataStructure(*(bc->dispBcPtr));
779 
780  MOFEM_LOG("BcMngWorld", Sev::verbose)
781  << "Found block DISPLACEMENTSET id = " << m->getMeshsetId();
782  MOFEM_LOG("BcMngWorld", Sev::verbose) << *bc->dispBcPtr;
783 
784  MOFEM_LOG("BcMngSync", Sev::noisy)
785  << "Found block DISPLACEMENTSET id = " << m->getMeshsetId()
786  << " nb. of entities " << bc->bcEnts.size()
787  << " highest dim of entities "
788  << BcManagerImplTools::get_dim(bc->bcEnts);
789  MOFEM_LOG("BcMngSync", Sev::noisy) << *bc->dispBcPtr;
790  MOFEM_LOG_SEVERITY_SYNC(m_field.get_comm(), Sev::noisy);
791 
792  if (problem_name.size()) {
793 
794  if (bc->dispBcPtr->data.flag1)
795  CHKERR prb_mng->modifyMarkDofs(problem_name, ROW, field_name, 0, 0,
796  ProblemsManager::MarkOP::OR, 1,
797  bc->bcMarkers);
798  if (bc->dispBcPtr->data.flag2)
799  CHKERR prb_mng->modifyMarkDofs(problem_name, ROW, field_name, 1, 1,
800  ProblemsManager::MarkOP::OR, 1,
801  bc->bcMarkers);
802  if (bc->dispBcPtr->data.flag3)
803  CHKERR prb_mng->modifyMarkDofs(problem_name, ROW, field_name, 2, 2,
804  ProblemsManager::MarkOP::OR, 1,
805  bc->bcMarkers);
806  if (bc->dispBcPtr->data.flag4) {
807 
808  CHKERR prb_mng->modifyMarkDofs(problem_name, ROW, field_name, 1, 1,
809  ProblemsManager::MarkOP::OR, 1,
810  bc->bcMarkers);
811  CHKERR prb_mng->modifyMarkDofs(problem_name, ROW, field_name, 2, 2,
812  ProblemsManager::MarkOP::OR, 1,
813  bc->bcMarkers);
814  }
815  if (bc->dispBcPtr->data.flag5) {
816 
817  CHKERR prb_mng->modifyMarkDofs(problem_name, ROW, field_name, 0, 0,
818  ProblemsManager::MarkOP::OR, 1,
819  bc->bcMarkers);
820  CHKERR prb_mng->modifyMarkDofs(problem_name, ROW, field_name, 2, 2,
821  ProblemsManager::MarkOP::OR, 1,
822  bc->bcMarkers);
823  }
824  if (bc->dispBcPtr->data.flag6) {
825 
826  CHKERR prb_mng->modifyMarkDofs(problem_name, ROW, field_name, 0, 0,
827  ProblemsManager::MarkOP::OR, 1,
828  bc->bcMarkers);
829  CHKERR prb_mng->modifyMarkDofs(problem_name, ROW, field_name, 1, 1,
830  ProblemsManager::MarkOP::OR, 1,
831  bc->bcMarkers);
832  }
833  }
834 
835  if (get_low_dim_ents) {
836  auto low_dim_ents =
837  BcManagerImplTools::get_adj_ents(m_field.get_moab(), bc->bcEnts);
838  bc->bcEnts.swap(low_dim_ents);
839  }
840 
841  CHKERR m_field.getInterface<CommInterface>()->synchroniseEntities(
842  bc->bcEnts);
843  if (problem_name.size())
844  CHKERR prb_mng->markDofs(problem_name, ROW, ProblemsManager::AND,
845  bc->bcEnts, bc->bcMarkers);
846 
847  const std::string bc_id =
848  problem_name + "_" + field_name + "_DISPLACEMENTSET" +
849  boost::lexical_cast<std::string>(m->getMeshsetId());
850  bcMapByBlockName[bc_id] = bc;
851  }
853  };
854 
855  CHKERR iterate_meshsets(
856 
859 
860  );
861 
863  };
864 
865  CHKERR fix_disp();
866 
868 }
869 
870 template <>
871 MoFEMErrorCode BcManager::pushMarkDOFsOnEntities<BcMeshsetType<TEMPERATURESET>>(
872  const std::string problem_name, const std::string field_name,
873  bool get_low_dim_ents, bool block_name_field_prefix) {
874  Interface &m_field = cOre;
875  auto prb_mng = m_field.getInterface<ProblemsManager>();
877 
878  if (block_name_field_prefix)
879  MOFEM_LOG("BcMngWorld", Sev::warning)
880  << "Argument block_name_field_prefix=true has no effect";
881 
882  auto fix_temp = [&]() {
884 
885  auto iterate_meshsets = [&](auto &&meshset_vec_ptr) {
887  for (auto m : meshset_vec_ptr) {
888  auto bc = boost::make_shared<BCs>();
889  CHKERR m_field.get_moab().get_entities_by_handle(m->getMeshset(),
890  bc->bcEnts, true);
891  bc->tempBcPtr = boost::make_shared<TemperatureCubitBcData>();
892  CHKERR m->getBcDataStructure(*(bc->tempBcPtr));
893 
894  MOFEM_LOG("BcMngWorld", Sev::verbose)
895  << "Found block TEMPERATURESET id = " << m->getMeshsetId();
896  MOFEM_LOG("BcMngWorld", Sev::verbose) << *bc->tempBcPtr;
897 
898  CHKERR prb_mng->modifyMarkDofs(
899  problem_name, ROW, field_name, 0, MAX_DOFS_ON_ENTITY,
900  ProblemsManager::MarkOP::OR, 1, bc->bcMarkers);
901 
902  if (get_low_dim_ents) {
903  auto low_dim_ents =
904  BcManagerImplTools::get_adj_ents(m_field.get_moab(), bc->bcEnts);
905  bc->bcEnts.swap(low_dim_ents);
906  }
907 
908  CHKERR m_field.getInterface<CommInterface>()->synchroniseEntities(
909  bc->bcEnts);
910  CHKERR prb_mng->markDofs(problem_name, ROW, ProblemsManager::AND,
911  bc->bcEnts, bc->bcMarkers);
912 
913  const std::string bc_id =
914  problem_name + "_" + field_name + "_TEMPERATURESET" +
915  boost::lexical_cast<std::string>(m->getMeshsetId());
916  bcMapByBlockName[bc_id] = bc;
917  }
919  };
920 
921  CHKERR iterate_meshsets(
922 
925 
926  );
927 
929  };
930 
931  CHKERR fix_temp();
932 
934 }
935 
936 template <>
937 MoFEMErrorCode BcManager::pushMarkDOFsOnEntities<BcMeshsetType<HEATFLUXSET>>(
938  const std::string problem_name, const std::string field_name,
939  bool get_low_dim_ents, bool block_name_field_prefix) {
940  Interface &m_field = cOre;
941  auto prb_mng = m_field.getInterface<ProblemsManager>();
943 
944  if (block_name_field_prefix)
945  MOFEM_LOG("BcMngWorld", Sev::warning)
946  << "Argument block_name_field_prefix=true has no effect";
947 
948  auto fix_disp = [&]() {
950 
951  auto iterate_meshsets = [&](auto &&meshset_vec_ptr) {
953  for (auto m : meshset_vec_ptr) {
954  auto bc = boost::make_shared<BCs>();
955  CHKERR m_field.get_moab().get_entities_by_handle(m->getMeshset(),
956  bc->bcEnts, true);
957  bc->heatFluxBcPtr = boost::make_shared<HeatFluxCubitBcData>();
958  CHKERR m->getBcDataStructure(*(bc->heatFluxBcPtr));
959 
960  CHKERR prb_mng->modifyMarkDofs(
961  problem_name, ROW, field_name, 0, MAX_DOFS_ON_ENTITY,
962  ProblemsManager::MarkOP::OR, 1, bc->bcMarkers);
963 
964  MOFEM_LOG("BcMngWorld", Sev::inform)
965  << "Found block HEATFLUX id = " << m->getMeshsetId();
966  MOFEM_LOG("BcMngWorld", Sev::inform) << *bc->heatFluxBcPtr;
967 
968  if (get_low_dim_ents) {
969  auto low_dim_ents =
970  BcManagerImplTools::get_adj_ents(m_field.get_moab(), bc->bcEnts);
971  bc->bcEnts.swap(low_dim_ents);
972  }
973 
974  CHKERR m_field.getInterface<CommInterface>()->synchroniseEntities(
975  bc->bcEnts);
976  CHKERR prb_mng->markDofs(problem_name, ROW, ProblemsManager::AND,
977  bc->bcEnts, bc->bcMarkers);
978 
979  const std::string bc_id =
980  problem_name + "_" + field_name + "_HEATFLUXSET" +
981  boost::lexical_cast<std::string>(m->getMeshsetId());
982  bcMapByBlockName[bc_id] = bc;
983  }
985  };
986 
987  CHKERR iterate_meshsets(
988 
990  HEATFLUXSET)
991 
992  );
993 
995  };
996 
997  CHKERR fix_disp();
998 
1000 }
1001 
1002 template <>
1004 BcManager::pushMarkDOFsOnEntities<BcDisplacementMeshsetType<BLOCKSET>>(
1005  const std::string problem_name, const std::string field_name,
1006  bool get_low_dim_ents, bool block_name_field_prefix) {
1008 
1009  auto mark_dofs = [&](const string block_name, const int &idx_0,
1010  const int &idx_1) {
1012  if (block_name_field_prefix) {
1013  const string field_block = field_name + "_" + block_name;
1014  CHKERR pushMarkDOFsOnEntities(problem_name, field_block, field_name,
1015  idx_0, idx_1, get_low_dim_ents);
1016  } else {
1017 
1018  CHKERR pushMarkDOFsOnEntities(problem_name, block_name, field_name, idx_0,
1019  idx_1, get_low_dim_ents);
1020  }
1022  };
1023 
1024  // displacement
1025  CHKERR mark_dofs("FIX_X", 0, 0);
1026  CHKERR mark_dofs("FIX_Y", 1, 1);
1027  CHKERR mark_dofs("FIX_Z", 2, 2);
1028  CHKERR mark_dofs("FIX_ALL", 0, MAX_DOFS_ON_ENTITY);
1029 
1030  // rotation
1031  CHKERR mark_dofs("ROTATE_X", 1, 1);
1032  CHKERR mark_dofs("ROTATE_X", 2, 2);
1033  CHKERR mark_dofs("ROTATE_Y", 0, 0);
1034  CHKERR mark_dofs("ROTATE_Y", 2, 2);
1035  CHKERR mark_dofs("ROTATE_Z", 0, 0);
1036  CHKERR mark_dofs("ROTATE_Z", 1, 1);
1037  CHKERR mark_dofs("ROTATE_ALL", 0, MAX_DOFS_ON_ENTITY);
1038 
1039  std::string regex_str;
1040  if (block_name_field_prefix) {
1041  regex_str = (boost::format("%s_%s_%s_(.*)") % problem_name % field_name %
1042  field_name)
1043  .str();
1044  } else {
1045  regex_str = (boost::format("%s_%s_(.*)") % problem_name % field_name).str();
1046  }
1047 
1048  for (auto &m : bcMapByBlockName) {
1049  auto &bc_id = m.first;
1050  if (std::regex_match(bc_id, std::regex(regex_str))) {
1051  auto &bc = m.second;
1052  if (std::regex_match(bc_id, std::regex("(.*)_FIX_X(.*)"))) {
1053  bc->dispBcPtr = boost::make_shared<DisplacementCubitBcData>();
1054  bc->dispBcPtr->data.flag1 = 1;
1055  if (bc->bcAttributes.empty()) {
1056  bc->dispBcPtr->data.value1 = 0;
1057  MOFEM_LOG("BcMngWorld", Sev::warning)
1058  << "Expected one attribute on block but have "
1059  << bc->bcAttributes.size();
1060  } else if (bc->bcAttributes.size() >= 1) {
1061  bc->dispBcPtr->data.value1 = bc->bcAttributes[0];
1062  }
1063  MOFEM_LOG("BcMngWorld", Sev::inform) << "Add X " << bc_id;
1064  MOFEM_LOG("BcMngWorld", Sev::inform) << *bc->dispBcPtr;
1065  } else if (std::regex_match(bc_id, std::regex("(.*)_FIX_Y(.*)"))) {
1066  bc->dispBcPtr = boost::make_shared<DisplacementCubitBcData>();
1067  bc->dispBcPtr->data.flag2 = 1;
1068  if (bc->bcAttributes.empty()) {
1069  bc->dispBcPtr->data.value2 = 0;
1070  MOFEM_LOG("BcMngWorld", Sev::warning)
1071  << "Expected one attribute on block but have "
1072  << bc->bcAttributes.size();
1073  } else if (bc->bcAttributes.size() == 1) {
1074  bc->dispBcPtr->data.value2 = bc->bcAttributes[0];
1075  } else if (bc->bcAttributes.size() >= 2) {
1076  bc->dispBcPtr->data.value2 = bc->bcAttributes[1];
1077  }
1078  MOFEM_LOG("BcMngWorld", Sev::inform) << "Add Y " << bc_id;
1079  MOFEM_LOG("BcMngWorld", Sev::inform) << *(bc->dispBcPtr);
1080  } else if (std::regex_match(bc_id, std::regex("(.*)_FIX_Z(.*)"))) {
1081  bc->dispBcPtr = boost::make_shared<DisplacementCubitBcData>();
1082  bc->dispBcPtr->data.flag3 = 1;
1083  if (bc->bcAttributes.empty()) {
1084  bc->dispBcPtr->data.value3 = 0;
1085  MOFEM_LOG("BcMngWorld", Sev::warning)
1086  << "Expected one attribute on block but have "
1087  << bc->bcAttributes.size();
1088  } else if (bc->bcAttributes.size() == 1) {
1089  bc->dispBcPtr->data.value3 = bc->bcAttributes[0];
1090  } else if (bc->bcAttributes.size() == 3) {
1091  bc->dispBcPtr->data.value3 = bc->bcAttributes[2];
1092  }
1093  MOFEM_LOG("BcMngWorld", Sev::inform) << "Add Z " << bc_id;
1094  MOFEM_LOG("BcMngWorld", Sev::inform) << *(bc->dispBcPtr);
1095  } else if (std::regex_match(bc_id, std::regex("(.*)_FIX_ALL(.*)"))) {
1096  bc->dispBcPtr = boost::make_shared<DisplacementCubitBcData>();
1097  bc->dispBcPtr->data.flag1 = 1;
1098  bc->dispBcPtr->data.flag2 = 1;
1099  bc->dispBcPtr->data.flag3 = 1;
1100  if (bc->bcAttributes.size() >= 1) {
1101  bc->dispBcPtr->data.value1 = bc->bcAttributes[0];
1102  }
1103  if (bc->bcAttributes.size() >= 2) {
1104  bc->dispBcPtr->data.value2 = bc->bcAttributes[1];
1105  }
1106  if (bc->bcAttributes.size() >= 3) {
1107  bc->dispBcPtr->data.value3 = bc->bcAttributes[2];
1108  }
1109  MOFEM_LOG("BcMngWorld", Sev::inform) << "Add ALL " << bc_id;
1110  MOFEM_LOG("BcMngWorld", Sev::inform) << *(bc->dispBcPtr);
1111  } else if (std::regex_match(bc_id, std::regex("(.*)_ROTATE_X(.*)"))) {
1112  bc->dispBcPtr =
1113  boost::make_shared<DisplacementCubitBcDataWithRotation>();
1114  bc->dispBcPtr->data.flag4 = 1;
1115  bc->dispBcPtr->data.flag5 = 0;
1116  bc->dispBcPtr->data.flag6 = 0;
1117  // for the ROTATE_X block the angles can be specified with either one or
1118  // three attributes, e.g. 1, coords or 1,0,0,coords
1119  if (bc->bcAttributes.empty()) {
1120  bc->dispBcPtr->data.value4 = 0;
1121  MOFEM_LOG("BcMngWorld", Sev::warning)
1122  << "Expected one attribute on block on block (angle (1 or 3), "
1123  "center coords(3) but have "
1124  << bc->bcAttributes.size();
1125  } else if (bc->bcAttributes.size() >= 1) {
1126  bc->dispBcPtr->data.value4 = bc->bcAttributes[0];
1127  }
1128  MOFEM_LOG("BcMngWorld", Sev::inform) << "Add X " << bc_id;
1129  MOFEM_LOG("BcMngWorld", Sev::inform) << *bc->dispBcPtr;
1130  if (bc->bcAttributes.size() == 4 || bc->bcAttributes.size() == 6) {
1131  if (auto ext_disp_bc =
1132  dynamic_cast<DisplacementCubitBcDataWithRotation *>(
1133  bc->dispBcPtr.get())) {
1134  auto &o = ext_disp_bc->rotOffset;
1135  for (int a = 0; a != 3; ++a)
1136  o[a] = bc->bcAttributes[bc->bcAttributes.size() - 3 + a];
1137  MOFEM_LOG("BcMngWorld", Sev::inform)
1138  << "Add Rotate X Center: " << o[0] << " " << o[1] << " "
1139  << o[2];
1140  }
1141  }
1142  } else if (std::regex_match(bc_id, std::regex("(.*)_ROTATE_Y(.*)"))) {
1143  bc->dispBcPtr =
1144  boost::make_shared<DisplacementCubitBcDataWithRotation>();
1145  bc->dispBcPtr->data.flag4 = 0;
1146  bc->dispBcPtr->data.flag5 = 1;
1147  bc->dispBcPtr->data.flag6 = 0;
1148  // for the ROTATE_Y block the angles can be specified with either one or
1149  // three attributes, e.g. 1, coords or 0,1,0,coords
1150  if (bc->bcAttributes.empty()) {
1151  bc->dispBcPtr->data.value5 = 0;
1152  MOFEM_LOG("BcMngWorld", Sev::warning)
1153  << "Expected one attribute on block on block (angle (1 or 3), "
1154  "center coords(3) but have "
1155  << bc->bcAttributes.size();
1156  } else if (bc->bcAttributes.size() == 1 ||
1157  bc->bcAttributes.size() == 4) {
1158  bc->dispBcPtr->data.value5 = bc->bcAttributes[0];
1159  } else if (bc->bcAttributes.size() == 6) {
1160  bc->dispBcPtr->data.value5 = bc->bcAttributes[1];
1161  }
1162  MOFEM_LOG("BcMngWorld", Sev::inform) << "Add Y " << bc_id;
1163  MOFEM_LOG("BcMngWorld", Sev::inform) << *(bc->dispBcPtr);
1164  if (bc->bcAttributes.size() == 4 || bc->bcAttributes.size() == 6) {
1165  if (auto ext_disp_bc =
1166  dynamic_cast<DisplacementCubitBcDataWithRotation *>(
1167  bc->dispBcPtr.get())) {
1168  auto &o = ext_disp_bc->rotOffset;
1169  for (int a = 0; a != 3; ++a)
1170  o[a] = bc->bcAttributes[bc->bcAttributes.size() - 3 + a];
1171  MOFEM_LOG("BcMngWorld", Sev::inform)
1172  << "Add Rotate Y Center: " << o[0] << " " << o[1] << " "
1173  << o[2];
1174  }
1175  }
1176  } else if (std::regex_match(bc_id, std::regex("(.*)_ROTATE_Z(.*)"))) {
1177  bc->dispBcPtr =
1178  boost::make_shared<DisplacementCubitBcDataWithRotation>();
1179  bc->dispBcPtr->data.flag4 = 0;
1180  bc->dispBcPtr->data.flag5 = 0;
1181  bc->dispBcPtr->data.flag6 = 1;
1182  // for the ROTATE_Z block the angles can be specified with either one or
1183  // three attributes, e.g. 1, coords or 0,0,1,coords
1184  if (bc->bcAttributes.empty()) {
1185  bc->dispBcPtr->data.value6 = 0;
1186  MOFEM_LOG("BcMngWorld", Sev::warning)
1187  << "Expected one attribute on block (angle (1 or 3), center "
1188  "coords(3) but have "
1189  << bc->bcAttributes.size();
1190  } else if (bc->bcAttributes.size() == 1 ||
1191  bc->bcAttributes.size() == 4) {
1192  bc->dispBcPtr->data.value6 = bc->bcAttributes[0];
1193  } else if (bc->bcAttributes.size() == 3 ||
1194  bc->bcAttributes.size() == 6) {
1195  bc->dispBcPtr->data.value6 = bc->bcAttributes[2];
1196  }
1197  MOFEM_LOG("BcMngWorld", Sev::inform) << "Add Z " << bc_id;
1198  MOFEM_LOG("BcMngWorld", Sev::inform) << *(bc->dispBcPtr);
1199  if (bc->bcAttributes.size() == 4 || bc->bcAttributes.size() == 6) {
1200  if (auto ext_disp_bc =
1201  dynamic_cast<DisplacementCubitBcDataWithRotation *>(
1202  bc->dispBcPtr.get())) {
1203  auto &o = ext_disp_bc->rotOffset;
1204  for (int a = 0; a != 3; ++a)
1205  o[a] = bc->bcAttributes[bc->bcAttributes.size() - 3 + a];
1206  MOFEM_LOG("BcMngWorld", Sev::inform)
1207  << "Add Rotate Z Center: " << o[0] << " " << o[1] << " "
1208  << o[2];
1209  }
1210  }
1211  } else if (std::regex_match(bc_id, std::regex("(.*)_ROTATE_ALL(.*)"))) {
1212  bc->dispBcPtr =
1213  boost::make_shared<DisplacementCubitBcDataWithRotation>();
1214  bc->dispBcPtr->data.flag4 = 1;
1215  bc->dispBcPtr->data.flag5 = 1;
1216  bc->dispBcPtr->data.flag6 = 1;
1217  if (bc->bcAttributes.size() >= 1) {
1218  bc->dispBcPtr->data.value4 = bc->bcAttributes[0];
1219  }
1220  if (bc->bcAttributes.size() >= 2) {
1221  bc->dispBcPtr->data.value5 = bc->bcAttributes[1];
1222  }
1223  if (bc->bcAttributes.size() >= 3) {
1224  bc->dispBcPtr->data.value6 = bc->bcAttributes[2];
1225  }
1226  MOFEM_LOG("BcMngWorld", Sev::inform) << "Add ALL " << bc_id;
1227  MOFEM_LOG("BcMngWorld", Sev::inform) << *(bc->dispBcPtr);
1228  if (bc->bcAttributes.size() > 3) {
1229  if (auto ext_disp_bc =
1230  dynamic_cast<DisplacementCubitBcDataWithRotation *>(
1231  bc->dispBcPtr.get())) {
1232  auto &o = ext_disp_bc->rotOffset;
1233  for (int a = 0; a != 3; ++a)
1234  o[a] = bc->bcAttributes[3 + a];
1235  MOFEM_LOG("BcMngWorld", Sev::inform)
1236  << "Add Rotate ALL Center: " << o[0] << " " << o[1] << " "
1237  << o[2];
1238  }
1239  }
1240  }
1241  }
1242  }
1243 
1245 }
1246 
1247 template <>
1248 MoFEMErrorCode BcManager::pushMarkDOFsOnEntities<BcScalarMeshsetType<BLOCKSET>>(
1249  const std::string problem_name, const std::string block_name,
1250  const std::string field_name, bool get_low_dim_ents) {
1252 
1253  CHKERR pushMarkDOFsOnEntities(problem_name, block_name, field_name, 0,
1254  MAX_DOFS_ON_ENTITY, get_low_dim_ents);
1255 
1256  auto regex_str =
1257  (boost::format("%s_%s_%s(.*)") % problem_name % field_name % block_name)
1258  .str();
1259 
1260  for (auto &m : bcMapByBlockName) {
1261 
1262  auto &bc_id = m.first;
1263 
1264  if (std::regex_match(bc_id, std::regex(regex_str))) {
1265 
1266  auto &bc = m.second;
1267  bc->tempBcPtr = boost::make_shared<TemperatureCubitBcData>();
1268  bc->tempBcPtr->data.flag1 = 1;
1269  if (bc->bcAttributes.empty()) {
1270  bc->tempBcPtr->data.value1 = 0;
1271  MOFEM_LOG("BcMngWorld", Sev::warning)
1272  << "Expected one attribute on block but have "
1273  << bc->bcAttributes.size();
1274  } else if (bc->bcAttributes.size() >= 1) {
1275  bc->tempBcPtr->data.value1 = bc->bcAttributes[0];
1276  }
1277  }
1278  }
1279 
1281 }
1282 
1283 template <>
1284 MoFEMErrorCode BcManager::pushMarkDOFsOnEntities<DisplacementCubitBcData>(
1285  const std::string problem_name, const std::string field_name,
1286  bool get_low_dim_ents, bool block_name_field_prefix) {
1288  // that marks DOFs and create data when are set by cubit nodesets.
1289  CHKERR pushMarkDOFsOnEntities<BcMeshsetType<DISPLACEMENTSET>>(
1290  problem_name, field_name, get_low_dim_ents, block_name_field_prefix);
1291  // that marks DOFs and create data when are set by blocsket.
1292  CHKERR pushMarkDOFsOnEntities<BcDisplacementMeshsetType<BLOCKSET>>(
1293  problem_name, field_name, get_low_dim_ents, block_name_field_prefix);
1295 }
1296 
1297 template <>
1298 MoFEMErrorCode BcManager::removeBlockDOFsOnEntities<DisplacementCubitBcData>(
1299  const std::string problem_name, const std::string field_name,
1300  bool get_low_dim_ents, bool block_name_field_prefix,
1301  bool is_distributed_mesh) {
1303  // that remove DOFs when are set by cubit nodesets.
1304  CHKERR removeBlockDOFsOnEntities<BcMeshsetType<DISPLACEMENTSET>>(
1305  problem_name, field_name, get_low_dim_ents, block_name_field_prefix,
1306  is_distributed_mesh);
1307  // that remove DOFs when are by blocksets
1308  CHKERR removeBlockDOFsOnEntities<BcDisplacementMeshsetType<BLOCKSET>>(
1309  problem_name, field_name, get_low_dim_ents, block_name_field_prefix,
1310  is_distributed_mesh);
1311  // add more ways to remove bcs when appropiate
1313 }
1314 
1315 template <>
1316 MoFEMErrorCode BcManager::pushMarkDOFsOnEntities<TemperatureCubitBcData>(
1317  const std::string problem_name, const std::string field_name,
1318  bool get_low_dim_ents, bool block_name_field_prefix) {
1320  CHKERR pushMarkDOFsOnEntities<BcMeshsetType<TEMPERATURESET>>(
1321  problem_name, field_name, get_low_dim_ents, block_name_field_prefix);
1322 
1323  auto get_block_name = [&]() {
1324  if (block_name_field_prefix)
1325  return (boost::format("%s_FIX_SCALAR") % field_name).str();
1326  else
1327  return field_name;
1328  };
1329 
1330  CHKERR pushMarkDOFsOnEntities<BcScalarMeshsetType<BLOCKSET>>(
1331  problem_name, get_block_name(), field_name, get_low_dim_ents);
1333 }
1334 
1335 template <>
1336 MoFEMErrorCode BcManager::removeBlockDOFsOnEntities<TemperatureCubitBcData>(
1337  const std::string problem_name, const std::string field_name,
1338  bool get_low_dim_ents, bool block_name_field_prefix,
1339  bool is_distributed_mesh) {
1341  CHKERR removeBlockDOFsOnEntities<BcMeshsetType<TEMPERATURESET>>(
1342  problem_name, field_name, get_low_dim_ents, block_name_field_prefix,
1343  is_distributed_mesh);
1344 
1345  auto get_block_name = [&]() {
1346  if (block_name_field_prefix)
1347  return (boost::format("%s_FIX_SCALAR") % field_name).str();
1348  else
1349  return field_name;
1350  };
1351 
1352  CHKERR removeBlockDOFsOnEntities<BcScalarMeshsetType<BLOCKSET>>(
1353  problem_name, get_block_name(), field_name, get_low_dim_ents,
1354  is_distributed_mesh);
1356 }
1357 
1358 template <>
1359 MoFEMErrorCode BcManager::pushMarkDOFsOnEntities<HeatFluxCubitBcData>(
1360  const std::string problem_name, const std::string field_name,
1361  bool get_low_dim_ents, bool block_name_field_prefix) {
1363  CHKERR pushMarkDOFsOnEntities<BcMeshsetType<HEATFLUXSET>>(
1364  problem_name, field_name, get_low_dim_ents, block_name_field_prefix);
1366 }
1367 
1368 template <>
1369 MoFEMErrorCode BcManager::removeBlockDOFsOnEntities<HeatFluxCubitBcData>(
1370  const std::string problem_name, const std::string field_name,
1371  bool get_low_dim_ents, bool block_name_field_prefix,
1372  bool is_distributed_mesh) {
1374  CHKERR removeBlockDOFsOnEntities<BcMeshsetType<HEATFLUXSET>>(
1375  problem_name, field_name, get_low_dim_ents, block_name_field_prefix,
1376  is_distributed_mesh);
1378 }
1379 
1380 std::pair<std::string, std::string>
1381 BcManager::extractStringFromBlockId(const std::string block_id,
1382  const std::string prb_name) {
1383 
1384  // Assumes that field name is consist with letters and numbers.
1385  // No special characters.
1386  auto field_rgx_str =
1387  (boost::format("%s_([a-zA-Z0-9]*)_(.*)") % prb_name).str();
1388  std::regex field_rgx(field_rgx_str);
1389  std::smatch match_field_name;
1390  std::string field_name;
1391  std::string block_name;
1392 
1393  if (std::regex_search(block_id, match_field_name, field_rgx)) {
1394  field_name = match_field_name[1];
1395  block_name = match_field_name[2];
1396  } else {
1398  "Field name and block name can not be resolved");
1399  }
1400 
1401  return std::make_pair(field_name, block_name);
1402 }
1403 
1404 template <>
1405 MoFEMErrorCode BcManager::pushMarkDOFsOnEntities<BcMeshsetType<FORCESET>>(
1406  const std::string problem_name, const std::string field_name,
1407  bool get_low_dim_ents, bool block_name_field_prefix) {
1408  Interface &m_field = cOre;
1409  auto prb_mng = m_field.getInterface<ProblemsManager>();
1411 
1412  if (problem_name.empty())
1413  MOFEM_LOG("BcMngWorld", Sev::warning) << "Argument problem_name is empty";
1414 
1415  if (block_name_field_prefix)
1416  MOFEM_LOG("BcMngWorld", Sev::warning)
1417  << "Argument block_name_field_prefix=true has no effect";
1418 
1419  auto fix_force = [&]() {
1421 
1422  auto iterate_meshsets = [&](auto &&meshset_vec_ptr) {
1424  for (auto m : meshset_vec_ptr) {
1425  auto bc = boost::make_shared<BCs>();
1426  CHKERR m_field.get_moab().get_entities_by_handle(m->getMeshset(),
1427  bc->bcEnts, true);
1428  bc->forceBcPtr = boost::make_shared<ForceCubitBcData>();
1429  CHKERR m->getBcDataStructure(*(bc->forceBcPtr));
1430 
1431  MOFEM_LOG("BcMngWorld", Sev::verbose)
1432  << "Found block FORCESET id = " << m->getMeshsetId();
1433  MOFEM_LOG("BcMngWorld", Sev::verbose) << *bc->forceBcPtr;
1434 
1435  MOFEM_LOG("BcMngSync", Sev::noisy)
1436  << "Found block FORCESET id = " << m->getMeshsetId()
1437  << " nb. of entities " << bc->bcEnts.size()
1438  << " highest dim of entities "
1439  << BcManagerImplTools::get_dim(bc->bcEnts);
1440  MOFEM_LOG("BcMngSync", Sev::noisy) << *bc->forceBcPtr;
1441  MOFEM_LOG_SEVERITY_SYNC(m_field.get_comm(), Sev::noisy);
1442 
1443  if (problem_name.size()) {
1444 
1445  if (bc->forceBcPtr->data.value2 > 0)
1446  CHKERR prb_mng->modifyMarkDofs(problem_name, ROW, field_name, 0, 0,
1447  ProblemsManager::MarkOP::OR, 1,
1448  bc->bcMarkers);
1449  if (bc->forceBcPtr->data.value3 > 0)
1450  CHKERR prb_mng->modifyMarkDofs(problem_name, ROW, field_name, 1, 1,
1451  ProblemsManager::MarkOP::OR, 1,
1452  bc->bcMarkers);
1453  if (bc->forceBcPtr->data.value4 > 0)
1454  CHKERR prb_mng->modifyMarkDofs(problem_name, ROW, field_name, 2, 2,
1455  ProblemsManager::MarkOP::OR, 1,
1456  bc->bcMarkers);
1457 
1458  if (bc->forceBcPtr->data.value5 > 0) {
1459 
1460  CHKERR prb_mng->modifyMarkDofs(problem_name, ROW, field_name, 1, 1,
1461  ProblemsManager::MarkOP::OR, 1,
1462  bc->bcMarkers);
1463  CHKERR prb_mng->modifyMarkDofs(problem_name, ROW, field_name, 2, 2,
1464  ProblemsManager::MarkOP::OR, 1,
1465  bc->bcMarkers);
1466  }
1467  if (bc->forceBcPtr->data.value5) {
1468 
1469  CHKERR prb_mng->modifyMarkDofs(problem_name, ROW, field_name, 0, 0,
1470  ProblemsManager::MarkOP::OR, 1,
1471  bc->bcMarkers);
1472  CHKERR prb_mng->modifyMarkDofs(problem_name, ROW, field_name, 2, 2,
1473  ProblemsManager::MarkOP::OR, 1,
1474  bc->bcMarkers);
1475  }
1476  if (bc->forceBcPtr->data.value6) {
1477 
1478  CHKERR prb_mng->modifyMarkDofs(problem_name, ROW, field_name, 0, 0,
1479  ProblemsManager::MarkOP::OR, 1,
1480  bc->bcMarkers);
1481  CHKERR prb_mng->modifyMarkDofs(problem_name, ROW, field_name, 1, 1,
1482  ProblemsManager::MarkOP::OR, 1,
1483  bc->bcMarkers);
1484  }
1485  }
1486 
1487  if (get_low_dim_ents) {
1488  auto low_dim_ents =
1489  BcManagerImplTools::get_adj_ents(m_field.get_moab(), bc->bcEnts);
1490  bc->bcEnts.swap(low_dim_ents);
1491  }
1492 
1493  CHKERR m_field.getInterface<CommInterface>()->synchroniseEntities(
1494  bc->bcEnts);
1495  if (problem_name.size())
1496  CHKERR prb_mng->markDofs(problem_name, ROW, ProblemsManager::AND,
1497  bc->bcEnts, bc->bcMarkers);
1498 
1499  const std::string bc_id =
1500  problem_name + "_" + field_name + "_FORCESET" +
1501  boost::lexical_cast<std::string>(m->getMeshsetId());
1502  bcMapByBlockName[bc_id] = bc;
1503  }
1505  };
1506 
1507  CHKERR iterate_meshsets(
1508 
1510  FORCESET)
1511 
1512  );
1513 
1515  };
1516 
1517  CHKERR fix_force();
1518 
1520 }
1521 
1522 template <>
1523 MoFEMErrorCode BcManager::pushMarkDOFsOnEntities<BcForceMeshsetType<BLOCKSET>>(
1524  const std::string problem_name, const std::string field_name,
1525  bool get_low_dim_ents, bool block_name_field_prefix) {
1526  Interface &m_field = cOre;
1527  auto prb_mng = m_field.getInterface<ProblemsManager>();
1529 
1530  if (problem_name.size() == 0)
1531  MOFEM_LOG("BcMngWorld", Sev::warning) << "Argument problem_name is empty";
1532 
1533  auto get_force_block = [&](auto block_name) {
1535 
1536  for (auto m :
1537  m_field.getInterface<MeshsetsManager>()->getCubitMeshsetPtr(std::regex(
1538 
1539  (boost::format("%s(.*)") % block_name).str()
1540 
1541  ))
1542 
1543  ) {
1544 
1545  const auto block_name = m->getName();
1546 
1547  MOFEM_LOG("BcMngWorld", Sev::inform)
1548  << "Found force block " << block_name;
1549 
1550  auto bc = boost::make_shared<BCs>();
1551  CHKERR m_field.get_moab().get_entities_by_handle(m->getMeshset(),
1552  bc->bcEnts, true);
1553 
1554  CHKERR m->getAttributes(bc->bcAttributes);
1555  if (bc->bcAttributes.size() != 3) {
1556  SETERRQ(m_field.get_comm(), MOFEM_DATA_INCONSISTENCY,
1557  "Expect three block attributes for force block");
1558  }
1559 
1560  bc->forceBcPtr = boost::make_shared<ForceCubitBcData>();
1561  // For details look at ForceCubitBcData in
1562  // mofem/src/multi_indices/BCData.hpp
1563  bc->forceBcPtr->data.value1 = 1;
1564  bc->forceBcPtr->data.value3 = bc->bcAttributes[0];
1565  bc->forceBcPtr->data.value4 = bc->bcAttributes[1];
1566  bc->forceBcPtr->data.value5 = bc->bcAttributes[2];
1567 
1568  MOFEM_LOG("BcMngWorld", Sev::inform) << *bc->forceBcPtr;
1569  MOFEM_LOG("BcMngSync", Sev::noisy)
1570  << "Found block FORCESET id = " << m->getMeshsetId()
1571  << " nb. of entities " << bc->bcEnts.size()
1572  << " highest dim of entities "
1573  << BcManagerImplTools::get_dim(bc->bcEnts);
1574  MOFEM_LOG("BcMngSync", Sev::noisy) << *bc->forceBcPtr;
1575  MOFEM_LOG_SEVERITY_SYNC(m_field.get_comm(), Sev::noisy);
1576 
1577  if (problem_name.size()) {
1578 
1579  if (bc->forceBcPtr->data.value2 > 0)
1580  CHKERR prb_mng->modifyMarkDofs(problem_name, ROW, field_name, 0, 0,
1581  ProblemsManager::MarkOP::OR, 1,
1582  bc->bcMarkers);
1583  if (bc->forceBcPtr->data.value3 > 0)
1584  CHKERR prb_mng->modifyMarkDofs(problem_name, ROW, field_name, 1, 1,
1585  ProblemsManager::MarkOP::OR, 1,
1586  bc->bcMarkers);
1587  if (bc->forceBcPtr->data.value4 > 0)
1588  CHKERR prb_mng->modifyMarkDofs(problem_name, ROW, field_name, 2, 2,
1589  ProblemsManager::MarkOP::OR, 1,
1590  bc->bcMarkers);
1591  }
1592 
1593  if (get_low_dim_ents) {
1594  auto low_dim_ents =
1595  BcManagerImplTools::get_adj_ents(m_field.get_moab(), bc->bcEnts);
1596  bc->bcEnts.swap(low_dim_ents);
1597  }
1598 
1599  CHKERR m_field.getInterface<CommInterface>()->synchroniseEntities(
1600  bc->bcEnts);
1601  if (problem_name.size())
1602  CHKERR prb_mng->markDofs(problem_name, ROW, ProblemsManager::AND,
1603  bc->bcEnts, bc->bcMarkers);
1604 
1605  const std::string bc_id =
1606  problem_name + "_" + field_name + "_" + block_name;
1607  bcMapByBlockName[bc_id] = bc;
1608  }
1610  };
1611 
1612  CHKERR get_force_block("FORCE");
1613 
1615 }
1616 
1617 template <>
1618 MoFEMErrorCode BcManager::pushMarkDOFsOnEntities<ForceCubitBcData>(
1619  const std::string problem_name, const std::string field_name,
1620  bool get_low_dim_ents, bool block_name_field_prefix) {
1622 
1623  CHKERR pushMarkDOFsOnEntities<BcMeshsetType<FORCESET>>(
1624  problem_name, field_name, get_low_dim_ents, block_name_field_prefix);
1625 
1626  CHKERR pushMarkDOFsOnEntities<BcForceMeshsetType<BLOCKSET>>(
1627  problem_name, field_name, get_low_dim_ents, block_name_field_prefix);
1628 
1630 }
1631 
1632 template <>
1633 MoFEMErrorCode BcManager::removeBlockDOFsOnEntities<BcMeshsetType<FORCESET>>(
1634  const std::string problem_name, const std::string field_name,
1635  bool get_low_dim_ents, bool block_name_field_prefix,
1636  bool is_distributed_mesh) {
1637  Interface &m_field = cOre;
1638  auto prb_mng = m_field.getInterface<ProblemsManager>();
1640 
1641  CHKERR pushMarkDOFsOnEntities<BcMeshsetType<FORCESET>>(
1642  problem_name, field_name, get_low_dim_ents);
1643 
1644  std::array<Range, 3> ents_to_remove;
1645 
1646  for (auto m :
1647 
1649  FORCESET)) {
1650 
1651  const auto block_name = m->getName();
1652  std::string bc_id = problem_name + "_" + field_name + "_" + block_name;
1653 
1654  auto str = boost::format("%s_%s_%s(.*)")
1655 
1656  % problem_name % field_name % block_name;
1657 
1658  if (std::regex_match(bc_id, std::regex(str.str()))) {
1659 
1660  auto bc = bcMapByBlockName.at(bc_id);
1661 
1662  if (auto force_bc = bc->forceBcPtr) {
1663  if (force_bc->data.value3 > 0) {
1664  ents_to_remove[0].merge(bc->bcEnts);
1665  }
1666  if (force_bc->data.value4 > 0) {
1667  ents_to_remove[1].merge(bc->bcEnts);
1668  }
1669  if (force_bc->data.value5 > 0) {
1670  ents_to_remove[2].merge(bc->bcEnts);
1671  }
1672  } else {
1673  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
1674  "BC type not implemented");
1675  }
1676  }
1677  }
1678 
1679  auto remove_dofs_on_ents = [&](const Range &ents, const int lo,
1680  const int hi) {
1681  if (is_distributed_mesh)
1682  return prb_mng->removeDofsOnEntities(problem_name, field_name, ents, lo,
1683  hi);
1684  else
1685  return prb_mng->removeDofsOnEntitiesNotDistributed(
1686  problem_name, field_name, ents, lo, hi);
1687  };
1688 
1689  CHKERR remove_dofs_on_ents(ents_to_remove[0], 0, 0);
1690  CHKERR remove_dofs_on_ents(ents_to_remove[1], 1, 1);
1691  CHKERR remove_dofs_on_ents(ents_to_remove[2], 2, 2);
1692 
1694 }
1695 
1696 template <>
1698 BcManager::removeBlockDOFsOnEntities<BcForceMeshsetType<BLOCKSET>>(
1699  const std::string problem_name, const std::string field_name,
1700  bool get_low_dim_ents, bool block_name_field_prefix,
1701  bool is_distributed_mesh) {
1702  Interface &m_field = cOre;
1703  auto prb_mng = m_field.getInterface<ProblemsManager>();
1705 
1706  CHKERR pushMarkDOFsOnEntities<BcForceMeshsetType<BLOCKSET>>(
1707  problem_name, field_name, get_low_dim_ents);
1708 
1709  std::array<Range, 3> ents_to_remove;
1710 
1711  for (auto m :
1712 
1714  BLOCKSET | UNKNOWNNAME)) {
1715 
1716  const auto block_name = m->getName();
1717  std::string bc_id = problem_name + "_" + field_name + "_" + block_name;
1718 
1719  auto str = boost::format("%s_%s_%s(.*)")
1720 
1721  % problem_name % field_name % block_name;
1722 
1723  if (std::regex_match(bc_id, std::regex(str.str()))) {
1724 
1725  auto bc = bcMapByBlockName.at(bc_id);
1726 
1727  if (auto force_bc = bc->forceBcPtr) {
1728  if (force_bc->data.value3 > 0) {
1729  ents_to_remove[0].merge(bc->bcEnts);
1730  }
1731  if (force_bc->data.value4 > 0) {
1732  ents_to_remove[1].merge(bc->bcEnts);
1733  }
1734  if (force_bc->data.value5 > 0) {
1735  ents_to_remove[2].merge(bc->bcEnts);
1736  }
1737  } else {
1738  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
1739  "BC type not implemented");
1740  }
1741  }
1742  }
1743 
1744  auto remove_dofs_on_ents = [&](const Range &ents, const int lo,
1745  const int hi) {
1746  if (is_distributed_mesh)
1747  return prb_mng->removeDofsOnEntities(problem_name, field_name, ents, lo,
1748  hi);
1749  else
1750  return prb_mng->removeDofsOnEntitiesNotDistributed(
1751  problem_name, field_name, ents, lo, hi);
1752  };
1753 
1754  CHKERR remove_dofs_on_ents(ents_to_remove[0], 0, 0);
1755  CHKERR remove_dofs_on_ents(ents_to_remove[1], 1, 1);
1756  CHKERR remove_dofs_on_ents(ents_to_remove[2], 2, 2);
1757 
1759 }
1760 
1761 template <>
1762 MoFEMErrorCode BcManager::removeBlockDOFsOnEntities<ForceCubitBcData>(
1763  const std::string problem_name, const std::string field_name,
1764  bool get_low_dim_ents, bool block_name_field_prefix,
1765  bool is_distributed_mesh) {
1767 
1768  CHKERR removeBlockDOFsOnEntities<BcMeshsetType<FORCESET>>(
1769  problem_name, field_name, get_low_dim_ents, block_name_field_prefix,
1770  is_distributed_mesh);
1771 
1772  CHKERR removeBlockDOFsOnEntities<BcForceMeshsetType<BLOCKSET>>(
1773  problem_name, field_name, get_low_dim_ents, block_name_field_prefix,
1774  is_distributed_mesh);
1775 
1777 };
1778 
1780 
1781  const std::string problem_name, const std::string block_name,
1782  const std::string field_name, int bridge_dim, int lo, int hi
1783 
1784 ) {
1785  Interface &m_field = cOre;
1787 
1788  if (problem_name.empty())
1790 
1791  auto iterate_meshsets = [&](auto &&meshset_vec_ptr) {
1792  auto prb_mng = m_field.getInterface<ProblemsManager>();
1794  for (auto m : meshset_vec_ptr) {
1795  auto bc = boost::make_shared<BCs>();
1796  CHKERR m_field.get_moab().get_entities_by_handle(m->getMeshset(),
1797  bc->bcEnts, true);
1798  CHKERR m->getAttributes(bc->bcAttributes);
1799 
1800  bc->dofsViewPtr = boost::make_shared<BCs::DofsView>();
1801 
1802  CHKERR m_field.getInterface<CommInterface>()->synchroniseEntities(
1803  bc->bcEnts);
1804  CHKERR prb_mng->getSideDofsOnBrokenSpaceEntities(
1805  *(bc->dofsViewPtr), problem_name, ROW, field_name, bc->bcEnts,
1806  bridge_dim, lo, hi);
1807  CHKERR prb_mng->markDofs(problem_name, ROW, *(bc->dofsViewPtr),
1808  ProblemsManager::OR, bc->bcMarkers);
1809 
1810  MOFEM_LOG("BcMngWorld", Sev::inform)
1811  << "Found block " << m->getName() << " number of attributes "
1812  << bc->bcAttributes.size() << " number of entities "
1813  << bc->bcEnts.size();
1814 
1815  const std::string bc_id =
1816  problem_name + "_" + field_name + "_" + m->getName();
1817  bcMapByBlockName[bc_id] = bc;
1818  }
1820  };
1821 
1822  CHKERR iterate_meshsets(
1823 
1824  m_field.getInterface<MeshsetsManager>()->getCubitMeshsetPtr(std::regex(
1825 
1826  (boost::format("%s(.*)") % block_name).str()
1827 
1828  ))
1829 
1830  );
1831 
1833 }
1834 
1835 MoFEMErrorCode BcManager::removeSideDOFs(const std::string problem_name,
1836  const std::string block_name,
1837  const std::string field_name,
1838  int bridge_dim, int lo, int hi,
1839  bool is_distributed_mesh) {
1840  Interface &m_field = cOre;
1842 
1843  CHKERR pushMarkSideDofs(problem_name, block_name, field_name, bridge_dim, lo,
1844  hi);
1845 
1846  auto iterate_meshsets = [&](auto &&meshset_vec_ptr) {
1848  auto prb_mng = m_field.getInterface<ProblemsManager>();
1849  for (auto m : meshset_vec_ptr) {
1850  const std::string bc_id =
1851  problem_name + "_" + field_name + "_" + m->getName();
1852  auto &bc = bcMapByBlockName.at(bc_id);
1853  CHKERR prb_mng->removeDofs(problem_name, ROW, *(bc->dofsViewPtr), lo, hi);
1854  CHKERR prb_mng->removeDofs(problem_name, COL, *(bc->dofsViewPtr), lo, hi);
1855  }
1857  };
1858 
1859  if (is_distributed_mesh) {
1860 
1861  CHKERR iterate_meshsets(
1862 
1863  m_field.getInterface<MeshsetsManager>()->getCubitMeshsetPtr(std::regex(
1864 
1865  (boost::format("%s(.*)") % block_name).str()
1866 
1867  ))
1868 
1869  );
1870  } else {
1871  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "Not implemented");
1872  }
1873 
1875 }
1876 
1877 } // namespace MoFEM
MoFEMFunctionReturnHot
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:460
MoFEM::UnknownInterface::getInterface
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
Definition: UnknownInterface.hpp:93
CHK_MOAB_THROW
#define CHK_MOAB_THROW(err, msg)
Check error code of MoAB function and throw MoFEM exception.
Definition: definitions.h:589
MoFEM::BcManager::extractStringFromBlockId
static std::pair< std::string, std::string > extractStringFromBlockId(const std::string block_id, const std::string prb_name)
Extract block name and block name form block id.
Definition: BcManager.cpp:1381
SIDESET
@ SIDESET
Definition: definitions.h:160
MoFEM::LogManager::checkIfChannelExist
static bool checkIfChannelExist(const std::string channel)
Check if channel exist.
Definition: LogManager.cpp:404
MoFEM::CoreTmp< 0 >
Core (interface) class.
Definition: Core.hpp:82
MoFEM::BcManager::popMarkDOFsOnEntities
boost::shared_ptr< BCs > popMarkDOFsOnEntities(const std::string block_name)
Get bc data and remove element.
Definition: BcManager.cpp:349
MOFEM_LOG_SEVERITY_SYNC
#define MOFEM_LOG_SEVERITY_SYNC(comm, severity)
Synchronise "SYNC" on curtain severity level.
Definition: LogManager.hpp:352
MoFEM::BcManager::query_interface
MoFEMErrorCode query_interface(boost::typeindex::type_index type_index, UnknownInterface **iface) const
Definition: BcManager.cpp:34
MoFEM::ProblemsManager
Problem manager is used to build and partition problems.
Definition: ProblemsManager.hpp:21
MoFEM::BcManager::pushMarkSideDofs
MoFEMErrorCode pushMarkSideDofs(const std::string problem_name, const std::string block_name, const std::string field_name, int bridge_dim, int lo, int hi)
Mark side DOFs.
Definition: BcManager.cpp:1779
MoFEM::BcManager::BcManager
BcManager(const MoFEM::Core &core)
Definition: BcManager.cpp:41
MoFEM::CoreInterface::get_comm
virtual MPI_Comm & get_comm() const =0
MoFEM::BcManager::cOre
MoFEM::Core & cOre
Definition: BcManager.hpp:368
CHK_THROW_MESSAGE
#define CHK_THROW_MESSAGE(err, msg)
Check and throw MoFEM exception.
Definition: definitions.h:609
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::BcManager::getOptions
MoFEMErrorCode getOptions()
get options
Definition: BcManager.cpp:65
MoFEM::BcManager::removeBlockDOFsOnEntities
MoFEMErrorCode removeBlockDOFsOnEntities(const std::string problem_name, const std::string block_name, const std::string field_name, int lo, int hi, bool get_low_dim_ents=true, bool is_distributed_mesh=true)
Remove DOFs from problem.
Definition: BcManager.cpp:73
MoFEM::BcManager::getMergedBlocksMarker
BcMarkerPtr getMergedBlocksMarker(std::vector< std::regex > bc_regex_vec)
Get the Merged Boundary Marker object.
Definition: BcManager.cpp:374
MoFEM::ProblemsManager::OR
@ OR
Definition: ProblemsManager.hpp:416
MoFEM::BcManager::removeSideDOFs
MoFEMErrorCode removeSideDOFs(const std::string problem_name, const std::string block_name, const std::string field_name, int bridge_dim, int lo, int hi, bool is_distributed_mesh=true)
Remove side DOFs.
Definition: BcManager.cpp:1835
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::DeprecatedCoreInterface
Deprecated interface functions.
Definition: DeprecatedCoreInterface.hpp:16
UNKNOWNNAME
@ UNKNOWNNAME
Definition: definitions.h:171
ROW
@ ROW
Definition: definitions.h:136
MoFEM::Interface
DeprecatedCoreInterface Interface
Definition: Interface.hpp:2010
MoFEM::BcManager::pushMarkDOFsOnEntities
MoFEMErrorCode pushMarkDOFsOnEntities(const std::string problem_name, const std::string block_name, const std::string field_name, int lo, int hi, bool get_low_dim_ents=true)
Mark block DOFs.
Definition: BcManager.cpp:111
NODESET
@ NODESET
Definition: definitions.h:159
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:548
MoFEM::CoreInterface::get_moab
virtual moab::Interface & get_moab()=0
MoFEM
implementation of Data Operators for Forces and Sources
Definition: Common.hpp:10
MoFEM::ISManager
Section manager is used to create indexes and sections.
Definition: ISManager.hpp:23
a
constexpr double a
Definition: approx_sphere.cpp:30
MoFEM::BcManager
Simple interface for fast problem set-up.
Definition: BcManager.hpp:25
FORCESET
@ FORCESET
Definition: definitions.h:164
MoFEM::BcManagerImplTools::get_dim
auto get_dim(const Range &ents)
Definition: BcManager.cpp:10
MoFEM::BcManager::getBcMapByBlockName
BcMapByBlockName & getBcMapByBlockName()
Get the bc map.
Definition: BcManager.hpp:243
DISPLACEMENTSET
@ DISPLACEMENTSET
Definition: definitions.h:163
MoFEM::LogManager::getStrmSync
static boost::shared_ptr< std::ostream > getStrmSync()
Get the strm sync object.
Definition: LogManager.cpp:348
COL
@ COL
Definition: definitions.h:136
MoFEM::BcManager::getBlockIS
SmartPetscObj< IS > getBlockIS(const std::string block_prefix, const std::string block_name, const std::string field_name, const std::string problem_name, int lo, int hi, SmartPetscObj< IS > is_expand=SmartPetscObj< IS >())
Get block IS.
Definition: BcManager.cpp:405
MoFEM::LogManager::getStrmWorld
static boost::shared_ptr< std::ostream > getStrmWorld()
Get the strm world object.
Definition: LogManager.cpp:344
MoFEM::LogManager::getStrmSelf
static boost::shared_ptr< std::ostream > getStrmSelf()
Get the strm self object.
Definition: LogManager.cpp:340
MoFEM::BcManager::getMergedBlocksRange
Range getMergedBlocksRange(std::vector< std::regex > bc_regex_vec)
Merge block ranges.
Definition: BcManager.cpp:359
MoFEM::BcManager::BcMarkerPtr
boost::shared_ptr< std::vector< char unsigned > > BcMarkerPtr
Definition: BcManager.hpp:227
MOFEM_LOG_TAG
#define MOFEM_LOG_TAG(channel, tag)
Tag channel.
Definition: LogManager.hpp:339
i
FTensor::Index< 'i', SPACE_DIM > i
Definition: hcurl_divergence_operator_2d.cpp:27
field_name
constexpr auto field_name
Definition: poisson_2d_homogeneous.cpp:13
MoFEM::ProblemsManager::AND
@ AND
Definition: ProblemsManager.hpp:416
MoFEM::BcManager::addBlockDOFsToMPCs
MoFEMErrorCode addBlockDOFsToMPCs(const std::string problem_name, const std::string field_name, bool get_low_dim_ents=false, bool block_name_field_prefix=false, bool is_distributed_mesh=false)
Definition: BcManager.cpp:185
MoFEM::UnknownInterface
base class for all interface classes
Definition: UnknownInterface.hpp:34
MAX_DOFS_ON_ENTITY
#define MAX_DOFS_ON_ENTITY
Maximal number of DOFs on entity.
Definition: definitions.h:249
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_LOG
#define MOFEM_LOG(channel, severity)
Log.
Definition: LogManager.hpp:308
MoFEM::CommInterface
Managing BitRefLevels.
Definition: CommInterface.hpp:21
BLOCKSET
@ BLOCKSET
Definition: definitions.h:161
MoFEM::BcManager::bcMapByBlockName
BcMapByBlockName bcMapByBlockName
Definition: BcManager.hpp:370
MoFEM::Exceptions::ierr
static MoFEMErrorCodeGeneric< PetscErrorCode > ierr
Definition: Exceptions.hpp:76
MOFEM_DATA_INCONSISTENCY
@ MOFEM_DATA_INCONSISTENCY
Definition: definitions.h:31
TEMPERATURESET
@ TEMPERATURESET
Definition: definitions.h:168
MoFEM::MeshsetsManager
Interface for managing meshsets containing materials and boundary conditions.
Definition: MeshsetsManager.hpp:104
sdf_hertz_2d_axisymm_plane.d
float d
Definition: sdf_hertz_2d_axisymm_plane.py:4
m
FTensor::Index< 'm', 3 > m
Definition: shallow_wave.cpp:80
MoFEMFunctionBeginHot
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:453
HEATFLUXSET
@ HEATFLUXSET
Definition: definitions.h:169
MoFEM::DisplacementCubitBcDataWithRotation
A specialized version of DisplacementCubitBcData that includes an additional rotation offset.
Definition: EssentialDisplacementCubitBcData.hpp:24
MoFEM::MeshsetsManager::getCubitMeshsetPtr
MoFEMErrorCode getCubitMeshsetPtr(const int ms_id, const CubitBCType cubit_bc_type, const CubitMeshSets **cubit_meshset_ptr) const
get cubit meshset
Definition: MeshsetsManager.cpp:578
MoFEM::LogManager::setLog
static LoggerType & setLog(const std::string channel)
Set ans resset chanel logger.
Definition: LogManager.cpp:389
MoFEM::SmartPetscObj< IS >
MoFEM::MeshsetsManager::checkMeshset
bool checkMeshset(const int ms_id, const CubitBCType cubit_bc_type) const
check for CUBIT Id and CUBIT type
Definition: MeshsetsManager.cpp:360
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:429
CHKERRG
#define CHKERRG(n)
Check error code of MoFEM/MOAB/PETSc function.
Definition: definitions.h:496
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:359
l
FTensor::Index< 'l', 3 > l
Definition: matrix_function.cpp:21
MoFEM::BcManagerImplTools::get_adj_ents
auto get_adj_ents(moab::Interface &moab, const Range &ents)
Definition: BcManager.cpp:17
MoFEM::MPC::COUPLING
@ COUPLING