xref: /petsc/include/petscdmda.h (revision 3923b477fd0dced8a2d147b4fb4519fe3af97d3f)
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 physical 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, you can put values into them and then apply a stencil that uses those ghost locations),
43    DMDA_BOUNDARY_MIRROR (not yet implemented), or DMDA_BOUNDARY_PERIODIC
44    (ghost nodes filled by the opposite edge of the domain).
45 
46    Note: This is information for the boundary of the __PHYSICAL__ domain. It has nothing to do with boundaries between
47      processes, that width is always determined by the stencil width, see DMDASetStencilWidth().
48 
49 .seealso: DMDASetBoundaryType(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDACreate()
50 E*/
51 typedef enum { DMDA_BOUNDARY_NONE, DMDA_BOUNDARY_GHOSTED, DMDA_BOUNDARY_MIRROR, DMDA_BOUNDARY_PERIODIC } DMDABoundaryType;
52 
53 PETSC_EXTERN const char *const DMDABoundaryTypes[];
54 
55 /*E
56     DMDAInterpolationType - Defines the type of interpolation that will be returned by
57        DMCreateInterpolation.
58 
59    Level: beginner
60 
61 .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMCreateInterpolation(), DMDASetInterpolationType(), DMDACreate()
62 E*/
63 typedef enum { DMDA_Q0, DMDA_Q1 } DMDAInterpolationType;
64 
65 PETSC_EXTERN PetscErrorCode DMDASetInterpolationType(DM,DMDAInterpolationType);
66 PETSC_EXTERN PetscErrorCode DMDAGetInterpolationType(DM,DMDAInterpolationType*);
67 
68 /*E
69     DMDAElementType - Defines the type of elements that will be returned by
70        DMDAGetElements()
71 
72    Level: beginner
73 
74 .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMCreateInterpolation(), DMDASetInterpolationType(),
75           DMDASetElementType(), DMDAGetElements(), DMDARestoreElements(), DMDACreate()
76 E*/
77 typedef enum { DMDA_ELEMENT_P1, DMDA_ELEMENT_Q1 } DMDAElementType;
78 
79 PETSC_EXTERN PetscErrorCode DMDASetElementType(DM,DMDAElementType);
80 PETSC_EXTERN PetscErrorCode DMDAGetElementType(DM,DMDAElementType*);
81 PETSC_EXTERN PetscErrorCode DMDAGetElements(DM,PetscInt *,PetscInt *,const PetscInt*[]);
82 PETSC_EXTERN PetscErrorCode DMDARestoreElements(DM,PetscInt *,PetscInt *,const PetscInt*[]);
83 
84 typedef enum { DMDA_X,DMDA_Y,DMDA_Z } DMDADirection;
85 
86 #define MATSEQUSFFT        "sequsfft"
87 
88 PETSC_EXTERN PetscErrorCode DMDACreate(MPI_Comm,DM*);
89 PETSC_EXTERN PetscErrorCode DMDASetDim(DM,PetscInt);
90 PETSC_EXTERN PetscErrorCode DMDASetSizes(DM,PetscInt,PetscInt,PetscInt);
91 PETSC_EXTERN PetscErrorCode DMDACreate1d(MPI_Comm,DMDABoundaryType,PetscInt,PetscInt,PetscInt,const PetscInt[],DM *);
92 PETSC_EXTERN PetscErrorCode DMDACreate2d(MPI_Comm,DMDABoundaryType,DMDABoundaryType,DMDAStencilType,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],DM*);
93 PETSC_EXTERN PetscErrorCode DMDACreate3d(MPI_Comm,DMDABoundaryType,DMDABoundaryType,DMDABoundaryType,DMDAStencilType,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscInt[],DM*);
94 
95 PETSC_EXTERN PetscErrorCode DMDAGlobalToNaturalBegin(DM,Vec,InsertMode,Vec);
96 PETSC_EXTERN PetscErrorCode DMDAGlobalToNaturalEnd(DM,Vec,InsertMode,Vec);
97 PETSC_EXTERN PetscErrorCode DMDANaturalToGlobalBegin(DM,Vec,InsertMode,Vec);
98 PETSC_EXTERN PetscErrorCode DMDANaturalToGlobalEnd(DM,Vec,InsertMode,Vec);
99 PETSC_EXTERN PetscErrorCode DMDALocalToLocalBegin(DM,Vec,InsertMode,Vec);
100 PETSC_EXTERN PetscErrorCode DMDALocalToLocalEnd(DM,Vec,InsertMode,Vec);
101 PETSC_EXTERN PetscErrorCode DMDACreateNaturalVector(DM,Vec *);
102 
103 PETSC_EXTERN PetscErrorCode DMDAGetCorners(DM,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
104 PETSC_EXTERN PetscErrorCode DMDAGetGhostCorners(DM,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
105 PETSC_EXTERN PetscErrorCode DMDAGetInfo(DM,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,DMDABoundaryType*,DMDABoundaryType*,DMDABoundaryType*,DMDAStencilType*);
106 PETSC_EXTERN PetscErrorCode DMDAGetProcessorSubset(DM,DMDADirection,PetscInt,MPI_Comm*);
107 PETSC_EXTERN PetscErrorCode DMDAGetProcessorSubsets(DM,DMDADirection,MPI_Comm*);
108 
109 PETSC_EXTERN PetscErrorCode DMDAGlobalToNaturalAllCreate(DM,VecScatter*);
110 PETSC_EXTERN PetscErrorCode DMDANaturalAllToGlobalCreate(DM,VecScatter*);
111 
112 PETSC_EXTERN PetscErrorCode DMDAGetGlobalIndices(DM,PetscInt*,PetscInt**);
113 
114 PETSC_EXTERN PetscErrorCode DMDAGetScatter(DM,VecScatter*,VecScatter*,VecScatter*);
115 PETSC_EXTERN PetscErrorCode DMDAGetNeighbors(DM,const PetscMPIInt**);
116 
117 PETSC_EXTERN PetscErrorCode DMDAGetAO(DM,AO*);
118 PETSC_EXTERN PetscErrorCode DMDASetUniformCoordinates(DM,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal);
119 PETSC_EXTERN PetscErrorCode DMDAGetBoundingBox(DM,PetscReal[],PetscReal[]);
120 PETSC_EXTERN PetscErrorCode DMDAGetLocalBoundingBox(DM,PetscReal[],PetscReal[]);
121 /* function to wrap coordinates around boundary */
122 PETSC_EXTERN PetscErrorCode DMDAMapCoordsToPeriodicDomain(DM,PetscScalar*,PetscScalar*);
123 
124 PETSC_EXTERN PetscErrorCode DMDAGetReducedDA(DM,PetscInt,DM*);
125 
126 PETSC_EXTERN PetscErrorCode DMDASetFieldName(DM,PetscInt,const char[]);
127 PETSC_EXTERN PetscErrorCode DMDAGetFieldName(DM,PetscInt,const char**);
128 
129 PETSC_EXTERN PetscErrorCode DMDASetBoundaryType(DM,DMDABoundaryType,DMDABoundaryType,DMDABoundaryType);
130 PETSC_EXTERN PetscErrorCode DMDASetDof(DM, PetscInt);
131 PETSC_EXTERN PetscErrorCode DMDASetOverlap(DM, PetscInt);
132 PETSC_EXTERN PetscErrorCode DMDASetStencilWidth(DM, PetscInt);
133 PETSC_EXTERN PetscErrorCode DMDASetOwnershipRanges(DM,const PetscInt[],const PetscInt[],const PetscInt[]);
134 PETSC_EXTERN PetscErrorCode DMDAGetOwnershipRanges(DM,const PetscInt**,const PetscInt**,const PetscInt**);
135 PETSC_EXTERN PetscErrorCode DMDASetNumProcs(DM, PetscInt, PetscInt, PetscInt);
136 PETSC_EXTERN PetscErrorCode DMDASetStencilType(DM, DMDAStencilType);
137 
138 PETSC_EXTERN PetscErrorCode DMDAVecGetArray(DM,Vec,void *);
139 PETSC_EXTERN PetscErrorCode DMDAVecRestoreArray(DM,Vec,void *);
140 
141 PETSC_EXTERN PetscErrorCode DMDAVecGetArrayDOF(DM,Vec,void *);
142 PETSC_EXTERN PetscErrorCode DMDAVecRestoreArrayDOF(DM,Vec,void *);
143 
144 PETSC_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 point 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       DMGetCoordinates(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(), DMGetCoordinates(), 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       DMGetCoordinates(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(), DMGetCoordinates(), DMDAGetGhostCoordinates()
255 M*/
256 typedef struct {PetscScalar x,y,z;} DMDACoor3d;
257 
258 PETSC_EXTERN PetscErrorCode DMDAGetLocalInfo(DM,DMDALocalInfo*);
259 PETSC_EXTERN PetscErrorCode DMDAGetLocalBlockInfo(DM,DMDALocalInfo*);
260 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDALocalFunction1)(DMDALocalInfo*,void*,void*,void*);
261 PETSC_EXTERN PetscErrorCode DMDAComputeFunctionLocal(DM, DMDALocalFunction1, Vec, Vec, void *);
262 PETSC_EXTERN PetscErrorCode DMDAComputeFunctionLocalGhost(DM, DMDALocalFunction1, Vec, Vec, void *);
263 PETSC_EXTERN PetscErrorCode DMDAFormJacobianLocal(DM, DMDALocalFunction1, Vec, Mat, void *);
264 PETSC_EXTERN PetscErrorCode DMDAComputeFunction1(DM,Vec,Vec,void*);
265 PETSC_EXTERN PetscErrorCode DMDAComputeFunction(DM,PetscErrorCode (*)(void),Vec,Vec,void*);
266 PETSC_EXTERN PetscErrorCode DMDAComputeFunctioni1(DM,PetscInt,Vec,PetscScalar*,void*);
267 PETSC_EXTERN PetscErrorCode DMDAComputeFunctionib1(DM,PetscInt,Vec,PetscScalar*,void*);
268 PETSC_EXTERN PetscErrorCode DMDAComputeJacobian1WithAdic(DM,Vec,Mat,void*);
269 PETSC_EXTERN PetscErrorCode DMDAComputeJacobian1WithAdifor(DM,Vec,Mat,void*);
270 PETSC_EXTERN PetscErrorCode DMDAMultiplyByJacobian1WithAdic(DM,Vec,Vec,Vec,void*);
271 PETSC_EXTERN PetscErrorCode DMDAMultiplyByJacobian1WithAdifor(DM,Vec,Vec,Vec,void*);
272 PETSC_EXTERN PetscErrorCode DMDAMultiplyByJacobian1WithAD(DM,Vec,Vec,Vec,void*);
273 PETSC_EXTERN PetscErrorCode DMDAComputeJacobian1(DM,Vec,Mat,void*);
274 PETSC_EXTERN PetscErrorCode DMDAGetLocalFunction(DM,DMDALocalFunction1*);
275 PETSC_EXTERN PetscErrorCode DMDASetLocalFunction(DM,DMDALocalFunction1);
276 PETSC_EXTERN PetscErrorCode DMDASetLocalFunctioni(DM,PetscErrorCode (*)(DMDALocalInfo*,MatStencil*,void*,PetscScalar*,void*));
277 PETSC_EXTERN PetscErrorCode DMDASetLocalFunctionib(DM,PetscErrorCode (*)(DMDALocalInfo*,MatStencil*,void*,PetscScalar*,void*));
278 PETSC_EXTERN PetscErrorCode DMDAGetLocalJacobian(DM,DMDALocalFunction1*);
279 PETSC_EXTERN PetscErrorCode DMDASetLocalJacobian(DM,DMDALocalFunction1);
280 PETSC_EXTERN PetscErrorCode DMDASetLocalAdicFunction_Private(DM,DMDALocalFunction1);
281 
282 PETSC_EXTERN PetscErrorCode MatSetDM(Mat,DM);
283 PETSC_EXTERN PetscErrorCode MatRegisterDAAD(void);
284 PETSC_EXTERN PetscErrorCode MatCreateDAAD(DM,Mat*);
285 PETSC_EXTERN PetscErrorCode MatCreateSeqUSFFT(Vec,DM,Mat*);
286 
287 /*MC
288        DMDASetLocalAdicFunction - Caches in a DM a local function computed by ADIC/ADIFOR
289 
290    Synopsis:
291    PetscErrorCode DMDASetLocalAdicFunction(DM da,DMDALocalFunction1 ad_lf)
292 
293    Logically Collective on DM
294 
295    Input Parameter:
296 +  da - initial distributed array
297 -  ad_lf - the local function as computed by ADIC/ADIFOR
298 
299    Level: intermediate
300 
301 .keywords:  distributed array, refine
302 
303 .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDAGetLocalFunction(), DMDASetLocalFunction(),
304           DMDASetLocalJacobian()
305 M*/
306 #if defined(PETSC_HAVE_ADIC)
307 #  define DMDASetLocalAdicFunction(a,d) DMDASetLocalAdicFunction_Private(a,(DMDALocalFunction1)d)
308 #else
309 #  define DMDASetLocalAdicFunction(a,d) DMDASetLocalAdicFunction_Private(a,0)
310 #endif
311 
312 PETSC_EXTERN PetscErrorCode DMDASetLocalAdicMFFunction_Private(DM,DMDALocalFunction1);
313 #if defined(PETSC_HAVE_ADIC)
314 #  define DMDASetLocalAdicMFFunction(a,d) DMDASetLocalAdicMFFunction_Private(a,(DMDALocalFunction1)d)
315 #else
316 #  define DMDASetLocalAdicMFFunction(a,d) DMDASetLocalAdicMFFunction_Private(a,0)
317 #endif
318 PETSC_EXTERN PetscErrorCode DMDASetLocalAdicFunctioni_Private(DM,PetscErrorCode (*)(DMDALocalInfo*,MatStencil*,void*,void*,void*));
319 #if defined(PETSC_HAVE_ADIC)
320 #  define DMDASetLocalAdicFunctioni(a,d) DMDASetLocalAdicFunctioni_Private(a,(PetscErrorCode (*)(DMDALocalInfo*,MatStencil*,void*,void*,void*))d)
321 #else
322 #  define DMDASetLocalAdicFunctioni(a,d) DMDASetLocalAdicFunctioni_Private(a,0)
323 #endif
324 PETSC_EXTERN PetscErrorCode DMDASetLocalAdicMFFunctioni_Private(DM,PetscErrorCode (*)(DMDALocalInfo*,MatStencil*,void*,void*,void*));
325 #if defined(PETSC_HAVE_ADIC)
326 #  define DMDASetLocalAdicMFFunctioni(a,d) DMDASetLocalAdicMFFunctioni_Private(a,(PetscErrorCode (*)(DMDALocalInfo*,MatStencil*,void*,void*,void*))d)
327 #else
328 #  define DMDASetLocalAdicMFFunctioni(a,d) DMDASetLocalAdicMFFunctioni_Private(a,0)
329 #endif
330 
331 PETSC_EXTERN PetscErrorCode DMDASetLocalAdicFunctionib_Private(DM,PetscErrorCode (*)(DMDALocalInfo*,MatStencil*,void*,void*,void*));
332 #if defined(PETSC_HAVE_ADIC)
333 #  define DMDASetLocalAdicFunctionib(a,d) DMDASetLocalAdicFunctionib_Private(a,(PetscErrorCode (*)(DMDALocalInfo*,MatStencil*,void*,void*,void*))d)
334 #else
335 #  define DMDASetLocalAdicFunctionib(a,d) DMDASetLocalAdicFunctionib_Private(a,0)
336 #endif
337 PETSC_EXTERN PetscErrorCode DMDASetLocalAdicMFFunctionib_Private(DM,PetscErrorCode (*)(DMDALocalInfo*,MatStencil*,void*,void*,void*));
338 #if defined(PETSC_HAVE_ADIC)
339 #  define DMDASetLocalAdicMFFunctionib(a,d) DMDASetLocalAdicMFFunctionib_Private(a,(PetscErrorCode (*)(DMDALocalInfo*,MatStencil*,void*,void*,void*))d)
340 #else
341 #  define DMDASetLocalAdicMFFunctionib(a,d) DMDASetLocalAdicMFFunctionib_Private(a,0)
342 #endif
343 
344 PETSC_EXTERN PetscErrorCode DMDAComputeFunctioniTest1(DM,void*);
345 PETSC_EXTERN PetscErrorCode DMDASetGetMatrix(DM,PetscErrorCode (*)(DM, MatType,Mat *));
346 PETSC_EXTERN PetscErrorCode DMDASetBlockFills(DM,PetscInt*,PetscInt*);
347 PETSC_EXTERN PetscErrorCode DMDASetRefinementFactor(DM,PetscInt,PetscInt,PetscInt);
348 PETSC_EXTERN PetscErrorCode DMDAGetRefinementFactor(DM,PetscInt*,PetscInt*,PetscInt*);
349 
350 PETSC_EXTERN PetscErrorCode DMDAGetAdicArray(DM,PetscBool ,void*,void*,PetscInt*);
351 PETSC_EXTERN PetscErrorCode DMDARestoreAdicArray(DM,PetscBool ,void*,void*,PetscInt*);
352 PETSC_EXTERN PetscErrorCode DMDAGetAdicMFArray(DM,PetscBool ,void*,void*,PetscInt*);
353 PETSC_EXTERN PetscErrorCode DMDAGetAdicMFArray4(DM,PetscBool ,void*,void*,PetscInt*);
354 PETSC_EXTERN PetscErrorCode DMDAGetAdicMFArray9(DM,PetscBool ,void*,void*,PetscInt*);
355 PETSC_EXTERN PetscErrorCode DMDAGetAdicMFArrayb(DM,PetscBool ,void*,void*,PetscInt*);
356 PETSC_EXTERN PetscErrorCode DMDARestoreAdicMFArray(DM,PetscBool ,void*,void*,PetscInt*);
357 PETSC_EXTERN PetscErrorCode DMDAGetArray(DM,PetscBool ,void*);
358 PETSC_EXTERN PetscErrorCode DMDARestoreArray(DM,PetscBool ,void*);
359 PETSC_EXTERN PetscErrorCode ad_DAGetArray(DM,PetscBool ,void*);
360 PETSC_EXTERN PetscErrorCode ad_DARestoreArray(DM,PetscBool ,void*);
361 PETSC_EXTERN PetscErrorCode admf_DAGetArray(DM,PetscBool ,void*);
362 PETSC_EXTERN PetscErrorCode admf_DARestoreArray(DM,PetscBool ,void*);
363 
364 PETSC_EXTERN PetscErrorCode DMDACreatePF(DM,PF*);
365 
366 PETSC_EXTERN PetscErrorCode DMDAGetNumCells(DM, PetscInt *);
367 PETSC_EXTERN PetscErrorCode DMDAGetNumVertices(DM, PetscInt *, PetscInt *, PetscInt *, PetscInt *);
368 PETSC_EXTERN PetscErrorCode DMDAGetNumFaces(DM, PetscInt *, PetscInt *, PetscInt *, PetscInt *, PetscInt *, PetscInt *);
369 PETSC_EXTERN PetscErrorCode DMDAGetHeightStratum(DM, PetscInt, PetscInt *, PetscInt *);
370 PETSC_EXTERN PetscErrorCode DMDACreateSection(DM, PetscInt[], PetscInt[], PetscInt[], PetscInt[]);
371 PETSC_EXTERN PetscErrorCode DMDAComputeCellGeometry(DM, PetscInt, PetscQuadrature *, PetscReal [], PetscReal [], PetscReal [], PetscReal []);
372 PETSC_EXTERN PetscErrorCode DMDAVecGetClosure(DM, PetscSection, Vec, PetscInt, const PetscScalar **);
373 PETSC_EXTERN PetscErrorCode DMDAVecSetClosure(DM, PetscSection, Vec, PetscInt, const PetscScalar *, InsertMode);
374 PETSC_EXTERN PetscErrorCode DMDAGetClosure(DM,PetscSection,PetscInt,PetscInt*,const PetscInt**);
375 PETSC_EXTERN PetscErrorCode DMDARestoreClosure(DM,PetscSection,PetscInt,PetscInt*,const PetscInt**);
376 PETSC_EXTERN PetscErrorCode DMDAGetClosureScalar(DM,PetscSection,PetscInt,PetscScalar*,const PetscScalar**);
377 PETSC_EXTERN PetscErrorCode DMDARestoreClosureScalar(DM,PetscSection,PetscInt,PetscScalar*,const PetscScalar**);
378 PETSC_EXTERN PetscErrorCode DMDASetClosureScalar(DM,PetscSection,PetscInt,PetscScalar*,const PetscScalar*,InsertMode);
379 PETSC_EXTERN PetscErrorCode DMDAConvertToCell(DM, MatStencil, PetscInt *);
380 
381 #define DMDA_FILE_CLASSID 1211220
382 #endif
383