xref: /petsc/include/petscdmda.h (revision 3c48a1e8da19189ff2402a4e41a2fc082d52c349)
1 #if !defined(__PETSCDMDA_H)
2 #define __PETSCDMDA_H
3 
4 #include "petscdm.h"
5 #include "petscpf.h"
6 
7 /*E
8     DMDAStencilType - Determines if the stencil extends only along the coordinate directions, or also
9       to the northeast, northwest etc
10 
11    Level: beginner
12 
13 .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDACreate()
14 E*/
15 typedef enum { DMDA_STENCIL_STAR,DMDA_STENCIL_BOX } DMDAStencilType;
16 
17 /*MC
18      DMDA_STENCIL_STAR - "Star"-type stencil. In logical grid coordinates, only (i,j,k), (i+s,j,k), (i,j+s,k),
19                        (i,j,k+s) are in the stencil  NOT, for example, (i+s,j+s,k)
20 
21      Level: beginner
22 
23 .seealso: DMDA_STENCIL_BOX, DMDAStencilType
24 M*/
25 
26 /*MC
27      DMDA_STENCIL_BOX - "Box"-type stencil. In logical grid coordinates, any of (i,j,k), (i+s,j+r,k+t) may
28                       be in the stencil.
29 
30      Level: beginner
31 
32 .seealso: DMDA_STENCIL_STAR, DMDAStencilType
33 M*/
34 
35 /*E
36     DMDABoundaryType - Describes the choice for fill of ghost cells on domain boundaries.
37 
38    Level: beginner
39 
40    A boundary may be of type DMDA_BOUNDARY_NONE (no ghost nodes), DMDA_BOUNDARY_GHOST (ghost nodes
41    exist but aren't filled), DMDA_BOUNDARY_MIRROR (not yet implemented), or DMDA_BOUNDARY_PERIODIC
42    (ghost nodes filled by the opposite edge of the domain).
43 
44 .seealso: DMDASetBoundaryType(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDACreate()
45 E*/
46 typedef enum { DMDA_BOUNDARY_NONE, DMDA_BOUNDARY_GHOSTED, DMDA_BOUNDARY_MIRROR, DMDA_BOUNDARY_PERIODIC } DMDABoundaryType;
47 
48 extern const char *DMDABoundaryTypes[];
49 
50 /*E
51     DMDAInterpolationType - Defines the type of interpolation that will be returned by
52        DMGetInterpolation.
53 
54    Level: beginner
55 
56 .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMGetInterpolation(), DMDASetInterpolationType(), DMDACreate()
57 E*/
58 typedef enum { DMDA_Q0, DMDA_Q1 } DMDAInterpolationType;
59 
60 extern PetscErrorCode   DMDASetInterpolationType(DM,DMDAInterpolationType);
61 
62 /*E
63     DMDAElementType - Defines the type of elements that will be returned by
64        DMGetElements()
65 
66    Level: beginner
67 
68 .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMGetInterpolation(), DMDASetInterpolationType(),
69           DMDASetElementType(), DMGetElements(), DMRestoreElements(), DMDACreate()
70 E*/
71 typedef enum { DMDA_ELEMENT_P1, DMDA_ELEMENT_Q1 } DMDAElementType;
72 
73 extern PetscErrorCode   DMDASetElementType(DM,DMDAElementType);
74 extern PetscErrorCode   DMDAGetElementType(DM,DMDAElementType*);
75 
76 typedef enum { DMDA_X,DMDA_Y,DMDA_Z } DMDADirection;
77 
78 #define MATSEQUSFFT        "sequsfft"
79 
80 extern PetscErrorCode  DMDACreate(MPI_Comm,DM*);
81 extern PetscErrorCode  DMDASetDim(DM,PetscInt);
82 extern PetscErrorCode  DMDASetSizes(DM,PetscInt,PetscInt,PetscInt);
83 extern PetscErrorCode     DMDACreate1d(MPI_Comm,DMDABoundaryType,PetscInt,PetscInt,PetscInt,const PetscInt[],DM *);
84 extern PetscErrorCode     DMDACreate2d(MPI_Comm,DMDABoundaryType,DMDABoundaryType,DMDAStencilType,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],DM*);
85 extern PetscErrorCode     DMDACreate3d(MPI_Comm,DMDABoundaryType,DMDABoundaryType,DMDABoundaryType,DMDAStencilType,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscInt[],DM*);
86 extern PetscErrorCode  DMSetOptionsPrefix(DM,const char []);
87 extern PetscErrorCode  DMSetVecType(DM,const VecType);
88 
89 extern PetscErrorCode     DMDAGlobalToNaturalBegin(DM,Vec,InsertMode,Vec);
90 extern PetscErrorCode     DMDAGlobalToNaturalEnd(DM,Vec,InsertMode,Vec);
91 extern PetscErrorCode     DMDANaturalToGlobalBegin(DM,Vec,InsertMode,Vec);
92 extern PetscErrorCode     DMDANaturalToGlobalEnd(DM,Vec,InsertMode,Vec);
93 extern PetscErrorCode     DMDALocalToLocalBegin(DM,Vec,InsertMode,Vec);
94 extern PetscErrorCode     DMDALocalToLocalEnd(DM,Vec,InsertMode,Vec);
95 extern PetscErrorCode     DMDACreateNaturalVector(DM,Vec *);
96 
97 extern PetscErrorCode     DMDALoad(PetscViewer,PetscInt,PetscInt,PetscInt,DM *);
98 extern PetscErrorCode     DMDAGetCorners(DM,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
99 extern PetscErrorCode     DMDAGetGhostCorners(DM,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
100 extern PetscErrorCode     DMDAGetInfo(DM,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,DMDABoundaryType*,DMDABoundaryType*,DMDABoundaryType*,DMDAStencilType*);
101 extern PetscErrorCode     DMDAGetProcessorSubset(DM,DMDADirection,PetscInt,MPI_Comm*);
102 extern PetscErrorCode     DMDAGetProcessorSubsets(DM,DMDADirection,MPI_Comm*);
103 
104 extern PetscErrorCode     DMDAGlobalToNaturalAllCreate(DM,VecScatter*);
105 extern PetscErrorCode     DMDANaturalAllToGlobalCreate(DM,VecScatter*);
106 
107 extern PetscErrorCode     DMDAGetGlobalIndices(DM,PetscInt*,PetscInt**);
108 
109 extern PetscErrorCode     DMDAGetScatter(DM,VecScatter*,VecScatter*,VecScatter*);
110 extern PetscErrorCode     DMDAGetNeighbors(DM,const PetscMPIInt**);
111 
112 extern PetscErrorCode     DMDAGetAO(DM,AO*);
113 extern PetscErrorCode     DMDASetCoordinates(DM,Vec);
114 extern PetscErrorCode     DMDASetGhostedCoordinates(DM,Vec);
115 extern PetscErrorCode     DMDAGetCoordinates(DM,Vec *);
116 extern PetscErrorCode     DMDAGetGhostedCoordinates(DM,Vec *);
117 extern PetscErrorCode     DMDAGetCoordinateDA(DM,DM *);
118 extern PetscErrorCode     DMDASetUniformCoordinates(DM,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal);
119 extern PetscErrorCode     DMDAGetBoundingBox(DM,PetscReal[],PetscReal[]);
120 extern PetscErrorCode     DMDAGetLocalBoundingBox(DM,PetscReal[],PetscReal[]);
121 extern PetscErrorCode     DMDASetFieldName(DM,PetscInt,const char[]);
122 extern PetscErrorCode     DMDAGetFieldName(DM,PetscInt,const char**);
123 
124 extern PetscErrorCode  DMDASetBoundaryType(DM,DMDABoundaryType,DMDABoundaryType,DMDABoundaryType);
125 extern PetscErrorCode  DMDASetDof(DM, int);
126 extern PetscErrorCode  DMDASetStencilWidth(DM, PetscInt);
127 extern PetscErrorCode  DMDASetOwnershipRanges(DM,const PetscInt[],const PetscInt[],const PetscInt[]);
128 extern PetscErrorCode  DMDAGetOwnershipRanges(DM,const PetscInt**,const PetscInt**,const PetscInt**);
129 extern PetscErrorCode  DMDASetNumProcs(DM, PetscInt, PetscInt, PetscInt);
130 extern PetscErrorCode  DMDASetStencilType(DM, DMDAStencilType);
131 
132 extern PetscErrorCode     DMDAVecGetArray(DM,Vec,void *);
133 extern PetscErrorCode     DMDAVecRestoreArray(DM,Vec,void *);
134 
135 extern PetscErrorCode     DMDAVecGetArrayDOF(DM,Vec,void *);
136 extern PetscErrorCode     DMDAVecRestoreArrayDOF(DM,Vec,void *);
137 
138 extern PetscErrorCode     DMDASplitComm2d(MPI_Comm,PetscInt,PetscInt,PetscInt,MPI_Comm*);
139 
140 /*S
141      DMDALocalInfo - C struct that contains information about a structured grid and a processors logical
142               location in it.
143 
144    Level: beginner
145 
146   Concepts: distributed array
147 
148   Developer note: Then entries in this struct are int instead of PetscInt so that the elements may
149                   be extracted in Fortran as if from an integer array
150 
151 .seealso:  DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DM, DMDAGetLocalInfo(), DMDAGetInfo()
152 S*/
153 typedef struct {
154   PetscInt       dim,dof,sw;
155   PetscInt       mx,my,mz;    /* global number of grid points in each direction */
156   PetscInt       xs,ys,zs;    /* starting pointd of this processor, excluding ghosts */
157   PetscInt       xm,ym,zm;    /* number of grid points on this processor, excluding ghosts */
158   PetscInt       gxs,gys,gzs;    /* starting point of this processor including ghosts */
159   PetscInt       gxm,gym,gzm;    /* number of grid points on this processor including ghosts */
160   DMDABoundaryType bx,by,bz; /* type of ghost nodes at boundary */
161   DMDAStencilType  st;
162   DM             da;
163 } DMDALocalInfo;
164 
165 /*MC
166       DMDAForEachPointBegin2d - Starts a loop over the local part of a two dimensional DMDA
167 
168    Synopsis:
169    void  DMDAForEachPointBegin2d(DALocalInfo *info,PetscInt i,PetscInt j);
170 
171    Not Collective
172 
173    Level: intermediate
174 
175 .seealso: DMDAForEachPointEnd2d(), DMDAVecGetArray()
176 M*/
177 #define DMDAForEachPointBegin2d(info,i,j) {\
178   PetscInt _xints = info->xs,_xinte = info->xs+info->xm,_yints = info->ys,_yinte = info->ys+info->ym;\
179   for (j=_yints; j<_yinte; j++) {\
180     for (i=_xints; i<_xinte; i++) {\
181 
182 /*MC
183       DMDAForEachPointEnd2d - Ends a loop over the local part of a two dimensional DMDA
184 
185    Synopsis:
186    void  DMDAForEachPointEnd2d;
187 
188    Not Collective
189 
190    Level: intermediate
191 
192 .seealso: DMDAForEachPointBegin2d(), DMDAVecGetArray()
193 M*/
194 #define DMDAForEachPointEnd2d }}}
195 
196 /*MC
197       DMDACoor2d - Structure for holding 2d (x and y) coordinates.
198 
199     Level: intermediate
200 
201     Sample Usage:
202       DMDACoor2d **coors;
203       Vec      vcoors;
204       DM       cda;
205 
206       DMDAGetCoordinates(da,&vcoors);
207       DMDAGetCoordinateDA(da,&cda);
208       DMDAVecGetArray(cda,vcoors,&coors);
209       DMDAGetCorners(cda,&mstart,&nstart,0,&m,&n,0)
210       for (i=mstart; i<mstart+m; i++) {
211         for (j=nstart; j<nstart+n; j++) {
212           x = coors[j][i].x;
213           y = coors[j][i].y;
214           ......
215         }
216       }
217       DMDAVecRestoreArray(dac,vcoors,&coors);
218 
219 .seealso: DMDACoor3d, DMDAForEachPointBegin(), DMDAGetCoordinateDA(), DMDAGetCoordinates(), DMDAGetGhostCoordinates()
220 M*/
221 typedef struct {PetscScalar x,y;} DMDACoor2d;
222 
223 /*MC
224       DMDACoor3d - Structure for holding 3d (x, y and z) coordinates.
225 
226     Level: intermediate
227 
228     Sample Usage:
229       DMDACoor3d ***coors;
230       Vec      vcoors;
231       DM       cda;
232 
233       DMDAGetCoordinates(da,&vcoors);
234       DMDAGetCoordinateDA(da,&cda);
235       DMDAVecGetArray(cda,vcoors,&coors);
236       DMDAGetCorners(cda,&mstart,&nstart,&pstart,&m,&n,&p)
237       for (i=mstart; i<mstart+m; i++) {
238         for (j=nstart; j<nstart+n; j++) {
239           for (k=pstart; k<pstart+p; k++) {
240             x = coors[k][j][i].x;
241             y = coors[k][j][i].y;
242             z = coors[k][j][i].z;
243           ......
244         }
245       }
246       DMDAVecRestoreArray(dac,vcoors,&coors);
247 
248 .seealso: DMDACoor2d, DMDAForEachPointBegin(), DMDAGetCoordinateDA(), DMDAGetCoordinates(), DMDAGetGhostCoordinates()
249 M*/
250 typedef struct {PetscScalar x,y,z;} DMDACoor3d;
251 
252 extern PetscErrorCode   DMDAGetLocalInfo(DM,DMDALocalInfo*);
253 typedef PetscErrorCode (*DMDALocalFunction1)(DMDALocalInfo*,void*,void*,void*);
254 extern PetscErrorCode   DMDAFormFunctionLocal(DM, DMDALocalFunction1, Vec, Vec, void *);
255 extern PetscErrorCode   DMDAFormFunctionLocalGhost(DM, DMDALocalFunction1, Vec, Vec, void *);
256 extern PetscErrorCode   DMDAFormJacobianLocal(DM, DMDALocalFunction1, Vec, Mat, void *);
257 extern PetscErrorCode   DMDAFormFunction1(DM,Vec,Vec,void*);
258 extern PetscErrorCode   DMDAFormFunction(DM,PetscErrorCode (*)(void),Vec,Vec,void*);
259 extern PetscErrorCode   DMDAFormFunctioni1(DM,PetscInt,Vec,PetscScalar*,void*);
260 extern PetscErrorCode   DMDAFormFunctionib1(DM,PetscInt,Vec,PetscScalar*,void*);
261 extern PetscErrorCode   DMDAComputeJacobian1WithAdic(DM,Vec,Mat,void*);
262 extern PetscErrorCode   DMDAComputeJacobian1WithAdifor(DM,Vec,Mat,void*);
263 extern PetscErrorCode   DMDAMultiplyByJacobian1WithAdic(DM,Vec,Vec,Vec,void*);
264 extern PetscErrorCode   DMDAMultiplyByJacobian1WithAdifor(DM,Vec,Vec,Vec,void*);
265 extern PetscErrorCode   DMDAMultiplyByJacobian1WithAD(DM,Vec,Vec,Vec,void*);
266 extern PetscErrorCode   DMDAComputeJacobian1(DM,Vec,Mat,void*);
267 extern PetscErrorCode   DMDAGetLocalFunction(DM,DMDALocalFunction1*);
268 extern PetscErrorCode   DMDASetLocalFunction(DM,DMDALocalFunction1);
269 extern PetscErrorCode   DMDASetLocalFunctioni(DM,PetscErrorCode (*)(DMDALocalInfo*,MatStencil*,void*,PetscScalar*,void*));
270 extern PetscErrorCode   DMDASetLocalFunctionib(DM,PetscErrorCode (*)(DMDALocalInfo*,MatStencil*,void*,PetscScalar*,void*));
271 extern PetscErrorCode   DMDAGetLocalJacobian(DM,DMDALocalFunction1*);
272 extern PetscErrorCode   DMDASetLocalJacobian(DM,DMDALocalFunction1);
273 extern PetscErrorCode   DMDASetLocalAdicFunction_Private(DM,DMDALocalFunction1);
274 
275 extern PetscErrorCode MatSetDA(Mat,DM);
276 extern PetscErrorCode MatRegisterDAAD(void);
277 extern PetscErrorCode MatCreateDAAD(DM,Mat*);
278 extern PetscErrorCode MatCreateSeqUSFFT(Vec,DM,Mat*);
279 
280 /*MC
281        DMDASetLocalAdicFunction - Caches in a DM a local function computed by ADIC/ADIFOR
282 
283    Synopsis:
284    PetscErrorCode DMDASetLocalAdicFunction(DM da,DMDALocalFunction1 ad_lf)
285 
286    Logically Collective on DM
287 
288    Input Parameter:
289 +  da - initial distributed array
290 -  ad_lf - the local function as computed by ADIC/ADIFOR
291 
292    Level: intermediate
293 
294 .keywords:  distributed array, refine
295 
296 .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDAGetLocalFunction(), DMDASetLocalFunction(),
297           DMDASetLocalJacobian()
298 M*/
299 #if defined(PETSC_HAVE_ADIC)
300 #  define DMDASetLocalAdicFunction(a,d) DMDASetLocalAdicFunction_Private(a,(DMDALocalFunction1)d)
301 #else
302 #  define DMDASetLocalAdicFunction(a,d) DMDASetLocalAdicFunction_Private(a,0)
303 #endif
304 
305 extern PetscErrorCode   DMDASetLocalAdicMFFunction_Private(DM,DMDALocalFunction1);
306 #if defined(PETSC_HAVE_ADIC)
307 #  define DMDASetLocalAdicMFFunction(a,d) DMDASetLocalAdicMFFunction_Private(a,(DMDALocalFunction1)d)
308 #else
309 #  define DMDASetLocalAdicMFFunction(a,d) DMDASetLocalAdicMFFunction_Private(a,0)
310 #endif
311 extern PetscErrorCode   DMDASetLocalAdicFunctioni_Private(DM,PetscErrorCode (*)(DMDALocalInfo*,MatStencil*,void*,void*,void*));
312 #if defined(PETSC_HAVE_ADIC)
313 #  define DMDASetLocalAdicFunctioni(a,d) DMDASetLocalAdicFunctioni_Private(a,(PetscErrorCode (*)(DMDALocalInfo*,MatStencil*,void*,void*,void*))d)
314 #else
315 #  define DMDASetLocalAdicFunctioni(a,d) DMDASetLocalAdicFunctioni_Private(a,0)
316 #endif
317 extern PetscErrorCode   DMDASetLocalAdicMFFunctioni_Private(DM,PetscErrorCode (*)(DMDALocalInfo*,MatStencil*,void*,void*,void*));
318 #if defined(PETSC_HAVE_ADIC)
319 #  define DMDASetLocalAdicMFFunctioni(a,d) DMDASetLocalAdicMFFunctioni_Private(a,(PetscErrorCode (*)(DMDALocalInfo*,MatStencil*,void*,void*,void*))d)
320 #else
321 #  define DMDASetLocalAdicMFFunctioni(a,d) DMDASetLocalAdicMFFunctioni_Private(a,0)
322 #endif
323 
324 extern PetscErrorCode   DMDASetLocalAdicFunctionib_Private(DM,PetscErrorCode (*)(DMDALocalInfo*,MatStencil*,void*,void*,void*));
325 #if defined(PETSC_HAVE_ADIC)
326 #  define DMDASetLocalAdicFunctionib(a,d) DMDASetLocalAdicFunctionib_Private(a,(PetscErrorCode (*)(DMDALocalInfo*,MatStencil*,void*,void*,void*))d)
327 #else
328 #  define DMDASetLocalAdicFunctionib(a,d) DMDASetLocalAdicFunctionib_Private(a,0)
329 #endif
330 extern PetscErrorCode   DMDASetLocalAdicMFFunctionib_Private(DM,PetscErrorCode (*)(DMDALocalInfo*,MatStencil*,void*,void*,void*));
331 #if defined(PETSC_HAVE_ADIC)
332 #  define DMDASetLocalAdicMFFunctionib(a,d) DMDASetLocalAdicMFFunctionib_Private(a,(PetscErrorCode (*)(DMDALocalInfo*,MatStencil*,void*,void*,void*))d)
333 #else
334 #  define DMDASetLocalAdicMFFunctionib(a,d) DMDASetLocalAdicMFFunctionib_Private(a,0)
335 #endif
336 
337 extern PetscErrorCode   DMDAFormFunctioniTest1(DM,void*);
338 extern PetscErrorCode   DMDASetGetMatrix(DM,PetscErrorCode (*)(DM, const MatType,Mat *));
339 extern PetscErrorCode   DMDASetBlockFills(DM,PetscInt*,PetscInt*);
340 extern PetscErrorCode   DMDASetMatPreallocateOnly(DM,PetscBool );
341 extern PetscErrorCode   DMDASetRefinementFactor(DM,PetscInt,PetscInt,PetscInt);
342 extern PetscErrorCode   DMDAGetRefinementFactor(DM,PetscInt*,PetscInt*,PetscInt*);
343 
344 extern PetscErrorCode   DMDAGetAdicArray(DM,PetscBool ,void*,void*,PetscInt*);
345 extern PetscErrorCode   DMDARestoreAdicArray(DM,PetscBool ,void*,void*,PetscInt*);
346 extern PetscErrorCode   DMDAGetAdicMFArray(DM,PetscBool ,void*,void*,PetscInt*);
347 extern PetscErrorCode   DMDAGetAdicMFArray4(DM,PetscBool ,void*,void*,PetscInt*);
348 extern PetscErrorCode   DMDAGetAdicMFArray9(DM,PetscBool ,void*,void*,PetscInt*);
349 extern PetscErrorCode   DMDAGetAdicMFArrayb(DM,PetscBool ,void*,void*,PetscInt*);
350 extern PetscErrorCode   DMDARestoreAdicMFArray(DM,PetscBool ,void*,void*,PetscInt*);
351 extern PetscErrorCode   DMDAGetArray(DM,PetscBool ,void*);
352 extern PetscErrorCode   DMDARestoreArray(DM,PetscBool ,void*);
353 extern PetscErrorCode   ad_DAGetArray(DM,PetscBool ,void*);
354 extern PetscErrorCode   ad_DARestoreArray(DM,PetscBool ,void*);
355 extern PetscErrorCode   admf_DAGetArray(DM,PetscBool ,void*);
356 extern PetscErrorCode   admf_DARestoreArray(DM,PetscBool ,void*);
357 
358 extern PetscErrorCode   DMDACreatePF(DM,PF*);
359 
360 #endif
361