xref: /petsc/include/petscdm.h (revision 9bb5e83dfd8d08918f3df3058fee52935e3d9cc1)
1 /*
2       Objects to manage the interactions between the mesh data structures and the algebraic objects
3 */
4 #if !defined(__PETSCDM_H)
5 #define __PETSCDM_H
6 #include <petscmat.h>
7 #include <petscdmtypes.h>
8 
9 PETSC_EXTERN PetscErrorCode DMInitializePackage(void);
10 
11 PETSC_EXTERN PetscClassId DM_CLASSID;
12 
13 /*J
14     DMType - String with the name of a PETSc DM
15 
16    Level: beginner
17 
18 .seealso: DMSetType(), DM
19 J*/
20 typedef const char* DMType;
21 #define DMDA        "da"
22 #define DMADDA      "adda"
23 #define DMCOMPOSITE "composite"
24 #define DMSLICED    "sliced"
25 #define DMSHELL     "shell"
26 #define DMMESH      "mesh"
27 #define DMPLEX      "plex"
28 #define DMCARTESIAN "cartesian"
29 #define DMREDUNDANT "redundant"
30 #define DMPATCH     "patch"
31 #define DMMOAB      "moab"
32 
33 PETSC_EXTERN PetscFunctionList DMList;
34 PETSC_EXTERN PetscBool         DMRegisterAllCalled;
35 PETSC_EXTERN PetscErrorCode DMCreate(MPI_Comm,DM*);
36 PETSC_EXTERN PetscErrorCode DMSetType(DM, DMType);
37 PETSC_EXTERN PetscErrorCode DMGetType(DM, DMType *);
38 PETSC_EXTERN PetscErrorCode DMRegister(const char[],PetscErrorCode (*)(DM));
39 PETSC_EXTERN PetscErrorCode DMRegisterAll(void);
40 PETSC_EXTERN PetscErrorCode DMRegisterDestroy(void);
41 
42 PETSC_EXTERN PetscErrorCode DMView(DM,PetscViewer);
43 PETSC_EXTERN PetscErrorCode DMLoad(DM,PetscViewer);
44 PETSC_EXTERN PetscErrorCode DMDestroy(DM*);
45 PETSC_EXTERN PetscErrorCode DMCreateGlobalVector(DM,Vec*);
46 PETSC_EXTERN PetscErrorCode DMCreateLocalVector(DM,Vec*);
47 PETSC_EXTERN PetscErrorCode DMGetLocalVector(DM,Vec *);
48 PETSC_EXTERN PetscErrorCode DMRestoreLocalVector(DM,Vec *);
49 PETSC_EXTERN PetscErrorCode DMGetGlobalVector(DM,Vec *);
50 PETSC_EXTERN PetscErrorCode DMRestoreGlobalVector(DM,Vec *);
51 PETSC_EXTERN PetscErrorCode DMClearGlobalVectors(DM);
52 PETSC_EXTERN PetscErrorCode DMGetNamedGlobalVector(DM,const char*,Vec*);
53 PETSC_EXTERN PetscErrorCode DMRestoreNamedGlobalVector(DM,const char*,Vec*);
54 PETSC_EXTERN PetscErrorCode DMGetNamedLocalVector(DM,const char*,Vec*);
55 PETSC_EXTERN PetscErrorCode DMRestoreNamedLocalVector(DM,const char*,Vec*);
56 PETSC_EXTERN PetscErrorCode DMGetLocalToGlobalMapping(DM,ISLocalToGlobalMapping*);
57 PETSC_EXTERN PetscErrorCode DMGetLocalToGlobalMappingBlock(DM,ISLocalToGlobalMapping*);
58 PETSC_EXTERN PetscErrorCode DMCreateFieldIS(DM,PetscInt*,char***,IS**);
59 PETSC_EXTERN PetscErrorCode DMGetBlockSize(DM,PetscInt*);
60 PETSC_EXTERN PetscErrorCode DMCreateColoring(DM,ISColoringType,MatType,ISColoring*);
61 PETSC_EXTERN PetscErrorCode DMCreateMatrix(DM,MatType,Mat*);
62 PETSC_EXTERN PetscErrorCode DMSetMatrixPreallocateOnly(DM,PetscBool);
63 PETSC_EXTERN PetscErrorCode DMCreateInterpolation(DM,DM,Mat*,Vec*);
64 PETSC_EXTERN PetscErrorCode DMCreateInjection(DM,DM,VecScatter*);
65 PETSC_EXTERN PetscErrorCode DMGetWorkArray(DM,PetscInt,PetscDataType,void*);
66 PETSC_EXTERN PetscErrorCode DMRestoreWorkArray(DM,PetscInt,PetscDataType,void*);
67 PETSC_EXTERN PetscErrorCode DMRefine(DM,MPI_Comm,DM*);
68 PETSC_EXTERN PetscErrorCode DMCoarsen(DM,MPI_Comm,DM*);
69 PETSC_EXTERN PetscErrorCode DMRefineHierarchy(DM,PetscInt,DM[]);
70 PETSC_EXTERN PetscErrorCode DMCoarsenHierarchy(DM,PetscInt,DM[]);
71 PETSC_EXTERN PetscErrorCode DMCoarsenHookAdd(DM,PetscErrorCode (*)(DM,DM,void*),PetscErrorCode (*)(DM,Mat,Vec,Mat,DM,void*),void*);
72 PETSC_EXTERN PetscErrorCode DMRefineHookAdd(DM,PetscErrorCode (*)(DM,DM,void*),PetscErrorCode (*)(DM,Mat,DM,void*),void*);
73 PETSC_EXTERN PetscErrorCode DMRestrict(DM,Mat,Vec,Mat,DM);
74 PETSC_EXTERN PetscErrorCode DMInterpolate(DM,Mat,DM);
75 PETSC_EXTERN PetscErrorCode DMSetFromOptions(DM);
76 PETSC_EXTERN PetscErrorCode DMViewFromOptions(DM,const char[], const char[]);
77 
78 PETSC_EXTERN PetscErrorCode DMSetUp(DM);
79 PETSC_EXTERN PetscErrorCode DMCreateInterpolationScale(DM,DM,Mat,Vec*);
80 PETSC_EXTERN PetscErrorCode DMCreateAggregates(DM,DM,Mat*);
81 PETSC_EXTERN PetscErrorCode DMGlobalToLocalHookAdd(DM,PetscErrorCode (*)(DM,Vec,InsertMode,Vec,void*),PetscErrorCode (*)(DM,Vec,InsertMode,Vec,void*),void*);
82 PETSC_EXTERN PetscErrorCode DMGlobalToLocalBegin(DM,Vec,InsertMode,Vec);
83 PETSC_EXTERN PetscErrorCode DMGlobalToLocalEnd(DM,Vec,InsertMode,Vec);
84 PETSC_EXTERN PetscErrorCode DMLocalToGlobalBegin(DM,Vec,InsertMode,Vec);
85 PETSC_EXTERN PetscErrorCode DMLocalToGlobalEnd(DM,Vec,InsertMode,Vec);
86 PETSC_EXTERN PetscErrorCode DMConvert(DM,DMType,DM*);
87 
88 PETSC_EXTERN PetscErrorCode DMGetCoordinateDM(DM,DM*);
89 PETSC_EXTERN PetscErrorCode DMGetCoordinates(DM,Vec*);
90 PETSC_EXTERN PetscErrorCode DMSetCoordinates(DM,Vec);
91 PETSC_EXTERN PetscErrorCode DMGetCoordinatesLocal(DM,Vec*);
92 PETSC_EXTERN PetscErrorCode DMSetCoordinatesLocal(DM,Vec);
93 PETSC_EXTERN PetscErrorCode DMLocatePoints(DM,Vec,IS*);
94 
95 /* block hook interface */
96 PETSC_EXTERN PetscErrorCode DMSubDomainHookAdd(DM,PetscErrorCode (*)(DM,DM,void*),PetscErrorCode (*)(DM,VecScatter,VecScatter,DM,void*),void*);
97 PETSC_EXTERN PetscErrorCode DMSubDomainRestrict(DM,VecScatter,VecScatter,DM);
98 
99 PETSC_EXTERN PetscErrorCode DMSetOptionsPrefix(DM,const char []);
100 PETSC_EXTERN PetscErrorCode DMSetVecType(DM,VecType);
101 PETSC_EXTERN PetscErrorCode DMSetMatType(DM,MatType);
102 PETSC_EXTERN PetscErrorCode DMSetApplicationContext(DM,void*);
103 PETSC_EXTERN PetscErrorCode DMSetApplicationContextDestroy(DM,PetscErrorCode (*)(void**));
104 PETSC_EXTERN PetscErrorCode DMGetApplicationContext(DM,void*);
105 PETSC_EXTERN PetscErrorCode DMSetVariableBounds(DM,PetscErrorCode (*)(DM,Vec,Vec));
106 PETSC_EXTERN PetscErrorCode DMHasVariableBounds(DM,PetscBool *);
107 PETSC_EXTERN PetscErrorCode DMHasColoring(DM,PetscBool *);
108 PETSC_EXTERN PetscErrorCode DMComputeVariableBounds(DM,Vec,Vec);
109 
110 PETSC_EXTERN PetscErrorCode DMCreateSubDM(DM, PetscInt, PetscInt[], IS *, DM *);
111 PETSC_EXTERN PetscErrorCode DMCreateFieldDecomposition(DM,PetscInt*,char***,IS**,DM**);
112 PETSC_EXTERN PetscErrorCode DMCreateDomainDecomposition(DM,PetscInt*,char***,IS**,IS**,DM**);
113 PETSC_EXTERN PetscErrorCode DMCreateDomainDecompositionScatters(DM,PetscInt,DM*,VecScatter**,VecScatter**,VecScatter**);
114 
115 PETSC_EXTERN PetscErrorCode DMGetRefineLevel(DM,PetscInt*);
116 PETSC_EXTERN PetscErrorCode DMGetCoarsenLevel(DM,PetscInt*);
117 PETSC_EXTERN PetscErrorCode DMFinalizePackage(void);
118 
119 PETSC_EXTERN PetscErrorCode VecGetDM(Vec, DM*);
120 PETSC_EXTERN PetscErrorCode VecSetDM(Vec, DM);
121 PETSC_EXTERN PetscErrorCode MatGetDM(Mat, DM*);
122 PETSC_EXTERN PetscErrorCode MatSetDM(Mat, DM);
123 
124 typedef struct NLF_DAAD* NLF;
125 
126 #define DM_FILE_CLASSID 1211221
127 
128 /* FEM support */
129 PETSC_EXTERN PetscErrorCode DMPrintCellVector(PetscInt, const char [], PetscInt, const PetscScalar []);
130 PETSC_EXTERN PetscErrorCode DMPrintCellMatrix(PetscInt, const char [], PetscInt, PetscInt, const PetscScalar []);
131 
132 PETSC_EXTERN PetscErrorCode DMGetDefaultSection(DM, PetscSection *);
133 PETSC_EXTERN PetscErrorCode DMSetDefaultSection(DM, PetscSection);
134 PETSC_EXTERN PetscErrorCode DMGetDefaultGlobalSection(DM, PetscSection *);
135 PETSC_EXTERN PetscErrorCode DMSetDefaultGlobalSection(DM, PetscSection);
136 PETSC_EXTERN PetscErrorCode DMGetDefaultSF(DM, PetscSF *);
137 PETSC_EXTERN PetscErrorCode DMSetDefaultSF(DM, PetscSF);
138 PETSC_EXTERN PetscErrorCode DMCreateDefaultSF(DM, PetscSection, PetscSection);
139 PETSC_EXTERN PetscErrorCode DMGetPointSF(DM, PetscSF *);
140 PETSC_EXTERN PetscErrorCode DMSetPointSF(DM, PetscSF);
141 
142 PETSC_EXTERN PetscErrorCode DMGetNumFields(DM, PetscInt *);
143 PETSC_EXTERN PetscErrorCode DMSetNumFields(DM, PetscInt);
144 PETSC_EXTERN PetscErrorCode DMGetField(DM, PetscInt, PetscObject *);
145 
146 typedef struct {
147   PetscInt         numQuadPoints; /* The number of quadrature points on an element */
148   const PetscReal *quadPoints;    /* The quadrature point coordinates */
149   const PetscReal *quadWeights;   /* The quadrature weights */
150   PetscInt         numBasisFuncs; /* The number of finite element basis functions on an element */
151   PetscInt         numComponents; /* The number of components for each basis function */
152   const PetscReal *basis;         /* The basis functions tabulated at the quadrature points */
153   const PetscReal *basisDer;      /* The basis function derivatives tabulated at the quadrature points */
154 } PetscQuadrature;
155 
156 typedef struct {
157   PetscQuadrature *quad;
158   void (**f0Funcs)(const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]); /* The f_0 functions for each field */
159   void (**f1Funcs)(const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]); /* The f_1 functions for each field */
160   void (**g0Funcs)(const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]); /* The g_0 functions for each field pair */
161   void (**g1Funcs)(const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]); /* The g_1 functions for each field pair */
162   void (**g2Funcs)(const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]); /* The g_2 functions for each field pair */
163   void (**g3Funcs)(const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]); /* The g_3 functions for each field pair */
164   void (**bcFuncs)(const PetscReal[], PetscScalar *); /* The boundary condition function for each field component */
165   PetscQuadrature *quadBd;
166   void (**f0BdFuncs)(const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]); /* The f_0 functions for each field */
167   void (**f1BdFuncs)(const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]); /* The f_1 functions for each field */
168 } PetscFEM;
169 
170 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;
171 
172 struct _DMInterpolationInfo {
173   MPI_Comm   comm;
174   PetscInt   dim;    /*1 The spatial dimension of points */
175   PetscInt   nInput; /* The number of input points */
176   PetscReal *points; /* The input point coordinates */
177   PetscInt  *cells;  /* The cell containing each point */
178   PetscInt   n;      /* The number of local points */
179   Vec        coords; /* The point coordinates */
180   PetscInt   dof;    /* The number of components to interpolate */
181 };
182 typedef struct _DMInterpolationInfo *DMInterpolationInfo;
183 
184 PetscErrorCode DMInterpolationCreate(MPI_Comm, DMInterpolationInfo *);
185 PetscErrorCode DMInterpolationSetDim(DMInterpolationInfo, PetscInt);
186 PetscErrorCode DMInterpolationGetDim(DMInterpolationInfo, PetscInt *);
187 PetscErrorCode DMInterpolationSetDof(DMInterpolationInfo, PetscInt);
188 PetscErrorCode DMInterpolationGetDof(DMInterpolationInfo, PetscInt *);
189 PetscErrorCode DMInterpolationAddPoints(DMInterpolationInfo, PetscInt, PetscReal[]);
190 PetscErrorCode DMInterpolationSetUp(DMInterpolationInfo, DM, PetscBool);
191 PetscErrorCode DMInterpolationGetCoordinates(DMInterpolationInfo, Vec *);
192 PetscErrorCode DMInterpolationGetVector(DMInterpolationInfo, Vec *);
193 PetscErrorCode DMInterpolationRestoreVector(DMInterpolationInfo, Vec *);
194 PetscErrorCode DMInterpolationEvaluate(DMInterpolationInfo, DM, Vec, Vec);
195 PetscErrorCode DMInterpolationDestroy(DMInterpolationInfo *);
196 #endif
197