libMesh
Functions
petscdmlibmesh.h File Reference

Go to the source code of this file.

Functions

PETSC_EXTERN PetscErrorCode DMlibMeshSetSystem (DM, libMesh::NonlinearImplicitSystem &)
 Any functional implementation of the DMlibMesh API must compose the following functions with the DM object. More...
 
PETSC_EXTERN PetscErrorCode DMlibMeshGetSystem (DM, libMesh::NonlinearImplicitSystem *&)
 
EXTERN_C_BEGIN PETSC_EXTERN PetscErrorCode DMCreate_libMesh (DM)
 

Function Documentation

◆ DMCreate_libMesh()

EXTERN_C_BEGIN PETSC_EXTERN PetscErrorCode DMCreate_libMesh ( DM  )

Definition at line 944 of file petscdmlibmeshimpl.C.

References DM_libMesh::blockids, DM_libMesh::blocknames, DM_libMesh::decomposition, DM_libMesh::decomposition_type, DMCreateDomainDecomposition_libMesh(), DMCreateFieldDecomposition_libMesh(), DMCreateGlobalVector_libMesh(), DMCreateMatrix_libMesh(), DMDestroy_libMesh(), DMlibMeshGetSystem_libMesh(), DMlibMeshSetSystem_libMesh(), DMSetUp_libMesh(), DMView_libMesh(), libMesh::LibmeshPetscCallQ(), libMesh::PetscFunctionReturn(), DM_libMesh::varids, and DM_libMesh::varnames.

945 {
946  DM_libMesh * dlm;
947 
948  PetscFunctionBegin;
949  PetscValidHeaderSpecific(dm,DM_CLASSID,1);
950  LibmeshPetscCallQ(PetscNew(&dlm));
951  dm->data = dlm;
952 
953  /* DMlibMesh impl */
954  dlm->varids = new(std::map<std::string, unsigned int>);
955  dlm->blockids = new(std::map<std::string, unsigned int>);
956  dlm->varnames = new(std::map<unsigned int, std::string>);
957  dlm->blocknames = new(std::map<unsigned int, std::string>);
958  dlm->decomposition = LIBMESH_PETSC_NULLPTR;
959  dlm->decomposition_type = DMLIBMESH_NO_DECOMPOSITION;
960 
961  /* DM API */
962  dm->ops->createglobalvector = DMCreateGlobalVector_libMesh;
963  dm->ops->createlocalvector = 0; // DMCreateLocalVector_libMesh;
964  dm->ops->getcoloring = 0; // DMGetColoring_libMesh;
965  dm->ops->creatematrix = DMCreateMatrix_libMesh;
966  dm->ops->createinterpolation= 0; // DMCreateInterpolation_libMesh;
967 
968  dm->ops->refine = 0; // DMRefine_libMesh;
969  dm->ops->coarsen = 0; // DMCoarsen_libMesh;
970 
971  // * dm->ops->getinjection was renamed to dm->ops->createinjection in PETSc 5a84ad338 (5 Jul 2019)
972  // * dm->ops-getaggregates was removed in PETSc 97779f9a (5 Jul 2019)
973  // * Both changes were merged into PETSc master in 94aad3ce (7 Jul 2019).
974 #if PETSC_VERSION_LESS_THAN(3,12,0)
975  dm->ops->getinjection = 0; // DMGetInjection_libMesh;
976  dm->ops->getaggregates = 0; // DMGetAggregates_libMesh;
977 #else
978  dm->ops->createinjection = 0;
979 #endif
980 
981 
982  dm->ops->createfielddecomposition = DMCreateFieldDecomposition_libMesh;
983  dm->ops->createdomaindecomposition = DMCreateDomainDecomposition_libMesh;
984 
985  dm->ops->destroy = DMDestroy_libMesh;
986  dm->ops->view = DMView_libMesh;
987  dm->ops->setfromoptions = 0; // DMSetFromOptions_libMesh;
988  dm->ops->setup = DMSetUp_libMesh;
989 
990  /* DMlibMesh API */
991  LibmeshPetscCallQ(PetscObjectComposeFunction((PetscObject)dm,"DMlibMeshSetSystem_C",DMlibMeshSetSystem_libMesh));
992  LibmeshPetscCallQ(PetscObjectComposeFunction((PetscObject)dm,"DMlibMeshGetSystem_C",DMlibMeshGetSystem_libMesh));
993 
994  PetscFunctionReturn(LIBMESH_PETSC_SUCCESS);
995 }
static PetscErrorCode DMCreateDomainDecomposition_libMesh(DM dm, PetscInt *len, char ***namelist, IS **innerislist, IS **outerislist, DM **dmlist)
static PetscErrorCode DMCreateGlobalVector_libMesh(DM dm, Vec *x)
std::map< std::string, unsigned int > * blockids
PetscErrorCode DMlibMeshSetSystem_libMesh(DM dm, NonlinearImplicitSystem &sys)
std::map< std::string, unsigned int > * varids
static PetscErrorCode DMCreateMatrix_libMesh(DM dm, Mat *A)
static PetscErrorCode DMCreateFieldDecomposition_libMesh(DM dm, PetscInt *len, char ***namelist, IS **islist, DM **dmlist)
std::vector< std::set< unsigned int > > * decomposition
static PetscErrorCode DMSetUp_libMesh(DM dm)
PetscErrorCode DMlibMeshGetSystem_libMesh(DM dm, NonlinearImplicitSystem *&sys)
std::map< unsigned int, std::string > * varnames
std::map< unsigned int, std::string > * blocknames
LibmeshPetscCallQ(DMShellGetContext(dm, &ctx))
PetscFunctionReturn(LIBMESH_PETSC_SUCCESS)
static PetscErrorCode DMView_libMesh(DM dm, PetscViewer viewer)
unsigned int decomposition_type
static PetscErrorCode DMDestroy_libMesh(DM dm)

◆ DMlibMeshGetSystem()

PETSC_EXTERN PetscErrorCode DMlibMeshGetSystem ( DM  ,
libMesh::NonlinearImplicitSystem *&   
)

Definition at line 48 of file petscdmlibmesh.C.

References libMesh::LibmeshPetscCallQ(), and libMesh::PetscFunctionReturn().

Referenced by DMlibMeshFunction(), DMlibMeshJacobian(), and DMVariableBounds_libMesh().

49 {
50  PetscErrorCode (*f)(DM,libMesh::NonlinearImplicitSystem *&) = nullptr;
51 
52  PetscFunctionBegin;
53  PetscValidHeaderSpecific(dm,DM_CLASSID,1);
54  LibmeshPetscCallQ(PetscObjectQueryFunction((PetscObject)dm,"DMlibMeshGetSystem_C",&f));
55  if (!f) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP, "DM has no implementation for DMlibMeshGetSystem");
56  LibmeshPetscCallQ((*f)(dm,sys));
57  PetscFunctionReturn(LIBMESH_PETSC_SUCCESS);
58 }
Manages consistently variables, degrees of freedom, coefficient vectors, matrices and non-linear solv...
LibmeshPetscCallQ(DMShellGetContext(dm, &ctx))
PetscFunctionReturn(LIBMESH_PETSC_SUCCESS)

◆ DMlibMeshSetSystem()

PETSC_EXTERN PetscErrorCode DMlibMeshSetSystem ( DM  ,
libMesh::NonlinearImplicitSystem  
)

Any functional implementation of the DMlibMesh API must compose the following functions with the DM object.

(See PETSc documentation on PetscObjectComposeFunction(), a polymorphism mechanism.) The following functions are called in PetscNonlinear Solver (others can be called by users): DMlibMeshSetSystem(), DMlibMeshGetSystem()

Any implementation needs to register its creation routine, DMCreate_libMesh, with PETSc using DMRegister().

Definition at line 34 of file petscdmlibmesh.C.

References libMesh::LibmeshPetscCallQ(), and libMesh::PetscFunctionReturn().

Referenced by libMesh::PetscNonlinearSolver< Number >::init().

35 {
36  PetscErrorCode (*f)(DM,libMesh::NonlinearImplicitSystem &) = nullptr;
37 
38  PetscFunctionBegin;
39  PetscValidHeaderSpecific(dm,DM_CLASSID,1);
40  LibmeshPetscCallQ(PetscObjectQueryFunction((PetscObject)dm,"DMlibMeshSetSystem_C",&f));
41  if (!f) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP, "DM has no implementation for DMlibMeshSetSystem");
42  LibmeshPetscCallQ((*f)(dm,sys));
43  PetscFunctionReturn(LIBMESH_PETSC_SUCCESS);
44 }
Manages consistently variables, degrees of freedom, coefficient vectors, matrices and non-linear solv...
LibmeshPetscCallQ(DMShellGetContext(dm, &ctx))
PetscFunctionReturn(LIBMESH_PETSC_SUCCESS)