v0.14.0
Loading...
Searching...
No Matches
Distributed mesh manager

Implementation of PETSc DM, managing interactions between mesh data structures and vectors and matrices. More...

Collaboration diagram for Distributed mesh manager:

Classes

struct  MoFEM::DMCtx
 PETSc Discrete Manager data structure. More...
 
struct  MoFEM::DMMGViaApproxOrdersCtx
 Structure for DM for multi-grid via approximation orders. More...
 

Functions

PetscErrorCode MoFEM::DMRegister_MoFEM (const char sname[])
 Register MoFEM problem.
 
PetscErrorCode MoFEM::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 MoFEM::DMMoFEMCreateSubDM (DM subdm, DM dm, const char problem_name[])
 Must be called by user to set Sub DM MoFEM data structures.
 
PetscErrorCode MoFEM::DMoFEMGetInterfacePtr (DM dm, MoFEM::Interface **m_field_ptr)
 Get pointer to MoFEM::Interface.
 
PetscErrorCode MoFEM::DMMoFEMGetProblemPtr (DM dm, const MoFEM::Problem **problem_ptr)
 Get pointer to problem data structure.
 
PetscErrorCode MoFEM::DMMoFEMSetSquareProblem (DM dm, PetscBool square_problem)
 set squared problem
 
PetscErrorCode MoFEM::DMMoFEMGetSquareProblem (DM dm, PetscBool *square_problem)
 get squared problem
 
PetscErrorCode MoFEM::DMMoFEMResolveSharedFiniteElements (DM dm, std::string fe_name)
 Resolve shared entities.
 
PetscErrorCode MoFEM::DMMoFEMGetProblemFiniteElementLayout (DM dm, std::string fe_name, PetscLayout *layout)
 Get finite elements layout in the problem.
 
PetscErrorCode MoFEM::DMMoFEMAddElement (DM dm, std::string fe_name)
 add element to dm
 
PetscErrorCode MoFEM::DMMoFEMAddElement (DM dm, std::vector< std::string > fe_name)
 add element to dm
 
PetscErrorCode MoFEM::DMMoFEMUnSetElement (DM dm, std::string fe_name)
 unset element from dm
 
PetscErrorCode MoFEM::DMoFEMMeshToLocalVector (DM dm, Vec l, InsertMode mode, ScatterMode scatter_mode)
 set local (or ghosted) vector values on mesh for partition only
 
PetscErrorCode MoFEM::DMoFEMMeshToGlobalVector (DM dm, Vec g, InsertMode mode, ScatterMode scatter_mode)
 set ghosted vector values on all existing mesh entities
 
PetscErrorCode MoFEM::DMoFEMPreProcessFiniteElements (DM dm, MoFEM::FEMethod *method)
 execute finite element method for each element in dm (problem)
 
PetscErrorCode MoFEM::DMoFEMPostProcessFiniteElements (DM dm, MoFEM::FEMethod *method)
 execute finite element method for each element in dm (problem)
 
PetscErrorCode MoFEM::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 MoFEM::DMoFEMLoopFiniteElementsUpAndLowRank (DM dm, const std::string fe_name, boost::shared_ptr< MoFEM::FEMethod > method, int low_rank, int up_rank, CacheTupleWeakPtr cache_ptr=CacheTupleSharedPtr())
 Executes FEMethod for finite elements in DM.
 
PetscErrorCode MoFEM::DMoFEMLoopFiniteElements (DM dm, const char fe_name[], MoFEM::FEMethod *method, CacheTupleWeakPtr cache_ptr=CacheTupleSharedPtr())
 Executes FEMethod for finite elements in DM.
 
PetscErrorCode MoFEM::DMoFEMLoopFiniteElements (DM dm, const std::string fe_name, boost::shared_ptr< MoFEM::FEMethod > method, CacheTupleWeakPtr cache_ptr=CacheTupleSharedPtr())
 Executes FEMethod for finite elements in DM.
 
PetscErrorCode MoFEM::DMoFEMLoopDofs (DM dm, const char field_name[], MoFEM::DofMethod *method)
 execute method for dofs on field in problem
 
PetscErrorCode MoFEM::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 MoFEM::DMMoFEMKSPSetComputeRHS (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 KSP right hand side evaluation function
 
PetscErrorCode MoFEM::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 MoFEM::DMMoFEMKSPSetComputeOperators (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 KSP operators and push mofem finite element methods.
 
PetscErrorCode MoFEM::DMMoFEMSNESSetFunction (DM dm, const char fe_name[], MoFEM::FEMethod *method, MoFEM::BasicMethod *pre_only, MoFEM::BasicMethod *post_only)
 set SNES residual evaluation function
 
PetscErrorCode MoFEM::DMMoFEMSNESSetFunction (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 SNES residual evaluation function
 
PetscErrorCode MoFEM::DMMoFEMSNESSetJacobian (DM dm, const char fe_name[], MoFEM::FEMethod *method, MoFEM::BasicMethod *pre_only, MoFEM::BasicMethod *post_only)
 set SNES Jacobian evaluation function
 
PetscErrorCode MoFEM::DMMoFEMSNESSetJacobian (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 SNES Jacobian evaluation function
 
PetscErrorCode MoFEM::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 MoFEM::DMMoFEMTSSetIFunction (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 MoFEM::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
 
PetscErrorCode MoFEM::DMMoFEMTSSetIJacobian (DM dm, const char fe_name[], MoFEM::FEMethod *method, MoFEM::BasicMethod *pre_only, MoFEM::BasicMethod *post_only)
 set TS Jacobian evaluation function
 
PetscErrorCode MoFEM::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 MoFEM::DMMoFEMTSSetI2Function (DM dm, const char fe_name[], MoFEM::FEMethod *method, MoFEM::BasicMethod *pre_only, MoFEM::BasicMethod *post_only)
 set TS implicit function evaluation function
 
PetscErrorCode MoFEM::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 MoFEM::DMMoFEMTSSetI2Jacobian (DM dm, const char fe_name[], MoFEM::FEMethod *method, MoFEM::BasicMethod *pre_only, MoFEM::BasicMethod *post_only)
 set TS Jacobian evaluation function
 
PetscErrorCode MoFEM::DMMoFEMGetKspCtx (DM dm, MoFEM::KspCtx **ksp_ctx)
 get MoFEM::KspCtx data structure
 
PetscErrorCode MoFEM::DMMoFEMGetKspCtx (DM dm, const boost::shared_ptr< MoFEM::KspCtx > &ksp_ctx)
 get MoFEM::KspCtx data structure
 
PetscErrorCode MoFEM::DMMoFEMSetKspCtx (DM dm, boost::shared_ptr< MoFEM::KspCtx > ksp_ctx)
 set MoFEM::KspCtx data structure
 
PetscErrorCode MoFEM::DMMoFEMGetSnesCtx (DM dm, MoFEM::SnesCtx **snes_ctx)
 get MoFEM::SnesCtx data structure
 
PetscErrorCode MoFEM::DMMoFEMGetSnesCtx (DM dm, const boost::shared_ptr< MoFEM::SnesCtx > &snes_ctx)
 get MoFEM::SnesCtx data structure
 
PetscErrorCode MoFEM::DMMoFEMSetSnesCtx (DM dm, boost::shared_ptr< MoFEM::SnesCtx > snes_ctx)
 Set MoFEM::SnesCtx data structure.
 
PetscErrorCode MoFEM::DMMoFEMGetTsCtx (DM dm, MoFEM::TsCtx **ts_ctx)
 get MoFEM::TsCtx data structure
 
PetscErrorCode MoFEM::DMMoFEMGetTsCtx (DM dm, const boost::shared_ptr< MoFEM::TsCtx > &ts_ctx)
 get MoFEM::TsCtx data structure
 
PetscErrorCode MoFEM::DMMoFEMSetTsCtx (DM dm, boost::shared_ptr< MoFEM::TsCtx > ts_ctx)
 Set MoFEM::TsCtx data structure.
 
PetscErrorCode MoFEM::DMMoFEMSetIsPartitioned (DM dm, PetscBool is_partitioned)
 
PetscErrorCode MoFEM::DMMoFEMGetIsPartitioned (DM dm, PetscBool *is_partitioned)
 
PetscErrorCode MoFEM::DMSetOperators_MoFEM (DM dm)
 Set operators for MoFEM dm.
 
PetscErrorCode MoFEM::DMCreate_MoFEM (DM dm)
 Create dm data structure with MoFEM data structure.
 
PetscErrorCode MoFEM::DMDestroy_MoFEM (DM dm)
 Destroys dm with MoFEM data structure.
 
PetscErrorCode MoFEM::DMCreateGlobalVector_MoFEM (DM dm, Vec *g)
 DMShellSetCreateGlobalVector.
 
PetscErrorCode MoFEM::DMCreateGlobalVector_MoFEM (DM dm, SmartPetscObj< Vec > &g_ptr)
 DMShellSetCreateGlobalVector.
 
PetscErrorCode MoFEM::DMCreateLocalVector_MoFEM (DM dm, Vec *l)
 DMShellSetCreateLocalVector.
 
PetscErrorCode MoFEM::DMCreateMatrix_MoFEM (DM dm, Mat *M)
 
PetscErrorCode MoFEM::DMCreateMatrix_MoFEM (DM dm, SmartPetscObj< Mat > &M)
 
PetscErrorCode MoFEM::DMSetFromOptions_MoFEM (DM dm)
 
PetscErrorCode MoFEM::DMSetUp_MoFEM (DM dm)
 
PetscErrorCode MoFEM::DMSubDMSetUp_MoFEM (DM subdm)
 
PetscErrorCode MoFEM::DMMoFEMAddSubFieldRow (DM dm, const char field_name[])
 
PetscErrorCode MoFEM::DMMoFEMAddSubFieldRow (DM dm, const char field_name[], boost::shared_ptr< Range > r_ptr)
 Add field to sub dm problem on rows.
 
PetscErrorCode MoFEM::DMMoFEMAddSubFieldCol (DM dm, const char field_name[])
 
PetscErrorCode MoFEM::DMMoFEMGetIsSubDM (DM dm, PetscBool *is_sub_dm)
 
PetscErrorCode MoFEM::DMMoFEMAddRowCompositeProblem (DM dm, const char prb_name[])
 Add problem to composite DM on row.
 
PetscErrorCode MoFEM::DMMoFEMAddColCompositeProblem (DM dm, const char prb_name[])
 Add problem to composite DM on col.
 
PetscErrorCode MoFEM::DMMoFEMGetIsCompDM (DM dm, PetscBool *is_comp_dm)
 Get if this DM is composite DM.
 
PetscErrorCode MoFEM::DMGlobalToLocalBegin_MoFEM (DM dm, Vec, InsertMode, Vec)
 
PetscErrorCode MoFEM::DMGlobalToLocalEnd_MoFEM (DM dm, Vec, InsertMode, Vec)
 
PetscErrorCode MoFEM::DMLocalToGlobalBegin_MoFEM (DM, Vec, InsertMode, Vec)
 
PetscErrorCode MoFEM::DMLocalToGlobalEnd_MoFEM (DM, Vec, InsertMode, Vec)
 
PetscErrorCode MoFEM::DMCreateFieldIS_MoFEM (DM dm, PetscInt *numFields, char ***fieldNames, IS **fields)
 
PetscErrorCode MoFEM::DMMoFEMGetFieldIS (DM dm, RowColData rc, const char field_name[], IS *is)
 get field is in the problem
 
MoFEMErrorCode MoFEM::DMMoFEMCreateBlockMat (DM dm, Mat *mat)
 Create block matrix.
 
MoFEMErrorCode MoFEM::DMMoFEMCreateBlockMat (DM dm, SmartPetscObj< Mat > &mat)
 Create block matrix.
 
auto MoFEM::createDMMatrix (DM dm)
 Get smart matrix from DM.
 
auto MoFEM::createDMVector (DM dm)
 Get smart vector from DM.
 
MoFEMErrorCode MoFEM::DMMGViaApproxOrdersPushBackCoarseningIS (DM, IS is, Mat A, bool create_sub_matrix, bool shell_sub_a)
 Push back coarsening level to MG via approximation orders.
 
MoFEMErrorCode MoFEM::DMMGViaApproxOrdersPopBackCoarseningIS (DM)
 Pop IS form MG via approximation orders.
 
MoFEMErrorCode MoFEM::DMMGViaApproxOrdersClearCoarseningIS (DM)
 Clear approximation orders.
 
MoFEMErrorCode MoFEM::DMRegister_MGViaApproxOrders (const char sname[])
 Register DM for Multi-Grid via approximation orders.
 
MoFEMErrorCode MoFEM::DMCreateMatrix_MGViaApproxOrders (DM dm, Mat *M)
 Create matrix for Multi-Grid via approximation orders.
 
MoFEMErrorCode MoFEM::DMCoarsen_MGViaApproxOrders (DM dm, MPI_Comm comm, DM *dmc)
 Coarsen DM.
 

Detailed Description

Implementation of PETSc DM, managing interactions between mesh data structures and vectors and matrices.

DM objects are used to manage communication between the algebraic structures in PETSc (Vec and Mat) and mesh data structures in PDE-based (or other) simulations.

DM is abstract interface, here is it particular implementation for MoFEM code.

Function Documentation

◆ createDMMatrix()

auto MoFEM::createDMMatrix ( DM dm)
inline

Get smart matrix from DM.

Examples
free_surface.cpp, level_set.cpp, and navier_stokes.cpp.

Definition at line 1056 of file DMMoFEM.hpp.

1056 {
1057 SmartPetscObj<Mat> a;
1059 CHKERRABORT(getCommFromPetscObject(reinterpret_cast<PetscObject>(dm)), ierr);
1060 return a;
1061};
constexpr double a
static PetscErrorCode ierr
PetscErrorCode DMCreateMatrix_MoFEM(DM dm, Mat *M)
Definition DMMoFEM.cpp:1197
MPI_Comm getCommFromPetscObject(PetscObject obj)
Get the Comm From Petsc Object object.

◆ createDMVector()

auto MoFEM::createDMVector ( DM dm)
inline

Get smart vector from DM.

Examples
Electrostatics.cpp, adolc_plasticity.cpp, approx_sphere.cpp, child_and_parent.cpp, dg_projection.cpp, dynamic_first_order_con_law.cpp, eigen_elastic.cpp, free_surface.cpp, hanging_node_approx.cpp, heat_method.cpp, helmholtz.cpp, higher_derivatives.cpp, level_set.cpp, mixed_poisson.cpp, navier_stokes.cpp, phase.cpp, plate.cpp, poisson_2d_dis_galerkin.cpp, poisson_2d_homogeneous.cpp, and shallow_wave.cpp.

Definition at line 1099 of file DMMoFEM.hpp.

1099 {
1100 SmartPetscObj<Vec> v;
1102 CHKERRABORT(getCommFromPetscObject(reinterpret_cast<PetscObject>(dm)), ierr);
1103 return v;
1104};
PetscErrorCode DMCreateGlobalVector_MoFEM(DM dm, Vec *g)
DMShellSetCreateGlobalVector.
Definition DMMoFEM.cpp:1167
const double v
phase velocity of light in medium (cm/ns)

◆ DMCoarsen_MGViaApproxOrders()

MoFEMErrorCode MoFEM::DMCoarsen_MGViaApproxOrders ( DM dm,
MPI_Comm comm,
DM * dmc )

Coarsen DM.

Not used directly by user

Parameters
dmDistributed mesh data structure
commCommunicator
dmcCoarse distributed mesh data structure
Returns
Error code

Definition at line 376 of file PCMGSetUpViaApproxOrders.cpp.

376 {
377 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
379 GET_DM_FIELD(dm);
380 CHKERR PetscObjectGetComm((PetscObject)dm, &comm);
381 CHKERR DMCreate(comm, dmc);
382 (*dmc)->data = dm->data;
383 DMType type;
384 CHKERR DMGetType(dm, &type);
385 CHKERR DMSetType(*dmc, type);
386 dm_field->coarsenDM = SmartPetscObj<DM>(*dmc, false);
387 PetscInfo1(dm, "Coarsen DMMGViaApproxOrders leveldown = %d\n", dm->leveldown);
389}
#define GET_DM_FIELD(DM)
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define CHKERR
Inline error check.

◆ DMCreate_MoFEM()

PetscErrorCode MoFEM::DMCreate_MoFEM ( DM dm)

Create dm data structure with MoFEM data structure.

Definition at line 71 of file DMMoFEM.cpp.

71 {
72 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
74 dm->data = new DMCtxImpl();
77}
PetscErrorCode DMSetOperators_MoFEM(DM dm)
Set operators for MoFEM dm.
Definition DMMoFEM.cpp:49

◆ DMCreateFieldIS_MoFEM()

PetscErrorCode MoFEM::DMCreateFieldIS_MoFEM ( DM dm,
PetscInt * numFields,
char *** fieldNames,
IS ** fields )

Creates a set of IS objects with the global indices of dofs for each field

Parameters
dmThe number of fields (or NULL if not requested)

Output:

Parameters
numFieldsThe number of fields (or NULL if not requested)
fieldNamesThe name for each field (or NULL if not requested)
fieldsThe global indices for each field (or NULL if not requested)
Returns
error code
Note
The user is responsible for freeing all requested arrays. In particular, every entry of names should be freed with PetscFree(), every entry of fields should be destroyed with ISDestroy(), and both arrays should be freed with PetscFree().

Definition at line 1469 of file DMMoFEM.cpp.

1470 {
1471 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1473
1474 if (numFields) {
1475 *numFields = 0;
1476 }
1477 if (fieldNames) {
1478 *fieldNames = NULL;
1479 }
1480 if (fields) {
1481 *fields = NULL;
1482 }
1483
1484 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
1485 auto fields_ptr = dm_field->mField_ptr->get_fields();
1486 Field_multiIndex::iterator fit, hi_fit;
1487 fit = fields_ptr->begin();
1488 hi_fit = fields_ptr->end();
1489 *numFields = std::distance(fit, hi_fit);
1490
1491 if (fieldNames) {
1492 CHKERR PetscMalloc1(*numFields, fieldNames);
1493 }
1494 if (fields) {
1495 CHKERR PetscMalloc1(*numFields, fields);
1496 }
1497
1498 for (int f = 0; fit != hi_fit; fit++, f++) {
1499 if (fieldNames) {
1500 CHKERR PetscStrallocpy(fit->get()->getName().c_str(),
1501 (char **)&(*fieldNames)[f]);
1502 }
1503 if (fields) {
1504 CHKERR dm_field->mField_ptr->getInterface<ISManager>()
1505 ->isCreateProblemFieldAndRank(
1506 dm_field->problemPtr->getName(), ROW, fit->get()->getName(), 0,
1507 fit->get()->getNbOfCoeffs(), &(*fields)[f]);
1508 }
1509 }
1510
1512}
@ ROW

◆ DMCreateGlobalVector_MoFEM() [1/2]

PetscErrorCode MoFEM::DMCreateGlobalVector_MoFEM ( DM dm,
SmartPetscObj< Vec > & g_ptr )

DMShellSetCreateGlobalVector.

sets the routine to create a global vector associated with the shell DM

Definition at line 1177 of file DMMoFEM.cpp.

1177 {
1178 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1180 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
1181 CHKERR dm_field->mField_ptr->getInterface<VecManager>()->vecCreateGhost(
1182 dm_field->problemName, COL, g_ptr);
1183 CHKERR VecSetDM(g_ptr, dm);
1185}
@ COL

◆ DMCreateGlobalVector_MoFEM() [2/2]

PetscErrorCode MoFEM::DMCreateGlobalVector_MoFEM ( DM dm,
Vec * g )

DMShellSetCreateGlobalVector.

sets the routine to create a global vector associated with the shell DM

Definition at line 1167 of file DMMoFEM.cpp.

1167 {
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);
1173 CHKERR VecSetDM(*g, dm);
1175}
constexpr double g

◆ DMCreateLocalVector_MoFEM()

PetscErrorCode MoFEM::DMCreateLocalVector_MoFEM ( DM dm,
Vec * l )

DMShellSetCreateLocalVector.

sets the routine to create a local vector associated with the shell DM

Definition at line 1187 of file DMMoFEM.cpp.

1187 {
1188 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1190 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
1191 CHKERR dm_field->mField_ptr->getInterface<VecManager>()->vecCreateSeq(
1192 dm_field->problemName, COL, l);
1193 CHKERR VecSetDM(*l, dm);
1195}
FTensor::Index< 'l', 3 > l

◆ DMCreateMatrix_MGViaApproxOrders()

MoFEMErrorCode MoFEM::DMCreateMatrix_MGViaApproxOrders ( DM dm,
Mat * M )

Create matrix for Multi-Grid via approximation orders.

Not used directly by user

Parameters
dmDistributed mesh data structure
MMatrix
Returns
Error code

Definition at line 344 of file PCMGSetUpViaApproxOrders.cpp.

344 {
345
346 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
348 GET_DM_FIELD(dm);
349
350 int leveldown = dm->leveldown;
351
352 if (dm_field->kspOperators.empty()) {
354 } else {
355 MPI_Comm comm;
356 CHKERR PetscObjectGetComm((PetscObject)dm, &comm);
357 if (dm_field->kspOperators.empty()) {
358 SETERRQ(comm, MOFEM_DATA_INCONSISTENCY,
359 "data inconsistency, operator can not be set");
360 }
361 if (static_cast<int>(dm_field->kspOperators.size()) < leveldown) {
362 SETERRQ(comm, MOFEM_DATA_INCONSISTENCY,
363 "data inconsistency, no IS for that level");
364 }
365 *M = dm_field->kspOperators[dm_field->kspOperators.size() - 1 - leveldown];
366 CHKERR PetscObjectReference((PetscObject)*M);
367 CHKERR MatSetDM(*M, dm);
368 }
369
370 PetscInfo1(dm, "Create Matrix DMMGViaApproxOrders leveldown = %d\n",
371 leveldown);
372
374}
@ MOFEM_DATA_INCONSISTENCY
Definition definitions.h:31
FTensor::Index< 'M', 3 > M

◆ DMCreateMatrix_MoFEM() [1/2]

PetscErrorCode MoFEM::DMCreateMatrix_MoFEM ( DM dm,
Mat * M )

DMShellSetCreateMatrix

sets the routine to create a matrix associated with the shell DM

Examples
eigen_elastic.cpp.

Definition at line 1197 of file DMMoFEM.cpp.

1197 {
1198 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1200 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
1201
1202 if (strcmp(dm->mattype, MATSHELL) == 0) {
1203
1204 if (dm_field->blocMatDataPtr) {
1207 } else {
1208 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
1209 "Matrix shell data not set, or matrix type not implemented");
1210 }
1211
1212 } else if (strcmp(dm->mattype, MATMPIAIJ) == 0) {
1213 CHKERR dm_field->mField_ptr->getInterface<MatrixManager>()
1214 ->createMPIAIJWithArrays<PetscGlobalIdx_mi_tag>(dm_field->problemName,
1215 M);
1216 } else if (strcmp(dm->mattype, MATAIJ) == 0) {
1217 CHKERR dm_field->mField_ptr->getInterface<MatrixManager>()
1218 ->createSeqAIJWithArrays<PetscLocalIdx_mi_tag>(dm_field->problemName,
1219 M);
1220 } else if (strcmp(dm->mattype, MATAIJCUSPARSE) == 0) {
1221 CHKERR dm_field->mField_ptr->getInterface<MatrixManager>()
1222 ->createMPIAIJCUSPARSEWithArrays<PetscGlobalIdx_mi_tag>(
1223 dm_field->problemName, M);
1224 } else if (strcmp(dm->mattype, MATSEQAIJCUSPARSE) == 0) {
1225 CHKERR dm_field->mField_ptr->getInterface<MatrixManager>()
1226 ->createSeqAIJCUSPARSEWithArrays<PetscLocalIdx_mi_tag>(
1227 dm_field->problemName, M);
1228 } else {
1229 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
1230 "Matrix type not implemented");
1231 }
1232 CHKERR MatSetDM(*M, dm);
1234}
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
@ MOFEM_NOT_IMPLEMENTED
Definition definitions.h:32
MoFEMErrorCode DMMoFEMCreateBlockMat(DM dm, Mat *mat)
Create block matrix.
Definition DMMoFEM.cpp:1551

◆ DMCreateMatrix_MoFEM() [2/2]

PetscErrorCode MoFEM::DMCreateMatrix_MoFEM ( DM dm,
SmartPetscObj< Mat > & M )

DMShellSetCreateMatrix

sets the routine to create a matrix associated with the shell DM

Definition at line 1236 of file DMMoFEM.cpp.

1236 {
1237 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1239 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
1240
1241 if (strcmp(dm->mattype, MATSHELL) == 0) {
1242 if (dm_field->blocMatDataPtr) {
1243 Mat mat_raw;
1244 CHKERR DMMoFEMCreateBlockMat(dm, &mat_raw);
1245 M = SmartPetscObj<Mat>(mat_raw);
1247 } else {
1248 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
1249 "Matrix shell data not set, or matrix type not implemented");
1250 }
1251 } else if (strcmp(dm->mattype, MATMPIAIJ) == 0) {
1252 CHKERR dm_field->mField_ptr->getInterface<MatrixManager>()
1253 ->createMPIAIJWithArrays<PetscGlobalIdx_mi_tag>(dm_field->problemName,
1254 M);
1255 } else if (strcmp(dm->mattype, MATAIJ) == 0) {
1256 CHKERR dm_field->mField_ptr->getInterface<MatrixManager>()
1257 ->createSeqAIJWithArrays<PetscLocalIdx_mi_tag>(dm_field->problemName,
1258 M);
1259 } else if (strcmp(dm->mattype, MATAIJCUSPARSE) == 0) {
1260 CHKERR dm_field->mField_ptr->getInterface<MatrixManager>()
1261 ->createMPIAIJCUSPARSEWithArrays<PetscGlobalIdx_mi_tag>(
1262 dm_field->problemName, M);
1263 } else if (strcmp(dm->mattype, MATSEQAIJCUSPARSE) == 0) {
1264 CHKERR dm_field->mField_ptr->getInterface<MatrixManager>()
1265 ->createSeqAIJCUSPARSEWithArrays<PetscLocalIdx_mi_tag>(
1266 dm_field->problemName, M);
1267 } else {
1268 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
1269 "Matrix type not implemented");
1270 }
1271 CHKERR MatSetDM(M, dm);
1273}

◆ DMDestroy_MoFEM()

PetscErrorCode MoFEM::DMDestroy_MoFEM ( DM dm)

Destroys dm with MoFEM data structure.

destroy the MoFEM structure

Definition at line 79 of file DMMoFEM.cpp.

79 {
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}
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
#define MOFEM_LOG(channel, severity)
Log.

◆ DMGlobalToLocalBegin_MoFEM()

PetscErrorCode MoFEM::DMGlobalToLocalBegin_MoFEM ( DM dm,
Vec g,
InsertMode mode,
Vec l )

DMShellSetGlobalToLocal

the routine that begins the global to local scatter

Definition at line 1397 of file DMMoFEM.cpp.

1398 {
1399 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1401 CHKERR VecGhostUpdateBegin(g, INSERT_VALUES, SCATTER_FORWARD);
1403}

◆ DMGlobalToLocalEnd_MoFEM()

PetscErrorCode MoFEM::DMGlobalToLocalEnd_MoFEM ( DM dm,
Vec g,
InsertMode mode,
Vec l )

DMShellSetGlobalToLocal

the routine that begins the global to local scatter

Definition at line 1405 of file DMMoFEM.cpp.

1405 {
1407 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1409
1410 CHKERR VecGhostUpdateEnd(g, INSERT_VALUES, SCATTER_FORWARD);
1411
1412 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
1413 int nb_dofs = dm_field->problemPtr->getNbLocalDofsRow();
1414 int nb_ghost = dm_field->problemPtr->getNbGhostDofsRow();
1415
1416 double *array_loc, *array_glob;
1417 CHKERR VecGetArray(l, &array_loc);
1418 CHKERR VecGetArray(g, &array_glob);
1419 switch (mode) {
1420 case INSERT_VALUES:
1421 cblas_dcopy(nb_dofs + nb_ghost, array_glob, 1, array_loc, 1);
1422 break;
1423 case ADD_VALUES:
1424 cblas_daxpy(nb_dofs + nb_ghost, 1, array_glob, 1, array_loc, 1);
1425 break;
1426 default:
1427 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "not implemented");
1428 }
1429 CHKERR VecRestoreArray(l, &array_loc);
1430 CHKERR VecRestoreArray(g, &array_glob);
1432}

◆ DMLocalToGlobalBegin_MoFEM()

PetscErrorCode MoFEM::DMLocalToGlobalBegin_MoFEM ( DM dm,
Vec l,
InsertMode mode,
Vec g )

DMShellSetLocalToGlobal

the routine that begins the local to global scatter

DMShellSetLocalToGlobal

the routine that ends the local to global scatter

Definition at line 1434 of file DMMoFEM.cpp.

1435 {
1436
1437 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1439
1440 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
1441 int nb_dofs = dm_field->problemPtr->getNbLocalDofsRow();
1442 int nb_ghost = dm_field->problemPtr->getNbGhostDofsRow();
1443
1444 double *array_loc, *array_glob;
1445 CHKERR VecGetArray(l, &array_loc);
1446 CHKERR VecGetArray(g, &array_glob);
1447 switch (mode) {
1448 case INSERT_VALUES:
1449 cblas_dcopy(nb_dofs + nb_ghost, array_loc, 1, array_glob, 1);
1450 break;
1451 case ADD_VALUES:
1452 cblas_daxpy(nb_dofs + nb_ghost, 1, array_loc, 1, array_glob, 1);
1453 break;
1454 default:
1455 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "not implemented");
1456 }
1457 CHKERR VecRestoreArray(l, &array_loc);
1458 CHKERR VecRestoreArray(g, &array_glob);
1459
1461}

◆ DMLocalToGlobalEnd_MoFEM()

PetscErrorCode MoFEM::DMLocalToGlobalEnd_MoFEM ( DM ,
Vec l,
InsertMode mode,
Vec g )

DMShellSetLocalToGlobal

the routine that ends the local to global scatter

Definition at line 1463 of file DMMoFEM.cpp.

1463 {
1464 //
1467}

◆ DMMGViaApproxOrdersClearCoarseningIS()

MoFEMErrorCode MoFEM::DMMGViaApproxOrdersClearCoarseningIS ( DM )

Clear approximation orders.

Parameters
DMdm
Returns
Error code

Definition at line 293 of file PCMGSetUpViaApproxOrders.cpp.

293 {
294 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
296 GET_DM_FIELD(dm);
297 CHKERR dm_field->destroyCoarseningIS();
298 PetscInfo(dm, "Clear DMs data structures\n");
300}

◆ DMMGViaApproxOrdersPopBackCoarseningIS()

MoFEMErrorCode MoFEM::DMMGViaApproxOrdersPopBackCoarseningIS ( DM )

Pop IS form MG via approximation orders.

Parameters
DMdm
Returns
error code

Definition at line 281 of file PCMGSetUpViaApproxOrders.cpp.

281 {
282 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
284 GET_DM_FIELD(dm);
285 if (dm_field->coarseningIS.back()) {
286 dm_field->coarseningIS.pop_back();
287 }
288 dm_field->kspOperators.pop_back();
289 PetscInfo(dm, "Pop back IS to DMMGViaApproxOrders\n");
291}

◆ DMMGViaApproxOrdersPushBackCoarseningIS()

MoFEMErrorCode MoFEM::DMMGViaApproxOrdersPushBackCoarseningIS ( DM ,
IS is,
Mat A,
bool create_sub_matrix,
bool shell_sub_a )

Push back coarsening level to MG via approximation orders.

Parameters
DMdiscrete manager
isPush back IS used for coarsening
AGet sub-matrix of A using is (that sets operators for coarsening levels)
subAReturning pointer to created sub matrix
subAIf true create sub matrix, otherwise in subA has to be valid pointer to subA
Returns
Error code

Definition at line 227 of file PCMGSetUpViaApproxOrders.cpp.

229 {
230 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
232 GET_DM_FIELD(dm);
233 dm_field->coarseningIS.emplace_back(SmartPetscObj<IS>(is, true));
234 dm_field->shellMatrixCtxPtr.push_back(new PCMGSubMatrixCtx(A, is));
235 if (is) {
236 auto is2 = SmartPetscObj<IS>(is, true);
237 if (dm_field->aO) {
238 is2 = isDuplicate(is);
239 CHKERR ISCopy(is, is2);
240 CHKERR AOApplicationToPetscIS(dm_field->aO, is2);
241 }
242 if (create_sub_matrix) {
243 if (shell_sub_a) {
244 int n, N;
245 CHKERR ISGetSize(is, &N);
246 CHKERR ISGetLocalSize(is, &n);
247 MPI_Comm comm;
248 CHKERR PetscObjectGetComm((PetscObject)A, &comm);
249 Mat sub_a_raw;
250 CHKERR MatCreateShell(comm, n, n, N, N,
251 &(dm_field->shellMatrixCtxPtr.back()),
252 &sub_a_raw);
253 CHKERR MatShellSetOperation(sub_a_raw, MATOP_MULT,
254 (void (*)(void))sub_mat_mult);
255 CHKERR MatShellSetOperation(sub_a_raw, MATOP_MULT_ADD,
256 (void (*)(void))sub_mat_mult_add);
257 CHKERR MatShellSetOperation(sub_a_raw, MATOP_SOR,
258 (void (*)(void))sub_mat_sor);
259 dm_field->kspOperators.emplace_back(
260 SmartPetscObj<Mat>(sub_a_raw, false));
261 } else {
262 Mat sub_a_raw;
263#if PETSC_VERSION_GE(3, 8, 0)
264 CHKERR MatCreateSubMatrix(A, is2, is2, MAT_INITIAL_MATRIX, &sub_a_raw);
265#else
266 CHKERR MatGetSubMatrix(A, is2, is2, MAT_INITIAL_MATRIX, &sub_a_raw);
267#endif
268 dm_field->kspOperators.emplace_back(
269 SmartPetscObj<Mat>(sub_a_raw, false));
270 }
271 } else {
272 dm_field->kspOperators.emplace_back(SmartPetscObj<Mat>(A, true));
273 }
274 } else {
275 SETERRQ(PETSC_COMM_WORLD, MOFEM_DATA_INCONSISTENCY, "data inconsistency");
276 }
277 PetscInfo(dm, "Push back IS to DMMGViaApproxOrders\n");
279}
const double n
refractive index of diffusive medium
auto isDuplicate(IS is)
const int N
Definition speed_test.cpp:3

◆ DMMoFEMAddColCompositeProblem()

PetscErrorCode MoFEM::DMMoFEMAddColCompositeProblem ( DM dm,
const char prb_name[] )

Add problem to composite DM on col.

This create block on col with DOFs from problem of given name

Parameters
dmthe DM object
prb_nameadd problem name
Returns
error code

Definition at line 385 of file DMMoFEM.cpp.

385 {
386 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
388 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
389 if (!dm->data) {
390 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
391 "data structure for MoFEM not yet created");
392 }
393 if (!dm_field->isCompDM) {
394 dm_field->isCompDM = PETSC_TRUE;
395 }
396 if (dm_field->isSquareMatrix) {
397 SETERRQ(PETSC_COMM_SELF, MOFEM_INVALID_DATA,
398 "No need to add problem on column when problem block structurally "
399 "symmetric");
400 }
401 dm_field->colCompPrb.push_back(prb_name);
403}
@ MOFEM_INVALID_DATA
Definition definitions.h:36

◆ DMMoFEMAddElement() [1/2]

PetscErrorCode MoFEM::DMMoFEMAddElement ( DM dm,
std::string fe_name )

add element to dm

Note
add_file is a collective, should be executed on all processors. Otherwise could lead to deadlock.
Examples
Electrostatics.cpp, free_surface.cpp, level_set.cpp, navier_stokes.cpp, and photon_diffusion.cpp.

Definition at line 497 of file DMMoFEM.cpp.

497 {
498 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
500 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
501 CHKERR dm_field->mField_ptr->modify_problem_add_finite_element(
502 dm_field->problemName, fe_name);
504}

◆ DMMoFEMAddElement() [2/2]

PetscErrorCode MoFEM::DMMoFEMAddElement ( DM dm,
std::vector< std::string > fe_name )

add element to dm

Note
add_file is a collective, should be executed on all processors. Otherwise could lead to deadlock.

Definition at line 506 of file DMMoFEM.cpp.

506 {
508 for (auto fe : fe_name) {
510 }
512}
PetscErrorCode DMMoFEMAddElement(DM dm, std::string fe_name)
add element to dm
Definition DMMoFEM.cpp:497

◆ DMMoFEMAddRowCompositeProblem()

PetscErrorCode MoFEM::DMMoFEMAddRowCompositeProblem ( DM dm,
const char prb_name[] )

Add problem to composite DM on row.

This create block on row with DOFs from problem of given name

Parameters
dmthe DM object
prb_nameadd problem name
Returns
error code

Definition at line 367 of file DMMoFEM.cpp.

367 {
368 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
370 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
371 if (!dm->data) {
372 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
373 "data structure for MoFEM not yet created");
374 }
375 if (!dm_field->isCompDM) {
376 dm_field->isCompDM = PETSC_TRUE;
377 }
378 dm_field->rowCompPrb.push_back(prb_name);
379 if (dm_field->isSquareMatrix) {
380 dm_field->colCompPrb.push_back(prb_name);
381 }
383}

◆ DMMoFEMAddSubFieldCol()

PetscErrorCode MoFEM::DMMoFEMAddSubFieldCol ( DM dm,
const char field_name[] )

Add field to sub dm problem on columns

Examples
free_surface.cpp.

Definition at line 280 of file DMMoFEM.cpp.

280 {
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}
constexpr auto field_name

◆ DMMoFEMAddSubFieldRow() [1/2]

PetscErrorCode MoFEM::DMMoFEMAddSubFieldRow ( DM dm,
const char field_name[] )

Add field to sub dm problem on rows

Examples
free_surface.cpp, and level_set.cpp.

Definition at line 238 of file DMMoFEM.cpp.

238 {
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}

◆ DMMoFEMAddSubFieldRow() [2/2]

PetscErrorCode MoFEM::DMMoFEMAddSubFieldRow ( DM dm,
const char field_name[],
boost::shared_ptr< Range > r_ptr )

Add field to sub dm problem on rows.

Parameters
dm
field_name
m
Returns
PetscErrorCode

Definition at line 258 of file DMMoFEM.cpp.

259 {
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}

◆ DMMoFEMCreateBlockMat() [1/2]

MoFEMErrorCode MoFEM::DMMoFEMCreateBlockMat ( DM dm,
Mat * mat )

Create block matrix.

Parameters
dm
mat
Returns
MoFEMErrorCode

Definition at line 1551 of file DMMoFEM.cpp.

1551 {
1552 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1554 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
1555 auto mat_data = createBlockMat(dm, dm_field->blocMatDataPtr);
1556 *mat = mat_data.first;
1557 CHKERR PetscObjectReference((PetscObject)(*mat));
1559}
SchurShellMatData createBlockMat(DM dm, boost::shared_ptr< BlockStructure > data)
Create a Schur Mat object.
Definition Schur.cpp:1379

◆ DMMoFEMCreateBlockMat() [2/2]

MoFEMErrorCode MoFEM::DMMoFEMCreateBlockMat ( DM dm,
SmartPetscObj< Mat > & mat )

Create block matrix.

Parameters
dm
matsmart pointer
Returns
MoFEMErrorCode

◆ DMMoFEMCreateMoFEM()

PetscErrorCode MoFEM::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.

Note
If problem exist function create DM for it. If you set bit levels, those bits are to existing bits. Thus if you do not like to change bit ref level for existing problem, set bits to zero.
Examples
free_surface.cpp, and navier_stokes.cpp.

Definition at line 114 of file DMMoFEM.cpp.

117 {
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 }
140 CHKERR dm_field->mField_ptr->modify_problem_ref_level_add_bit(
141 dm_field->problemName, bit_level);
142 CHKERR dm_field->mField_ptr->modify_problem_mask_ref_level_add_bit(
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}
@ MF_EXCL
virtual bool check_problem(const std::string name)=0
check if problem exist
virtual MPI_Comm & get_comm() const =0

◆ DMMoFEMCreateSubDM()

PetscErrorCode MoFEM::DMMoFEMCreateSubDM ( DM subdm,
DM dm,
const char problem_name[] )

Must be called by user to set Sub DM MoFEM data structures.

Examples
free_surface.cpp, and level_set.cpp.

Definition at line 215 of file DMMoFEM.cpp.

215 {
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}
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 DMSubDMSetUp_MoFEM(DM subdm)
Definition DMMoFEM.cpp:1350

◆ DMMoFEMGetFieldIS()

PetscErrorCode MoFEM::DMMoFEMGetFieldIS ( DM dm,
RowColData rc,
const char field_name[],
IS * is )

get field is in the problem

Parameters
dmthe DM objects
rcROW or COL (no difference is problem is squared)
field_namename of the field
isreturned the IS object
Returns
error code
IS is;
ierr = DMMoFEMGetFieldIS(dm,ROW,"DISP",&is_disp); CHKERRG(ierr);
#define CHKERRG(n)
Check error code of MoFEM/MOAB/PETSc function.
PetscErrorCode DMMoFEMGetFieldIS(DM dm, RowColData rc, const char field_name[], IS *is)
get field is in the problem
Definition DMMoFEM.cpp:1514
static MoFEMErrorCodeGeneric< PetscErrorCode > ierr

Definition at line 1514 of file DMMoFEM.cpp.

1515 {
1516 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1518 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
1519 CHKERR dm_field->mField_ptr->getInterface<ISManager>()
1520 ->isCreateProblemFieldAndRank(dm_field->problemPtr->getName(), ROW,
1521 field_name, 0, 1000, is);
1523}

◆ DMMoFEMGetIsCompDM()

PetscErrorCode MoFEM::DMMoFEMGetIsCompDM ( DM dm,
PetscBool * is_comp_dm )

Get if this DM is composite DM.

Parameters
dmthe DM object
is_comp_dmreturn true if composite problem here
Returns
error code

Definition at line 405 of file DMMoFEM.cpp.

405 {
407 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
409 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
410 *is_comp_dm = dm_field->isCompDM;
412}

◆ DMMoFEMGetIsPartitioned()

PetscErrorCode MoFEM::DMMoFEMGetIsPartitioned ( DM dm,
PetscBool * is_partitioned )

get if read mesh is partitioned

Definition at line 1134 of file DMMoFEM.cpp.

1134 {
1135 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1137 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
1138 *is_partitioned = dm_field->isPartitioned;
1140}

◆ DMMoFEMGetIsSubDM()

PetscErrorCode MoFEM::DMMoFEMGetIsSubDM ( DM dm,
PetscBool * is_sub_dm )

Return true if this DM is sub problem

Parameters
dmthe DM object
is_subproblemtrue if subproblem
Returns
error code

Definition at line 322 of file DMMoFEM.cpp.

322 {
324 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
326 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
327 *is_sub_dm = dm_field->isSubDM;
329}

◆ DMMoFEMGetKspCtx() [1/2]

PetscErrorCode MoFEM::DMMoFEMGetKspCtx ( DM dm,
const boost::shared_ptr< MoFEM::KspCtx > & ksp_ctx )

get MoFEM::KspCtx data structure

Definition at line 1077 of file DMMoFEM.cpp.

1077 {
1078 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1080 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
1081 const_cast<boost::shared_ptr<MoFEM::KspCtx> &>(ksp_ctx) = dm_field->kspCtx;
1083}

◆ DMMoFEMGetKspCtx() [2/2]

PetscErrorCode MoFEM::DMMoFEMGetKspCtx ( DM dm,
MoFEM::KspCtx ** ksp_ctx )

get MoFEM::KspCtx data structure

Definition at line 1068 of file DMMoFEM.cpp.

1068 {
1069 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1071 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
1072 *ksp_ctx = dm_field->kspCtx.get();
1074}

◆ DMMoFEMGetProblemFiniteElementLayout()

PetscErrorCode MoFEM::DMMoFEMGetProblemFiniteElementLayout ( DM dm,
std::string fe_name,
PetscLayout * layout )

Get finite elements layout in the problem.

In layout is stored information how many elements is on each processor, for more information look int petsc documentation http://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/IS/PetscLayoutCreate.html#PetscLayoutCreate

Parameters
dmdiscrete manager for this problem
fe_namefinite element name
layoutpointer to layout, for created layout user takes responsibility for destroying it.
Returns
error code

Definition at line 473 of file DMMoFEM.cpp.

474 {
476 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
478 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
479
480 MPI_Comm comm;
481 CHKERR PetscObjectGetComm((PetscObject)dm, &comm);
482 CHKERR dm_field->problemPtr->getNumberOfElementsByNameAndPart(comm, fe_name,
483 layout);
485}

◆ DMMoFEMGetProblemPtr()

PetscErrorCode MoFEM::DMMoFEMGetProblemPtr ( DM dm,
const MoFEM::Problem ** problem_ptr )

Get pointer to problem data structure.

Examples
poisson_2d_dis_galerkin.cpp.

Definition at line 426 of file DMMoFEM.cpp.

426 {
427 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
429 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
430 if (!dm->data) {
431 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
432 "data structure for MoFEM not yet created");
433 }
434 *problem_ptr = dm_field->problemPtr;
436}

◆ DMMoFEMGetSnesCtx() [1/2]

PetscErrorCode MoFEM::DMMoFEMGetSnesCtx ( DM dm,
const boost::shared_ptr< MoFEM::SnesCtx > & snes_ctx )

get MoFEM::SnesCtx data structure

Definition at line 1103 of file DMMoFEM.cpp.

1103 {
1104 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1106 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
1107 const_cast<boost::shared_ptr<MoFEM::SnesCtx> &>(snes_ctx) = dm_field->snesCtx;
1109}

◆ DMMoFEMGetSnesCtx() [2/2]

PetscErrorCode MoFEM::DMMoFEMGetSnesCtx ( DM dm,
MoFEM::SnesCtx ** snes_ctx )

get MoFEM::SnesCtx data structure

Examples
navier_stokes.cpp.

Definition at line 1094 of file DMMoFEM.cpp.

1094 {
1095 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1097 DMCtxImpl *dm_field = (DMCtxImpl *)dm->data;
1098 *snes_ctx = dm_field->snesCtx.get();
1100}

◆ DMMoFEMGetSquareProblem()

PetscErrorCode MoFEM::DMMoFEMGetSquareProblem ( DM dm,
PetscBool * square_problem )

get squared problem

It if true is assumed that matrix has the same indexing on rows and columns. This reduces interprocessor communication.

Definition at line 487 of file DMMoFEM.cpp.

487 {
490 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
492 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
493 *square_problem = dm_field->isSquareMatrix;
495}

◆ DMMoFEMGetTsCtx() [1/2]

PetscErrorCode MoFEM::DMMoFEMGetTsCtx ( DM dm,
const boost::shared_ptr< MoFEM::TsCtx > & ts_ctx )

get MoFEM::TsCtx data structure

Definition at line 1150 of file DMMoFEM.cpp.

1151 {
1152 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1154 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
1155 const_cast<boost::shared_ptr<MoFEM::TsCtx> &>(ts_ctx) = dm_field->tsCtx;
1157}
MoFEM::TsCtx * ts_ctx

◆ DMMoFEMGetTsCtx() [2/2]

PetscErrorCode MoFEM::DMMoFEMGetTsCtx ( DM dm,
MoFEM::TsCtx ** ts_ctx )

get MoFEM::TsCtx data structure

Definition at line 1142 of file DMMoFEM.cpp.

1142 {
1143 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1145 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
1146 *ts_ctx = dm_field->tsCtx.get();
1148}

◆ DMMoFEMKSPSetComputeOperators() [1/2]

PetscErrorCode MoFEM::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.

Parameters
dmDM
fe_namefinite element name
methodmethod on the element (executed for each element in the problem which given name)
pre_onlymethod for pre-process before element method
post_onlymethod for post-process after element method
Returns
error code
Examples
level_set.cpp.

Definition at line 678 of file DMMoFEM.cpp.

681 {
682 return DMMoFEMKSPSetComputeOperators<const char *, MoFEM::FEMethod *,
685 dm, fe_name, method, pre_only, post_only);
686}
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:678
Data structure to exchange data between mofem and User Loop Methods.
structure for User Loop Methods on finite elements

◆ DMMoFEMKSPSetComputeOperators() [2/2]

PetscErrorCode MoFEM::DMMoFEMKSPSetComputeOperators ( 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 KSP operators and push mofem finite element methods.

Parameters
dmDM
fe_namefinite element name
methodmethod on the element (executed for each element in the problem which given name)
pre_onlymethod for pre-process before element method
post_onlymethod for post-process after element method
Returns
error code

Definition at line 689 of file DMMoFEM.cpp.

692 {
693 return DMMoFEMKSPSetComputeOperators<const std::string,
694 boost::shared_ptr<MoFEM::FEMethod>>(
695 dm, fe_name, method, pre_only, post_only);
696}

◆ DMMoFEMKSPSetComputeRHS() [1/2]

PetscErrorCode MoFEM::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.

Parameters
dmDM
Returns
error code

set KSP right hand side evaluation function

Examples
Electrostatics.cpp, and level_set.cpp.

Definition at line 637 of file DMMoFEM.cpp.

640 {
641 return DMMoFEMKSPSetComputeRHS<const char *, MoFEM::FEMethod *,
643 dm, fe_name, method, pre_only, post_only);
644}
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:637

◆ DMMoFEMKSPSetComputeRHS() [2/2]

PetscErrorCode MoFEM::DMMoFEMKSPSetComputeRHS ( 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 KSP right hand side evaluation function

Definition at line 647 of file DMMoFEM.cpp.

650 {
651 return DMMoFEMKSPSetComputeRHS<const std::string,
652 boost::shared_ptr<MoFEM::FEMethod>,
653 boost::shared_ptr<MoFEM::BasicMethod>,
654 boost::shared_ptr<MoFEM::BasicMethod>>(
655 dm, fe_name, method, pre_only, post_only);
656}

◆ DMMoFEMResolveSharedFiniteElements()

PetscErrorCode MoFEM::DMMoFEMResolveSharedFiniteElements ( DM dm,
std::string fe_name )

Resolve shared entities.

Parameters
dmdm
fe_namefinite element for which shared entities are resolved
Returns
error code
Note
This function is valid for parallel algebra and serial mesh. It should be run collectively, i.e. on all processors.

This allows for tag reduction or tag exchange, f.e.

CHKERR mField.get_moab().tag_get_handle("ADAPT_ORDER",th);
ParallelComm* pcomm =
ParallelComm::get_pcomm(&mField.get_moab(),MYPCOMM_INDEX);
// CHKERR pcomm->reduce_tags(th,MPI_SUM,prisms);
CHKERR pcomm->exchange_tags(th,prisms);
#define MYPCOMM_INDEX
default communicator number PCOMM
PetscErrorCode DMMoFEMResolveSharedFiniteElements(DM dm, std::string fe_name)
Resolve shared entities.
Definition DMMoFEM.cpp:464

Definition at line 464 of file DMMoFEM.cpp.

464 {
465 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
467 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
468 CHKERR dm_field->mField_ptr->getInterface<CommInterface>()
469 ->resolveSharedFiniteElements(dm_field->problemPtr, fe_name);
471}

◆ DMMoFEMSetIsPartitioned()

PetscErrorCode MoFEM::DMMoFEMSetIsPartitioned ( DM dm,
PetscBool is_partitioned )

sets if read mesh is partitioned

get if read mesh is partitioned

Examples
navier_stokes.cpp.

Definition at line 1123 of file DMMoFEM.cpp.

1123 {
1124 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1126 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
1127 dm_field->isPartitioned = is_partitioned;
1129}

◆ DMMoFEMSetKspCtx()

PetscErrorCode MoFEM::DMMoFEMSetKspCtx ( DM dm,
boost::shared_ptr< MoFEM::KspCtx > ksp_ctx )

set MoFEM::KspCtx data structure

Definition at line 1085 of file DMMoFEM.cpp.

1086 {
1087 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1089 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
1090 dm_field->kspCtx = ksp_ctx;
1092}

◆ DMMoFEMSetSnesCtx()

PetscErrorCode MoFEM::DMMoFEMSetSnesCtx ( DM dm,
boost::shared_ptr< MoFEM::SnesCtx > snes_ctx )

Set MoFEM::SnesCtx data structure.

Definition at line 1111 of file DMMoFEM.cpp.

1112 {
1113 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1115 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
1116 dm_field->snesCtx = snes_ctx;
1118}

◆ DMMoFEMSetSquareProblem()

PetscErrorCode MoFEM::DMMoFEMSetSquareProblem ( DM dm,
PetscBool square_problem )

set squared problem

It if true is assumed that matrix has the same indexing on rows and columns. This reduces interprocessor communication.

Examples
free_surface.cpp, and level_set.cpp.

Definition at line 456 of file DMMoFEM.cpp.

456 {
457 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
459 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
460 dm_field->isSquareMatrix = square_problem;
462}

◆ DMMoFEMSetTsCtx()

PetscErrorCode MoFEM::DMMoFEMSetTsCtx ( DM dm,
boost::shared_ptr< MoFEM::TsCtx > ts_ctx )

Set MoFEM::TsCtx data structure.

It take over pointer, do not delete it, DM will destroy pointer when is destroyed.

Definition at line 1159 of file DMMoFEM.cpp.

1159 {
1160 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1162 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
1163 dm_field->tsCtx = ts_ctx;
1165}

◆ DMMoFEMSNESSetFunction() [1/2]

PetscErrorCode MoFEM::DMMoFEMSNESSetFunction ( DM dm,
const char fe_name[],
MoFEM::FEMethod * method,
MoFEM::BasicMethod * pre_only,
MoFEM::BasicMethod * post_only )

set SNES residual evaluation function

Examples
navier_stokes.cpp.

Definition at line 718 of file DMMoFEM.cpp.

721 {
722 return DMMoFEMSNESSetFunction<const char *, MoFEM::FEMethod *,
724 dm, fe_name, method, pre_only, post_only);
725}
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:718

◆ DMMoFEMSNESSetFunction() [2/2]

PetscErrorCode MoFEM::DMMoFEMSNESSetFunction ( 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 SNES residual evaluation function

Definition at line 728 of file DMMoFEM.cpp.

731 {
732 return DMMoFEMSNESSetFunction<const std::string,
733 boost::shared_ptr<MoFEM::FEMethod>,
734 boost::shared_ptr<MoFEM::BasicMethod>,
735 boost::shared_ptr<MoFEM::BasicMethod>>(
736 dm, fe_name, method, pre_only, post_only);
737}

◆ DMMoFEMSNESSetJacobian() [1/2]

PetscErrorCode MoFEM::DMMoFEMSNESSetJacobian ( DM dm,
const char fe_name[],
MoFEM::FEMethod * method,
MoFEM::BasicMethod * pre_only,
MoFEM::BasicMethod * post_only )

set SNES Jacobian evaluation function

Examples
navier_stokes.cpp.

Definition at line 759 of file DMMoFEM.cpp.

762 {
763 return DMMoFEMSNESSetJacobian<const char *, MoFEM::FEMethod *,
765 dm, fe_name, method, pre_only, post_only);
766}
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:759

◆ DMMoFEMSNESSetJacobian() [2/2]

PetscErrorCode MoFEM::DMMoFEMSNESSetJacobian ( 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 SNES Jacobian evaluation function

Definition at line 769 of file DMMoFEM.cpp.

772 {
773 return DMMoFEMSNESSetJacobian<const std::string,
774 boost::shared_ptr<MoFEM::FEMethod>,
775 boost::shared_ptr<MoFEM::BasicMethod>,
776 boost::shared_ptr<MoFEM::BasicMethod>>(
777 dm, fe_name, method, pre_only, post_only);
778}

◆ DMMoFEMTSSetI2Function() [1/2]

PetscErrorCode MoFEM::DMMoFEMTSSetI2Function ( DM dm,
const char fe_name[],
MoFEM::FEMethod * method,
MoFEM::BasicMethod * pre_only,
MoFEM::BasicMethod * post_only )

set TS implicit function evaluation function

Definition at line 964 of file DMMoFEM.cpp.

967 {
968 return DMMoFEMTSSetI2Function<const char *, MoFEM::FEMethod *,
970 dm, fe_name, method, pre_only, post_only);
972}
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:975

◆ DMMoFEMTSSetI2Function() [2/2]

PetscErrorCode MoFEM::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 at line 975 of file DMMoFEM.cpp.

978 {
979 return DMMoFEMTSSetI2Function<const std::string,
980 boost::shared_ptr<MoFEM::FEMethod>,
981 boost::shared_ptr<MoFEM::BasicMethod>,
982 boost::shared_ptr<MoFEM::BasicMethod>>(
983 dm, fe_name, method, pre_only, post_only);
985}

◆ DMMoFEMTSSetI2Jacobian() [1/2]

PetscErrorCode MoFEM::DMMoFEMTSSetI2Jacobian ( DM dm,
const char fe_name[],
MoFEM::FEMethod * method,
MoFEM::BasicMethod * pre_only,
MoFEM::BasicMethod * post_only )

set TS Jacobian evaluation function

Definition at line 1007 of file DMMoFEM.cpp.

1010 {
1011 return DMMoFEMTSSetI2Jacobian<const char *, FEMethod *, MoFEM::BasicMethod *,
1012 MoFEM::BasicMethod *>(dm, fe_name, method,
1013 pre_only, post_only);
1014}
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:1017

◆ DMMoFEMTSSetI2Jacobian() [2/2]

PetscErrorCode MoFEM::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 at line 1017 of file DMMoFEM.cpp.

1020 {
1021 return DMMoFEMTSSetI2Jacobian<const std::string,
1022 boost::shared_ptr<MoFEM::FEMethod>,
1023 boost::shared_ptr<MoFEM::BasicMethod>,
1024 boost::shared_ptr<MoFEM::BasicMethod>>(
1025 dm, fe_name, method, pre_only, post_only);
1026}

◆ DMMoFEMTSSetIFunction() [1/2]

PetscErrorCode MoFEM::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 at line 800 of file DMMoFEM.cpp.

803 {
804 return DMMoFEMTSSetIFunction<const char *, MoFEM::FEMethod *,
806 dm, fe_name, method, pre_only, post_only);
808}
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:800

◆ DMMoFEMTSSetIFunction() [2/2]

PetscErrorCode MoFEM::DMMoFEMTSSetIFunction ( 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 at line 811 of file DMMoFEM.cpp.

814 {
815 return DMMoFEMTSSetIFunction<const std::string,
816 boost::shared_ptr<MoFEM::FEMethod>,
817 boost::shared_ptr<MoFEM::BasicMethod>,
818 boost::shared_ptr<MoFEM::BasicMethod>>(
819 dm, fe_name, method, pre_only, post_only);
821}

◆ DMMoFEMTSSetIJacobian() [1/2]

PetscErrorCode MoFEM::DMMoFEMTSSetIJacobian ( DM dm,
const char fe_name[],
MoFEM::FEMethod * method,
MoFEM::BasicMethod * pre_only,
MoFEM::BasicMethod * post_only )

set TS Jacobian evaluation function

Definition at line 843 of file DMMoFEM.cpp.

846 {
847 return DMMoFEMTSSetIJacobian<const char *, FEMethod *, MoFEM::BasicMethod *,
848 MoFEM::BasicMethod *>(dm, fe_name, method,
849 pre_only, post_only);
850}
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:853

◆ DMMoFEMTSSetIJacobian() [2/2]

PetscErrorCode MoFEM::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 at line 853 of file DMMoFEM.cpp.

856 {
857 return DMMoFEMTSSetIJacobian<const std::string,
858 boost::shared_ptr<MoFEM::FEMethod>,
859 boost::shared_ptr<MoFEM::BasicMethod>,
860 boost::shared_ptr<MoFEM::BasicMethod>>(
861 dm, fe_name, method, pre_only, post_only);
862}

◆ DMMoFEMUnSetElement()

PetscErrorCode MoFEM::DMMoFEMUnSetElement ( DM dm,
std::string fe_name )

unset element from dm

Definition at line 514 of file DMMoFEM.cpp.

514 {
515 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
517 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
518 CHKERR dm_field->mField_ptr->modify_problem_unset_finite_element(
519 dm_field->problemName, fe_name);
521}

◆ DMoFEMGetInterfacePtr()

PetscErrorCode MoFEM::DMoFEMGetInterfacePtr ( DM dm,
MoFEM::Interface ** m_field_ptr )

Get pointer to MoFEM::Interface.

Parameters
dmDistributed mesh manager
m_field_ptrPointer to pointer of field interface
Returns
Error code
Examples
adolc_plasticity.cpp, and free_surface.cpp.

Definition at line 414 of file DMMoFEM.cpp.

414 {
415 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
417 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
418 if (!dm->data) {
419 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
420 "data structure for MoFEM not yet created");
421 }
422 *m_field_ptr = dm_field->mField_ptr;
424}

◆ DMoFEMLoopDofs()

PetscErrorCode MoFEM::DMoFEMLoopDofs ( DM dm,
const char field_name[],
MoFEM::DofMethod * method )

execute method for dofs on field in problem

Definition at line 605 of file DMMoFEM.cpp.

606 {
607 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
609 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
610 ierr =
611 dm_field->mField_ptr->loop_dofs(dm_field->problemPtr, field_name, COL,
612 *method, dm_field->rAnk, dm_field->rAnk);
613 CHKERRG(ierr);
615}

◆ DMoFEMLoopFiniteElements() [1/2]

PetscErrorCode MoFEM::DMoFEMLoopFiniteElements ( DM dm,
const char fe_name[],
MoFEM::FEMethod * method,
CacheTupleWeakPtr cache_ptr = CacheTupleSharedPtr() )

Executes FEMethod for finite elements in DM.

Parameters
dmMoFEM discrete manager
fe_namename of element
methodpointer to MOFEM::FEMethod
Returns
Error code
Examples
adolc_plasticity.cpp, dynamic_first_order_con_law.cpp, free_surface.cpp, heat_equation.cpp, level_set.cpp, navier_stokes.cpp, photon_diffusion.cpp, plastic.cpp, poisson_2d_homogeneous.cpp, shallow_wave.cpp, and wave_equation.cpp.

Definition at line 586 of file DMMoFEM.cpp.

588 {
589 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
591 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
593 dm, fe_name, method, dm_field->rAnk, dm_field->rAnk, cache_ptr);
594 CHKERRG(ierr);
596}
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:567

◆ DMoFEMLoopFiniteElements() [2/2]

PetscErrorCode MoFEM::DMoFEMLoopFiniteElements ( DM dm,
const std::string fe_name,
boost::shared_ptr< MoFEM::FEMethod > method,
CacheTupleWeakPtr cache_ptr = CacheTupleSharedPtr() )

Executes FEMethod for finite elements in DM.

Parameters
dmMoFEM discrete manager
fe_namename of element
methodpointer to MOFEM::FEMethod
Returns
Error code

Definition at line 599 of file DMMoFEM.cpp.

601 {
602 return DMoFEMLoopFiniteElements(dm, fe_name.c_str(), method.get(), cache_ptr);
603}
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:586

◆ DMoFEMLoopFiniteElementsUpAndLowRank() [1/2]

PetscErrorCode MoFEM::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.

Parameters
dmMoFEM discrete manager
fe_namename of finite element
methodpointer to MoFEM::FEMethod
low_ranklowest rank of processor
up_rankupper run of processor
Returns
Error code
Examples
Electrostatics.cpp.

Definition at line 567 of file DMMoFEM.cpp.

569 {
571 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
572 ierr = dm_field->mField_ptr->loop_finite_elements(
573 dm_field->problemPtr, fe_name, *method, low_rank, up_rank, nullptr,
574 MF_EXIST, cache_ptr);
575 CHKERRG(ierr);
577}
@ MF_EXIST

◆ DMoFEMLoopFiniteElementsUpAndLowRank() [2/2]

PetscErrorCode MoFEM::DMoFEMLoopFiniteElementsUpAndLowRank ( DM dm,
const std::string fe_name,
boost::shared_ptr< MoFEM::FEMethod > method,
int low_rank,
int up_rank,
CacheTupleWeakPtr cache_ptr = CacheTupleSharedPtr() )

Executes FEMethod for finite elements in DM.

Parameters
dmMoFEM discrete manager
fe_namename of finite element
methodpointer to MoFEM::FEMethod
low_ranklowest rank of processor
up_rankupper run of processor
Returns
Error code

Definition at line 579 of file DMMoFEM.cpp.

581 {
582 return DMoFEMLoopFiniteElementsUpAndLowRank(dm, fe_name.c_str(), method.get(),
583 low_rank, up_rank, cache_ptr);
584}

◆ DMoFEMMeshToGlobalVector()

PetscErrorCode MoFEM::DMoFEMMeshToGlobalVector ( DM dm,
Vec g,
InsertMode mode,
ScatterMode scatter_mode )

set ghosted vector values on all existing mesh entities

Parameters
gvector
modesee petsc manual for VecSetValue (ADD_VALUES or INSERT_VALUES)
scatter_modesee petsc manual for ScatterMode (The available modes are: SCATTER_FORWARD or SCATTER_REVERSE)

SCATTER_REVERSE set data to field entities from V vector.

SCATTER_FORWARD set vector V from data field entities

Examples
navier_stokes.cpp.

Definition at line 535 of file DMMoFEM.cpp.

536 {
537 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
539 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
540 ierr = dm_field->mField_ptr->getInterface<VecManager>()->setGlobalGhostVector(
541 dm_field->problemPtr, COL, g, mode, scatter_mode);
542 CHKERRG(ierr);
544}

◆ DMoFEMMeshToLocalVector()

PetscErrorCode MoFEM::DMoFEMMeshToLocalVector ( DM dm,
Vec l,
InsertMode mode,
ScatterMode scatter_mode )

set local (or ghosted) vector values on mesh for partition only

Parameters
lvector
modesee petsc manual for VecSetValue (ADD_VALUES or INSERT_VALUES)
scatter_modesee petsc manual for ScatterMode (The available modes are: SCATTER_FORWARD or SCATTER_REVERSE)

SCATTER_REVERSE set data to field entities from V vector.

SCATTER_FORWARD set vector V from data field entities

Examples
Electrostatics.cpp, adolc_plasticity.cpp, approx_sphere.cpp, child_and_parent.cpp, dg_projection.cpp, dynamic_first_order_con_law.cpp, free_surface.cpp, hanging_node_approx.cpp, heat_method.cpp, helmholtz.cpp, higher_derivatives.cpp, level_set.cpp, mixed_poisson.cpp, phase.cpp, plate.cpp, poisson_2d_dis_galerkin.cpp, poisson_2d_homogeneous.cpp, and shallow_wave.cpp.

Definition at line 523 of file DMMoFEM.cpp.

524 {
526 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
528 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
529 ierr = dm_field->mField_ptr->getInterface<VecManager>()->setLocalGhostVector(
530 dm_field->problemPtr, COL, l, mode, scatter_mode);
531 CHKERRG(ierr);
533}

◆ DMoFEMPostProcessFiniteElements()

PetscErrorCode MoFEM::DMoFEMPostProcessFiniteElements ( DM dm,
MoFEM::FEMethod * method )

execute finite element method for each element in dm (problem)

Examples
free_surface.cpp.

Definition at line 556 of file DMMoFEM.cpp.

556 {
557 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
559 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
560 ierr = dm_field->mField_ptr->problem_basic_method_postProcess(
561 dm_field->problemPtr, *method);
562 CHKERRG(ierr);
564}

◆ DMoFEMPreProcessFiniteElements()

PetscErrorCode MoFEM::DMoFEMPreProcessFiniteElements ( DM dm,
MoFEM::FEMethod * method )

execute finite element method for each element in dm (problem)

Examples
free_surface.cpp, and navier_stokes.cpp.

Definition at line 546 of file DMMoFEM.cpp.

546 {
547 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
549 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
550 ierr = dm_field->mField_ptr->problem_basic_method_preProcess(
551 dm_field->problemPtr, *method);
552 CHKERRG(ierr);
554}

◆ DMRegister_MGViaApproxOrders()

MoFEMErrorCode MoFEM::DMRegister_MGViaApproxOrders ( const char sname[])

Register DM for Multi-Grid via approximation orders.

Parameters
snameproblem/dm registered name
Returns
error code

Definition at line 302 of file PCMGSetUpViaApproxOrders.cpp.

302 {
304 CHKERR DMRegister(sname, DMCreate_MGViaApproxOrders);
306}

◆ DMRegister_MoFEM()

PetscErrorCode MoFEM::DMRegister_MoFEM ( const char sname[])

Register MoFEM problem.

Examples
Electrostatics.cpp, mixed_poisson.cpp, navier_stokes.cpp, and photon_diffusion.cpp.

Definition at line 43 of file DMMoFEM.cpp.

43 {
45 CHKERR DMRegister(sname, DMCreate_MoFEM);
47}

◆ DMSetFromOptions_MoFEM()

PetscErrorCode MoFEM::DMSetFromOptions_MoFEM ( DM dm)

Set options for MoFEM DM

Definition at line 1281 of file DMMoFEM.cpp.

1281 {
1282#endif
1283
1284 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1286 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
1287#if PETSC_VERSION_GE(3, 5, 3)
1288 ierr = PetscOptionsHead(PetscOptionsObject, "DMMoFEM Options");
1289 CHKERRG(ierr);
1290#else
1291 ierr = PetscOptionsHead("DMMoFEM Options");
1292 CHKERRG(ierr);
1293#endif
1294 ierr = PetscOptionsBool("-dm_is_partitioned",
1295 "set if mesh is partitioned (works which native MOAB "
1296 "file format, i.e. h5m",
1297 "DMSetUp", dm_field->isPartitioned,
1298 &dm_field->isPartitioned, NULL);
1299 CHKERRG(ierr);
1301}

◆ DMSetOperators_MoFEM()

PetscErrorCode MoFEM::DMSetOperators_MoFEM ( DM dm)

Set operators for MoFEM dm.

Parameters
dm
Returns
error code

Definition at line 49 of file DMMoFEM.cpp.

49 {
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}
PetscErrorCode DMGlobalToLocalBegin_MoFEM(DM dm, Vec, InsertMode, Vec)
Definition DMMoFEM.cpp:1397
PetscErrorCode DMGlobalToLocalEnd_MoFEM(DM dm, Vec, InsertMode, Vec)
Definition DMMoFEM.cpp:1405
PetscErrorCode DMCreateLocalVector_MoFEM(DM dm, Vec *l)
DMShellSetCreateLocalVector.
Definition DMMoFEM.cpp:1187
PetscErrorCode DMDestroy_MoFEM(DM dm)
Destroys dm with MoFEM data structure.
Definition DMMoFEM.cpp:79
PetscErrorCode DMLocalToGlobalBegin_MoFEM(DM, Vec, InsertMode, Vec)
Definition DMMoFEM.cpp:1434
PetscErrorCode DMSetFromOptions_MoFEM(DM dm)
Definition DMMoFEM.cpp:1281
PetscErrorCode DMCreateFieldIS_MoFEM(DM dm, PetscInt *numFields, char ***fieldNames, IS **fields)
Definition DMMoFEM.cpp:1469
PetscErrorCode DMSetUp_MoFEM(DM dm)
Definition DMMoFEM.cpp:1303
PetscErrorCode DMLocalToGlobalEnd_MoFEM(DM, Vec, InsertMode, Vec)
Definition DMMoFEM.cpp:1463

◆ DMSetUp_MoFEM()

PetscErrorCode MoFEM::DMSetUp_MoFEM ( DM dm)

sets up the MoFEM structures inside a DM object

Examples
mixed_poisson.cpp.

Definition at line 1303 of file DMMoFEM.cpp.

1303 {
1304 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1305 ProblemsManager *prb_mng_ptr;
1307 DMCtxImpl *dm_field = static_cast<DMCtxImpl *>(dm->data);
1308 CHKERR dm_field->mField_ptr->getInterface(prb_mng_ptr);
1309
1310 if (dm_field->isCompDM) {
1311 // It is composite probelm
1312 CHKERR prb_mng_ptr->buildComposedProblem(
1313 dm_field->problemName, dm_field->rowCompPrb, dm_field->colCompPrb,
1314 dm_field->isSquareMatrix == PETSC_TRUE, dm_field->verbosity);
1315 } else {
1316 if (dm_field->isPartitioned) {
1317 CHKERR prb_mng_ptr->buildProblemOnDistributedMesh(
1318 dm_field->problemName, dm_field->isSquareMatrix == PETSC_TRUE,
1319 dm_field->verbosity);
1320 } else {
1321 CHKERR prb_mng_ptr->buildProblem(dm_field->problemName,
1322 dm_field->isSquareMatrix == PETSC_TRUE,
1323 dm_field->verbosity);
1324 CHKERR prb_mng_ptr->partitionProblem(dm_field->problemName,
1325 dm_field->verbosity);
1326 }
1327 }
1328
1329 // Partition finite elements
1330 if (dm_field->isPartitioned) {
1331 CHKERR prb_mng_ptr->partitionFiniteElements(
1332 dm_field->problemName, true, 0, dm_field->sIze, dm_field->verbosity);
1333 CHKERR prb_mng_ptr->partitionGhostDofsOnDistributedMesh(
1334 dm_field->problemName, dm_field->verbosity);
1335 } else {
1336 // partition finite elemnets
1337 CHKERR prb_mng_ptr->partitionFiniteElements(dm_field->problemName, false,
1338 -1, -1, dm_field->verbosity);
1339 // Get ghost DOFs
1340 CHKERR prb_mng_ptr->partitionGhostDofs(dm_field->problemName,
1341 dm_field->verbosity);
1342 }
1343
1344 // Set flag that problem is build and partitioned
1345 dm_field->isProblemBuild = PETSC_TRUE;
1346
1348}

◆ DMSubDMSetUp_MoFEM()

PetscErrorCode MoFEM::DMSubDMSetUp_MoFEM ( DM subdm)

Sets up the MoFEM structures inside a DM object for sub dm

Examples
level_set.cpp.

Definition at line 1350 of file DMMoFEM.cpp.

1350 {
1351 PetscValidHeaderSpecific(subdm, DM_CLASSID, 1);
1352 ProblemsManager *prb_mng_ptr;
1354
1355 DMCtxImpl *subdm_field = static_cast<DMCtxImpl *>(subdm->data);
1356
1357 // build sub dm problem
1358 CHKERR subdm_field->mField_ptr->getInterface(prb_mng_ptr);
1359
1360 map<std::string, boost::shared_ptr<Range>> *entity_map_row = nullptr;
1361 map<std::string, boost::shared_ptr<Range>> *entity_map_col = nullptr;
1362
1363 if (subdm_field->mapTypeRow.size())
1364 entity_map_row = &subdm_field->mapTypeRow;
1365 if (subdm_field->mapTypeCol.size())
1366 entity_map_col = &subdm_field->mapTypeCol;
1367
1368 CHKERR prb_mng_ptr->buildSubProblem(
1369 subdm_field->problemName, subdm_field->rowSubFields,
1370 subdm_field->colSubFields, subdm_field->problemMainOfSubPtr->getName(),
1371 subdm_field->isSquareMatrix == PETSC_TRUE, entity_map_row, entity_map_col,
1372 subdm_field->verbosity);
1373
1374 // partition problem
1375 subdm_field->isPartitioned = subdm_field->isPartitioned;
1376 if (subdm_field->isPartitioned) {
1377 CHKERR prb_mng_ptr->partitionFiniteElements(subdm_field->problemName, true,
1378 0, subdm_field->sIze,
1379 subdm_field->verbosity);
1380 // set ghost nodes
1381 CHKERR prb_mng_ptr->partitionGhostDofsOnDistributedMesh(
1382 subdm_field->problemName, subdm_field->verbosity);
1383 } else {
1384 // partition finite elements
1385 CHKERR prb_mng_ptr->partitionFiniteElements(subdm_field->problemName, false,
1386 -1, -1, subdm_field->verbosity);
1387 // set ghost nodes
1388 CHKERR prb_mng_ptr->partitionGhostDofs(subdm_field->problemName,
1389 subdm_field->verbosity);
1390 }
1391
1392 subdm_field->isProblemBuild = PETSC_TRUE;
1393
1395}