v0.14.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_NULL);
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 IFACE> MoFEMErrorCode Core::regSubInterface() {
152 CHKERR registerInterface<IFACE>(true);
153 IFACE *ptr = new IFACE(*this);
154
155 // If sub interface has function getSubInterfaceOptions run
156 // it after construction. getSubInterfaceOptions is used to
157 // get parameters from command line.
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
161 auto get_sub_iface_options = [](auto *const ptr) {
162 return get_sub_iface_options_imp(ptr, 0);
163 };
164 CHKERR get_sub_iface_options(ptr);
165
166 auto type_idx = boost::typeindex::type_id<IFACE>();
167 iFaces.insert(type_idx, ptr);
169}
170
172 MPI_Comm comm, const int verbose) {
174
175 // This is deprecated ONE should use MoFEM::Core::Initialize
176 if (!isGloballyInitialised)
177 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
178 "MoFEM globally is not initialised, call MoFEM::Core::Initialize");
179
180 // Create duplicate communicator
181 wrapMPIMOABComm = boost::make_shared<WrapMPIComm>(comm, false);
182
183 MPI_Comm_size(mofemComm, &sIze);
184 MPI_Comm_rank(mofemComm, &rAnk);
185
186 // CHeck if moab has set communicator if not set communicator interbally
187 ParallelComm *pComm = ParallelComm::get_pcomm(&moab, MYPCOMM_INDEX);
188 if (pComm == NULL)
189 pComm = new ParallelComm(&moab, wrapMPIMOABComm->get_comm());
190
191 // Register interfaces for this implementation
192 CHKERR registerInterface<UnknownInterface>();
193 CHKERR registerInterface<CoreInterface>();
194 CHKERR registerInterface<DeprecatedCoreInterface>();
195
196 // Register MOFEM events in PETSc
197 PetscLogEventRegister("FE_preProcess", 0, &MOFEM_EVENT_preProcess);
198 PetscLogEventRegister("FE_operator", 0, &MOFEM_EVENT_operator);
199 PetscLogEventRegister("FE_postProcess", 0, &MOFEM_EVENT_postProcess);
200 PetscLogEventRegister("MoFEMCreateMat", 0, &MOFEM_EVENT_createMat);
201
202 MOFEM_LOG_CHANNEL("WORLD");
203 MOFEM_LOG_CHANNEL("SELF");
204 MOFEM_LOG_CHANNEL("SYNC");
205
207}
208
209Core::CoreTmp(moab::Interface &moab, ///< MoAB interface
210 MPI_Comm comm, ///< MPI communicator
211 const int verbose ///< Verbosity level
212
213 )
214 : CoreTmp(moab, comm, verbose, CoreValue<0>()) {
215
216 // Register sub-interfaces
217 ierr = this->registerSubInterfaces();
218 CHKERRABORT(comm, ierr);
219 ierr = this->clearMap();
220 CHKERRABORT(comm, ierr);
221 ierr = this->getTags();
222 CHKERRABORT(comm, ierr);
223 ierr = this->getOptions(verbose);
224 CHKERRABORT(comm, ierr);
225
226 this->basicEntityDataPtr = boost::make_shared<BasicEntityData>(moab);
227 setRefEntBasicDataPtr(*this, this->basicEntityDataPtr);
228
229 ierr = this->initialiseDatabaseFromMesh(verbose);
230 CHKERRABORT(comm, ierr);
231}
232
233CoreTmp<-1>::CoreTmp(moab::Interface &moab, ///< MoAB interface
234 MPI_Comm comm, ///< MPI communicator
235 const int verbose ///< Verbosity level
236
237 )
238 : CoreTmp<0>(moab, comm, verbose, CoreValue<-1>()) {
239
240 // Register sub-interfaces
241 ierr = this->registerSubInterfaces();
242 CHKERRABORT(comm, ierr);
243 ierr = this->clearMap();
244 CHKERRABORT(comm, ierr);
245 ierr = this->getTags();
246 CHKERRABORT(comm, ierr);
247 ierr = this->getOptions(verbose);
248 CHKERRABORT(comm, ierr);
249
250 this->basicEntityDataPtr = boost::make_shared<BasicEntityData>(moab);
251 setRefEntBasicDataPtr(*this, this->basicEntityDataPtr);
252
253 ierr = this->initialiseDatabaseFromMesh(verbose);
254 CHKERRABORT(comm, ierr);
255}
256
258 PetscBool is_finalized;
259 PetscFinalized(&is_finalized);
260 // Destroy interfaces
261 iFaces.clear();
262 // This is deprecated ONE should use MoFEM::Core::Initialize
263 if (isGloballyInitialised && is_finalized) {
264 isGloballyInitialised = false;
265 }
266}
267
269 MOFEM_LOG_CHANNEL("WORLD");
271 if (verb == -1)
272 verb = verbose;
273
274 Range ref_elems_to_add;
275
276 // Initialize database
277 Range meshsets;
278 CHKERR get_moab().get_entities_by_type(0, MBENTITYSET, meshsets, false);
279 Range special_meshsets;
280 for (auto mit : meshsets) {
281 BitFieldId field_id;
282 // Get bit id form field tag
283 CHKERR get_moab().tag_get_data(th_FieldId, &mit, 1, &field_id);
284 // Check if meshset if field meshset
285 if (field_id != 0) {
286
287 const void *tag_name;
288 int tag_name_size;
289 CHKERR get_moab().tag_get_by_ptr(
290 th_FieldName, &mit, 1, (const void **)&tag_name, &tag_name_size);
291
292 if (verb > QUIET)
293 MOFEM_LOG("WORLD", Sev::verbose)
294 << "Read field "
295 << boost::string_ref((char *)tag_name, tag_name_size);
296
297 auto p = fIelds.insert(boost::make_shared<Field>(moab, mit));
298
299 if (!p.second) {
300 // Field meshset exists, remove duplicate meshsets from other
301 // processors.
302 Range ents;
303 CHKERR get_moab().get_entities_by_handle(mit, ents, true);
304 CHKERR get_moab().add_entities((*p.first)->getMeshset(), ents);
305 CHKERR get_moab().delete_entities(&mit, 1);
306 } else {
307 special_meshsets.insert(mit);
308 }
309 }
310 // Check for finite elements
311 BitFieldId fe_id;
312 // Get bit id from fe tag
313 CHKERR get_moab().tag_get_data(th_FEId, &mit, 1, &fe_id);
314 // check if meshset is finite element meshset
315 if (fe_id != 0) {
316 std::pair<FiniteElement_multiIndex::iterator, bool> p =
317 finiteElements.insert(
318 boost::shared_ptr<FiniteElement>(new FiniteElement(moab, mit)));
319 if (verb > QUIET)
320 MOFEM_LOG("WORLD", Sev::verbose) << "Read finite element " << **p.first;
321
322 Range ents;
323 CHKERR get_moab().get_entities_by_type(mit, MBENTITYSET, ents, false);
324 CHKERR get_moab().get_entities_by_handle(mit, ents, true);
325 ref_elems_to_add.merge(ents);
326 if (!p.second) {
327 // Finite element mesh set exist, could be created on other processor.
328 // Remove duplicate.
329 CHKERR get_moab().add_entities((*p.first)->getMeshset(), ents);
330 CHKERR get_moab().delete_entities(&mit, 1);
331 } else {
332 special_meshsets.insert(mit);
333 }
334 }
335 BitProblemId problem_id;
336 // get bit id form problem tag
337 CHKERR get_moab().tag_get_data(th_ProblemId, &mit, 1, &problem_id);
338 // check if meshset if problem meshset
339 if (problem_id != 0) {
340 std::pair<Problem_multiIndex::iterator, bool> p =
341 pRoblems.insert(Problem(moab, mit));
342 if (verb > QUIET) {
343 MOFEM_LOG("WORLD", Sev::verbose) << "Read problem " << *p.first;
344 MOFEM_LOG("WORLD", Sev::noisy)
345 << "\tBitRef " << p.first->getBitRefLevel() << " BitMask "
346 << p.first->getBitRefLevelMask();
347 }
348
349 if (!p.second) {
350 // Problem meshset exists, could be created on other processor.
351 // Remove duplicate.
352 Range ents;
353 CHKERR get_moab().get_entities_by_handle(mit, ents, true);
354 CHKERR get_moab().get_entities_by_type(mit, MBENTITYSET, ents, true);
355 CHKERR get_moab().add_entities(p.first->meshset, ents);
356 CHKERR get_moab().delete_entities(&mit, 1);
357 } else {
358 special_meshsets.insert(mit);
359 }
360 }
361 }
362
363 // Add entities to database
364 Range bit_ref_ents;
365 CHKERR get_moab().get_entities_by_handle(0, bit_ref_ents, false);
366 bit_ref_ents = subtract(bit_ref_ents, special_meshsets);
367 CHKERR getInterface<BitRefManager>()->filterEntitiesByRefLevel(
368 BitRefLevel().set(), BitRefLevel().set(), bit_ref_ents);
369 CHKERR getInterface<BitRefManager>()->setEntitiesBitRefLevel(bit_ref_ents);
370 CHKERR getInterface<BitRefManager>()->setElementsBitRefLevel(
371 ref_elems_to_add);
372
373 // Build field entities
374 for (auto field : fIelds) {
375 if (field->getSpace() != NOSPACE) {
376 Range ents_of_id_meshset;
377 CHKERR get_moab().get_entities_by_handle(field->getMeshset(),
378 ents_of_id_meshset, false);
379 CHKERR set_field_order(ents_of_id_meshset, field->getId(), -1, verb);
380 }
381 }
382
383 if (initaliseAndBuildField || initaliseAndBuildFiniteElements) {
384 CHKERR build_fields(verb);
385 if (initaliseAndBuildFiniteElements) {
386 CHKERR build_finite_elements(verb);
387 }
388 }
389
390 if (verb > VERY_NOISY) {
391 list_fields();
392 list_finite_elements();
393 list_problem();
394 }
395
396 // Initialize interfaces
397 MeshsetsManager *m_manger_ptr;
398 CHKERR getInterface(m_manger_ptr);
399 CHKERR m_manger_ptr->initialiseDatabaseFromMesh(verb);
400 SeriesRecorder *series_recorder_ptr;
401 CHKERR getInterface(series_recorder_ptr);
402 CHKERR series_recorder_ptr->initialiseDatabaseFromMesh(verb);
403
405}
406
407MoFEMErrorCode Core::setMoabInterface(moab::Interface &new_moab, int verb) {
409 if (verb == -1)
410 verb = verbose;
411
412 // clear moab database
413 CHKERR clearMap();
414
415 // set new reference
416 moab = std::ref(new_moab);
417
418 // check if moab has set communicator if not set communicator internally
419 ParallelComm *pComm = ParallelComm::get_pcomm(&new_moab, MYPCOMM_INDEX);
420 if (pComm == NULL) {
421 pComm = new ParallelComm(&new_moab, wrapMPIMOABComm->get_comm());
422 }
423
424 // create MoFEM tags
425 CHKERR getTags();
426
427 // Create basic entity data struture
428 basicEntityDataPtr = boost::make_shared<BasicEntityData>(moab);
429 setRefEntBasicDataPtr(*this, this->basicEntityDataPtr);
430
431 // Initalise database
432 CHKERR this->initialiseDatabaseFromMesh(verb);
433
435};
436
439
440 iFaces.clear();
441
442 // Register sub interfaces
443 CHKERR regSubInterface<LogManager>();
444 CHKERR regSubInterface<Simple>();
445 CHKERR regSubInterface<OperatorsTester>();
446 CHKERR regSubInterface<PipelineManager>();
447 CHKERR regSubInterface<ProblemsManager>();
448 CHKERR regSubInterface<MatrixManager>();
449 CHKERR regSubInterface<ISManager>();
450 CHKERR regSubInterface<VecManager>();
451 CHKERR regSubInterface<FieldBlas>();
452 CHKERR regSubInterface<BitRefManager>();
453 CHKERR regSubInterface<Tools>();
454 CHKERR regSubInterface<CommInterface>();
455 CHKERR regSubInterface<MeshsetsManager>();
456 CHKERR regSubInterface<NodeMergerInterface>();
457 CHKERR regSubInterface<PrismsFromSurfaceInterface>();
458 CHKERR regSubInterface<MeshRefinement>();
459 CHKERR regSubInterface<PrismInterface>();
460 CHKERR regSubInterface<CutMeshInterface>();
461 CHKERR regSubInterface<SeriesRecorder>();
462#ifdef WITH_TETGEN
463 CHKERR regSubInterface<TetGenInterface>();
464#endif
465#ifdef WITH_MED
466 CHKERR regSubInterface<MedInterface>();
467#endif
468 CHKERR regSubInterface<FieldEvaluatorInterface>();
469 CHKERR regSubInterface<BcManager>();
470
472};
473
476 // Cleaning databases in interfaces
477 CHKERR getInterface<SeriesRecorder>()->clearMap();
478 CHKERR getInterface<MeshsetsManager>()->clearMap();
479 CHKERR getInterface<CutMeshInterface>()->clearMap();
480 // Cleaning databases
481 refinedEntities.clear();
482 refinedFiniteElements.clear();
483 fIelds.clear();
484 entsFields.clear();
485 dofsField.clear();
486 finiteElements.clear();
487 entsFiniteElements.clear();
488 entFEAdjacencies.clear();
489 pRoblems.clear();
491}
492
495 if (verb == -1)
496 verb = verbose;
497 std::pair<RefEntity_multiIndex::iterator, bool> p_ent;
498 p_ent = refinedEntities.insert(
499 boost::make_shared<RefEntity>(basicEntityDataPtr, prism));
500 if (p_ent.second) {
501 std::pair<RefElement_multiIndex::iterator, bool> p;
502 p = refinedFiniteElements.insert(
503 boost::shared_ptr<RefElement>(new RefElement_PRISM(*p_ent.first)));
504 int num_nodes;
505 const EntityHandle *conn;
506 CHKERR get_moab().get_connectivity(prism, conn, num_nodes, true);
507 Range face_side3, face_side4;
508 CHKERR get_moab().get_adjacencies(conn, 3, 2, false, face_side3);
509 CHKERR get_moab().get_adjacencies(&conn[3], 3, 2, false, face_side4);
510 if (face_side3.size() != 1)
511 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
512 "prism don't have side face 3");
513 if (face_side4.size() != 1)
514 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
515 "prims don't have side face 4");
516 p.first->get()->getSideNumberPtr(*face_side3.begin());
517 p.first->get()->getSideNumberPtr(*face_side4.begin());
518 }
520}
521
524
525 const EntityHandle root_meshset = get_moab().get_root_set();
526 if (root_meshset) {
527 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
528 "Root meshset should be 0");
529 }
530
531 // Set version
532 {
533 Version version;
534 CHKERR getFileVersion(moab, version);
535 }
536
537 // Global Variables
538 {
539
540 auto check_tag_allocated = [](auto &rval) {
542 if (rval == MB_ALREADY_ALLOCATED)
543 rval = MB_SUCCESS;
544 else
545 CHKERRG(rval);
547 };
548
549 // Safety nets
550 int def_bool = 0;
551 rval = get_moab().tag_get_handle("_MoFEMBuild", 1, MB_TYPE_INTEGER,
552 th_MoFEMBuild, MB_TAG_CREAT | MB_TAG_MESH,
553 &def_bool);
554 CHKERR check_tag_allocated(rval);
555
556 CHKERR get_moab().tag_get_by_ptr(th_MoFEMBuild, &root_meshset, 1,
557 (const void **)&buildMoFEM);
558 }
559
560 // Tags saved in vtk-files
561 {
562 const int def_part = -1;
563 CHKERR get_moab().tag_get_handle("PARTITION", 1, MB_TYPE_INTEGER, th_Part,
564 MB_TAG_CREAT | MB_TAG_SPARSE, &def_part);
565 }
566
567 // Tags Ref
568 {
569
570 // Fix size of bir ref level tags
572
573 const int def_part = -1;
574 CHKERR get_moab().tag_get_handle("_MeshsetPartition", 1, MB_TYPE_INTEGER,
575 th_Part, MB_TAG_CREAT | MB_TAG_SPARSE,
576 &def_part);
577 EntityHandle def_handle = 0;
578 CHKERR get_moab().tag_get_handle("_RefParentHandle", 1, MB_TYPE_HANDLE,
579 th_RefParentHandle,
580 MB_TAG_CREAT | MB_TAG_SPARSE, &def_handle);
581 BitRefLevel def_bit_level = 0;
582 CHKERR get_moab().tag_get_handle(
583 "_RefBitLevel", sizeof(BitRefLevel), MB_TYPE_OPAQUE, th_RefBitLevel,
584 MB_TAG_CREAT | MB_TAG_BYTES | MB_TAG_SPARSE, &def_bit_level);
585 BitRefLevel def_bit_level_mask = BitRefLevel().set();
586 CHKERR get_moab().tag_get_handle(
587 "_RefBitLevelMask", sizeof(BitRefLevel), MB_TYPE_OPAQUE,
588 th_RefBitLevel_Mask, MB_TAG_CREAT | MB_TAG_BYTES | MB_TAG_SPARSE,
589 &def_bit_level_mask);
590 BitRefEdges def_bit_edge = 0;
591 CHKERR get_moab().tag_get_handle(
592 "_RefBitEdge", sizeof(BitRefEdges), MB_TYPE_OPAQUE, th_RefBitEdge,
593 MB_TAG_CREAT | MB_TAG_SPARSE | MB_TAG_BYTES, &def_bit_edge);
594 }
595
596 // Tags Field
597 {
598 const unsigned long int def_id = 0;
599 CHKERR get_moab().tag_get_handle(
600 "_FieldId", sizeof(BitFieldId), MB_TYPE_OPAQUE, th_FieldId,
601 MB_TAG_CREAT | MB_TAG_BYTES | MB_TAG_SPARSE, &def_id);
602 FieldSpace def_space = LASTSPACE;
603 CHKERR get_moab().tag_get_handle(
604 "_FieldSpace", sizeof(FieldSpace), MB_TYPE_OPAQUE, th_FieldSpace,
605 MB_TAG_CREAT | MB_TAG_BYTES | MB_TAG_SPARSE, &def_space);
607 CHKERR get_moab().tag_get_handle(
608 "_FieldBase", sizeof(FieldApproximationBase), MB_TYPE_OPAQUE,
609 th_FieldBase, MB_TAG_CREAT | MB_TAG_BYTES | MB_TAG_SPARSE, &def_base);
610 const int def_val_len = 0;
611 CHKERR get_moab().tag_get_handle(
612 "_FieldName", def_val_len, MB_TYPE_OPAQUE, th_FieldName,
613 MB_TAG_CREAT | MB_TAG_BYTES | MB_TAG_VARLEN | MB_TAG_SPARSE, NULL);
614 CHKERR get_moab().tag_get_handle(
615 "_FieldName_DataNamePrefix", def_val_len, MB_TYPE_OPAQUE,
616 th_FieldName_DataNamePrefix,
617 MB_TAG_CREAT | MB_TAG_BYTES | MB_TAG_VARLEN | MB_TAG_SPARSE, NULL);
618 }
619
620 // Tags FE
621 {
622 const unsigned long int def_id = 0;
623 const int def_val_len = 0;
624 CHKERR get_moab().tag_get_handle(
625 "_FEId", sizeof(BitFEId), MB_TYPE_OPAQUE, th_FEId,
626 MB_TAG_CREAT | MB_TAG_BYTES | MB_TAG_SPARSE, &def_id);
627 CHKERR get_moab().tag_get_handle(
628 "_FEName", def_val_len, MB_TYPE_OPAQUE, th_FEName,
629 MB_TAG_CREAT | MB_TAG_BYTES | MB_TAG_VARLEN | MB_TAG_SPARSE, NULL);
630 CHKERR get_moab().tag_get_handle(
631 "_FEIdCol", sizeof(BitFieldId), MB_TYPE_OPAQUE, th_FEIdCol,
632 MB_TAG_CREAT | MB_TAG_BYTES | MB_TAG_SPARSE, &def_id);
633 CHKERR get_moab().tag_get_handle(
634 "_FEIdRow", sizeof(BitFieldId), MB_TYPE_OPAQUE, th_FEIdRow,
635 MB_TAG_CREAT | MB_TAG_BYTES | MB_TAG_SPARSE, &def_id);
636 CHKERR get_moab().tag_get_handle(
637 "_FEIdData", sizeof(BitFieldId), MB_TYPE_OPAQUE, th_FEIdData,
638 MB_TAG_CREAT | MB_TAG_BYTES | MB_TAG_SPARSE, &def_id);
639 }
640
641 // Tags Problem
642 {
643 const unsigned long int def_id = 0;
644 const int def_val_len = 0;
645 CHKERR get_moab().tag_get_handle(
646 "_ProblemId", sizeof(BitProblemId), MB_TYPE_OPAQUE, th_ProblemId,
647 MB_TAG_CREAT | MB_TAG_BYTES | MB_TAG_SPARSE, &def_id);
648 CHKERR get_moab().tag_get_handle(
649 "_ProblemFEId", sizeof(BitFEId), MB_TYPE_OPAQUE, th_ProblemFEId,
650 MB_TAG_CREAT | MB_TAG_BYTES | MB_TAG_SPARSE, &def_id);
651 CHKERR get_moab().tag_get_handle(
652 "_ProblemName", def_val_len, MB_TYPE_OPAQUE, th_ProblemName,
653 MB_TAG_CREAT | MB_TAG_BYTES | MB_TAG_VARLEN | MB_TAG_SPARSE, NULL);
654 DofIdx def_nbdofs = 0;
655 CHKERR get_moab().tag_get_handle(
656 "_ProblemNbDofsRow", sizeof(DofIdx), MB_TYPE_OPAQUE,
657 th_ProblemNbDofsRow, MB_TAG_CREAT | MB_TAG_BYTES | MB_TAG_SPARSE,
658 &def_nbdofs);
659 CHKERR get_moab().tag_get_handle(
660 "_ProblemNbDofsCol", sizeof(DofIdx), MB_TYPE_OPAQUE,
661 th_ProblemNbDofsCol, MB_TAG_CREAT | MB_TAG_BYTES | MB_TAG_SPARSE,
662 &def_nbdofs);
663 CHKERR get_moab().tag_get_handle(
664 "_ProblemLocalNbDofsRow", sizeof(DofIdx), MB_TYPE_OPAQUE,
665 th_ProblemLocalNbDofRow, MB_TAG_CREAT | MB_TAG_BYTES | MB_TAG_SPARSE,
666 &def_nbdofs);
667 CHKERR get_moab().tag_get_handle(
668 "_ProblemGhostNbDofsRow", sizeof(DofIdx), MB_TYPE_OPAQUE,
669 th_ProblemGhostNbDofRow, MB_TAG_CREAT | MB_TAG_BYTES | MB_TAG_SPARSE,
670 &def_nbdofs);
671 CHKERR get_moab().tag_get_handle(
672 "_ProblemLocalNbDofsCol", sizeof(DofIdx), MB_TYPE_OPAQUE,
673 th_ProblemLocalNbDofCol, MB_TAG_CREAT | MB_TAG_BYTES | MB_TAG_SPARSE,
674 &def_nbdofs);
675 CHKERR get_moab().tag_get_handle(
676 "_ProblemGhostNbDofsCol", sizeof(DofIdx), MB_TYPE_OPAQUE,
677 th_ProblemGhostNbDofCol, MB_TAG_CREAT | MB_TAG_BYTES | MB_TAG_SPARSE,
678 &def_nbdofs);
679 }
680
681 // Meshsets with boundary conditions and material sets
682 MeshsetsManager *meshsets_manager_ptr;
683 CHKERR getInterface(meshsets_manager_ptr);
684 CHKERR meshsets_manager_ptr->getTags(verb);
685
686 // Series recorder
687 SeriesRecorder *series_recorder_ptr;
688 CHKERR getInterface(series_recorder_ptr);
689 CHKERR series_recorder_ptr->getTags(verb);
690
692}
693
696 if (verb == -1)
697 verb = verbose;
698 CHKERR clearMap();
700}
701
704 if (verb == -1)
705 verb = verbose;
706 CHKERR this->clearMap();
707 CHKERR this->getTags(verb);
708 CHKERR this->initialiseDatabaseFromMesh(verb);
710}
711
712MoFEMErrorCode Core::set_moab_interface(moab::Interface &new_moab, int verb) {
713 return this->setMoabInterface(new_moab, verb);
714};
715
717 int verb) {
718 return this->setMoabInterface(new_moab, verb);
719};
720
723 if (verb == -1)
724 verb = verbose;
725
726 CHKERR PetscOptionsBegin(mofemComm, optionsPrefix.c_str(), "Mesh cut options",
727 "See MoFEM documentation");
728
729 CHKERR PetscOptionsBool(
730 "-mofem_init_fields", "Initialise fields on construction", "",
731 initaliseAndBuildField, &initaliseAndBuildField, NULL);
732
733 CHKERR PetscOptionsBool(
734 "-mofem_init_fields", "Initialise fields on construction", "",
735 initaliseAndBuildFiniteElements, &initaliseAndBuildFiniteElements, NULL);
736
737 // TODO: Add read verbosity level
738 // TODO: Add option to initalise problems ??? - DO WE REALLY NEED THAT
739
740 ierr = PetscOptionsEnd();
741 CHKERRG(ierr);
742
744}
745
746// cubit meshsets
747
750 *fields_ptr = &fIelds;
752}
753
755Core::get_ref_ents(const RefEntity_multiIndex **refined_entities_ptr) const {
757 *refined_entities_ptr = &refinedEntities;
759}
761 const RefElement_multiIndex **refined_finite_elements_ptr) const {
763 *refined_finite_elements_ptr = &refinedFiniteElements;
765}
766
767MoFEMErrorCode Core::get_problem(const std::string &problem_name,
768 const Problem **problem_ptr) const {
770 typedef Problem_multiIndex::index<Problem_mi_tag>::type ProblemsByName;
771 const ProblemsByName &problems = pRoblems.get<Problem_mi_tag>();
772 ProblemsByName::iterator p_miit = problems.find(problem_name);
773 if (p_miit == problems.end()) {
774 SETERRQ1(PETSC_COMM_SELF, MOFEM_OPERATION_UNSUCCESSFUL,
775 "problem < %s > not found, (top tip: check spelling)",
776 problem_name.c_str());
777 }
778 *problem_ptr = &*p_miit;
780}
781
783Core::get_problems(const Problem_multiIndex **problems_ptr) const {
785 *problems_ptr = &pRoblems;
787}
788
792 *field_ents = &entsFields;
794}
797 *dofs_ptr = &dofsField;
799}
800
804 *fe_ptr = &finiteElements;
806}
807
809 const EntFiniteElement_multiIndex **fe_ent_ptr) const {
811 *fe_ent_ptr = &entsFiniteElements;
813}
814
816 MeshsetsManager *meshsets_manager_ptr;
817 getInterface(meshsets_manager_ptr);
818 return meshsets_manager_ptr;
819}
820
822 MeshsetsManager *meshsets_manager_ptr;
823 getInterface(meshsets_manager_ptr);
824 return meshsets_manager_ptr;
825}
826
829 *dofs_elements_adjacency) const {
831 *dofs_elements_adjacency = &entFEAdjacencies;
833}
834
837 return &entFEAdjacencies;
838}
839
840const Field_multiIndex *Core::get_fields() const { return &fIelds; }
842 return &refinedEntities;
843}
845 return &refinedFiniteElements;
846}
848 return &finiteElements;
849}
851 return &entsFiniteElements;
852}
854 return &entsFields;
855}
856const DofEntity_multiIndex *Core::get_dofs() const { return &dofsField; }
857const Problem *Core::get_problem(const std::string problem_name) const {
858 const Problem *prb;
859 CHK_THROW_MESSAGE(get_problem(problem_name, &prb),
860 "Problem of given name not found");
861 return prb;
862}
863const Problem_multiIndex *Core::get_problems() const { return &pRoblems; }
864
865template <int V, typename std::enable_if<(V >= 0), int>::type * = nullptr>
866void set_ref_ent_basic_data_ptr_impl(boost::shared_ptr<BasicEntityData> &ptr) {
868};
869
870template <int V, typename std::enable_if<(V < 0), int>::type * = nullptr>
871void set_ref_ent_basic_data_ptr_impl(boost::shared_ptr<BasicEntityData> &ptr) {
872 return;
873};
874
875void Core::setRefEntBasicDataPtr(MoFEM::Interface &m_field,
876 boost::shared_ptr<BasicEntityData> &ptr) {
877
878 switch (m_field.getValue()) {
879 case -1:
880 set_ref_ent_basic_data_ptr_impl<-1>(ptr);
881 break;
882 case 0:
883 set_ref_ent_basic_data_ptr_impl<0>(ptr);
884 break;
885 case 1:
886 set_ref_ent_basic_data_ptr_impl<1>(ptr);
887 break;
888 default:
889 THROW_MESSAGE("Core index can vary from -1 to MAX_CORE_TMP");
890 }
891};
892
893boost::shared_ptr<RefEntityTmp<0>>
894Core::makeSharedRefEntity(MoFEM::Interface &m_field, const EntityHandle ent) {
895
896 boost::shared_ptr<RefEntityTmp<0>> ref_ent_ptr;
897
898 switch (m_field.getValue()) {
899 case -1:
900 ref_ent_ptr = boost::shared_ptr<RefEntityTmp<0>>(
901
902 new RefEntityTmp<-1>(m_field.get_basic_entity_data_ptr(), ent)
903
904 );
905 break;
906 case 0:
907 ref_ent_ptr = boost::shared_ptr<RefEntityTmp<0>>(
908
909 new RefEntityTmp<0>(m_field.get_basic_entity_data_ptr(), ent)
910
911 );
912 break;
913 case 1:
914 ref_ent_ptr = boost::shared_ptr<RefEntityTmp<0>>(
915
916 new RefEntityTmp<1>(m_field.get_basic_entity_data_ptr(), ent)
917
918 );
919 break;
920 default:
921 THROW_MESSAGE("Core index can vary from -1 to MAX_CORE_TMP");
922 }
923
924 return ref_ent_ptr;
925}
926
927boost::shared_ptr<RefEntityTmp<0>>
928Core::make_shared_ref_entity(const EntityHandle ent) {
929 return this->makeSharedRefEntity(*this, ent);
930}
931
932boost::shared_ptr<RefEntityTmp<0>>
934 return this->makeSharedRefEntity(*this, ent);
935}
936
937} // 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.
static Index< 'p', 3 > p
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
#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
static auto get_sub_iface_options_imp(T *const ptr, int) -> decltype(ptr->getSubInterfaceOptions())
Definition Core.cpp:138
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
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
void set_ref_ent_basic_data_ptr_impl(boost::shared_ptr< BasicEntityData > &ptr)
Definition Core.cpp:866
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 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:474
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 subinterfac in core interface.
Definition Core.cpp:150
static bool isGloballyInitialised
Core base globally initialized.
Definition Core.hpp:1010
const RefElement_multiIndex * get_ref_finite_elements() const
Get the ref finite elements object.
Definition Core.cpp:844
const FieldEntity_multiIndex * get_field_ents() const
Get the field ents object.
Definition Core.cpp:853
MoFEMErrorCode initialiseDatabaseFromMesh(int verb=DEFAULT_VERBOSITY)
Initialize database getting information on mesh.
Definition Core.cpp:268
MoFEMErrorCode registerSubInterfaces()
Register insterfaces.
Definition Core.cpp:437
const Field_multiIndex * get_fields() const
Get the fields object.
Definition Core.cpp:840
const Problem_multiIndex * get_problems() const
Get the problems object.
Definition Core.cpp:863
MoFEMErrorCode coreGenericConstructor(moab::Interface &moab, MPI_Comm comm, const int verbose)
Definition Core.cpp:171
const DofEntity_multiIndex * get_dofs() const
Get the dofs object.
Definition Core.cpp:856
MoFEMErrorCode getOptions(int verb=DEFAULT_VERBOSITY)
Get core options from command line.
Definition Core.cpp:721
static PetscBool isInitialized
petsc was initialised by other agent
Definition Core.hpp:1012
const RefEntity_multiIndex * get_ref_ents() const
Get the ref ents object.
Definition Core.cpp:841
const FieldEntityEntFiniteElementAdjacencyMap_multiIndex * get_ents_elements_adjacency() const
Get the dofs elements adjacency object.
Definition Core.cpp:836
const EntFiniteElement_multiIndex * get_ents_finite_elements() const
Get the ents finite elements object.
Definition Core.cpp:850
const FiniteElement_multiIndex * get_finite_elements() const
Get the finite elements object.
Definition Core.cpp:847
MoFEMErrorCode addPrismToDatabase(const EntityHandle prism, int verb=DEFAULT_VERBOSITY)
add prim element
Definition Core.cpp:493
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:712
MoFEMErrorCode setMoabInterface(moab::Interface &new_moab, int verb=VERBOSE)
Definition Core.cpp:407
MoFEMErrorCode get_problem(const std::string &problem_name, const Problem **problem_ptr) const
Get problem database (data structure)
Definition Core.cpp:767
MoFEMErrorCode getTags(int verb=DEFAULT_VERBOSITY)
Get tag handles.
Definition Core.cpp:522
CoreTmp(moab::Interface &moab, MPI_Comm comm=PETSC_COMM_WORLD, const int verbose=VERBOSE)
Definition Core.cpp:209
static int mpiInitialised
mpi was initialised by other agent
Definition Core.hpp:1011
MeshsetsManager * get_meshsets_manager_ptr()
get MeshsetsManager pointer
Definition Core.cpp:815
MoFEMErrorCode rebuild_database(int verb=DEFAULT_VERBOSITY)
Clear database and initialize it once again.
Definition Core.cpp:702
MoFEMErrorCode clear_database(int verb=DEFAULT_VERBOSITY)
Clear database.
Definition Core.cpp:694
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 initialiseDatabaseFromMesh(int verb=DEFAULT_VERBOSITY)
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
MoFEMErrorCode initialiseDatabaseFromMesh(int verb=0)
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