xref: /petsc/include/petscdmda.h (revision ffeef943c8ee50edff320d8a3135bb0c94853e4c)
1a4963045SJacob Faibussowitsch #pragma once
23c48a1e8SJed Brown 
32c8e378dSBarry Smith #include <petscdm.h>
41e25c274SJed Brown #include <petscdmdatypes.h>
52c8e378dSBarry Smith #include <petscpf.h>
62c8e378dSBarry Smith #include <petscao.h>
7bb619f44SMatthew G. Knepley #include <petscfe.h>
83c48a1e8SJed Brown 
9ac09b921SBarry Smith /* SUBMANSEC = DMDA */
10ac09b921SBarry Smith 
113c48a1e8SJed Brown /*MC
123c48a1e8SJed Brown      DMDA_STENCIL_STAR - "Star"-type stencil. In logical grid coordinates, only (i,j,k), (i+s,j,k), (i,j+s,k),
133c48a1e8SJed Brown                          (i,j,k+s) are in the stencil  NOT, for example, (i+s,j+s,k)
143c48a1e8SJed Brown 
153c48a1e8SJed Brown      Level: beginner
163c48a1e8SJed Brown 
1716a05f60SBarry Smith      Note:
1895bd0b28SBarry Smith      Determines what ghost point values are brought over to each process in `DMGlobalToLocalBegin()`/ `DMGlobalToLocalEnd()`; in this case the "corner" values are not
191957e957SBarry Smith      brought over and hence should not be accessed locally
201957e957SBarry Smith 
21af27ebaaSBarry Smith .seealso: [](ch_dmbase), `DMDA`, `DMDA_STENCIL_BOX`, `DMDAStencilType`, `DMDASetStencilType()`
223c48a1e8SJed Brown M*/
233c48a1e8SJed Brown 
243c48a1e8SJed Brown /*MC
253c48a1e8SJed Brown      DMDA_STENCIL_BOX - "Box"-type stencil. In logical grid coordinates, any of (i,j,k), (i+s,j+r,k+t) may
263c48a1e8SJed Brown                         be in the stencil.
273c48a1e8SJed Brown 
283c48a1e8SJed Brown      Level: beginner
293c48a1e8SJed Brown 
3095bd0b28SBarry Smith      Note:
3195bd0b28SBarry Smith      Determines what ghost point values are brought over to each process in `DMGlobalToLocalBegin()`/ `DMGlobalToLocalEnd()`
3295bd0b28SBarry Smith 
33af27ebaaSBarry Smith .seealso: [](ch_dmbase), `DMDA`, `DMDA_STENCIL_STAR`, `DMDAStencilType`, `DMDASetStencilType()`
343c48a1e8SJed Brown M*/
353c48a1e8SJed Brown 
36014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDASetInterpolationType(DM, DMDAInterpolationType);
37014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGetInterpolationType(DM, DMDAInterpolationType *);
3897779f9aSLisandro Dalcin PETSC_EXTERN PetscErrorCode DMDACreateAggregates(DM, DM, Mat *);
393c48a1e8SJed Brown 
40d4a6ed37SStefano Zampini /* FEM */
41014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDASetElementType(DM, DMDAElementType);
42014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGetElementType(DM, DMDAElementType *);
43014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGetElements(DM, PetscInt *, PetscInt *, const PetscInt *[]);
44014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDARestoreElements(DM, PetscInt *, PetscInt *, const PetscInt *[]);
45d4a6ed37SStefano Zampini PETSC_EXTERN PetscErrorCode DMDAGetElementsSizes(DM, PetscInt *, PetscInt *, PetscInt *);
46d4a6ed37SStefano Zampini PETSC_EXTERN PetscErrorCode DMDAGetElementsCorners(DM, PetscInt *, PetscInt *, PetscInt *);
47d4a6ed37SStefano Zampini PETSC_EXTERN PetscErrorCode DMDAGetSubdomainCornersIS(DM, IS *);
48d4a6ed37SStefano Zampini PETSC_EXTERN PetscErrorCode DMDARestoreSubdomainCornersIS(DM, IS *);
493c48a1e8SJed Brown 
503c48a1e8SJed Brown #define MATSEQUSFFT "sequsfft"
513c48a1e8SJed Brown 
52014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDACreate(MPI_Comm, DM *);
53014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDASetSizes(DM, PetscInt, PetscInt, PetscInt);
54bff4a2f0SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMDACreate1d(MPI_Comm, DMBoundaryType, PetscInt, PetscInt, PetscInt, const PetscInt[], DM *);
55bff4a2f0SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMDACreate2d(MPI_Comm, DMBoundaryType, DMBoundaryType, DMDAStencilType, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], DM *);
56bff4a2f0SMatthew 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 *);
573c48a1e8SJed Brown 
58014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGlobalToNaturalBegin(DM, Vec, InsertMode, Vec);
59014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGlobalToNaturalEnd(DM, Vec, InsertMode, Vec);
60014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDANaturalToGlobalBegin(DM, Vec, InsertMode, Vec);
61014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDANaturalToGlobalEnd(DM, Vec, InsertMode, Vec);
62edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 5, 0, "DMLocalToLocalBegin()", ) static inline PetscErrorCode DMDALocalToLocalBegin(DM dm, Vec g, InsertMode mode, Vec l)
63d71ae5a4SJacob Faibussowitsch {
649371c9d4SSatish Balay   return DMLocalToLocalBegin(dm, g, mode, l);
659371c9d4SSatish Balay }
66edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 5, 0, "DMLocalToLocalEnd()", ) static inline PetscErrorCode DMDALocalToLocalEnd(DM dm, Vec g, InsertMode mode, Vec l)
67d71ae5a4SJacob Faibussowitsch {
689371c9d4SSatish Balay   return DMLocalToLocalEnd(dm, g, mode, l);
699371c9d4SSatish Balay }
70014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDACreateNaturalVector(DM, Vec *);
713c48a1e8SJed Brown 
72014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGetCorners(DM, PetscInt *, PetscInt *, PetscInt *, PetscInt *, PetscInt *, PetscInt *);
73014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGetGhostCorners(DM, PetscInt *, PetscInt *, PetscInt *, PetscInt *, PetscInt *, PetscInt *);
74bff4a2f0SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMDAGetInfo(DM, PetscInt *, PetscInt *, PetscInt *, PetscInt *, PetscInt *, PetscInt *, PetscInt *, PetscInt *, PetscInt *, DMBoundaryType *, DMBoundaryType *, DMBoundaryType *, DMDAStencilType *);
753ee9839eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMDAGetProcessorSubset(DM, DMDirection, PetscInt, MPI_Comm *);
763ee9839eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMDAGetProcessorSubsets(DM, DMDirection, MPI_Comm *);
773ee9839eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMDAGetRay(DM, DMDirection, PetscInt, Vec *, VecScatter *);
783c48a1e8SJed Brown 
79014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGlobalToNaturalAllCreate(DM, VecScatter *);
80014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDANaturalAllToGlobalCreate(DM, VecScatter *);
813c48a1e8SJed Brown 
82bd1fc5aeSBarry Smith PETSC_EXTERN PetscErrorCode DMDAGetScatter(DM, VecScatter *, VecScatter *);
83014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGetNeighbors(DM, const PetscMPIInt **);
843c48a1e8SJed Brown 
859db3d8bcSStefano Zampini PETSC_EXTERN PetscErrorCode DMDASetAOType(DM, AOType);
86014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGetAO(DM, AO *);
87014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDASetUniformCoordinates(DM, PetscReal, PetscReal, PetscReal, PetscReal, PetscReal, PetscReal);
888272889dSSatish Balay PETSC_EXTERN PetscErrorCode DMDASetGLLCoordinates(DM, PetscInt, PetscReal *);
89c593f006SBarry Smith PETSC_EXTERN PetscErrorCode DMDAGetCoordinateArray(DM, void *);
90c593f006SBarry Smith PETSC_EXTERN PetscErrorCode DMDARestoreCoordinateArray(DM, void *);
91a74ba6f7SBarry Smith PETSC_EXTERN PetscErrorCode DMDAGetLogicalCoordinate(DM, PetscScalar, PetscScalar, PetscScalar, PetscInt *, PetscInt *, PetscInt *, PetscScalar *, PetscScalar *, PetscScalar *);
92f0aa4865SJed Brown /* function to wrap coordinates around boundary */
93014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAMapCoordsToPeriodicDomain(DM, PetscScalar *, PetscScalar *);
94f0aa4865SJed Brown 
95211d3afbSPatrick Sanan PETSC_EXTERN PetscErrorCode DMDACreateCompatibleDMDA(DM, PetscInt, DM *);
96edd03b47SJacob Faibussowitsch PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 10, 0, "DMDACreateCompatibleDMDA()", ) PetscErrorCode DMDAGetReducedDMDA(DM, PetscInt, DM *);
973c48a1e8SJed Brown 
98014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDASetFieldName(DM, PetscInt, const char[]);
99014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGetFieldName(DM, PetscInt, const char **);
100c629b14aSBarry Smith PETSC_EXTERN PetscErrorCode DMDASetFieldNames(DM, const char *const *);
101c629b14aSBarry Smith PETSC_EXTERN PetscErrorCode DMDAGetFieldNames(DM, const char *const **);
102109c9344SBarry Smith PETSC_EXTERN PetscErrorCode DMDASetCoordinateName(DM, PetscInt, const char[]);
103109c9344SBarry Smith PETSC_EXTERN PetscErrorCode DMDAGetCoordinateName(DM, PetscInt, const char **);
1043c48a1e8SJed Brown 
105bff4a2f0SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMDASetBoundaryType(DM, DMBoundaryType, DMBoundaryType, DMBoundaryType);
106fe58071bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMDAGetBoundaryType(DM, DMBoundaryType *, DMBoundaryType *, DMBoundaryType *);
107014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDASetDof(DM, PetscInt);
108fb6725baSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMDAGetDof(DM, PetscInt *);
1097ddda789SPeter Brune PETSC_EXTERN PetscErrorCode DMDASetOverlap(DM, PetscInt, PetscInt, PetscInt);
1107ddda789SPeter Brune PETSC_EXTERN PetscErrorCode DMDAGetOverlap(DM, PetscInt *, PetscInt *, PetscInt *);
1113e7870d2SPeter Brune PETSC_EXTERN PetscErrorCode DMDASetNumLocalSubDomains(DM, PetscInt);
1123e7870d2SPeter Brune PETSC_EXTERN PetscErrorCode DMDAGetNumLocalSubDomains(DM, PetscInt *);
11395c13181SPeter Brune PETSC_EXTERN PetscErrorCode DMDAGetOffset(DM, PetscInt *, PetscInt *, PetscInt *, PetscInt *, PetscInt *, PetscInt *);
11495c13181SPeter Brune PETSC_EXTERN PetscErrorCode DMDASetOffset(DM, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt);
11540234c92SPeter Brune PETSC_EXTERN PetscErrorCode DMDAGetNonOverlappingRegion(DM, PetscInt *, PetscInt *, PetscInt *, PetscInt *, PetscInt *, PetscInt *);
11640234c92SPeter Brune PETSC_EXTERN PetscErrorCode DMDASetNonOverlappingRegion(DM, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt);
117014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDASetStencilWidth(DM, PetscInt);
118fb6725baSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMDAGetStencilWidth(DM, PetscInt *);
11938fb4e8eSJunchao Zhang PETSC_EXTERN PetscErrorCode DMDAMapMatStencilToGlobal(DM, PetscInt, const MatStencil[], PetscInt[]);
120014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDASetOwnershipRanges(DM, const PetscInt[], const PetscInt[], const PetscInt[]);
121014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGetOwnershipRanges(DM, const PetscInt **, const PetscInt **, const PetscInt **);
122014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDASetNumProcs(DM, PetscInt, PetscInt, PetscInt);
123014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDASetStencilType(DM, DMDAStencilType);
124fb6725baSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMDAGetStencilType(DM, DMDAStencilType *);
1253c48a1e8SJed Brown 
126014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAVecGetArray(DM, Vec, void *);
127014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAVecRestoreArray(DM, Vec, void *);
128fdc842d1SBarry Smith PETSC_EXTERN PetscErrorCode DMDAVecGetArrayWrite(DM, Vec, void *);
129fdc842d1SBarry Smith PETSC_EXTERN PetscErrorCode DMDAVecRestoreArrayWrite(DM, Vec, void *);
1303c48a1e8SJed Brown 
131014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAVecGetArrayDOF(DM, Vec, void *);
132014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAVecRestoreArrayDOF(DM, Vec, void *);
1333c48a1e8SJed Brown 
1345edff71fSBarry Smith PETSC_EXTERN PetscErrorCode DMDAVecGetArrayRead(DM, Vec, void *);
1355edff71fSBarry Smith PETSC_EXTERN PetscErrorCode DMDAVecRestoreArrayRead(DM, Vec, void *);
1365edff71fSBarry Smith 
1375edff71fSBarry Smith PETSC_EXTERN PetscErrorCode DMDAVecGetArrayDOFRead(DM, Vec, void *);
1385edff71fSBarry Smith PETSC_EXTERN PetscErrorCode DMDAVecRestoreArrayDOFRead(DM, Vec, void *);
1395edff71fSBarry Smith 
1401e5d2365SBarry Smith PETSC_EXTERN PetscErrorCode DMDAVecGetArrayDOFWrite(DM, Vec, void *);
1411e5d2365SBarry Smith PETSC_EXTERN PetscErrorCode DMDAVecRestoreArrayDOFWrite(DM, Vec, void *);
1421e5d2365SBarry Smith 
1433b5e53cdSSajid Ali PETSC_EXTERN PetscErrorCode DMDACreatePatchIS(DM, MatStencil *, MatStencil *, IS *, PetscBool);
144ff9846d9SPeter Brune 
1453c48a1e8SJed Brown /*MC
14687497f52SBarry Smith       DMDACoor2d - Structure for holding 2d (x and y) coordinates when working with `DMDA`
1473c48a1e8SJed Brown 
1488e897e54SSatish Balay     Synopsis:
14987497f52SBarry Smith .vb
1503c48a1e8SJed Brown       DMDACoor2d **coors;
1513c48a1e8SJed Brown       Vec      vcoors;
1523c48a1e8SJed Brown       DM       cda;
1532150357eSBarry Smith       DMGetCoordinates(da,&vcoors);
1541957e957SBarry Smith       DMGetCoordinateDM(da,&cda);
1553c48a1e8SJed Brown       DMDAVecGetArray(cda,vcoors,&coors);
1563c48a1e8SJed Brown       DMDAGetCorners(cda,&mstart,&nstart,0,&m,&n,0)
1573c48a1e8SJed Brown       for (i=mstart; i<mstart+m; i++) {
1583c48a1e8SJed Brown         for (j=nstart; j<nstart+n; j++) {
1593c48a1e8SJed Brown           x = coors[j][i].x;
1603c48a1e8SJed Brown           y = coors[j][i].y;
1613c48a1e8SJed Brown           ......
1623c48a1e8SJed Brown         }
1633c48a1e8SJed Brown       }
1643c48a1e8SJed Brown       DMDAVecRestoreArray(dac,vcoors,&coors);
16587497f52SBarry Smith .ve
1663c48a1e8SJed Brown 
16716a05f60SBarry Smith     Level: intermediate
16816a05f60SBarry Smith 
169af27ebaaSBarry Smith .seealso: [](ch_dmbase), `DMDA`, `DMDACoor3d`, `DMDAVecRestoreArray()`, `DMDAVecGetArray()`, `DMGetCoordinateDM()`, `DMGetCoordinates()`
1703c48a1e8SJed Brown M*/
1719371c9d4SSatish Balay typedef struct {
1729371c9d4SSatish Balay   PetscScalar x, y;
1739371c9d4SSatish Balay } DMDACoor2d;
1743c48a1e8SJed Brown 
1753c48a1e8SJed Brown /*MC
17687497f52SBarry Smith       DMDACoor3d - Structure for holding 3d (x, y and z) coordinates  coordinates when working with `DMDA`
1773c48a1e8SJed Brown 
1788e897e54SSatish Balay     Synopsis:
17987497f52SBarry Smith .vb
1803c48a1e8SJed Brown       DMDACoor3d ***coors;
1813c48a1e8SJed Brown       Vec      vcoors;
1823c48a1e8SJed Brown       DM       cda;
1832150357eSBarry Smith       DMGetCoordinates(da,&vcoors);
1841957e957SBarry Smith       DMGetCoordinateDM(da,&cda);
1853c48a1e8SJed Brown       DMDAVecGetArray(cda,vcoors,&coors);
1863c48a1e8SJed Brown       DMDAGetCorners(cda,&mstart,&nstart,&pstart,&m,&n,&p)
1873c48a1e8SJed Brown       for (i=mstart; i<mstart+m; i++) {
1883c48a1e8SJed Brown         for (j=nstart; j<nstart+n; j++) {
1893c48a1e8SJed Brown           for (k=pstart; k<pstart+p; k++) {
1903c48a1e8SJed Brown             x = coors[k][j][i].x;
1913c48a1e8SJed Brown             y = coors[k][j][i].y;
1923c48a1e8SJed Brown             z = coors[k][j][i].z;
1933c48a1e8SJed Brown           ......
1943c48a1e8SJed Brown         }
1953c48a1e8SJed Brown       }
1963c48a1e8SJed Brown       DMDAVecRestoreArray(dac,vcoors,&coors);
19787497f52SBarry Smith .ve
1983c48a1e8SJed Brown 
19916a05f60SBarry Smith     Level: intermediate
20016a05f60SBarry Smith 
201af27ebaaSBarry Smith .seealso: [](ch_dmbase), `DMDA`, `DMDACoor2d`, `DMDAVecRestoreArray()`, `DMDAVecGetArray()`, `DMGetCoordinateDM()`, `DMGetCoordinates()`
2023c48a1e8SJed Brown M*/
2039371c9d4SSatish Balay typedef struct {
2049371c9d4SSatish Balay   PetscScalar x, y, z;
2059371c9d4SSatish Balay } DMDACoor3d;
2063c48a1e8SJed Brown 
207014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGetLocalInfo(DM, DMDALocalInfo *);
2083c48a1e8SJed Brown 
209014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRegisterDAAD(void);
210014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqUSFFT(Vec, DM, Mat *);
2113c48a1e8SJed Brown 
212b412c318SBarry Smith PETSC_EXTERN PetscErrorCode DMDASetGetMatrix(DM, PetscErrorCode (*)(DM, Mat *));
213ce308e1dSBarry Smith PETSC_EXTERN PetscErrorCode DMDASetBlockFills(DM, const PetscInt *, const PetscInt *);
21409e28618SBarry Smith PETSC_EXTERN PetscErrorCode DMDASetBlockFillsSparse(DM, const PetscInt *, const PetscInt *);
215014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDASetRefinementFactor(DM, PetscInt, PetscInt, PetscInt);
216014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGetRefinementFactor(DM, PetscInt *, PetscInt *, PetscInt *);
2173c48a1e8SJed Brown 
218014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGetArray(DM, PetscBool, void *);
219014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDARestoreArray(DM, PetscBool, void *);
2203c48a1e8SJed Brown 
221014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDACreatePF(DM, PF *);
2223c48a1e8SJed Brown 
2233582350cSMatthew G. Knepley /* Emulation of DMPlex */
2243582350cSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMDAGetNumCells(DM, PetscInt *, PetscInt *, PetscInt *, PetscInt *);
22548710be7SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMDAGetCellPoint(DM, PetscInt, PetscInt, PetscInt, PetscInt *);
22657459e9aSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMDAGetNumVertices(DM, PetscInt *, PetscInt *, PetscInt *, PetscInt *);
22757459e9aSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMDAGetNumFaces(DM, PetscInt *, PetscInt *, PetscInt *, PetscInt *, PetscInt *, PetscInt *);
22857459e9aSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMDAGetHeightStratum(DM, PetscInt, PetscInt *, PetscInt *);
229793f3fe5SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMDAGetDepthStratum(DM, PetscInt, PetscInt *, PetscInt *);
230f0e3914cSSatish Balay PETSC_EXTERN PetscErrorCode DMDAConvertToCell(DM, MatStencil, PetscInt *);
2313385ff7eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMDASetVertexCoordinates(DM, PetscReal, PetscReal, PetscReal, PetscReal, PetscReal, PetscReal);
2328ca2f13cSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMDASetPreallocationCenterDimension(DM, PetscInt);
2338ca2f13cSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMDAGetPreallocationCenterDimension(DM, PetscInt *);
234*ffeef943SBarry Smith 
235*ffeef943SBarry Smith PETSC_EXTERN PetscErrorCode DMDAVTKWriteAll(PetscObject, PetscViewer);
236