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