v0.8.23
Classes | Functions
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  DMMGViaApproxOrdersCtx
 Structure for DM for multi-grid via approximation orders. More...
 

Functions

PetscErrorCode MoFEM::DMRegister_MoFEM (const char sname[])
 Register MoFEM problem. More...
 
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. More...
 
PetscErrorCode MoFEM::DMMoFEMCreateSubDM (DM subdm, DM dm, const char problem_name[])
 Must be called by user to set Sub DM MoFEM data structures. More...
 
PetscErrorCode MoFEM::DMoFEMGetInterfacePtr (DM dm, MoFEM::Interface **m_field_ptr)
 Get pointer to MoFEM::Interface. More...
 
PetscErrorCode MoFEM::DMMoFEMGetProblemPtr (DM dm, const MoFEM::Problem **problem_ptr)
 Get pointer to problem data structure. More...
 
PetscErrorCode MoFEM::DMMoFEMSetSquareProblem (DM dm, PetscBool square_problem)
 set squared problemIt if true is assumed that matrix has the same indexing on rows and columns. This reduces interprocessor communication. More...
 
PetscErrorCode MoFEM::DMMoFEMGetSquareProblem (DM dm, PetscBool *square_problem)
 get squared problemIt if true is assumed that matrix has the same indexing on rows and columns. This reduces interprocessor communication. More...
 
PetscErrorCode MoFEM::DMMoFEMResolveSharedEntities (DM dm, const char fe_name[])
 Resolve shared entities. More...
 
PetscErrorCode MoFEM::DMMoFEMGetProblemFiniteElementLayout (DM dm, const char fe_name[], PetscLayout *layout)
 Get finite elements layout in the problem. More...
 
PetscErrorCode MoFEM::DMMoFEMAddElement (DM dm, const char fe_name[])
 add element to dm More...
 
PetscErrorCode MoFEM::DMMoFEMUnSetElement (DM dm, const char fe_name[])
 unset element from dm More...
 
PetscErrorCode MoFEM::DMoFEMMeshToLocalVector (DM dm, Vec l, InsertMode mode, ScatterMode scatter_mode)
 set local (or ghosted) vector values on mesh for partition only More...
 
PetscErrorCode MoFEM::DMoFEMMeshToGlobalVector (DM dm, Vec g, InsertMode mode, ScatterMode scatter_mode)
 set ghosted vector values on all existing mesh entities More...
 
PetscErrorCode MoFEM::DMoFEMPreProcessFiniteElements (DM dm, MoFEM::FEMethod *method)
 execute finite element method for each element in dm (problem) More...
 
PetscErrorCode MoFEM::DMoFEMPostProcessFiniteElements (DM dm, MoFEM::FEMethod *method)
 execute finite element method for each element in dm (problem) More...
 
PetscErrorCode MoFEM::DMoFEMLoopFiniteElementsUpAndLowRank (DM dm, const char fe_name[], MoFEM::FEMethod *method, int low_rank, int up_rank)
 Executes FEMethod for finite elements in DM. More...
 
PetscErrorCode MoFEM::DMoFEMLoopFiniteElementsUpAndLowRank (DM dm, const std::string fe_name, boost::shared_ptr< MoFEM::FEMethod > method, int low_rank, int up_rank)
 Executes FEMethod for finite elements in DM. More...
 
PetscErrorCode MoFEM::DMoFEMLoopFiniteElements (DM dm, const char fe_name[], MoFEM::FEMethod *method)
 Executes FEMethod for finite elements in DM. More...
 
PetscErrorCode MoFEM::DMoFEMLoopFiniteElements (DM dm, const std::string fe_name, boost::shared_ptr< MoFEM::FEMethod > method)
 Executes FEMethod for finite elements in DM. More...
 
PetscErrorCode MoFEM::DMoFEMLoopDofs (DM dm, const char field_name[], MoFEM::DofMethod *method)
 execute method for dofs on field in problem More...
 
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. More...
 
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 More...
 
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. More...
 
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. More...
 
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 More...
 
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 More...
 
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 More...
 
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 More...
 
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 More...
 
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 More...
 
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 More...
 
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 More...
 
PetscErrorCode MoFEM::DMMoFEMGetKspCtx (DM dm, MoFEM::KspCtx **ksp_ctx)
 get MoFEM::KspCtx data structure More...
 
PetscErrorCode MoFEM::DMMoFEMGetKspCtx (DM dm, const boost::shared_ptr< MoFEM::KspCtx > &ksp_ctx)
 get MoFEM::KspCtx data structure More...
 
PetscErrorCode MoFEM::DMMoFEMSetKspCtx (DM dm, boost::shared_ptr< MoFEM::KspCtx > &ksp_ctx)
 set MoFEM::KspCtx data structure More...
 
PetscErrorCode MoFEM::DMMoFEMGetSnesCtx (DM dm, MoFEM::SnesCtx **snes_ctx)
 get MoFEM::SnesCtx data structure More...
 
PetscErrorCode MoFEM::DMMoFEMGetSnesCtx (DM dm, const boost::shared_ptr< MoFEM::SnesCtx > &snes_ctx)
 get MoFEM::SnesCtx data structure More...
 
PetscErrorCode MoFEM::DMMoFEMSetSnesCtx (DM dm, boost::shared_ptr< MoFEM::SnesCtx > &snes_ctx)
 Set MoFEM::SnesCtx data structure. More...
 
PetscErrorCode MoFEM::DMMoFEMGetTsCtx (DM dm, MoFEM::TsCtx **ts_ctx)
 get MoFEM::TsCtx data structure More...
 
PetscErrorCode MoFEM::DMMoFEMGetTsCtx (DM dm, const boost::shared_ptr< MoFEM::TsCtx > &ts_ctx)
 get MoFEM::TsCtx data structure More...
 
PetscErrorCode MoFEM::DMMoFEMSetTsCtx (DM dm, boost::shared_ptr< MoFEM::TsCtx > &ts_ctx)
 Set MoFEM::TsCtx data structureIt take over pointer, do not delete it, DM will destroy pointer when is destroyed. More...
 
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. More...
 
PetscErrorCode MoFEM::DMCreate_MoFEM (DM dm)
 Create dm data structure with MoFEM data structure. More...
 
PetscErrorCode MoFEM::DMDestroy_MoFEM (DM dm)
 Destroys dm with MoFEM data structure. More...
 
PetscErrorCode MoFEM::DMCreateGlobalVector_MoFEM (DM dm, Vec *g)
 DMShellSetCreateGlobalVectorsets the routine to create a global vector associated with the shell DM. More...
 
PetscErrorCode MoFEM::DMCreateGlobalVector_MoFEM (DM dm, SmartPetscObj< Vec > &g_ptr)
 DMShellSetCreateGlobalVectorsets the routine to create a global vector associated with the shell DM. More...
 
PetscErrorCode MoFEM::DMCreateLocalVector_MoFEM (DM dm, Vec *l)
 DMShellSetCreateLocalVectorsets the routine to create a local vector associated with the shell DM. More...
 
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[], EntityType lo_type=MBVERTEX, EntityType hi_type=MBMAXTYPE)
 
PetscErrorCode MoFEM::DMMoFEMAddSubFieldCol (DM dm, const char field_name[], EntityType lo_type=MBVERTEX, EntityType hi_type=MBMAXTYPE)
 
PetscErrorCode MoFEM::DMMoFEMGetIsSubDM (DM dm, PetscBool *is_sub_dm)
 
PetscErrorCode MoFEM::DMMoFEMAddRowCompositeProblem (DM dm, const char prb_name[])
 Add problem to composite DM on rowThis create block on row with DOFs from problem of given name. More...
 
PetscErrorCode MoFEM::DMMoFEMAddColCompositeProblem (DM dm, const char prb_name[])
 Add problem to composite DM on colThis create block on col with DOFs from problem of given name. More...
 
PetscErrorCode MoFEM::DMMoFEMGetIsCompDM (DM dm, PetscBool *is_comp_dm)
 Get if this DM is composite DM. More...
 
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 More...
 
MoFEMErrorCode DMMGViaApproxOrdersPushBackCoarseningIS (DM, IS is, Mat A, Mat *subA, bool create_sub_matrix, bool shell_sub_a)
 Push back coarsening level to MG via approximation orders. More...
 
MoFEMErrorCode DMMGViaApproxOrdersPopBackCoarseningIS (DM)
 Pop is form MG via approximation orders. More...
 
MoFEMErrorCode DMMGViaApproxOrdersClearCoarseningIS (DM)
 Clear approximation orders. More...
 
MoFEMErrorCode DMRegister_MGViaApproxOrders (const char sname[])
 Register DM for Multi-Grid via approximation orders. More...
 
MoFEMErrorCode DMCreateMatrix_MGViaApproxOrders (DM dm, Mat *M)
 Create matrix for Multi-Grid via approximation orders. More...
 
MoFEMErrorCode DMCoarsen_MGViaApproxOrders (DM dm, MPI_Comm comm, DM *dmc)
 Coarsen DM. More...
 

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

◆ DMCoarsen_MGViaApproxOrders()

MoFEMErrorCode 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 477 of file PCMGSetUpViaApproxOrders.cpp.

477  {
478  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
480  GET_DM_FIELD(dm);
481  CHKERR PetscObjectGetComm((PetscObject)dm, &comm);
482  CHKERR DMCreate(comm, dmc);
483  (*dmc)->data = dm->data;
484  DMType type;
485  CHKERR DMGetType(dm, &type);
486  CHKERR DMSetType(*dmc, type);
487  CHKERR PetscObjectReference((PetscObject)(*dmc));
488  PetscInfo1(dm, "Coarsen DMMGViaApproxOrders leveldown = %d\n", dm->leveldown);
490 }
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:476
#define GET_DM_FIELD(DM)
#define CHKERR
Inline error check.
Definition: definitions.h:595
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:406

◆ DMCreate_MoFEM()

PetscErrorCode MoFEM::DMCreate_MoFEM ( DM  dm)

Create dm data structure with MoFEM data structure.

Definition at line 81 of file DMMMoFEM.cpp.

81  {
82  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
84  dm->data = new DMCtx();
87 }
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:476
PetscErrorCode DMSetOperators_MoFEM(DM dm)
Set operators for MoFEM dm.
Definition: DMMMoFEM.cpp:59
#define CHKERR
Inline error check.
Definition: definitions.h:595
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:406

◆ 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 1246 of file DMMMoFEM.cpp.

1247  {
1248  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1250 
1251  if (numFields) {
1252  *numFields = 0;
1253  }
1254  if (fieldNames) {
1255  *fieldNames = NULL;
1256  }
1257  if (fields) {
1258  *fields = NULL;
1259  }
1260 
1261  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
1262  const Field_multiIndex *fields_ptr;
1263  CHKERR dm_field->mField_ptr->get_fields(&fields_ptr);
1264  Field_multiIndex::iterator fit, hi_fit;
1265  fit = fields_ptr->begin();
1266  hi_fit = fields_ptr->end();
1267  *numFields = std::distance(fit, hi_fit);
1268 
1269  if (fieldNames) {
1270  CHKERR PetscMalloc1(*numFields, fieldNames);
1271  }
1272  if (fields) {
1273  CHKERR PetscMalloc1(*numFields, fields);
1274  }
1275 
1276  for (int f = 0; fit != hi_fit; fit++, f++) {
1277  if (fieldNames) {
1278  CHKERR PetscStrallocpy(fit->get()->getName().c_str(),
1279  (char **)&(*fieldNames)[f]);
1280  }
1281  if (fields) {
1282  CHKERR dm_field->mField_ptr->getInterface<ISManager>()
1283  ->isCreateProblemFieldAndRank(
1284  dm_field->problemPtr->getName(), ROW, fit->get()->getName(), 0,
1285  fit->get()->getNbOfCoeffs(), &(*fields)[f]);
1286  }
1287  }
1288 
1290 }
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:476
multi_index_container< boost::shared_ptr< Field >, indexed_by< hashed_unique< tag< BitFieldId_mi_tag >, const_mem_fun< Field, const BitFieldId &, &Field::getId >, HashBit< BitFieldId >, EqBit< BitFieldId > >, ordered_unique< tag< Meshset_mi_tag >, member< Field, EntityHandle, &Field::meshSet > >, ordered_unique< tag< FieldName_mi_tag >, const_mem_fun< Field, boost::string_ref, &Field::getNameRef > >, ordered_non_unique< tag< BitFieldId_space_mi_tag >, const_mem_fun< Field, FieldSpace, &Field::getSpace > > > > Field_multiIndex
Field_multiIndex for Field.
#define CHKERR
Inline error check.
Definition: definitions.h:595
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:406

◆ DMCreateGlobalVector_MoFEM() [1/2]

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

DMShellSetCreateGlobalVectorsets the routine to create a global vector associated with the shell DM.

Examples
elasticity_mixed_formulation.cpp, and reaction_diffusion_equation.cpp.

Definition at line 987 of file DMMMoFEM.cpp.

987  {
988  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
990  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
991  CHKERR dm_field->mField_ptr->getInterface<VecManager>()->vecCreateGhost(
992  dm_field->problemName, COL, g);
994 }
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:500
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:507
#define CHKERR
Inline error check.
Definition: definitions.h:595

◆ DMCreateGlobalVector_MoFEM() [2/2]

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

DMShellSetCreateGlobalVectorsets the routine to create a global vector associated with the shell DM.

Definition at line 996 of file DMMMoFEM.cpp.

996  {
997  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
999  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
1000  CHKERR dm_field->mField_ptr->getInterface<VecManager>()->vecCreateGhost(
1001  dm_field->problemName, COL, g_ptr);
1003 }
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:500
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:507
#define CHKERR
Inline error check.
Definition: definitions.h:595

◆ DMCreateLocalVector_MoFEM()

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

DMShellSetCreateLocalVectorsets the routine to create a local vector associated with the shell DM.

Definition at line 1005 of file DMMMoFEM.cpp.

1005  {
1006  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1008  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
1009  CHKERR dm_field->mField_ptr->getInterface<VecManager>()->vecCreateSeq(
1010  dm_field->problemName, COL, l);
1012 }
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:476
#define CHKERR
Inline error check.
Definition: definitions.h:595
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:406

◆ DMCreateMatrix_MGViaApproxOrders()

MoFEMErrorCode 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 446 of file PCMGSetUpViaApproxOrders.cpp.

446  {
447 
448  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
450  GET_DM_FIELD(dm);
451 
452  int leveldown = dm->leveldown;
453 
454  if (dm_field->kspOperators.empty()) {
456  } else {
457  MPI_Comm comm;
458  CHKERR PetscObjectGetComm((PetscObject)dm, &comm);
459  if (dm_field->kspOperators.empty()) {
460  SETERRQ(comm, MOFEM_DATA_INCONSISTENCY,
461  "data inconsistency, operator can not be set");
462  }
463  if (static_cast<int>(dm_field->kspOperators.size()) < leveldown) {
464  SETERRQ(comm, MOFEM_DATA_INCONSISTENCY,
465  "data inconsistency, no IS for that level");
466  }
467  *M = dm_field->kspOperators[dm_field->kspOperators.size() - 1 - leveldown];
468  CHKERR PetscObjectReference((PetscObject)*M);
469  }
470 
471  PetscInfo1(dm, "Create Matrix DMMGViaApproxOrders leveldown = %d\n",
472  leveldown);
473 
475 }
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:476
PetscErrorCode DMCreateMatrix_MoFEM(DM dm, Mat *M)
Definition: DMMMoFEM.cpp:1014
#define GET_DM_FIELD(DM)
#define CHKERR
Inline error check.
Definition: definitions.h:595
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:406

◆ 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
elasticity_mixed_formulation.cpp, and reaction_diffusion_equation.cpp.

Definition at line 1014 of file DMMMoFEM.cpp.

1014  {
1015  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1017  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
1018  if (strcmp(dm->mattype, MATMPIAIJ) == 0) {
1019  CHKERR dm_field->mField_ptr->getInterface<MatrixManager>()
1020  ->createMPIAIJWithArrays<PetscGlobalIdx_mi_tag>(dm_field->problemName,
1021  M);
1022  } else if (strcmp(dm->mattype, MATAIJ) == 0) {
1023  CHKERR dm_field->mField_ptr->getInterface<MatrixManager>()
1024  ->createSeqAIJWithArrays<PetscLocalIdx_mi_tag>(dm_field->problemName,
1025  M);
1026  } else {
1027  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
1028  "Matrix type not implemented");
1029  }
1031 }
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:476
#define CHKERR
Inline error check.
Definition: definitions.h:595
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:406

◆ 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 1033 of file DMMMoFEM.cpp.

1033  {
1034  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1036  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
1037  if (strcmp(dm->mattype, MATMPIAIJ) == 0) {
1038  CHKERR dm_field->mField_ptr->getInterface<MatrixManager>()
1039  ->createMPIAIJWithArrays<PetscGlobalIdx_mi_tag>(dm_field->problemName,
1040  M);
1041  } else if (strcmp(dm->mattype, MATAIJ) == 0) {
1042  CHKERR dm_field->mField_ptr->getInterface<MatrixManager>()
1043  ->createSeqAIJWithArrays<PetscLocalIdx_mi_tag>(dm_field->problemName,
1044  M);
1045  } else {
1046  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
1047  "Matrix type not implemented");
1048  }
1050 }
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:476
#define CHKERR
Inline error check.
Definition: definitions.h:595
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:406

◆ DMDestroy_MoFEM()

PetscErrorCode MoFEM::DMDestroy_MoFEM ( DM  dm)

Destroys dm with MoFEM data structure.

destroy the MoFEM structure

Definition at line 89 of file DMMMoFEM.cpp.

89  {
90  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
91  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
93  if (((DMCtx *)dm->data)->referenceNumber == 0) {
94  if (dm_field->destroyProblem) {
95  if (dm_field->mField_ptr->check_problem(dm_field->problemName)) {
96  dm_field->mField_ptr->delete_problem(dm_field->problemName);
97  } // else problem has to be deleted by the user
98  }
99  // cerr << "Destroy " << dm_field->problemName << endl;
100  delete ((DMCtx *)dm->data);
101  } else {
102  // cerr << "Dereference " << dm_field->problemName << " " <<
103  // ((DMCtx*)dm->data)->referenceNumber << endl;
104  (((DMCtx *)dm->data)->referenceNumber)--;
105  }
107 }
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:500
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:507

◆ 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 1174 of file DMMMoFEM.cpp.

1175  {
1176  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1178  CHKERR VecGhostUpdateBegin(g, INSERT_VALUES, SCATTER_FORWARD);
1180 }
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:476
#define CHKERR
Inline error check.
Definition: definitions.h:595
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:406

◆ 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 1182 of file DMMMoFEM.cpp.

1182  {
1184  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1186 
1187  CHKERR VecGhostUpdateEnd(g, INSERT_VALUES, SCATTER_FORWARD);
1188 
1189  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
1190  int nb_dofs = dm_field->problemPtr->getNbLocalDofsRow();
1191  int nb_ghost = dm_field->problemPtr->getNbGhostDofsRow();
1192 
1193  double *array_loc, *array_glob;
1194  CHKERR VecGetArray(l, &array_loc);
1195  CHKERR VecGetArray(g, &array_glob);
1196  switch (mode) {
1197  case INSERT_VALUES:
1198  cblas_dcopy(nb_dofs + nb_ghost, array_glob, 1, array_loc, 1);
1199  break;
1200  case ADD_VALUES:
1201  cblas_daxpy(nb_dofs + nb_ghost, 1, array_glob, 1, array_loc, 1);
1202  break;
1203  default:
1204  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "not implemented");
1205  }
1206  CHKERR VecRestoreArray(l, &array_loc);
1207  CHKERR VecRestoreArray(g, &array_glob);
1209 }
void cblas_daxpy(const int N, const double alpha, const double *X, const int incX, double *Y, const int incY)
Definition: cblas_daxpy.c:11
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:500
void cblas_dcopy(const int N, const double *X, const int incX, double *Y, const int incY)
Definition: cblas_dcopy.c:11
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:476
#define CHKERR
Inline error check.
Definition: definitions.h:595
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:406

◆ 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 1211 of file DMMMoFEM.cpp.

1212  {
1213 
1214  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1216 
1217  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
1218  int nb_dofs = dm_field->problemPtr->getNbLocalDofsRow();
1219  int nb_ghost = dm_field->problemPtr->getNbGhostDofsRow();
1220 
1221  double *array_loc, *array_glob;
1222  CHKERR VecGetArray(l, &array_loc);
1223  CHKERR VecGetArray(g, &array_glob);
1224  switch (mode) {
1225  case INSERT_VALUES:
1226  cblas_dcopy(nb_dofs + nb_ghost, array_loc, 1, array_glob, 1);
1227  break;
1228  case ADD_VALUES:
1229  cblas_daxpy(nb_dofs + nb_ghost, 1, array_loc, 1, array_glob, 1);
1230  break;
1231  default:
1232  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "not implemented");
1233  }
1234  CHKERR VecRestoreArray(l, &array_loc);
1235  CHKERR VecRestoreArray(g, &array_glob);
1236 
1238 }
void cblas_daxpy(const int N, const double alpha, const double *X, const int incX, double *Y, const int incY)
Definition: cblas_daxpy.c:11
void cblas_dcopy(const int N, const double *X, const int incX, double *Y, const int incY)
Definition: cblas_dcopy.c:11
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:476
#define CHKERR
Inline error check.
Definition: definitions.h:595
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:406

◆ 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 1240 of file DMMMoFEM.cpp.

1240  {
1241  //
1244 }
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:500
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:507

◆ DMMGViaApproxOrdersClearCoarseningIS()

MoFEMErrorCode DMMGViaApproxOrdersClearCoarseningIS ( DM  )

Clear approximation orders.

Parameters
DMdm
Returns
Error code

Definition at line 296 of file PCMGSetUpViaApproxOrders.cpp.

296  {
297  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
299  GET_DM_FIELD(dm);
300  CHKERR dm_field->destroyCoarseningIS();
301  PetscInfo(dm, "Clear DMs data structures\n");
303 }
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:476
#define GET_DM_FIELD(DM)
#define CHKERR
Inline error check.
Definition: definitions.h:595
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:406

◆ DMMGViaApproxOrdersPopBackCoarseningIS()

MoFEMErrorCode DMMGViaApproxOrdersPopBackCoarseningIS ( DM  )

Pop is form MG via approximation orders.

Parameters
DMdm
ispop back IS
Returns
error code

Definition at line 280 of file PCMGSetUpViaApproxOrders.cpp.

280  {
281  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
283  GET_DM_FIELD(dm);
284  if (dm_field->coarseningIS.back()) {
285  CHKERR ISDestroy(&dm_field->coarseningIS.back());
286  dm_field->coarseningIS.pop_back();
287  }
288  if (dm_field->kspOperators.back()) {
289  CHKERR MatDestroy(&dm_field->kspOperators.back());
290  }
291  dm_field->kspOperators.pop_back();
292  PetscInfo(dm, "Pop back IS to DMMGViaApproxOrders\n");
294 }
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:476
#define GET_DM_FIELD(DM)
#define CHKERR
Inline error check.
Definition: definitions.h:595
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:406

◆ DMMGViaApproxOrdersPushBackCoarseningIS()

MoFEMErrorCode DMMGViaApproxOrdersPushBackCoarseningIS ( DM  ,
IS  is,
Mat  A,
Mat *  subA,
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 230 of file PCMGSetUpViaApproxOrders.cpp.

233  {
234  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
236  GET_DM_FIELD(dm);
237  dm_field->coarseningIS.push_back(is);
238  dm_field->shellMatrixCtxPtr.push_back(new PCMGSubMatrixCtx_private(A, is));
239  if (is) {
240  CHKERR PetscObjectReference((PetscObject)is);
241  }
242  if (is) {
243  IS is2 = is;
244  if (dm_field->aO) {
245  CHKERR ISDuplicate(is, &is2);
246  CHKERR ISCopy(is, is2);
247  CHKERR AOApplicationToPetscIS(dm_field->aO, is2);
248  }
249  if (create_sub_matrix) {
250  if (shell_sub_a) {
251  int n, N;
252  CHKERR ISGetSize(is, &N);
253  CHKERR ISGetLocalSize(is, &n);
254  MPI_Comm comm;
255  CHKERR PetscObjectGetComm((PetscObject)A, &comm);
256  CHKERR MatCreateShell(comm, n, n, N, N,
257  &(dm_field->shellMatrixCtxPtr.back()), subA);
258  CHKERR MatShellSetOperation(*subA, MATOP_MULT,
259  (void (*)(void))sub_mat_mult);
260  CHKERR MatShellSetOperation(*subA, MATOP_MULT_ADD,
261  (void (*)(void))sub_mat_mult_add);
262  CHKERR MatShellSetOperation(*subA, MATOP_SOR,
263  (void (*)(void))sub_mat_sor);
264  } else {
265  CHKERR MatGetSubMatrix(A, is2, is2, MAT_INITIAL_MATRIX, subA);
266  }
267  }
268  if (dm_field->aO) {
269  CHKERR ISDestroy(&is2);
270  }
271  dm_field->kspOperators.push_back(*subA);
272  CHKERR PetscObjectReference((PetscObject)(*subA));
273  } else {
274  SETERRQ(PETSC_COMM_WORLD, MOFEM_DATA_INCONSISTENCY, "data inconsistency");
275  }
276  PetscInfo(dm, "Push back IS to DMMGViaApproxOrders\n");
278 }
MoFEMErrorCode sub_mat_mult_add(Mat a, Vec x, Vec f)
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:476
MoFEMErrorCode sub_mat_sor(Mat mat, Vec b, PetscReal omega, MatSORType flag, PetscReal shift, PetscInt its, PetscInt lits, Vec x)
#define GET_DM_FIELD(DM)
#define CHKERR
Inline error check.
Definition: definitions.h:595
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:406
MoFEMErrorCode sub_mat_mult(Mat a, Vec x, Vec f)
const int N
Definition: speed_test.cpp:3

◆ DMMoFEMAddColCompositeProblem()

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

Add problem to composite DM on colThis 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 296 of file DMMMoFEM.cpp.

296  {
297  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
299  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
300  if (!dm->data) {
301  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
302  "data structure for MoFEM not yet created");
303  }
304  if (!dm_field->isCompDM) {
305  dm_field->isCompDM = PETSC_TRUE;
306  }
307  if (dm_field->isSquareMatrix) {
308  SETERRQ(PETSC_COMM_SELF, MOFEM_INVALID_DATA,
309  "No need to add problem on column when problem block structurally "
310  "symmetric");
311  }
312  dm_field->colCompPrb.push_back(prb_name);
314 }
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:476
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:406

◆ DMMoFEMAddElement()

PetscErrorCode MoFEM::DMMoFEMAddElement ( DM  dm,
const char  fe_name[] 
)

add element to dm

Note
add_file is a collective, should be executed on all processors. Otherwise could lead to deadlock.
Examples
analytical_poisson_field_split.cpp, cell_forces.cpp, dm_build_partitioned_mesh.cpp, dm_create_subdm.cpp, elasticity.cpp, elasticity_mixed_formulation.cpp, EshelbianPlasticity.cpp, lorentz_force.cpp, MagneticElement.hpp, minimal_surface_area.cpp, Remodeling.cpp, simple_elasticity.cpp, and UnsaturatedFlow.hpp.

Definition at line 410 of file DMMMoFEM.cpp.

410  {
411  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
413  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
414  ierr = dm_field->mField_ptr->modify_problem_add_finite_element(
415  dm_field->problemName, fe_name);
416  CHKERRG(ierr);
418 }
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:500
#define CHKERRG(n)
Check error code of MoFEM/MOAB/PETSc function.
Definition: definitions.h:543
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:507

◆ DMMoFEMAddRowCompositeProblem()

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

Add problem to composite DM on rowThis 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 278 of file DMMMoFEM.cpp.

278  {
279  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
281  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
282  if (!dm->data) {
283  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
284  "data structure for MoFEM not yet created");
285  }
286  if (!dm_field->isCompDM) {
287  dm_field->isCompDM = PETSC_TRUE;
288  }
289  dm_field->rowCompPrb.push_back(prb_name);
290  if (dm_field->isSquareMatrix) {
291  dm_field->colCompPrb.push_back(prb_name);
292  }
294 }
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:476
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:406

◆ DMMoFEMAddSubFieldCol()

PetscErrorCode MoFEM::DMMoFEMAddSubFieldCol ( DM  dm,
const char  field_name[],
EntityType  lo_type = MBVERTEX,
EntityType  hi_type = MBMAXTYPE 
)

Add field to sub dm problem on columns

Examples
analytical_poisson_field_split.cpp, cell_forces.cpp, and dm_create_subdm.cpp.

Definition at line 210 of file DMMMoFEM.cpp.

211  {
212  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
214  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
215  if (!dm->data) {
216  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
217  "data structure for MoFEM not yet created");
218  }
219  if (!dm_field->isSubDM) {
220  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "this is not sub-dm");
221  }
222  dm_field->colFields.push_back(field_name);
223  if (lo_type != MBVERTEX || hi_type != MBMAXTYPE) {
224  if (!dm_field->mapTypeCol)
225  dm_field->mapTypeCol = boost::make_shared<
226  std::map<std::string, std::pair<EntityType, EntityType>>>();
227  (*dm_field->mapTypeCol)[field_name] =
228  std::pair<EntityType, EntityType>(lo_type, hi_type);
229  }
231 }
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:500
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:507

◆ DMMoFEMAddSubFieldRow()

PetscErrorCode MoFEM::DMMoFEMAddSubFieldRow ( DM  dm,
const char  field_name[],
EntityType  lo_type = MBVERTEX,
EntityType  hi_type = MBMAXTYPE 
)

Add field to sub dm problem on rows

Examples
analytical_poisson_field_split.cpp, cell_forces.cpp, dm_create_subdm.cpp, and EshelbianPlasticity.cpp.

Definition at line 187 of file DMMMoFEM.cpp.

188  {
189  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
191  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
192  if (!dm->data) {
193  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
194  "data structure for MoFEM not yet created");
195  }
196  if (!dm_field->isSubDM) {
197  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "this is not sub-dm");
198  }
199  dm_field->rowFields.push_back(field_name);
200  if (lo_type != MBVERTEX || hi_type != MBMAXTYPE) {
201  if (!dm_field->mapTypeRow)
202  dm_field->mapTypeRow = boost::make_shared<
203  std::map<std::string, std::pair<EntityType, EntityType>>>();
204  (*dm_field->mapTypeRow)[field_name] =
205  std::pair<EntityType, EntityType>(lo_type, hi_type);
206  }
208 }
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:500
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:507

◆ 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.

Examples
cell_forces.cpp, dm_build_partitioned_mesh.cpp, dm_create_subdm.cpp, elasticity.cpp, elasticity_mixed_formulation.cpp, EshelbianPlasticity.cpp, lorentz_force.cpp, MagneticElement.hpp, minimal_surface_area.cpp, Remodeling.cpp, and UnsaturatedFlow.hpp.

Definition at line 109 of file DMMMoFEM.cpp.

112  {
114 
115  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
116  if (!dm->data) {
117  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
118  "data structure for MoFEM not yet created");
119  }
120  if (!m_field_ptr) {
121  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
122  "DM function not implemented into MoFEM");
123  }
124  dm_field->mField_ptr = m_field_ptr;
125  dm_field->problemName = problem_name;
126  if (!m_field_ptr->check_problem(dm_field->problemName)) {
127  // problem is not defined, declare problem internally set bool to
128  // destroyProblem problem with DM
129  dm_field->destroyProblem = PETSC_TRUE;
130  CHKERR dm_field->mField_ptr->add_problem(dm_field->problemName, MF_EXCL,
131  dm_field->verbosity);
132  } else {
133  dm_field->destroyProblem = PETSC_FALSE;
134  }
135  CHKERR dm_field->mField_ptr->modify_problem_ref_level_add_bit(
136  dm_field->problemName, bit_level);
137  CHKERR dm_field->mField_ptr->modify_problem_mask_ref_level_add_bit(
138  dm_field->problemName, bit_mask);
139  dm_field->kspCtx =
140  boost::shared_ptr<KspCtx>(new KspCtx(*m_field_ptr, problem_name));
141  dm_field->snesCtx =
142  boost::shared_ptr<SnesCtx>(new SnesCtx(*m_field_ptr, problem_name));
143  dm_field->tsCtx =
144  boost::shared_ptr<TsCtx>(new TsCtx(*m_field_ptr, problem_name));
145 
146  MPI_Comm comm;
147  CHKERR PetscObjectGetComm((PetscObject)dm, &comm);
148  int result = 0;
149  MPI_Comm_compare(comm, m_field_ptr->get_comm(), &result);
150  // std::cerr << result << " " << MPI_IDENT << " " << MPI_CONGRUENT << " " <<
151  // MPI_SIMILAR << " " << MPI_UNEQUAL << std::endl;
152  if (result > MPI_CONGRUENT) {
153  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
154  "MoFEM and DM using different communicators");
155  }
156  MPI_Comm_size(comm, &dm_field->sIze);
157  MPI_Comm_rank(comm, &dm_field->rAnk);
158 
159  // problem structure
160  CHKERR dm_field->mField_ptr->get_problem(dm_field->problemName,
161  &dm_field->problemPtr);
162 
164 }
virtual bool check_problem(const std::string name)=0
check if problem exist
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:476
#define CHKERR
Inline error check.
Definition: definitions.h:595
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:406
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
analytical_poisson_field_split.cpp, cell_forces.cpp, dm_create_subdm.cpp, and EshelbianPlasticity.cpp.

Definition at line 166 of file DMMMoFEM.cpp.

166  {
168 
169  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
170  if (!dm->data) {
171  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
172  "data structure for MoFEM not yet created");
173  }
174  CHKERR DMMoFEMCreateMoFEM(subdm, dm_field->mField_ptr, problem_name,
175  dm_field->problemPtr->getBitRefLevel());
176 
177  DMCtx *subdm_field = (DMCtx *)subdm->data;
178  subdm_field->isSubDM = PETSC_TRUE;
179  subdm_field->problemMainOfSubPtr = dm_field->problemPtr;
180  subdm_field->isPartitioned = dm_field->isPartitioned;
181  subdm_field->isSquareMatrix = PETSC_FALSE;
182  subdm->ops->setup = DMSubDMSetUp_MoFEM;
183 
185 }
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:476
PetscErrorCode DMSubDMSetUp_MoFEM(DM subdm)
Definition: DMMMoFEM.cpp:1127
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: DMMMoFEM.cpp:109
#define CHKERR
Inline error check.
Definition: definitions.h:595
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:406

◆ 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);

Definition at line 1292 of file DMMMoFEM.cpp.

1293  {
1294  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1296  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
1297  CHKERR dm_field->mField_ptr->getInterface<ISManager>()
1298  ->isCreateProblemFieldAndRank(dm_field->problemPtr->getName(), ROW,
1299  field_name, 0, 1000, is);
1301 }
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:476
#define CHKERR
Inline error check.
Definition: definitions.h:595
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:406

◆ 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 316 of file DMMMoFEM.cpp.

316  {
318  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
320  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
321  *is_comp_dm = dm_field->isCompDM;
323 }
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:500
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:507

◆ DMMoFEMGetIsPartitioned()

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

get if read mesh is partitioned

Definition at line 954 of file DMMMoFEM.cpp.

954  {
955  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
957  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
958  *is_partitioned = dm_field->isPartitioned;
960 }
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:500
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:507

◆ 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 233 of file DMMMoFEM.cpp.

233  {
235  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
237  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
238  *is_sub_dm = dm_field->isSubDM;
240 }
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:500
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:507

◆ DMMoFEMGetKspCtx() [1/2]

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

get MoFEM::KspCtx data structure

Definition at line 888 of file DMMMoFEM.cpp.

888  {
889  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
891  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
892  *ksp_ctx = dm_field->kspCtx.get();
894 }
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:500
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:507

◆ DMMoFEMGetKspCtx() [2/2]

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

get MoFEM::KspCtx data structure

Definition at line 897 of file DMMMoFEM.cpp.

897  {
898  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
900  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
901  const_cast<boost::shared_ptr<MoFEM::KspCtx> &>(ksp_ctx) = dm_field->kspCtx;
903 }
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:500
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:507

◆ DMMoFEMGetProblemFiniteElementLayout()

PetscErrorCode MoFEM::DMMoFEMGetProblemFiniteElementLayout ( DM  dm,
const char  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 386 of file DMMMoFEM.cpp.

387  {
389  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
391  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
392 
393  MPI_Comm comm;
394  CHKERR PetscObjectGetComm((PetscObject)dm, &comm);
395  CHKERR dm_field->problemPtr->getNumberOfElementsByNameAndPart(comm, fe_name,
396  layout);
398 }
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:500
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:476
#define CHKERR
Inline error check.
Definition: definitions.h:595
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:406

◆ DMMoFEMGetProblemPtr()

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

Get pointer to problem data structure.

Examples
EshelbianPlasticity.cpp, field_evaluator.cpp, loop_entities.cpp, lorentz_force.cpp, MagneticElement.hpp, and Remodeling.cpp.

Definition at line 337 of file DMMMoFEM.cpp.

337  {
338  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
340  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
341  if (!dm->data) {
342  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
343  "data structure for MoFEM not yet created");
344  }
345  *problem_ptr = dm_field->problemPtr;
347 }
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:500
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:507

◆ DMMoFEMGetSnesCtx() [1/2]

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

get MoFEM::SnesCtx data structure

Examples
minimal_surface_area.cpp, testing_jacobian_of_hook_element.cpp, and testing_jacobian_of_hook_scaled_with_density_element.cpp.

Definition at line 914 of file DMMMoFEM.cpp.

914  {
915  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
917  DMCtx *dm_field = (DMCtx *)dm->data;
918  *snes_ctx = dm_field->snesCtx.get();
920 }
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:500
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:507

◆ DMMoFEMGetSnesCtx() [2/2]

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

get MoFEM::SnesCtx data structure

Definition at line 923 of file DMMMoFEM.cpp.

923  {
924  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
926  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
927  const_cast<boost::shared_ptr<MoFEM::SnesCtx> &>(snes_ctx) = dm_field->snesCtx;
929 }
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:500
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:507

◆ DMMoFEMGetSquareProblem()

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

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

Definition at line 400 of file DMMMoFEM.cpp.

400  {
403  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
405  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
406  *square_problem = dm_field->isSquareMatrix;
408 }
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:500
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:507

◆ DMMoFEMGetTsCtx() [1/2]

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

get MoFEM::TsCtx data structure

Examples
EshelbianPlasticity.cpp, Remodeling.cpp, and UnsaturatedFlow.hpp.

Definition at line 962 of file DMMMoFEM.cpp.

962  {
963  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
965  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
966  *ts_ctx = dm_field->tsCtx.get();
968 }
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:500
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:507

◆ DMMoFEMGetTsCtx() [2/2]

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

get MoFEM::TsCtx data structure

Definition at line 970 of file DMMMoFEM.cpp.

971  {
972  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
974  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
975  const_cast<boost::shared_ptr<MoFEM::TsCtx> &>(ts_ctx) = dm_field->tsCtx;
977 }
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:500
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:507

◆ 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
analytical_poisson.cpp, elasticity_mixed_formulation.cpp, and simple_elasticity.cpp.

Definition at line 582 of file DMMMoFEM.cpp.

585  {
586  return DMMoFEMKSPSetComputeOperators<const char *, MoFEM::FEMethod *,
589  dm, fe_name, method, pre_only, post_only);
590 }
structure for User Loop Methods on finite elementsIt can be used to calculate stiffness matrices,...
static PetscErrorCode DMMoFEMKSPSetComputeOperators(DM dm, S fe_name, T0 method, T1 pre_only, T2 post_only)
Definition: DMMMoFEM.cpp:563
Data structure to exchange data between mofem and User Loop Methods.It allows to exchange data betwee...

◆ 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 593 of file DMMMoFEM.cpp.

596  {
597  return DMMoFEMKSPSetComputeOperators<const std::string,
598  boost::shared_ptr<MoFEM::FEMethod>>(
599  dm, fe_name, method, pre_only, post_only);
600 }
static PetscErrorCode DMMoFEMKSPSetComputeOperators(DM dm, S fe_name, T0 method, T1 pre_only, T2 post_only)
Definition: DMMMoFEM.cpp:563

◆ 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 codeset KSP right hand side evaluation function
Examples
analytical_poisson.cpp, and simple_elasticity.cpp.

Definition at line 541 of file DMMMoFEM.cpp.

544  {
545  return DMMoFEMKSPSetComputeRHS<const char *, MoFEM::FEMethod *,
547  dm, fe_name, method, pre_only, post_only);
548 }
structure for User Loop Methods on finite elementsIt can be used to calculate stiffness matrices,...
static PetscErrorCode DMMoFEMKSPSetComputeRHS(DM dm, S fe_name, T0 method, T1 pre_only, T2 post_only)
Definition: DMMMoFEM.cpp:522
Data structure to exchange data between mofem and User Loop Methods.It allows to exchange data betwee...

◆ 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 551 of file DMMMoFEM.cpp.

554  {
555  return DMMoFEMKSPSetComputeRHS<const std::string,
556  boost::shared_ptr<MoFEM::FEMethod>,
557  boost::shared_ptr<MoFEM::BasicMethod>,
558  boost::shared_ptr<MoFEM::BasicMethod>>(
559  dm, fe_name, method, pre_only, post_only);
560 }
static PetscErrorCode DMMoFEMKSPSetComputeRHS(DM dm, S fe_name, T0 method, T1 pre_only, T2 post_only)
Definition: DMMMoFEM.cpp:522

◆ DMMoFEMResolveSharedEntities()

PetscErrorCode MoFEM::DMMoFEMResolveSharedEntities ( DM  dm,
const char  fe_name[] 
)

Resolve shared entities.

Parameters
dmdm
fe_namefinite element for which shared entities are resolved
Returns
error code

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

ierr = DMMoFEMGetSquareProblem(dm,"SHELL_ELEMENT"); CHKERRG(ierr);
Tag th;
rval = mField.get_moab().tag_get_handle("ADAPT_ORDER",th); CHKERRQ_MOAB(rval);
ParallelComm* pcomm =
ParallelComm::get_pcomm(&mField.get_moab(),MYPCOMM_INDEX);
// rval = pcomm->reduce_tags(th,MPI_SUM,prisms);
rval = pcomm->exchange_tags(th,prisms);

Definition at line 376 of file DMMMoFEM.cpp.

376  {
378  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
380  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
381  CHKERR dm_field->mField_ptr->resolve_shared_ents(dm_field->problemPtr,
382  fe_name);
384 }
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:500
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:476
#define CHKERR
Inline error check.
Definition: definitions.h:595
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:406

◆ DMMoFEMSetIsPartitioned()

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

sets if read mesh is partitioned

get if read mesh is partitioned

Examples
cell_forces.cpp, elasticity.cpp, elasticity_mixed_formulation.cpp, EshelbianPlasticity.cpp, Remodeling.cpp, simple_elasticity.cpp, and UnsaturatedFlow.hpp.

Definition at line 943 of file DMMMoFEM.cpp.

943  {
944  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
946  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
947  dm_field->isPartitioned = is_partitioned;
949 }
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:500
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:507

◆ DMMoFEMSetKspCtx()

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

set MoFEM::KspCtx data structure

Definition at line 905 of file DMMMoFEM.cpp.

906  {
907  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
909  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
910  dm_field->kspCtx = ksp_ctx;
912 }
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:500
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:507

◆ DMMoFEMSetSnesCtx()

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

Set MoFEM::SnesCtx data structure.

Definition at line 931 of file DMMMoFEM.cpp.

932  {
933  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
935  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
936  dm_field->snesCtx = snes_ctx;
938 }
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:500
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:507

◆ DMMoFEMSetSquareProblem()

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

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

Examples
analytical_poisson_field_split.cpp, cell_forces.cpp, dm_build_partitioned_mesh.cpp, dm_create_subdm.cpp, and EshelbianPlasticity.cpp.

Definition at line 367 of file DMMMoFEM.cpp.

367  {
369  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
371  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
372  dm_field->isSquareMatrix = square_problem;
374 }
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:500
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:507

◆ DMMoFEMSetTsCtx()

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

Set MoFEM::TsCtx data structureIt take over pointer, do not delete it, DM will destroy pointer when is destroyed.

Definition at line 979 of file DMMMoFEM.cpp.

979  {
980  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
982  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
983  dm_field->tsCtx = ts_ctx;
985 }
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:500
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:507

◆ 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
analytical_nonlinear_poisson.cpp, mesh_smoothing.cpp, minimal_surface_area.cpp, testing_jacobian_of_hook_element.cpp, and testing_jacobian_of_hook_scaled_with_density_element.cpp.

Definition at line 622 of file DMMMoFEM.cpp.

625  {
626  return DMMoFEMSNESSetFunction<const char *, MoFEM::FEMethod *,
628  dm, fe_name, method, pre_only, post_only);
629 }
structure for User Loop Methods on finite elementsIt can be used to calculate stiffness matrices,...
static PetscErrorCode DMMoFEMSNESSetFunction(DM dm, S fe_name, T0 method, T1 pre_only, T2 post_only)
Definition: DMMMoFEM.cpp:603
Data structure to exchange data between mofem and User Loop Methods.It allows to exchange data betwee...

◆ 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 632 of file DMMMoFEM.cpp.

635  {
636  return DMMoFEMSNESSetFunction<const std::string,
637  boost::shared_ptr<MoFEM::FEMethod>,
638  boost::shared_ptr<MoFEM::BasicMethod>,
639  boost::shared_ptr<MoFEM::BasicMethod>>(
640  dm, fe_name, method, pre_only, post_only);
641 }
static PetscErrorCode DMMoFEMSNESSetFunction(DM dm, S fe_name, T0 method, T1 pre_only, T2 post_only)
Definition: DMMMoFEM.cpp:603

◆ 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
analytical_nonlinear_poisson.cpp, mesh_smoothing.cpp, minimal_surface_area.cpp, testing_jacobian_of_hook_element.cpp, and testing_jacobian_of_hook_scaled_with_density_element.cpp.

Definition at line 663 of file DMMMoFEM.cpp.

666  {
667  return DMMoFEMSNESSetJacobian<const char *, MoFEM::FEMethod *,
669  dm, fe_name, method, pre_only, post_only);
670 }
structure for User Loop Methods on finite elementsIt can be used to calculate stiffness matrices,...
static PetscErrorCode DMMoFEMSNESSetJacobian(DM dm, S fe_name, T0 method, T1 pre_only, T2 post_only)
Definition: DMMMoFEM.cpp:644
Data structure to exchange data between mofem and User Loop Methods.It allows to exchange data betwee...

◆ 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 673 of file DMMMoFEM.cpp.

676  {
677  return DMMoFEMSNESSetJacobian<const std::string,
678  boost::shared_ptr<MoFEM::FEMethod>,
679  boost::shared_ptr<MoFEM::BasicMethod>,
680  boost::shared_ptr<MoFEM::BasicMethod>>(
681  dm, fe_name, method, pre_only, post_only);
682 }
static PetscErrorCode DMMoFEMSNESSetJacobian(DM dm, S fe_name, T0 method, T1 pre_only, T2 post_only)
Definition: DMMMoFEM.cpp:644

◆ 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

Examples
bone_adaptation.cpp, EshelbianPlasticity.cpp, reaction_diffusion_equation.cpp, and UnsaturatedFlow.hpp.

Definition at line 704 of file DMMMoFEM.cpp.

707  {
708  return DMMoFEMTSSetIFunction<const char *, MoFEM::FEMethod *,
710  dm, fe_name, method, pre_only, post_only);
712 }
structure for User Loop Methods on finite elementsIt can be used to calculate stiffness matrices,...
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:507
static PetscErrorCode DMMoFEMTSSetIFunction(DM dm, S fe_name, T0 method, T1 pre_only, T2 post_only)
Definition: DMMMoFEM.cpp:685
Data structure to exchange data between mofem and User Loop Methods.It allows to exchange data betwee...

◆ 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 715 of file DMMMoFEM.cpp.

718  {
719  return DMMoFEMTSSetIFunction<const std::string,
720  boost::shared_ptr<MoFEM::FEMethod>,
721  boost::shared_ptr<MoFEM::BasicMethod>,
722  boost::shared_ptr<MoFEM::BasicMethod>>(
723  dm, fe_name, method, pre_only, post_only);
725 }
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:507
static PetscErrorCode DMMoFEMTSSetIFunction(DM dm, S fe_name, T0 method, T1 pre_only, T2 post_only)
Definition: DMMMoFEM.cpp:685

◆ DMMoFEMTSSetIJacobian() [1/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

Examples
bone_adaptation.cpp, EshelbianPlasticity.cpp, reaction_diffusion_equation.cpp, and UnsaturatedFlow.hpp.

Definition at line 757 of file DMMMoFEM.cpp.

760  {
761  return DMMoFEMTSSetIJacobian<const std::string,
762  boost::shared_ptr<MoFEM::FEMethod>,
763  boost::shared_ptr<MoFEM::BasicMethod>,
764  boost::shared_ptr<MoFEM::BasicMethod>>(
765  dm, fe_name, method, pre_only, post_only);
766 }
static PetscErrorCode DMMoFEMTSSetIJacobian(DM dm, S fe_name, T0 method, T1 pre_only, T2 post_only)
Definition: DMMMoFEM.cpp:728

◆ DMMoFEMTSSetIJacobian() [2/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 747 of file DMMMoFEM.cpp.

750  {
751  return DMMoFEMTSSetIJacobian<const char *, FEMethod *, MoFEM::BasicMethod *,
752  MoFEM::BasicMethod *>(dm, fe_name, method,
753  pre_only, post_only);
754 }
static PetscErrorCode DMMoFEMTSSetIJacobian(DM dm, S fe_name, T0 method, T1 pre_only, T2 post_only)
Definition: DMMMoFEM.cpp:728
Data structure to exchange data between mofem and User Loop Methods.It allows to exchange data betwee...

◆ DMMoFEMUnSetElement()

PetscErrorCode MoFEM::DMMoFEMUnSetElement ( DM  dm,
const char  fe_name[] 
)

unset element from dm

Definition at line 420 of file DMMMoFEM.cpp.

420  {
421  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
423  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
424  ierr = dm_field->mField_ptr->modify_problem_unset_finite_element(
425  dm_field->problemName, fe_name);
426  CHKERRG(ierr);
428 }
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:500
#define CHKERRG(n)
Check error code of MoFEM/MOAB/PETSc function.
Definition: definitions.h:543
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:507

◆ 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
HookeElement.cpp, and mesh_smoothing.cpp.

Definition at line 325 of file DMMMoFEM.cpp.

325  {
326  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
328  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
329  if (!dm->data) {
330  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
331  "data structure for MoFEM not yet created");
332  }
333  *m_field_ptr = dm_field->mField_ptr;
335 }
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:500
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:507

◆ DMoFEMLoopDofs()

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

execute method for dofs on field in problem

Definition at line 509 of file DMMMoFEM.cpp.

510  {
511  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
513  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
514  ierr =
515  dm_field->mField_ptr->loop_dofs(dm_field->problemPtr, field_name, COL,
516  *method, dm_field->rAnk, dm_field->rAnk);
517  CHKERRG(ierr);
519 }
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:500
#define CHKERRG(n)
Check error code of MoFEM/MOAB/PETSc function.
Definition: definitions.h:543
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:507

◆ DMoFEMLoopFiniteElements() [1/2]

PetscErrorCode MoFEM::DMoFEMLoopFiniteElements ( DM  dm,
const char  fe_name[],
MoFEM::FEMethod method 
)

Executes FEMethod for finite elements in DM.

Parameters
dmMoFEM discrete manager
fe_namename of element
methodpointer to MOFEM::FEMethod
Returns
Error code
Examples
analytical_nonlinear_poisson.cpp, analytical_poisson.cpp, analytical_poisson_field_split.cpp, cell_forces.cpp, elasticity.cpp, elasticity_mixed_formulation.cpp, EshelbianPlasticity.cpp, hello_world.cpp, HookeElement.cpp, MagneticElement.hpp, mesh_smoothing.cpp, minimal_surface_area.cpp, reaction_diffusion_equation.cpp, Remodeling.cpp, simple_elasticity.cpp, simple_interface.cpp, and UnsaturatedFlow.hpp.

Definition at line 492 of file DMMMoFEM.cpp.

493  {
494  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
496  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
497  ierr = DMoFEMLoopFiniteElementsUpAndLowRank(dm, fe_name, method,
498  dm_field->rAnk, dm_field->rAnk);
499  CHKERRG(ierr);
501 }
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:500
#define CHKERRG(n)
Check error code of MoFEM/MOAB/PETSc function.
Definition: definitions.h:543
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:507
PetscErrorCode DMoFEMLoopFiniteElementsUpAndLowRank(DM dm, const char fe_name[], MoFEM::FEMethod *method, int low_rank, int up_rank)
Executes FEMethod for finite elements in DM.
Definition: DMMMoFEM.cpp:473

◆ DMoFEMLoopFiniteElements() [2/2]

PetscErrorCode MoFEM::DMoFEMLoopFiniteElements ( DM  dm,
const std::string  fe_name,
boost::shared_ptr< MoFEM::FEMethod method 
)

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 504 of file DMMMoFEM.cpp.

505  {
506  return DMoFEMLoopFiniteElements(dm, fe_name.c_str(), method.get());
507 }
PetscErrorCode DMoFEMLoopFiniteElements(DM dm, const char fe_name[], MoFEM::FEMethod *method)
Executes FEMethod for finite elements in DM.
Definition: DMMMoFEM.cpp:492

◆ DMoFEMLoopFiniteElementsUpAndLowRank() [1/2]

PetscErrorCode MoFEM::DMoFEMLoopFiniteElementsUpAndLowRank ( DM  dm,
const char  fe_name[],
MoFEM::FEMethod method,
int  low_rank,
int  up_rank 
)

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 473 of file DMMMoFEM.cpp.

475  {
477  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
478  ierr = dm_field->mField_ptr->loop_finite_elements(
479  dm_field->problemPtr, fe_name, *method, low_rank, up_rank);
480  CHKERRG(ierr);
482 }
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:500
#define CHKERRG(n)
Check error code of MoFEM/MOAB/PETSc function.
Definition: definitions.h:543
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:507

◆ 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 
)

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 485 of file DMMMoFEM.cpp.

487  {
488  return DMoFEMLoopFiniteElementsUpAndLowRank(dm, fe_name.c_str(), method.get(),
489  low_rank, up_rank);
490 }
PetscErrorCode DMoFEMLoopFiniteElementsUpAndLowRank(DM dm, const char fe_name[], MoFEM::FEMethod *method, int low_rank, int up_rank)
Executes FEMethod for finite elements in DM.
Definition: DMMMoFEM.cpp:473

◆ 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
analytical_nonlinear_poisson.cpp, analytical_poisson.cpp, analytical_poisson_field_split.cpp, elasticity_mixed_formulation.cpp, mesh_smoothing.cpp, minimal_surface_area.cpp, and simple_elasticity.cpp.

Definition at line 442 of file DMMMoFEM.cpp.

443  {
444  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
446  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
447  ierr = dm_field->mField_ptr->getInterface<VecManager>()->setGlobalGhostVector(
448  dm_field->problemPtr, COL, g, mode, scatter_mode);
449  CHKERRG(ierr);
451 }
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:500
#define CHKERRG(n)
Check error code of MoFEM/MOAB/PETSc function.
Definition: definitions.h:543
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:507

◆ 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
bone_adaptation.cpp, cell_forces.cpp, elasticity.cpp, elasticity_mixed_formulation.cpp, EshelbianPlasticity.cpp, MagneticElement.hpp, mesh_smoothing.cpp, minimal_surface_area.cpp, reaction_diffusion_equation.cpp, testing_jacobian_of_hook_element.cpp, testing_jacobian_of_hook_scaled_with_density_element.cpp, and UnsaturatedFlow.hpp.

Definition at line 430 of file DMMMoFEM.cpp.

431  {
433  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
435  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
436  ierr = dm_field->mField_ptr->getInterface<VecManager>()->setLocalGhostVector(
437  dm_field->problemPtr, COL, l, mode, scatter_mode);
438  CHKERRG(ierr);
440 }
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:500
#define CHKERRG(n)
Check error code of MoFEM/MOAB/PETSc function.
Definition: definitions.h:543
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:507

◆ DMoFEMPostProcessFiniteElements()

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

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

Examples
cell_forces.cpp, elasticity.cpp, elasticity_mixed_formulation.cpp, and simple_elasticity.cpp.

Definition at line 463 of file DMMMoFEM.cpp.

463  {
464  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
466  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
467  ierr = dm_field->mField_ptr->problem_basic_method_postProcess(
468  dm_field->problemPtr, *method);
469  CHKERRG(ierr);
471 }
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:500
#define CHKERRG(n)
Check error code of MoFEM/MOAB/PETSc function.
Definition: definitions.h:543
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:507

◆ DMoFEMPreProcessFiniteElements()

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

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

Examples
cell_forces.cpp, elasticity.cpp, and elasticity_mixed_formulation.cpp.

Definition at line 453 of file DMMMoFEM.cpp.

453  {
454  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
456  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
457  ierr = dm_field->mField_ptr->problem_basic_method_preProcess(
458  dm_field->problemPtr, *method);
459  CHKERRG(ierr);
461 }
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:500
#define CHKERRG(n)
Check error code of MoFEM/MOAB/PETSc function.
Definition: definitions.h:543
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:507

◆ DMRegister_MGViaApproxOrders()

MoFEMErrorCode DMRegister_MGViaApproxOrders ( const char  sname[])

Register DM for Multi-Grid via approximation orders.

Parameters
snameproblem/dm registered name
Returns
error code

Definition at line 392 of file PCMGSetUpViaApproxOrders.cpp.

392  {
394  CHKERR DMRegister(sname, DMCreate_MGViaApproxOrders);
396 }
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:476
MoFEMErrorCode DMCreate_MGViaApproxOrders(DM dm)
Create DM data structure for Multi-Grid via approximation orders.
#define CHKERR
Inline error check.
Definition: definitions.h:595
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:406

◆ DMRegister_MoFEM()

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

Register MoFEM problem.

Examples
analytical_nonlinear_poisson.cpp, analytical_poisson.cpp, analytical_poisson_field_split.cpp, cell_forces.cpp, dm_build_partitioned_mesh.cpp, dm_create_subdm.cpp, elasticity_mixed_formulation.cpp, ep.cpp, field_blas_set_vertex_dofs.cpp, field_evaluator.cpp, hello_world.cpp, loop_entities.cpp, lorentz_force.cpp, MagneticElement.hpp, mesh_smoothing.cpp, minimal_surface_area.cpp, reaction_diffusion_equation.cpp, Remodeling.cpp, simple_elasticity.cpp, simple_interface.cpp, testing_jacobian_of_hook_element.cpp, testing_jacobian_of_hook_scaled_with_density_element.cpp, and unsaturated_transport.cpp.

Definition at line 53 of file DMMMoFEM.cpp.

53  {
55  CHKERR DMRegister(sname, DMCreate_MoFEM);
57 }
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:476
PetscErrorCode DMCreate_MoFEM(DM dm)
Create dm data structure with MoFEM data structure.
Definition: DMMMoFEM.cpp:81
#define CHKERR
Inline error check.
Definition: definitions.h:595
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:406

◆ DMSetFromOptions_MoFEM()

PetscErrorCode MoFEM::DMSetFromOptions_MoFEM ( DM  dm)

Set options for MoFEM DM

Definition at line 1058 of file DMMMoFEM.cpp.

1058  {
1059 #endif
1060 
1061  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1063  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
1064 #if PETSC_VERSION_GE(3, 5, 3)
1065  ierr = PetscOptionsHead(PetscOptionsObject, "DMMoFEM Options");
1066  CHKERRG(ierr);
1067 #else
1068  ierr = PetscOptionsHead("DMMoFEM Options");
1069  CHKERRG(ierr);
1070 #endif
1071  ierr = PetscOptionsBool("-dm_is_partitioned",
1072  "set if mesh is partitioned (works which native MOAB "
1073  "file format, i.e. h5m",
1074  "DMSetUp", dm_field->isPartitioned,
1075  &dm_field->isPartitioned, NULL);
1076  CHKERRG(ierr);
1078 }
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:500
#define CHKERRG(n)
Check error code of MoFEM/MOAB/PETSc function.
Definition: definitions.h:543
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:507

◆ DMSetOperators_MoFEM()

PetscErrorCode MoFEM::DMSetOperators_MoFEM ( DM  dm)

Set operators for MoFEM dm.

Parameters
dm
Returns
error code

Definition at line 59 of file DMMMoFEM.cpp.

59  {
60  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
62 
63  dm->ops->createglobalvector = DMCreateGlobalVector_MoFEM;
64  dm->ops->createlocalvector = DMCreateLocalVector_MoFEM;
65  dm->ops->creatematrix = DMCreateMatrix_MoFEM;
66  dm->ops->setup = DMSetUp_MoFEM;
67  dm->ops->destroy = DMDestroy_MoFEM;
68  dm->ops->setfromoptions = DMSetFromOptions_MoFEM;
69  dm->ops->globaltolocalbegin = DMGlobalToLocalBegin_MoFEM;
70  dm->ops->globaltolocalend = DMGlobalToLocalEnd_MoFEM;
71  dm->ops->localtoglobalbegin = DMLocalToGlobalBegin_MoFEM;
72  dm->ops->localtoglobalend = DMLocalToGlobalEnd_MoFEM;
73  dm->ops->createfieldis = DMCreateFieldIS_MoFEM;
74 
75  // Default matrix type
76  CHKERR DMSetMatType(dm, MATMPIAIJ);
77 
79 }
PetscErrorCode DMGlobalToLocalBegin_MoFEM(DM dm, Vec, InsertMode, Vec)
Definition: DMMMoFEM.cpp:1174
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:476
PetscErrorCode DMSetUp_MoFEM(DM dm)
Definition: DMMMoFEM.cpp:1080
PetscErrorCode DMLocalToGlobalBegin_MoFEM(DM, Vec, InsertMode, Vec)
Definition: DMMMoFEM.cpp:1211
PetscErrorCode DMCreateMatrix_MoFEM(DM dm, Mat *M)
Definition: DMMMoFEM.cpp:1014
PetscErrorCode DMDestroy_MoFEM(DM dm)
Destroys dm with MoFEM data structure.
Definition: DMMMoFEM.cpp:89
PetscErrorCode DMCreateGlobalVector_MoFEM(DM dm, Vec *g)
DMShellSetCreateGlobalVectorsets the routine to create a global vector associated with the shell DM.
Definition: DMMMoFEM.cpp:987
PetscErrorCode DMGlobalToLocalEnd_MoFEM(DM dm, Vec, InsertMode, Vec)
Definition: DMMMoFEM.cpp:1182
#define CHKERR
Inline error check.
Definition: definitions.h:595
PetscErrorCode DMCreateFieldIS_MoFEM(DM dm, PetscInt *numFields, char ***fieldNames, IS **fields)
Definition: DMMMoFEM.cpp:1246
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:406
PetscErrorCode DMLocalToGlobalEnd_MoFEM(DM, Vec, InsertMode, Vec)
Definition: DMMMoFEM.cpp:1240
PetscErrorCode DMSetFromOptions_MoFEM(DM dm)
Definition: DMMMoFEM.cpp:1058
PetscErrorCode DMCreateLocalVector_MoFEM(DM dm, Vec *l)
DMShellSetCreateLocalVectorsets the routine to create a local vector associated with the shell DM.
Definition: DMMMoFEM.cpp:1005

◆ DMSetUp_MoFEM()

PetscErrorCode MoFEM::DMSetUp_MoFEM ( DM  dm)

sets up the MoFEM structures inside a DM object

Definition at line 1080 of file DMMMoFEM.cpp.

1080  {
1081  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1082  ProblemsManager *prb_mng_ptr;
1084  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
1085  CHKERR dm_field->mField_ptr->getInterface(prb_mng_ptr);
1086 
1087  if (dm_field->isCompDM) {
1088  // It is composite probelm
1089  CHKERR prb_mng_ptr->buildCompsedProblem(
1090  dm_field->problemName, dm_field->rowCompPrb, dm_field->colCompPrb,
1091  dm_field->isSquareMatrix == PETSC_TRUE, dm_field->verbosity);
1092  } else {
1093  if (dm_field->isPartitioned) {
1094  CHKERR prb_mng_ptr->buildProblemOnDistributedMesh(
1095  dm_field->problemName, dm_field->isSquareMatrix == PETSC_TRUE,
1096  dm_field->verbosity);
1097  } else {
1098  CHKERR prb_mng_ptr->buildProblem(dm_field->problemName,
1099  dm_field->isSquareMatrix == PETSC_TRUE,
1100  dm_field->verbosity);
1101  CHKERR prb_mng_ptr->partitionProblem(dm_field->problemName,
1102  dm_field->verbosity);
1103  }
1104  }
1105 
1106  // Partition finite elements
1107  if (dm_field->isPartitioned) {
1108  CHKERR prb_mng_ptr->partitionFiniteElements(
1109  dm_field->problemName, true, 0, dm_field->sIze, dm_field->verbosity);
1110  CHKERR prb_mng_ptr->partitionGhostDofsOnDistributedMesh(
1111  dm_field->problemName, dm_field->verbosity);
1112  } else {
1113  // partition finite elemnets
1114  CHKERR prb_mng_ptr->partitionFiniteElements(dm_field->problemName, false,
1115  -1, -1, dm_field->verbosity);
1116  // Get ghost DOFs
1117  CHKERR prb_mng_ptr->partitionGhostDofs(dm_field->problemName,
1118  dm_field->verbosity);
1119  }
1120 
1121  // Set flag that problem is build and partitioned
1122  dm_field->isProblemBuild = PETSC_TRUE;
1123 
1125 }
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:476
#define CHKERR
Inline error check.
Definition: definitions.h:595
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:406

◆ DMSubDMSetUp_MoFEM()

PetscErrorCode MoFEM::DMSubDMSetUp_MoFEM ( DM  subdm)

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

Definition at line 1127 of file DMMMoFEM.cpp.

1127  {
1128  PetscValidHeaderSpecific(subdm, DM_CLASSID, 1);
1129  ProblemsManager *prb_mng_ptr;
1131 
1132  DMCtx *subdm_field = static_cast<DMCtx *>(subdm->data);
1133 
1134  // build sub dm problem
1135  CHKERR subdm_field->mField_ptr->getInterface(prb_mng_ptr);
1136 
1137  map<std::string, std::pair<EntityType, EntityType>> *entity_map_row = nullptr;
1138  map<std::string, std::pair<EntityType, EntityType>> *entity_map_col = nullptr;
1139 
1140  if (subdm_field->mapTypeRow)
1141  entity_map_row = subdm_field->mapTypeRow.get();
1142  if (subdm_field->mapTypeCol)
1143  entity_map_row = subdm_field->mapTypeCol.get();
1144 
1145  CHKERR prb_mng_ptr->buildSubProblem(
1146  subdm_field->problemName, subdm_field->rowFields, subdm_field->colFields,
1147  subdm_field->problemMainOfSubPtr->getName(),
1148  subdm_field->isSquareMatrix == PETSC_TRUE, entity_map_row, entity_map_col,
1149  subdm_field->verbosity);
1150 
1151  // partition problem
1152  subdm_field->isPartitioned = subdm_field->isPartitioned;
1153  if (subdm_field->isPartitioned) {
1154  CHKERR prb_mng_ptr->partitionFiniteElements(subdm_field->problemName, true,
1155  0, subdm_field->sIze,
1156  subdm_field->verbosity);
1157  // set ghost nodes
1158  CHKERR prb_mng_ptr->partitionGhostDofsOnDistributedMesh(
1159  subdm_field->problemName, subdm_field->verbosity);
1160  } else {
1161  // partition finite elements
1162  CHKERR prb_mng_ptr->partitionFiniteElements(subdm_field->problemName, false,
1163  -1, -1, subdm_field->verbosity);
1164  // set ghost nodes
1165  CHKERR prb_mng_ptr->partitionGhostDofs(subdm_field->problemName,
1166  subdm_field->verbosity);
1167  }
1168 
1169  subdm_field->isProblemBuild = PETSC_TRUE;
1170 
1172 }
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:476
#define CHKERR
Inline error check.
Definition: definitions.h:595
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:406