6#if PETSC_VERSION_GE(3, 6, 0)
7#include <petsc/private/dmimpl.h>
10#include <petsc-private/dmimpl.h>
11#include <petsc-private/vecimpl.h>
21 *iface =
const_cast<DMCtx *
>(
this);
27 auto core_log = logging::core::get();
50 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
66 CHKERR DMSetMatType(dm, MATMPIAIJ);
72 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
80 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
85 CHKERR PetscObjectGetComm((PetscObject)dm, &comm);
88 MPI_Comm_compare(comm, PETSC_COMM_SELF, &result);
89 if (result == MPI_IDENT)
91 <<
"MoFEM DM destroy for problem " << dm_field->
problemName
95 <<
"MoFEM DM destroy for problem " << dm_field->
problemName
106 delete static_cast<DMCtxImpl *
>(dm->data);
115 const char problem_name[],
123 "data structure for MoFEM not yet created");
127 "DM function not implemented into MoFEM");
145 boost::shared_ptr<KspCtx>(
new KspCtx(*m_field_ptr, problem_name));
147 boost::shared_ptr<SnesCtx>(
new SnesCtx(*m_field_ptr, problem_name));
149 boost::shared_ptr<TsCtx>(
new TsCtx(*m_field_ptr, problem_name));
152 CHKERR PetscObjectGetComm((PetscObject)dm, &comm);
154 MPI_Comm_compare(comm, m_field_ptr->
get_comm(), &result);
155 if (result > MPI_CONGRUENT) {
157 "MoFEM and DM using different communicators");
159 MPI_Comm_size(comm, &dm_field->
sIze);
160 MPI_Comm_rank(comm, &dm_field->
rAnk);
166 MPI_Comm_compare(comm, PETSC_COMM_SELF, &result);
167 if (result == MPI_IDENT) {
169 <<
"MoFEM DM created for problem " << dm_field->
problemName;
173 <<
"MoFEM DM created for problem " << dm_field->
problemName;
185 "data structure for MoFEM not yet created");
188 delete static_cast<DMCtxImpl *
>(dm_duplicate->data);
190 dm_duplicate->data = dm->data;
200 "data structure for MoFEM not yet created on dm");
203 "data structure for MoFEM not yet created on swap dm");
205 auto *dm_field =
static_cast<DMCtxImpl *
>(dm->data);
206 auto *dm_field_swap =
static_cast<DMCtxImpl *
>(dm_swap->data);
208 auto tmp_field = dm_field;
209 dm_field = dm_field_swap;
210 dm_field_swap = tmp_field;
220 "data structure for MoFEM not yet created");
229 subdm_field->
isSubDM = PETSC_TRUE;
239 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
244 "data structure for MoFEM not yet created");
259 boost::shared_ptr<Range> r_ptr) {
260 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
265 "data structure for MoFEM not yet created");
276 boost::shared_ptr<Range> r_ptr) {
281 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
286 "data structure for MoFEM not yet created");
301 boost::shared_ptr<Range> r_ptr) {
302 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
307 "data structure for MoFEM not yet created");
318 boost::shared_ptr<Range> r_ptr) {
323 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
326 *is_sub_dm = dm_field->
isSubDM;
331 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
334 if (dm_field->
isSubDM != PETSC_TRUE) {
336 "This DM is not created as a SubDM");
341 boost::shared_ptr<Problem::SubProblemData> sub_data =
343 CHKERR sub_data->getRowIs(is);
348 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
351 if (dm_field->
isSubDM != PETSC_TRUE) {
353 "This DM is not created as a SubDM");
358 boost::shared_ptr<Problem::SubProblemData> sub_data =
360 CHKERR sub_data->getColIs(is);
365 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
370 "data structure for MoFEM not yet created");
383 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
388 "data structure for MoFEM not yet created");
395 "No need to add problem on column when problem block structurally "
403 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
411 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
416 "data structure for MoFEM not yet created");
423 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
428 "data structure for MoFEM not yet created");
435 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
443 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
451 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
459 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
463 ->resolveSharedFiniteElements(dm_field->
problemPtr, fe_name);
468 PetscLayout *layout) {
469 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
474 CHKERR PetscObjectGetComm((PetscObject)dm, &comm);
481 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
489 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
499 for (
auto fe : fe_name) {
506 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
515 ScatterMode scatter_mode) {
516 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
526 ScatterMode scatter_mode) {
527 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
537 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
547 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
563 dm_field->
problemPtr, fe_name, *method, low_rank, up_rank,
nullptr,
570 DM dm,
const std::string fe_name, boost::shared_ptr<MoFEM::FEMethod> method,
573 low_rank, up_rank, cache_ptr);
579 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
583 dm, fe_name, method, dm_field->
rAnk, dm_field->
rAnk, cache_ptr);
590 boost::shared_ptr<MoFEM::FEMethod> method,
597 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
602 *method, dm_field->
rAnk, dm_field->
rAnk);
607template <
class S,
class T0,
class T1,
class T2>
609 T1 pre_only, T2 post_only) {
610 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
614 dm_field->
kspCtx->getPreProcComputeRhs().push_back(pre_only);
617 dm_field->
kspCtx->getComputeRhs().push_back(
621 dm_field->
kspCtx->getPostProcComputeRhs().push_back(post_only);
633 dm, fe_name, method, pre_only, post_only);
638 boost::shared_ptr<MoFEM::FEMethod> method,
639 boost::shared_ptr<MoFEM::BasicMethod> pre_only,
640 boost::shared_ptr<MoFEM::BasicMethod> post_only) {
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);
648template <
class S,
class T0,
class T1,
class T2>
650 T1 pre_only, T2 post_only) {
651 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
655 dm_field->
kspCtx->getPreProcSetOperators().push_back(pre_only);
658 dm_field->
kspCtx->getSetOperators().push_back(
662 dm_field->
kspCtx->getPostProcSetOperators().push_back(post_only);
675 dm, fe_name, method, pre_only, post_only);
680 boost::shared_ptr<MoFEM::FEMethod> method,
681 boost::shared_ptr<MoFEM::BasicMethod> pre_only,
682 boost::shared_ptr<MoFEM::BasicMethod> post_only) {
684 boost::shared_ptr<MoFEM::FEMethod>>(
685 dm, fe_name, method, pre_only, post_only);
688template <
class S,
class T0,
class T1,
class T2>
690 T1 pre_only, T2 post_only) {
691 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
695 dm_field->
snesCtx->getPreProcComputeRhs().push_back(pre_only);
698 dm_field->
snesCtx->getComputeRhs().push_back(
702 dm_field->
snesCtx->getPostProcComputeRhs().push_back(post_only);
714 dm, fe_name, method, pre_only, post_only);
719 boost::shared_ptr<MoFEM::FEMethod> method,
720 boost::shared_ptr<MoFEM::BasicMethod> pre_only,
721 boost::shared_ptr<MoFEM::BasicMethod> post_only) {
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);
729template <
class S,
class T0,
class T1,
class T2>
731 T1 pre_only, T2 post_only) {
732 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
736 dm_field->
snesCtx->getPreProcSetOperators().push_back(pre_only);
739 dm_field->
snesCtx->getSetOperators().push_back(
743 dm_field->
snesCtx->getPostProcSetOperators().push_back(post_only);
755 dm, fe_name, method, pre_only, post_only);
760 boost::shared_ptr<MoFEM::FEMethod> method,
761 boost::shared_ptr<MoFEM::BasicMethod> pre_only,
762 boost::shared_ptr<MoFEM::BasicMethod> post_only) {
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);
770template <
class S,
class T0,
class T1,
class T2>
772 T1 pre_only, T2 post_only) {
773 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
777 dm_field->
tsCtx->getPreProcessIFunction().push_back(pre_only);
780 dm_field->
tsCtx->getLoopsIFunction().push_back(
784 dm_field->
tsCtx->getPostProcessIFunction().push_back(post_only);
796 dm, fe_name, method, pre_only, post_only);
802 boost::shared_ptr<MoFEM::FEMethod> method,
803 boost::shared_ptr<MoFEM::BasicMethod> pre_only,
804 boost::shared_ptr<MoFEM::BasicMethod> post_only) {
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);
813template <
class S,
class T0,
class T1,
class T2>
815 T1 pre_only, T2 post_only) {
816 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
820 dm_field->
tsCtx->getPreProcessIJacobian().push_back(pre_only);
823 dm_field->
tsCtx->getLoopsIJacobian().push_back(
827 dm_field->
tsCtx->getPostProcessIJacobian().push_back(post_only);
839 pre_only, post_only);
844 boost::shared_ptr<MoFEM::FEMethod> method,
845 boost::shared_ptr<MoFEM::BasicMethod> pre_only,
846 boost::shared_ptr<MoFEM::BasicMethod> post_only) {
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);
854template <
class S,
class T0,
class T1,
class T2>
856 T1 pre_only, T2 post_only) {
857 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
861 dm_field->
tsCtx->getPreProcessRHSFunction().push_back(pre_only);
863 dm_field->
tsCtx->getLoopsRHSFunction().push_back(
866 dm_field->
tsCtx->getPostProcessRHSFunction().push_back(post_only);
873 boost::shared_ptr<MoFEM::FEMethod> method,
874 boost::shared_ptr<MoFEM::BasicMethod> pre_only,
875 boost::shared_ptr<MoFEM::BasicMethod> post_only) {
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);
890 dm, fe_name, method, pre_only, post_only);
894template <
class S,
class T0,
class T1,
class T2>
896 T1 pre_only, T2 post_only) {
897 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
901 dm_field->
tsCtx->getPreProcessRHSFunction().push_back(pre_only);
903 dm_field->
tsCtx->getLoopsRHSFunction().push_back(
906 dm_field->
tsCtx->getPostProcessRHSFunction().push_back(post_only);
913 boost::shared_ptr<MoFEM::FEMethod> method,
914 boost::shared_ptr<MoFEM::BasicMethod> pre_only,
915 boost::shared_ptr<MoFEM::BasicMethod> post_only) {
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);
930 dm, fe_name, method, pre_only, post_only);
934template <
class S,
class T0,
class T1,
class T2>
936 T1 pre_only, T2 post_only) {
937 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
941 dm_field->
tsCtx->getPreProcessIFunction().push_back(pre_only);
944 dm_field->
tsCtx->getLoopsIFunction().push_back(
948 dm_field->
tsCtx->getPostProcessIFunction().push_back(post_only);
960 dm, fe_name, method, pre_only, post_only);
966 boost::shared_ptr<MoFEM::FEMethod> method,
967 boost::shared_ptr<MoFEM::BasicMethod> pre_only,
968 boost::shared_ptr<MoFEM::BasicMethod> post_only) {
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);
977template <
class S,
class T0,
class T1,
class T2>
979 T1 pre_only, T2 post_only) {
980 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
984 dm_field->
tsCtx->getPreProcessIJacobian().push_back(pre_only);
987 dm_field->
tsCtx->getLoopsIJacobian().push_back(
991 dm_field->
tsCtx->getPostProcessIJacobian().push_back(post_only);
1003 pre_only, post_only);
1008 boost::shared_ptr<MoFEM::FEMethod> method,
1009 boost::shared_ptr<MoFEM::BasicMethod> pre_only,
1010 boost::shared_ptr<MoFEM::BasicMethod> post_only) {
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);
1018template <
class S,
class T0,
class T1,
class T2>
1020 T1 pre_only, T2 post_only) {
1021 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1025 dm_field->
tsCtx->getPreProcessMonitor().push_back(pre_only);
1027 dm_field->
tsCtx->getLoopsMonitor().push_back(
1030 dm_field->
tsCtx->getPostProcessMonitor().push_back(post_only);
1041 dm, ts, fe_name, method, pre_only, post_only);
1047 boost::shared_ptr<MoFEM::FEMethod> method,
1048 boost::shared_ptr<MoFEM::BasicMethod> pre_only,
1049 boost::shared_ptr<MoFEM::BasicMethod> post_only) {
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);
1059 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1062 *ksp_ctx = dm_field->
kspCtx.get();
1068 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1071 const_cast<boost::shared_ptr<MoFEM::KspCtx> &
>(ksp_ctx) = dm_field->
kspCtx;
1076 boost::shared_ptr<MoFEM::KspCtx> ksp_ctx) {
1077 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1080 dm_field->
kspCtx = ksp_ctx;
1085 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1088 *snes_ctx = dm_field->
snesCtx.get();
1094 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1097 const_cast<boost::shared_ptr<MoFEM::SnesCtx> &
>(snes_ctx) = dm_field->
snesCtx;
1102 boost::shared_ptr<MoFEM::SnesCtx> snes_ctx) {
1103 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1114 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1125 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1133 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1141 const boost::shared_ptr<MoFEM::TsCtx> &
ts_ctx) {
1142 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1145 const_cast<boost::shared_ptr<MoFEM::TsCtx> &
>(
ts_ctx) = dm_field->
tsCtx;
1150 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1158 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1168 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1173 CHKERR VecSetDM(g_ptr, dm);
1178 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1188 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1192 if (strcmp(dm->mattype, MATSHELL) == 0) {
1200 "Matrix shell data not set, or matrix type not implemented");
1203 }
else if (strcmp(dm->mattype, MATMPIAIJ) == 0) {
1205 ->createMPIAIJWithArrays<PetscGlobalIdx_mi_tag>(dm_field->
problemName,
1207 }
else if (strcmp(dm->mattype, MATAIJ) == 0) {
1209 ->createSeqAIJWithArrays<PetscLocalIdx_mi_tag>(dm_field->
problemName,
1211 }
else if (strcmp(dm->mattype, MATAIJCUSPARSE) == 0) {
1213 ->createMPIAIJCUSPARSEWithArrays<PetscGlobalIdx_mi_tag>(
1215 }
else if (strcmp(dm->mattype, MATSEQAIJCUSPARSE) == 0) {
1217 ->createSeqAIJCUSPARSEWithArrays<PetscLocalIdx_mi_tag>(
1221 "Matrix type not implemented");
1228 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1232 if (strcmp(dm->mattype, MATSHELL) == 0) {
1240 "Matrix shell data not set, or matrix type not implemented");
1242 }
else if (strcmp(dm->mattype, MATMPIAIJ) == 0) {
1244 ->createMPIAIJWithArrays<PetscGlobalIdx_mi_tag>(dm_field->
problemName,
1246 }
else if (strcmp(dm->mattype, MATAIJ) == 0) {
1248 ->createSeqAIJWithArrays<PetscLocalIdx_mi_tag>(dm_field->
problemName,
1250 }
else if (strcmp(dm->mattype, MATAIJCUSPARSE) == 0) {
1252 ->createMPIAIJCUSPARSEWithArrays<PetscGlobalIdx_mi_tag>(
1254 }
else if (strcmp(dm->mattype, MATSEQAIJCUSPARSE) == 0) {
1256 ->createSeqAIJCUSPARSEWithArrays<PetscLocalIdx_mi_tag>(
1260 "Matrix type not implemented");
1266#if PETSC_VERSION_GE(3, 17, 0)
1268 PetscOptionItems *PetscOptionsObject) {
1269#elif PETSC_VERSION_GE(3, 7, 0)
1272#elif PETSC_VERSION_GE(3, 5, 3)
1278 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
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");
1286 PetscOptionsHead(
"DMMoFEM Options");
1288 ierr = PetscOptionsBool(
"-dm_is_partitioned",
1289 "set if mesh is partitioned (works which native MOAB "
1290 "file format, i.e. h5m",
1298 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1345 PetscValidHeaderSpecific(subdm, DM_CLASSID, 1);
1354 map<std::string, boost::shared_ptr<Range>> *entity_map_row =
nullptr;
1355 map<std::string, boost::shared_ptr<Range>> *entity_map_col =
nullptr;
1365 subdm_field->
isSquareMatrix == PETSC_TRUE, entity_map_row, entity_map_col,
1372 0, subdm_field->
sIze,
1393 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1395 CHKERR VecGhostUpdateBegin(
g, INSERT_VALUES, SCATTER_FORWARD);
1400 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1403 CHKERR VecGhostUpdateEnd(
g, INSERT_VALUES, SCATTER_FORWARD);
1409 double *array_loc, *array_glob;
1410 CHKERR VecGetArray(
l, &array_loc);
1411 CHKERR VecGetArray(
g, &array_glob);
1414 cblas_dcopy(nb_dofs + nb_ghost, array_glob, 1, array_loc, 1);
1417 cblas_daxpy(nb_dofs + nb_ghost, 1, array_glob, 1, array_loc, 1);
1422 CHKERR VecRestoreArray(
l, &array_loc);
1423 CHKERR VecRestoreArray(
g, &array_glob);
1430 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1437 double *array_loc, *array_glob;
1438 CHKERR VecGetArray(
l, &array_loc);
1439 CHKERR VecGetArray(
g, &array_glob);
1442 cblas_dcopy(nb_dofs + nb_ghost, array_loc, 1, array_glob, 1);
1445 cblas_daxpy(nb_dofs + nb_ghost, 1, array_loc, 1, array_glob, 1);
1450 CHKERR VecRestoreArray(
l, &array_loc);
1451 CHKERR VecRestoreArray(
g, &array_glob);
1462 char ***fieldNames, IS **fields) {
1463 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
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);
1484 CHKERR PetscMalloc1(*numFields, fieldNames);
1487 CHKERR PetscMalloc1(*numFields, fields);
1490 for (
int f = 0; fit != hi_fit; fit++, f++) {
1492 CHKERR PetscStrallocpy(fit->get()->getName().c_str(),
1493 (
char **)&(*fieldNames)[f]);
1497 ->isCreateProblemFieldAndRank(
1499 fit->get()->getNbOfCoeffs(), &(*fields)[f]);
1508 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1518 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1526 boost::shared_ptr<BlockStructure> data) {
1527 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1535 boost::shared_ptr<BlockStructure> &data) {
1536 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1544 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1548 *mat = mat_data.first;
1549 CHKERR PetscObjectReference((PetscObject)(*mat));
1555 boost::shared_ptr<NestSchurData> data) {
1556 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1565 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1569 *mat = mat_data.first;
1570 CHKERR PetscObjectReference((PetscObject)(*mat));
1575 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
Implementation of DM context. You should not use it directly.
Discrete manager interface for MoFEM.
#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
#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)
PetscErrorCode DMGlobalToLocalEnd_MoFEM(DM dm, Vec, InsertMode, Vec)
PetscErrorCode DMMoFEMSetIsPartitioned(DM dm, PetscBool is_partitioned)
PetscErrorCode DMMoFEMCreateSubDM(DM subdm, DM dm, const char problem_name[])
Must be called by user to set Sub DM MoFEM data structures.
PetscErrorCode DMMoFEMGetIsPartitioned(DM dm, PetscBool *is_partitioned)
PetscErrorCode DMMoFEMAddElement(DM dm, std::string fe_name)
add element to dm
PetscErrorCode DMMoFEMSNESSetFunction(DM dm, const char fe_name[], MoFEM::FEMethod *method, MoFEM::BasicMethod *pre_only, MoFEM::BasicMethod *post_only)
set SNES residual evaluation function
PetscErrorCode DMMoFEMGetProblemFiniteElementLayout(DM dm, std::string fe_name, PetscLayout *layout)
Get finite elements layout in the problem.
PetscErrorCode DMCreateLocalVector_MoFEM(DM dm, Vec *l)
DMShellSetCreateLocalVector.
PetscErrorCode DMMoFEMSetSquareProblem(DM dm, PetscBool square_problem)
set squared problem
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
PetscErrorCode DMDestroy_MoFEM(DM dm)
Destroys dm with MoFEM data structure.
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.
PetscErrorCode DMoFEMPostProcessFiniteElements(DM dm, MoFEM::FEMethod *method)
execute finite element method for each element in dm (problem)
PetscErrorCode DMoFEMMeshToLocalVector(DM dm, Vec l, InsertMode mode, ScatterMode scatter_mode)
set local (or ghosted) vector values on mesh for partition only
PetscErrorCode DMMoFEMAddColCompositeProblem(DM dm, const char prb_name[])
Add problem to composite DM on col.
PetscErrorCode DMLocalToGlobalBegin_MoFEM(DM, Vec, InsertMode, Vec)
PetscErrorCode DMMoFEMGetProblemPtr(DM dm, const MoFEM::Problem **problem_ptr)
Get pointer to problem data structure.
PetscErrorCode DMCreate_MoFEM(DM dm)
Create dm data structure with MoFEM data structure.
PetscErrorCode DMMoFEMSetSnesCtx(DM dm, boost::shared_ptr< MoFEM::SnesCtx > snes_ctx)
Set MoFEM::SnesCtx data structure.
PetscErrorCode DMMoFEMSetKspCtx(DM dm, boost::shared_ptr< MoFEM::KspCtx > ksp_ctx)
set MoFEM::KspCtx data structure
PetscErrorCode DMMoFEMAddSubFieldRow(DM dm, const char field_name[])
PetscErrorCode DMMoFEMGetIsSubDM(DM dm, PetscBool *is_sub_dm)
PetscErrorCode DMMoFEMAddRowCompositeProblem(DM dm, const char prb_name[])
Add problem to composite DM on row.
PetscErrorCode DMMoFEMSNESSetJacobian(DM dm, const char fe_name[], MoFEM::FEMethod *method, MoFEM::BasicMethod *pre_only, MoFEM::BasicMethod *post_only)
set SNES Jacobian evaluation function
PetscErrorCode DMCreateMatrix_MoFEM(DM dm, Mat *M)
PetscErrorCode DMSetFromOptions_MoFEM(DM dm)
PetscErrorCode DMMoFEMGetTsCtx(DM dm, MoFEM::TsCtx **ts_ctx)
get MoFEM::TsCtx data structure
PetscErrorCode DMRegister_MoFEM(const char sname[])
Register MoFEM problem.
PetscErrorCode DMMoFEMGetKspCtx(DM dm, MoFEM::KspCtx **ksp_ctx)
get MoFEM::KspCtx data structure
PetscErrorCode DMMoFEMGetSquareProblem(DM dm, PetscBool *square_problem)
get squared problem
PetscErrorCode DMoFEMLoopDofs(DM dm, const char field_name[], MoFEM::DofMethod *method)
execute method for dofs on field in problem
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.
PetscErrorCode DMMoFEMResolveSharedFiniteElements(DM dm, std::string fe_name)
Resolve shared entities.
PetscErrorCode DMMoFEMSetTsCtx(DM dm, boost::shared_ptr< MoFEM::TsCtx > ts_ctx)
Set MoFEM::TsCtx data structure.
PetscErrorCode DMoFEMLoopFiniteElements(DM dm, const char fe_name[], MoFEM::FEMethod *method, CacheTupleWeakPtr cache_ptr=CacheTupleSharedPtr())
Executes FEMethod for finite elements in DM.
PetscErrorCode DMMoFEMGetFieldIS(DM dm, RowColData rc, const char field_name[], IS *is)
get field is in the problem
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
MoFEMErrorCode DMMoFEMCreateBlockMat(DM dm, Mat *mat)
Create block matrix.
PetscErrorCode DMoFEMGetInterfacePtr(DM dm, MoFEM::Interface **m_field_ptr)
Get pointer to MoFEM::Interface.
PetscErrorCode DMMoFEMAddSubFieldCol(DM dm, const char field_name[])
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
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
PetscErrorCode DMCreateGlobalVector_MoFEM(DM dm, Vec *g)
DMShellSetCreateGlobalVector.
PetscErrorCode DMMoFEMGetSnesCtx(DM dm, MoFEM::SnesCtx **snes_ctx)
get MoFEM::SnesCtx data structure
PetscErrorCode DMoFEMMeshToGlobalVector(DM dm, Vec g, InsertMode mode, ScatterMode scatter_mode)
set ghosted vector values on all existing mesh entities
PetscErrorCode DMMoFEMGetIsCompDM(DM dm, PetscBool *is_comp_dm)
Get if this DM is composite DM.
PetscErrorCode DMSubDMSetUp_MoFEM(DM subdm)
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.
PetscErrorCode DMCreateFieldIS_MoFEM(DM dm, PetscInt *numFields, char ***fieldNames, IS **fields)
PetscErrorCode DMMoFEMUnSetElement(DM dm, std::string fe_name)
unset element from dm
PetscErrorCode DMSetUp_MoFEM(DM dm)
PetscErrorCode DMLocalToGlobalEnd_MoFEM(DM, Vec, InsertMode, Vec)
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.
PetscErrorCode DMoFEMPreProcessFiniteElements(DM dm, MoFEM::FEMethod *method)
execute finite element method for each element in dm (problem)
PetscErrorCode DMSetOperators_MoFEM(DM dm)
Set operators for MoFEM dm.
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
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.
implementation of Data Operators for Forces and Sources
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.
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.
PetscErrorCode TsMonitorSet(TS ts, PetscInt step, PetscReal t, Vec u, void *ctx)
Set monitor for TS solver.
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
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.
PetscErrorCode DMMoFEMSetDestroyProblem(DM dm, PetscBool destroy_problem)
MoFEMErrorCode DMMoFEMCreateHybridL2Mat(DM dm, SmartPetscObj< Mat > &mat)
Create matrix for hybridised system.
PetscErrorCode DMMoFEMGetDestroyProblem(DM dm, PetscBool *destroy_problem)
MoFEMErrorCode DMMoFEMGetBlocMatData(DM dm, boost::shared_ptr< BlockStructure > &)
Get data for block mat.
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.
PetscErrorCode TsSetRHSFunction(TS ts, PetscReal t, Vec u, Vec F, void *ctx)
TS solver function.
std::pair< SmartPetscObj< Mat >, boost::shared_ptr< NestSchurData > > createSchurNestedMatrix(boost::shared_ptr< NestSchurData > schur_net_data_ptr)
Create a Mat Diag Blocks object.
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.
PetscErrorCode DMMoFEMGetSubRowIS(DM dm, IS *is)
get sub problem is
MoFEMErrorCode DMMoFEMCreateNestSchurMat(DM dm, Mat *mat)
Create nest schur matrix.
PetscErrorCode TsSetIFunction(TS ts, PetscReal t, Vec u, Vec u_t, Vec F, void *ctx)
Set IFunction for TS solver.
PetscErrorCode DMMoFEMGetSubColIS(DM dm, IS *is)
get sub problem is
PetscErrorCode KspRhs(KSP ksp, Vec f, void *ctx)
Run over elements in the lists.
MoFEMErrorCode DMMoFEMSetBlocMatData(DM dm, boost::shared_ptr< BlockStructure >)
Set data for block mat.
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.
PetscErrorCode KspMat(KSP ksp, Mat A, Mat B, void *ctx)
Run over elements in the list.
PetscErrorCode TsSetRHSJacobian(TS ts, PetscReal t, Vec u, Mat A, Mat B, void *ctx)
TS solver function.
MoFEMErrorCode DMMoFEMSetNestSchurData(DM dm, boost::shared_ptr< NestSchurData >)
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
PetscErrorCode DMMoFEMSetVerbosity(DM dm, const int verb)
Set verbosity level.
PetscErrorCode DMMoFEMSwapDMCtx(DM dm, DM dm_swap)
Swap internal data struture.
PetscErrorCode DMMoFEMDuplicateDMCtx(DM dm, DM dm_duplicate)
Duplicate internal data struture.
SchurShellMatData createBlockMat(DM dm, boost::shared_ptr< BlockStructure > data)
Create a Schur Mat object.
constexpr auto field_name
Data structure to exchange data between mofem and User Loop Methods.
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
std::vector< std::string > rowCompPrb
boost::shared_ptr< NestSchurData > nestedSchurDataPtr
PetscBool destroyProblem
If true destroy problem with DM.
const Problem * problemMainOfSubPtr
pointer to main problem to sub-problem
std::string problemName
Problem name.
PetscBool isPartitioned
true if read mesh is on parts
std::map< std::string, boost::shared_ptr< Range > > mapTypeCol
const Problem * problemPtr
pointer to problem data structure
std::map< std::string, boost::shared_ptr< Range > > mapTypeRow
std::vector< std::string > colSubFields
boost::shared_ptr< BlockStructure > blocMatDataPtr
Interface * mField_ptr
MoFEM interface.
std::vector< std::string > rowSubFields
PetscBool isProblemBuild
True if problem is build.
std::vector< std::string > colCompPrb
PETSc Discrete Manager data structure.
boost::shared_ptr< SnesCtx > snesCtx
data structure SNES
boost::shared_ptr< TsCtx > tsCtx
data structure for TS solver
boost::shared_ptr< KspCtx > kspCtx
data structure KSP
MoFEMErrorCode query_interface(boost::typeindex::type_index type_index, UnknownInterface **iface) const
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.
Interface for linear (KSP) solver.
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.
Interface for Time Stepping (TS) solver.
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.