xref: /petsc/include/petscdmda.h (revision d9822059a04b3bd78629edfbf3ebfea74b8f6e57)
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 
87 extern PetscErrorCode     DMDAGlobalToNaturalBegin(DM,Vec,InsertMode,Vec);
88 extern PetscErrorCode     DMDAGlobalToNaturalEnd(DM,Vec,InsertMode,Vec);
89 extern PetscErrorCode     DMDANaturalToGlobalBegin(DM,Vec,InsertMode,Vec);
90 extern PetscErrorCode     DMDANaturalToGlobalEnd(DM,Vec,InsertMode,Vec);
91 extern PetscErrorCode     DMDALocalToLocalBegin(DM,Vec,InsertMode,Vec);
92 extern PetscErrorCode     DMDALocalToLocalEnd(DM,Vec,InsertMode,Vec);
93 extern PetscErrorCode     DMDACreateNaturalVector(DM,Vec *);
94 
95 extern PetscErrorCode     DMDALoad(PetscViewer,PetscInt,PetscInt,PetscInt,DM *);
96 extern PetscErrorCode     DMDAGetCorners(DM,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
97 extern PetscErrorCode     DMDAGetGhostCorners(DM,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
98 extern PetscErrorCode     DMDAGetInfo(DM,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,DMDABoundaryType*,DMDABoundaryType*,DMDABoundaryType*,DMDAStencilType*);
99 extern PetscErrorCode     DMDAGetProcessorSubset(DM,DMDADirection,PetscInt,MPI_Comm*);
100 extern PetscErrorCode     DMDAGetProcessorSubsets(DM,DMDADirection,MPI_Comm*);
101 
102 extern PetscErrorCode     DMDAGlobalToNaturalAllCreate(DM,VecScatter*);
103 extern PetscErrorCode     DMDANaturalAllToGlobalCreate(DM,VecScatter*);
104 
105 extern PetscErrorCode     DMDAGetGlobalIndices(DM,PetscInt*,PetscInt**);
106 
107 extern PetscErrorCode     DMDAGetScatter(DM,VecScatter*,VecScatter*,VecScatter*);
108 extern PetscErrorCode     DMDAGetNeighbors(DM,const PetscMPIInt**);
109 
110 extern PetscErrorCode     DMDAGetAO(DM,AO*);
111 extern PetscErrorCode     DMDASetCoordinates(DM,Vec);
112 extern PetscErrorCode     DMDASetGhostedCoordinates(DM,Vec);
113 extern PetscErrorCode     DMDAGetCoordinates(DM,Vec *);
114 extern PetscErrorCode     DMDAGetGhostedCoordinates(DM,Vec *);
115 extern PetscErrorCode     DMDAGetCoordinateDA(DM,DM *);
116 extern PetscErrorCode     DMDASetUniformCoordinates(DM,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal);
117 extern PetscErrorCode     DMDAGetBoundingBox(DM,PetscReal[],PetscReal[]);
118 extern PetscErrorCode     DMDAGetLocalBoundingBox(DM,PetscReal[],PetscReal[]);
119 extern PetscErrorCode     DMDASetFieldName(DM,PetscInt,const char[]);
120 extern PetscErrorCode     DMDAGetFieldName(DM,PetscInt,const char**);
121 
122 extern PetscErrorCode  DMDASetBoundaryType(DM,DMDABoundaryType,DMDABoundaryType,DMDABoundaryType);
123 extern PetscErrorCode  DMDASetDof(DM, int);
124 extern PetscErrorCode  DMDASetStencilWidth(DM, PetscInt);
125 extern PetscErrorCode  DMDASetOwnershipRanges(DM,const PetscInt[],const PetscInt[],const PetscInt[]);
126 extern PetscErrorCode  DMDAGetOwnershipRanges(DM,const PetscInt**,const PetscInt**,const PetscInt**);
127 extern PetscErrorCode  DMDASetNumProcs(DM, PetscInt, PetscInt, PetscInt);
128 extern PetscErrorCode  DMDASetStencilType(DM, DMDAStencilType);
129 
130 extern PetscErrorCode     DMDAVecGetArray(DM,Vec,void *);
131 extern PetscErrorCode     DMDAVecRestoreArray(DM,Vec,void *);
132 
133 extern PetscErrorCode     DMDAVecGetArrayDOF(DM,Vec,void *);
134 extern PetscErrorCode     DMDAVecRestoreArrayDOF(DM,Vec,void *);
135 
136 extern PetscErrorCode     DMDASplitComm2d(MPI_Comm,PetscInt,PetscInt,PetscInt,MPI_Comm*);
137 
138 /*S
139      DMDALocalInfo - C struct that contains information about a structured grid and a processors logical
140               location in it.
141 
142    Level: beginner
143 
144   Concepts: distributed array
145 
146   Developer note: Then entries in this struct are int instead of PetscInt so that the elements may
147                   be extracted in Fortran as if from an integer array
148 
149 .seealso:  DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DM, DMDAGetLocalInfo(), DMDAGetInfo()
150 S*/
151 typedef struct {
152   PetscInt       dim,dof,sw;
153   PetscInt       mx,my,mz;    /* global number of grid points in each direction */
154   PetscInt       xs,ys,zs;    /* starting pointd of this processor, excluding ghosts */
155   PetscInt       xm,ym,zm;    /* number of grid points on this processor, excluding ghosts */
156   PetscInt       gxs,gys,gzs;    /* starting point of this processor including ghosts */
157   PetscInt       gxm,gym,gzm;    /* number of grid points on this processor including ghosts */
158   DMDABoundaryType bx,by,bz; /* type of ghost nodes at boundary */
159   DMDAStencilType  st;
160   DM             da;
161 } DMDALocalInfo;
162 
163 /*MC
164       DMDAForEachPointBegin2d - Starts a loop over the local part of a two dimensional DMDA
165 
166    Synopsis:
167    void  DMDAForEachPointBegin2d(DALocalInfo *info,PetscInt i,PetscInt j);
168 
169    Not Collective
170 
171    Level: intermediate
172 
173 .seealso: DMDAForEachPointEnd2d(), DMDAVecGetArray()
174 M*/
175 #define DMDAForEachPointBegin2d(info,i,j) {\
176   PetscInt _xints = info->xs,_xinte = info->xs+info->xm,_yints = info->ys,_yinte = info->ys+info->ym;\
177   for (j=_yints; j<_yinte; j++) {\
178     for (i=_xints; i<_xinte; i++) {\
179 
180 /*MC
181       DMDAForEachPointEnd2d - Ends a loop over the local part of a two dimensional DMDA
182 
183    Synopsis:
184    void  DMDAForEachPointEnd2d;
185 
186    Not Collective
187 
188    Level: intermediate
189 
190 .seealso: DMDAForEachPointBegin2d(), DMDAVecGetArray()
191 M*/
192 #define DMDAForEachPointEnd2d }}}
193 
194 /*MC
195       DMDACoor2d - Structure for holding 2d (x and y) coordinates.
196 
197     Level: intermediate
198 
199     Sample Usage:
200       DMDACoor2d **coors;
201       Vec      vcoors;
202       DM       cda;
203 
204       DMDAGetCoordinates(da,&vcoors);
205       DMDAGetCoordinateDA(da,&cda);
206       DMDAVecGetArray(cda,vcoors,&coors);
207       DMDAGetCorners(cda,&mstart,&nstart,0,&m,&n,0)
208       for (i=mstart; i<mstart+m; i++) {
209         for (j=nstart; j<nstart+n; j++) {
210           x = coors[j][i].x;
211           y = coors[j][i].y;
212           ......
213         }
214       }
215       DMDAVecRestoreArray(dac,vcoors,&coors);
216 
217 .seealso: DMDACoor3d, DMDAForEachPointBegin(), DMDAGetCoordinateDA(), DMDAGetCoordinates(), DMDAGetGhostCoordinates()
218 M*/
219 typedef struct {PetscScalar x,y;} DMDACoor2d;
220 
221 /*MC
222       DMDACoor3d - Structure for holding 3d (x, y and z) coordinates.
223 
224     Level: intermediate
225 
226     Sample Usage:
227       DMDACoor3d ***coors;
228       Vec      vcoors;
229       DM       cda;
230 
231       DMDAGetCoordinates(da,&vcoors);
232       DMDAGetCoordinateDA(da,&cda);
233       DMDAVecGetArray(cda,vcoors,&coors);
234       DMDAGetCorners(cda,&mstart,&nstart,&pstart,&m,&n,&p)
235       for (i=mstart; i<mstart+m; i++) {
236         for (j=nstart; j<nstart+n; j++) {
237           for (k=pstart; k<pstart+p; k++) {
238             x = coors[k][j][i].x;
239             y = coors[k][j][i].y;
240             z = coors[k][j][i].z;
241           ......
242         }
243       }
244       DMDAVecRestoreArray(dac,vcoors,&coors);
245 
246 .seealso: DMDACoor2d, DMDAForEachPointBegin(), DMDAGetCoordinateDA(), DMDAGetCoordinates(), DMDAGetGhostCoordinates()
247 M*/
248 typedef struct {PetscScalar x,y,z;} DMDACoor3d;
249 
250 extern PetscErrorCode   DMDAGetLocalInfo(DM,DMDALocalInfo*);
251 typedef PetscErrorCode (*DMDALocalFunction1)(DMDALocalInfo*,void*,void*,void*);
252 extern PetscErrorCode   DMDAFormFunctionLocal(DM, DMDALocalFunction1, Vec, Vec, void *);
253 extern PetscErrorCode   DMDAFormFunctionLocalGhost(DM, DMDALocalFunction1, Vec, Vec, void *);
254 extern PetscErrorCode   DMDAFormJacobianLocal(DM, DMDALocalFunction1, Vec, Mat, void *);
255 extern PetscErrorCode   DMDAFormFunction1(DM,Vec,Vec,void*);
256 extern PetscErrorCode   DMDAFormFunction(DM,PetscErrorCode (*)(void),Vec,Vec,void*);
257 extern PetscErrorCode   DMDAFormFunctioni1(DM,PetscInt,Vec,PetscScalar*,void*);
258 extern PetscErrorCode   DMDAFormFunctionib1(DM,PetscInt,Vec,PetscScalar*,void*);
259 extern PetscErrorCode   DMDAComputeJacobian1WithAdic(DM,Vec,Mat,void*);
260 extern PetscErrorCode   DMDAComputeJacobian1WithAdifor(DM,Vec,Mat,void*);
261 extern PetscErrorCode   DMDAMultiplyByJacobian1WithAdic(DM,Vec,Vec,Vec,void*);
262 extern PetscErrorCode   DMDAMultiplyByJacobian1WithAdifor(DM,Vec,Vec,Vec,void*);
263 extern PetscErrorCode   DMDAMultiplyByJacobian1WithAD(DM,Vec,Vec,Vec,void*);
264 extern PetscErrorCode   DMDAComputeJacobian1(DM,Vec,Mat,void*);
265 extern PetscErrorCode   DMDAGetLocalFunction(DM,DMDALocalFunction1*);
266 extern PetscErrorCode   DMDASetLocalFunction(DM,DMDALocalFunction1);
267 extern PetscErrorCode   DMDASetLocalFunctioni(DM,PetscErrorCode (*)(DMDALocalInfo*,MatStencil*,void*,PetscScalar*,void*));
268 extern PetscErrorCode   DMDASetLocalFunctionib(DM,PetscErrorCode (*)(DMDALocalInfo*,MatStencil*,void*,PetscScalar*,void*));
269 extern PetscErrorCode   DMDAGetLocalJacobian(DM,DMDALocalFunction1*);
270 extern PetscErrorCode   DMDASetLocalJacobian(DM,DMDALocalFunction1);
271 extern PetscErrorCode   DMDASetLocalAdicFunction_Private(DM,DMDALocalFunction1);
272 
273 extern PetscErrorCode MatSetDA(Mat,DM);
274 extern PetscErrorCode MatRegisterDAAD(void);
275 extern PetscErrorCode MatCreateDAAD(DM,Mat*);
276 extern PetscErrorCode MatCreateSeqUSFFT(Vec,DM,Mat*);
277 
278 /*MC
279        DMDASetLocalAdicFunction - Caches in a DM a local function computed by ADIC/ADIFOR
280 
281    Synopsis:
282    PetscErrorCode DMDASetLocalAdicFunction(DM da,DMDALocalFunction1 ad_lf)
283 
284    Logically Collective on DM
285 
286    Input Parameter:
287 +  da - initial distributed array
288 -  ad_lf - the local function as computed by ADIC/ADIFOR
289 
290    Level: intermediate
291 
292 .keywords:  distributed array, refine
293 
294 .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDAGetLocalFunction(), DMDASetLocalFunction(),
295           DMDASetLocalJacobian()
296 M*/
297 #if defined(PETSC_HAVE_ADIC)
298 #  define DMDASetLocalAdicFunction(a,d) DMDASetLocalAdicFunction_Private(a,(DMDALocalFunction1)d)
299 #else
300 #  define DMDASetLocalAdicFunction(a,d) DMDASetLocalAdicFunction_Private(a,0)
301 #endif
302 
303 extern PetscErrorCode   DMDASetLocalAdicMFFunction_Private(DM,DMDALocalFunction1);
304 #if defined(PETSC_HAVE_ADIC)
305 #  define DMDASetLocalAdicMFFunction(a,d) DMDASetLocalAdicMFFunction_Private(a,(DMDALocalFunction1)d)
306 #else
307 #  define DMDASetLocalAdicMFFunction(a,d) DMDASetLocalAdicMFFunction_Private(a,0)
308 #endif
309 extern PetscErrorCode   DMDASetLocalAdicFunctioni_Private(DM,PetscErrorCode (*)(DMDALocalInfo*,MatStencil*,void*,void*,void*));
310 #if defined(PETSC_HAVE_ADIC)
311 #  define DMDASetLocalAdicFunctioni(a,d) DMDASetLocalAdicFunctioni_Private(a,(PetscErrorCode (*)(DMDALocalInfo*,MatStencil*,void*,void*,void*))d)
312 #else
313 #  define DMDASetLocalAdicFunctioni(a,d) DMDASetLocalAdicFunctioni_Private(a,0)
314 #endif
315 extern PetscErrorCode   DMDASetLocalAdicMFFunctioni_Private(DM,PetscErrorCode (*)(DMDALocalInfo*,MatStencil*,void*,void*,void*));
316 #if defined(PETSC_HAVE_ADIC)
317 #  define DMDASetLocalAdicMFFunctioni(a,d) DMDASetLocalAdicMFFunctioni_Private(a,(PetscErrorCode (*)(DMDALocalInfo*,MatStencil*,void*,void*,void*))d)
318 #else
319 #  define DMDASetLocalAdicMFFunctioni(a,d) DMDASetLocalAdicMFFunctioni_Private(a,0)
320 #endif
321 
322 extern PetscErrorCode   DMDASetLocalAdicFunctionib_Private(DM,PetscErrorCode (*)(DMDALocalInfo*,MatStencil*,void*,void*,void*));
323 #if defined(PETSC_HAVE_ADIC)
324 #  define DMDASetLocalAdicFunctionib(a,d) DMDASetLocalAdicFunctionib_Private(a,(PetscErrorCode (*)(DMDALocalInfo*,MatStencil*,void*,void*,void*))d)
325 #else
326 #  define DMDASetLocalAdicFunctionib(a,d) DMDASetLocalAdicFunctionib_Private(a,0)
327 #endif
328 extern PetscErrorCode   DMDASetLocalAdicMFFunctionib_Private(DM,PetscErrorCode (*)(DMDALocalInfo*,MatStencil*,void*,void*,void*));
329 #if defined(PETSC_HAVE_ADIC)
330 #  define DMDASetLocalAdicMFFunctionib(a,d) DMDASetLocalAdicMFFunctionib_Private(a,(PetscErrorCode (*)(DMDALocalInfo*,MatStencil*,void*,void*,void*))d)
331 #else
332 #  define DMDASetLocalAdicMFFunctionib(a,d) DMDASetLocalAdicMFFunctionib_Private(a,0)
333 #endif
334 
335 extern PetscErrorCode   DMDAFormFunctioniTest1(DM,void*);
336 extern PetscErrorCode   DMDASetGetMatrix(DM,PetscErrorCode (*)(DM, const MatType,Mat *));
337 extern PetscErrorCode   DMDASetBlockFills(DM,PetscInt*,PetscInt*);
338 extern PetscErrorCode   DMDASetMatPreallocateOnly(DM,PetscBool );
339 extern PetscErrorCode   DMDASetRefinementFactor(DM,PetscInt,PetscInt,PetscInt);
340 extern PetscErrorCode   DMDAGetRefinementFactor(DM,PetscInt*,PetscInt*,PetscInt*);
341 
342 extern PetscErrorCode   DMDAGetAdicArray(DM,PetscBool ,void*,void*,PetscInt*);
343 extern PetscErrorCode   DMDARestoreAdicArray(DM,PetscBool ,void*,void*,PetscInt*);
344 extern PetscErrorCode   DMDAGetAdicMFArray(DM,PetscBool ,void*,void*,PetscInt*);
345 extern PetscErrorCode   DMDAGetAdicMFArray4(DM,PetscBool ,void*,void*,PetscInt*);
346 extern PetscErrorCode   DMDAGetAdicMFArray9(DM,PetscBool ,void*,void*,PetscInt*);
347 extern PetscErrorCode   DMDAGetAdicMFArrayb(DM,PetscBool ,void*,void*,PetscInt*);
348 extern PetscErrorCode   DMDARestoreAdicMFArray(DM,PetscBool ,void*,void*,PetscInt*);
349 extern PetscErrorCode   DMDAGetArray(DM,PetscBool ,void*);
350 extern PetscErrorCode   DMDARestoreArray(DM,PetscBool ,void*);
351 extern PetscErrorCode   ad_DAGetArray(DM,PetscBool ,void*);
352 extern PetscErrorCode   ad_DARestoreArray(DM,PetscBool ,void*);
353 extern PetscErrorCode   admf_DAGetArray(DM,PetscBool ,void*);
354 extern PetscErrorCode   admf_DARestoreArray(DM,PetscBool ,void*);
355 
356 extern PetscErrorCode   DMDACreatePF(DM,PF*);
357 
358 #endif
359