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