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