xref: /petsc/include/petscdmda.h (revision 95c1318154c5a85ba7f10eb2e425da809e409add)
13c48a1e8SJed Brown #if !defined(__PETSCDMDA_H)
23c48a1e8SJed Brown #define __PETSCDMDA_H
33c48a1e8SJed Brown 
42c8e378dSBarry Smith #include <petscdm.h>
52c8e378dSBarry Smith #include <petscpf.h>
62c8e378dSBarry Smith #include <petscao.h>
73c48a1e8SJed Brown 
83c48a1e8SJed Brown /*E
93c48a1e8SJed Brown     DMDAStencilType - Determines if the stencil extends only along the coordinate directions, or also
103c48a1e8SJed Brown       to the northeast, northwest etc
113c48a1e8SJed Brown 
123c48a1e8SJed Brown    Level: beginner
133c48a1e8SJed Brown 
14db05f41bSBarry Smith .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDACreate(), DMDASetStencilType()
153c48a1e8SJed Brown E*/
163c48a1e8SJed Brown typedef enum { DMDA_STENCIL_STAR,DMDA_STENCIL_BOX } DMDAStencilType;
173c48a1e8SJed Brown 
183c48a1e8SJed Brown /*MC
193c48a1e8SJed Brown      DMDA_STENCIL_STAR - "Star"-type stencil. In logical grid coordinates, only (i,j,k), (i+s,j,k), (i,j+s,k),
203c48a1e8SJed Brown                        (i,j,k+s) are in the stencil  NOT, for example, (i+s,j+s,k)
213c48a1e8SJed Brown 
223c48a1e8SJed Brown      Level: beginner
233c48a1e8SJed Brown 
24db05f41bSBarry Smith .seealso: DMDA_STENCIL_BOX, DMDAStencilType, DMDASetStencilType()
253c48a1e8SJed Brown M*/
263c48a1e8SJed Brown 
273c48a1e8SJed Brown /*MC
283c48a1e8SJed Brown      DMDA_STENCIL_BOX - "Box"-type stencil. In logical grid coordinates, any of (i,j,k), (i+s,j+r,k+t) may
293c48a1e8SJed Brown                       be in the stencil.
303c48a1e8SJed Brown 
313c48a1e8SJed Brown      Level: beginner
323c48a1e8SJed Brown 
33db05f41bSBarry Smith .seealso: DMDA_STENCIL_STAR, DMDAStencilType, DMDASetStencilType()
343c48a1e8SJed Brown M*/
353c48a1e8SJed Brown 
363c48a1e8SJed Brown /*E
375d26d725SBarry Smith     DMDABoundaryType - Describes the choice for fill of ghost cells on physical domain boundaries.
383c48a1e8SJed Brown 
393c48a1e8SJed Brown    Level: beginner
403c48a1e8SJed Brown 
413c48a1e8SJed Brown    A boundary may be of type DMDA_BOUNDARY_NONE (no ghost nodes), DMDA_BOUNDARY_GHOST (ghost nodes
425d26d725SBarry Smith    exist but aren't filled, you can put values into them and then apply a stencil that uses those ghost locations),
43db05f41bSBarry Smith    DMDA_BOUNDARY_MIRROR (not yet implemented for 3d), or DMDA_BOUNDARY_PERIODIC
443c48a1e8SJed Brown    (ghost nodes filled by the opposite edge of the domain).
453c48a1e8SJed Brown 
465d26d725SBarry Smith    Note: This is information for the boundary of the __PHYSICAL__ domain. It has nothing to do with boundaries between
475d26d725SBarry Smith      processes, that width is always determined by the stencil width, see DMDASetStencilWidth().
485d26d725SBarry Smith 
493c48a1e8SJed Brown .seealso: DMDASetBoundaryType(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDACreate()
503c48a1e8SJed Brown E*/
513c48a1e8SJed Brown typedef enum { DMDA_BOUNDARY_NONE, DMDA_BOUNDARY_GHOSTED, DMDA_BOUNDARY_MIRROR, DMDA_BOUNDARY_PERIODIC } DMDABoundaryType;
523c48a1e8SJed Brown 
536a6fc655SJed Brown PETSC_EXTERN const char *const DMDABoundaryTypes[];
543c48a1e8SJed Brown 
553c48a1e8SJed Brown /*E
563c48a1e8SJed Brown     DMDAInterpolationType - Defines the type of interpolation that will be returned by
57e727c939SJed Brown        DMCreateInterpolation.
583c48a1e8SJed Brown 
593c48a1e8SJed Brown    Level: beginner
603c48a1e8SJed Brown 
61e727c939SJed Brown .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMCreateInterpolation(), DMDASetInterpolationType(), DMDACreate()
623c48a1e8SJed Brown E*/
633c48a1e8SJed Brown typedef enum { DMDA_Q0, DMDA_Q1 } DMDAInterpolationType;
643c48a1e8SJed Brown 
65014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDASetInterpolationType(DM,DMDAInterpolationType);
66014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGetInterpolationType(DM,DMDAInterpolationType*);
673c48a1e8SJed Brown 
683c48a1e8SJed Brown /*E
693c48a1e8SJed Brown     DMDAElementType - Defines the type of elements that will be returned by
702dde6fd4SLisandro Dalcin        DMDAGetElements()
713c48a1e8SJed Brown 
723c48a1e8SJed Brown    Level: beginner
733c48a1e8SJed Brown 
74e727c939SJed Brown .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMCreateInterpolation(), DMDASetInterpolationType(),
752dde6fd4SLisandro Dalcin           DMDASetElementType(), DMDAGetElements(), DMDARestoreElements(), DMDACreate()
763c48a1e8SJed Brown E*/
773c48a1e8SJed Brown typedef enum { DMDA_ELEMENT_P1, DMDA_ELEMENT_Q1 } DMDAElementType;
783c48a1e8SJed Brown 
79014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDASetElementType(DM,DMDAElementType);
80014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGetElementType(DM,DMDAElementType*);
81014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGetElements(DM,PetscInt *,PetscInt *,const PetscInt*[]);
82014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDARestoreElements(DM,PetscInt *,PetscInt *,const PetscInt*[]);
833c48a1e8SJed Brown 
843c48a1e8SJed Brown typedef enum { DMDA_X,DMDA_Y,DMDA_Z } DMDADirection;
853c48a1e8SJed Brown 
863c48a1e8SJed Brown #define MATSEQUSFFT        "sequsfft"
873c48a1e8SJed Brown 
88014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDACreate(MPI_Comm,DM*);
89014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDASetDim(DM,PetscInt);
90014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDASetSizes(DM,PetscInt,PetscInt,PetscInt);
91014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDACreate1d(MPI_Comm,DMDABoundaryType,PetscInt,PetscInt,PetscInt,const PetscInt[],DM *);
92014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDACreate2d(MPI_Comm,DMDABoundaryType,DMDABoundaryType,DMDAStencilType,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],DM*);
93014dd563SJed Brown 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*);
943c48a1e8SJed Brown 
95014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGlobalToNaturalBegin(DM,Vec,InsertMode,Vec);
96014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGlobalToNaturalEnd(DM,Vec,InsertMode,Vec);
97014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDANaturalToGlobalBegin(DM,Vec,InsertMode,Vec);
98014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDANaturalToGlobalEnd(DM,Vec,InsertMode,Vec);
99014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDALocalToLocalBegin(DM,Vec,InsertMode,Vec);
100014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDALocalToLocalEnd(DM,Vec,InsertMode,Vec);
101014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDACreateNaturalVector(DM,Vec *);
1023c48a1e8SJed Brown 
103014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGetCorners(DM,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
104014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGetGhostCorners(DM,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
105014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGetInfo(DM,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,DMDABoundaryType*,DMDABoundaryType*,DMDABoundaryType*,DMDAStencilType*);
106014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGetProcessorSubset(DM,DMDADirection,PetscInt,MPI_Comm*);
107014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGetProcessorSubsets(DM,DMDADirection,MPI_Comm*);
108db05f41bSBarry Smith PETSC_EXTERN PetscErrorCode DMDAGetRay(DM,DMDADirection,PetscInt,Vec*,VecScatter*);
1093c48a1e8SJed Brown 
110014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGlobalToNaturalAllCreate(DM,VecScatter*);
111014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDANaturalAllToGlobalCreate(DM,VecScatter*);
1123c48a1e8SJed Brown 
113014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGetGlobalIndices(DM,PetscInt*,PetscInt**);
1143c48a1e8SJed Brown 
115014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGetScatter(DM,VecScatter*,VecScatter*,VecScatter*);
116014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGetNeighbors(DM,const PetscMPIInt**);
1173c48a1e8SJed Brown 
118014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGetAO(DM,AO*);
119014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDASetUniformCoordinates(DM,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal);
120014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGetBoundingBox(DM,PetscReal[],PetscReal[]);
121014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGetLocalBoundingBox(DM,PetscReal[],PetscReal[]);
122f0aa4865SJed Brown /* function to wrap coordinates around boundary */
123014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAMapCoordsToPeriodicDomain(DM,PetscScalar*,PetscScalar*);
124f0aa4865SJed Brown 
125db05f41bSBarry Smith PETSC_EXTERN PetscErrorCode DMDAGetReducedDMDA(DM,PetscInt,DM*);
1263c48a1e8SJed Brown 
127014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDASetFieldName(DM,PetscInt,const char[]);
128014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGetFieldName(DM,PetscInt,const char**);
129109c9344SBarry Smith PETSC_EXTERN PetscErrorCode DMDASetCoordinateName(DM,PetscInt,const char[]);
130109c9344SBarry Smith PETSC_EXTERN PetscErrorCode DMDAGetCoordinateName(DM,PetscInt,const char**);
1313c48a1e8SJed Brown 
132014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDASetBoundaryType(DM,DMDABoundaryType,DMDABoundaryType,DMDABoundaryType);
133014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDASetDof(DM, PetscInt);
13488661749SPeter Brune PETSC_EXTERN PetscErrorCode DMDASetOverlap(DM, PetscInt);
135*95c13181SPeter Brune PETSC_EXTERN PetscErrorCode DMDAGetOffset(DM,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
136*95c13181SPeter Brune PETSC_EXTERN PetscErrorCode DMDASetOffset(DM,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt);
137014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDASetStencilWidth(DM, PetscInt);
138014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDASetOwnershipRanges(DM,const PetscInt[],const PetscInt[],const PetscInt[]);
139014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGetOwnershipRanges(DM,const PetscInt**,const PetscInt**,const PetscInt**);
140014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDASetNumProcs(DM, PetscInt, PetscInt, PetscInt);
141014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDASetStencilType(DM, DMDAStencilType);
1423c48a1e8SJed Brown 
143014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAVecGetArray(DM,Vec,void *);
144014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAVecRestoreArray(DM,Vec,void *);
1453c48a1e8SJed Brown 
146014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAVecGetArrayDOF(DM,Vec,void *);
147014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAVecRestoreArrayDOF(DM,Vec,void *);
1483c48a1e8SJed Brown 
149014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDASplitComm2d(MPI_Comm,PetscInt,PetscInt,PetscInt,MPI_Comm*);
1503c48a1e8SJed Brown 
151ff9846d9SPeter Brune PETSC_EXTERN PetscErrorCode DMDACreatePatchIS(DM,MatStencil*,MatStencil*,IS*);
152ff9846d9SPeter Brune 
1533c48a1e8SJed Brown /*S
1543c48a1e8SJed Brown      DMDALocalInfo - C struct that contains information about a structured grid and a processors logical
1553c48a1e8SJed Brown               location in it.
1563c48a1e8SJed Brown 
1573c48a1e8SJed Brown    Level: beginner
1583c48a1e8SJed Brown 
1593c48a1e8SJed Brown   Concepts: distributed array
1603c48a1e8SJed Brown 
1613c48a1e8SJed Brown   Developer note: Then entries in this struct are int instead of PetscInt so that the elements may
1623c48a1e8SJed Brown                   be extracted in Fortran as if from an integer array
1633c48a1e8SJed Brown 
1643c48a1e8SJed Brown .seealso:  DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DM, DMDAGetLocalInfo(), DMDAGetInfo()
1653c48a1e8SJed Brown S*/
1663c48a1e8SJed Brown typedef struct {
1673c48a1e8SJed Brown   PetscInt         dim,dof,sw;
1683c48a1e8SJed Brown   PetscInt         mx,my,mz;    /* global number of grid points in each direction */
1696ee1aa1aSDmitry Karpeev   PetscInt         xs,ys,zs;    /* starting point of this processor, excluding ghosts */
1703c48a1e8SJed Brown   PetscInt         xm,ym,zm;    /* number of grid points on this processor, excluding ghosts */
1713c48a1e8SJed Brown   PetscInt         gxs,gys,gzs;    /* starting point of this processor including ghosts */
1723c48a1e8SJed Brown   PetscInt         gxm,gym,gzm;    /* number of grid points on this processor including ghosts */
1733c48a1e8SJed Brown   DMDABoundaryType bx,by,bz; /* type of ghost nodes at boundary */
1743c48a1e8SJed Brown   DMDAStencilType  st;
1753c48a1e8SJed Brown   DM               da;
1763c48a1e8SJed Brown } DMDALocalInfo;
1773c48a1e8SJed Brown 
1783c48a1e8SJed Brown /*MC
1793c48a1e8SJed Brown       DMDAForEachPointBegin2d - Starts a loop over the local part of a two dimensional DMDA
1803c48a1e8SJed Brown 
1813c48a1e8SJed Brown    Synopsis:
182f2ba6396SBarry Smith    #include "petscdm.h"
1833c48a1e8SJed Brown    void  DMDAForEachPointBegin2d(DALocalInfo *info,PetscInt i,PetscInt j);
1843c48a1e8SJed Brown 
1853c48a1e8SJed Brown    Not Collective
1863c48a1e8SJed Brown 
1873c48a1e8SJed Brown    Level: intermediate
1883c48a1e8SJed Brown 
1893c48a1e8SJed Brown .seealso: DMDAForEachPointEnd2d(), DMDAVecGetArray()
1903c48a1e8SJed Brown M*/
1913c48a1e8SJed Brown #define DMDAForEachPointBegin2d(info,i,j) {\
1923c48a1e8SJed Brown   PetscInt _xints = info->xs,_xinte = info->xs+info->xm,_yints = info->ys,_yinte = info->ys+info->ym;\
1933c48a1e8SJed Brown   for (j=_yints; j<_yinte; j++) {\
1943c48a1e8SJed Brown     for (i=_xints; i<_xinte; i++) {\
1953c48a1e8SJed Brown 
1963c48a1e8SJed Brown /*MC
1973c48a1e8SJed Brown       DMDAForEachPointEnd2d - Ends a loop over the local part of a two dimensional DMDA
1983c48a1e8SJed Brown 
1993c48a1e8SJed Brown    Synopsis:
200f2ba6396SBarry Smith    #include "petscdm.h"
2013c48a1e8SJed Brown    void  DMDAForEachPointEnd2d;
2023c48a1e8SJed Brown 
2033c48a1e8SJed Brown    Not Collective
2043c48a1e8SJed Brown 
2053c48a1e8SJed Brown    Level: intermediate
2063c48a1e8SJed Brown 
2073c48a1e8SJed Brown .seealso: DMDAForEachPointBegin2d(), DMDAVecGetArray()
2083c48a1e8SJed Brown M*/
2093c48a1e8SJed Brown #define DMDAForEachPointEnd2d }}}
2103c48a1e8SJed Brown 
2113c48a1e8SJed Brown /*MC
2123c48a1e8SJed Brown       DMDACoor2d - Structure for holding 2d (x and y) coordinates.
2133c48a1e8SJed Brown 
2143c48a1e8SJed Brown     Level: intermediate
2153c48a1e8SJed Brown 
2163c48a1e8SJed Brown     Sample Usage:
2173c48a1e8SJed Brown       DMDACoor2d **coors;
2183c48a1e8SJed Brown       Vec      vcoors;
2193c48a1e8SJed Brown       DM       cda;
2203c48a1e8SJed Brown 
2212150357eSBarry Smith       DMGetCoordinates(da,&vcoors);
2223c48a1e8SJed Brown       DMDAGetCoordinateDA(da,&cda);
2233c48a1e8SJed Brown       DMDAVecGetArray(cda,vcoors,&coors);
2243c48a1e8SJed Brown       DMDAGetCorners(cda,&mstart,&nstart,0,&m,&n,0)
2253c48a1e8SJed Brown       for (i=mstart; i<mstart+m; i++) {
2263c48a1e8SJed Brown         for (j=nstart; j<nstart+n; j++) {
2273c48a1e8SJed Brown           x = coors[j][i].x;
2283c48a1e8SJed Brown           y = coors[j][i].y;
2293c48a1e8SJed Brown           ......
2303c48a1e8SJed Brown         }
2313c48a1e8SJed Brown       }
2323c48a1e8SJed Brown       DMDAVecRestoreArray(dac,vcoors,&coors);
2333c48a1e8SJed Brown 
2342150357eSBarry Smith .seealso: DMDACoor3d, DMDAForEachPointBegin(), DMDAGetCoordinateDA(), DMGetCoordinates(), DMDAGetGhostCoordinates()
2353c48a1e8SJed Brown M*/
2363c48a1e8SJed Brown typedef struct {PetscScalar x,y;} DMDACoor2d;
2373c48a1e8SJed Brown 
2383c48a1e8SJed Brown /*MC
2393c48a1e8SJed Brown       DMDACoor3d - Structure for holding 3d (x, y and z) coordinates.
2403c48a1e8SJed Brown 
2413c48a1e8SJed Brown     Level: intermediate
2423c48a1e8SJed Brown 
2433c48a1e8SJed Brown     Sample Usage:
2443c48a1e8SJed Brown       DMDACoor3d ***coors;
2453c48a1e8SJed Brown       Vec      vcoors;
2463c48a1e8SJed Brown       DM       cda;
2473c48a1e8SJed Brown 
2482150357eSBarry Smith       DMGetCoordinates(da,&vcoors);
2493c48a1e8SJed Brown       DMDAGetCoordinateDA(da,&cda);
2503c48a1e8SJed Brown       DMDAVecGetArray(cda,vcoors,&coors);
2513c48a1e8SJed Brown       DMDAGetCorners(cda,&mstart,&nstart,&pstart,&m,&n,&p)
2523c48a1e8SJed Brown       for (i=mstart; i<mstart+m; i++) {
2533c48a1e8SJed Brown         for (j=nstart; j<nstart+n; j++) {
2543c48a1e8SJed Brown           for (k=pstart; k<pstart+p; k++) {
2553c48a1e8SJed Brown             x = coors[k][j][i].x;
2563c48a1e8SJed Brown             y = coors[k][j][i].y;
2573c48a1e8SJed Brown             z = coors[k][j][i].z;
2583c48a1e8SJed Brown           ......
2593c48a1e8SJed Brown         }
2603c48a1e8SJed Brown       }
2613c48a1e8SJed Brown       DMDAVecRestoreArray(dac,vcoors,&coors);
2623c48a1e8SJed Brown 
2632150357eSBarry Smith .seealso: DMDACoor2d, DMDAForEachPointBegin(), DMDAGetCoordinateDA(), DMGetCoordinates(), DMDAGetGhostCoordinates()
2643c48a1e8SJed Brown M*/
2653c48a1e8SJed Brown typedef struct {PetscScalar x,y,z;} DMDACoor3d;
2663c48a1e8SJed Brown 
267014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGetLocalInfo(DM,DMDALocalInfo*);
2683c48a1e8SJed Brown 
269014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRegisterDAAD(void);
270014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateDAAD(DM,Mat*);
271014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqUSFFT(Vec,DM,Mat*);
2723c48a1e8SJed Brown 
27319fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode DMDASetGetMatrix(DM,PetscErrorCode (*)(DM, MatType,Mat *));
274ce308e1dSBarry Smith PETSC_EXTERN PetscErrorCode DMDASetBlockFills(DM,const PetscInt*,const PetscInt*);
275014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDASetRefinementFactor(DM,PetscInt,PetscInt,PetscInt);
276014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGetRefinementFactor(DM,PetscInt*,PetscInt*,PetscInt*);
2773c48a1e8SJed Brown 
278014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGetArray(DM,PetscBool ,void*);
279014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDARestoreArray(DM,PetscBool ,void*);
2803c48a1e8SJed Brown 
281014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDACreatePF(DM,PF*);
2823c48a1e8SJed Brown 
28357459e9aSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMDAGetNumCells(DM, PetscInt *);
28457459e9aSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMDAGetNumVertices(DM, PetscInt *, PetscInt *, PetscInt *, PetscInt *);
28557459e9aSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMDAGetNumFaces(DM, PetscInt *, PetscInt *, PetscInt *, PetscInt *, PetscInt *, PetscInt *);
28657459e9aSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMDAGetHeightStratum(DM, PetscInt, PetscInt *, PetscInt *);
287014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDACreateSection(DM, PetscInt[], PetscInt[], PetscInt[], PetscInt[]);
28857459e9aSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMDAComputeCellGeometry(DM, PetscInt, PetscQuadrature *, PetscReal [], PetscReal [], PetscReal [], PetscReal []);
28957459e9aSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMDAVecGetClosure(DM, PetscSection, Vec, PetscInt, const PetscScalar **);
29057459e9aSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMDAVecSetClosure(DM, PetscSection, Vec, PetscInt, const PetscScalar *, InsertMode);
291aa1993deSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMDAGetClosure(DM,PetscSection,PetscInt,PetscInt*,const PetscInt**);
292aa1993deSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMDARestoreClosure(DM,PetscSection,PetscInt,PetscInt*,const PetscInt**);
293aa1993deSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMDAGetClosureScalar(DM,PetscSection,PetscInt,PetscScalar*,const PetscScalar**);
294aa1993deSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMDARestoreClosureScalar(DM,PetscSection,PetscInt,PetscScalar*,const PetscScalar**);
295aa1993deSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMDASetClosureScalar(DM,PetscSection,PetscInt,PetscScalar*,const PetscScalar*,InsertMode);
296f0e3914cSSatish Balay PETSC_EXTERN PetscErrorCode DMDAConvertToCell(DM, MatStencil, PetscInt *);
29780800b1aSMatthew G Knepley 
2983c48a1e8SJed Brown #endif
299