xref: /petsc/include/petscdm.h (revision 09fe277da1f73a743f13ea3855cd2bdc4d5de01f)
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 #include <petscfetypes.h>
9 #include <petscdstypes.h>
10 #include <petscdmlabel.h>
11 
12 PETSC_EXTERN PetscErrorCode DMInitializePackage(void);
13 
14 PETSC_EXTERN PetscClassId DM_CLASSID;
15 
16 #define DMLOCATEPOINT_POINT_NOT_FOUND -367
17 
18 /*J
19     DMType - String with the name of a PETSc DM
20 
21    Level: beginner
22 
23 .seealso: DMSetType(), DM
24 J*/
25 typedef const char* DMType;
26 #define DMDA        "da"
27 #define DMCOMPOSITE "composite"
28 #define DMSLICED    "sliced"
29 #define DMSHELL     "shell"
30 #define DMPLEX      "plex"
31 #define DMCARTESIAN "cartesian"
32 #define DMREDUNDANT "redundant"
33 #define DMPATCH     "patch"
34 #define DMMOAB      "moab"
35 #define DMNETWORK   "network"
36 #define DMFOREST    "forest"
37 #define DMP4EST     "p4est"
38 #define DMP8EST     "p8est"
39 #define DMSWARM     "swarm"
40 
41 PETSC_EXTERN const char *const DMBoundaryTypes[];
42 PETSC_EXTERN PetscFunctionList DMList;
43 PETSC_EXTERN PetscErrorCode DMCreate(MPI_Comm,DM*);
44 PETSC_EXTERN PetscErrorCode DMClone(DM,DM*);
45 PETSC_EXTERN PetscErrorCode DMSetType(DM, DMType);
46 PETSC_EXTERN PetscErrorCode DMGetType(DM, DMType *);
47 PETSC_EXTERN PetscErrorCode DMRegister(const char[],PetscErrorCode (*)(DM));
48 PETSC_EXTERN PetscErrorCode DMRegisterDestroy(void);
49 
50 PETSC_EXTERN PetscErrorCode DMView(DM,PetscViewer);
51 PETSC_EXTERN PetscErrorCode DMLoad(DM,PetscViewer);
52 PETSC_EXTERN PetscErrorCode DMDestroy(DM*);
53 PETSC_EXTERN PetscErrorCode DMCreateGlobalVector(DM,Vec*);
54 PETSC_EXTERN PetscErrorCode DMCreateLocalVector(DM,Vec*);
55 PETSC_EXTERN PetscErrorCode DMGetLocalVector(DM,Vec *);
56 PETSC_EXTERN PetscErrorCode DMRestoreLocalVector(DM,Vec *);
57 PETSC_EXTERN PetscErrorCode DMGetGlobalVector(DM,Vec *);
58 PETSC_EXTERN PetscErrorCode DMRestoreGlobalVector(DM,Vec *);
59 PETSC_EXTERN PetscErrorCode DMClearGlobalVectors(DM);
60 PETSC_EXTERN PetscErrorCode DMClearLocalVectors(DM);
61 PETSC_EXTERN PetscErrorCode DMHasNamedGlobalVector(DM,const char*,PetscBool*);
62 PETSC_EXTERN PetscErrorCode DMGetNamedGlobalVector(DM,const char*,Vec*);
63 PETSC_EXTERN PetscErrorCode DMRestoreNamedGlobalVector(DM,const char*,Vec*);
64 PETSC_EXTERN PetscErrorCode DMHasNamedLocalVector(DM,const char*,PetscBool*);
65 PETSC_EXTERN PetscErrorCode DMGetNamedLocalVector(DM,const char*,Vec*);
66 PETSC_EXTERN PetscErrorCode DMRestoreNamedLocalVector(DM,const char*,Vec*);
67 PETSC_EXTERN PetscErrorCode DMGetLocalToGlobalMapping(DM,ISLocalToGlobalMapping*);
68 PETSC_EXTERN PetscErrorCode DMCreateFieldIS(DM,PetscInt*,char***,IS**);
69 PETSC_EXTERN PetscErrorCode DMGetBlockSize(DM,PetscInt*);
70 PETSC_EXTERN PetscErrorCode DMCreateColoring(DM,ISColoringType,ISColoring*);
71 PETSC_EXTERN PetscErrorCode DMCreateMatrix(DM,Mat*);
72 PETSC_EXTERN PetscErrorCode DMSetMatrixPreallocateOnly(DM,PetscBool);
73 PETSC_EXTERN PetscErrorCode DMCreateInterpolation(DM,DM,Mat*,Vec*);
74 PETSC_EXTERN PetscErrorCode DMCreateRestriction(DM,DM,Mat*);
75 PETSC_EXTERN PetscErrorCode DMCreateInjection(DM,DM,Mat*);
76 PETSC_EXTERN PetscErrorCode DMGetWorkArray(DM,PetscInt,PetscDataType,void*);
77 PETSC_EXTERN PetscErrorCode DMRestoreWorkArray(DM,PetscInt,PetscDataType,void*);
78 PETSC_EXTERN PetscErrorCode DMRefine(DM,MPI_Comm,DM*);
79 PETSC_EXTERN PetscErrorCode DMCoarsen(DM,MPI_Comm,DM*);
80 PETSC_EXTERN PetscErrorCode DMGetCoarseDM(DM,DM*);
81 PETSC_EXTERN PetscErrorCode DMSetCoarseDM(DM,DM);
82 PETSC_EXTERN PetscErrorCode DMGetFineDM(DM,DM*);
83 PETSC_EXTERN PetscErrorCode DMSetFineDM(DM,DM);
84 PETSC_EXTERN PetscErrorCode DMRefineHierarchy(DM,PetscInt,DM[]);
85 PETSC_EXTERN PetscErrorCode DMCoarsenHierarchy(DM,PetscInt,DM[]);
86 PETSC_EXTERN PetscErrorCode DMCoarsenHookAdd(DM,PetscErrorCode (*)(DM,DM,void*),PetscErrorCode (*)(DM,Mat,Vec,Mat,DM,void*),void*);
87 PETSC_EXTERN PetscErrorCode DMRefineHookAdd(DM,PetscErrorCode (*)(DM,DM,void*),PetscErrorCode (*)(DM,Mat,DM,void*),void*);
88 PETSC_EXTERN PetscErrorCode DMRestrict(DM,Mat,Vec,Mat,DM);
89 PETSC_EXTERN PetscErrorCode DMInterpolate(DM,Mat,DM);
90 PETSC_EXTERN PetscErrorCode DMSetFromOptions(DM);
91 PETSC_STATIC_INLINE PetscErrorCode DMViewFromOptions(DM A,PetscObject obj,const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,obj,name);}
92 
93 PETSC_EXTERN PetscErrorCode DMSetUp(DM);
94 PETSC_EXTERN PetscErrorCode DMCreateInterpolationScale(DM,DM,Mat,Vec*);
95 PETSC_EXTERN PetscErrorCode DMCreateAggregates(DM,DM,Mat*);
96 PETSC_EXTERN PetscErrorCode DMGlobalToLocalHookAdd(DM,PetscErrorCode (*)(DM,Vec,InsertMode,Vec,void*),PetscErrorCode (*)(DM,Vec,InsertMode,Vec,void*),void*);
97 PETSC_EXTERN PetscErrorCode DMLocalToGlobalHookAdd(DM,PetscErrorCode (*)(DM,Vec,InsertMode,Vec,void*),PetscErrorCode (*)(DM,Vec,InsertMode,Vec,void*),void*);
98 PETSC_EXTERN PetscErrorCode DMGlobalToLocalBegin(DM,Vec,InsertMode,Vec);
99 PETSC_EXTERN PetscErrorCode DMGlobalToLocalEnd(DM,Vec,InsertMode,Vec);
100 PETSC_EXTERN PetscErrorCode DMLocalToGlobalBegin(DM,Vec,InsertMode,Vec);
101 PETSC_EXTERN PetscErrorCode DMLocalToGlobalEnd(DM,Vec,InsertMode,Vec);
102 PETSC_EXTERN PetscErrorCode DMLocalToLocalBegin(DM,Vec,InsertMode,Vec);
103 PETSC_EXTERN PetscErrorCode DMLocalToLocalEnd(DM,Vec,InsertMode,Vec);
104 PETSC_EXTERN PetscErrorCode DMConvert(DM,DMType,DM*);
105 
106 /* Topology support */
107 PETSC_EXTERN PetscErrorCode DMGetDimension(DM,PetscInt*);
108 PETSC_EXTERN PetscErrorCode DMSetDimension(DM,PetscInt);
109 PETSC_EXTERN PetscErrorCode DMGetDimPoints(DM,PetscInt,PetscInt*,PetscInt*);
110 PETSC_EXTERN PetscErrorCode DMGetUseNatural(DM,PetscBool*);
111 PETSC_EXTERN PetscErrorCode DMSetUseNatural(DM,PetscBool);
112 
113 /* Coordinate support */
114 PETSC_EXTERN PetscErrorCode DMGetCoordinateDM(DM,DM*);
115 PETSC_EXTERN PetscErrorCode DMSetCoordinateDM(DM,DM);
116 PETSC_EXTERN PetscErrorCode DMGetCoordinateDim(DM,PetscInt*);
117 PETSC_EXTERN PetscErrorCode DMSetCoordinateDim(DM,PetscInt);
118 PETSC_EXTERN PetscErrorCode DMGetCoordinateSection(DM,PetscSection*);
119 PETSC_EXTERN PetscErrorCode DMSetCoordinateSection(DM,PetscInt,PetscSection);
120 PETSC_EXTERN PetscErrorCode DMGetCoordinates(DM,Vec*);
121 PETSC_EXTERN PetscErrorCode DMSetCoordinates(DM,Vec);
122 PETSC_EXTERN PetscErrorCode DMGetCoordinatesLocal(DM,Vec*);
123 PETSC_EXTERN PetscErrorCode DMSetCoordinatesLocal(DM,Vec);
124 PETSC_EXTERN PetscErrorCode DMLocatePoints(DM,Vec,DMPointLocationType,PetscSF*);
125 PETSC_EXTERN PetscErrorCode DMGetPeriodicity(DM,const PetscReal**,const PetscReal**,const DMBoundaryType**);
126 PETSC_EXTERN PetscErrorCode DMSetPeriodicity(DM,const PetscReal[],const PetscReal[],const DMBoundaryType[]);
127 PETSC_EXTERN PetscErrorCode DMLocalizeCoordinate(DM, const PetscScalar[], PetscScalar[]);
128 PETSC_EXTERN PetscErrorCode DMLocalizeCoordinates(DM);
129 PETSC_EXTERN PetscErrorCode DMGetCoordinatesLocalized(DM,PetscBool*);
130 PETSC_EXTERN PetscErrorCode DMGetNeighbors(DM,PetscInt*,const PetscMPIInt**);
131 
132 /* block hook interface */
133 PETSC_EXTERN PetscErrorCode DMSubDomainHookAdd(DM,PetscErrorCode (*)(DM,DM,void*),PetscErrorCode (*)(DM,VecScatter,VecScatter,DM,void*),void*);
134 PETSC_EXTERN PetscErrorCode DMSubDomainRestrict(DM,VecScatter,VecScatter,DM);
135 
136 PETSC_EXTERN PetscErrorCode DMSetOptionsPrefix(DM,const char []);
137 PETSC_EXTERN PetscErrorCode DMAppendOptionsPrefix(DM,const char []);
138 PETSC_EXTERN PetscErrorCode DMGetOptionsPrefix(DM,const char*[]);
139 PETSC_EXTERN PetscErrorCode DMSetVecType(DM,VecType);
140 PETSC_EXTERN PetscErrorCode DMGetVecType(DM,VecType*);
141 PETSC_EXTERN PetscErrorCode DMSetMatType(DM,MatType);
142 PETSC_EXTERN PetscErrorCode DMGetMatType(DM,MatType*);
143 PETSC_EXTERN PetscErrorCode DMSetApplicationContext(DM,void*);
144 PETSC_EXTERN PetscErrorCode DMSetApplicationContextDestroy(DM,PetscErrorCode (*)(void**));
145 PETSC_EXTERN PetscErrorCode DMGetApplicationContext(DM,void*);
146 PETSC_EXTERN PetscErrorCode DMSetVariableBounds(DM,PetscErrorCode (*)(DM,Vec,Vec));
147 PETSC_EXTERN PetscErrorCode DMHasVariableBounds(DM,PetscBool *);
148 PETSC_EXTERN PetscErrorCode DMHasColoring(DM,PetscBool *);
149 PETSC_EXTERN PetscErrorCode DMHasCreateRestriction(DM,PetscBool *);
150 PETSC_EXTERN PetscErrorCode DMComputeVariableBounds(DM,Vec,Vec);
151 
152 PETSC_EXTERN PetscErrorCode DMCreateSubDM(DM, PetscInt, PetscInt[], IS *, DM *);
153 PETSC_EXTERN PetscErrorCode DMCreateFieldDecomposition(DM,PetscInt*,char***,IS**,DM**);
154 PETSC_EXTERN PetscErrorCode DMCreateDomainDecomposition(DM,PetscInt*,char***,IS**,IS**,DM**);
155 PETSC_EXTERN PetscErrorCode DMCreateDomainDecompositionScatters(DM,PetscInt,DM*,VecScatter**,VecScatter**,VecScatter**);
156 
157 PETSC_EXTERN PetscErrorCode DMGetRefineLevel(DM,PetscInt*);
158 PETSC_EXTERN PetscErrorCode DMSetRefineLevel(DM,PetscInt);
159 PETSC_EXTERN PetscErrorCode DMGetCoarsenLevel(DM,PetscInt*);
160 PETSC_EXTERN PetscErrorCode DMFinalizePackage(void);
161 
162 PETSC_EXTERN PetscErrorCode VecGetDM(Vec, DM*);
163 PETSC_EXTERN PetscErrorCode VecSetDM(Vec, DM);
164 PETSC_EXTERN PetscErrorCode MatGetDM(Mat, DM*);
165 PETSC_EXTERN PetscErrorCode MatSetDM(Mat, DM);
166 
167 typedef struct NLF_DAAD* NLF;
168 
169 #define DM_FILE_CLASSID 1211221
170 
171 /* FEM support */
172 PETSC_EXTERN PetscErrorCode DMPrintCellVector(PetscInt, const char [], PetscInt, const PetscScalar []);
173 PETSC_EXTERN PetscErrorCode DMPrintCellMatrix(PetscInt, const char [], PetscInt, PetscInt, const PetscScalar []);
174 PETSC_EXTERN PetscErrorCode DMPrintLocalVec(DM, const char [], PetscReal, Vec);
175 
176 PETSC_EXTERN PetscErrorCode DMGetDefaultSection(DM, PetscSection *);
177 PETSC_EXTERN PetscErrorCode DMSetDefaultSection(DM, PetscSection);
178 PETSC_EXTERN PetscErrorCode DMGetDefaultConstraints(DM, PetscSection *, Mat *);
179 PETSC_EXTERN PetscErrorCode DMSetDefaultConstraints(DM, PetscSection, Mat);
180 PETSC_EXTERN PetscErrorCode DMGetDefaultGlobalSection(DM, PetscSection *);
181 PETSC_EXTERN PetscErrorCode DMSetDefaultGlobalSection(DM, PetscSection);
182 PETSC_EXTERN PetscErrorCode DMGetDefaultSF(DM, PetscSF *);
183 PETSC_EXTERN PetscErrorCode DMSetDefaultSF(DM, PetscSF);
184 PETSC_EXTERN PetscErrorCode DMCreateDefaultSF(DM, PetscSection, PetscSection);
185 PETSC_EXTERN PetscErrorCode DMGetPointSF(DM, PetscSF *);
186 PETSC_EXTERN PetscErrorCode DMSetPointSF(DM, PetscSF);
187 
188 PETSC_EXTERN PetscErrorCode DMGetOutputDM(DM, DM *);
189 PETSC_EXTERN PetscErrorCode DMGetOutputSequenceNumber(DM, PetscInt *, PetscReal *);
190 PETSC_EXTERN PetscErrorCode DMSetOutputSequenceNumber(DM, PetscInt, PetscReal);
191 PETSC_EXTERN PetscErrorCode DMOutputSequenceLoad(DM, PetscViewer, const char *, PetscInt, PetscReal *);
192 
193 PETSC_EXTERN PetscErrorCode DMGetDS(DM, PetscDS *);
194 PETSC_EXTERN PetscErrorCode DMSetDS(DM, PetscDS);
195 PETSC_EXTERN PetscErrorCode DMGetNumFields(DM, PetscInt *);
196 PETSC_EXTERN PetscErrorCode DMSetNumFields(DM, PetscInt);
197 PETSC_EXTERN PetscErrorCode DMGetField(DM, PetscInt, PetscObject *);
198 PETSC_EXTERN PetscErrorCode DMSetField(DM, PetscInt, PetscObject);
199 
200 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;
201 
202 struct _DMInterpolationInfo {
203   MPI_Comm   comm;
204   PetscInt   dim;    /*1 The spatial dimension of points */
205   PetscInt   nInput; /* The number of input points */
206   PetscReal *points; /* The input point coordinates */
207   PetscInt  *cells;  /* The cell containing each point */
208   PetscInt   n;      /* The number of local points */
209   Vec        coords; /* The point coordinates */
210   PetscInt   dof;    /* The number of components to interpolate */
211 };
212 typedef struct _DMInterpolationInfo *DMInterpolationInfo;
213 
214 PETSC_EXTERN PetscErrorCode DMInterpolationCreate(MPI_Comm, DMInterpolationInfo *);
215 PETSC_EXTERN PetscErrorCode DMInterpolationSetDim(DMInterpolationInfo, PetscInt);
216 PETSC_EXTERN PetscErrorCode DMInterpolationGetDim(DMInterpolationInfo, PetscInt *);
217 PETSC_EXTERN PetscErrorCode DMInterpolationSetDof(DMInterpolationInfo, PetscInt);
218 PETSC_EXTERN PetscErrorCode DMInterpolationGetDof(DMInterpolationInfo, PetscInt *);
219 PETSC_EXTERN PetscErrorCode DMInterpolationAddPoints(DMInterpolationInfo, PetscInt, PetscReal[]);
220 PETSC_EXTERN PetscErrorCode DMInterpolationSetUp(DMInterpolationInfo, DM, PetscBool);
221 PETSC_EXTERN PetscErrorCode DMInterpolationGetCoordinates(DMInterpolationInfo, Vec *);
222 PETSC_EXTERN PetscErrorCode DMInterpolationGetVector(DMInterpolationInfo, Vec *);
223 PETSC_EXTERN PetscErrorCode DMInterpolationRestoreVector(DMInterpolationInfo, Vec *);
224 PETSC_EXTERN PetscErrorCode DMInterpolationEvaluate(DMInterpolationInfo, DM, Vec, Vec);
225 PETSC_EXTERN PetscErrorCode DMInterpolationDestroy(DMInterpolationInfo *);
226 
227 PETSC_EXTERN PetscErrorCode DMCreateLabel(DM, const char []);
228 PETSC_EXTERN PetscErrorCode DMGetLabelValue(DM, const char[], PetscInt, PetscInt *);
229 PETSC_EXTERN PetscErrorCode DMSetLabelValue(DM, const char[], PetscInt, PetscInt);
230 PETSC_EXTERN PetscErrorCode DMClearLabelValue(DM, const char[], PetscInt, PetscInt);
231 PETSC_EXTERN PetscErrorCode DMGetLabelSize(DM, const char[], PetscInt *);
232 PETSC_EXTERN PetscErrorCode DMGetLabelIdIS(DM, const char[], IS *);
233 PETSC_EXTERN PetscErrorCode DMGetStratumSize(DM, const char [], PetscInt, PetscInt *);
234 PETSC_EXTERN PetscErrorCode DMGetStratumIS(DM, const char [], PetscInt, IS *);
235 PETSC_EXTERN PetscErrorCode DMClearLabelStratum(DM, const char[], PetscInt);
236 PETSC_EXTERN PetscErrorCode DMGetLabelOutput(DM, const char[], PetscBool *);
237 PETSC_EXTERN PetscErrorCode DMSetLabelOutput(DM, const char[], PetscBool);
238 
239 PETSC_EXTERN PetscErrorCode DMGetNumLabels(DM, PetscInt *);
240 PETSC_EXTERN PetscErrorCode DMGetLabelName(DM, PetscInt, const char **);
241 PETSC_EXTERN PetscErrorCode DMHasLabel(DM, const char [], PetscBool *);
242 PETSC_EXTERN PetscErrorCode DMGetLabel(DM, const char *, DMLabel *);
243 PETSC_EXTERN PetscErrorCode DMGetLabelByNum(DM, PetscInt, DMLabel *);
244 PETSC_EXTERN PetscErrorCode DMAddLabel(DM, DMLabel);
245 PETSC_EXTERN PetscErrorCode DMRemoveLabel(DM, const char [], DMLabel *);
246 PETSC_EXTERN PetscErrorCode DMCopyLabels(DM, DM);
247 
248 PETSC_EXTERN PetscErrorCode DMAddBoundary(DM, PetscBool, const char[], const char[], PetscInt, PetscInt, const PetscInt *, void (*)(), PetscInt, const PetscInt *, void *);
249 PETSC_EXTERN PetscErrorCode DMGetNumBoundary(DM, PetscInt *);
250 PETSC_EXTERN PetscErrorCode DMGetBoundary(DM, PetscInt, PetscBool *, const char **, const char **, PetscInt *, PetscInt *, const PetscInt **, void (**)(), PetscInt *, const PetscInt **, void **);
251 PETSC_EXTERN PetscErrorCode DMIsBoundaryPoint(DM, PetscInt, PetscBool *);
252 PETSC_EXTERN PetscErrorCode DMCopyBoundary(DM, DM);
253 
254 PETSC_EXTERN PetscErrorCode DMProjectFunction(DM,PetscReal,PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal[],PetscInt,PetscScalar *,void *),void**,InsertMode,Vec);
255 PETSC_EXTERN PetscErrorCode DMProjectFunctionLocal(DM,PetscReal,PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal[],PetscInt,PetscScalar *,void *),void**,InsertMode,Vec);
256 PETSC_EXTERN PetscErrorCode DMProjectFunctionLabelLocal(DM,PetscReal,DMLabel,PetscInt,const PetscInt[],PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal[],PetscInt,PetscScalar *,void *),void **,InsertMode,Vec);
257 PETSC_EXTERN PetscErrorCode DMProjectFieldLocal(DM,Vec,void (**)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,const PetscReal[],PetscScalar[]),InsertMode,Vec);
258 PETSC_EXTERN PetscErrorCode DMComputeL2Diff(DM,PetscReal,PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal[],PetscInt,PetscScalar *,void *),void **,Vec,PetscReal *);
259 PETSC_EXTERN PetscErrorCode DMComputeL2GradientDiff(DM, PetscReal, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal [], const PetscReal [], PetscInt, PetscScalar *, void *), void **, Vec, const PetscReal [], PetscReal *);
260 PETSC_EXTERN PetscErrorCode DMComputeL2FieldDiff(DM,PetscReal,PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal[],PetscInt,PetscScalar *,void *),void **,Vec,PetscReal *);
261 #endif
262