v0.15.0
Loading...
Searching...
No Matches
DMMoFEM.cpp
Go to the documentation of this file.
1
2
3// #undef PETSC_VERSION_RELEASE
4// #define PETSC_VERSION_RELEASE 1
5
6#if PETSC_VERSION_GE(3, 6, 0)
7#include <petsc/private/dmimpl.h> /*I "petscdm.h" I*/
8// #include <petsc/private/vecimpl.h> /*I "petscdm.h" I*/
9#else
10#include <petsc-private/dmimpl.h> /*I "petscdm.h" I*/
11#include <petsc-private/vecimpl.h> /*I "petscdm.h" I*/
12#endif
13
14#include <DMMoFEM.hpp>
15#include <DMCtxImpl.hpp>
16
17namespace MoFEM {
18
19MoFEMErrorCode DMCtx::query_interface(boost::typeindex::type_index type_index,
20 UnknownInterface **iface) const {
21 *iface = const_cast<DMCtx *>(this);
22 return 0;
23}
24
26 if (!LogManager::checkIfChannelExist("DMWORLD")) {
27 auto core_log = logging::core::get();
28 core_log->add_sink(
30 core_log->add_sink(
32 core_log->add_sink(
34 LogManager::setLog("DMWORLD");
35 LogManager::setLog("DMSYNC");
36 LogManager::setLog("DMSELF");
37 MOFEM_LOG_TAG("DMWORLD", "DM");
38 MOFEM_LOG_TAG("DMSYNC", "DM");
39 MOFEM_LOG_TAG("DMSELF", "DM");
40 }
41}
42
43PetscErrorCode DMRegister_MoFEM(const char sname[]) {
45 CHKERR DMRegister(sname, DMCreate_MoFEM);
47}
48
49PetscErrorCode DMSetOperators_MoFEM(DM dm) {
50 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
52
53 dm->ops->createglobalvector = DMCreateGlobalVector_MoFEM;
54 dm->ops->createlocalvector = DMCreateLocalVector_MoFEM;
55 dm->ops->creatematrix = DMCreateMatrix_MoFEM;
56 dm->ops->setup = DMSetUp_MoFEM;
57 dm->ops->destroy = DMDestroy_MoFEM;
58 dm->ops->setfromoptions = DMSetFromOptions_MoFEM;
59 dm->ops->globaltolocalbegin = DMGlobalToLocalBegin_MoFEM;
60 dm->ops->globaltolocalend = DMGlobalToLocalEnd_MoFEM;
61 dm->ops->localtoglobalbegin = DMLocalToGlobalBegin_MoFEM;
62 dm->ops->localtoglobalend = DMLocalToGlobalEnd_MoFEM;
63 dm->ops->createfieldis = DMCreateFieldIS_MoFEM;
64
65 // Default matrix type
66 CHKERR DMSetMatType(dm, MATMPIAIJ);
67
69}
70
71PetscErrorCode DMCreate_MoFEM(DM dm) {
72 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
74 dm->data = new DMCtxImpl();
77}
78
79PetscErrorCode DMDestroy_MoFEM(DM dm) {
80 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
81 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
83
84 MPI_Comm comm;
85 CHKERR PetscObjectGetComm((PetscObject)dm, &comm);
86
87 int result;
88 MPI_Comm_compare(comm, PETSC_COMM_SELF, &result);
89 if (result == MPI_IDENT)
90 MOFEM_LOG("DMSELF", Sev::noisy)
91 << "MoFEM DM destroy for problem " << dm_field->problemName
92 << " referenceNumber " << dm_field->referenceNumber;
93 else
94 MOFEM_LOG("DMWORLD", Sev::noisy)
95 << "MoFEM DM destroy for problem " << dm_field->problemName
96 << " referenceNumber " << dm_field->referenceNumber;
97
98 if (dm_field->referenceNumber == 0) {
99 if (dm_field->destroyProblem) {
100
101 if (dm_field->mField_ptr->check_problem(dm_field->problemName)) {
102 dm_field->mField_ptr->delete_problem(dm_field->problemName);
103 } // else problem has to be deleted by the user
104 }
105
106 delete static_cast<DMCtxImpl *>(dm->data);
107
108 } else
109 --dm_field->referenceNumber;
110
112}
113
114PetscErrorCode DMMoFEMCreateMoFEM(DM dm, MoFEM::Interface *m_field_ptr,
115 const char problem_name[],
116 const MoFEM::BitRefLevel bit_level,
117 const MoFEM::BitRefLevel bit_mask) {
119
120 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
121 if (!dm->data) {
122 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
123 "data structure for MoFEM not yet created");
124 }
125 if (!m_field_ptr) {
126 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
127 "DM function not implemented into MoFEM");
128 }
129 dm_field->mField_ptr = m_field_ptr;
130 dm_field->problemName = problem_name;
131 if (!m_field_ptr->check_problem(dm_field->problemName)) {
132 // problem is not defined, declare problem internally set bool to
133 // destroyProblem problem with DM
134 dm_field->destroyProblem = PETSC_TRUE;
135 CHKERR dm_field->mField_ptr->add_problem(dm_field->problemName, MF_EXCL,
136 dm_field->verbosity);
137 } else {
138 dm_field->destroyProblem = PETSC_FALSE;
139 }
141 dm_field->problemName, bit_level);
143 dm_field->problemName, bit_mask);
144 dm_field->kspCtx =
145 boost::shared_ptr<KspCtx>(new KspCtx(*m_field_ptr, problem_name));
146 dm_field->snesCtx =
147 boost::shared_ptr<SnesCtx>(new SnesCtx(*m_field_ptr, problem_name));
148 dm_field->tsCtx =
149 boost::shared_ptr<TsCtx>(new TsCtx(*m_field_ptr, problem_name));
150
151 MPI_Comm comm;
152 CHKERR PetscObjectGetComm((PetscObject)dm, &comm);
153 int result = 0;
154 MPI_Comm_compare(comm, m_field_ptr->get_comm(), &result);
155 if (result > MPI_CONGRUENT) {
156 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
157 "MoFEM and DM using different communicators");
158 }
159 MPI_Comm_size(comm, &dm_field->sIze);
160 MPI_Comm_rank(comm, &dm_field->rAnk);
161
162 // problem structure
163 CHKERR dm_field->mField_ptr->get_problem(dm_field->problemName,
164 &dm_field->problemPtr);
165
166 MPI_Comm_compare(comm, PETSC_COMM_SELF, &result);
167 if (result == MPI_IDENT) {
168 MOFEM_LOG("DMSELF", Sev::noisy)
169 << "MoFEM DM created for problem " << dm_field->problemName;
170 MOFEM_LOG("DMSELF", Sev::noisy) << *dm_field->problemPtr;
171 } else {
172 MOFEM_LOG("DMWORLD", Sev::noisy)
173 << "MoFEM DM created for problem " << dm_field->problemName;
174 MOFEM_LOG("DMWORLD", Sev::noisy) << *dm_field->problemPtr;
175 }
176
178}
179
180PetscErrorCode DMMoFEMDuplicateDMCtx(DM dm, DM dm_duplicate) {
182
183 if (!dm->data)
184 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
185 "data structure for MoFEM not yet created");
186
187 if (static_cast<DMCtxImpl *>(dm_duplicate->data)->referenceNumber == 0)
188 delete static_cast<DMCtxImpl *>(dm_duplicate->data);
189
190 dm_duplicate->data = dm->data;
191 ++(static_cast<DMCtxImpl *>(dm_duplicate->data)->referenceNumber);
192
194}
195
196PetscErrorCode DMMoFEMSwapDMCtx(DM dm, DM dm_swap) {
198 if (!dm->data)
199 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
200 "data structure for MoFEM not yet created on dm");
201 if (!dm_swap->data)
202 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
203 "data structure for MoFEM not yet created on swap dm");
204
205 auto *dm_field = static_cast<DMCtxImpl *>(dm->data);
206 auto *dm_field_swap = static_cast<DMCtxImpl *>(dm_swap->data);
207
208 auto tmp_field = dm_field;
209 dm_field = dm_field_swap;
210 dm_field_swap = tmp_field;
211
213}
214
215PetscErrorCode DMMoFEMCreateSubDM(DM subdm, DM dm, const char problem_name[]) {
217
218 if (!dm->data) {
219 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
220 "data structure for MoFEM not yet created");
221 }
222 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
223
224 CHKERR DMMoFEMCreateMoFEM(subdm, dm_field->mField_ptr, problem_name,
225 dm_field->problemPtr->getBitRefLevel(),
226 dm_field->problemPtr->getBitRefLevelMask());
227
228 DMCtxImpl *subdm_field = (DMCtxImpl *)subdm->data;
229 subdm_field->isSubDM = PETSC_TRUE;
230 subdm_field->problemMainOfSubPtr = dm_field->problemPtr;
231 subdm_field->isPartitioned = dm_field->isPartitioned;
232 subdm_field->isSquareMatrix = PETSC_FALSE;
233 subdm->ops->setup = DMSubDMSetUp_MoFEM;
234
236}
237
238PetscErrorCode DMMoFEMAddSubFieldRow(DM dm, const char field_name[]) {
239 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
241 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
242 if (!dm->data) {
243 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
244 "data structure for MoFEM not yet created");
245 }
246 if (!dm_field->isSubDM) {
247 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "this is not sub-dm");
248 }
249 dm_field->rowSubFields.push_back(field_name);
250 dm_field->mapTypeRow.erase(field_name);
252}
253
254PetscErrorCode DMMoFEMAddSubFieldRow(DM dm, std::string field_name) {
255 return DMMoFEMAddSubFieldRow(dm, field_name.c_str());
256}
257
258PetscErrorCode DMMoFEMAddSubFieldRow(DM dm, const char field_name[],
259 boost::shared_ptr<Range> r_ptr) {
260 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
262 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
263 if (!dm->data) {
264 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
265 "data structure for MoFEM not yet created");
266 }
267 if (!dm_field->isSubDM) {
268 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "this is not sub-dm");
269 }
270 dm_field->rowSubFields.push_back(field_name);
271 dm_field->mapTypeRow[field_name] = r_ptr;
273}
274
275PetscErrorCode DMMoFEMAddSubFieldRow(DM dm, std::string field_name,
276 boost::shared_ptr<Range> r_ptr) {
277 return DMMoFEMAddSubFieldRow(dm, field_name.c_str(), r_ptr);
278}
279
280PetscErrorCode DMMoFEMAddSubFieldCol(DM dm, const char field_name[]) {
281 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
283 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
284 if (!dm->data) {
285 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
286 "data structure for MoFEM not yet created");
287 }
288 if (!dm_field->isSubDM) {
289 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "this is not sub-dm");
290 }
291 dm_field->colSubFields.push_back(field_name);
292 dm_field->mapTypeCol.erase(field_name);
294}
295
296PetscErrorCode DMMoFEMAddSubFieldCol(DM dm, std::string field_name) {
297 return DMMoFEMAddSubFieldCol(dm, field_name.c_str());
298}
299
300PetscErrorCode DMMoFEMAddSubFieldCol(DM dm, const char field_name[],
301 boost::shared_ptr<Range> r_ptr) {
302 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
304 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
305 if (!dm->data) {
306 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
307 "data structure for MoFEM not yet created");
308 }
309 if (!dm_field->isSubDM) {
310 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "this is not sub-dm");
311 }
312 dm_field->colSubFields.push_back(field_name);
313 dm_field->mapTypeCol[field_name] = r_ptr;
315}
316
317PetscErrorCode DMMoFEMAddSubFieldCol(DM dm, std::string field_name,
318 boost::shared_ptr<Range> r_ptr) {
319 return DMMoFEMAddSubFieldCol(dm, field_name.c_str(), r_ptr);
320}
321
322PetscErrorCode DMMoFEMGetIsSubDM(DM dm, PetscBool *is_sub_dm) {
323 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
325 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
326 *is_sub_dm = dm_field->isSubDM;
328}
329
330PetscErrorCode DMMoFEMGetSubRowIS(DM dm, IS *is) {
331 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
333 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
334 if (dm_field->isSubDM != PETSC_TRUE) {
335 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
336 "This DM is not created as a SubDM");
337 }
338 if (dm_field->isProblemBuild != PETSC_TRUE) {
339 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "Problem is not build");
340 }
341 boost::shared_ptr<Problem::SubProblemData> sub_data =
342 dm_field->problemPtr->getSubData();
343 CHKERR sub_data->getRowIs(is);
345}
346
347PetscErrorCode DMMoFEMGetSubColIS(DM dm, IS *is) {
348 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
350 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
351 if (dm_field->isSubDM != PETSC_TRUE) {
352 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
353 "This DM is not created as a SubDM");
354 }
355 if (dm_field->isProblemBuild != PETSC_TRUE) {
356 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "Problem is not build");
357 }
358 boost::shared_ptr<Problem::SubProblemData> sub_data =
359 dm_field->problemPtr->getSubData();
360 CHKERR sub_data->getColIs(is);
362}
363
364PetscErrorCode DMMoFEMAddRowCompositeProblem(DM dm, const char prb_name[]) {
365 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
367 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
368 if (!dm->data) {
369 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
370 "data structure for MoFEM not yet created");
371 }
372 if (!dm_field->isCompDM) {
373 dm_field->isCompDM = PETSC_TRUE;
374 }
375 dm_field->rowCompPrb.push_back(prb_name);
376 if (dm_field->isSquareMatrix) {
377 dm_field->colCompPrb.push_back(prb_name);
378 }
380}
381
382PetscErrorCode DMMoFEMAddColCompositeProblem(DM dm, const char prb_name[]) {
383 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
385 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
386 if (!dm->data) {
387 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
388 "data structure for MoFEM not yet created");
389 }
390 if (!dm_field->isCompDM) {
391 dm_field->isCompDM = PETSC_TRUE;
392 }
393 if (dm_field->isSquareMatrix) {
394 SETERRQ(PETSC_COMM_SELF, MOFEM_INVALID_DATA,
395 "No need to add problem on column when problem block structurally "
396 "symmetric");
397 }
398 dm_field->colCompPrb.push_back(prb_name);
400}
401
402PetscErrorCode DMMoFEMGetIsCompDM(DM dm, PetscBool *is_comp_dm) {
403 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
405 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
406 *is_comp_dm = dm_field->isCompDM;
408}
409
410PetscErrorCode DMoFEMGetInterfacePtr(DM dm, MoFEM::Interface **m_field_ptr) {
411 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
413 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
414 if (!dm->data) {
415 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
416 "data structure for MoFEM not yet created");
417 }
418 *m_field_ptr = dm_field->mField_ptr;
420}
421
422PetscErrorCode DMMoFEMGetProblemPtr(DM dm, const MoFEM::Problem **problem_ptr) {
423 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
425 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
426 if (!dm->data) {
427 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
428 "data structure for MoFEM not yet created");
429 }
430 *problem_ptr = dm_field->problemPtr;
432}
433
434PetscErrorCode DMMoFEMSetDestroyProblem(DM dm, PetscBool destroy_problem) {
435 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
437 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
438 dm_field->destroyProblem = destroy_problem;
440}
441
442PetscErrorCode DMMoFEMGetDestroyProblem(DM dm, PetscBool *destroy_problem) {
443 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
445 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
446 *destroy_problem = dm_field->destroyProblem;
448}
449
450PetscErrorCode DMMoFEMSetSquareProblem(DM dm, PetscBool square_problem) {
451 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
453 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
454 dm_field->isSquareMatrix = square_problem;
456}
457
458PetscErrorCode DMMoFEMResolveSharedFiniteElements(DM dm, std::string fe_name) {
459 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
461 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
463 ->resolveSharedFiniteElements(dm_field->problemPtr, fe_name);
465}
466
467PetscErrorCode DMMoFEMGetProblemFiniteElementLayout(DM dm, std::string fe_name,
468 PetscLayout *layout) {
469 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
471 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
472
473 MPI_Comm comm;
474 CHKERR PetscObjectGetComm((PetscObject)dm, &comm);
475 CHKERR dm_field->problemPtr->getNumberOfElementsByNameAndPart(comm, fe_name,
476 layout);
478}
479
480PetscErrorCode DMMoFEMGetSquareProblem(DM dm, PetscBool *square_problem) {
481 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
483 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
484 *square_problem = dm_field->isSquareMatrix;
486}
487
488PetscErrorCode DMMoFEMAddElement(DM dm, std::string fe_name) {
489 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
491 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
493 dm_field->problemName, fe_name);
495}
496
497PetscErrorCode DMMoFEMAddElement(DM dm, std::vector<std::string> fe_name) {
499 for (auto fe : fe_name) {
501 }
503}
504
505PetscErrorCode DMMoFEMUnSetElement(DM dm, std::string fe_name) {
506 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
508 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
510 dm_field->problemName, fe_name);
512}
513
514PetscErrorCode DMoFEMMeshToLocalVector(DM dm, Vec l, InsertMode mode,
515 ScatterMode scatter_mode) {
516 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
518 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
519 ierr = dm_field->mField_ptr->getInterface<VecManager>()->setLocalGhostVector(
520 dm_field->problemPtr, COL, l, mode, scatter_mode);
521 CHKERRG(ierr);
523}
524
525PetscErrorCode DMoFEMMeshToGlobalVector(DM dm, Vec g, InsertMode mode,
526 ScatterMode scatter_mode) {
527 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
529 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
530 ierr = dm_field->mField_ptr->getInterface<VecManager>()->setGlobalGhostVector(
531 dm_field->problemPtr, COL, g, mode, scatter_mode);
532 CHKERRG(ierr);
534}
535
536PetscErrorCode DMoFEMPreProcessFiniteElements(DM dm, MoFEM::FEMethod *method) {
537 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
539 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
541 dm_field->problemPtr, *method);
542 CHKERRG(ierr);
544}
545
546PetscErrorCode DMoFEMPostProcessFiniteElements(DM dm, MoFEM::FEMethod *method) {
547 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
549 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
551 dm_field->problemPtr, *method);
552 CHKERRG(ierr);
554}
555
556PetscErrorCode
557DMoFEMLoopFiniteElementsUpAndLowRank(DM dm, const char fe_name[],
558 MoFEM::FEMethod *method, int low_rank,
559 int up_rank, CacheTupleWeakPtr cache_ptr) {
561 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
563 dm_field->problemPtr, fe_name, *method, low_rank, up_rank, nullptr,
564 MF_EXIST, cache_ptr);
565 CHKERRG(ierr);
567}
568
570 DM dm, const std::string fe_name, boost::shared_ptr<MoFEM::FEMethod> method,
571 int low_rank, int up_rank, CacheTupleWeakPtr cache_ptr) {
572 return DMoFEMLoopFiniteElementsUpAndLowRank(dm, fe_name.c_str(), method.get(),
573 low_rank, up_rank, cache_ptr);
574}
575
576PetscErrorCode DMoFEMLoopFiniteElements(DM dm, const char fe_name[],
577 MoFEM::FEMethod *method,
578 CacheTupleWeakPtr cache_ptr) {
579 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
581 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
583 dm, fe_name, method, dm_field->rAnk, dm_field->rAnk, cache_ptr);
584 CHKERRG(ierr);
586}
587
588PetscErrorCode
589DMoFEMLoopFiniteElements(DM dm, const std::string fe_name,
590 boost::shared_ptr<MoFEM::FEMethod> method,
591 CacheTupleWeakPtr cache_ptr) {
592 return DMoFEMLoopFiniteElements(dm, fe_name.c_str(), method.get(), cache_ptr);
593}
594
595PetscErrorCode DMoFEMLoopDofs(DM dm, const char field_name[],
596 MoFEM::DofMethod *method) {
597 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
599 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
600 ierr =
601 dm_field->mField_ptr->loop_dofs(dm_field->problemPtr, field_name, COL,
602 *method, dm_field->rAnk, dm_field->rAnk);
603 CHKERRG(ierr);
605}
606
607template <class S, class T0, class T1, class T2>
608static PetscErrorCode DMMoFEMKSPSetComputeRHS(DM dm, S fe_name, T0 method,
609 T1 pre_only, T2 post_only) {
610 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
612 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
613 if (pre_only) {
614 dm_field->kspCtx->getPreProcComputeRhs().push_back(pre_only);
615 }
616 if (method) {
617 dm_field->kspCtx->getComputeRhs().push_back(
618 PairNameFEMethodPtr(fe_name, method));
619 }
620 if (post_only) {
621 dm_field->kspCtx->getPostProcComputeRhs().push_back(post_only);
622 }
623 CHKERR DMKSPSetComputeRHS(dm, KspRhs, dm_field->kspCtx.get());
625}
626
627PetscErrorCode DMMoFEMKSPSetComputeRHS(DM dm, const char fe_name[],
628 MoFEM::FEMethod *method,
629 MoFEM::BasicMethod *pre_only,
630 MoFEM::BasicMethod *post_only) {
631 return DMMoFEMKSPSetComputeRHS<const char *, MoFEM::FEMethod *,
633 dm, fe_name, method, pre_only, post_only);
634}
635
636PetscErrorCode
637DMMoFEMKSPSetComputeRHS(DM dm, const std::string fe_name,
638 boost::shared_ptr<MoFEM::FEMethod> method,
639 boost::shared_ptr<MoFEM::BasicMethod> pre_only,
640 boost::shared_ptr<MoFEM::BasicMethod> post_only) {
641 return DMMoFEMKSPSetComputeRHS<const std::string,
642 boost::shared_ptr<MoFEM::FEMethod>,
643 boost::shared_ptr<MoFEM::BasicMethod>,
644 boost::shared_ptr<MoFEM::BasicMethod>>(
645 dm, fe_name, method, pre_only, post_only);
646}
647
648template <class S, class T0, class T1, class T2>
649static PetscErrorCode DMMoFEMKSPSetComputeOperators(DM dm, S fe_name, T0 method,
650 T1 pre_only, T2 post_only) {
651 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
653 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
654 if (pre_only) {
655 dm_field->kspCtx->getPreProcSetOperators().push_back(pre_only);
656 }
657 if (method) {
658 dm_field->kspCtx->getSetOperators().push_back(
659 PairNameFEMethodPtr(fe_name, method));
660 }
661 if (post_only) {
662 dm_field->kspCtx->getPostProcSetOperators().push_back(post_only);
663 }
664 CHKERR DMKSPSetComputeOperators(dm, KspMat, dm_field->kspCtx.get());
666}
667
668PetscErrorCode DMMoFEMKSPSetComputeOperators(DM dm, const char fe_name[],
669 MoFEM::FEMethod *method,
670 MoFEM::BasicMethod *pre_only,
671 MoFEM::BasicMethod *post_only) {
672 return DMMoFEMKSPSetComputeOperators<const char *, MoFEM::FEMethod *,
675 dm, fe_name, method, pre_only, post_only);
676}
677
678PetscErrorCode
679DMMoFEMKSPSetComputeOperators(DM dm, const std::string fe_name,
680 boost::shared_ptr<MoFEM::FEMethod> method,
681 boost::shared_ptr<MoFEM::BasicMethod> pre_only,
682 boost::shared_ptr<MoFEM::BasicMethod> post_only) {
683 return DMMoFEMKSPSetComputeOperators<const std::string,
684 boost::shared_ptr<MoFEM::FEMethod>>(
685 dm, fe_name, method, pre_only, post_only);
686}
687
688template <class S, class T0, class T1, class T2>
689static PetscErrorCode DMMoFEMSNESSetFunction(DM dm, S fe_name, T0 method,
690 T1 pre_only, T2 post_only) {
691 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
693 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
694 if (pre_only) {
695 dm_field->snesCtx->getPreProcComputeRhs().push_back(pre_only);
696 }
697 if (method) {
698 dm_field->snesCtx->getComputeRhs().push_back(
699 PairNameFEMethodPtr(fe_name, method));
700 }
701 if (post_only) {
702 dm_field->snesCtx->getPostProcComputeRhs().push_back(post_only);
703 }
704 CHKERR DMSNESSetFunction(dm, SnesRhs, dm_field->snesCtx.get());
706}
707
708PetscErrorCode DMMoFEMSNESSetFunction(DM dm, const char fe_name[],
709 MoFEM::FEMethod *method,
710 MoFEM::BasicMethod *pre_only,
711 MoFEM::BasicMethod *post_only) {
712 return DMMoFEMSNESSetFunction<const char *, MoFEM::FEMethod *,
714 dm, fe_name, method, pre_only, post_only);
715}
716
717PetscErrorCode
718DMMoFEMSNESSetFunction(DM dm, const std::string fe_name,
719 boost::shared_ptr<MoFEM::FEMethod> method,
720 boost::shared_ptr<MoFEM::BasicMethod> pre_only,
721 boost::shared_ptr<MoFEM::BasicMethod> post_only) {
722 return DMMoFEMSNESSetFunction<const std::string,
723 boost::shared_ptr<MoFEM::FEMethod>,
724 boost::shared_ptr<MoFEM::BasicMethod>,
725 boost::shared_ptr<MoFEM::BasicMethod>>(
726 dm, fe_name, method, pre_only, post_only);
727}
728
729template <class S, class T0, class T1, class T2>
730static PetscErrorCode DMMoFEMSNESSetJacobian(DM dm, S fe_name, T0 method,
731 T1 pre_only, T2 post_only) {
732 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
734 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
735 if (pre_only) {
736 dm_field->snesCtx->getPreProcSetOperators().push_back(pre_only);
737 }
738 if (method) {
739 dm_field->snesCtx->getSetOperators().push_back(
740 PairNameFEMethodPtr(fe_name, method));
741 }
742 if (post_only) {
743 dm_field->snesCtx->getPostProcSetOperators().push_back(post_only);
744 }
745 CHKERR DMSNESSetJacobian(dm, SnesMat, dm_field->snesCtx.get());
747}
748
749PetscErrorCode DMMoFEMSNESSetJacobian(DM dm, const char fe_name[],
750 MoFEM::FEMethod *method,
751 MoFEM::BasicMethod *pre_only,
752 MoFEM::BasicMethod *post_only) {
753 return DMMoFEMSNESSetJacobian<const char *, MoFEM::FEMethod *,
755 dm, fe_name, method, pre_only, post_only);
756}
757
758PetscErrorCode
759DMMoFEMSNESSetJacobian(DM dm, const std::string fe_name,
760 boost::shared_ptr<MoFEM::FEMethod> method,
761 boost::shared_ptr<MoFEM::BasicMethod> pre_only,
762 boost::shared_ptr<MoFEM::BasicMethod> post_only) {
763 return DMMoFEMSNESSetJacobian<const std::string,
764 boost::shared_ptr<MoFEM::FEMethod>,
765 boost::shared_ptr<MoFEM::BasicMethod>,
766 boost::shared_ptr<MoFEM::BasicMethod>>(
767 dm, fe_name, method, pre_only, post_only);
768}
769
770template <class S, class T0, class T1, class T2>
771static PetscErrorCode DMMoFEMTSSetIFunction(DM dm, S fe_name, T0 method,
772 T1 pre_only, T2 post_only) {
773 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
775 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
776 if (pre_only) {
777 dm_field->tsCtx->getPreProcessIFunction().push_back(pre_only);
778 }
779 if (method) {
780 dm_field->tsCtx->getLoopsIFunction().push_back(
781 PairNameFEMethodPtr(fe_name, method));
782 }
783 if (post_only) {
784 dm_field->tsCtx->getPostProcessIFunction().push_back(post_only);
785 }
786 CHKERR DMTSSetIFunction(dm, TsSetIFunction, dm_field->tsCtx.get());
788}
789
790PetscErrorCode DMMoFEMTSSetIFunction(DM dm, const char fe_name[],
791 MoFEM::FEMethod *method,
792 MoFEM::BasicMethod *pre_only,
793 MoFEM::BasicMethod *post_only) {
794 return DMMoFEMTSSetIFunction<const char *, MoFEM::FEMethod *,
796 dm, fe_name, method, pre_only, post_only);
798}
799
800PetscErrorCode
801DMMoFEMTSSetIFunction(DM dm, const std::string fe_name,
802 boost::shared_ptr<MoFEM::FEMethod> method,
803 boost::shared_ptr<MoFEM::BasicMethod> pre_only,
804 boost::shared_ptr<MoFEM::BasicMethod> post_only) {
805 return DMMoFEMTSSetIFunction<const std::string,
806 boost::shared_ptr<MoFEM::FEMethod>,
807 boost::shared_ptr<MoFEM::BasicMethod>,
808 boost::shared_ptr<MoFEM::BasicMethod>>(
809 dm, fe_name, method, pre_only, post_only);
811}
812
813template <class S, class T0, class T1, class T2>
814static PetscErrorCode DMMoFEMTSSetIJacobian(DM dm, S fe_name, T0 method,
815 T1 pre_only, T2 post_only) {
816 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
818 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
819 if (pre_only) {
820 dm_field->tsCtx->getPreProcessIJacobian().push_back(pre_only);
821 }
822 if (method) {
823 dm_field->tsCtx->getLoopsIJacobian().push_back(
824 PairNameFEMethodPtr(fe_name, method));
825 }
826 if (post_only) {
827 dm_field->tsCtx->getPostProcessIJacobian().push_back(post_only);
828 }
829 CHKERR DMTSSetIJacobian(dm, TsSetIJacobian, dm_field->tsCtx.get());
831}
832
833PetscErrorCode DMMoFEMTSSetIJacobian(DM dm, const char fe_name[],
834 MoFEM::FEMethod *method,
835 MoFEM::BasicMethod *pre_only,
836 MoFEM::BasicMethod *post_only) {
837 return DMMoFEMTSSetIJacobian<const char *, FEMethod *, MoFEM::BasicMethod *,
838 MoFEM::BasicMethod *>(dm, fe_name, method,
839 pre_only, post_only);
840}
841
842PetscErrorCode
843DMMoFEMTSSetIJacobian(DM dm, const std::string fe_name,
844 boost::shared_ptr<MoFEM::FEMethod> method,
845 boost::shared_ptr<MoFEM::BasicMethod> pre_only,
846 boost::shared_ptr<MoFEM::BasicMethod> post_only) {
847 return DMMoFEMTSSetIJacobian<const std::string,
848 boost::shared_ptr<MoFEM::FEMethod>,
849 boost::shared_ptr<MoFEM::BasicMethod>,
850 boost::shared_ptr<MoFEM::BasicMethod>>(
851 dm, fe_name, method, pre_only, post_only);
852}
853
854template <class S, class T0, class T1, class T2>
855static PetscErrorCode DMMoFEMTSSetRHSFunction(DM dm, S fe_name, T0 method,
856 T1 pre_only, T2 post_only) {
857 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
859 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
860 if (pre_only)
861 dm_field->tsCtx->getPreProcessRHSFunction().push_back(pre_only);
862 if (method)
863 dm_field->tsCtx->getLoopsRHSFunction().push_back(
864 PairNameFEMethodPtr(fe_name, method));
865 if (post_only)
866 dm_field->tsCtx->getPostProcessRHSFunction().push_back(post_only);
867 CHKERR DMTSSetRHSFunction(dm, TsSetRHSFunction, dm_field->tsCtx.get());
869}
870
871PetscErrorCode
872DMMoFEMTSSetRHSFunction(DM dm, const std::string fe_name,
873 boost::shared_ptr<MoFEM::FEMethod> method,
874 boost::shared_ptr<MoFEM::BasicMethod> pre_only,
875 boost::shared_ptr<MoFEM::BasicMethod> post_only) {
876 return DMMoFEMTSSetRHSFunction<const std::string,
877 boost::shared_ptr<MoFEM::FEMethod>,
878 boost::shared_ptr<MoFEM::BasicMethod>,
879 boost::shared_ptr<MoFEM::BasicMethod>>(
880 dm, fe_name, method, pre_only, post_only);
882}
883
884PetscErrorCode DMMoFEMTSSetRHSFunction(DM dm, const char fe_name[],
885 MoFEM::FEMethod *method,
886 MoFEM::BasicMethod *pre_only,
887 MoFEM::BasicMethod *post_only) {
888 return DMMoFEMTSSetRHSFunction<const char *, MoFEM::FEMethod *,
890 dm, fe_name, method, pre_only, post_only);
892}
893
894template <class S, class T0, class T1, class T2>
895static PetscErrorCode DMMoFEMTSSetRHSJacobian(DM dm, S fe_name, T0 method,
896 T1 pre_only, T2 post_only) {
897 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
899 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
900 if (pre_only)
901 dm_field->tsCtx->getPreProcessRHSFunction().push_back(pre_only);
902 if (method)
903 dm_field->tsCtx->getLoopsRHSFunction().push_back(
904 PairNameFEMethodPtr(fe_name, method));
905 if (post_only)
906 dm_field->tsCtx->getPostProcessRHSFunction().push_back(post_only);
907 CHKERR DMTSSetRHSJacobian(dm, TsSetRHSJacobian, dm_field->tsCtx.get());
909}
910
911PetscErrorCode
912DMMoFEMTSSetRHSJacobian(DM dm, const std::string fe_name,
913 boost::shared_ptr<MoFEM::FEMethod> method,
914 boost::shared_ptr<MoFEM::BasicMethod> pre_only,
915 boost::shared_ptr<MoFEM::BasicMethod> post_only) {
916 return DMMoFEMTSSetRHSJacobian<const std::string,
917 boost::shared_ptr<MoFEM::FEMethod>,
918 boost::shared_ptr<MoFEM::BasicMethod>,
919 boost::shared_ptr<MoFEM::BasicMethod>>(
920 dm, fe_name, method, pre_only, post_only);
922}
923
924PetscErrorCode DMMoFEMTSSetRHSJacobian(DM dm, const char fe_name[],
925 MoFEM::FEMethod *method,
926 MoFEM::BasicMethod *pre_only,
927 MoFEM::BasicMethod *post_only) {
928 return DMMoFEMTSSetRHSJacobian<const char *, MoFEM::FEMethod *,
930 dm, fe_name, method, pre_only, post_only);
932}
933
934template <class S, class T0, class T1, class T2>
935static PetscErrorCode DMMoFEMTSSetI2Function(DM dm, S fe_name, T0 method,
936 T1 pre_only, T2 post_only) {
937 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
939 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
940 if (pre_only) {
941 dm_field->tsCtx->getPreProcessIFunction().push_back(pre_only);
942 }
943 if (method) {
944 dm_field->tsCtx->getLoopsIFunction().push_back(
945 PairNameFEMethodPtr(fe_name, method));
946 }
947 if (post_only) {
948 dm_field->tsCtx->getPostProcessIFunction().push_back(post_only);
949 }
950 CHKERR DMTSSetI2Function(dm, TsSetI2Function, dm_field->tsCtx.get());
952}
953
954PetscErrorCode DMMoFEMTSSetI2Function(DM dm, const char fe_name[],
955 MoFEM::FEMethod *method,
956 MoFEM::BasicMethod *pre_only,
957 MoFEM::BasicMethod *post_only) {
958 return DMMoFEMTSSetI2Function<const char *, MoFEM::FEMethod *,
960 dm, fe_name, method, pre_only, post_only);
962}
963
964PetscErrorCode
965DMMoFEMTSSetI2Function(DM dm, const std::string fe_name,
966 boost::shared_ptr<MoFEM::FEMethod> method,
967 boost::shared_ptr<MoFEM::BasicMethod> pre_only,
968 boost::shared_ptr<MoFEM::BasicMethod> post_only) {
969 return DMMoFEMTSSetI2Function<const std::string,
970 boost::shared_ptr<MoFEM::FEMethod>,
971 boost::shared_ptr<MoFEM::BasicMethod>,
972 boost::shared_ptr<MoFEM::BasicMethod>>(
973 dm, fe_name, method, pre_only, post_only);
975}
976
977template <class S, class T0, class T1, class T2>
978static PetscErrorCode DMMoFEMTSSetI2Jacobian(DM dm, S fe_name, T0 method,
979 T1 pre_only, T2 post_only) {
980 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
982 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
983 if (pre_only) {
984 dm_field->tsCtx->getPreProcessIJacobian().push_back(pre_only);
985 }
986 if (method) {
987 dm_field->tsCtx->getLoopsIJacobian().push_back(
988 PairNameFEMethodPtr(fe_name, method));
989 }
990 if (post_only) {
991 dm_field->tsCtx->getPostProcessIJacobian().push_back(post_only);
992 }
993 CHKERR DMTSSetI2Jacobian(dm, TsSetI2Jacobian, dm_field->tsCtx.get());
995}
996
997PetscErrorCode DMMoFEMTSSetI2Jacobian(DM dm, const char fe_name[],
998 MoFEM::FEMethod *method,
999 MoFEM::BasicMethod *pre_only,
1000 MoFEM::BasicMethod *post_only) {
1001 return DMMoFEMTSSetI2Jacobian<const char *, FEMethod *, MoFEM::BasicMethod *,
1002 MoFEM::BasicMethod *>(dm, fe_name, method,
1003 pre_only, post_only);
1004}
1005
1006PetscErrorCode
1007DMMoFEMTSSetI2Jacobian(DM dm, const std::string fe_name,
1008 boost::shared_ptr<MoFEM::FEMethod> method,
1009 boost::shared_ptr<MoFEM::BasicMethod> pre_only,
1010 boost::shared_ptr<MoFEM::BasicMethod> post_only) {
1011 return DMMoFEMTSSetI2Jacobian<const std::string,
1012 boost::shared_ptr<MoFEM::FEMethod>,
1013 boost::shared_ptr<MoFEM::BasicMethod>,
1014 boost::shared_ptr<MoFEM::BasicMethod>>(
1015 dm, fe_name, method, pre_only, post_only);
1016}
1017
1018template <class S, class T0, class T1, class T2>
1019static PetscErrorCode DMMoFEMTSSetMonitor(DM dm, TS ts, S fe_name, T0 method,
1020 T1 pre_only, T2 post_only) {
1021 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1023 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
1024 if (pre_only)
1025 dm_field->tsCtx->getPreProcessMonitor().push_back(pre_only);
1026 if (method)
1027 dm_field->tsCtx->getLoopsMonitor().push_back(
1028 PairNameFEMethodPtr(fe_name, method));
1029 if (post_only)
1030 dm_field->tsCtx->getPostProcessMonitor().push_back(post_only);
1031 CHKERR TSMonitorSet(ts, TsMonitorSet, dm_field->tsCtx.get(), PETSC_NULLPTR);
1033}
1034
1035PetscErrorCode DMMoFEMTSSetMonitor(DM dm, TS ts, const char fe_name[],
1036 MoFEM::FEMethod *method,
1037 MoFEM::BasicMethod *pre_only,
1038 MoFEM::BasicMethod *post_only) {
1039 return DMMoFEMTSSetMonitor<const char *, MoFEM::FEMethod *,
1041 dm, ts, fe_name, method, pre_only, post_only);
1043}
1044
1045PetscErrorCode
1046DMMoFEMTSSetMonitor(DM dm, TS ts, const std::string fe_name,
1047 boost::shared_ptr<MoFEM::FEMethod> method,
1048 boost::shared_ptr<MoFEM::BasicMethod> pre_only,
1049 boost::shared_ptr<MoFEM::BasicMethod> post_only) {
1050 return DMMoFEMTSSetMonitor<const std::string,
1051 boost::shared_ptr<MoFEM::FEMethod>,
1052 boost::shared_ptr<MoFEM::BasicMethod>,
1053 boost::shared_ptr<MoFEM::BasicMethod>>(
1054 dm, ts, fe_name, method, pre_only, post_only);
1056}
1057
1058PetscErrorCode DMMoFEMGetKspCtx(DM dm, MoFEM::KspCtx **ksp_ctx) {
1059 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1061 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
1062 *ksp_ctx = dm_field->kspCtx.get();
1064}
1065
1066PetscErrorCode
1067DMMoFEMGetKspCtx(DM dm, const boost::shared_ptr<MoFEM::KspCtx> &ksp_ctx) {
1068 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1070 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
1071 const_cast<boost::shared_ptr<MoFEM::KspCtx> &>(ksp_ctx) = dm_field->kspCtx;
1073}
1074
1075PetscErrorCode DMMoFEMSetKspCtx(DM dm,
1076 boost::shared_ptr<MoFEM::KspCtx> ksp_ctx) {
1077 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1079 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
1080 dm_field->kspCtx = ksp_ctx;
1082}
1083
1084PetscErrorCode DMMoFEMGetSnesCtx(DM dm, MoFEM::SnesCtx **snes_ctx) {
1085 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1087 DMCtxImpl *dm_field = (DMCtxImpl *)dm->data;
1088 *snes_ctx = dm_field->snesCtx.get();
1090}
1091
1092PetscErrorCode
1093DMMoFEMGetSnesCtx(DM dm, const boost::shared_ptr<MoFEM::SnesCtx> &snes_ctx) {
1094 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1096 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
1097 const_cast<boost::shared_ptr<MoFEM::SnesCtx> &>(snes_ctx) = dm_field->snesCtx;
1099}
1100
1101PetscErrorCode DMMoFEMSetSnesCtx(DM dm,
1102 boost::shared_ptr<MoFEM::SnesCtx> snes_ctx) {
1103 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1105 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
1106 dm_field->snesCtx = snes_ctx;
1108}
1109
1110/** get if read mesh is partitioned
1111 * \ingroup dm
1112 */
1113PetscErrorCode DMMoFEMSetIsPartitioned(DM dm, PetscBool is_partitioned) {
1114 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1116 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
1117 dm_field->isPartitioned = is_partitioned;
1119}
1120
1121/** get if read mesh is partitioned
1122 * \ingroup dm
1123 */
1124PetscErrorCode DMMoFEMGetIsPartitioned(DM dm, PetscBool *is_partitioned) {
1125 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1127 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
1128 *is_partitioned = dm_field->isPartitioned;
1130}
1131
1132PetscErrorCode DMMoFEMGetTsCtx(DM dm, MoFEM::TsCtx **ts_ctx) {
1133 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1135 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
1136 *ts_ctx = dm_field->tsCtx.get();
1138}
1139
1140PetscErrorCode DMMoFEMGetTsCtx(DM dm,
1141 const boost::shared_ptr<MoFEM::TsCtx> &ts_ctx) {
1142 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1144 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
1145 const_cast<boost::shared_ptr<MoFEM::TsCtx> &>(ts_ctx) = dm_field->tsCtx;
1147}
1148
1149PetscErrorCode DMMoFEMSetTsCtx(DM dm, boost::shared_ptr<MoFEM::TsCtx> ts_ctx) {
1150 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1152 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
1153 dm_field->tsCtx = ts_ctx;
1155}
1156
1157PetscErrorCode DMCreateGlobalVector_MoFEM(DM dm, Vec *g) {
1158 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1160 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
1161 CHKERR dm_field->mField_ptr->getInterface<VecManager>()->vecCreateGhost(
1162 dm_field->problemName, COL, g);
1163 CHKERR VecSetDM(*g, dm);
1165}
1166
1167PetscErrorCode DMCreateGlobalVector_MoFEM(DM dm, SmartPetscObj<Vec> &g_ptr) {
1168 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1170 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
1171 CHKERR dm_field->mField_ptr->getInterface<VecManager>()->vecCreateGhost(
1172 dm_field->problemName, COL, g_ptr);
1173 CHKERR VecSetDM(g_ptr, dm);
1175}
1176
1177PetscErrorCode DMCreateLocalVector_MoFEM(DM dm, Vec *l) {
1178 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1180 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
1181 CHKERR dm_field->mField_ptr->getInterface<VecManager>()->vecCreateSeq(
1182 dm_field->problemName, COL, l);
1183 CHKERR VecSetDM(*l, dm);
1185}
1186
1187PetscErrorCode DMCreateMatrix_MoFEM(DM dm, Mat *M) {
1188 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1190 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
1191
1192 if (strcmp(dm->mattype, MATSHELL) == 0) {
1193
1194 if (dm_field->blocMatDataPtr) {
1196 CHKERR MatSetDM(*M, dm);
1198 } else {
1199 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
1200 "Matrix shell data not set, or matrix type not implemented");
1201 }
1202
1203 } else if (strcmp(dm->mattype, MATMPIAIJ) == 0) {
1205 ->createMPIAIJWithArrays<PetscGlobalIdx_mi_tag>(dm_field->problemName,
1206 M);
1207 } else if (strcmp(dm->mattype, MATAIJ) == 0) {
1209 ->createSeqAIJWithArrays<PetscLocalIdx_mi_tag>(dm_field->problemName,
1210 M);
1211 } else if (strcmp(dm->mattype, MATAIJCUSPARSE) == 0) {
1213 ->createMPIAIJCUSPARSEWithArrays<PetscGlobalIdx_mi_tag>(
1214 dm_field->problemName, M);
1215 } else if (strcmp(dm->mattype, MATSEQAIJCUSPARSE) == 0) {
1217 ->createSeqAIJCUSPARSEWithArrays<PetscLocalIdx_mi_tag>(
1218 dm_field->problemName, M);
1219 } else {
1220 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
1221 "Matrix type not implemented");
1222 }
1223 CHKERR MatSetDM(*M, dm);
1225}
1226
1227PetscErrorCode DMCreateMatrix_MoFEM(DM dm, SmartPetscObj<Mat> &M) {
1228 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1230 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
1231
1232 if (strcmp(dm->mattype, MATSHELL) == 0) {
1233 if (dm_field->blocMatDataPtr) {
1234 Mat mat_raw;
1235 CHKERR DMMoFEMCreateBlockMat(dm, &mat_raw);
1236 M = SmartPetscObj<Mat>(mat_raw);
1238 } else {
1239 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
1240 "Matrix shell data not set, or matrix type not implemented");
1241 }
1242 } else if (strcmp(dm->mattype, MATMPIAIJ) == 0) {
1244 ->createMPIAIJWithArrays<PetscGlobalIdx_mi_tag>(dm_field->problemName,
1245 M);
1246 } else if (strcmp(dm->mattype, MATAIJ) == 0) {
1248 ->createSeqAIJWithArrays<PetscLocalIdx_mi_tag>(dm_field->problemName,
1249 M);
1250 } else if (strcmp(dm->mattype, MATAIJCUSPARSE) == 0) {
1252 ->createMPIAIJCUSPARSEWithArrays<PetscGlobalIdx_mi_tag>(
1253 dm_field->problemName, M);
1254 } else if (strcmp(dm->mattype, MATSEQAIJCUSPARSE) == 0) {
1256 ->createSeqAIJCUSPARSEWithArrays<PetscLocalIdx_mi_tag>(
1257 dm_field->problemName, M);
1258 } else {
1259 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
1260 "Matrix type not implemented");
1261 }
1262 CHKERR MatSetDM(M, dm);
1264}
1265
1266#if PETSC_VERSION_GE(3, 17, 0)
1267PetscErrorCode DMSetFromOptions_MoFEM(DM dm,
1268 PetscOptionItems *PetscOptionsObject) {
1269#elif PETSC_VERSION_GE(3, 7, 0)
1270PetscErrorCode DMSetFromOptions_MoFEM(PetscOptionItems *PetscOptionsObject,
1271 DM dm) {
1272#elif PETSC_VERSION_GE(3, 5, 3)
1273PetscErrorCode DMSetFromOptions_MoFEM(PetscOptions *PetscOptionsObject, DM dm) {
1274#else
1275PetscErrorCode DMSetFromOptions_MoFEM(DM dm) {
1276#endif
1277
1278 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1280 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
1281#if PETSC_VERSION_GE(3, 18, 0)
1282 PetscOptionsHeadBegin(PetscOptionsObject, "DMMoFEM Options");
1283#elif PETSC_VERSION_GE(3, 5, 3)
1284 PetscOptionsHead(PetscOptionsObject, "DMMoFEM Options");
1285#else
1286 PetscOptionsHead("DMMoFEM Options");
1287#endif
1288 ierr = PetscOptionsBool("-dm_is_partitioned",
1289 "set if mesh is partitioned (works which native MOAB "
1290 "file format, i.e. h5m",
1291 "DMSetUp", dm_field->isPartitioned,
1292 &dm_field->isPartitioned, NULL);
1293 CHKERRG(ierr);
1295}
1296
1297PetscErrorCode DMSetUp_MoFEM(DM dm) {
1298 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1299 ProblemsManager *prb_mng_ptr;
1301 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
1302 CHKERR dm_field->mField_ptr->getInterface(prb_mng_ptr);
1303
1304 if (dm_field->isCompDM) {
1305 // It is composite probelm
1306 CHKERR prb_mng_ptr->buildComposedProblem(
1307 dm_field->problemName, dm_field->rowCompPrb, dm_field->colCompPrb,
1308 dm_field->isSquareMatrix == PETSC_TRUE, dm_field->verbosity);
1309 } else {
1310 if (dm_field->isPartitioned) {
1312 dm_field->problemName, dm_field->isSquareMatrix == PETSC_TRUE,
1313 dm_field->verbosity);
1314 } else {
1315 CHKERR prb_mng_ptr->buildProblem(dm_field->problemName,
1316 dm_field->isSquareMatrix == PETSC_TRUE,
1317 dm_field->verbosity);
1318 CHKERR prb_mng_ptr->partitionProblem(dm_field->problemName,
1319 dm_field->verbosity);
1320 }
1321 }
1322
1323 // Partition finite elements
1324 if (dm_field->isPartitioned) {
1325 CHKERR prb_mng_ptr->partitionFiniteElements(
1326 dm_field->problemName, true, 0, dm_field->sIze, dm_field->verbosity);
1328 dm_field->problemName, dm_field->verbosity);
1329 } else {
1330 // partition finite elemnets
1331 CHKERR prb_mng_ptr->partitionFiniteElements(dm_field->problemName, false,
1332 -1, -1, dm_field->verbosity);
1333 // Get ghost DOFs
1334 CHKERR prb_mng_ptr->partitionGhostDofs(dm_field->problemName,
1335 dm_field->verbosity);
1336 }
1337
1338 // Set flag that problem is build and partitioned
1339 dm_field->isProblemBuild = PETSC_TRUE;
1340
1342}
1343
1344PetscErrorCode DMSubDMSetUp_MoFEM(DM subdm) {
1345 PetscValidHeaderSpecific(subdm, DM_CLASSID, 1);
1346 ProblemsManager *prb_mng_ptr;
1348
1349 DMCtxImpl *subdm_field = static_cast<DMCtxImpl *>(subdm->data);
1350
1351 // build sub dm problem
1352 CHKERR subdm_field->mField_ptr->getInterface(prb_mng_ptr);
1353
1354 map<std::string, boost::shared_ptr<Range>> *entity_map_row = nullptr;
1355 map<std::string, boost::shared_ptr<Range>> *entity_map_col = nullptr;
1356
1357 if (subdm_field->mapTypeRow.size())
1358 entity_map_row = &subdm_field->mapTypeRow;
1359 if (subdm_field->mapTypeCol.size())
1360 entity_map_col = &subdm_field->mapTypeCol;
1361
1362 CHKERR prb_mng_ptr->buildSubProblem(
1363 subdm_field->problemName, subdm_field->rowSubFields,
1364 subdm_field->colSubFields, subdm_field->problemMainOfSubPtr->getName(),
1365 subdm_field->isSquareMatrix == PETSC_TRUE, entity_map_row, entity_map_col,
1366 subdm_field->verbosity);
1367
1368 // partition problem
1369 subdm_field->isPartitioned = subdm_field->isPartitioned;
1370 if (subdm_field->isPartitioned) {
1371 CHKERR prb_mng_ptr->partitionFiniteElements(subdm_field->problemName, true,
1372 0, subdm_field->sIze,
1373 subdm_field->verbosity);
1374 // set ghost nodes
1376 subdm_field->problemName, subdm_field->verbosity);
1377 } else {
1378 // partition finite elements
1379 CHKERR prb_mng_ptr->partitionFiniteElements(subdm_field->problemName, false,
1380 -1, -1, subdm_field->verbosity);
1381 // set ghost nodes
1382 CHKERR prb_mng_ptr->partitionGhostDofs(subdm_field->problemName,
1383 subdm_field->verbosity);
1384 }
1385
1386 subdm_field->isProblemBuild = PETSC_TRUE;
1387
1389}
1390
1391PetscErrorCode DMGlobalToLocalBegin_MoFEM(DM dm, Vec g, InsertMode mode,
1392 Vec l) {
1393 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1395 CHKERR VecGhostUpdateBegin(g, INSERT_VALUES, SCATTER_FORWARD);
1397}
1398
1399PetscErrorCode DMGlobalToLocalEnd_MoFEM(DM dm, Vec g, InsertMode mode, Vec l) {
1400 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1402
1403 CHKERR VecGhostUpdateEnd(g, INSERT_VALUES, SCATTER_FORWARD);
1404
1405 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
1406 int nb_dofs = dm_field->problemPtr->getNbLocalDofsRow();
1407 int nb_ghost = dm_field->problemPtr->getNbGhostDofsRow();
1408
1409 double *array_loc, *array_glob;
1410 CHKERR VecGetArray(l, &array_loc);
1411 CHKERR VecGetArray(g, &array_glob);
1412 switch (mode) {
1413 case INSERT_VALUES:
1414 cblas_dcopy(nb_dofs + nb_ghost, array_glob, 1, array_loc, 1);
1415 break;
1416 case ADD_VALUES:
1417 cblas_daxpy(nb_dofs + nb_ghost, 1, array_glob, 1, array_loc, 1);
1418 break;
1419 default:
1420 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "not implemented");
1421 }
1422 CHKERR VecRestoreArray(l, &array_loc);
1423 CHKERR VecRestoreArray(g, &array_glob);
1425}
1426
1427PetscErrorCode DMLocalToGlobalBegin_MoFEM(DM dm, Vec l, InsertMode mode,
1428 Vec g) {
1429
1430 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1432
1433 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
1434 int nb_dofs = dm_field->problemPtr->getNbLocalDofsRow();
1435 int nb_ghost = dm_field->problemPtr->getNbGhostDofsRow();
1436
1437 double *array_loc, *array_glob;
1438 CHKERR VecGetArray(l, &array_loc);
1439 CHKERR VecGetArray(g, &array_glob);
1440 switch (mode) {
1441 case INSERT_VALUES:
1442 cblas_dcopy(nb_dofs + nb_ghost, array_loc, 1, array_glob, 1);
1443 break;
1444 case ADD_VALUES:
1445 cblas_daxpy(nb_dofs + nb_ghost, 1, array_loc, 1, array_glob, 1);
1446 break;
1447 default:
1448 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "not implemented");
1449 }
1450 CHKERR VecRestoreArray(l, &array_loc);
1451 CHKERR VecRestoreArray(g, &array_glob);
1452
1454}
1455
1456PetscErrorCode DMLocalToGlobalEnd_MoFEM(DM, Vec l, InsertMode mode, Vec g) {
1459}
1460
1461PetscErrorCode DMCreateFieldIS_MoFEM(DM dm, PetscInt *numFields,
1462 char ***fieldNames, IS **fields) {
1463 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1465
1466 if (numFields) {
1467 *numFields = 0;
1468 }
1469 if (fieldNames) {
1470 *fieldNames = NULL;
1471 }
1472 if (fields) {
1473 *fields = NULL;
1474 }
1475
1476 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
1477 auto fields_ptr = dm_field->mField_ptr->get_fields();
1478 Field_multiIndex::iterator fit, hi_fit;
1479 fit = fields_ptr->begin();
1480 hi_fit = fields_ptr->end();
1481 *numFields = std::distance(fit, hi_fit);
1482
1483 if (fieldNames) {
1484 CHKERR PetscMalloc1(*numFields, fieldNames);
1485 }
1486 if (fields) {
1487 CHKERR PetscMalloc1(*numFields, fields);
1488 }
1489
1490 for (int f = 0; fit != hi_fit; fit++, f++) {
1491 if (fieldNames) {
1492 CHKERR PetscStrallocpy(fit->get()->getName().c_str(),
1493 (char **)&(*fieldNames)[f]);
1494 }
1495 if (fields) {
1497 ->isCreateProblemFieldAndRank(
1498 dm_field->problemPtr->getName(), ROW, fit->get()->getName(), 0,
1499 fit->get()->getNbOfCoeffs(), &(*fields)[f]);
1500 }
1501 }
1502
1504}
1505
1506PetscErrorCode DMMoFEMGetFieldIS(DM dm, RowColData rc, const char field_name[],
1507 IS *is) {
1508 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1510 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
1512 ->isCreateProblemFieldAndRank(dm_field->problemPtr->getName(), ROW,
1513 field_name, 0, 1000, is);
1515}
1516
1517PetscErrorCode DMMoFEMSetVerbosity(DM dm, const int verb) {
1518 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1520 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
1521 dm_field->verbosity = verb;
1523}
1524
1526 boost::shared_ptr<BlockStructure> data) {
1527 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1529 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
1530 dm_field->blocMatDataPtr = data;
1532}
1533
1535 boost::shared_ptr<BlockStructure> &data) {
1536 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1538 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
1539 data = dm_field->blocMatDataPtr;
1541}
1542
1544 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1546 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
1547 auto mat_data = createBlockMat(dm, dm_field->blocMatDataPtr);
1548 *mat = mat_data.first;
1549 CHKERR PetscObjectReference((PetscObject)(*mat));
1551}
1552
1553template <>
1555 boost::shared_ptr<NestSchurData> data) {
1556 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1558 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
1559 dm_field->nestedSchurDataPtr = data;
1560 dm_field->blocMatDataPtr = get<2>(*(data));
1562}
1563
1565 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1567 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
1568 auto mat_data = createSchurNestedMatrix(dm_field->nestedSchurDataPtr);
1569 *mat = mat_data.first;
1570 CHKERR PetscObjectReference((PetscObject)(*mat));
1572}
1573
1575 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1577 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
1578 auto matrix_mng = dm_field->mField_ptr->getInterface<MatrixManager>();
1580 dm_field->problemName, mat);
1582}
1583
1584} // namespace MoFEM
Implementation of DM context. You should not use it directly.
Discrete manager interface for MoFEM.
RowColData
RowColData.
@ COL
@ ROW
@ MF_EXIST
@ MF_EXCL
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#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_DATA_INCONSISTENCY
Definition definitions.h:31
@ MOFEM_INVALID_DATA
Definition definitions.h:36
@ MOFEM_NOT_IMPLEMENTED
Definition definitions.h:32
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define CHKERR
Inline error check.
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
PetscErrorCode DMGlobalToLocalBegin_MoFEM(DM dm, Vec, InsertMode, Vec)
Definition DMMoFEM.cpp:1391
PetscErrorCode DMGlobalToLocalEnd_MoFEM(DM dm, Vec, InsertMode, Vec)
Definition DMMoFEM.cpp:1399
PetscErrorCode DMMoFEMSetIsPartitioned(DM dm, PetscBool is_partitioned)
Definition DMMoFEM.cpp:1113
PetscErrorCode DMMoFEMCreateSubDM(DM subdm, DM dm, const char problem_name[])
Must be called by user to set Sub DM MoFEM data structures.
Definition DMMoFEM.cpp:215
PetscErrorCode DMMoFEMGetIsPartitioned(DM dm, PetscBool *is_partitioned)
Definition DMMoFEM.cpp:1124
PetscErrorCode DMMoFEMAddElement(DM dm, std::string fe_name)
add element to dm
Definition DMMoFEM.cpp:488
PetscErrorCode DMMoFEMSNESSetFunction(DM dm, const char fe_name[], MoFEM::FEMethod *method, MoFEM::BasicMethod *pre_only, MoFEM::BasicMethod *post_only)
set SNES residual evaluation function
Definition DMMoFEM.cpp:708
PetscErrorCode DMMoFEMGetProblemFiniteElementLayout(DM dm, std::string fe_name, PetscLayout *layout)
Get finite elements layout in the problem.
Definition DMMoFEM.cpp:467
PetscErrorCode DMCreateLocalVector_MoFEM(DM dm, Vec *l)
DMShellSetCreateLocalVector.
Definition DMMoFEM.cpp:1177
PetscErrorCode DMMoFEMSetSquareProblem(DM dm, PetscBool square_problem)
set squared problem
Definition DMMoFEM.cpp:450
PetscErrorCode DMMoFEMTSSetIFunction(DM dm, const char fe_name[], MoFEM::FEMethod *method, MoFEM::BasicMethod *pre_only, MoFEM::BasicMethod *post_only)
set TS implicit function evaluation function
Definition DMMoFEM.cpp:790
PetscErrorCode DMDestroy_MoFEM(DM dm)
Destroys dm with MoFEM data structure.
Definition DMMoFEM.cpp:79
PetscErrorCode DMMoFEMCreateMoFEM(DM dm, MoFEM::Interface *m_field_ptr, const char problem_name[], const MoFEM::BitRefLevel bit_level, const MoFEM::BitRefLevel bit_mask=MoFEM::BitRefLevel().set())
Must be called by user to set MoFEM data structures.
Definition DMMoFEM.cpp:114
PetscErrorCode DMoFEMPostProcessFiniteElements(DM dm, MoFEM::FEMethod *method)
execute finite element method for each element in dm (problem)
Definition DMMoFEM.cpp:546
PetscErrorCode DMoFEMMeshToLocalVector(DM dm, Vec l, InsertMode mode, ScatterMode scatter_mode)
set local (or ghosted) vector values on mesh for partition only
Definition DMMoFEM.cpp:514
PetscErrorCode DMMoFEMAddColCompositeProblem(DM dm, const char prb_name[])
Add problem to composite DM on col.
Definition DMMoFEM.cpp:382
PetscErrorCode DMLocalToGlobalBegin_MoFEM(DM, Vec, InsertMode, Vec)
Definition DMMoFEM.cpp:1427
PetscErrorCode DMMoFEMGetProblemPtr(DM dm, const MoFEM::Problem **problem_ptr)
Get pointer to problem data structure.
Definition DMMoFEM.cpp:422
PetscErrorCode DMCreate_MoFEM(DM dm)
Create dm data structure with MoFEM data structure.
Definition DMMoFEM.cpp:71
PetscErrorCode DMMoFEMSetSnesCtx(DM dm, boost::shared_ptr< MoFEM::SnesCtx > snes_ctx)
Set MoFEM::SnesCtx data structure.
Definition DMMoFEM.cpp:1101
PetscErrorCode DMMoFEMSetKspCtx(DM dm, boost::shared_ptr< MoFEM::KspCtx > ksp_ctx)
set MoFEM::KspCtx data structure
Definition DMMoFEM.cpp:1075
PetscErrorCode DMMoFEMAddSubFieldRow(DM dm, const char field_name[])
Definition DMMoFEM.cpp:238
PetscErrorCode DMMoFEMGetIsSubDM(DM dm, PetscBool *is_sub_dm)
Definition DMMoFEM.cpp:322
PetscErrorCode DMMoFEMAddRowCompositeProblem(DM dm, const char prb_name[])
Add problem to composite DM on row.
Definition DMMoFEM.cpp:364
PetscErrorCode DMMoFEMSNESSetJacobian(DM dm, const char fe_name[], MoFEM::FEMethod *method, MoFEM::BasicMethod *pre_only, MoFEM::BasicMethod *post_only)
set SNES Jacobian evaluation function
Definition DMMoFEM.cpp:749
PetscErrorCode DMCreateMatrix_MoFEM(DM dm, Mat *M)
Definition DMMoFEM.cpp:1187
PetscErrorCode DMSetFromOptions_MoFEM(DM dm)
Definition DMMoFEM.cpp:1275
PetscErrorCode DMMoFEMGetTsCtx(DM dm, MoFEM::TsCtx **ts_ctx)
get MoFEM::TsCtx data structure
Definition DMMoFEM.cpp:1132
PetscErrorCode DMRegister_MoFEM(const char sname[])
Register MoFEM problem.
Definition DMMoFEM.cpp:43
PetscErrorCode DMMoFEMGetKspCtx(DM dm, MoFEM::KspCtx **ksp_ctx)
get MoFEM::KspCtx data structure
Definition DMMoFEM.cpp:1058
PetscErrorCode DMMoFEMGetSquareProblem(DM dm, PetscBool *square_problem)
get squared problem
Definition DMMoFEM.cpp:480
PetscErrorCode DMoFEMLoopDofs(DM dm, const char field_name[], MoFEM::DofMethod *method)
execute method for dofs on field in problem
Definition DMMoFEM.cpp:595
PetscErrorCode DMMoFEMKSPSetComputeRHS(DM dm, const char fe_name[], MoFEM::FEMethod *method, MoFEM::BasicMethod *pre_only, MoFEM::BasicMethod *post_only)
Set compute operator for KSP solver via sub-matrix and IS.
Definition DMMoFEM.cpp:627
PetscErrorCode DMMoFEMResolveSharedFiniteElements(DM dm, std::string fe_name)
Resolve shared entities.
Definition DMMoFEM.cpp:458
PetscErrorCode DMMoFEMSetTsCtx(DM dm, boost::shared_ptr< MoFEM::TsCtx > ts_ctx)
Set MoFEM::TsCtx data structure.
Definition DMMoFEM.cpp:1149
PetscErrorCode DMoFEMLoopFiniteElements(DM dm, const char fe_name[], MoFEM::FEMethod *method, CacheTupleWeakPtr cache_ptr=CacheTupleSharedPtr())
Executes FEMethod for finite elements in DM.
Definition DMMoFEM.cpp:576
PetscErrorCode DMMoFEMGetFieldIS(DM dm, RowColData rc, const char field_name[], IS *is)
get field is in the problem
Definition DMMoFEM.cpp:1506
PetscErrorCode DMMoFEMTSSetIJacobian(DM dm, const std::string fe_name, boost::shared_ptr< MoFEM::FEMethod > method, boost::shared_ptr< MoFEM::BasicMethod > pre_only, boost::shared_ptr< MoFEM::BasicMethod > post_only)
set TS Jacobian evaluation function
Definition DMMoFEM.cpp:843
MoFEMErrorCode DMMoFEMCreateBlockMat(DM dm, Mat *mat)
Create block matrix.
Definition DMMoFEM.cpp:1543
PetscErrorCode DMoFEMGetInterfacePtr(DM dm, MoFEM::Interface **m_field_ptr)
Get pointer to MoFEM::Interface.
Definition DMMoFEM.cpp:410
PetscErrorCode DMMoFEMAddSubFieldCol(DM dm, const char field_name[])
Definition DMMoFEM.cpp:280
PetscErrorCode DMMoFEMTSSetI2Jacobian(DM dm, const std::string fe_name, boost::shared_ptr< MoFEM::FEMethod > method, boost::shared_ptr< MoFEM::BasicMethod > pre_only, boost::shared_ptr< MoFEM::BasicMethod > post_only)
set TS Jacobian evaluation function
Definition DMMoFEM.cpp:1007
PetscErrorCode DMMoFEMTSSetI2Function(DM dm, const std::string fe_name, boost::shared_ptr< MoFEM::FEMethod > method, boost::shared_ptr< MoFEM::BasicMethod > pre_only, boost::shared_ptr< MoFEM::BasicMethod > post_only)
set TS implicit function evaluation function
Definition DMMoFEM.cpp:965
PetscErrorCode DMCreateGlobalVector_MoFEM(DM dm, Vec *g)
DMShellSetCreateGlobalVector.
Definition DMMoFEM.cpp:1157
PetscErrorCode DMMoFEMGetSnesCtx(DM dm, MoFEM::SnesCtx **snes_ctx)
get MoFEM::SnesCtx data structure
Definition DMMoFEM.cpp:1084
PetscErrorCode DMoFEMMeshToGlobalVector(DM dm, Vec g, InsertMode mode, ScatterMode scatter_mode)
set ghosted vector values on all existing mesh entities
Definition DMMoFEM.cpp:525
PetscErrorCode DMMoFEMGetIsCompDM(DM dm, PetscBool *is_comp_dm)
Get if this DM is composite DM.
Definition DMMoFEM.cpp:402
PetscErrorCode DMSubDMSetUp_MoFEM(DM subdm)
Definition DMMoFEM.cpp:1344
PetscErrorCode DMMoFEMKSPSetComputeOperators(DM dm, const char fe_name[], MoFEM::FEMethod *method, MoFEM::BasicMethod *pre_only, MoFEM::BasicMethod *post_only)
Set KSP operators and push mofem finite element methods.
Definition DMMoFEM.cpp:668
PetscErrorCode DMCreateFieldIS_MoFEM(DM dm, PetscInt *numFields, char ***fieldNames, IS **fields)
Definition DMMoFEM.cpp:1461
PetscErrorCode DMMoFEMUnSetElement(DM dm, std::string fe_name)
unset element from dm
Definition DMMoFEM.cpp:505
PetscErrorCode DMSetUp_MoFEM(DM dm)
Definition DMMoFEM.cpp:1297
PetscErrorCode DMLocalToGlobalEnd_MoFEM(DM, Vec, InsertMode, Vec)
Definition DMMoFEM.cpp:1456
PetscErrorCode DMoFEMLoopFiniteElementsUpAndLowRank(DM dm, const char fe_name[], MoFEM::FEMethod *method, int low_rank, int up_rank, CacheTupleWeakPtr cache_ptr=CacheTupleSharedPtr())
Executes FEMethod for finite elements in DM.
Definition DMMoFEM.cpp:557
PetscErrorCode DMoFEMPreProcessFiniteElements(DM dm, MoFEM::FEMethod *method)
execute finite element method for each element in dm (problem)
Definition DMMoFEM.cpp:536
PetscErrorCode DMSetOperators_MoFEM(DM dm)
Set operators for MoFEM dm.
Definition DMMoFEM.cpp:49
virtual const Problem * get_problem(const std::string problem_name) const =0
Get the problem object.
virtual const Field_multiIndex * get_fields() const =0
Get the fields object.
static LoggerType & setLog(const std::string channel)
Set ans resset chanel logger.
#define MOFEM_LOG(channel, severity)
Log.
#define MOFEM_LOG_TAG(channel, tag)
Tag channel.
virtual MoFEMErrorCode loop_dofs(const Problem *problem_ptr, const std::string &field_name, RowColData rc, DofMethod &method, int lower_rank, int upper_rank, int verb=DEFAULT_VERBOSITY)=0
Make a loop over dofs.
virtual MoFEMErrorCode problem_basic_method_postProcess(const Problem *problem_ptr, BasicMethod &method, int verb=DEFAULT_VERBOSITY)=0
Set data for BasicMethod.
virtual MoFEMErrorCode loop_finite_elements(const std::string problem_name, const std::string &fe_name, FEMethod &method, boost::shared_ptr< NumeredEntFiniteElement_multiIndex > fe_ptr=nullptr, MoFEMTypes bh=MF_EXIST, CacheTupleWeakPtr cache_ptr=CacheTupleSharedPtr(), int verb=DEFAULT_VERBOSITY)=0
Make a loop over finite elements.
MoFEMErrorCode buildProblemOnDistributedMesh(const std::string name, const bool square_matrix, int verb=VERBOSE)
build problem data structures, assuming that mesh is distributed (collective)
MoFEMErrorCode partitionGhostDofs(const std::string name, int verb=VERBOSE)
determine ghost nodes
MoFEMErrorCode buildComposedProblem(const std::string out_name, const std::vector< std::string > add_row_problems, const std::vector< std::string > add_col_problems, const bool square_matrix=true, int verb=1)
build composite problem
MoFEMErrorCode buildProblem(const std::string name, const bool square_matrix, int verb=VERBOSE)
build problem data structures
MoFEMErrorCode buildSubProblem(const std::string out_name, const std::vector< std::string > &fields_row, const std::vector< std::string > &fields_col, const std::string main_problem, const bool square_matrix=true, const map< std::string, boost::shared_ptr< Range > > *entityMapRow=nullptr, const map< std::string, boost::shared_ptr< Range > > *entityMapCol=nullptr, int verb=VERBOSE)
build sub problem
MoFEMErrorCode partitionProblem(const std::string name, int verb=VERBOSE)
partition problem dofs (collective)
MoFEMErrorCode partitionGhostDofsOnDistributedMesh(const std::string name, int verb=VERBOSE)
determine ghost nodes on distributed meshes
MoFEMErrorCode partitionFiniteElements(const std::string name, bool part_from_moab=false, int low_proc=-1, int hi_proc=-1, int verb=VERBOSE)
partition finite elements
virtual MoFEMErrorCode add_problem(const std::string &name, enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
Add problem.
virtual MoFEMErrorCode modify_problem_unset_finite_element(const std::string name_problem, const std::string &fe_name)=0
unset finite element from problem, this remove entities assigned to finite element to a particular pr...
virtual MoFEMErrorCode delete_problem(const std::string name)=0
Delete problem.
virtual MoFEMErrorCode modify_problem_ref_level_add_bit(const std::string &name_problem, const BitRefLevel &bit)=0
add ref level to problem
virtual bool check_problem(const std::string name)=0
check if problem exist
virtual MoFEMErrorCode modify_problem_mask_ref_level_add_bit(const std::string &name_problem, const BitRefLevel &bit)=0
set dof mask ref level for problem
virtual MoFEMErrorCode modify_problem_add_finite_element(const std::string name_problem, const std::string &fe_name)=0
add finite element to problem, this add entities assigned to finite element to a particular problem
MoFEM::TsCtx * ts_ctx
FTensor::Index< 'l', 3 > l
static MoFEMErrorCodeGeneric< PetscErrorCode > ierr
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
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
boost::weak_ptr< CacheTuple > CacheTupleWeakPtr
PetscErrorCode DMMoFEMTSSetMonitor(DM dm, TS ts, const std::string fe_name, boost::shared_ptr< MoFEM::FEMethod > method, boost::shared_ptr< MoFEM::BasicMethod > pre_only, boost::shared_ptr< MoFEM::BasicMethod > post_only)
Set Monitor To TS solver.
Definition DMMoFEM.cpp:1046
PetscErrorCode TsSetIJacobian(TS ts, PetscReal t, Vec u, Vec u_t, PetscReal a, Mat A, Mat B, void *ctx)
Set function evaluating jacobian in TS solver.
Definition TsCtx.cpp:169
PetscErrorCode TsMonitorSet(TS ts, PetscInt step, PetscReal t, Vec u, void *ctx)
Set monitor for TS solver.
Definition TsCtx.cpp:263
PetscErrorCode DMMoFEMTSSetRHSFunction(DM dm, const std::string fe_name, boost::shared_ptr< MoFEM::FEMethod > method, boost::shared_ptr< MoFEM::BasicMethod > pre_only, boost::shared_ptr< MoFEM::BasicMethod > post_only)
set TS the right hand side function
Definition DMMoFEM.cpp:872
PetscErrorCode SnesMat(SNES snes, Vec x, Mat A, Mat B, void *ctx)
This is MoFEM implementation for the left hand side (tangent matrix) evaluation in SNES solver.
Definition SnesCtx.cpp:142
PetscErrorCode DMMoFEMSetDestroyProblem(DM dm, PetscBool destroy_problem)
Definition DMMoFEM.cpp:434
MoFEMErrorCode DMMoFEMCreateHybridL2Mat(DM dm, SmartPetscObj< Mat > &mat)
Create matrix for hybridised system.
Definition DMMoFEM.cpp:1574
PetscErrorCode DMMoFEMGetDestroyProblem(DM dm, PetscBool *destroy_problem)
Definition DMMoFEM.cpp:442
MoFEMErrorCode DMMoFEMGetBlocMatData(DM dm, boost::shared_ptr< BlockStructure > &)
Get data for block mat.
Definition DMMoFEM.cpp:1534
PetscErrorCode TsSetI2Function(TS ts, PetscReal t, Vec u, Vec u_t, Vec u_tt, Vec F, void *ctx)
Calculation the right hand side for second order PDE in time.
Definition TsCtx.cpp:620
PetscErrorCode TsSetRHSFunction(TS ts, PetscReal t, Vec u, Vec F, void *ctx)
TS solver function.
Definition TsCtx.cpp:327
std::pair< SmartPetscObj< Mat >, boost::shared_ptr< NestSchurData > > createSchurNestedMatrix(boost::shared_ptr< NestSchurData > schur_net_data_ptr)
Create a Mat Diag Blocks object.
Definition Schur.cpp:2552
PetscErrorCode SnesRhs(SNES snes, Vec x, Vec f, void *ctx)
This is MoFEM implementation for the right hand side (residual vector) evaluation in SNES solver.
Definition SnesCtx.cpp:27
PetscErrorCode DMMoFEMGetSubRowIS(DM dm, IS *is)
get sub problem is
Definition DMMoFEM.cpp:330
MoFEMErrorCode DMMoFEMCreateNestSchurMat(DM dm, Mat *mat)
Create nest schur matrix.
Definition DMMoFEM.cpp:1564
PetscErrorCode TsSetIFunction(TS ts, PetscReal t, Vec u, Vec u_t, Vec F, void *ctx)
Set IFunction for TS solver.
Definition TsCtx.cpp:56
PetscErrorCode DMMoFEMGetSubColIS(DM dm, IS *is)
get sub problem is
Definition DMMoFEM.cpp:347
PetscErrorCode KspRhs(KSP ksp, Vec f, void *ctx)
Run over elements in the lists.
Definition KspCtx.cpp:21
MoFEMErrorCode DMMoFEMSetBlocMatData(DM dm, boost::shared_ptr< BlockStructure >)
Set data for block mat.
Definition DMMoFEM.cpp:1525
PetscErrorCode TsSetI2Jacobian(TS ts, PetscReal t, Vec u, Vec u_t, Vec u_tt, PetscReal a, PetscReal aa, Mat A, Mat B, void *ctx)
Calculation Jacobian for second order PDE in time.
Definition TsCtx.cpp:519
PetscErrorCode KspMat(KSP ksp, Mat A, Mat B, void *ctx)
Run over elements in the list.
Definition KspCtx.cpp:86
PetscErrorCode TsSetRHSJacobian(TS ts, PetscReal t, Vec u, Mat A, Mat B, void *ctx)
TS solver function.
Definition TsCtx.cpp:430
MoFEMErrorCode DMMoFEMSetNestSchurData(DM dm, boost::shared_ptr< NestSchurData >)
Definition DMMoFEM.cpp:1554
PetscErrorCode DMMoFEMTSSetRHSJacobian(DM dm, const std::string fe_name, boost::shared_ptr< MoFEM::FEMethod > method, boost::shared_ptr< MoFEM::BasicMethod > pre_only, boost::shared_ptr< MoFEM::BasicMethod > post_only)
set TS the right hand side jacobian
Definition DMMoFEM.cpp:912
PetscErrorCode DMMoFEMSetVerbosity(DM dm, const int verb)
Set verbosity level.
Definition DMMoFEM.cpp:1517
PetscErrorCode DMMoFEMSwapDMCtx(DM dm, DM dm_swap)
Swap internal data struture.
Definition DMMoFEM.cpp:196
PetscErrorCode DMMoFEMDuplicateDMCtx(DM dm, DM dm_duplicate)
Duplicate internal data struture.
Definition DMMoFEM.cpp:180
SchurShellMatData createBlockMat(DM dm, boost::shared_ptr< BlockStructure > data)
Create a Schur Mat object.
Definition Schur.cpp:1452
constexpr auto field_name
constexpr double g
Data structure to exchange data between mofem and User Loop Methods.
Managing BitRefLevels.
virtual MoFEMErrorCode problem_basic_method_preProcess(const Problem *problem_ptr, BasicMethod &method, int verb=DEFAULT_VERBOSITY)=0
Set data for BasicMethod.
virtual MPI_Comm & get_comm() const =0
PetscBool isSquareMatrix
true if rows equals to cols
Definition DMCtxImpl.hpp:37
std::vector< std::string > rowCompPrb
Definition DMCtxImpl.hpp:30
boost::shared_ptr< NestSchurData > nestedSchurDataPtr
Definition DMCtxImpl.hpp:49
PetscBool destroyProblem
If true destroy problem with DM.
Definition DMCtxImpl.hpp:38
const Problem * problemMainOfSubPtr
pointer to main problem to sub-problem
Definition DMCtxImpl.hpp:27
std::string problemName
Problem name.
Definition DMCtxImpl.hpp:45
PetscBool isSubDM
Definition DMCtxImpl.hpp:24
PetscBool isPartitioned
true if read mesh is on parts
Definition DMCtxImpl.hpp:36
std::map< std::string, boost::shared_ptr< Range > > mapTypeCol
Definition DMCtxImpl.hpp:33
const Problem * problemPtr
pointer to problem data structure
Definition DMCtxImpl.hpp:44
std::map< std::string, boost::shared_ptr< Range > > mapTypeRow
Definition DMCtxImpl.hpp:32
std::vector< std::string > colSubFields
Definition DMCtxImpl.hpp:26
PetscBool isCompDM
Definition DMCtxImpl.hpp:29
boost::shared_ptr< BlockStructure > blocMatDataPtr
Definition DMCtxImpl.hpp:48
Interface * mField_ptr
MoFEM interface.
Definition DMCtxImpl.hpp:41
int verbosity
verbosity
Definition DMCtxImpl.hpp:20
std::vector< std::string > rowSubFields
Definition DMCtxImpl.hpp:25
PetscBool isProblemBuild
True if problem is build.
Definition DMCtxImpl.hpp:39
std::vector< std::string > colCompPrb
Definition DMCtxImpl.hpp:31
PETSc Discrete Manager data structure.
Definition DMMoFEM.hpp:1015
boost::shared_ptr< SnesCtx > snesCtx
data structure SNES
Definition DMMoFEM.hpp:1023
boost::shared_ptr< TsCtx > tsCtx
data structure for TS solver
Definition DMMoFEM.hpp:1024
boost::shared_ptr< KspCtx > kspCtx
data structure KSP
Definition DMMoFEM.hpp:1022
MoFEMErrorCode query_interface(boost::typeindex::type_index type_index, UnknownInterface **iface) const
Definition DMMoFEM.cpp:19
Deprecated interface functions.
Data structure to exchange data between mofem and User Loop Methods on entities.
structure for User Loop Methods on finite elements
Section manager is used to create indexes and sections.
Definition ISManager.hpp:23
Interface for linear (KSP) solver.
Definition KspCtx.hpp:13
static boost::shared_ptr< SinkType > createSink(boost::shared_ptr< std::ostream > stream_ptr, std::string comm_filter)
Create a sink object.
static boost::shared_ptr< std::ostream > getStrmWorld()
Get the strm world object.
static boost::shared_ptr< std::ostream > getStrmSync()
Get the strm sync object.
static bool checkIfChannelExist(const std::string channel)
Check if channel exist.
static boost::shared_ptr< std::ostream > getStrmSelf()
Get the strm self object.
Matrix manager is used to build and partition problems.
MoFEMErrorCode createHybridL2MPIAIJ(const std::string problem_name, SmartPetscObj< Mat > &aij_ptr, int verb=QUIET)
Create a Hybrid MPIAij object.
keeps basic data about problem
DofIdx getNbLocalDofsRow() const
MoFEMErrorCode getNumberOfElementsByNameAndPart(MPI_Comm comm, const std::string name, PetscLayout *layout) const
Get number of finite elements by name on processors.
BitRefLevel getBitRefLevel() const
BitRefLevel getBitRefLevelMask() const
boost::shared_ptr< SubProblemData > & getSubData() const
Get main problem of sub-problem is.
DofIdx getNbGhostDofsRow() const
Problem manager is used to build and partition problems.
intrusive_ptr for managing petsc objects
Interface for nonlinear (SNES) solver.
Definition SnesCtx.hpp:13
Interface for Time Stepping (TS) solver.
Definition TsCtx.hpp:17
base class for all interface classes
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
Vector manager is used to create vectors \mofem_vectors.