v0.15.5
Loading...
Searching...
No Matches
Core.cpp
Go to the documentation of this file.
1/** \file Core.cpp
2 * \brief Multi-index containers, data structures and other low-level functions
3 */
4
5
6#include <MoFEM.hpp>
7
9
10extern "C" {
12}
13
14namespace MoFEM {
15
16WrapMPIComm::WrapMPIComm(MPI_Comm comm, bool petsc)
17 : comm(comm), isPetscComm(petsc) {
18 if (isPetscComm) {
19 ierr = PetscCommDuplicate(comm, &duplicatedComm, NULL);
20 CHKERRABORT(comm, ierr);
21 } else {
22 int ierr = MPI_Comm_dup(comm, &duplicatedComm);
23 if (ierr) {
24 THROW_MESSAGE("MPI_Comm_dup not working");
25 }
26 }
27}
29 if (isPetscComm) {
30 ierr = PetscCommDestroy(&duplicatedComm);
31 CHKERRABORT(comm, ierr);
32 } else {
33 int ierr = MPI_Comm_free(&duplicatedComm);
34 if (ierr) {
35 CHKERRABORT(comm, MOFEM_DATA_INCONSISTENCY);
36 }
37 }
38}
39
40constexpr const int CoreTmp<0>::value;
41constexpr const int CoreTmp<-1>::value;
42
43MoFEMErrorCode Core::query_interface(boost::typeindex::type_index type_index,
44 UnknownInterface **iface) const {
46 *iface = NULL;
47 if (type_index == boost::typeindex::type_id<CoreInterface>()) {
48 *iface = static_cast<CoreInterface *>(const_cast<Core *>(this));
50 } else if (type_index ==
51 boost::typeindex::type_id<DeprecatedCoreInterface>()) {
52 *iface = static_cast<DeprecatedCoreInterface *>(const_cast<Core *>(this));
54 }
55
56 // Get sub-interface
57 auto it = iFaces.find(type_index);
58 if (it != iFaces.end()) {
59 *iface = it->second;
61 }
62
63 *iface = NULL;
64 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "unknown interface");
66}
67
70PetscBool Core::isInitialized;
71
72MoFEMErrorCode Core::Initialize(int *argc, char ***args, const char file[],
73 const char help[]) {
74
75 MPI_Initialized(&mpiInitialised);
76 if (!mpiInitialised)
77 MPI_Init(argc, args);
78
79 PetscInitialized(&isInitialized);
80 if (isInitialized == PETSC_FALSE) {
81 for (int i = 0; i < *argc - 1; ++i) {
82 if (std::strcmp((*args)[i], "-param_file") == 0) {
83 file = (*args)[i + 1];
84 break;
85 }
86 }
87 PetscInitialize(argc, args, file, help);
88 PetscPushErrorHandler(mofem_error_handler, PETSC_NULLPTR);
89 }
90
91 LogManager::createDefaultSinks(MPI_COMM_WORLD);
92 PetscVFPrintf = LogManager::logPetscFPrintf;
94 isGloballyInitialised = true;
95
96 MOFEM_LOG_CHANNEL("WORLD");
97 char petsc_version[255];
98 CHKERR PetscGetVersion(petsc_version, 255);
99 MOFEM_LOG_C("WORLD", Sev::inform, "MoFEM version %d.%d.%d (%s %s)",
100 MoFEM_VERSION_MAJOR, MoFEM_VERSION_MINOR, MoFEM_VERSION_BUILD,
101 MOAB_VERSION_STRING, petsc_version);
102 MOFEM_LOG_C("WORLD", Sev::inform, "git commit id %s", GIT_SHA1_NAME);
103
104 auto log_time = [&](const auto perefix, auto time) {
105 MOFEM_LOG("WORLD", Sev::inform)
106 << perefix << time.date().year() << "-" << time.date().month() << "-"
107 << time.date().day() << " " << time.time_of_day().hours() << ":"
108 << time.time_of_day().minutes() << ":" << time.time_of_day().seconds();
109 };
110
111 // Get current system time
112 log_time("Local time: ", boost::posix_time::second_clock::local_time());
113 log_time("UTC time: ", boost::posix_time::second_clock::universal_time());
114
115 return MOFEM_SUCCESS;
116}
117
119 if (isGloballyInitialised) {
120 PetscPopErrorHandler();
121 isGloballyInitialised = false;
122
123 if (isInitialized == PETSC_FALSE) {
124 PetscBool is_finalized;
125 PetscFinalized(&is_finalized);
126 if (!is_finalized)
127 PetscFinalize();
128 }
129
130 if (!mpiInitialised) {
131 int mpi_finalized;
132 MPI_Finalized(&mpi_finalized);
133 if (!mpi_finalized)
134 MPI_Finalize();
135 }
136 }
137
138 return 0;
139}
140
141// Use SFINAE to decide which template should be run,
142// if exist getSubInterfaceOptions run this one.
143template <class T>
144static auto get_sub_iface_options_imp(T *const ptr, int)
145 -> decltype(ptr->getSubInterfaceOptions()) {
146 return ptr->getSubInterfaceOptions();
147};
148
149// Use SFINAE to decide which template should be run,
150// if getSubInterfaceOptions not exist run this one.
151template <class T>
152static auto get_sub_iface_options_imp(T *const ptr, long) -> MoFEMErrorCode {
153 return 0;
154};
155
156template <class T>
157static auto get_event_options_imp(T *const ptr, int)
158 -> decltype(ptr->getEventOptions()) {
159 return ptr->getEventptions();
160};
161
162// Use SFINAE to decide which template should be run,
163// if getSubInterfaceOptions not exist run this one.
164// See SFINAE:
165// https://stackoverflow.com/questions/257288/is-it-possible-to-write-a-template-to-check-for-a-functions-existence
166// https://en.wikipedia.org/wiki/Substitution_failure_is_not_an_error
167template <class T>
168static auto get_event_options_imp(T *const ptr, long) -> MoFEMErrorCode {
169 return 0;
170};
171
172template <class IFACE> MoFEMErrorCode Core::regSubInterface() {
174 CHKERR registerInterface<IFACE>(true);
175 IFACE *ptr = new IFACE(*this);
176
177 // If sub interface has function getSubInterfaceOptions run
178 // it after construction. getSubInterfaceOptions is used to
179 // get parameters from command line.
180 auto get_sub_iface_options = [](auto *const ptr) {
181 return get_sub_iface_options_imp(ptr, 0);
182 };
183 CHKERR get_sub_iface_options(ptr);
184
185 auto type_idx = boost::typeindex::type_id<IFACE>();
186 iFaces.insert(type_idx, ptr);
188}
189
190template <class IFACE> MoFEMErrorCode Core::regEvents() {
192 auto ptr = boost::make_shared<IFACE>();
193 // See SFINAE:
194 // https://stackoverflow.com/questions/257288/is-it-possible-to-write-a-template-to-check-for-a-functions-existence
195 // https://en.wikipedia.org/wiki/Substitution_failure_is_not_an_error
196 auto get_event_options = [](auto *const ptr) {
197 return get_event_options_imp(ptr, 0);
198 };
199 CHKERR get_event_options(ptr.get());
201}
202
204 MPI_Comm comm, const int verbose) {
206
207 // This is deprecated ONE should use MoFEM::Core::Initialize
208 if (!isGloballyInitialised)
209 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
210 "MoFEM globally is not initialised, call MoFEM::Core::Initialize");
211
212 // Create duplicate communicator
213 wrapMPIMOABComm = boost::make_shared<WrapMPIComm>(comm, false);
214
215 MPI_Comm_size(mofemComm, &sIze);
216 MPI_Comm_rank(mofemComm, &rAnk);
217
218 // CHeck if moab has set communicator if not set communicator internally
219 ParallelComm *pComm = ParallelComm::get_pcomm(&moab, MYPCOMM_INDEX);
220 if (pComm == NULL)
221 pComm = new ParallelComm(&moab, wrapMPIMOABComm->get_comm());
222
223 // Register interfaces for this implementation
224 CHKERR registerInterface<UnknownInterface>();
225 CHKERR registerInterface<CoreInterface>();
226 CHKERR registerInterface<DeprecatedCoreInterface>();
227
228 // Register MOFEM events in PETSc
229 PetscLogEventRegister("FE_preProcess", 0, &MOFEM_EVENT_preProcess);
230 PetscLogEventRegister("FE_operator", 0, &MOFEM_EVENT_operator);
231 PetscLogEventRegister("FE_postProcess", 0, &MOFEM_EVENT_postProcess);
232 PetscLogEventRegister("MoFEMCreateMat", 0, &MOFEM_EVENT_createMat);
233
234 MOFEM_LOG_CHANNEL("WORLD");
235 MOFEM_LOG_CHANNEL("SELF");
236 MOFEM_LOG_CHANNEL("SYNC");
237
239}
240
241Core::CoreTmp(moab::Interface &moab, ///< MoAB interface
242 MPI_Comm comm, ///< MPI communicator
243 const int verbose ///< Verbosity level
244
245 )
246 : CoreTmp(moab, comm, verbose, CoreValue<0>()) {
247
248 // Register sub-interfaces
249 ierr = this->registerSubInterfaces();
250 CHKERRABORT(comm, ierr);
251 ierr = this->clearMap();
252 CHKERRABORT(comm, ierr);
253 ierr = this->getTags();
254 CHKERRABORT(comm, ierr);
255 ierr = this->getOptions(verbose);
256 CHKERRABORT(comm, ierr);
257
258 this->basicEntityDataPtr = boost::make_shared<BasicEntityData>(moab);
259 setRefEntBasicDataPtr(*this, this->basicEntityDataPtr);
260
261 ierr = this->initialiseDatabaseFromMesh(verbose);
262 CHKERRABORT(comm, ierr);
263}
264
265CoreTmp<-1>::CoreTmp(moab::Interface &moab, ///< MoAB interface
266 MPI_Comm comm, ///< MPI communicator
267 const int verbose ///< Verbosity level
268
269 )
270 : CoreTmp<0>(moab, comm, verbose, CoreValue<-1>()) {
271
272 // Register sub-interfaces
273 ierr = this->registerSubInterfaces();
274 CHKERRABORT(comm, ierr);
275 ierr = this->clearMap();
276 CHKERRABORT(comm, ierr);
277 ierr = this->getTags();
278 CHKERRABORT(comm, ierr);
279 ierr = this->getOptions(verbose);
280 CHKERRABORT(comm, ierr);
281
282 this->basicEntityDataPtr = boost::make_shared<BasicEntityData>(moab);
283 setRefEntBasicDataPtr(*this, this->basicEntityDataPtr);
284
285 ierr = this->initialiseDatabaseFromMesh(verbose);
286 CHKERRABORT(comm, ierr);
287}
288
290 PetscBool is_finalized = PETSC_FALSE;
291 PetscFinalized(&is_finalized);
292 if (!is_finalized && !iFaces.empty())
293 clearMap();
294 // Destroy interfaces
295 iFaces.clear();
296 // This is deprecated ONE should use MoFEM::Core::Initialize
297 if (isGloballyInitialised && is_finalized) {
298 isGloballyInitialised = false;
299 }
300}
301
303 BOOST_LOG_SCOPED_THREAD_ATTR("Timeline", attrs::timer());
304 MOFEM_LOG_CHANNEL("WORLD");
305
307 if (verb == -1)
308 verb = verbose;
309
310 Range ref_elems_to_add;
311
312 MOFEM_LOG("WORLD", Sev::verbose) << "Get MoFEM meshsets";
313 // Initialize database
314 Range meshsets;
315 CHKERR get_moab().get_entities_by_type(0, MBENTITYSET, meshsets, false);
316 Range special_meshsets;
317 for (auto mit : meshsets) {
318 BitFieldId field_id;
319 // Get bit id form field tag
320 CHKERR get_moab().tag_get_data(th_FieldId, &mit, 1, &field_id);
321 // Check if meshset if field meshset
322 if (field_id != 0) {
323
324 const void *tag_name;
325 int tag_name_size;
326 CHKERR get_moab().tag_get_by_ptr(
327 th_FieldName, &mit, 1, (const void **)&tag_name, &tag_name_size);
328
329 if (verb > QUIET)
330 MOFEM_LOG("WORLD", Sev::verbose)
331 << "Read field "
332 << boost::string_ref((char *)tag_name, tag_name_size);
333
334 auto p = fIelds.insert(boost::make_shared<Field>(moab, mit));
335
336 if (!p.second) {
337 // Field meshset exists, remove duplicate meshsets from other
338 // processors.
339 Range ents;
340 CHKERR get_moab().get_entities_by_handle(mit, ents, true);
341 CHKERR get_moab().add_entities((*p.first)->getMeshset(), ents);
342 CHKERR get_moab().delete_entities(&mit, 1);
343 } else {
344 special_meshsets.insert(mit);
345 }
346 }
347 // Check for finite elements
348 BitFieldId fe_id;
349 // Get bit id from fe tag
350 CHKERR get_moab().tag_get_data(th_FEId, &mit, 1, &fe_id);
351 // check if meshset is finite element meshset
352 if (fe_id != 0) {
353 std::pair<FiniteElement_multiIndex::iterator, bool> p =
354 finiteElements.insert(
355 boost::shared_ptr<FiniteElement>(new FiniteElement(moab, mit)));
356 if (verb > QUIET)
357 MOFEM_LOG("WORLD", Sev::verbose) << "Read finite element " << **p.first;
358
359 Range ents;
360 CHKERR get_moab().get_entities_by_type(mit, MBENTITYSET, ents, false);
361 CHKERR get_moab().get_entities_by_handle(mit, ents, true);
362 ref_elems_to_add.merge(ents);
363 if (!p.second) {
364 // Finite element mesh set exist, could be created on other processor.
365 // Remove duplicate.
366 CHKERR get_moab().add_entities((*p.first)->getMeshset(), ents);
367 CHKERR get_moab().delete_entities(&mit, 1);
368 } else {
369 special_meshsets.insert(mit);
370 }
371 }
372 BitProblemId problem_id;
373 // get bit id form problem tag
374 CHKERR get_moab().tag_get_data(th_ProblemId, &mit, 1, &problem_id);
375 // check if meshset if problem meshset
376 if (problem_id != 0) {
377 std::pair<Problem_multiIndex::iterator, bool> p =
378 pRoblems.insert(Problem(moab, mit));
379 if (verb > QUIET) {
380 MOFEM_LOG("WORLD", Sev::verbose) << "Read problem " << *p.first;
381 MOFEM_LOG("WORLD", Sev::noisy)
382 << "\tBitRef " << p.first->getBitRefLevel() << " BitMask "
383 << p.first->getBitRefLevelMask();
384 }
385
386 if (!p.second) {
387 // Problem meshset exists, could be created on other processor.
388 // Remove duplicate.
389 Range ents;
390 CHKERR get_moab().get_entities_by_handle(mit, ents, true);
391 CHKERR get_moab().get_entities_by_type(mit, MBENTITYSET, ents, true);
392 CHKERR get_moab().add_entities(p.first->meshset, ents);
393 CHKERR get_moab().delete_entities(&mit, 1);
394 } else {
395 special_meshsets.insert(mit);
396 }
397 }
398 }
399 MOFEM_LOG("WORLD", Sev::verbose) << "Get MoFEM meshsets <- done";
400
401 // Add entities to database
402 MOFEM_LOG("WORLD", Sev::verbose) << "Add entities to database";
403 Range bit_ref_ents;
404 CHKERR get_moab().get_entities_by_handle(0, bit_ref_ents, false);
405 bit_ref_ents = subtract(bit_ref_ents, special_meshsets);
406 CHKERR getInterface<BitRefManager>()->filterEntitiesByRefLevel(
407 BitRefLevel().set(), BitRefLevel().set(), bit_ref_ents);
408 CHKERR getInterface<BitRefManager>()->setEntitiesBitRefLevel(bit_ref_ents);
409 CHKERR getInterface<BitRefManager>()->setElementsBitRefLevel(
410 ref_elems_to_add);
411 MOFEM_LOG("WORLD", Sev::verbose) << "Add entities to database <- done";
412
413 // Build field entities
414 MOFEM_LOG("WORLD", Sev::verbose) << "Add field to database";
415 for (auto field : fIelds) {
416 if (field->getSpace() != NOSPACE) {
417 Range ents_of_id_meshset;
418 CHKERR get_moab().get_entities_by_handle(field->getMeshset(),
419 ents_of_id_meshset, false);
420 CHKERR set_field_order(ents_of_id_meshset, field->getId(), -1, verb);
421 }
422 }
423 MOFEM_LOG("WORLD", Sev::verbose) << "Add field to database <- done";
424
425 if (initaliseAndBuildField || initaliseAndBuildFiniteElements) {
426 MOFEM_LOG("WORLD", Sev::verbose) << "Build fields elements";
427 CHKERR build_fields(verb);
428 MOFEM_LOG("WORLD", Sev::verbose) << "Build fields elements <- done";
429 if (initaliseAndBuildFiniteElements) {
430 MOFEM_LOG("WORLD", Sev::verbose) << "Build finite elements";
431 CHKERR build_finite_elements(verb);
432 MOFEM_LOG("WORLD", Sev::verbose) << "Build finite elements <- done";
433 }
434 }
435
436 if (verb > VERY_NOISY) {
437 list_fields();
438 list_finite_elements();
439 list_problem();
440 }
441
442 // Initialize interfaces
443 MOFEM_LOG("WORLD", Sev::verbose) << "Initialise interfaces from mesh";
444
445 MOFEM_LOG("WORLD", Sev::verbose) << "Initialise MeshsetManager";
446 CHKERR getInterface<MeshsetsManager>() -> initialiseDatabaseFromMesh(verb);
447 MOFEM_LOG("WORLD", Sev::verbose) << "Initialise MeshsetManager <- done";
448 MOFEM_LOG("WORLD", Sev::verbose) << "Initialise SeriesRecorder";
449 CHKERR getInterface<SeriesRecorder>() -> initialiseDatabaseFromMesh(verb);
450 MOFEM_LOG("WORLD", Sev::verbose) << "Initialise SeriesRecorder <- done";
451
452 MOFEM_LOG("WORLD", Sev::verbose) << "Initialise interfaces from mesh <- done";
453
455}
456
457MoFEMErrorCode Core::setMoabInterface(moab::Interface &new_moab, int verb) {
459 if (verb == -1)
460 verb = verbose;
461
462 // clear moab database
463 CHKERR clearMap();
464
465 // set new reference
466 moab = std::ref(new_moab);
467
468 // check if moab has set communicator if not set communicator internally
469 ParallelComm *pComm = ParallelComm::get_pcomm(&new_moab, MYPCOMM_INDEX);
470 if (pComm == NULL) {
471 pComm = new ParallelComm(&new_moab, wrapMPIMOABComm->get_comm());
472 }
473
474 // create MoFEM tags
475 CHKERR getTags();
476
477 // Create basic entity data struture
478 basicEntityDataPtr = boost::make_shared<BasicEntityData>(moab);
479 setRefEntBasicDataPtr(*this, this->basicEntityDataPtr);
480
481 // Initalise database
482 CHKERR this->initialiseDatabaseFromMesh(verb);
483
485};
486
489
490 iFaces.clear();
491
492 // Register sub interfaces
493 CHKERR regSubInterface<LogManager>();
494 CHKERR regSubInterface<Simple>();
495 CHKERR regSubInterface<OperatorsTester>();
496 CHKERR regSubInterface<PipelineManager>();
497 CHKERR regSubInterface<PipelineGraph>();
498 CHKERR regSubInterface<ProblemsManager>();
499 CHKERR regSubInterface<MatrixManager>();
500 CHKERR regSubInterface<ISManager>();
501 CHKERR regSubInterface<VecManager>();
502 CHKERR regSubInterface<FieldBlas>();
503 CHKERR regSubInterface<BitRefManager>();
504 CHKERR regSubInterface<Tools>();
505 CHKERR regSubInterface<CommInterface>();
506 CHKERR regSubInterface<MeshsetsManager>();
507 CHKERR regSubInterface<NodeMergerInterface>();
508 CHKERR regSubInterface<PrismsFromSurfaceInterface>();
509 CHKERR regSubInterface<MeshRefinement>();
510 CHKERR regSubInterface<PrismInterface>();
511 CHKERR regSubInterface<CutMeshInterface>();
512 CHKERR regSubInterface<SeriesRecorder>();
513#ifdef WITH_TETGEN
514 CHKERR regSubInterface<TetGenInterface>();
515#endif
516#ifdef WITH_MED
517 CHKERR regSubInterface<MedInterface>();
518#endif
519 CHKERR regSubInterface<FieldEvaluatorInterface>();
520 CHKERR regSubInterface<BcManager>();
521
522 // Register events
523 CHKERR regEvents<SchurEvents>();
524
526};
527
530 // Cleaning databases in interfaces
531 CHKERR getInterface<SeriesRecorder>()->clearMap();
532 CHKERR getInterface<MeshsetsManager>()->clearMap();
533 CHKERR getInterface<CutMeshInterface>()->clearMap();
534 // Cleaning databases
535 refinedEntities.clear();
536 refinedFiniteElements.clear();
537 fIelds.clear();
538 entsFields.clear();
539 dofsField.clear();
540 finiteElements.clear();
541 entsFiniteElements.clear();
542 entFEAdjacencies.clear();
543 pRoblems.clear();
545}
546
549 if (verb == -1)
550 verb = verbose;
551 std::pair<RefEntity_multiIndex::iterator, bool> p_ent;
552 p_ent = refinedEntities.insert(
553 boost::make_shared<RefEntity>(basicEntityDataPtr, prism));
554 if (p_ent.second) {
555 std::pair<RefElement_multiIndex::iterator, bool> p;
556 p = refinedFiniteElements.insert(
557 boost::shared_ptr<RefElement>(new RefElement_PRISM(*p_ent.first)));
558 int num_nodes;
559 const EntityHandle *conn;
560 CHKERR get_moab().get_connectivity(prism, conn, num_nodes, true);
561 Range face_side3, face_side4;
562 CHKERR get_moab().get_adjacencies(conn, 3, 2, false, face_side3);
563 CHKERR get_moab().get_adjacencies(&conn[3], 3, 2, false, face_side4);
564 if (face_side3.size() != 1)
565 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
566 "prism don't have side face 3");
567 if (face_side4.size() != 1)
568 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
569 "prims don't have side face 4");
570 p.first->get()->getSideNumberPtr(*face_side3.begin());
571 p.first->get()->getSideNumberPtr(*face_side4.begin());
572 }
574}
575
578
579 const EntityHandle root_meshset = get_moab().get_root_set();
580 if (root_meshset) {
581 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
582 "Root meshset should be 0");
583 }
584
585 // Set version
586 {
587 Version version;
588 CHKERR getFileVersion(moab, version);
589 }
590
591 // Global Variables
592 {
593
594 auto check_tag_allocated = [](auto &rval) {
596 if (rval == MB_ALREADY_ALLOCATED)
597 rval = MB_SUCCESS;
598 else
599 CHKERRG(rval);
601 };
602
603 // Safety nets
604 int def_bool = 0;
605 rval = get_moab().tag_get_handle("_MoFEMBuild", 1, MB_TYPE_INTEGER,
606 th_MoFEMBuild, MB_TAG_CREAT | MB_TAG_MESH,
607 &def_bool);
608 CHKERR check_tag_allocated(rval);
609
610 CHKERR get_moab().tag_get_by_ptr(th_MoFEMBuild, &root_meshset, 1,
611 (const void **)&buildMoFEM);
612 }
613
614 // Tags saved in vtk-files
615 {
616 const int def_part = -1;
617 CHKERR get_moab().tag_get_handle("PARTITION", 1, MB_TYPE_INTEGER, th_Part,
618 MB_TAG_CREAT | MB_TAG_SPARSE, &def_part);
619 }
620
621 // Tags Ref
622 {
623
624 // Fix size of bir ref level tags
626
627 const int def_part = -1;
628 CHKERR get_moab().tag_get_handle("_MeshsetPartition", 1, MB_TYPE_INTEGER,
629 th_Part, MB_TAG_CREAT | MB_TAG_SPARSE,
630 &def_part);
631 EntityHandle def_handle = 0;
632 CHKERR get_moab().tag_get_handle("_RefParentHandle", 1, MB_TYPE_HANDLE,
633 th_RefParentHandle,
634 MB_TAG_CREAT | MB_TAG_SPARSE, &def_handle);
635 BitRefLevel def_bit_level = 0;
636 CHKERR get_moab().tag_get_handle(
637 "_RefBitLevel", sizeof(BitRefLevel), MB_TYPE_OPAQUE, th_RefBitLevel,
638 MB_TAG_CREAT | MB_TAG_BYTES | MB_TAG_SPARSE, &def_bit_level);
639 BitRefLevel def_bit_level_mask = BitRefLevel().set();
640 CHKERR get_moab().tag_get_handle(
641 "_RefBitLevelMask", sizeof(BitRefLevel), MB_TYPE_OPAQUE,
642 th_RefBitLevel_Mask, MB_TAG_CREAT | MB_TAG_BYTES | MB_TAG_SPARSE,
643 &def_bit_level_mask);
644 BitRefEdges def_bit_edge = 0;
645 CHKERR get_moab().tag_get_handle(
646 "_RefBitEdge", sizeof(BitRefEdges), MB_TYPE_OPAQUE, th_RefBitEdge,
647 MB_TAG_CREAT | MB_TAG_SPARSE | MB_TAG_BYTES, &def_bit_edge);
648 }
649
650 // Tags Field
651 {
652 const unsigned long int def_id = 0;
653 CHKERR get_moab().tag_get_handle(
654 "_FieldId", sizeof(BitFieldId), MB_TYPE_OPAQUE, th_FieldId,
655 MB_TAG_CREAT | MB_TAG_BYTES | MB_TAG_SPARSE, &def_id);
656 FieldSpace def_space = LASTSPACE;
657 CHKERR get_moab().tag_get_handle(
658 "_FieldSpace", sizeof(FieldSpace), MB_TYPE_OPAQUE, th_FieldSpace,
659 MB_TAG_CREAT | MB_TAG_BYTES | MB_TAG_SPARSE, &def_space);
660 FieldContinuity def_continuity = LASTCONTINUITY;
661 CHKERR get_moab().tag_get_handle(
662 "_FieldContinuity", sizeof(FieldContinuity), MB_TYPE_OPAQUE,
663 th_FieldContinuity, MB_TAG_CREAT | MB_TAG_BYTES | MB_TAG_SPARSE,
664 &def_continuity);
666 CHKERR get_moab().tag_get_handle(
667 "_FieldBase", sizeof(FieldApproximationBase), MB_TYPE_OPAQUE,
668 th_FieldBase, MB_TAG_CREAT | MB_TAG_BYTES | MB_TAG_SPARSE, &def_base);
669 const int def_val_len = 0;
670 CHKERR get_moab().tag_get_handle(
671 "_FieldName", def_val_len, MB_TYPE_OPAQUE, th_FieldName,
672 MB_TAG_CREAT | MB_TAG_BYTES | MB_TAG_VARLEN | MB_TAG_SPARSE, NULL);
673 CHKERR get_moab().tag_get_handle(
674 "_FieldName_DataNamePrefix", def_val_len, MB_TYPE_OPAQUE,
675 th_FieldName_DataNamePrefix,
676 MB_TAG_CREAT | MB_TAG_BYTES | MB_TAG_VARLEN | MB_TAG_SPARSE, NULL);
677 }
678
679 // Tags FE
680 {
681 const unsigned long int def_id = 0;
682 const int def_val_len = 0;
683 CHKERR get_moab().tag_get_handle(
684 "_FEId", sizeof(BitFEId), MB_TYPE_OPAQUE, th_FEId,
685 MB_TAG_CREAT | MB_TAG_BYTES | MB_TAG_SPARSE, &def_id);
686 CHKERR get_moab().tag_get_handle(
687 "_FEName", def_val_len, MB_TYPE_OPAQUE, th_FEName,
688 MB_TAG_CREAT | MB_TAG_BYTES | MB_TAG_VARLEN | MB_TAG_SPARSE, NULL);
689 CHKERR get_moab().tag_get_handle(
690 "_FEIdCol", sizeof(BitFieldId), MB_TYPE_OPAQUE, th_FEIdCol,
691 MB_TAG_CREAT | MB_TAG_BYTES | MB_TAG_SPARSE, &def_id);
692 CHKERR get_moab().tag_get_handle(
693 "_FEIdRow", sizeof(BitFieldId), MB_TYPE_OPAQUE, th_FEIdRow,
694 MB_TAG_CREAT | MB_TAG_BYTES | MB_TAG_SPARSE, &def_id);
695 CHKERR get_moab().tag_get_handle(
696 "_FEIdData", sizeof(BitFieldId), MB_TYPE_OPAQUE, th_FEIdData,
697 MB_TAG_CREAT | MB_TAG_BYTES | MB_TAG_SPARSE, &def_id);
698 }
699
700 // Tags Problem
701 {
702 const unsigned long int def_id = 0;
703 const int def_val_len = 0;
704 CHKERR get_moab().tag_get_handle(
705 "_ProblemId", sizeof(BitProblemId), MB_TYPE_OPAQUE, th_ProblemId,
706 MB_TAG_CREAT | MB_TAG_BYTES | MB_TAG_SPARSE, &def_id);
707 CHKERR get_moab().tag_get_handle(
708 "_ProblemFEId", sizeof(BitFEId), MB_TYPE_OPAQUE, th_ProblemFEId,
709 MB_TAG_CREAT | MB_TAG_BYTES | MB_TAG_SPARSE, &def_id);
710 CHKERR get_moab().tag_get_handle(
711 "_ProblemName", def_val_len, MB_TYPE_OPAQUE, th_ProblemName,
712 MB_TAG_CREAT | MB_TAG_BYTES | MB_TAG_VARLEN | MB_TAG_SPARSE, NULL);
713 DofIdx def_nbdofs = 0;
714 CHKERR get_moab().tag_get_handle(
715 "_ProblemNbDofsRow", sizeof(DofIdx), MB_TYPE_OPAQUE,
716 th_ProblemNbDofsRow, MB_TAG_CREAT | MB_TAG_BYTES | MB_TAG_SPARSE,
717 &def_nbdofs);
718 CHKERR get_moab().tag_get_handle(
719 "_ProblemNbDofsCol", sizeof(DofIdx), MB_TYPE_OPAQUE,
720 th_ProblemNbDofsCol, MB_TAG_CREAT | MB_TAG_BYTES | MB_TAG_SPARSE,
721 &def_nbdofs);
722 CHKERR get_moab().tag_get_handle(
723 "_ProblemLocalNbDofsRow", sizeof(DofIdx), MB_TYPE_OPAQUE,
724 th_ProblemLocalNbDofRow, MB_TAG_CREAT | MB_TAG_BYTES | MB_TAG_SPARSE,
725 &def_nbdofs);
726 CHKERR get_moab().tag_get_handle(
727 "_ProblemGhostNbDofsRow", sizeof(DofIdx), MB_TYPE_OPAQUE,
728 th_ProblemGhostNbDofRow, MB_TAG_CREAT | MB_TAG_BYTES | MB_TAG_SPARSE,
729 &def_nbdofs);
730 CHKERR get_moab().tag_get_handle(
731 "_ProblemLocalNbDofsCol", sizeof(DofIdx), MB_TYPE_OPAQUE,
732 th_ProblemLocalNbDofCol, MB_TAG_CREAT | MB_TAG_BYTES | MB_TAG_SPARSE,
733 &def_nbdofs);
734 CHKERR get_moab().tag_get_handle(
735 "_ProblemGhostNbDofsCol", sizeof(DofIdx), MB_TYPE_OPAQUE,
736 th_ProblemGhostNbDofCol, MB_TAG_CREAT | MB_TAG_BYTES | MB_TAG_SPARSE,
737 &def_nbdofs);
738 }
739
740 // Meshsets with boundary conditions and material sets
741 MeshsetsManager *meshsets_manager_ptr;
742 CHKERR getInterface(meshsets_manager_ptr);
743 CHKERR meshsets_manager_ptr->getTags(verb);
744
745 // Series recorder
746 SeriesRecorder *series_recorder_ptr;
747 CHKERR getInterface(series_recorder_ptr);
748 CHKERR series_recorder_ptr->getTags(verb);
749
751}
752
755 if (verb == -1)
756 verb = verbose;
757 CHKERR clearMap();
759}
760
763 if (verb == -1)
764 verb = verbose;
765 CHKERR this->clearMap();
766 CHKERR this->getTags(verb);
767 CHKERR this->initialiseDatabaseFromMesh(verb);
769}
770
771MoFEMErrorCode Core::set_moab_interface(moab::Interface &new_moab, int verb) {
772 return this->setMoabInterface(new_moab, verb);
773};
774
776 int verb) {
777 return this->setMoabInterface(new_moab, verb);
778};
779
782 if (verb == -1)
783 verb = verbose;
784
785 PetscOptionsBegin(mofemComm, optionsPrefix.c_str(), "Mesh cut options",
786 "See MoFEM documentation");
787
788 CHKERR PetscOptionsBool(
789 "-mofem_init_fields", "Initialise fields on construction", "",
790 initaliseAndBuildField, &initaliseAndBuildField, NULL);
791
792 CHKERR PetscOptionsBool(
793 "-mofem_init_fields", "Initialise fields on construction", "",
794 initaliseAndBuildFiniteElements, &initaliseAndBuildFiniteElements, NULL);
795
796 // TODO: Add read verbosity level
797 // TODO: Add option to initalise problems ??? - DO WE REALLY NEED THAT
798
799 PetscOptionsEnd();
800
802}
803
804// cubit meshsets
805
808 *fields_ptr = &fIelds;
810}
811
813Core::get_ref_ents(const RefEntity_multiIndex **refined_entities_ptr) const {
815 *refined_entities_ptr = &refinedEntities;
817}
819 const RefElement_multiIndex **refined_finite_elements_ptr) const {
821 *refined_finite_elements_ptr = &refinedFiniteElements;
823}
824
825MoFEMErrorCode Core::get_problem(const std::string &problem_name,
826 const Problem **problem_ptr) const {
828 typedef Problem_multiIndex::index<Problem_mi_tag>::type ProblemsByName;
829 const ProblemsByName &problems = pRoblems.get<Problem_mi_tag>();
830 ProblemsByName::iterator p_miit = problems.find(problem_name);
831 if (p_miit == problems.end()) {
832 SETERRQ(PETSC_COMM_SELF, MOFEM_OPERATION_UNSUCCESSFUL,
833 "problem < %s > not found, (top tip: check spelling)",
834 problem_name.c_str());
835 }
836 *problem_ptr = &*p_miit;
838}
839
841Core::get_problems(const Problem_multiIndex **problems_ptr) const {
843 *problems_ptr = &pRoblems;
845}
846
850 *field_ents = &entsFields;
852}
855 *dofs_ptr = &dofsField;
857}
858
862 *fe_ptr = &finiteElements;
864}
865
867 const EntFiniteElement_multiIndex **fe_ent_ptr) const {
869 *fe_ent_ptr = &entsFiniteElements;
871}
872
874 MeshsetsManager *meshsets_manager_ptr;
875 getInterface(meshsets_manager_ptr);
876 return meshsets_manager_ptr;
877}
878
880 MeshsetsManager *meshsets_manager_ptr;
881 getInterface(meshsets_manager_ptr);
882 return meshsets_manager_ptr;
883}
884
887 *dofs_elements_adjacency) const {
889 *dofs_elements_adjacency = &entFEAdjacencies;
891}
892
895 return &entFEAdjacencies;
896}
897
898const Field_multiIndex *Core::get_fields() const { return &fIelds; }
900 return &refinedEntities;
901}
903 return &refinedFiniteElements;
904}
906 return &finiteElements;
907}
909 return &entsFiniteElements;
910}
912 return &entsFields;
913}
914const DofEntity_multiIndex *Core::get_dofs() const { return &dofsField; }
915const Problem *Core::get_problem(const std::string problem_name) const {
916 const Problem *prb;
917 CHK_THROW_MESSAGE(get_problem(problem_name, &prb),
918 "Problem of given name not found");
919 return prb;
920}
921const Problem_multiIndex *Core::get_problems() const { return &pRoblems; }
922
923template <int V, typename std::enable_if<(V >= 0), int>::type * = nullptr>
924void set_ref_ent_basic_data_ptr_impl(boost::shared_ptr<BasicEntityData> &ptr) {
926};
927
928template <int V, typename std::enable_if<(V < 0), int>::type * = nullptr>
929void set_ref_ent_basic_data_ptr_impl(boost::shared_ptr<BasicEntityData> &ptr) {
930 return;
931};
932
933void Core::setRefEntBasicDataPtr(MoFEM::Interface &m_field,
934 boost::shared_ptr<BasicEntityData> &ptr) {
935
936 switch (m_field.getValue()) {
937 case -1:
938 set_ref_ent_basic_data_ptr_impl<-1>(ptr);
939 break;
940 case 0:
941 set_ref_ent_basic_data_ptr_impl<0>(ptr);
942 break;
943 case 1:
944 set_ref_ent_basic_data_ptr_impl<1>(ptr);
945 break;
946 default:
947 THROW_MESSAGE("Core index can vary from -1 to MAX_CORE_TMP");
948 }
949};
950
951boost::shared_ptr<RefEntityTmp<0>>
952Core::makeSharedRefEntity(MoFEM::Interface &m_field, const EntityHandle ent) {
953
954 boost::shared_ptr<RefEntityTmp<0>> ref_ent_ptr;
955
956 switch (m_field.getValue()) {
957 case -1:
958 ref_ent_ptr = boost::shared_ptr<RefEntityTmp<0>>(
959
960 new RefEntityTmp<-1>(m_field.get_basic_entity_data_ptr(), ent)
961
962 );
963 break;
964 case 0:
965 ref_ent_ptr = boost::shared_ptr<RefEntityTmp<0>>(
966
967 new RefEntityTmp<0>(m_field.get_basic_entity_data_ptr(), ent)
968
969 );
970 break;
971 case 1:
972 ref_ent_ptr = boost::shared_ptr<RefEntityTmp<0>>(
973
974 new RefEntityTmp<1>(m_field.get_basic_entity_data_ptr(), ent)
975
976 );
977 break;
978 default:
979 THROW_MESSAGE("Core index can vary from -1 to MAX_CORE_TMP");
980 }
981
982 return ref_ent_ptr;
983}
984
985boost::shared_ptr<RefEntityTmp<0>>
986Core::make_shared_ref_entity(const EntityHandle ent) {
987 return this->makeSharedRefEntity(*this, ent);
988}
989
990boost::shared_ptr<RefEntityTmp<0>>
992 return this->makeSharedRefEntity(*this, ent);
993}
994
995} // namespace MoFEM
multi_index_container< FieldEntityEntFiniteElementAdjacencyMap, indexed_by< ordered_unique< tag< Composite_Unique_mi_tag >, composite_key< FieldEntityEntFiniteElementAdjacencyMap, const_mem_fun< FieldEntityEntFiniteElementAdjacencyMap, UId, &FieldEntityEntFiniteElementAdjacencyMap::getEntUniqueId >, const_mem_fun< FieldEntityEntFiniteElementAdjacencyMap, UId, &FieldEntityEntFiniteElementAdjacencyMap::getFeUniqueId > > >, ordered_non_unique< tag< Unique_mi_tag >, const_mem_fun< FieldEntityEntFiniteElementAdjacencyMap, UId, &FieldEntityEntFiniteElementAdjacencyMap::getEntUniqueId > >, ordered_non_unique< tag< FE_Unique_mi_tag >, const_mem_fun< FieldEntityEntFiniteElementAdjacencyMap, UId, &FieldEntityEntFiniteElementAdjacencyMap::getFeUniqueId > >, ordered_non_unique< tag< FEEnt_mi_tag >, const_mem_fun< FieldEntityEntFiniteElementAdjacencyMap, EntityHandle, &FieldEntityEntFiniteElementAdjacencyMap::getFeHandle > >, ordered_non_unique< tag< Ent_mi_tag >, const_mem_fun< FieldEntityEntFiniteElementAdjacencyMap, EntityHandle, &FieldEntityEntFiniteElementAdjacencyMap::getEntHandle > > > > FieldEntityEntFiniteElementAdjacencyMap_multiIndex
MultiIndex container keeps Adjacencies Element and dof entities adjacencies and vice versa.
void macro_is_deprecated_using_deprecated_function()
Is used to indicate that macro is deprecated Do nothing just triggers error at the compilation.
Definition Core.cpp:11
static PetscErrorCode mofem_error_handler(MPI_Comm comm, int line, const char *fun, const char *file, PetscErrorCode n, PetscErrorType p, const char *mess, void *ctx)
multi_index_container< boost::shared_ptr< Field >, indexed_by< hashed_unique< tag< BitFieldId_mi_tag >, const_mem_fun< Field, const BitFieldId &, &Field::getId >, HashBit< BitFieldId >, EqBit< BitFieldId > >, ordered_unique< tag< Meshset_mi_tag >, member< Field, EntityHandle, &Field::meshSet > >, ordered_unique< tag< FieldName_mi_tag >, const_mem_fun< Field, boost::string_ref, &Field::getNameRef > >, ordered_non_unique< tag< BitFieldId_space_mi_tag >, const_mem_fun< Field, FieldSpace, &Field::getSpace > > > > Field_multiIndex
Multi-index container for field storage and retrieval.
#define MOFEM_LOG_C(channel, severity, format,...)
static char help[]
@ QUIET
@ VERY_NOISY
FieldApproximationBase
approximation base
Definition definitions.h:58
@ LASTBASE
Definition definitions.h:69
#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
@ LASTSPACE
FieldSpace in [ 0, LASTSPACE )
Definition definitions.h:89
@ NOSPACE
Definition definitions.h:83
#define MYPCOMM_INDEX
default communicator number PCOMM
FieldContinuity
Field continuity.
Definition definitions.h:99
@ LASTCONTINUITY
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
#define CHKERRG(n)
Check error code of MoFEM/MOAB/PETSc function.
@ MOFEM_OPERATION_UNSUCCESSFUL
Definition definitions.h:34
@ MOFEM_DATA_INCONSISTENCY
Definition definitions.h:31
@ MOFEM_SUCCESS
Definition definitions.h:30
#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.
multi_index_container< boost::shared_ptr< DofEntity >, indexed_by< ordered_unique< tag< Unique_mi_tag >, const_mem_fun< DofEntity, UId, &DofEntity::getLocalUniqueId > >, ordered_non_unique< tag< Ent_mi_tag >, const_mem_fun< DofEntity, EntityHandle, &DofEntity::getEnt > > > > DofEntity_multiIndex
MultiIndex container keeps DofEntity.
multi_index_container< boost::shared_ptr< EntFiniteElement >, indexed_by< ordered_unique< tag< Unique_mi_tag >, const_mem_fun< EntFiniteElement, UId, &EntFiniteElement::getLocalUniqueId > >, ordered_non_unique< tag< Ent_mi_tag >, const_mem_fun< EntFiniteElement::interface_type_RefEntity, EntityHandle, &EntFiniteElement::getEnt > > > > EntFiniteElement_multiIndex
MultiIndex container for EntFiniteElement.
multi_index_container< Problem, indexed_by< ordered_unique< tag< Meshset_mi_tag >, member< Problem, EntityHandle, &Problem::meshset > >, hashed_unique< tag< BitProblemId_mi_tag >, const_mem_fun< Problem, BitProblemId, &Problem::getId >, HashBit< BitProblemId >, EqBit< BitProblemId > >, hashed_unique< tag< Problem_mi_tag >, const_mem_fun< Problem, std::string, &Problem::getName > > > > Problem_multiIndex
MultiIndex for entities for Problem.
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.
#define MOFEM_LOG(channel, severity)
Log.
#define MOFEM_LOG_CHANNEL(channel)
Set and reset channel.
MoFEMErrorCode getTags(int verb=-1)
Get tag handles used on meshsets.
FTensor::Index< 'i', SPACE_DIM > i
static MoFEMErrorCodeGeneric< PetscErrorCode > ierr
static MoFEMErrorCodeGeneric< moab::ErrorCode > rval
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
std::bitset< BITFEID_SIZE > BitFEId
Finite element Id.
Definition Types.hpp:43
int DofIdx
Index of DOF.
Definition Types.hpp:18
std::bitset< BITPROBLEMID_SIZE > BitProblemId
Problem Id.
Definition Types.hpp:44
std::bitset< BITFIELDID_SIZE > BitFieldId
Field Id.
Definition Types.hpp:42
std::bitset< BITREFEDGES_SIZE > BitRefEdges
Definition Types.hpp:34
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
multi_index_container< boost::shared_ptr< RefEntity >, indexed_by< ordered_unique< tag< Ent_mi_tag >, const_mem_fun< RefEntity, EntityHandle, &RefEntity::getEnt > >, ordered_non_unique< tag< Ent_Ent_mi_tag >, const_mem_fun< RefEntity, EntityHandle, &RefEntity::getParentEnt > >, ordered_non_unique< tag< Composite_EntType_and_ParentEntType_mi_tag >, composite_key< RefEntity, const_mem_fun< RefEntity, EntityType, &RefEntity::getEntType >, const_mem_fun< RefEntity, EntityType, &RefEntity::getParentEntType > > >, ordered_non_unique< tag< Composite_ParentEnt_And_EntType_mi_tag >, composite_key< RefEntity, const_mem_fun< RefEntity, EntityType, &RefEntity::getEntType >, const_mem_fun< RefEntity, EntityHandle, &RefEntity::getParentEnt > > > > > RefEntity_multiIndex
static auto get_sub_iface_options_imp(T *const ptr, int) -> decltype(ptr->getSubInterfaceOptions())
Definition Core.cpp:144
void set_ref_ent_basic_data_ptr_impl(boost::shared_ptr< BasicEntityData > &ptr)
Definition Core.cpp:924
multi_index_container< boost::shared_ptr< RefElement >, indexed_by< ordered_unique< tag< Ent_mi_tag >, const_mem_fun< RefElement::interface_type_RefEntity, EntityHandle, &RefElement::getEnt > > > > RefElement_multiIndex
static auto get_event_options_imp(T *const ptr, int) -> decltype(ptr->getEventOptions())
Definition Core.cpp:157
multi_index_container< boost::shared_ptr< FieldEntity >, indexed_by< ordered_unique< tag< Unique_mi_tag >, member< FieldEntity, UId, &FieldEntity::localUId > >, ordered_non_unique< tag< Ent_mi_tag >, const_mem_fun< FieldEntity::interface_type_RefEntity, EntityHandle, &FieldEntity::getEnt > > > > FieldEntity_multiIndex
static MoFEMErrorCode fixTagSize(moab::Interface &moab, bool *changed=nullptr)
Fix tag size when BITREFLEVEL_SIZE of core library is different than file BITREFLEVEL_SIZE.
Core (interface) class.
Definition Core.hpp:82
MoFEMErrorCode clearMap()
Cleaning database.
Definition Core.cpp:528
MoFEMErrorCode query_interface(boost::typeindex::type_index type_index, UnknownInterface **iface) const
Getting interface of core database.
Definition Core.cpp:43
static MoFEMErrorCode Initialize(int *argc, char ***args, const char file[], const char help[])
Initializes the MoFEM database PETSc, MOAB and MPI.
Definition Core.cpp:72
MoFEMErrorCode regSubInterface()
Register sub-interfaces in core interface.
Definition Core.cpp:172
static bool isGloballyInitialised
Core base globally initialized.
Definition Core.hpp:1056
const RefElement_multiIndex * get_ref_finite_elements() const
Get the ref finite elements object.
Definition Core.cpp:902
const FieldEntity_multiIndex * get_field_ents() const
Get the field ents object.
Definition Core.cpp:911
MoFEMErrorCode initialiseDatabaseFromMesh(int verb=DEFAULT_VERBOSITY)
Initialize database getting information on mesh.
Definition Core.cpp:302
MoFEMErrorCode registerSubInterfaces()
Register insterfaces.
Definition Core.cpp:487
const Field_multiIndex * get_fields() const
Get the fields object.
Definition Core.cpp:898
const Problem_multiIndex * get_problems() const
Get the problems object.
Definition Core.cpp:921
MoFEMErrorCode coreGenericConstructor(moab::Interface &moab, MPI_Comm comm, const int verbose)
Definition Core.cpp:203
const DofEntity_multiIndex * get_dofs() const
Get the dofs object.
Definition Core.cpp:914
MoFEMErrorCode getOptions(int verb=DEFAULT_VERBOSITY)
Get core options from command line.
Definition Core.cpp:780
static PetscBool isInitialized
petsc was initialised by other agent
Definition Core.hpp:1058
const RefEntity_multiIndex * get_ref_ents() const
Get the ref ents object.
Definition Core.cpp:899
const FieldEntityEntFiniteElementAdjacencyMap_multiIndex * get_ents_elements_adjacency() const
Get the dofs elements adjacency object.
Definition Core.cpp:894
const EntFiniteElement_multiIndex * get_ents_finite_elements() const
Get the ents finite elements object.
Definition Core.cpp:908
const FiniteElement_multiIndex * get_finite_elements() const
Get the finite elements object.
Definition Core.cpp:905
MoFEMErrorCode addPrismToDatabase(const EntityHandle prism, int verb=DEFAULT_VERBOSITY)
add prim element
Definition Core.cpp:547
MoFEMErrorCode regEvents()
Register petsc events.
Definition Core.cpp:190
virtual ~CoreTmp()
Definition Core.cpp:289
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
Definition Core.cpp:118
MoFEMErrorCode set_moab_interface(moab::Interface &new_moab, int verb=VERBOSE)
Set the moab interface object.
Definition Core.cpp:771
MoFEMErrorCode setMoabInterface(moab::Interface &new_moab, int verb=VERBOSE)
Definition Core.cpp:457
MoFEMErrorCode get_problem(const std::string &problem_name, const Problem **problem_ptr) const
Get problem database (data structure)
Definition Core.cpp:825
MoFEMErrorCode getTags(int verb=DEFAULT_VERBOSITY)
Get tag handles.
Definition Core.cpp:576
CoreTmp(moab::Interface &moab, MPI_Comm comm=PETSC_COMM_WORLD, const int verbose=VERBOSE)
Definition Core.cpp:241
static int mpiInitialised
mpi was initialised by other agent
Definition Core.hpp:1057
MeshsetsManager * get_meshsets_manager_ptr()
get MeshsetsManager pointer
Definition Core.cpp:873
MoFEMErrorCode rebuild_database(int verb=DEFAULT_VERBOSITY)
Clear database and initialize it once again.
Definition Core.cpp:761
MoFEMErrorCode clear_database(int verb=DEFAULT_VERBOSITY)
Clear database.
Definition Core.cpp:753
MoFEMErrorCode set_moab_interface(moab::Interface &new_moab, int verb)
CoreTmp(moab::Interface &moab, MPI_Comm comm=PETSC_COMM_WORLD, const int verbose=VERBOSE)
Deprecated interface functions.
Finite element definition.
static MoFEMErrorCode getOptions()
Get logger option.
static void createDefaultSinks(MPI_Comm comm)
Create default sinks.
static PetscErrorCode logPetscFPrintf(FILE *fd, const char format[], va_list Argp)
Use to handle PETSc output.
Interface for managing meshsets containing materials and boundary conditions.
keeps basic data about problem
keeps data about abstract PRISM finite element
MoFEMErrorCode getTags(int verb=-1)
get tags handlers used on meshsets
base class for all interface classes
MPI_Comm duplicatedComm
Definition Core.hpp:27
WrapMPIComm(MPI_Comm comm, bool petsc)
Definition Core.cpp:16
MPI_Comm comm
Definition Core.hpp:26