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