v0.15.0
Loading...
Searching...
No Matches
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
6namespace MoFEM {
7
8MoFEMErrorCode Simple::query_interface(boost::typeindex::type_index type_index,
9 UnknownInterface **iface) const {
10 *iface = const_cast<Simple *>(this);
11 return 0;
12}
13
14template <int DIM>
16 static_assert(DIM == 2 || DIM == 3, "not implemented");
18}
19
20template <>
21MoFEMErrorCode 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
34template <>
35MoFEMErrorCode 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
51template <>
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");
70}
71
72template <> 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
103template <> 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
143MoFEMErrorCode 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
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("SimpSetUp", 0, &MOFEM_EVENT_SimpleSetUP);
171 PetscLogEventRegister("SimpLoadMesh", 0, &MOFEM_EVENT_SimpleLoadMesh);
172 PetscLogEventRegister("SimpBuildFields", 0, &MOFEM_EVENT_SimpleBuildFields);
173 PetscLogEventRegister("SimpBuildFEs", 0,
175 PetscLogEventRegister("SimpSetUp", 0, &MOFEM_EVENT_SimpleBuildProblem);
176 PetscLogEventRegister("SimpKSPSolve", 0, &MOFEM_EVENT_SimpleKSPSolve);
177 strcpy(meshFileName, "mesh.h5m");
178}
179
181 PetscBool flg = PETSC_TRUE;
183 PetscOptionsBegin(PETSC_COMM_WORLD, "", "Simple interface options",
184 "none");
185 ierr = PetscOptionsString("-file_name", "file name", "", "mesh.h5m",
186 meshFileName, 255, &flg);
187 PetscOptionsEnd();
189}
190
191MoFEMErrorCode Simple::loadFile(const std::string options,
192 const std::string mesh_file_name,
193 LoadFileFunc loadFunc) {
194 Interface &m_field = cOre;
196 PetscLogEventBegin(MOFEM_EVENT_SimpleLoadMesh, 0, 0, 0, 0);
197
198 if (!mesh_file_name.empty())
199 strcpy(meshFileName, mesh_file_name.c_str());
200
201 // Invoke the boost function, passing in required arguments
202 CHKERR loadFunc(m_field, meshFileName, options.c_str());
203 CHKERR m_field.rebuild_database();
204
205 // determine problem dimension
206 if (dIm == -1) {
207 int nb_ents_3d;
208 CHKERR m_field.get_moab().get_number_entities_by_dimension(
209 meshSet, 3, nb_ents_3d, true);
210 if (nb_ents_3d > 0) {
211 dIm = 3;
212 } else {
213 int nb_ents_2d;
214 CHKERR m_field.get_moab().get_number_entities_by_dimension(
215 meshSet, 2, nb_ents_2d, true);
216 if (nb_ents_2d > 0) {
217 dIm = 2;
218 } else {
219 dIm = 1;
220 }
221 }
222 }
223
224 if (!boundaryMeshset)
226 if (!skeletonMeshset)
228 if (addSkeletonFE)
230
231 if (bitLevel.any()) {
232 Range ents;
233 CHKERR m_field.get_moab().get_entities_by_dimension(meshSet, dIm, ents,
234 true);
235 CHKERR m_field.getInterface<BitRefManager>()->setBitRefLevel(ents, bitLevel,
236 false);
237 } else {
238 MOFEM_LOG("WORLD", Sev::warning) << "BitRefLevel is none and not set";
239 CHKERR m_field.getInterface<BitRefManager>()->addToDatabaseBitRefLevelByDim(
240 dIm, BitRefLevel().set(), BitRefLevel().set());
241 }
242
243 PetscLogEventEnd(MOFEM_EVENT_SimpleLoadMesh, 0, 0, 0, 0);
245}
246
249 Interface &m_field = cOre;
250 if (m_field.get_comm_size() == 1)
252 else
253 CHKERR loadFile("PARALLEL=READ_PART;"
254 "PARALLEL_RESOLVE_SHARED_ENTS;"
255 "PARTITION=PARALLEL_PARTITION;",
258}
259
261Simple::addDomainField(const std::string name, const FieldSpace space,
262 const FieldApproximationBase base,
263 const FieldCoefficientsNumber nb_of_coefficients,
264 const TagType tag_type, const enum MoFEMTypes bh,
265 int verb) {
266
267 Interface &m_field = cOre;
269 CHKERR m_field.add_field(name, space, base, nb_of_coefficients, tag_type, bh,
270 verb);
271 domainFields.push_back(name);
272 if (space == NOFIELD)
273 SETERRQ(
274 PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
275 "NOFIELD space for boundary field not implemented in Simple interface");
276
277
279}
280
282Simple::addDomainBrokenField(const std::string name, const FieldSpace space,
283 const FieldApproximationBase base,
284 const FieldCoefficientsNumber nb_of_coefficients,
285 const TagType tag_type, const enum MoFEMTypes bh,
286 int verb) {
287
288 Interface &m_field = cOre;
290
291 auto get_side_map_hcurl = [&]() {
292 return std::vector<
293
294 std::pair<EntityType,
296
297 >>{
298
299 {MBTRI,
300 [&](BaseFunction::DofsSideMap &dofs_side_map) -> MoFEMErrorCode {
302 dofs_side_map);
303 }}
304
305 };
306 };
307
308 auto get_side_map_hdiv = [&]() {
309 return std::vector<
310
311 std::pair<EntityType,
313
314 >>{
315
316 {MBTET,
317 [&](BaseFunction::DofsSideMap &dofs_side_map) -> MoFEMErrorCode {
319 dofs_side_map);
320 }}
321
322 };
323 };
324
325 auto get_side_map = [&](auto space) {
326 switch (space) {
327 case HCURL:
328 return get_side_map_hcurl();
329 case HDIV:
330 return get_side_map_hdiv();
331 default:
332 MOFEM_LOG_CHANNEL("WORLD");
333 MOFEM_LOG("WORLD", Sev::warning)
334 << "Side dof map for broken space <" << space << "> not implemented";
335 }
336 return std::vector<
337
338 std::pair<EntityType,
340
341 >>{};
342 };
343
344 CHKERR m_field.add_broken_field(name, space, base, nb_of_coefficients,
345 get_side_map(space), tag_type, bh, verb);
346 domainFields.push_back(name);
347 if (space == NOFIELD)
348 SETERRQ(
349 PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
350 "NOFIELD space for boundary field not implemented in Simple interface");
352}
353
355Simple::addBoundaryField(const std::string name, const FieldSpace space,
356 const FieldApproximationBase base,
357 const FieldCoefficientsNumber nb_of_coefficients,
358 const TagType tag_type, const enum MoFEMTypes bh,
359 int verb) {
360 Interface &m_field = cOre;
362 CHKERR m_field.add_field(name, space, base, nb_of_coefficients, tag_type, bh,
363 verb);
364 boundaryFields.push_back(name);
365 if (space == NOFIELD)
366 SETERRQ(
367 PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
368 "NOFIELD space for boundary field not implemented in Simple interface");
370}
371
373Simple::addSkeletonField(const std::string name, const FieldSpace space,
374 const FieldApproximationBase base,
375 const FieldCoefficientsNumber nb_of_coefficients,
376 const TagType tag_type, const enum MoFEMTypes bh,
377 int verb) {
378
379 Interface &m_field = cOre;
381 CHKERR m_field.add_field(name, space, base, nb_of_coefficients, tag_type, bh,
382 verb);
383 skeletonFields.push_back(name);
384 if (space == NOFIELD)
385 SETERRQ(
386 PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
387 "NOFIELD space for boundary field not implemented in Simple interface");
388
390}
391
393Simple::addDataField(const std::string name, const FieldSpace space,
394 const FieldApproximationBase base,
395 const FieldCoefficientsNumber nb_of_coefficients,
396 const TagType tag_type, const enum MoFEMTypes bh,
397 int verb) {
398
399 Interface &m_field = cOre;
401 CHKERR m_field.add_field(name, space, base, nb_of_coefficients, tag_type, bh,
402 verb);
403 if (space == NOFIELD)
404 noFieldDataFields.push_back(name);
405 else
406 dataFields.push_back(name);
408}
409
411Simple::addMeshsetField(const std::string name, const FieldSpace space,
412 const FieldApproximationBase base,
413 const FieldCoefficientsNumber nb_of_coefficients,
414 const TagType tag_type, const enum MoFEMTypes bh,
415 int verb) {
416
417 Interface &m_field = cOre;
419 CHKERR m_field.add_field(name, space, base, nb_of_coefficients, tag_type, bh,
420 verb);
421 if (space != NOFIELD) {
422 meshsetFields.push_back(name);
423 } else {
424 noFieldFields.push_back(name);
425 }
427}
428
431
432 auto remove_field_from_list = [&](auto &vec) {
433 auto it = std::find(vec.begin(), vec.end(), name);
434 if (it != vec.end())
435 vec.erase(it);
436 };
437
438 remove_field_from_list(domainFields);
439
441}
442
445
446 auto remove_field_from_list = [&](auto &vec) {
447 auto it = std::find(vec.begin(), vec.end(), name);
448 if (it != vec.end())
449 vec.erase(it);
450 };
451
452 remove_field_from_list(boundaryFields);
453
455}
456
459
460 auto remove_field_from_list = [&](auto &vec) {
461 auto it = std::find(vec.begin(), vec.end(), name);
462 if (it != vec.end())
463 vec.erase(it);
464 };
465
466 remove_field_from_list(skeletonFields);
467
469}
470
472 Interface &m_field = cOre;
474
475 auto clear_rows_and_cols = [&](auto &fe_name) {
477 auto fe_ptr = m_field.get_finite_elements();
478 auto &fe_by_name = const_cast<FiniteElement_multiIndex *>(fe_ptr)
479 ->get<FiniteElement_name_mi_tag>();
480 auto it_fe = fe_by_name.find(fe_name);
481 if (it_fe != fe_by_name.end()) {
482
483 if (!fe_by_name.modify(it_fe, FiniteElement_row_change_bit_reset()))
484 SETERRQ(m_field.get_comm(), MOFEM_OPERATION_UNSUCCESSFUL,
485 "modification unsuccessful");
486
487 if (!fe_by_name.modify(it_fe, FiniteElement_col_change_bit_reset()))
488 SETERRQ(m_field.get_comm(), MOFEM_OPERATION_UNSUCCESSFUL,
489 "modification unsuccessful");
490 }
492 };
493 CHKERR clear_rows_and_cols(domainFE);
494 CHKERR clear_rows_and_cols(boundaryFE);
495 CHKERR clear_rows_and_cols(skeletonFE);
496
497 // Define finite elements
499
500 auto add_fields = [&](auto &fe_name, auto &fields) {
502 for (auto &field : fields) {
503 CHKERR m_field.modify_finite_element_add_field_row(fe_name, field);
504 CHKERR m_field.modify_finite_element_add_field_col(fe_name, field);
505 CHKERR m_field.modify_finite_element_add_field_data(fe_name, field);
506 }
508 };
509
510 auto add_data_fields = [&](auto &fe_name, auto &fields) {
512 for (auto &field : fields)
513 CHKERR m_field.modify_finite_element_add_field_data(fe_name, field);
515 };
516
517 CHKERR add_fields(domainFE, domainFields);
518 CHKERR add_data_fields(domainFE, dataFields);
519 CHKERR add_fields(domainFE, noFieldFields);
520 CHKERR add_data_fields(domainFE, noFieldDataFields);
521
522 if (addBoundaryFE || !boundaryFields.empty()) {
524 CHKERR add_fields(boundaryFE, domainFields);
525 if (!boundaryFields.empty())
526 CHKERR add_fields(boundaryFE, boundaryFields);
527 CHKERR add_data_fields(boundaryFE, dataFields);
528 CHKERR add_data_fields(boundaryFE, noFieldDataFields);
529 CHKERR add_fields(boundaryFE, noFieldFields);
530 }
531 if (addSkeletonFE || !skeletonFields.empty()) {
533 CHKERR add_fields(skeletonFE, domainFields);
534 if (!skeletonFields.empty())
535 CHKERR add_fields(skeletonFE, skeletonFields);
536 CHKERR add_data_fields(skeletonFE, dataFields);
537 CHKERR add_data_fields(skeletonFE, noFieldDataFields);
538 CHKERR add_fields(skeletonFE, noFieldFields);
539 }
540
541 if (noFieldFields.size() || meshsetFields.size()) {
543
544 auto create_meshset = [&]() {
545 EntityHandle fe_meshset;
546 auto create = [&]() {
548 CHKERR m_field.get_moab().create_meshset(MESHSET_SET, fe_meshset);
549 {
550 EntityHandle meshset;
551 CHKERR m_field.get_moab().create_meshset(MESHSET_SET, meshset);
552 CHKERR m_field.getInterface<BitRefManager>()->setBitLevelToMeshset(
553 meshset, BitRefLevel().set());
554 ParallelComm *pcomm =
555 ParallelComm::get_pcomm(&m_field.get_moab(), MYPCOMM_INDEX);
556 Tag part_tag = pcomm->part_tag();
557 int rank = pcomm->rank();
558 CHKERR m_field.get_moab().tag_set_data(part_tag, &meshset, 1, &rank);
559 }
561 };
562 CHK_THROW_MESSAGE(create(), "Failed to create meshset");
563 return fe_meshset;
564 };
565
566 CHKERR m_field.add_ents_to_finite_element_by_MESHSET(create_meshset(),
567 meshsetFE, false);
568
569 CHKERR add_fields(meshsetFE, meshsetFields);
570 CHKERR add_fields(meshsetFE, noFieldFields);
571 CHKERR add_data_fields(meshsetFE, noFieldDataFields);
572 }
573
575}
576
577MoFEMErrorCode Simple::defineProblem(const PetscBool is_partitioned) {
578 Interface &m_field = cOre;
580 // Create dm instance
581 dM = createDM(m_field.get_comm(), "DMMOFEM");
582 // set dm data structure which created mofem data structures
585 CHKERR DMSetFromOptions(dM);
587 if (addBoundaryFE || !boundaryFields.empty()) {
589 }
590 if (addSkeletonFE || !skeletonFields.empty()) {
592 }
593 if (noFieldFields.size() || meshsetFields.size()) {
595 }
597 CHKERR DMMoFEMSetIsPartitioned(dM, is_partitioned);
599}
600
602 const int order, const Range *ents) {
604 fieldsOrder.emplace_back(field_name, order,
605 ents == NULL ? Range() : Range(*ents),
606 ents == NULL ? false : true);
608}
609
611 Interface &m_field = cOre;
613 PetscLogEventBegin(MOFEM_EVENT_SimpleBuildFields, 0, 0, 0, 0);
614
615 auto comm_interface_ptr = m_field.getInterface<CommInterface>();
616 auto bit_ref_ptr = m_field.getInterface<BitRefManager>();
617
618 // Add entities to the fields
619 auto add_ents_to_field = [&](auto meshset, auto dim, auto &fields) {
621 for (auto &field : fields) {
622 CHKERR m_field.add_ents_to_field_by_dim(meshset, dim, field);
623 CHKERR comm_interface_ptr->synchroniseFieldEntities(field, 0);
624 }
626 };
627
628 auto make_no_field_ents = [&](auto &fields) {
630 for (auto &field : fields) {
631 std::array<double, 6> coords = {0, 0, 0, 0, 0, 0};
632 CHKERR m_field.create_vertices_and_add_to_field(field, coords.data(), 1);
633 CHKERR comm_interface_ptr->makeFieldEntitiesMultishared(field, 0);
634 CHKERR bit_ref_ptr->setFieldEntitiesBitRefLevel(field, bitLevel);
635 }
637 };
638
639 CHKERR add_ents_to_field(meshSet, dIm, domainFields);
640 CHKERR add_ents_to_field(meshSet, dIm, dataFields);
641 CHKERR add_ents_to_field(boundaryMeshset, dIm - 1, boundaryFields);
642 CHKERR add_ents_to_field(skeletonMeshset, dIm - 1, skeletonFields);
643
644 std::set<std::string> nofield_fields;
645 for (auto &f : noFieldFields)
646 nofield_fields.insert(f);
647 for (auto &f : noFieldDataFields)
648 nofield_fields.insert(f);
649
650 CHKERR make_no_field_ents(nofield_fields);
651
652 // Set order
653
654 auto get_field_ptr = [&](auto &f) { return m_field.get_field_structure(f); };
655 for (auto &t : fieldsOrder) {
656 const auto f = std::get<0>(t);
657 const auto order = std::get<1>(t);
658
659 MOFEM_TAG_AND_LOG("WORLD", Sev::inform, "Simple")
660 << "Set order to field " << f << " order " << order;
661 MOFEM_LOG_CHANNEL("SYNC");
662 if (std::get<3>(t)) {
663 MOFEM_TAG_AND_LOG("SYNC", Sev::verbose, "Simple")
664 << "To ents: " << std::endl
665 << std::get<2>(t) << std::endl;
666 }
667 MOFEM_LOG_SEVERITY_SYNC(m_field.get_comm(), Sev::verbose);
668
669 if (std::get<3>(t)) {
670
671 CHKERR m_field.set_field_order(std::get<2>(t), f, order);
672
673 } else {
674 auto f_ptr = get_field_ptr(f);
675
676 if (f_ptr->getSpace() == H1) {
677 if (f_ptr->getApproxBase() == AINSWORTH_BERNSTEIN_BEZIER_BASE) {
678 CHKERR m_field.set_field_order(meshSet, MBVERTEX, f, order);
679 } else {
680 CHKERR m_field.set_field_order(meshSet, MBVERTEX, f, 1);
681 }
682 }
683
684 for (auto d = 1; d <= dIm; ++d) {
685 for (EntityType t = CN::TypeDimensionMap[d].first;
686 t <= CN::TypeDimensionMap[d].second; ++t) {
687 CHKERR m_field.set_field_order(meshSet, t, f, order);
688 }
689 }
690 }
691 }
692 MOFEM_LOG_CHANNEL("WORLD");
693 // Build fields
694 CHKERR m_field.build_fields();
695 PetscLogEventEnd(MOFEM_EVENT_SimpleBuildFields, 0, 0, 0, 0);
697}
698
700 Interface &m_field = cOre;
702 PetscLogEventBegin(MOFEM_EVENT_SimpleBuildFiniteElements, 0, 0, 0, 0);
703 // Add finite elements
705 true);
707 if (addBoundaryFE || boundaryFields.size()) {
709 boundaryFE, true);
711 }
712 if (addSkeletonFE || skeletonFields.size()) {
714 skeletonFE, true);
716 }
717
718 if (noFieldFields.size() || meshsetFields.size()) {
719 auto add_fields_ents = [&](auto list) {
721 auto fe_meshset = m_field.get_finite_element_meshset(meshsetFE);
722 for (auto &fe_ents : meshsetFiniteElementEntities) {
723 EntityHandle meshset_entity_fe;
724 CHKERR m_field.get_moab().create_meshset(MESHSET_SET,
725 meshset_entity_fe);
726 CHKERR m_field.getInterface<BitRefManager>()->setBitLevelToMeshset(
727 meshset_entity_fe, BitRefLevel().set());
728 for (auto f : list) {
729 auto field_meshset = m_field.get_field_meshset(f);
730 Range ents;
731 CHKERR m_field.get_moab().get_entities_by_handle(field_meshset, ents,
732 true);
733 CHKERR m_field.get_moab().add_entities(meshset_entity_fe, ents);
734 }
735 CHKERR m_field.get_moab().add_entities(meshset_entity_fe, fe_ents);
736 CHKERR m_field.get_moab().add_entities(fe_meshset, &meshset_entity_fe,
737 1);
738 }
740 };
741 CHKERR add_fields_ents(noFieldFields);
743 }
744
745 for (std::vector<std::string>::iterator fit = otherFEs.begin();
746 fit != otherFEs.end(); ++fit) {
747 CHKERR m_field.build_finite_elements(*fit);
748 }
749 PetscLogEventEnd(MOFEM_EVENT_SimpleBuildFiniteElements, 0, 0, 0, 0);
751}
752
754 Interface &m_field = cOre;
756 PetscLogEventBegin(MOFEM_EVENT_SimpleBuildProblem, 0, 0, 0, 0);
758 // Set problem by the DOFs on the fields rather that by adding DOFs on the
759 // elements
760 m_field.getInterface<ProblemsManager>()->buildProblemFromFields = PETSC_TRUE;
762 m_field.getInterface<ProblemsManager>()->buildProblemFromFields = PETSC_FALSE;
763 PetscLogEventEnd(MOFEM_EVENT_SimpleBuildProblem, 0, 0, 0, 0);
765}
766
767MoFEMErrorCode Simple::setUp(const PetscBool is_partitioned) {
769
770 PetscLogEventBegin(MOFEM_EVENT_SimpleSetUP, 0, 0, 0, 0);
771
773
774 if (addSkeletonFE || !skeletonFields.empty()) {
776 if (addBoundaryFE || !boundaryFields.empty()) {
778 }
779 }
780
783
784 CHKERR defineProblem(is_partitioned);
788
789 PetscLogEventEnd(MOFEM_EVENT_SimpleSetUP, 0, 0, 0, 0);
791}
792
794 Interface &m_field = cOre;
796 PetscLogEventBegin(MOFEM_EVENT_SimpleBuildProblem, 0, 0, 0, 0);
797
798 if (!only_dm) {
801 if (addSkeletonFE || !skeletonFields.empty()) {
803 if (addBoundaryFE || !boundaryFields.empty()) {
805 }
806 }
810 }
811
813
814 const Problem *problem_ptr;
815 CHKERR DMMoFEMGetProblemPtr(dM, &problem_ptr);
816 const auto problem_name = problem_ptr->getName();
817 CHKERR m_field.modify_problem_ref_level_set_bit(problem_name, bitLevel);
820
821 // Set problem by the DOFs on the fields rather that by adding DOFs on the
822 // elements
823 m_field.getInterface<ProblemsManager>()->buildProblemFromFields = PETSC_TRUE;
825 m_field.getInterface<ProblemsManager>()->buildProblemFromFields = PETSC_FALSE;
826
827 PetscLogEventEnd(MOFEM_EVENT_SimpleBuildProblem, 0, 0, 0, 0);
829}
830
833 CHKERR PetscObjectReference(getPetscObject(dM.get()));
834 *dm = dM.get();
836}
837
838/**
839 * @brief Delete dm
840 *
841 * @return MoFEMErrorCode
842 */
848
849/**
850 * @brief Delete finite elements
851 *
852 * @return MoFEMErrorCode
853 */
855 Interface &m_field = cOre;
857 for (auto fe : {domainFE, boundaryFE, skeletonFE}) {
858 if (m_field.check_finite_element(fe)) {
859 CHKERR m_field.delete_finite_element(fe);
860 }
861 }
863}
864
866Simple::addFieldToEmptyFieldBlocks(const std::string row_field,
867 const std::string col_field) const {
868 Interface &m_field = cOre;
871 getProblemName(), row_field, col_field);
873}
874
876 Interface &m_field = cOre;
878 ParallelComm *pcomm =
879 ParallelComm::get_pcomm(&m_field.get_moab(), MYPCOMM_INDEX);
880
881 auto get_skin = [&](auto meshset) {
882 // filter not owned entities, those are not on boundary
883
884 Range domain_ents;
885 CHKERR m_field.get_moab().get_entities_by_dimension(meshset, dIm,
886 domain_ents, true);
887 CHKERR pcomm->filter_pstatus(domain_ents,
888 PSTATUS_SHARED | PSTATUS_MULTISHARED,
889 PSTATUS_NOT, -1, nullptr);
890
891 Skinner skin(&m_field.get_moab());
892 Range domain_skin;
893 CHKERR skin.find_skin(0, domain_ents, false, domain_skin);
894 CHKERR pcomm->filter_pstatus(domain_skin,
895 PSTATUS_SHARED | PSTATUS_MULTISHARED,
896 PSTATUS_NOT, -1, nullptr);
897 return domain_skin;
898 };
899
900 auto create_boundary_meshset = [&](auto &&domain_skin) {
902 // create boundary meshset
903 if (boundaryMeshset != 0) {
905 }
906 CHKERR m_field.get_moab().create_meshset(MESHSET_SET, boundaryMeshset);
907 CHKERR m_field.get_moab().add_entities(boundaryMeshset, domain_skin);
908 for (int dd = 0; dd != dIm - 1; dd++) {
909 Range adj;
910 CHKERR m_field.get_moab().get_adjacencies(domain_skin, dd, false, adj,
911 moab::Interface::UNION);
912 CHKERR m_field.get_moab().add_entities(boundaryMeshset, adj);
913 }
915 };
916
917 CHKERR create_boundary_meshset(get_skin(meshSet));
918
920}
921
923 Interface &m_field = cOre;
925
926 ParallelComm *pcomm =
927 ParallelComm::get_pcomm(&m_field.get_moab(), MYPCOMM_INDEX);
928
929 auto create_skeleton_meshset = [&](auto meshset) {
931 // create boundary meshset
932 if (skeletonMeshset != 0) {
934 }
935 Range boundary_ents, skeleton_ents;
936 CHKERR m_field.get_moab().get_entities_by_dimension(boundaryMeshset,
937 dIm - 1, boundary_ents);
938 Range domain_ents;
939 CHKERR m_field.get_moab().get_entities_by_dimension(meshset, dIm,
940 domain_ents, true);
941 CHKERR m_field.get_moab().get_adjacencies(
942 domain_ents, dIm - 1, true, skeleton_ents, moab::Interface::UNION);
943 skeleton_ents = subtract(skeleton_ents, boundary_ents);
944 CHKERR pcomm->filter_pstatus(skeleton_ents, PSTATUS_NOT_OWNED, PSTATUS_NOT,
945 -1, nullptr);
946 CHKERR m_field.get_moab().create_meshset(MESHSET_SET, skeletonMeshset);
947 CHKERR m_field.get_moab().add_entities(skeletonMeshset, skeleton_ents);
949 };
950
951 CHKERR create_skeleton_meshset(meshSet);
952
954}
955
957 Interface &m_field = cOre;
959 MOFEM_LOG("WORLD", Sev::verbose) << "Exchange ghost cells";
960
961 ParallelComm *pcomm =
962 ParallelComm::get_pcomm(&m_field.get_moab(), MYPCOMM_INDEX);
963 if (pcomm == NULL)
964 pcomm = new ParallelComm(&m_field.get_moab(), m_field.get_comm());
965
966 CHKERR pcomm->exchange_ghost_cells(getDim(), getDim() - 1, 1,
967 3 /**get all adjacent ghosted entities */,
968 true, true, meshSet ? &meshSet : nullptr);
969
970 Range shared;
971 CHKERR m_field.get_moab().get_entities_by_dimension(0, dIm, shared);
972 for (auto d = dIm - 1; d >= 0; --d) {
973 CHKERR m_field.get_moab().get_adjacencies(shared, d, false, shared,
974 moab::Interface::UNION);
975 }
976 CHKERR pcomm->filter_pstatus(shared, PSTATUS_SHARED | PSTATUS_MULTISHARED,
977 PSTATUS_OR, -1, &shared);
978 Tag part_tag = pcomm->part_tag();
979 CHKERR pcomm->exchange_tags(part_tag, shared);
980 CHKERR m_field.getInterface<CommInterface>()->resolveParentEntities(shared,
981 VERBOSE);
982
984}
985
986} // namespace MoFEM
static auto get_skin(MoFEM::Interface &m_field, Range body_ents)
#define MOFEM_LOG_SEVERITY_SYNC(comm, severity)
Synchronise "SYNC" on curtain severity level.
#define MOFEM_TAG_AND_LOG(channel, severity, tag)
Tag and log in channel.
@ VERBOSE
MoFEMTypes
Those types control how functions respond on arguments, f.e. error handling.
@ MF_ZERO
FieldApproximationBase
approximation base
Definition definitions.h:58
@ AINSWORTH_BERNSTEIN_BEZIER_BASE
Definition definitions.h:64
#define CHK_THROW_MESSAGE(err, msg)
Check and throw MoFEM exception.
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
FieldSpace
approximation spaces
Definition definitions.h:82
@ NOFIELD
scalar or vector of scalars describe (no true field)
Definition definitions.h:84
@ H1
continuous field
Definition definitions.h:85
@ HCURL
field with continuous tangents
Definition definitions.h:86
@ HDIV
field with continuous normal traction
Definition definitions.h:87
#define MYPCOMM_INDEX
default communicator number PCOMM
@ DISCONTINUOUS
Broken continuity (No effect on L2 space)
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
@ MOFEM_OPERATION_UNSUCCESSFUL
Definition definitions.h:34
@ MOFEM_NOT_IMPLEMENTED
Definition definitions.h:32
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define CHKERR
Inline error check.
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
#define THROW_MESSAGE(msg)
Throw MoFEM exception.
constexpr int order
PetscErrorCode DMMoFEMSetIsPartitioned(DM dm, PetscBool is_partitioned)
Definition DMMoFEM.cpp:1113
PetscErrorCode DMMoFEMAddElement(DM dm, std::string fe_name)
add element to dm
Definition DMMoFEM.cpp:488
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:114
PetscErrorCode DMMoFEMGetProblemPtr(DM dm, const MoFEM::Problem **problem_ptr)
Get pointer to problem data structure.
Definition DMMoFEM.cpp:422
PetscErrorCode DMSetUp_MoFEM(DM dm)
Definition DMMoFEM.cpp:1297
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.
virtual const FiniteElement_multiIndex * get_finite_elements() const =0
Get the finite elements object.
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
virtual MoFEMErrorCode add_finite_element(const std::string &fe_name, enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
add finite element
virtual EntityHandle get_finite_element_meshset(const std::string name) const =0
virtual MoFEMErrorCode build_finite_elements(int verb=DEFAULT_VERBOSITY)=0
Build finite elements.
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
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
virtual MoFEMErrorCode add_ents_to_finite_element_by_MESHSET(const EntityHandle meshset, const std::string &name, const bool recursive=false)=0
add MESHSET element to finite element database given by name
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
virtual MoFEMErrorCode modify_finite_element_add_field_data(const std::string &fe_name, const std::string name_field)=0
set finite element field data
virtual const Field * get_field_structure(const std::string &name, enum MoFEMTypes bh=MF_EXIST) const =0
get field structure
virtual MoFEMErrorCode build_fields(int verb=DEFAULT_VERBOSITY)=0
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.
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.
#define MOFEM_LOG(channel, severity)
Log.
#define MOFEM_LOG_CHANNEL(channel)
Set and reset channel.
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
virtual MoFEMErrorCode modify_problem_ref_level_set_bit(const std::string &name_problem, const BitRefLevel &bit)=0
set ref level for problem
char mesh_file_name[255]
static MoFEMErrorCodeGeneric< PetscErrorCode > ierr
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
int FieldCoefficientsNumber
Number of field coefficients.
Definition Types.hpp:27
std::bitset< BITREFLEVEL_SIZE > BitRefLevel
Bit structure attached to each entity identifying to what mesh entity is attached.
Definition Types.hpp:40
implementation of Data Operators for Forces and Sources
Definition Common.hpp:10
PetscObject getPetscObject(T obj)
auto createDM(MPI_Comm comm, const std::string dm_type_name)
Creates smart DM object.
constexpr double t
plate stiffness
Definition plate.cpp:58
constexpr auto field_name
multi_index_container< DofsSideMapData, indexed_by< ordered_non_unique< tag< TypeSide_mi_tag >, composite_key< DofsSideMapData, member< DofsSideMapData, EntityType, &DofsSideMapData::type >, member< DofsSideMapData, int, &DofsSideMapData::side > > >, ordered_unique< tag< EntDofIdx_mi_tag >, member< DofsSideMapData, int, &DofsSideMapData::dof > > > > DofsSideMap
Map entity stype and side to element/entity dof index.
Managing BitRefLevels.
Managing BitRefLevels.
virtual int get_comm_size() const =0
virtual moab::Interface & get_moab()=0
virtual EntityHandle get_field_meshset(const std::string name) const =0
get field meshset
virtual MoFEMErrorCode add_broken_field(const std::string name, const FieldSpace space, const FieldApproximationBase base, const FieldCoefficientsNumber nb_of_coefficients, const std::vector< std::pair< EntityType, std::function< MoFEMErrorCode(BaseFunction::DofsSideMap &)> > > list_dof_side_map, const TagType tag_type=MB_TAG_SPARSE, const enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
Add field.
virtual bool check_finite_element(const std::string &name) const =0
Check if finite element is in database.
virtual MoFEMErrorCode build_adjacencies(const Range &ents, int verb=DEFAULT_VERBOSITY)=0
build adjacencies
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.
virtual MPI_Comm & get_comm() const =0
virtual MoFEMErrorCode delete_finite_element(const std::string name, int verb=DEFAULT_VERBOSITY)=0
delete finite element from mofem database
virtual MoFEMErrorCode rebuild_database(int verb=DEFAULT_VERBOSITY)=0
Clear database and initialize it once again.
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.
Core (interface) class.
Definition Core.hpp:82
Deprecated interface functions.
keeps basic data about problem
Problem manager is used to build and partition problems.
Simple interface for fast problem set-up.
Definition Simple.hpp:27
boost::shared_ptr< ParentFiniteElementAdjacencyFunction< 2 > > parentAdjFunctionDim2
Definition Simple.hpp:618
EntityHandle boundaryMeshset
meshset with boundary
Definition Simple.hpp:566
MoFEMErrorCode buildProblem()
Build problem.
Definition Simple.cpp:753
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:261
std::vector< std::string > dataFields
Data fields.
Definition Simple.hpp:581
MoFEMErrorCode setParentAdjacency()
Definition Simple.cpp:67
std::list< std::tuple< std::string, int, Range, bool > > fieldsOrder
Definition Simple.hpp:587
std::vector< std::string > domainFields
domain fields
Definition Simple.hpp:578
MoFEMErrorCode defineFiniteElements()
Define finite elements.
Definition Simple.cpp:471
MoFEM::Core & cOre
Definition Simple.hpp:553
BitRefLevel bitAdjEnt
bit ref level for parent
Definition Simple.hpp:575
char meshFileName[255]
mesh file name
Definition Simple.hpp:600
EntityHandle meshSet
domain meshset
Definition Simple.hpp:565
int getDim() const
Get the problem dimension.
Definition Simple.hpp:331
MoFEMErrorCode loadFile(const std::string options, const std::string mesh_file_name, LoadFileFunc loadFunc=defaultLoadFileFunc)
Load mesh file.
Definition Simple.cpp:191
std::vector< std::string > otherFEs
Other finite elements.
Definition Simple.hpp:596
MoFEMErrorCode buildFiniteElements()
Build finite elements.
Definition Simple.cpp:699
MoFEMErrorCode addMeshsetField(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 meshset field.
Definition Simple.cpp:411
boost::shared_ptr< ParentFiniteElementAdjacencyFunctionSkeleton< 2 > > parentAdjSkeletonFunctionDim2
Definition Simple.hpp:623
MoFEMErrorCode exchangeGhostCells()
Definition Simple.cpp:956
PetscLogEvent MOFEM_EVENT_SimpleBuildFields
Definition Simple.hpp:560
PetscLogEvent MOFEM_EVENT_SimpleLoadMesh
Definition Simple.hpp:559
MoFEMErrorCode reSetUp(bool only_dm=false)
Rebuild internal MoFEM data structures.
Definition Simple.cpp:793
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:355
std::vector< std::string > noFieldFields
NOFIELD field name.
Definition Simple.hpp:583
std::string meshsetFE
meshset finite element
Definition Simple.hpp:594
bool addSkeletonFE
Add skeleton FE.
Definition Simple.hpp:569
BitRefLevel bitAdjParent
bit ref level for parent
Definition Simple.hpp:573
bool addBoundaryFE
Add boundary FE.
Definition Simple.hpp:570
MoFEMErrorCode addFieldToEmptyFieldBlocks(const std::string row_field, const std::string col_field) const
add empty block to problem
Definition Simple.cpp:866
std::string skeletonFE
skeleton finite element
Definition Simple.hpp:593
std::vector< std::string > noFieldDataFields
NOFIELD field name.
Definition Simple.hpp:584
std::vector< std::string > meshsetFields
meshset fields
Definition Simple.hpp:582
std::string boundaryFE
boundary finite element
Definition Simple.hpp:592
std::string nameOfProblem
problem name
Definition Simple.hpp:590
boost::shared_ptr< ParentFiniteElementAdjacencyFunction< 3 > > parentAdjFunctionDim3
Definition Simple.hpp:616
MoFEMErrorCode getOptions()
get options
Definition Simple.cpp:180
std::string domainFE
domain finite element
Definition Simple.hpp:591
std::vector< std::string > boundaryFields
boundary fields
Definition Simple.hpp:579
boost::shared_ptr< ParentFiniteElementAdjacencyFunctionSkeleton< 1 > > parentAdjSkeletonFunctionDim1
Definition Simple.hpp:625
MoFEMErrorCode removeSkeletonField(const std::string name)
Remove field form skeleton.
Definition Simple.cpp:457
BitRefLevel bitAdjEntMask
bit ref level for parent parent
Definition Simple.hpp:576
SmartPetscObj< DM > dM
Discrete manager (interface to PETSc/MoFEM functions)
Definition Simple.hpp:604
MoFEMErrorCode removeDomainField(const std::string name)
Remove field form domain.
Definition Simple.cpp:429
MoFEMErrorCode buildFields()
Build fields.
Definition Simple.cpp:610
SmartPetscObj< DM > getDM()
Return smart DM object.
Definition Simple.hpp:321
std::vector< Range > meshsetFiniteElementEntities
Meshset element entities.
Definition Simple.hpp:598
BitRefLevel bitLevel
BitRefLevel of the problem.
Definition Simple.hpp:555
MoFEMErrorCode setFieldOrder(const std::string field_name, const int order, const Range *ents=NULL)
Set field order.
Definition Simple.cpp:601
MoFEMErrorCode defineProblem(const PetscBool is_partitioned=PETSC_TRUE)
define problem
Definition Simple.cpp:577
std::vector< std::string > skeletonFields
fields on the skeleton
Definition Simple.hpp:580
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:393
MoFEMErrorCode query_interface(boost::typeindex::type_index type_index, UnknownInterface **iface) const
Definition Simple.cpp:8
MoFEMErrorCode setUp(const PetscBool is_partitioned=PETSC_TRUE)
Setup problem.
Definition Simple.cpp:767
boost::function< MoFEMErrorCode(Interface &, const char *, const char *)> LoadFileFunc
Definition Simple.hpp:43
MoFEMErrorCode removeBoundaryField(const std::string name)
Remove field form boundary.
Definition Simple.cpp:443
const std::string getProblemName() const
Get the Problem Name.
Definition Simple.hpp:408
int dIm
dimension of problem
Definition Simple.hpp:601
MoFEMErrorCode setSkeletonAdjacency(int dim=-1, std::string fe_name="")
Set the skeleton adjacency object.
Definition Simple.cpp:143
MoFEMErrorCode deleteDM()
Delete dm.
Definition Simple.cpp:843
PetscLogEvent MOFEM_EVENT_SimpleSetUP
Definition Simple.hpp:558
PetscLogEvent MOFEM_EVENT_SimpleKSPSolve
Definition Simple.hpp:563
MoFEMErrorCode createBoundaryMeshset()
Definition Simple.cpp:875
PetscLogEvent MOFEM_EVENT_SimpleBuildFiniteElements
Definition Simple.hpp:561
BitRefLevel bitLevelMask
BitRefLevel of the problem.
Definition Simple.hpp:556
BitRefLevel bitAdjParentMask
bit ref level for parent parent
Definition Simple.hpp:574
PetscLogEvent MOFEM_EVENT_SimpleBuildProblem
Definition Simple.hpp:562
boost::shared_ptr< ParentFiniteElementAdjacencyFunction< 1 > > parentAdjFunctionDim1
Definition Simple.hpp:620
MoFEMErrorCode deleteFiniteElements()
Delete finite elements.
Definition Simple.cpp:854
Simple(const MoFEM::Core &core)
Definition Simple.cpp:162
EntityHandle skeletonMeshset
skeleton meshset with boundary
Definition Simple.hpp:567
MoFEMErrorCode createSkeletonMeshset()
Definition Simple.cpp:922
bool addParentAdjacencies
If set to true parent adjacencies are build.
Definition Simple.hpp:571
MoFEMErrorCode addDomainBrokenField(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 broken field on domain.
Definition Simple.cpp:282
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:373
static MoFEMErrorCode setDofsSideMap(const FieldSpace space, const FieldContinuity continuity, const FieldApproximationBase base, DofsSideMap &)
Set map of dof to side number.
static MoFEMErrorCode setDofsSideMap(const FieldSpace space, const FieldContinuity continuity, const FieldApproximationBase base, DofsSideMap &)
Set map of dof to side number.
base class for all interface classes
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.