xref: /petsc/include/petscdm.h (revision 38221697d93bc7c9a6149445938a5de22e85bc83)
1e1589f56SBarry Smith /*
2e1589f56SBarry Smith       Objects to manage the interactions between the mesh data structures and the algebraic objects
3e1589f56SBarry Smith */
43c48a1e8SJed Brown #if !defined(__PETSCDM_H)
53c48a1e8SJed Brown #define __PETSCDM_H
62c8e378dSBarry Smith #include <petscmat.h>
71e25c274SJed Brown #include <petscdmtypes.h>
8e1589f56SBarry Smith 
9607a6623SBarry Smith PETSC_EXTERN PetscErrorCode DMInitializePackage(void);
10e1589f56SBarry Smith 
11014dd563SJed Brown PETSC_EXTERN PetscClassId DM_CLASSID;
12e1589f56SBarry Smith 
1376bdecfbSBarry Smith /*J
148f6c3df8SBarry Smith     DMType - String with the name of a PETSc DM
15e1589f56SBarry Smith 
16e1589f56SBarry Smith    Level: beginner
17e1589f56SBarry Smith 
18e1589f56SBarry Smith .seealso: DMSetType(), DM
1976bdecfbSBarry Smith J*/
2019fd82e9SBarry Smith typedef const char* DMType;
21e1589f56SBarry Smith #define DMDA        "da"
22e1589f56SBarry Smith #define DMADDA      "adda"
23e1589f56SBarry Smith #define DMCOMPOSITE "composite"
24e1589f56SBarry Smith #define DMSLICED    "sliced"
25fe1899a2SJed Brown #define DMSHELL     "shell"
26b30b9b2eSMatthew G Knepley #define DMMESH      "mesh"
27ab7f58a0SBarry Smith #define DMPLEX      "plex"
28b30b9b2eSMatthew G Knepley #define DMCARTESIAN "cartesian"
298ac4e037SJed Brown #define DMREDUNDANT "redundant"
303a19ef87SMatthew G Knepley #define DMPATCH     "patch"
311d72bce8STim Tautges #define DMMOAB      "moab"
32e1589f56SBarry Smith 
33140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList DMList;
34014dd563SJed Brown PETSC_EXTERN PetscBool         DMRegisterAllCalled;
35014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMCreate(MPI_Comm,DM*);
36*38221697SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMClone(DM,DM*);
3719fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode DMSetType(DM, DMType);
3819fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode DMGetType(DM, DMType *);
39bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode DMRegister(const char[],PetscErrorCode (*)(DM));
40607a6623SBarry Smith PETSC_EXTERN PetscErrorCode DMRegisterAll(void);
41014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMRegisterDestroy(void);
42e1589f56SBarry Smith 
43014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMView(DM,PetscViewer);
44014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMLoad(DM,PetscViewer);
45014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDestroy(DM*);
46014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMCreateGlobalVector(DM,Vec*);
47014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMCreateLocalVector(DM,Vec*);
48014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMGetLocalVector(DM,Vec *);
49014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMRestoreLocalVector(DM,Vec *);
50014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMGetGlobalVector(DM,Vec *);
51014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMRestoreGlobalVector(DM,Vec *);
52014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMClearGlobalVectors(DM);
53014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMGetNamedGlobalVector(DM,const char*,Vec*);
54014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMRestoreNamedGlobalVector(DM,const char*,Vec*);
552348bcf4SPeter Brune PETSC_EXTERN PetscErrorCode DMGetNamedLocalVector(DM,const char*,Vec*);
562348bcf4SPeter Brune PETSC_EXTERN PetscErrorCode DMRestoreNamedLocalVector(DM,const char*,Vec*);
57014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMGetLocalToGlobalMapping(DM,ISLocalToGlobalMapping*);
58014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMGetLocalToGlobalMappingBlock(DM,ISLocalToGlobalMapping*);
59014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMCreateFieldIS(DM,PetscInt*,char***,IS**);
60014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMGetBlockSize(DM,PetscInt*);
6119fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode DMCreateColoring(DM,ISColoringType,MatType,ISColoring*);
6219fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode DMCreateMatrix(DM,MatType,Mat*);
63014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMSetMatrixPreallocateOnly(DM,PetscBool);
64014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMCreateInterpolation(DM,DM,Mat*,Vec*);
65014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMCreateInjection(DM,DM,VecScatter*);
66aa1993deSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMGetWorkArray(DM,PetscInt,PetscDataType,void*);
67aa1993deSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMRestoreWorkArray(DM,PetscInt,PetscDataType,void*);
68014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMRefine(DM,MPI_Comm,DM*);
69014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMCoarsen(DM,MPI_Comm,DM*);
70014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMRefineHierarchy(DM,PetscInt,DM[]);
71014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMCoarsenHierarchy(DM,PetscInt,DM[]);
72014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMCoarsenHookAdd(DM,PetscErrorCode (*)(DM,DM,void*),PetscErrorCode (*)(DM,Mat,Vec,Mat,DM,void*),void*);
73014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMRefineHookAdd(DM,PetscErrorCode (*)(DM,DM,void*),PetscErrorCode (*)(DM,Mat,DM,void*),void*);
74014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMRestrict(DM,Mat,Vec,Mat,DM);
75014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMInterpolate(DM,Mat,DM);
76014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMSetFromOptions(DM);
77146574abSBarry Smith PETSC_EXTERN PetscErrorCode DMViewFromOptions(DM,const char[], const char[]);
78ca266f36SBarry Smith 
79014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMSetUp(DM);
80014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMCreateInterpolationScale(DM,DM,Mat,Vec*);
81014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMCreateAggregates(DM,DM,Mat*);
82baf369e7SPeter Brune PETSC_EXTERN PetscErrorCode DMGlobalToLocalHookAdd(DM,PetscErrorCode (*)(DM,Vec,InsertMode,Vec,void*),PetscErrorCode (*)(DM,Vec,InsertMode,Vec,void*),void*);
83014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMGlobalToLocalBegin(DM,Vec,InsertMode,Vec);
84014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMGlobalToLocalEnd(DM,Vec,InsertMode,Vec);
85014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMLocalToGlobalBegin(DM,Vec,InsertMode,Vec);
86014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMLocalToGlobalEnd(DM,Vec,InsertMode,Vec);
8719fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode DMConvert(DM,DMType,DM*);
88e1589f56SBarry Smith 
896636e97aSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMGetCoordinateDM(DM,DM*);
906636e97aSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMGetCoordinates(DM,Vec*);
916636e97aSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMSetCoordinates(DM,Vec);
926636e97aSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMGetCoordinatesLocal(DM,Vec*);
936636e97aSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMSetCoordinatesLocal(DM,Vec);
94e87bb0d3SMatthew G Knepley PETSC_EXTERN PetscErrorCode DMLocatePoints(DM,Vec,IS*);
956636e97aSMatthew G Knepley 
965dbd56e3SPeter Brune /* block hook interface */
97be081cd6SPeter Brune PETSC_EXTERN PetscErrorCode DMSubDomainHookAdd(DM,PetscErrorCode (*)(DM,DM,void*),PetscErrorCode (*)(DM,VecScatter,VecScatter,DM,void*),void*);
98be081cd6SPeter Brune PETSC_EXTERN PetscErrorCode DMSubDomainRestrict(DM,VecScatter,VecScatter,DM);
995dbd56e3SPeter Brune 
100014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMSetOptionsPrefix(DM,const char []);
10119fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode DMSetVecType(DM,VecType);
10219fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode DMSetMatType(DM,MatType);
103014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMSetApplicationContext(DM,void*);
104014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMSetApplicationContextDestroy(DM,PetscErrorCode (*)(void**));
105014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMGetApplicationContext(DM,void*);
106014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMSetVariableBounds(DM,PetscErrorCode (*)(DM,Vec,Vec));
107014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMHasVariableBounds(DM,PetscBool *);
108b0ae01b7SPeter Brune PETSC_EXTERN PetscErrorCode DMHasColoring(DM,PetscBool *);
109014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMComputeVariableBounds(DM,Vec,Vec);
11093d92d96SBarry Smith 
11181d26defSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMCreateSubDM(DM, PetscInt, PetscInt[], IS *, DM *);
11216621825SDmitry Karpeev PETSC_EXTERN PetscErrorCode DMCreateFieldDecomposition(DM,PetscInt*,char***,IS**,DM**);
1138d4ac253SDmitry Karpeev PETSC_EXTERN PetscErrorCode DMCreateDomainDecomposition(DM,PetscInt*,char***,IS**,IS**,DM**);
114e30e807fSPeter Brune PETSC_EXTERN PetscErrorCode DMCreateDomainDecompositionScatters(DM,PetscInt,DM*,VecScatter**,VecScatter**,VecScatter**);
115e7c4fc90SDmitry Karpeev 
116014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMGetRefineLevel(DM,PetscInt*);
117014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMGetCoarsenLevel(DM,PetscInt*);
118014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMFinalizePackage(void);
119e1589f56SBarry Smith 
1205f1ad066SMatthew G Knepley PETSC_EXTERN PetscErrorCode VecGetDM(Vec, DM*);
1215f1ad066SMatthew G Knepley PETSC_EXTERN PetscErrorCode VecSetDM(Vec, DM);
122c688c046SMatthew G Knepley PETSC_EXTERN PetscErrorCode MatGetDM(Mat, DM*);
123c688c046SMatthew G Knepley PETSC_EXTERN PetscErrorCode MatSetDM(Mat, DM);
1245f1ad066SMatthew G Knepley 
125e1589f56SBarry Smith typedef struct NLF_DAAD* NLF;
126e1589f56SBarry Smith 
127bc2bf880SBarry Smith #define DM_FILE_CLASSID 1211221
1287da65231SMatthew G Knepley 
1297da65231SMatthew G Knepley /* FEM support */
130014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMPrintCellVector(PetscInt, const char [], PetscInt, const PetscScalar []);
131014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMPrintCellMatrix(PetscInt, const char [], PetscInt, PetscInt, const PetscScalar []);
1327da65231SMatthew G Knepley 
133014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMGetDefaultSection(DM, PetscSection *);
134014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMSetDefaultSection(DM, PetscSection);
135014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMGetDefaultGlobalSection(DM, PetscSection *);
136eaf8d80aSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMSetDefaultGlobalSection(DM, PetscSection);
137014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMGetDefaultSF(DM, PetscSF *);
138014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMSetDefaultSF(DM, PetscSF);
139014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMCreateDefaultSF(DM, PetscSection, PetscSection);
140b21d0597SMatthew G Knepley PETSC_EXTERN PetscErrorCode DMGetPointSF(DM, PetscSF *);
141057b4bcdSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMSetPointSF(DM, PetscSF);
14288ed4aceSMatthew G Knepley 
143af122d2aSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMGetNumFields(DM, PetscInt *);
144af122d2aSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMSetNumFields(DM, PetscInt);
145af122d2aSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMGetField(DM, PetscInt, PetscObject *);
146af122d2aSMatthew G Knepley 
1477da65231SMatthew G Knepley typedef struct {
1487da65231SMatthew G Knepley   PetscInt         numQuadPoints; /* The number of quadrature points on an element */
1497da65231SMatthew G Knepley   const PetscReal *quadPoints;    /* The quadrature point coordinates */
1507da65231SMatthew G Knepley   const PetscReal *quadWeights;   /* The quadrature weights */
1517da65231SMatthew G Knepley   PetscInt         numBasisFuncs; /* The number of finite element basis functions on an element */
1527da65231SMatthew G Knepley   PetscInt         numComponents; /* The number of components for each basis function */
1537da65231SMatthew G Knepley   const PetscReal *basis;         /* The basis functions tabulated at the quadrature points */
1547da65231SMatthew G Knepley   const PetscReal *basisDer;      /* The basis function derivatives tabulated at the quadrature points */
1557da65231SMatthew G Knepley } PetscQuadrature;
15683e5549fSMatthew G Knepley 
15792fc3d87SMatthew G Knepley typedef struct {
15892fc3d87SMatthew G Knepley   PetscQuadrature *quad;
15992fc3d87SMatthew G Knepley   void (**f0Funcs)(const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]); /* The f_0 functions for each field */
16092fc3d87SMatthew G Knepley   void (**f1Funcs)(const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]); /* The f_1 functions for each field */
16192fc3d87SMatthew G Knepley   void (**g0Funcs)(const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]); /* The g_0 functions for each field pair */
16292fc3d87SMatthew G Knepley   void (**g1Funcs)(const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]); /* The g_1 functions for each field pair */
16392fc3d87SMatthew G Knepley   void (**g2Funcs)(const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]); /* The g_2 functions for each field pair */
16492fc3d87SMatthew G Knepley   void (**g3Funcs)(const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]); /* The g_3 functions for each field pair */
165bef01ed6SMatthew G. Knepley   void (**bcFuncs)(const PetscReal[], PetscScalar *); /* The boundary condition function for each field component */
166652b88e8SMatthew G. Knepley   PetscQuadrature *quadBd;
167652b88e8SMatthew G. Knepley   void (**f0BdFuncs)(const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]); /* The f_0 functions for each field */
168652b88e8SMatthew G. Knepley   void (**f1BdFuncs)(const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]); /* The f_1 functions for each field */
16992fc3d87SMatthew G Knepley } PetscFEM;
17092fc3d87SMatthew G Knepley 
17183e5549fSMatthew G Knepley typedef enum {PETSC_UNIT_LENGTH, PETSC_UNIT_MASS, PETSC_UNIT_TIME, PETSC_UNIT_CURRENT, PETSC_UNIT_TEMPERATURE, PETSC_UNIT_AMOUNT, PETSC_UNIT_LUMINOSITY, NUM_PETSC_UNITS} PetscUnit;
172e87bb0d3SMatthew G Knepley 
173e87bb0d3SMatthew G Knepley struct _DMInterpolationInfo {
174e87bb0d3SMatthew G Knepley   MPI_Comm   comm;
175e87bb0d3SMatthew G Knepley   PetscInt   dim;    /*1 The spatial dimension of points */
176e87bb0d3SMatthew G Knepley   PetscInt   nInput; /* The number of input points */
177e87bb0d3SMatthew G Knepley   PetscReal *points; /* The input point coordinates */
178e87bb0d3SMatthew G Knepley   PetscInt  *cells;  /* The cell containing each point */
179e87bb0d3SMatthew G Knepley   PetscInt   n;      /* The number of local points */
180e87bb0d3SMatthew G Knepley   Vec        coords; /* The point coordinates */
181e87bb0d3SMatthew G Knepley   PetscInt   dof;    /* The number of components to interpolate */
182e87bb0d3SMatthew G Knepley };
183e87bb0d3SMatthew G Knepley typedef struct _DMInterpolationInfo *DMInterpolationInfo;
184e87bb0d3SMatthew G Knepley 
185e87bb0d3SMatthew G Knepley PetscErrorCode DMInterpolationCreate(MPI_Comm, DMInterpolationInfo *);
186e87bb0d3SMatthew G Knepley PetscErrorCode DMInterpolationSetDim(DMInterpolationInfo, PetscInt);
187e87bb0d3SMatthew G Knepley PetscErrorCode DMInterpolationGetDim(DMInterpolationInfo, PetscInt *);
188e87bb0d3SMatthew G Knepley PetscErrorCode DMInterpolationSetDof(DMInterpolationInfo, PetscInt);
189e87bb0d3SMatthew G Knepley PetscErrorCode DMInterpolationGetDof(DMInterpolationInfo, PetscInt *);
190e87bb0d3SMatthew G Knepley PetscErrorCode DMInterpolationAddPoints(DMInterpolationInfo, PetscInt, PetscReal[]);
191e87bb0d3SMatthew G Knepley PetscErrorCode DMInterpolationSetUp(DMInterpolationInfo, DM, PetscBool);
192e87bb0d3SMatthew G Knepley PetscErrorCode DMInterpolationGetCoordinates(DMInterpolationInfo, Vec *);
193e87bb0d3SMatthew G Knepley PetscErrorCode DMInterpolationGetVector(DMInterpolationInfo, Vec *);
194e87bb0d3SMatthew G Knepley PetscErrorCode DMInterpolationRestoreVector(DMInterpolationInfo, Vec *);
195e87bb0d3SMatthew G Knepley PetscErrorCode DMInterpolationEvaluate(DMInterpolationInfo, DM, Vec, Vec);
196e87bb0d3SMatthew G Knepley PetscErrorCode DMInterpolationDestroy(DMInterpolationInfo *);
197e1589f56SBarry Smith #endif
198