https://mooseframework.inl.gov
Namespaces | Functions
PetscDMMoose.h File Reference

Go to the source code of this file.

Namespaces

 libMesh
 The following methods are specializations for using the libMesh::Parallel::packed_range_* routines for std::strings.
 

Functions

PetscErrorCode DMMooseRegisterAll ()
 
PetscErrorCode DMCreateMoose (MPI_Comm comm, NonlinearSystemBase &nl, const libMesh::DofMapBase &dof_map, const libMesh::System &system, const std::string &dm_name, DM *dm)
 Create a MOOSE DM. More...
 
PetscErrorCode DMMooseReset (DM)
 
PetscErrorCode DMMooseSetNonlinearSystem (DM, NonlinearSystemBase &)
 
PetscErrorCode DMMooseGetNonlinearSystem (DM, NonlinearSystemBase *&)
 
PetscErrorCode DMMooseSetDofMap (DM, const libMesh::DofMapBase &)
 
PetscErrorCode DMMooseGetBlocks (DM, std::vector< std::string > &)
 
PetscErrorCode DMMooseGetVariables (DM, std::vector< std::string > &)
 
PetscErrorCode DMMooseGetSides (DM, std::set< std::string > &)
 
PetscErrorCode DMMooseGetUnSides (DM, std::set< std::string > &)
 
PetscErrorCode DMMooseGetContacts (DM, std::vector< std::pair< std::string, std::string >> &, std::vector< bool > &)
 
PetscErrorCode DMMooseGetUnContacts (DM, std::vector< std::pair< std::string, std::string >> &, std::vector< bool > &)
 
PetscErrorCode DMMooseSetBlocks (DM, const std::vector< std::string > &)
 
PetscErrorCode DMMooseSetVariables (DM, const std::vector< std::string > &)
 
PetscErrorCode DMMooseSetSides (DM, const std::set< std::string > &)
 
PetscErrorCode DMMooseSetUnSides (DM, const std::set< std::string > &)
 
PetscErrorCode DMMooseSetContacts (DM, const std::vector< std::pair< std::string, std::string >> &, const std::vector< bool > &)
 
PetscErrorCode DMMooseSetUnContacts (DM, const std::vector< std::pair< std::string, std::string >> &, const std::vector< bool > &)
 
PetscErrorCode DMMooseSetSplitNames (DM, const std::vector< std::string > &)
 
PetscErrorCode DMMooseGetSplitNames (DM, const std::vector< std::string > &)
 
PetscErrorCode DMMooseSetSplitVars (DM, const std::string &, const std::set< std::string > &)
 
PetscErrorCode DMMooseGetSplitVars (DM, const std::string &, std::set< std::string > &)
 
PetscErrorCode DMMooseSetSplitBlocks (DM, const std::string &, const std::set< std::string > &)
 
PetscErrorCode DMMooseGetSplitBlocks (DM, const std::string &, std::set< std::string > &)
 
PetscErrorCode DMMooseSetSplitSides (DM, const std::string &, const std::set< std::string > &)
 
PetscErrorCode DMMooseGetSplitSides (DM, const std::string &, std::set< std::string > &)
 
PetscErrorCode SNESUpdateDMMoose (SNES snes, PetscInt iteration)
 

Function Documentation

◆ DMCreateMoose()

PetscErrorCode DMCreateMoose ( MPI_Comm  comm,
NonlinearSystemBase nl,
const libMesh::DofMapBase dof_map,
const libMesh::System system,
const std::string &  dm_name,
DM *  dm 
)

Create a MOOSE DM.

Parameters
commThe communicator that the DM should use
nlThe nonlinear system context that the DM is associated with
dof_mapA reference to the DoFMap, which you can get from the system
dm_nameA name to associate with the DM
dmA pointer to the PETSc DM

Definition at line 2061 of file PetscDMMoose.C.

Referenced by FieldSplitPreconditionerTempl< MoosePreconditioner >::createMooseDM(), and DMCreateFieldDecomposition_Moose().

2067 {
2069  LibmeshPetscCallQ(DMCreate(comm, dm));
2070  LibmeshPetscCallQ(DMSetType(*dm, DMMOOSE));
2074  LibmeshPetscCallQ(DMMooseSetName(*dm, dm_name));
2075  PetscFunctionReturn(PETSC_SUCCESS);
2076 }
PetscErrorCode PetscOptionItems *PetscErrorCode DM dm
PetscFunctionBegin
PetscErrorCode DMMooseSetNonlinearSystem(DM dm, NonlinearSystemBase &nl)
Definition: PetscDMMoose.C:225
LibmeshPetscCallQ(DMMooseValidityCheck(dm))
PetscFunctionReturn(PETSC_SUCCESS)
PetscErrorCode DMMooseSetName(DM dm, const std::string &dm_name)
Definition: PetscDMMoose.C:267
PetscErrorCode DMMooseSetDofMap(DM dm, const DofMapBase &dof_map)
Definition: PetscDMMoose.C:239
PetscErrorCode DMMooseSetSystem(DM dm, const System &system)
Definition: PetscDMMoose.C:253

◆ DMMooseGetBlocks()

PetscErrorCode DMMooseGetBlocks ( DM  ,
std::vector< std::string > &   
)

Definition at line 203 of file PetscDMMoose.C.

204 {
207  DM_Moose * dmm = (DM_Moose *)dm->data;
208  for (const auto & it : *(dmm->_block_ids))
209  block_names.push_back(it.first);
210  PetscFunctionReturn(PETSC_SUCCESS);
211 }
PetscErrorCode PetscOptionItems *PetscErrorCode DM dm
PetscFunctionBegin
LibmeshPetscCallQ(DMMooseValidityCheck(dm))
PetscFunctionReturn(PETSC_SUCCESS)
PetscErrorCode DMMooseValidityCheck(DM dm)
Definition: PetscDMMoose.C:132
std::map< std::string, subdomain_id_type > * _block_ids
Definition: PetscDMMoose.C:64
for(PetscInt i=0;i< nvars;++i)

◆ DMMooseGetContacts()

PetscErrorCode DMMooseGetContacts ( DM  ,
std::vector< std::pair< std::string, std::string >> &  ,
std::vector< bool > &   
)

◆ DMMooseGetNonlinearSystem()

PetscErrorCode DMMooseGetNonlinearSystem ( DM  ,
NonlinearSystemBase *&   
)

Definition at line 446 of file PetscDMMoose.C.

Referenced by DMMooseFunction(), DMMooseJacobian(), and DMVariableBounds_Moose().

447 {
450  DM_Moose * dmm = (DM_Moose *)(dm->data);
451  nl = dmm->_nl;
452  PetscFunctionReturn(PETSC_SUCCESS);
453 }
PetscErrorCode PetscOptionItems *PetscErrorCode DM dm
PetscFunctionBegin
LibmeshPetscCallQ(DMMooseValidityCheck(dm))
PetscFunctionReturn(PETSC_SUCCESS)
NonlinearSystemBase * _nl
Definition: PetscDMMoose.C:55
PetscErrorCode DMMooseValidityCheck(DM dm)
Definition: PetscDMMoose.C:132

◆ DMMooseGetSides()

PetscErrorCode DMMooseGetSides ( DM  ,
std::set< std::string > &   
)

◆ DMMooseGetSplitBlocks()

PetscErrorCode DMMooseGetSplitBlocks ( DM  ,
const std::string &  ,
std::set< std::string > &   
)

◆ DMMooseGetSplitNames()

PetscErrorCode DMMooseGetSplitNames ( DM  ,
const std::vector< std::string > &   
)

◆ DMMooseGetSplitSides()

PetscErrorCode DMMooseGetSplitSides ( DM  ,
const std::string &  ,
std::set< std::string > &   
)

◆ DMMooseGetSplitVars()

PetscErrorCode DMMooseGetSplitVars ( DM  ,
const std::string &  ,
std::set< std::string > &   
)

◆ DMMooseGetUnContacts()

PetscErrorCode DMMooseGetUnContacts ( DM  ,
std::vector< std::pair< std::string, std::string >> &  ,
std::vector< bool > &   
)

◆ DMMooseGetUnSides()

PetscErrorCode DMMooseGetUnSides ( DM  ,
std::set< std::string > &   
)

◆ DMMooseGetVariables()

PetscErrorCode DMMooseGetVariables ( DM  ,
std::vector< std::string > &   
)

Definition at line 214 of file PetscDMMoose.C.

215 {
218  DM_Moose * dmm = (DM_Moose *)(dm->data);
219  for (const auto & it : *(dmm->_var_ids))
220  var_names.push_back(it.first);
221  PetscFunctionReturn(PETSC_SUCCESS);
222 }
PetscErrorCode PetscOptionItems *PetscErrorCode DM dm
PetscFunctionBegin
LibmeshPetscCallQ(DMMooseValidityCheck(dm))
std::map< std::string, unsigned int > * _var_ids
Definition: PetscDMMoose.C:60
PetscFunctionReturn(PETSC_SUCCESS)
PetscErrorCode DMMooseValidityCheck(DM dm)
Definition: PetscDMMoose.C:132

◆ DMMooseRegisterAll()

PetscErrorCode DMMooseRegisterAll ( )

Definition at line 2180 of file PetscDMMoose.C.

Referenced by StaticCondensationFieldSplitPreconditioner::setupDM(), and FieldSplitPreconditioner::setupDM().

2181 {
2182  static PetscBool DMMooseRegisterAllCalled = PETSC_FALSE;
2183 
2185  if (!DMMooseRegisterAllCalled)
2186  {
2187  LibmeshPetscCallQ(DMRegister(DMMOOSE, DMCreate_Moose));
2188  DMMooseRegisterAllCalled = PETSC_TRUE;
2189  }
2190  PetscFunctionReturn(PETSC_SUCCESS);
2191 }
PetscFunctionBegin
LibmeshPetscCallQ(DMMooseValidityCheck(dm))
EXTERN_C_BEGIN PetscErrorCode DMCreate_Moose(DM dm)
PetscFunctionReturn(PETSC_SUCCESS)

◆ DMMooseReset()

PetscErrorCode DMMooseReset ( DM  )

Definition at line 1612 of file PetscDMMoose.C.

Referenced by SNESUpdateDMMoose().

1613 {
1614  PetscBool ismoose;
1615  DM_Moose * dmm = (DM_Moose *)(dm->data);
1616 
1618  LibmeshPetscCallQ(PetscObjectTypeCompare((PetscObject)dm, DMMOOSE, &ismoose));
1619  if (!ismoose)
1620  PetscFunctionReturn(PETSC_SUCCESS);
1621  if (!dmm->_nl)
1622  SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONGSTATE, "No Moose system set for DM_Moose");
1623  LibmeshPetscCallQ(ISDestroy(&dmm->_embedding));
1624  for (auto & it : *(dmm->_splits))
1625  {
1626  DM_Moose::SplitInfo & split = it.second;
1627  LibmeshPetscCallQ(ISDestroy(&split._rembedding));
1628  if (split._dm)
1630  }
1631  dm->setupcalled = PETSC_FALSE;
1632  PetscFunctionReturn(PETSC_SUCCESS);
1633 }
PetscErrorCode PetscOptionItems *PetscErrorCode DM dm
PetscFunctionBegin
PetscErrorCode DMMooseReset(DM dm)
PETSC_ERR_ARG_WRONGSTATE
LibmeshPetscCallQ(DMMooseValidityCheck(dm))
PetscFunctionReturn(PETSC_SUCCESS)
std::map< std::string, SplitInfo > * _splits
Definition: PetscDMMoose.C:99
tbb::split split
NonlinearSystemBase * _nl
Definition: PetscDMMoose.C:55

◆ DMMooseSetBlocks()

PetscErrorCode DMMooseSetBlocks ( DM  ,
const std::vector< std::string > &   
)

◆ DMMooseSetContacts()

PetscErrorCode DMMooseSetContacts ( DM  ,
const std::vector< std::pair< std::string, std::string >> &  ,
const std::vector< bool > &   
)

◆ DMMooseSetDofMap()

PetscErrorCode DMMooseSetDofMap ( DM  ,
const libMesh::DofMapBase  
)

Definition at line 239 of file PetscDMMoose.C.

Referenced by DMCreateMoose().

240 {
243  if (dm->setupcalled)
244  SETERRQ(((PetscObject)dm)->comm,
246  "Cannot reset the degree of freedom map after DM has been set up.");
247  DM_Moose * dmm = (DM_Moose *)(dm->data);
248  dmm->_dof_map = &dof_map;
249  PetscFunctionReturn(PETSC_SUCCESS);
250 }
const DofMapBase * _dof_map
Definition: PetscDMMoose.C:56
PetscErrorCode PetscOptionItems *PetscErrorCode DM dm
PetscFunctionBegin
PETSC_ERR_ARG_WRONGSTATE
LibmeshPetscCallQ(DMMooseValidityCheck(dm))
PetscFunctionReturn(PETSC_SUCCESS)
PetscErrorCode DMMooseValidityCheck(DM dm)
Definition: PetscDMMoose.C:132

◆ DMMooseSetNonlinearSystem()

PetscErrorCode DMMooseSetNonlinearSystem ( DM  ,
NonlinearSystemBase  
)

Definition at line 225 of file PetscDMMoose.C.

Referenced by DMCreateMoose().

226 {
229  if (dm->setupcalled)
230  SETERRQ(((PetscObject)dm)->comm,
232  "Cannot reset the NonlinearSystem after DM has been set up.");
233  DM_Moose * dmm = (DM_Moose *)(dm->data);
234  dmm->_nl = &nl;
235  PetscFunctionReturn(PETSC_SUCCESS);
236 }
PetscErrorCode PetscOptionItems *PetscErrorCode DM dm
PetscFunctionBegin
PETSC_ERR_ARG_WRONGSTATE
LibmeshPetscCallQ(DMMooseValidityCheck(dm))
PetscFunctionReturn(PETSC_SUCCESS)
NonlinearSystemBase * _nl
Definition: PetscDMMoose.C:55
PetscErrorCode DMMooseValidityCheck(DM dm)
Definition: PetscDMMoose.C:132

◆ DMMooseSetSides()

PetscErrorCode DMMooseSetSides ( DM  ,
const std::set< std::string > &   
)

Definition at line 340 of file PetscDMMoose.C.

341 {
342  DM_Moose * dmm = (DM_Moose *)dm->data;
343 
346  if (dm->setupcalled)
347  SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for an already setup DM");
348  if (dmm->_sides)
349  delete dmm->_sides;
350  dmm->_sides = new std::set<std::string>(sides);
351  PetscFunctionReturn(PETSC_SUCCESS);
352 }
PetscErrorCode PetscOptionItems *PetscErrorCode DM dm
PetscFunctionBegin
PETSC_ERR_ARG_WRONGSTATE
LibmeshPetscCallQ(DMMooseValidityCheck(dm))
std::set< std::string > * _sides
Definition: PetscDMMoose.C:67
char ** sides
PetscFunctionReturn(PETSC_SUCCESS)
PetscErrorCode DMMooseValidityCheck(DM dm)
Definition: PetscDMMoose.C:132

◆ DMMooseSetSplitBlocks()

PetscErrorCode DMMooseSetSplitBlocks ( DM  ,
const std::string &  ,
const std::set< std::string > &   
)

◆ DMMooseSetSplitNames()

PetscErrorCode DMMooseSetSplitNames ( DM  ,
const std::vector< std::string > &   
)

Definition at line 456 of file PetscDMMoose.C.

Referenced by if().

457 {
460  DM_Moose * dmm = (DM_Moose *)(dm->data);
461 
462  if (dmm->_splits)
463  {
464  for (auto & it : *(dmm->_splits))
465  {
466  LibmeshPetscCallQ(DMDestroy(&(it.second._dm)));
467  LibmeshPetscCallQ(ISDestroy(&(it.second._rembedding)));
468  }
469  delete dmm->_splits;
470  dmm->_splits = LIBMESH_PETSC_NULLPTR;
471  }
472  if (dmm->_splitlocs)
473  {
474  delete dmm->_splitlocs;
475  dmm->_splitlocs = LIBMESH_PETSC_NULLPTR;
476  }
477  dmm->_splits = new std::map<std::string, DM_Moose::SplitInfo>();
478  dmm->_splitlocs = new std::multimap<std::string, unsigned int>();
479  for (unsigned int i = 0; i < split_names.size(); ++i)
480  {
482  info._dm = LIBMESH_PETSC_NULLPTR;
483  info._rembedding = LIBMESH_PETSC_NULLPTR;
484  std::string name = split_names[i];
485  (*dmm->_splits)[name] = info;
486  dmm->_splitlocs->insert(std::make_pair(name, i));
487  }
488  PetscFunctionReturn(PETSC_SUCCESS);
489 }
std::string name(const ElemQuality q)
MPI_Info info
PetscErrorCode PetscOptionItems *PetscErrorCode DM dm
PetscFunctionBegin
LibmeshPetscCallQ(DMMooseValidityCheck(dm))
PetscFunctionReturn(PETSC_SUCCESS)
std::map< std::string, SplitInfo > * _splits
Definition: PetscDMMoose.C:99
PetscErrorCode DMMooseValidityCheck(DM dm)
Definition: PetscDMMoose.C:132
std::multimap< std::string, unsigned int > * _splitlocs
Definition: PetscDMMoose.C:93

◆ DMMooseSetSplitSides()

PetscErrorCode DMMooseSetSplitSides ( DM  ,
const std::string &  ,
const std::set< std::string > &   
)

◆ DMMooseSetSplitVars()

PetscErrorCode DMMooseSetSplitVars ( DM  ,
const std::string &  ,
const std::set< std::string > &   
)

◆ DMMooseSetUnContacts()

PetscErrorCode DMMooseSetUnContacts ( DM  ,
const std::vector< std::pair< std::string, std::string >> &  ,
const std::vector< bool > &   
)

◆ DMMooseSetUnSides()

PetscErrorCode DMMooseSetUnSides ( DM  ,
const std::set< std::string > &   
)

Definition at line 355 of file PetscDMMoose.C.

356 {
357  DM_Moose * dmm = (DM_Moose *)dm->data;
358 
361  if (dm->setupcalled)
362  SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for an already setup DM");
363  if (dmm->_unsides)
364  delete dmm->_unsides;
365  dmm->_unsides = new std::set<std::string>(unsides);
366  PetscFunctionReturn(PETSC_SUCCESS);
367 }
PetscErrorCode PetscOptionItems *PetscErrorCode DM dm
PetscFunctionBegin
PETSC_ERR_ARG_WRONGSTATE
LibmeshPetscCallQ(DMMooseValidityCheck(dm))
PetscFunctionReturn(PETSC_SUCCESS)
std::set< std::string > * _unsides
Definition: PetscDMMoose.C:70
PetscErrorCode DMMooseValidityCheck(DM dm)
Definition: PetscDMMoose.C:132

◆ DMMooseSetVariables()

PetscErrorCode DMMooseSetVariables ( DM  ,
const std::vector< std::string > &   
)

◆ SNESUpdateDMMoose()

PetscErrorCode SNESUpdateDMMoose ( SNES  snes,
PetscInt  iteration 
)

Definition at line 2141 of file PetscDMMoose.C.

2142 {
2143  /* This is called any time the structure of the problem changes in a way that affects the Jacobian
2144  sparsity pattern.
2145  For example, this may happen when NodeFaceConstraints change Jacobian's sparsity pattern based
2146  on newly-detected Penetration.
2147  In that case certain preconditioners (e.g., PCASM) will not work, unless we tell them that the
2148  sparsity pattern has changed.
2149  For now we are rebuilding the whole KSP, when necessary.
2150  */
2151  DM dm;
2152  KSP ksp;
2153  const char * prefix;
2154  MPI_Comm comm;
2155  PC pc;
2156 
2158  if (iteration)
2159  {
2160  /* TODO: limit this only to situations when displaced (un)contact splits are present, as is
2161  * DisplacedProblem(). */
2162  LibmeshPetscCallQ(SNESGetDM(snes, &dm));
2164  LibmeshPetscCallQ(DMSetUp(dm));
2165  LibmeshPetscCallQ(SNESGetKSP(snes, &ksp));
2166  /* Should we rebuild the whole KSP? */
2167  LibmeshPetscCallQ(PetscObjectGetOptionsPrefix((PetscObject)ksp, &prefix));
2168  LibmeshPetscCallQ(PetscObjectGetComm((PetscObject)ksp, &comm));
2169  LibmeshPetscCallQ(PCCreate(comm, &pc));
2170  LibmeshPetscCallQ(PCSetDM(pc, dm));
2171  LibmeshPetscCallQ(PCSetOptionsPrefix(pc, prefix));
2172  LibmeshPetscCallQ(PCSetFromOptions(pc));
2173  LibmeshPetscCallQ(KSPSetPC(ksp, pc));
2174  LibmeshPetscCallQ(PCDestroy(&pc));
2175  }
2176  PetscFunctionReturn(PETSC_SUCCESS);
2177 }
PetscErrorCode PetscOptionItems *PetscErrorCode DM dm
PetscFunctionBegin
PetscErrorCode DMMooseReset(DM dm)
LibmeshPetscCallQ(DMMooseValidityCheck(dm))
PetscFunctionReturn(PETSC_SUCCESS)