v0.14.0
Simple.cpp
Go to the documentation of this file.
1 /** \file Simple.cpp
2  * \brief Implementation of simple interface
3  * \ingroup mofem_simple_interface
4  */
5 
6 namespace MoFEM {
7 
8 MoFEMErrorCode Simple::query_interface(boost::typeindex::type_index type_index,
9  UnknownInterface **iface) const {
10  *iface = const_cast<Simple *>(this);
11  return 0;
12 }
13 
14 template <int DIM>
16  static_assert(DIM == 2 || DIM == 3, "not implemented");
17  return MOFEM_NOT_IMPLEMENTED;
18 }
19 
20 template <>
21 MoFEMErrorCode Simple::setSkeletonAdjacency<2>(std::string fe_name) {
22  Interface &m_field = cOre;
24 
26  boost::make_shared<ParentFiniteElementAdjacencyFunctionSkeleton<1>>(
29  fe_name, MBEDGE, *parentAdjSkeletonFunctionDim1);
30 
32 }
33 
34 template <>
35 MoFEMErrorCode Simple::setSkeletonAdjacency<3>(std::string fe_name) {
36  Interface &m_field = cOre;
38 
40  boost::make_shared<ParentFiniteElementAdjacencyFunctionSkeleton<2>>(
42 
44  fe_name, MBTRI, *parentAdjSkeletonFunctionDim2);
46  fe_name, MBQUAD, *parentAdjSkeletonFunctionDim2);
47 
49 }
50 
51 template <>
52 MoFEMErrorCode Simple::setSkeletonAdjacency<-1>(std::string fe_name) {
54 
55  switch (getDim()) {
56  case 1:
57  THROW_MESSAGE("Not implemented");
58  case 2:
59  return setSkeletonAdjacency<2>(fe_name);
60  case 3:
61  return setSkeletonAdjacency<3>(fe_name);
62  default:
63  THROW_MESSAGE("Not implemented");
64  }
66 }
68  static_assert(DIM == 2 || DIM == 3, "not implemented");
69  return MOFEM_NOT_IMPLEMENTED;
70 }
71 
72 template <> MoFEMErrorCode Simple::setParentAdjacency<3>() {
73  Interface &m_field = cOre;
75 
77  boost::make_shared<ParentFiniteElementAdjacencyFunction<3>>(
80  boost::make_shared<ParentFiniteElementAdjacencyFunction<2>>(
82 
87  if (addBoundaryFE || !boundaryFields.empty()) {
92  }
93  if (addSkeletonFE || !skeletonFields.empty()) {
98  }
99 
101 }
102 
103 template <> MoFEMErrorCode Simple::setParentAdjacency<2>() {
104  Interface &m_field = cOre;
106 
108  boost::make_shared<ParentFiniteElementAdjacencyFunction<2>>(
111  boost::make_shared<ParentFiniteElementAdjacencyFunction<1>>(
113 
118  if (addBoundaryFE || !boundaryFields.empty())
121  if (addSkeletonFE || !skeletonFields.empty())
124 
126 }
127 
130  switch (getDim()) {
131  case 1:
132  THROW_MESSAGE("Not implemented");
133  case 2:
134  return setParentAdjacency<2>();
135  case 3:
136  return setParentAdjacency<3>();
137  default:
138  THROW_MESSAGE("Not implemented");
139  }
141 }
142 
143 MoFEMErrorCode Simple::setSkeletonAdjacency(int dim, std::string fe_name) {
145  if (dim == -1)
146  dim = getDim();
147 
148  if (fe_name.empty())
149  fe_name = skeletonFE;
150 
151  switch (dim) {
152  case 2:
153  return setSkeletonAdjacency<2>(fe_name);
154  case 3:
155  return setSkeletonAdjacency<3>(fe_name);
156  default:
157  SETERRQ(PETSC_COMM_WORLD, MOFEM_NOT_IMPLEMENTED, "Not implemented");
158  }
160 }
161 
162 Simple::Simple(const Core &core)
163  : cOre(const_cast<Core &>(core)), bitLevel(BitRefLevel().set(0)),
164  bitLevelMask(BitRefLevel().set()), meshSet(0), boundaryMeshset(0),
165  skeletonMeshset(0), nameOfProblem("SimpleProblem"), domainFE("dFE"),
166  boundaryFE("bFE"), skeletonFE("sFE"), dIm(-1), addSkeletonFE(false),
167  addBoundaryFE(false), addParentAdjacencies(false),
168  bitAdjParent(BitRefLevel().set()), bitAdjParentMask(BitRefLevel().set()),
169  bitAdjEnt(BitRefLevel().set()), bitAdjEntMask(BitRefLevel().set()) {
170  PetscLogEventRegister("SimpleSetUp", 0, &MOFEM_EVENT_SimpleSetUP);
171  PetscLogEventRegister("SimpleLoadMesh", 0, &MOFEM_EVENT_SimpleLoadMesh);
172  PetscLogEventRegister("SimpleBuildFields", 0, &MOFEM_EVENT_SimpleBuildFields);
173  PetscLogEventRegister("SimpleBuildFiniteElements", 0,
175  PetscLogEventRegister("SimpleSetUp", 0, &MOFEM_EVENT_SimpleBuildProblem);
176  PetscLogEventRegister("SimpleKSPSolve", 0, &MOFEM_EVENT_SimpleKSPSolve);
177  strcpy(meshFileName, "mesh.h5m");
178 }
179 
181  PetscBool flg = PETSC_TRUE;
183  ierr = PetscOptionsBegin(PETSC_COMM_WORLD, "", "Simple interface options",
184  "none");
185  CHKERRG(ierr);
186  ierr = PetscOptionsString("-file_name", "file name", "", "mesh.h5m",
187  meshFileName, 255, &flg);
188  CHKERRG(ierr);
189  ierr = PetscOptionsEnd();
190  CHKERRG(ierr);
192 }
193 
194 MoFEMErrorCode Simple::loadFile(const std::string options,
195  const std::string mesh_file_name,
196  LoadFileFunc loadFunc) {
197  Interface &m_field = cOre;
199  PetscLogEventBegin(MOFEM_EVENT_SimpleLoadMesh, 0, 0, 0, 0);
200 
201  if (!mesh_file_name.empty())
202  strcpy(meshFileName, mesh_file_name.c_str());
203 
204  // Invoke the boost function, passing in required arguments
205  CHKERR loadFunc(m_field, meshFileName, options.c_str());
206  CHKERR m_field.rebuild_database();
207 
208  // determine problem dimension
209  if (dIm == -1) {
210  int nb_ents_3d;
211  CHKERR m_field.get_moab().get_number_entities_by_dimension(
212  meshSet, 3, nb_ents_3d, true);
213  if (nb_ents_3d > 0) {
214  dIm = 3;
215  } else {
216  int nb_ents_2d;
217  CHKERR m_field.get_moab().get_number_entities_by_dimension(
218  meshSet, 2, nb_ents_2d, true);
219  if (nb_ents_2d > 0) {
220  dIm = 2;
221  } else {
222  dIm = 1;
223  }
224  }
225  }
226 
227  if (!boundaryMeshset)
229  if (!skeletonMeshset)
231  if (addSkeletonFE)
233 
234  if (bitLevel.any()) {
235  Range ents;
236  CHKERR m_field.get_moab().get_entities_by_dimension(meshSet, dIm, ents,
237  true);
238  CHKERR m_field.getInterface<BitRefManager>()->setBitRefLevel(ents, bitLevel,
239  false);
240  } else {
241  MOFEM_LOG("WORLD", Sev::warning) << "BitRefLevel is none and not set";
242  CHKERR m_field.getInterface<BitRefManager>()->addToDatabaseBitRefLevelByDim(
243  dIm, BitRefLevel().set(), BitRefLevel().set());
244  }
245 
246  PetscLogEventEnd(MOFEM_EVENT_SimpleLoadMesh, 0, 0, 0, 0);
248 }
249 
252  Interface &m_field = cOre;
253  if (m_field.get_comm_size() == 1)
255  else
256  CHKERR loadFile("PARALLEL=READ_PART;"
257  "PARALLEL_RESOLVE_SHARED_ENTS;"
258  "PARTITION=PARALLEL_PARTITION;",
259  meshFileName);
261 }
262 
264 Simple::addDomainField(const std::string &name, const FieldSpace space,
265  const FieldApproximationBase base,
266  const FieldCoefficientsNumber nb_of_cooficients,
267  const TagType tag_type, const enum MoFEMTypes bh,
268  int verb) {
269 
270  Interface &m_field = cOre;
272  CHKERR m_field.add_field(name, space, base, nb_of_cooficients, tag_type, bh,
273  verb);
274  if (space == NOFIELD)
275  noFieldFields.push_back(name);
276  else
277  domainFields.push_back(name);
279 }
280 
282 Simple::addBoundaryField(const std::string &name, const FieldSpace space,
283  const FieldApproximationBase base,
284  const FieldCoefficientsNumber nb_of_cooficients,
285  const TagType tag_type, const enum MoFEMTypes bh,
286  int verb) {
287  Interface &m_field = cOre;
289  CHKERR m_field.add_field(name, space, base, nb_of_cooficients, tag_type, bh,
290  verb);
291  boundaryFields.push_back(name);
292  if (space == NOFIELD)
293  SETERRQ(
294  PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
295  "NOFIELD space for boundary filed not implemented in Simple interface");
297 }
298 
300 Simple::addSkeletonField(const std::string &name, const FieldSpace space,
301  const FieldApproximationBase base,
302  const FieldCoefficientsNumber nb_of_cooficients,
303  const TagType tag_type, const enum MoFEMTypes bh,
304  int verb) {
305 
306  Interface &m_field = cOre;
308  CHKERR m_field.add_field(name, space, base, nb_of_cooficients, tag_type, bh,
309  verb);
310  skeletonFields.push_back(name);
311  if (space == NOFIELD)
312  SETERRQ(
313  PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
314  "NOFIELD space for boundary filed not implemented in Simple interface");
315 
317 }
318 
320 Simple::addDataField(const std::string &name, const FieldSpace space,
321  const FieldApproximationBase base,
322  const FieldCoefficientsNumber nb_of_cooficients,
323  const TagType tag_type, const enum MoFEMTypes bh,
324  int verb) {
325 
326  Interface &m_field = cOre;
328  CHKERR m_field.add_field(name, space, base, nb_of_cooficients, tag_type, bh,
329  verb);
330  if (space == NOFIELD)
331  noFieldDataFields.push_back(name);
332  else
333  dataFields.push_back(name);
335 }
336 
337 MoFEMErrorCode Simple::removeDomainField(const std::string &name) {
339 
340  auto remove_field_from_list = [&](auto &vec) {
341  auto it = std::find(vec.begin(), vec.end(), name);
342  if (it != vec.end())
343  vec.erase(it);
344  };
345 
346  remove_field_from_list(noFieldFields);
347  remove_field_from_list(domainFields);
348 
350 }
351 
352 MoFEMErrorCode Simple::removeBoundaryField(const std::string &name) {
354 
355  auto remove_field_from_list = [&](auto &vec) {
356  auto it = std::find(vec.begin(), vec.end(), name);
357  if (it != vec.end())
358  vec.erase(it);
359  };
360 
361  remove_field_from_list(boundaryFields);
362 
364 }
365 
366 MoFEMErrorCode Simple::removeSkeletonField(const std::string &name) {
368 
369  auto remove_field_from_list = [&](auto &vec) {
370  auto it = std::find(vec.begin(), vec.end(), name);
371  if (it != vec.end())
372  vec.erase(it);
373  };
374 
375  remove_field_from_list(skeletonFields);
376 
378 }
379 
381  Interface &m_field = cOre;
383 
384  auto clear_rows_and_cols = [&](auto &fe_name) {
386  auto fe_ptr = m_field.get_finite_elements();
387  auto &fe_by_name = const_cast<FiniteElement_multiIndex *>(fe_ptr)
388  ->get<FiniteElement_name_mi_tag>();
389  auto it_fe = fe_by_name.find(fe_name);
390  if (it_fe != fe_by_name.end()) {
391 
392  if (!fe_by_name.modify(it_fe, FiniteElement_row_change_bit_reset()))
393  SETERRQ(m_field.get_comm(), MOFEM_OPERATION_UNSUCCESSFUL,
394  "modification unsuccessful");
395 
396  if (!fe_by_name.modify(it_fe, FiniteElement_col_change_bit_reset()))
397  SETERRQ(m_field.get_comm(), MOFEM_OPERATION_UNSUCCESSFUL,
398  "modification unsuccessful");
399  }
401  };
402  CHKERR clear_rows_and_cols(domainFE);
403  CHKERR clear_rows_and_cols(boundaryFE);
404  CHKERR clear_rows_and_cols(skeletonFE);
405 
406  // Define finite elements
408 
409  auto add_fields = [&](auto &fe_name, auto &fields) {
411  for (auto &field : fields) {
412  CHKERR m_field.modify_finite_element_add_field_row(fe_name, field);
413  CHKERR m_field.modify_finite_element_add_field_col(fe_name, field);
414  CHKERR m_field.modify_finite_element_add_field_data(fe_name, field);
415  }
417  };
418 
419  auto add_data_fields = [&](auto &fe_name, auto &fields) {
421  for (auto &field : fields)
422  CHKERR m_field.modify_finite_element_add_field_data(fe_name, field);
424  };
425 
426  CHKERR add_fields(domainFE, domainFields);
427  CHKERR add_data_fields(domainFE, dataFields);
428  CHKERR add_fields(domainFE, noFieldFields);
429  CHKERR add_data_fields(domainFE, noFieldDataFields);
430 
431  if (addBoundaryFE || !boundaryFields.empty()) {
433  CHKERR add_fields(boundaryFE, domainFields);
434  if (!boundaryFields.empty())
435  CHKERR add_fields(boundaryFE, boundaryFields);
436  CHKERR add_data_fields(boundaryFE, dataFields);
437  CHKERR add_data_fields(boundaryFE, noFieldDataFields);
438  CHKERR add_fields(boundaryFE, noFieldFields);
439  }
440  if (addSkeletonFE || !skeletonFields.empty()) {
442  CHKERR add_fields(skeletonFE, domainFields);
443  if (!skeletonFields.empty())
444  CHKERR add_fields(skeletonFE, skeletonFields);
445  CHKERR add_data_fields(skeletonFE, dataFields);
446  CHKERR add_data_fields(skeletonFE, noFieldDataFields);
447  CHKERR add_fields(skeletonFE, noFieldFields);
448  }
450 }
451 
452 MoFEMErrorCode Simple::defineProblem(const PetscBool is_partitioned) {
453  Interface &m_field = cOre;
455  // Create dm instance
456  dM = createDM(m_field.get_comm(), "DMMOFEM");
457  // set dm data structure which created mofem data structures
458  CHKERR DMMoFEMCreateMoFEM(dM, &m_field, nameOfProblem.c_str(), bitLevel,
459  bitLevelMask);
460  CHKERR DMSetFromOptions(dM);
462  if (addBoundaryFE || !boundaryFields.empty()) {
464  }
465  if (addSkeletonFE || !skeletonFields.empty()) {
467  }
469  CHKERR DMMoFEMSetIsPartitioned(dM, is_partitioned);
471 }
472 
474  const int order, const Range *ents) {
476  fieldsOrder.emplace_back(field_name, order,
477  ents == NULL ? Range() : Range(*ents),
478  ents == NULL ? false : true);
480 }
481 
483  Interface &m_field = cOre;
485  PetscLogEventBegin(MOFEM_EVENT_SimpleBuildFields, 0, 0, 0, 0);
486 
487  auto comm_interface_ptr = m_field.getInterface<CommInterface>();
488  auto bit_ref_ptr = m_field.getInterface<BitRefManager>();
489 
490  // Add entities to the fields
491  auto add_ents_to_field = [&](auto meshset, auto dim, auto &fields) {
493  for (auto &field : fields) {
494  CHKERR m_field.add_ents_to_field_by_dim(meshset, dim, field);
495  CHKERR comm_interface_ptr->synchroniseFieldEntities(field, 0);
496  }
498  };
499 
500  auto make_no_field_ents = [&](auto &fields) {
502  for (auto &field : fields) {
503  std::array<double, 6> coords = {0, 0, 0, 0, 0, 0};
504  CHKERR m_field.create_vertices_and_add_to_field(field, coords.data(), 1);
505  CHKERR comm_interface_ptr->makeFieldEntitiesMultishared(field, 0);
506  CHKERR bit_ref_ptr->setFieldEntitiesBitRefLevel(field, bitLevel);
507  }
509  };
510 
511  CHKERR add_ents_to_field(meshSet, dIm, domainFields);
512  CHKERR add_ents_to_field(meshSet, dIm, dataFields);
513  CHKERR add_ents_to_field(boundaryMeshset, dIm - 1, boundaryFields);
514  CHKERR add_ents_to_field(skeletonMeshset, dIm - 1, skeletonFields);
515 
516  std::set<std::string> nofield_fields;
517  for (auto &f : noFieldFields)
518  nofield_fields.insert(f);
519  for (auto &f : noFieldDataFields)
520  nofield_fields.insert(f);
521 
522  CHKERR make_no_field_ents(nofield_fields);
523 
524  // Set order
525 
526  auto get_field_ptr = [&](auto &f) { return m_field.get_field_structure(f); };
527  for (auto &t : fieldsOrder) {
528  const auto f = std::get<0>(t);
529  const auto order = std::get<1>(t);
530 
531  MOFEM_TAG_AND_LOG("WORLD", Sev::inform, "Simple")
532  << "Set order to field " << f << " order " << order;
533  MOFEM_LOG_CHANNEL("SYNC");
534  if (std::get<3>(t)) {
535  MOFEM_TAG_AND_LOG("SYNC", Sev::verbose, "Simple")
536  << "To ents: " << std::endl
537  << std::get<2>(t) << std::endl;
538  }
539  MOFEM_LOG_SEVERITY_SYNC(m_field.get_comm(), Sev::verbose);
540 
541  if (std::get<3>(t)) {
542 
543  CHKERR m_field.set_field_order(std::get<2>(t), f, order);
544 
545  } else {
546  auto f_ptr = get_field_ptr(f);
547 
548  if (f_ptr->getSpace() == H1) {
549  if (f_ptr->getApproxBase() == AINSWORTH_BERNSTEIN_BEZIER_BASE) {
550  CHKERR m_field.set_field_order(meshSet, MBVERTEX, f, order);
551  } else {
552  CHKERR m_field.set_field_order(meshSet, MBVERTEX, f, 1);
553  }
554  }
555 
556  for (auto d = 1; d <= dIm; ++d) {
557  for (EntityType t = CN::TypeDimensionMap[d].first;
558  t <= CN::TypeDimensionMap[d].second; ++t) {
559  CHKERR m_field.set_field_order(meshSet, t, f, order);
560  }
561  }
562  }
563  }
564  MOFEM_LOG_CHANNEL("WORLD");
565  // Build fields
566  CHKERR m_field.build_fields();
567  PetscLogEventEnd(MOFEM_EVENT_SimpleBuildFields, 0, 0, 0, 0);
569 }
570 
572  Interface &m_field = cOre;
574  PetscLogEventBegin(MOFEM_EVENT_SimpleBuildFiniteElements, 0, 0, 0, 0);
575  // Add finite elements
577  true);
579  if (addBoundaryFE || boundaryFields.size()) {
581  boundaryFE, true);
583  }
584  if (addSkeletonFE || skeletonFields.size()) {
586  skeletonFE, true);
588  }
589  for (std::vector<std::string>::iterator fit = otherFEs.begin();
590  fit != otherFEs.end(); ++fit) {
591  CHKERR m_field.build_finite_elements(*fit);
592  }
593  PetscLogEventEnd(MOFEM_EVENT_SimpleBuildFiniteElements, 0, 0, 0, 0);
595 }
596 
598  Interface &m_field = cOre;
600  PetscLogEventBegin(MOFEM_EVENT_SimpleBuildProblem, 0, 0, 0, 0);
602  // Set problem by the DOFs on the fields rather that by adding DOFs on the
603  // elements
604  m_field.getInterface<ProblemsManager>()->buildProblemFromFields = PETSC_TRUE;
606  m_field.getInterface<ProblemsManager>()->buildProblemFromFields = PETSC_FALSE;
607  PetscLogEventEnd(MOFEM_EVENT_SimpleBuildProblem, 0, 0, 0, 0);
609 }
610 
611 MoFEMErrorCode Simple::setUp(const PetscBool is_partitioned) {
613 
614  PetscLogEventBegin(MOFEM_EVENT_SimpleSetUP, 0, 0, 0, 0);
615 
617 
618  if (addSkeletonFE || !skeletonFields.empty())
620 
623 
624  CHKERR defineProblem(is_partitioned);
628 
629  PetscLogEventEnd(MOFEM_EVENT_SimpleSetUP, 0, 0, 0, 0);
631 }
632 
634  Interface &m_field = cOre;
636  PetscLogEventBegin(MOFEM_EVENT_SimpleBuildProblem, 0, 0, 0, 0);
637 
638  if (!only_dm) {
641  if (addSkeletonFE || !skeletonFields.empty())
646  }
647 
649 
650  const Problem *problem_ptr;
651  CHKERR DMMoFEMGetProblemPtr(dM, &problem_ptr);
652  const auto problem_name = problem_ptr->getName();
653  CHKERR m_field.modify_problem_ref_level_set_bit(problem_name, bitLevel);
654  CHKERR m_field.modify_problem_mask_ref_level_set_bit(problem_name,
655  bitLevelMask);
656 
657  // Set problem by the DOFs on the fields rather that by adding DOFs on the
658  // elements
659  m_field.getInterface<ProblemsManager>()->buildProblemFromFields = PETSC_TRUE;
661  m_field.getInterface<ProblemsManager>()->buildProblemFromFields = PETSC_FALSE;
662 
663  PetscLogEventEnd(MOFEM_EVENT_SimpleBuildProblem, 0, 0, 0, 0);
665 }
666 
669  CHKERR PetscObjectReference(getPetscObject(dM.get()));
670  *dm = dM.get();
672 }
673 
674 /**
675  * @brief Delete dm
676  *
677  * @return MoFEMErrorCode
678  */
681  dM.reset();
683 }
684 
685 /**
686  * @brief Delete finite elements
687  *
688  * @return MoFEMErrorCode
689  */
691  Interface &m_field = cOre;
693  for (auto fe : {domainFE, boundaryFE, skeletonFE}) {
694  if (m_field.check_finite_element(fe)) {
695  CHKERR m_field.delete_finite_element(fe);
696  }
697  }
699 }
700 
702 Simple::addFieldToEmptyFieldBlocks(const std::string row_field,
703  const std::string col_field) const {
704  Interface &m_field = cOre;
707  getProblemName(), row_field, col_field);
709 }
710 
712  Interface &m_field = cOre;
714  ParallelComm *pcomm =
715  ParallelComm::get_pcomm(&m_field.get_moab(), MYPCOMM_INDEX);
716 
717  auto get_skin = [&](auto meshset) {
718  // filter not owned entities, those are not on boundary
719 
720  Range domain_ents;
721  CHKERR m_field.get_moab().get_entities_by_dimension(meshset, dIm,
722  domain_ents, true);
723  CHKERR pcomm->filter_pstatus(domain_ents,
724  PSTATUS_SHARED | PSTATUS_MULTISHARED,
725  PSTATUS_NOT, -1, nullptr);
726 
727  Skinner skin(&m_field.get_moab());
728  Range domain_skin;
729  CHKERR skin.find_skin(0, domain_ents, false, domain_skin);
730  CHKERR pcomm->filter_pstatus(domain_skin,
731  PSTATUS_SHARED | PSTATUS_MULTISHARED,
732  PSTATUS_NOT, -1, nullptr);
733  return domain_skin;
734  };
735 
736  auto create_boundary_meshset = [&](auto &&domain_skin) {
738  // create boundary meshset
739  if (boundaryMeshset != 0) {
741  }
742  CHKERR m_field.get_moab().create_meshset(MESHSET_SET, boundaryMeshset);
743  CHKERR m_field.get_moab().add_entities(boundaryMeshset, domain_skin);
744  for (int dd = 0; dd != dIm - 1; dd++) {
745  Range adj;
746  CHKERR m_field.get_moab().get_adjacencies(domain_skin, dd, false, adj,
747  moab::Interface::UNION);
748  CHKERR m_field.get_moab().add_entities(boundaryMeshset, adj);
749  }
751  };
752 
753  CHKERR create_boundary_meshset(get_skin(meshSet));
754 
756 }
757 
759  Interface &m_field = cOre;
761 
762  ParallelComm *pcomm =
763  ParallelComm::get_pcomm(&m_field.get_moab(), MYPCOMM_INDEX);
764 
765  auto create_skeleton_meshset = [&](auto meshset) {
767  // create boundary meshset
768  if (skeletonMeshset != 0) {
770  }
771  Range boundary_ents, skeleton_ents;
772  CHKERR m_field.get_moab().get_entities_by_dimension(boundaryMeshset,
773  dIm - 1, boundary_ents);
774  Range domain_ents;
775  CHKERR m_field.get_moab().get_entities_by_dimension(meshset, dIm,
776  domain_ents, true);
777  CHKERR m_field.get_moab().get_adjacencies(
778  domain_ents, dIm - 1, false, skeleton_ents, moab::Interface::UNION);
779  skeleton_ents = subtract(skeleton_ents, boundary_ents);
780  CHKERR pcomm->filter_pstatus(skeleton_ents, PSTATUS_NOT_OWNED, PSTATUS_NOT,
781  -1, nullptr);
782  CHKERR m_field.get_moab().create_meshset(MESHSET_SET, skeletonMeshset);
783  CHKERR m_field.get_moab().add_entities(skeletonMeshset, skeleton_ents);
785  };
786 
787  CHKERR create_skeleton_meshset(meshSet);
788 
790 }
791 
793  Interface &m_field = cOre;
795  MOFEM_LOG("WORLD", Sev::verbose) << "Exchange ghost cells";
796 
797  ParallelComm *pcomm =
798  ParallelComm::get_pcomm(&m_field.get_moab(), MYPCOMM_INDEX);
799  if (pcomm == NULL)
800  pcomm = new ParallelComm(&m_field.get_moab(), m_field.get_comm());
801 
802  CHKERR pcomm->exchange_ghost_cells(getDim(), getDim() - 1, 1,
803  3 /**get all adjacent ghosted entities */,
804  true, true, meshSet ? &meshSet : nullptr);
805 
806  Range shared;
807  CHKERR m_field.get_moab().get_entities_by_dimension(0, dIm, shared);
808  for (auto d = dIm - 1; d >= 0; --d) {
809  CHKERR m_field.get_moab().get_adjacencies(shared, d, false, shared,
810  moab::Interface::UNION);
811  }
812  CHKERR pcomm->filter_pstatus(shared, PSTATUS_SHARED | PSTATUS_MULTISHARED,
813  PSTATUS_OR, -1, &shared);
814  Tag part_tag = pcomm->part_tag();
815  CHKERR pcomm->exchange_tags(part_tag, shared);
816  CHKERR m_field.getInterface<CommInterface>()->resolveParentEntities(shared,
817  VERBOSE);
818 
820 }
821 
822 } // namespace MoFEM
MoFEMFunctionReturnHot
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:447
MoFEM::UnknownInterface::getInterface
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface refernce to pointer of interface.
Definition: UnknownInterface.hpp:93
MoFEM::CoreInterface::modify_finite_element_adjacency_table
virtual MoFEMErrorCode modify_finite_element_adjacency_table(const std::string &fe_name, const EntityType type, ElementAdjacencyFunct function)=0
modify finite element table, only for advanced user
MoFEM::Simple::boundaryFields
std::vector< std::string > boundaryFields
boundary fields
Definition: Simple.hpp:518
MoFEM::Simple::parentAdjFunctionDim1
boost::shared_ptr< ParentFiniteElementAdjacencyFunction< 1 > > parentAdjFunctionDim1
Definition: Simple.hpp:555
MoFEM::Simple::setSkeletonAdjacency
MoFEMErrorCode setSkeletonAdjacency(int dim=-1, std::string fe_name="")
Set the skeleton adjacency object.
Definition: Simple.cpp:143
MYPCOMM_INDEX
#define MYPCOMM_INDEX
default communicator number PCOMM
Definition: definitions.h:215
MoFEM::Simple::otherFEs
std::vector< std::string > otherFEs
Other finite elements.
Definition: Simple.hpp:533
MoFEM::CoreTmp< 0 >
Core (interface) class.
Definition: Core.hpp:82
MoFEM::Simple::addParentAdjacencies
bool addParentAdjacencies
If set to true parent adjacencies are build.
Definition: Simple.hpp:510
H1
@ H1
continuous field
Definition: definitions.h:85
MoFEM::Simple::addBoundaryFE
bool addBoundaryFE
Add boundary FE.
Definition: Simple.hpp:509
MOFEM_LOG_SEVERITY_SYNC
#define MOFEM_LOG_SEVERITY_SYNC(comm, severity)
Synchronise "SYNC" on curtain severity level.
Definition: LogManager.hpp:352
MoFEM::Simple::exchangeGhostCells
MoFEMErrorCode exchangeGhostCells()
Definition: Simple.cpp:792
sdf_hertz.d
float d
Definition: sdf_hertz.py:5
MoFEM::Simple::fieldsOrder
std::list< std::tuple< std::string, int, Range, bool > > fieldsOrder
Definition: Simple.hpp:525
MoFEM::CoreInterface::modify_problem_ref_level_set_bit
virtual MoFEMErrorCode modify_problem_ref_level_set_bit(const std::string &name_problem, const BitRefLevel &bit)=0
set ref level for problem
MOFEM_LOG_CHANNEL
#define MOFEM_LOG_CHANNEL(channel)
Set and reset channel.
Definition: LogManager.hpp:284
MoFEM::ProblemsManager
Problem manager is used to build and partition problems.
Definition: ProblemsManager.hpp:21
MoFEM::Simple::Simple
Simple(const MoFEM::Core &core)
Definition: Simple.cpp:162
MoFEM::Simple::parentAdjFunctionDim2
boost::shared_ptr< ParentFiniteElementAdjacencyFunction< 2 > > parentAdjFunctionDim2
Definition: Simple.hpp:553
MoFEM::Simple::query_interface
MoFEMErrorCode query_interface(boost::typeindex::type_index type_index, UnknownInterface **iface) const
Definition: Simple.cpp:8
MoFEM::Simple::noFieldFields
std::vector< std::string > noFieldFields
NOFIELD field name.
Definition: Simple.hpp:521
MoFEM::Simple::meshSet
EntityHandle meshSet
domain meshset
Definition: Simple.hpp:504
MoFEM::FiniteElement_row_change_bit_reset
Reset field from row.
Definition: FEMultiIndices.hpp:954
MoFEM::CoreInterface::get_comm
virtual MPI_Comm & get_comm() const =0
MoFEM::CoreInterface::modify_finite_element_add_field_row
virtual MoFEMErrorCode modify_finite_element_add_field_row(const std::string &fe_name, const std::string name_row)=0
set field row which finite element use
MoFEM::Simple::getDM
SmartPetscObj< DM > getDM()
Return smart DM object.
Definition: Simple.hpp:275
MoFEM::Simple::createBoundaryMeshset
MoFEMErrorCode createBoundaryMeshset()
Definition: Simple.cpp:711
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::Simple::meshFileName
char meshFileName[255]
mesh file name
Definition: Simple.hpp:535
MoFEM::Simple::MOFEM_EVENT_SimpleKSPSolve
PetscLogEvent MOFEM_EVENT_SimpleKSPSolve
Definition: Simple.hpp:502
MoFEM::Simple::buildProblem
MoFEMErrorCode buildProblem()
Build problem.
Definition: Simple.cpp:597
MoFEM::Exceptions::MoFEMErrorCode
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:56
MoFEM::CoreInterface::create_vertices_and_add_to_field
virtual MoFEMErrorCode create_vertices_and_add_to_field(const std::string name, const double coords[], int size, int verb=DEFAULT_VERBOSITY)=0
Create a vertices and add to field object.
MoFEM::Simple::MOFEM_EVENT_SimpleBuildFields
PetscLogEvent MOFEM_EVENT_SimpleBuildFields
Definition: Simple.hpp:499
MoFEM::getPetscObject
PetscObject getPetscObject(T obj)
Definition: PetscSmartObj.hpp:9
MoFEM::Simple::loadFile
MoFEMErrorCode loadFile(const std::string options, const std::string mesh_file_name, LoadFileFunc loadFunc=defaultLoadFileFunc)
Load mesh file.
Definition: Simple.cpp:194
MoFEM::DMSetUp_MoFEM
PetscErrorCode DMSetUp_MoFEM(DM dm)
Definition: DMMoFEM.cpp:1285
MoFEM::Simple::dIm
int dIm
dimension of problem
Definition: Simple.hpp:536
MoFEM::CoreInterface::get_finite_elements
virtual const FiniteElement_multiIndex * get_finite_elements() const =0
Get the finite elements object.
MoFEM::Simple::addSkeletonFE
bool addSkeletonFE
Add skeleton FE.
Definition: Simple.hpp:508
MoFEM::Simple::nameOfProblem
std::string nameOfProblem
problem name
Definition: Simple.hpp:528
MoFEM::Simple::parentAdjSkeletonFunctionDim2
boost::shared_ptr< ParentFiniteElementAdjacencyFunctionSkeleton< 2 > > parentAdjSkeletonFunctionDim2
Definition: Simple.hpp:558
MoFEM::Simple
Simple interface for fast problem set-up.
Definition: Simple.hpp:27
MoFEM::CoreInterface::modify_problem_mask_ref_level_set_bit
virtual MoFEMErrorCode modify_problem_mask_ref_level_set_bit(const std::string &name_problem, const BitRefLevel &bit)=0
set dof mask ref level for problem
THROW_MESSAGE
#define THROW_MESSAGE(msg)
Throw MoFEM exception.
Definition: definitions.h:561
MoFEM::CoreInterface::add_ents_to_field_by_dim
virtual MoFEMErrorCode add_ents_to_field_by_dim(const Range &ents, const int dim, const std::string &name, int verb=DEFAULT_VERBOSITY)=0
Add entities to field meshset.
MoFEM::DMMoFEMAddElement
PetscErrorCode DMMoFEMAddElement(DM dm, std::string fe_name)
add element to dm
Definition: DMMoFEM.cpp:501
order
constexpr int order
Definition: dg_projection.cpp:18
MoFEM::DeprecatedCoreInterface
Deprecated interface functions.
Definition: DeprecatedCoreInterface.hpp:16
MoFEM::Simple::addSkeletonField
MoFEMErrorCode addSkeletonField(const std::string &name, const FieldSpace space, const FieldApproximationBase base, const FieldCoefficientsNumber nb_of_coefficients, const TagType tag_type=MB_TAG_SPARSE, const enum MoFEMTypes bh=MF_ZERO, int verb=-1)
Add field on skeleton.
Definition: Simple.cpp:300
MoFEM::Simple::removeDomainField
MoFEMErrorCode removeDomainField(const std::string &name)
Remove field form domain.
Definition: Simple.cpp:337
MoFEM::Simple::MOFEM_EVENT_SimpleLoadMesh
PetscLogEvent MOFEM_EVENT_SimpleLoadMesh
Definition: Simple.hpp:498
MoFEM::Simple::deleteDM
MoFEMErrorCode deleteDM()
Delete dm.
Definition: Simple.cpp:679
MoFEM::Simple::getOptions
MoFEMErrorCode getOptions()
get options
Definition: Simple.cpp:180
MoFEM::CoreInterface::rebuild_database
virtual MoFEMErrorCode rebuild_database(int verb=DEFAULT_VERBOSITY)=0
Clear database and initialize it once again.
MoFEM::Simple::bitLevelMask
BitRefLevel bitLevelMask
BitRefLevel of the problem.
Definition: Simple.hpp:495
MoFEM::Simple::dataFields
std::vector< std::string > dataFields
Data fields.
Definition: Simple.hpp:520
FieldSpace
FieldSpace
approximation spaces
Definition: definitions.h:82
VERBOSE
@ VERBOSE
Definition: definitions.h:209
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:535
MoFEM::Simple::bitAdjEntMask
BitRefLevel bitAdjEntMask
bit ref level for parent parent
Definition: Simple.hpp:515
MoFEM::CoreInterface::add_finite_element
virtual MoFEMErrorCode add_finite_element(const std::string &fe_name, enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
add finite element
MoFEM::CoreInterface::get_moab
virtual moab::Interface & get_moab()=0
MoFEM
implementation of Data Operators for Forces and Sources
Definition: Common.hpp:10
MoFEM::CoreInterface::modify_finite_element_add_field_col
virtual MoFEMErrorCode modify_finite_element_add_field_col(const std::string &fe_name, const std::string name_row)=0
set field col which finite element use
MoFEM::Simple::deleteFiniteElements
MoFEMErrorCode deleteFiniteElements()
Delete finite elements.
Definition: Simple.cpp:690
MoFEM::createDM
auto createDM(MPI_Comm comm, const std::string dm_type_name)
Creates smart DM object.
Definition: PetscSmartObj.hpp:137
MoFEM::Simple::domainFE
std::string domainFE
domain finite element
Definition: Simple.hpp:529
MoFEM::CoreInterface::build_finite_elements
virtual MoFEMErrorCode build_finite_elements(int verb=DEFAULT_VERBOSITY)=0
Build finite elements.
MoFEM::CoreInterface::add_field
virtual MoFEMErrorCode add_field(const std::string &name, const FieldSpace space, const FieldApproximationBase base, const FieldCoefficientsNumber nb_of_coefficients, const TagType tag_type=MB_TAG_SPARSE, const enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
Add field.
MoFEM::Types::FieldCoefficientsNumber
int FieldCoefficientsNumber
Number of field coefficients.
Definition: Types.hpp:27
MoFEM::Simple::MOFEM_EVENT_SimpleBuildFiniteElements
PetscLogEvent MOFEM_EVENT_SimpleBuildFiniteElements
Definition: Simple.hpp:500
MoFEM::Simple::LoadFileFunc
boost::function< MoFEMErrorCode(Interface &, const char *, const char *)> LoadFileFunc
Definition: Simple.hpp:41
MoFEM::Simple::createSkeletonMeshset
MoFEMErrorCode createSkeletonMeshset()
Definition: Simple.cpp:758
MoFEM::Simple::boundaryFE
std::string boundaryFE
boundary finite element
Definition: Simple.hpp:530
MoFEM::Simple::addDomainField
MoFEMErrorCode addDomainField(const std::string &name, const FieldSpace space, const FieldApproximationBase base, const FieldCoefficientsNumber nb_of_coefficients, const TagType tag_type=MB_TAG_SPARSE, const enum MoFEMTypes bh=MF_ZERO, int verb=-1)
Add field on domain.
Definition: Simple.cpp:264
MoFEM::CoreInterface::delete_finite_element
virtual MoFEMErrorCode delete_finite_element(const std::string name, int verb=DEFAULT_VERBOSITY)=0
delete finite element from mofem database
MoFEM::CoreInterface::get_comm_size
virtual int get_comm_size() const =0
MoFEM::Simple::parentAdjFunctionDim3
boost::shared_ptr< ParentFiniteElementAdjacencyFunction< 3 > > parentAdjFunctionDim3
Definition: Simple.hpp:551
MoFEM::Simple::bitAdjParentMask
BitRefLevel bitAdjParentMask
bit ref level for parent parent
Definition: Simple.hpp:513
MoFEM::Simple::skeletonFE
std::string skeletonFE
skeleton finite element
Definition: Simple.hpp:531
MoFEM::Simple::removeBoundaryField
MoFEMErrorCode removeBoundaryField(const std::string &name)
Remove field form boundary.
Definition: Simple.cpp:352
MOFEM_OPERATION_UNSUCCESSFUL
@ MOFEM_OPERATION_UNSUCCESSFUL
Definition: definitions.h:34
MoFEM::Simple::bitAdjEnt
BitRefLevel bitAdjEnt
bit ref level for parent
Definition: Simple.hpp:514
MoFEM::CoreInterface::add_ents_to_finite_element_by_dim
virtual MoFEMErrorCode add_ents_to_finite_element_by_dim(const EntityHandle entities, const int dim, const std::string &name, const bool recursive=true)=0
add entities to finite element
t
constexpr double t
plate stiffness
Definition: plate.cpp:59
MoFEM::DMMoFEMCreateMoFEM
PetscErrorCode DMMoFEMCreateMoFEM(DM dm, MoFEM::Interface *m_field_ptr, const char problem_name[], const MoFEM::BitRefLevel bit_level, const MoFEM::BitRefLevel bit_mask=MoFEM::BitRefLevel().set())
Must be called by user to set MoFEM data structures.
Definition: DMMoFEM.cpp:118
MoFEM::Simple::addDataField
MoFEMErrorCode addDataField(const std::string &name, const FieldSpace space, const FieldApproximationBase base, const FieldCoefficientsNumber nb_of_coefficients, const TagType tag_type=MB_TAG_SPARSE, const enum MoFEMTypes bh=MF_ZERO, int verb=-1)
Add field on domain.
Definition: Simple.cpp:320
field_name
constexpr auto field_name
Definition: poisson_2d_homogeneous.cpp:13
AINSWORTH_BERNSTEIN_BEZIER_BASE
@ AINSWORTH_BERNSTEIN_BEZIER_BASE
Definition: definitions.h:64
mesh_file_name
char mesh_file_name[255]
Definition: mesh_smoothing.cpp:23
MoFEM::Simple::setFieldOrder
MoFEMErrorCode setFieldOrder(const std::string field_name, const int order, const Range *ents=NULL)
Set field order.
Definition: Simple.cpp:473
MoFEM::Simple::buildFiniteElements
MoFEMErrorCode buildFiniteElements()
Build finite elements.
Definition: Simple.cpp:571
MoFEM::UnknownInterface
base class for all interface classes
Definition: UnknownInterface.hpp:34
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
MF_ZERO
@ MF_ZERO
Definition: definitions.h:98
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::Simple::domainFields
std::vector< std::string > domainFields
domain fields
Definition: Simple.hpp:517
MoFEM::CommInterface
Managing BitRefLevels.
Definition: CommInterface.hpp:21
HenckyOps::f
auto f
Definition: HenckyOps.hpp:15
MoFEM::Simple::skeletonFields
std::vector< std::string > skeletonFields
fields on the skeleton
Definition: Simple.hpp:519
MoFEM::Exceptions::ierr
static MoFEMErrorCodeGeneric< PetscErrorCode > ierr
Definition: Exceptions.hpp:76
MoFEM::DMMoFEMGetProblemPtr
PetscErrorCode DMMoFEMGetProblemPtr(DM dm, const MoFEM::Problem **problem_ptr)
Get pointer to problem data structure.
Definition: DMMoFEM.cpp:430
MoFEM::Simple::bitLevel
BitRefLevel bitLevel
BitRefLevel of the problem.
Definition: Simple.hpp:494
MoFEM::Simple::defineFiniteElements
MoFEMErrorCode defineFiniteElements()
Define finite elements.
Definition: Simple.cpp:380
MoFEM::Simple::getDim
int getDim() const
Get the problem dimension.
Definition: Simple.hpp:285
MoFEM::Simple::buildFields
MoFEMErrorCode buildFields()
Build fields.
Definition: Simple.cpp:482
MoFEM::CoreInterface::build_fields
virtual MoFEMErrorCode build_fields(int verb=DEFAULT_VERBOSITY)=0
MoFEM::Simple::bitAdjParent
BitRefLevel bitAdjParent
bit ref level for parent
Definition: Simple.hpp:512
MoFEM::Simple::noFieldDataFields
std::vector< std::string > noFieldDataFields
NOFIELD field name.
Definition: Simple.hpp:522
MoFEM::CoreInterface::modify_finite_element_add_field_data
virtual MoFEMErrorCode modify_finite_element_add_field_data(const std::string &fe_name, const std::string name_filed)=0
set finite element field data
FieldApproximationBase
FieldApproximationBase
approximation base
Definition: definitions.h:58
MoFEMTypes
MoFEMTypes
Those types control how functions respond on arguments, f.e. error handling.
Definition: definitions.h:97
MoFEM::Simple::addBoundaryField
MoFEMErrorCode addBoundaryField(const std::string &name, const FieldSpace space, const FieldApproximationBase base, const FieldCoefficientsNumber nb_of_coefficients, const TagType tag_type=MB_TAG_SPARSE, const enum MoFEMTypes bh=MF_ZERO, int verb=-1)
Add field on boundary.
Definition: Simple.cpp:282
MoFEM::Simple::reSetUp
MoFEMErrorCode reSetUp(bool only_dm=false)
Rebuild internal MoFEM data structures.
Definition: Simple.cpp:633
FiniteElement_multiIndex
multi_index_container< boost::shared_ptr< FiniteElement >, indexed_by< hashed_unique< tag< FiniteElement_Meshset_mi_tag >, member< FiniteElement, EntityHandle, &FiniteElement::meshset > >, hashed_unique< tag< BitFEId_mi_tag >, const_mem_fun< FiniteElement, BitFEId, &FiniteElement::getId >, HashBit< BitFEId >, EqBit< BitFEId > >, ordered_unique< tag< FiniteElement_name_mi_tag >, const_mem_fun< FiniteElement, boost::string_ref, &FiniteElement::getNameRef > > > > FiniteElement_multiIndex
MultiIndex for entities for FiniteElement.
Definition: FEMultiIndices.hpp:849
MoFEM::BitRefManager
Managing BitRefLevels.
Definition: BitRefManager.hpp:21
MoFEM::Simple::boundaryMeshset
EntityHandle boundaryMeshset
meshset with boundary
Definition: Simple.hpp:505
MoFEM::Simple::parentAdjSkeletonFunctionDim1
boost::shared_ptr< ParentFiniteElementAdjacencyFunctionSkeleton< 1 > > parentAdjSkeletonFunctionDim1
Definition: Simple.hpp:560
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
MoFEM::CoreInterface::build_adjacencies
virtual MoFEMErrorCode build_adjacencies(const Range &ents, int verb=DEFAULT_VERBOSITY)=0
build adjacencies
MoFEMFunctionBeginHot
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:440
MoFEM::Simple::defineProblem
MoFEMErrorCode defineProblem(const PetscBool is_partitioned=PETSC_TRUE)
define problem
Definition: Simple.cpp:452
MoFEM::Problem
keeps basic data about problem
Definition: ProblemsMultiIndices.hpp:54
MoFEM::CoreInterface::check_finite_element
virtual bool check_finite_element(const std::string &name) const =0
Check if finite element is in database.
MoFEM::CoreInterface::set_field_order
virtual MoFEMErrorCode set_field_order(const EntityHandle meshset, const EntityType type, const std::string &name, const ApproximationOrder order, int verb=DEFAULT_VERBOSITY)=0
Set order approximation of the entities in the field.
MoFEM::Simple::dM
SmartPetscObj< DM > dM
Discrete manager (interface to PETSc/MoFEM functions)
Definition: Simple.hpp:539
MoFEM::Simple::setParentAdjacency
MoFEMErrorCode setParentAdjacency()
Definition: Simple.cpp:67
MoFEM::Simple::setUp
MoFEMErrorCode setUp(const PetscBool is_partitioned=PETSC_TRUE)
Setup problem.
Definition: Simple.cpp:611
MoFEM::Problem::getName
auto getName() const
Definition: ProblemsMultiIndices.hpp:372
MoFEM::Simple::getProblemName
const std::string getProblemName() const
Get the Problem Name.
Definition: Simple.hpp:362
MoFEM::Simple::addFieldToEmptyFieldBlocks
MoFEMErrorCode addFieldToEmptyFieldBlocks(const std::string row_field, const std::string col_field) const
add empty block to problem
Definition: Simple.cpp:702
MoFEM::Simple::cOre
MoFEM::Core & cOre
Definition: Simple.hpp:492
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:416
MoFEM::Simple::skeletonMeshset
EntityHandle skeletonMeshset
skeleton meshset with boundary
Definition: Simple.hpp:506
MoFEM::Simple::MOFEM_EVENT_SimpleSetUP
PetscLogEvent MOFEM_EVENT_SimpleSetUP
Definition: Simple.hpp:497
MOFEM_NOT_IMPLEMENTED
@ MOFEM_NOT_IMPLEMENTED
Definition: definitions.h:32
MoFEM::FiniteElement_col_change_bit_reset
Reset field from column.
Definition: FEMultiIndices.hpp:944
CHKERRG
#define CHKERRG(n)
Check error code of MoFEM/MOAB/PETSc function.
Definition: definitions.h:483
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:346
NOFIELD
@ NOFIELD
scalar or vector of scalars describe (no true field)
Definition: definitions.h:84
MoFEM::DMMoFEMSetIsPartitioned
PetscErrorCode DMMoFEMSetIsPartitioned(DM dm, PetscBool is_partitioned)
Definition: DMMoFEM.cpp:1127
MoFEM::Simple::MOFEM_EVENT_SimpleBuildProblem
PetscLogEvent MOFEM_EVENT_SimpleBuildProblem
Definition: Simple.hpp:501
MoFEM::Simple::removeSkeletonField
MoFEMErrorCode removeSkeletonField(const std::string &name)
Remove field form skeleton.
Definition: Simple.cpp:366