xref: /petsc/include/petscdmda.h (revision 211d3afb70a2c85c8ea35e22f0c64a5e0dd01503)
13c48a1e8SJed Brown #if !defined(__PETSCDMDA_H)
23c48a1e8SJed Brown #define __PETSCDMDA_H
33c48a1e8SJed Brown 
42c8e378dSBarry Smith #include <petscdm.h>
51e25c274SJed Brown #include <petscdmdatypes.h>
62c8e378dSBarry Smith #include <petscpf.h>
72c8e378dSBarry Smith #include <petscao.h>
8bb619f44SMatthew G. Knepley #include <petscfe.h>
93c48a1e8SJed Brown 
103c48a1e8SJed Brown /*MC
113c48a1e8SJed Brown      DMDA_STENCIL_STAR - "Star"-type stencil. In logical grid coordinates, only (i,j,k), (i+s,j,k), (i,j+s,k),
123c48a1e8SJed Brown                        (i,j,k+s) are in the stencil  NOT, for example, (i+s,j+s,k)
133c48a1e8SJed Brown 
143c48a1e8SJed Brown      Level: beginner
153c48a1e8SJed Brown 
161957e957SBarry Smith      Determines what ghost point values are brought over to each process; in this case the "corner" values are not
171957e957SBarry Smith      brought over and hence should not be accessed locally
181957e957SBarry Smith 
19db05f41bSBarry Smith .seealso: DMDA_STENCIL_BOX, DMDAStencilType, DMDASetStencilType()
203c48a1e8SJed Brown M*/
213c48a1e8SJed Brown 
223c48a1e8SJed Brown /*MC
233c48a1e8SJed Brown      DMDA_STENCIL_BOX - "Box"-type stencil. In logical grid coordinates, any of (i,j,k), (i+s,j+r,k+t) may
243c48a1e8SJed Brown                       be in the stencil.
253c48a1e8SJed Brown 
263c48a1e8SJed Brown      Level: beginner
273c48a1e8SJed Brown 
28db05f41bSBarry Smith .seealso: DMDA_STENCIL_STAR, DMDAStencilType, DMDASetStencilType()
293c48a1e8SJed Brown M*/
303c48a1e8SJed Brown 
31014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDASetInterpolationType(DM,DMDAInterpolationType);
32014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGetInterpolationType(DM,DMDAInterpolationType*);
333c48a1e8SJed Brown 
34d4a6ed37SStefano Zampini /* FEM */
35014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDASetElementType(DM,DMDAElementType);
36014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGetElementType(DM,DMDAElementType*);
37014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGetElements(DM,PetscInt*,PetscInt*,const PetscInt*[]);
38014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDARestoreElements(DM,PetscInt*,PetscInt*,const PetscInt*[]);
39d4a6ed37SStefano Zampini PETSC_EXTERN PetscErrorCode DMDAGetElementsSizes(DM,PetscInt*,PetscInt*,PetscInt*);
40d4a6ed37SStefano Zampini PETSC_EXTERN PetscErrorCode DMDAGetElementsCorners(DM,PetscInt*,PetscInt*,PetscInt*);
41d4a6ed37SStefano Zampini PETSC_EXTERN PetscErrorCode DMDAGetSubdomainCornersIS(DM,IS*);
42d4a6ed37SStefano Zampini PETSC_EXTERN PetscErrorCode DMDARestoreSubdomainCornersIS(DM,IS*);
433c48a1e8SJed Brown 
443c48a1e8SJed Brown typedef enum { DMDA_X,DMDA_Y,DMDA_Z } DMDADirection;
453c48a1e8SJed Brown 
463c48a1e8SJed Brown #define MATSEQUSFFT        "sequsfft"
473c48a1e8SJed Brown 
48014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDACreate(MPI_Comm,DM*);
49014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDASetSizes(DM,PetscInt,PetscInt,PetscInt);
50bff4a2f0SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMDACreate1d(MPI_Comm,DMBoundaryType,PetscInt,PetscInt,PetscInt,const PetscInt[],DM *);
51bff4a2f0SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMDACreate2d(MPI_Comm,DMBoundaryType,DMBoundaryType,DMDAStencilType,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],DM*);
52bff4a2f0SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMDACreate3d(MPI_Comm,DMBoundaryType,DMBoundaryType,DMBoundaryType,DMDAStencilType,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscInt[],DM*);
533c48a1e8SJed Brown 
54014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGlobalToNaturalBegin(DM,Vec,InsertMode,Vec);
55014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGlobalToNaturalEnd(DM,Vec,InsertMode,Vec);
56014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDANaturalToGlobalBegin(DM,Vec,InsertMode,Vec);
57014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDANaturalToGlobalEnd(DM,Vec,InsertMode,Vec);
58d78e899eSRichard Tran Mills PETSC_DEPRECATED("Use DMLocalToLocalBegin()") PETSC_STATIC_INLINE PetscErrorCode DMDALocalToLocalBegin(DM dm,Vec g,InsertMode mode,Vec l) {return DMLocalToLocalBegin(dm,g,mode,l);}
59d78e899eSRichard Tran Mills PETSC_DEPRECATED("Use DMLocalToLocalEnd()") PETSC_STATIC_INLINE PetscErrorCode DMDALocalToLocalEnd(DM dm,Vec g,InsertMode mode,Vec l) {return DMLocalToLocalEnd(dm,g,mode,l);}
60014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDACreateNaturalVector(DM,Vec *);
613c48a1e8SJed Brown 
62014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGetCorners(DM,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
63014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGetGhostCorners(DM,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
64bff4a2f0SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMDAGetInfo(DM,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,DMBoundaryType*,DMBoundaryType*,DMBoundaryType*,DMDAStencilType*);
65014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGetProcessorSubset(DM,DMDADirection,PetscInt,MPI_Comm*);
66014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGetProcessorSubsets(DM,DMDADirection,MPI_Comm*);
67db05f41bSBarry Smith PETSC_EXTERN PetscErrorCode DMDAGetRay(DM,DMDADirection,PetscInt,Vec*,VecScatter*);
683c48a1e8SJed Brown 
69014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGlobalToNaturalAllCreate(DM,VecScatter*);
70014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDANaturalAllToGlobalCreate(DM,VecScatter*);
713c48a1e8SJed Brown 
72bd1fc5aeSBarry Smith PETSC_EXTERN PetscErrorCode DMDAGetScatter(DM,VecScatter*,VecScatter*);
73014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGetNeighbors(DM,const PetscMPIInt**);
743c48a1e8SJed Brown 
759db3d8bcSStefano Zampini PETSC_EXTERN PetscErrorCode DMDASetAOType(DM,AOType);
76014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGetAO(DM,AO*);
77014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDASetUniformCoordinates(DM,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal);
7862197512SBarry Smith #include <petscgll.h>
7962197512SBarry Smith PETSC_EXTERN PetscErrorCode DMDASetGLLCoordinates(DM,PetscGLL*);
80c593f006SBarry Smith PETSC_EXTERN PetscErrorCode DMDAGetCoordinateArray(DM,void*);
81c593f006SBarry Smith PETSC_EXTERN PetscErrorCode DMDARestoreCoordinateArray(DM,void*);
82014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGetBoundingBox(DM,PetscReal[],PetscReal[]);
83014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGetLocalBoundingBox(DM,PetscReal[],PetscReal[]);
84a74ba6f7SBarry Smith PETSC_EXTERN PetscErrorCode DMDAGetLogicalCoordinate(DM,PetscScalar,PetscScalar,PetscScalar,PetscInt*,PetscInt*,PetscInt*,PetscScalar*,PetscScalar*,PetscScalar*);
85f0aa4865SJed Brown /* function to wrap coordinates around boundary */
86014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAMapCoordsToPeriodicDomain(DM,PetscScalar*,PetscScalar*);
87f0aa4865SJed Brown 
88*211d3afbSPatrick Sanan PETSC_EXTERN PetscErrorCode DMDACreateCompatibleDMDA(DM,PetscInt,DM*);
893c48a1e8SJed Brown 
90014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDASetFieldName(DM,PetscInt,const char[]);
91014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGetFieldName(DM,PetscInt,const char**);
92c629b14aSBarry Smith PETSC_EXTERN PetscErrorCode DMDASetFieldNames(DM,const char * const *);
93c629b14aSBarry Smith PETSC_EXTERN PetscErrorCode DMDAGetFieldNames(DM,const char * const **);
94109c9344SBarry Smith PETSC_EXTERN PetscErrorCode DMDASetCoordinateName(DM,PetscInt,const char[]);
95109c9344SBarry Smith PETSC_EXTERN PetscErrorCode DMDAGetCoordinateName(DM,PetscInt,const char**);
963c48a1e8SJed Brown 
97bff4a2f0SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMDASetBoundaryType(DM,DMBoundaryType,DMBoundaryType,DMBoundaryType);
98014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDASetDof(DM, PetscInt);
99fb6725baSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMDAGetDof(DM, PetscInt*);
1007ddda789SPeter Brune PETSC_EXTERN PetscErrorCode DMDASetOverlap(DM,PetscInt,PetscInt,PetscInt);
1017ddda789SPeter Brune PETSC_EXTERN PetscErrorCode DMDAGetOverlap(DM,PetscInt*,PetscInt*,PetscInt*);
1023e7870d2SPeter Brune PETSC_EXTERN PetscErrorCode DMDASetNumLocalSubDomains(DM,PetscInt);
1033e7870d2SPeter Brune PETSC_EXTERN PetscErrorCode DMDAGetNumLocalSubDomains(DM,PetscInt*);
10495c13181SPeter Brune PETSC_EXTERN PetscErrorCode DMDAGetOffset(DM,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
10595c13181SPeter Brune PETSC_EXTERN PetscErrorCode DMDASetOffset(DM,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt);
10640234c92SPeter Brune PETSC_EXTERN PetscErrorCode DMDAGetNonOverlappingRegion(DM,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
10740234c92SPeter Brune PETSC_EXTERN PetscErrorCode DMDASetNonOverlappingRegion(DM,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt);
108014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDASetStencilWidth(DM, PetscInt);
109fb6725baSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMDAGetStencilWidth(DM, PetscInt*);
110014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDASetOwnershipRanges(DM,const PetscInt[],const PetscInt[],const PetscInt[]);
111014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGetOwnershipRanges(DM,const PetscInt**,const PetscInt**,const PetscInt**);
112014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDASetNumProcs(DM, PetscInt, PetscInt, PetscInt);
113014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDASetStencilType(DM, DMDAStencilType);
114fb6725baSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMDAGetStencilType(DM, DMDAStencilType*);
1153c48a1e8SJed Brown 
116014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAVecGetArray(DM,Vec,void *);
117014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAVecRestoreArray(DM,Vec,void *);
1183c48a1e8SJed Brown 
119014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAVecGetArrayDOF(DM,Vec,void *);
120014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAVecRestoreArrayDOF(DM,Vec,void *);
1213c48a1e8SJed Brown 
1225edff71fSBarry Smith PETSC_EXTERN PetscErrorCode DMDAVecGetArrayRead(DM,Vec,void *);
1235edff71fSBarry Smith PETSC_EXTERN PetscErrorCode DMDAVecRestoreArrayRead(DM,Vec,void *);
1245edff71fSBarry Smith 
1255edff71fSBarry Smith PETSC_EXTERN PetscErrorCode DMDAVecGetArrayDOFRead(DM,Vec,void *);
1265edff71fSBarry Smith PETSC_EXTERN PetscErrorCode DMDAVecRestoreArrayDOFRead(DM,Vec,void *);
1275edff71fSBarry Smith 
128014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDASplitComm2d(MPI_Comm,PetscInt,PetscInt,PetscInt,MPI_Comm*);
1293c48a1e8SJed Brown 
130ff9846d9SPeter Brune PETSC_EXTERN PetscErrorCode DMDACreatePatchIS(DM,MatStencil*,MatStencil*,IS*);
131ff9846d9SPeter Brune 
1323c48a1e8SJed Brown 
1333c48a1e8SJed Brown /*MC
1343c48a1e8SJed Brown       DMDACoor2d - Structure for holding 2d (x and y) coordinates.
1353c48a1e8SJed Brown 
1363c48a1e8SJed Brown     Level: intermediate
1373c48a1e8SJed Brown 
1388e897e54SSatish Balay     Synopsis:
1393c48a1e8SJed Brown       DMDACoor2d **coors;
1403c48a1e8SJed Brown       Vec      vcoors;
1413c48a1e8SJed Brown       DM       cda;
1422150357eSBarry Smith       DMGetCoordinates(da,&vcoors);
1431957e957SBarry Smith       DMGetCoordinateDM(da,&cda);
1443c48a1e8SJed Brown       DMDAVecGetArray(cda,vcoors,&coors);
1453c48a1e8SJed Brown       DMDAGetCorners(cda,&mstart,&nstart,0,&m,&n,0)
1463c48a1e8SJed Brown       for (i=mstart; i<mstart+m; i++) {
1473c48a1e8SJed Brown         for (j=nstart; j<nstart+n; j++) {
1483c48a1e8SJed Brown           x = coors[j][i].x;
1493c48a1e8SJed Brown           y = coors[j][i].y;
1503c48a1e8SJed Brown           ......
1513c48a1e8SJed Brown         }
1523c48a1e8SJed Brown       }
1533c48a1e8SJed Brown       DMDAVecRestoreArray(dac,vcoors,&coors);
1543c48a1e8SJed Brown 
1558e897e54SSatish Balay .seealso: DMDACoor3d, DMGetCoordinateDM(), DMGetCoordinates()
1563c48a1e8SJed Brown M*/
1573c48a1e8SJed Brown typedef struct {PetscScalar x,y;} DMDACoor2d;
1583c48a1e8SJed Brown 
1593c48a1e8SJed Brown /*MC
1603c48a1e8SJed Brown       DMDACoor3d - Structure for holding 3d (x, y and z) coordinates.
1613c48a1e8SJed Brown 
1623c48a1e8SJed Brown     Level: intermediate
1633c48a1e8SJed Brown 
1648e897e54SSatish Balay     Synopsis:
1653c48a1e8SJed Brown       DMDACoor3d ***coors;
1663c48a1e8SJed Brown       Vec      vcoors;
1673c48a1e8SJed Brown       DM       cda;
1682150357eSBarry Smith       DMGetCoordinates(da,&vcoors);
1691957e957SBarry Smith       DMGetCoordinateDM(da,&cda);
1703c48a1e8SJed Brown       DMDAVecGetArray(cda,vcoors,&coors);
1713c48a1e8SJed Brown       DMDAGetCorners(cda,&mstart,&nstart,&pstart,&m,&n,&p)
1723c48a1e8SJed Brown       for (i=mstart; i<mstart+m; i++) {
1733c48a1e8SJed Brown         for (j=nstart; j<nstart+n; j++) {
1743c48a1e8SJed Brown           for (k=pstart; k<pstart+p; k++) {
1753c48a1e8SJed Brown             x = coors[k][j][i].x;
1763c48a1e8SJed Brown             y = coors[k][j][i].y;
1773c48a1e8SJed Brown             z = coors[k][j][i].z;
1783c48a1e8SJed Brown           ......
1793c48a1e8SJed Brown         }
1803c48a1e8SJed Brown       }
1813c48a1e8SJed Brown       DMDAVecRestoreArray(dac,vcoors,&coors);
1823c48a1e8SJed Brown 
1838e897e54SSatish Balay .seealso: DMDACoor2d, DMGetCoordinateDM(), DMGetCoordinates()
1843c48a1e8SJed Brown M*/
1853c48a1e8SJed Brown typedef struct {PetscScalar x,y,z;} DMDACoor3d;
1863c48a1e8SJed Brown 
187014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGetLocalInfo(DM,DMDALocalInfo*);
1883c48a1e8SJed Brown 
189014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRegisterDAAD(void);
190014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateDAAD(DM,Mat*);
191014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqUSFFT(Vec,DM,Mat*);
1923c48a1e8SJed Brown 
193b412c318SBarry Smith PETSC_EXTERN PetscErrorCode DMDASetGetMatrix(DM,PetscErrorCode (*)(DM, Mat *));
194ce308e1dSBarry Smith PETSC_EXTERN PetscErrorCode DMDASetBlockFills(DM,const PetscInt*,const PetscInt*);
195014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDASetRefinementFactor(DM,PetscInt,PetscInt,PetscInt);
196014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGetRefinementFactor(DM,PetscInt*,PetscInt*,PetscInt*);
1973c48a1e8SJed Brown 
198014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGetArray(DM,PetscBool ,void*);
199014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDARestoreArray(DM,PetscBool ,void*);
2003c48a1e8SJed Brown 
201014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDACreatePF(DM,PF*);
2023c48a1e8SJed Brown 
2033582350cSMatthew G. Knepley /* Emulation of DMPlex */
2043582350cSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMDAGetNumCells(DM, PetscInt *, PetscInt *, PetscInt *, PetscInt *);
20548710be7SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMDAGetCellPoint(DM, PetscInt, PetscInt, PetscInt, PetscInt *);
20657459e9aSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMDAGetNumVertices(DM, PetscInt *, PetscInt *, PetscInt *, PetscInt *);
20757459e9aSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMDAGetNumFaces(DM, PetscInt *, PetscInt *, PetscInt *, PetscInt *, PetscInt *, PetscInt *);
20857459e9aSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMDAGetHeightStratum(DM, PetscInt, PetscInt *, PetscInt *);
209793f3fe5SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMDAGetDepthStratum(DM, PetscInt, PetscInt *, PetscInt *);
21044e1f9abSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMDACreateSection(DM, const PetscInt[], const PetscInt[], const PetscInt[], PetscSection *);
2118e0841e0SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMDAComputeCellGeometryFEM(DM, PetscInt, PetscQuadrature, PetscReal [], PetscReal [], PetscReal [], PetscReal []);
2126693a731SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMDAGetTransitiveClosure(DM, PetscInt, PetscBool, PetscInt *, PetscInt **);
2136693a731SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMDARestoreTransitiveClosure(DM, PetscInt, PetscBool, PetscInt *, PetscInt **);
214d7a12f93SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMDAVecGetClosure(DM, PetscSection, Vec, PetscInt, PetscInt *, PetscScalar **);
215d7a12f93SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMDAVecRestoreClosure(DM, PetscSection, Vec, PetscInt, PetscInt *, PetscScalar **);
21657459e9aSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMDAVecSetClosure(DM, PetscSection, Vec, PetscInt, const PetscScalar *, InsertMode);
217aa1993deSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMDAGetClosure(DM, PetscSection, PetscInt, PetscInt*, const PetscInt**);
218aa1993deSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMDARestoreClosure(DM, PetscSection, PetscInt, PetscInt*, const PetscInt**);
219d7a12f93SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMDAGetClosureScalar(DM, PetscSection, PetscInt, PetscScalar*, PetscInt*, PetscScalar**);
220d7a12f93SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMDARestoreClosureScalar(DM, PetscSection, PetscInt, PetscScalar*, PetscInt*, PetscScalar**);
221aa1993deSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMDASetClosureScalar(DM,PetscSection,PetscInt,PetscScalar*,const PetscScalar*,InsertMode);
222f0e3914cSSatish Balay PETSC_EXTERN PetscErrorCode DMDAConvertToCell(DM, MatStencil, PetscInt *);
2233385ff7eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMDASetVertexCoordinates(DM,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal);
2248ca2f13cSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMDASetPreallocationCenterDimension(DM, PetscInt);
2258ca2f13cSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMDAGetPreallocationCenterDimension(DM, PetscInt*);
22680800b1aSMatthew G Knepley 
2273c48a1e8SJed Brown #endif
228