xref: /petsc/include/petscdmda.h (revision 2150357e90fff80a116f8d4b85a252675c8fae1c)
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 
143c48a1e8SJed Brown .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDACreate()
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 
243c48a1e8SJed Brown .seealso: DMDA_STENCIL_BOX, DMDAStencilType
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 
333c48a1e8SJed Brown .seealso: DMDA_STENCIL_STAR, DMDAStencilType
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),
435d26d725SBarry Smith    DMDA_BOUNDARY_MIRROR (not yet implemented), 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*);
1083c48a1e8SJed Brown 
109014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGlobalToNaturalAllCreate(DM,VecScatter*);
110014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDANaturalAllToGlobalCreate(DM,VecScatter*);
1113c48a1e8SJed Brown 
112014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGetGlobalIndices(DM,PetscInt*,PetscInt**);
1133c48a1e8SJed Brown 
114014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGetScatter(DM,VecScatter*,VecScatter*,VecScatter*);
115014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGetNeighbors(DM,const PetscMPIInt**);
1163c48a1e8SJed Brown 
117014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGetAO(DM,AO*);
118014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDASetUniformCoordinates(DM,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal);
119014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGetBoundingBox(DM,PetscReal[],PetscReal[]);
120014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGetLocalBoundingBox(DM,PetscReal[],PetscReal[]);
121f0aa4865SJed Brown /* function to wrap coordinates around boundary */
122014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAMapCoordsToPeriodicDomain(DM,PetscScalar*,PetscScalar*);
123f0aa4865SJed Brown 
124014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGetReducedDA(DM,PetscInt,DM*);
1253c48a1e8SJed Brown 
126014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDASetFieldName(DM,PetscInt,const char[]);
127014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGetFieldName(DM,PetscInt,const char**);
1283c48a1e8SJed Brown 
129014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDASetBoundaryType(DM,DMDABoundaryType,DMDABoundaryType,DMDABoundaryType);
130014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDASetDof(DM, PetscInt);
13188661749SPeter Brune PETSC_EXTERN PetscErrorCode DMDASetOverlap(DM, PetscInt);
132014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDASetStencilWidth(DM, PetscInt);
133014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDASetOwnershipRanges(DM,const PetscInt[],const PetscInt[],const PetscInt[]);
134014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGetOwnershipRanges(DM,const PetscInt**,const PetscInt**,const PetscInt**);
135014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDASetNumProcs(DM, PetscInt, PetscInt, PetscInt);
136014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDASetStencilType(DM, DMDAStencilType);
1373c48a1e8SJed Brown 
138014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAVecGetArray(DM,Vec,void *);
139014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAVecRestoreArray(DM,Vec,void *);
1403c48a1e8SJed Brown 
141014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAVecGetArrayDOF(DM,Vec,void *);
142014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAVecRestoreArrayDOF(DM,Vec,void *);
1433c48a1e8SJed Brown 
144014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDASplitComm2d(MPI_Comm,PetscInt,PetscInt,PetscInt,MPI_Comm*);
1453c48a1e8SJed Brown 
1463c48a1e8SJed Brown /*S
1473c48a1e8SJed Brown      DMDALocalInfo - C struct that contains information about a structured grid and a processors logical
1483c48a1e8SJed Brown               location in it.
1493c48a1e8SJed Brown 
1503c48a1e8SJed Brown    Level: beginner
1513c48a1e8SJed Brown 
1523c48a1e8SJed Brown   Concepts: distributed array
1533c48a1e8SJed Brown 
1543c48a1e8SJed Brown   Developer note: Then entries in this struct are int instead of PetscInt so that the elements may
1553c48a1e8SJed Brown                   be extracted in Fortran as if from an integer array
1563c48a1e8SJed Brown 
1573c48a1e8SJed Brown .seealso:  DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DM, DMDAGetLocalInfo(), DMDAGetInfo()
1583c48a1e8SJed Brown S*/
1593c48a1e8SJed Brown typedef struct {
1603c48a1e8SJed Brown   PetscInt         dim,dof,sw;
1613c48a1e8SJed Brown   PetscInt         mx,my,mz;    /* global number of grid points in each direction */
1626ee1aa1aSDmitry Karpeev   PetscInt         xs,ys,zs;    /* starting point of this processor, excluding ghosts */
1633c48a1e8SJed Brown   PetscInt         xm,ym,zm;    /* number of grid points on this processor, excluding ghosts */
1643c48a1e8SJed Brown   PetscInt         gxs,gys,gzs;    /* starting point of this processor including ghosts */
1653c48a1e8SJed Brown   PetscInt         gxm,gym,gzm;    /* number of grid points on this processor including ghosts */
1663c48a1e8SJed Brown   DMDABoundaryType bx,by,bz; /* type of ghost nodes at boundary */
1673c48a1e8SJed Brown   DMDAStencilType  st;
1683c48a1e8SJed Brown   DM               da;
1693c48a1e8SJed Brown } DMDALocalInfo;
1703c48a1e8SJed Brown 
1713c48a1e8SJed Brown /*MC
1723c48a1e8SJed Brown       DMDAForEachPointBegin2d - Starts a loop over the local part of a two dimensional DMDA
1733c48a1e8SJed Brown 
1743c48a1e8SJed Brown    Synopsis:
1753c48a1e8SJed Brown    void  DMDAForEachPointBegin2d(DALocalInfo *info,PetscInt i,PetscInt j);
1763c48a1e8SJed Brown 
1773c48a1e8SJed Brown    Not Collective
1783c48a1e8SJed Brown 
1793c48a1e8SJed Brown    Level: intermediate
1803c48a1e8SJed Brown 
1813c48a1e8SJed Brown .seealso: DMDAForEachPointEnd2d(), DMDAVecGetArray()
1823c48a1e8SJed Brown M*/
1833c48a1e8SJed Brown #define DMDAForEachPointBegin2d(info,i,j) {\
1843c48a1e8SJed Brown   PetscInt _xints = info->xs,_xinte = info->xs+info->xm,_yints = info->ys,_yinte = info->ys+info->ym;\
1853c48a1e8SJed Brown   for (j=_yints; j<_yinte; j++) {\
1863c48a1e8SJed Brown     for (i=_xints; i<_xinte; i++) {\
1873c48a1e8SJed Brown 
1883c48a1e8SJed Brown /*MC
1893c48a1e8SJed Brown       DMDAForEachPointEnd2d - Ends a loop over the local part of a two dimensional DMDA
1903c48a1e8SJed Brown 
1913c48a1e8SJed Brown    Synopsis:
1923c48a1e8SJed Brown    void  DMDAForEachPointEnd2d;
1933c48a1e8SJed Brown 
1943c48a1e8SJed Brown    Not Collective
1953c48a1e8SJed Brown 
1963c48a1e8SJed Brown    Level: intermediate
1973c48a1e8SJed Brown 
1983c48a1e8SJed Brown .seealso: DMDAForEachPointBegin2d(), DMDAVecGetArray()
1993c48a1e8SJed Brown M*/
2003c48a1e8SJed Brown #define DMDAForEachPointEnd2d }}}
2013c48a1e8SJed Brown 
2023c48a1e8SJed Brown /*MC
2033c48a1e8SJed Brown       DMDACoor2d - Structure for holding 2d (x and y) coordinates.
2043c48a1e8SJed Brown 
2053c48a1e8SJed Brown     Level: intermediate
2063c48a1e8SJed Brown 
2073c48a1e8SJed Brown     Sample Usage:
2083c48a1e8SJed Brown       DMDACoor2d **coors;
2093c48a1e8SJed Brown       Vec      vcoors;
2103c48a1e8SJed Brown       DM       cda;
2113c48a1e8SJed Brown 
212*2150357eSBarry Smith       DMGetCoordinates(da,&vcoors);
2133c48a1e8SJed Brown       DMDAGetCoordinateDA(da,&cda);
2143c48a1e8SJed Brown       DMDAVecGetArray(cda,vcoors,&coors);
2153c48a1e8SJed Brown       DMDAGetCorners(cda,&mstart,&nstart,0,&m,&n,0)
2163c48a1e8SJed Brown       for (i=mstart; i<mstart+m; i++) {
2173c48a1e8SJed Brown         for (j=nstart; j<nstart+n; j++) {
2183c48a1e8SJed Brown           x = coors[j][i].x;
2193c48a1e8SJed Brown           y = coors[j][i].y;
2203c48a1e8SJed Brown           ......
2213c48a1e8SJed Brown         }
2223c48a1e8SJed Brown       }
2233c48a1e8SJed Brown       DMDAVecRestoreArray(dac,vcoors,&coors);
2243c48a1e8SJed Brown 
225*2150357eSBarry Smith .seealso: DMDACoor3d, DMDAForEachPointBegin(), DMDAGetCoordinateDA(), DMGetCoordinates(), DMDAGetGhostCoordinates()
2263c48a1e8SJed Brown M*/
2273c48a1e8SJed Brown typedef struct {PetscScalar x,y;} DMDACoor2d;
2283c48a1e8SJed Brown 
2293c48a1e8SJed Brown /*MC
2303c48a1e8SJed Brown       DMDACoor3d - Structure for holding 3d (x, y and z) coordinates.
2313c48a1e8SJed Brown 
2323c48a1e8SJed Brown     Level: intermediate
2333c48a1e8SJed Brown 
2343c48a1e8SJed Brown     Sample Usage:
2353c48a1e8SJed Brown       DMDACoor3d ***coors;
2363c48a1e8SJed Brown       Vec      vcoors;
2373c48a1e8SJed Brown       DM       cda;
2383c48a1e8SJed Brown 
239*2150357eSBarry Smith       DMGetCoordinates(da,&vcoors);
2403c48a1e8SJed Brown       DMDAGetCoordinateDA(da,&cda);
2413c48a1e8SJed Brown       DMDAVecGetArray(cda,vcoors,&coors);
2423c48a1e8SJed Brown       DMDAGetCorners(cda,&mstart,&nstart,&pstart,&m,&n,&p)
2433c48a1e8SJed Brown       for (i=mstart; i<mstart+m; i++) {
2443c48a1e8SJed Brown         for (j=nstart; j<nstart+n; j++) {
2453c48a1e8SJed Brown           for (k=pstart; k<pstart+p; k++) {
2463c48a1e8SJed Brown             x = coors[k][j][i].x;
2473c48a1e8SJed Brown             y = coors[k][j][i].y;
2483c48a1e8SJed Brown             z = coors[k][j][i].z;
2493c48a1e8SJed Brown           ......
2503c48a1e8SJed Brown         }
2513c48a1e8SJed Brown       }
2523c48a1e8SJed Brown       DMDAVecRestoreArray(dac,vcoors,&coors);
2533c48a1e8SJed Brown 
254*2150357eSBarry Smith .seealso: DMDACoor2d, DMDAForEachPointBegin(), DMDAGetCoordinateDA(), DMGetCoordinates(), DMDAGetGhostCoordinates()
2553c48a1e8SJed Brown M*/
2563c48a1e8SJed Brown typedef struct {PetscScalar x,y,z;} DMDACoor3d;
2573c48a1e8SJed Brown 
258014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGetLocalInfo(DM,DMDALocalInfo*);
2597aaf2d21SPeter Brune PETSC_EXTERN PetscErrorCode DMDAGetLocalBlockInfo(DM,DMDALocalInfo*);
260638cfed1SJed Brown PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDALocalFunction1)(DMDALocalInfo*,void*,void*,void*);
261014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAComputeFunctionLocal(DM, DMDALocalFunction1, Vec, Vec, void *);
262014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAComputeFunctionLocalGhost(DM, DMDALocalFunction1, Vec, Vec, void *);
263014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAFormJacobianLocal(DM, DMDALocalFunction1, Vec, Mat, void *);
264014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAComputeFunction1(DM,Vec,Vec,void*);
265014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAComputeFunction(DM,PetscErrorCode (*)(void),Vec,Vec,void*);
266014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAComputeFunctioni1(DM,PetscInt,Vec,PetscScalar*,void*);
267014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAComputeFunctionib1(DM,PetscInt,Vec,PetscScalar*,void*);
268014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAComputeJacobian1WithAdic(DM,Vec,Mat,void*);
269014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAComputeJacobian1WithAdifor(DM,Vec,Mat,void*);
270014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAMultiplyByJacobian1WithAdic(DM,Vec,Vec,Vec,void*);
271014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAMultiplyByJacobian1WithAdifor(DM,Vec,Vec,Vec,void*);
272014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAMultiplyByJacobian1WithAD(DM,Vec,Vec,Vec,void*);
273014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAComputeJacobian1(DM,Vec,Mat,void*);
274014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGetLocalFunction(DM,DMDALocalFunction1*);
275014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDASetLocalFunction(DM,DMDALocalFunction1);
276014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDASetLocalFunctioni(DM,PetscErrorCode (*)(DMDALocalInfo*,MatStencil*,void*,PetscScalar*,void*));
277014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDASetLocalFunctionib(DM,PetscErrorCode (*)(DMDALocalInfo*,MatStencil*,void*,PetscScalar*,void*));
278014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGetLocalJacobian(DM,DMDALocalFunction1*);
279014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDASetLocalJacobian(DM,DMDALocalFunction1);
280014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDASetLocalAdicFunction_Private(DM,DMDALocalFunction1);
2813c48a1e8SJed Brown 
282014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetDM(Mat,DM);
283014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRegisterDAAD(void);
284014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateDAAD(DM,Mat*);
285014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqUSFFT(Vec,DM,Mat*);
2863c48a1e8SJed Brown 
2873c48a1e8SJed Brown /*MC
2883c48a1e8SJed Brown        DMDASetLocalAdicFunction - Caches in a DM a local function computed by ADIC/ADIFOR
2893c48a1e8SJed Brown 
2903c48a1e8SJed Brown    Synopsis:
2913c48a1e8SJed Brown    PetscErrorCode DMDASetLocalAdicFunction(DM da,DMDALocalFunction1 ad_lf)
2923c48a1e8SJed Brown 
2933c48a1e8SJed Brown    Logically Collective on DM
2943c48a1e8SJed Brown 
2953c48a1e8SJed Brown    Input Parameter:
2963c48a1e8SJed Brown +  da - initial distributed array
2973c48a1e8SJed Brown -  ad_lf - the local function as computed by ADIC/ADIFOR
2983c48a1e8SJed Brown 
2993c48a1e8SJed Brown    Level: intermediate
3003c48a1e8SJed Brown 
3013c48a1e8SJed Brown .keywords:  distributed array, refine
3023c48a1e8SJed Brown 
3033c48a1e8SJed Brown .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDAGetLocalFunction(), DMDASetLocalFunction(),
3043c48a1e8SJed Brown           DMDASetLocalJacobian()
3053c48a1e8SJed Brown M*/
3063c48a1e8SJed Brown #if defined(PETSC_HAVE_ADIC)
3073c48a1e8SJed Brown #  define DMDASetLocalAdicFunction(a,d) DMDASetLocalAdicFunction_Private(a,(DMDALocalFunction1)d)
3083c48a1e8SJed Brown #else
3093c48a1e8SJed Brown #  define DMDASetLocalAdicFunction(a,d) DMDASetLocalAdicFunction_Private(a,0)
3103c48a1e8SJed Brown #endif
3113c48a1e8SJed Brown 
312014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDASetLocalAdicMFFunction_Private(DM,DMDALocalFunction1);
3133c48a1e8SJed Brown #if defined(PETSC_HAVE_ADIC)
3143c48a1e8SJed Brown #  define DMDASetLocalAdicMFFunction(a,d) DMDASetLocalAdicMFFunction_Private(a,(DMDALocalFunction1)d)
3153c48a1e8SJed Brown #else
3163c48a1e8SJed Brown #  define DMDASetLocalAdicMFFunction(a,d) DMDASetLocalAdicMFFunction_Private(a,0)
3173c48a1e8SJed Brown #endif
318014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDASetLocalAdicFunctioni_Private(DM,PetscErrorCode (*)(DMDALocalInfo*,MatStencil*,void*,void*,void*));
3193c48a1e8SJed Brown #if defined(PETSC_HAVE_ADIC)
3203c48a1e8SJed Brown #  define DMDASetLocalAdicFunctioni(a,d) DMDASetLocalAdicFunctioni_Private(a,(PetscErrorCode (*)(DMDALocalInfo*,MatStencil*,void*,void*,void*))d)
3213c48a1e8SJed Brown #else
3223c48a1e8SJed Brown #  define DMDASetLocalAdicFunctioni(a,d) DMDASetLocalAdicFunctioni_Private(a,0)
3233c48a1e8SJed Brown #endif
324014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDASetLocalAdicMFFunctioni_Private(DM,PetscErrorCode (*)(DMDALocalInfo*,MatStencil*,void*,void*,void*));
3253c48a1e8SJed Brown #if defined(PETSC_HAVE_ADIC)
3263c48a1e8SJed Brown #  define DMDASetLocalAdicMFFunctioni(a,d) DMDASetLocalAdicMFFunctioni_Private(a,(PetscErrorCode (*)(DMDALocalInfo*,MatStencil*,void*,void*,void*))d)
3273c48a1e8SJed Brown #else
3283c48a1e8SJed Brown #  define DMDASetLocalAdicMFFunctioni(a,d) DMDASetLocalAdicMFFunctioni_Private(a,0)
3293c48a1e8SJed Brown #endif
3303c48a1e8SJed Brown 
331014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDASetLocalAdicFunctionib_Private(DM,PetscErrorCode (*)(DMDALocalInfo*,MatStencil*,void*,void*,void*));
3323c48a1e8SJed Brown #if defined(PETSC_HAVE_ADIC)
3333c48a1e8SJed Brown #  define DMDASetLocalAdicFunctionib(a,d) DMDASetLocalAdicFunctionib_Private(a,(PetscErrorCode (*)(DMDALocalInfo*,MatStencil*,void*,void*,void*))d)
3343c48a1e8SJed Brown #else
3353c48a1e8SJed Brown #  define DMDASetLocalAdicFunctionib(a,d) DMDASetLocalAdicFunctionib_Private(a,0)
3363c48a1e8SJed Brown #endif
337014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDASetLocalAdicMFFunctionib_Private(DM,PetscErrorCode (*)(DMDALocalInfo*,MatStencil*,void*,void*,void*));
3383c48a1e8SJed Brown #if defined(PETSC_HAVE_ADIC)
3393c48a1e8SJed Brown #  define DMDASetLocalAdicMFFunctionib(a,d) DMDASetLocalAdicMFFunctionib_Private(a,(PetscErrorCode (*)(DMDALocalInfo*,MatStencil*,void*,void*,void*))d)
3403c48a1e8SJed Brown #else
3413c48a1e8SJed Brown #  define DMDASetLocalAdicMFFunctionib(a,d) DMDASetLocalAdicMFFunctionib_Private(a,0)
3423c48a1e8SJed Brown #endif
3433c48a1e8SJed Brown 
344014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAComputeFunctioniTest1(DM,void*);
34519fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode DMDASetGetMatrix(DM,PetscErrorCode (*)(DM, MatType,Mat *));
346014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDASetBlockFills(DM,PetscInt*,PetscInt*);
347014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDASetRefinementFactor(DM,PetscInt,PetscInt,PetscInt);
348014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGetRefinementFactor(DM,PetscInt*,PetscInt*,PetscInt*);
3493c48a1e8SJed Brown 
350014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGetAdicArray(DM,PetscBool ,void*,void*,PetscInt*);
351014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDARestoreAdicArray(DM,PetscBool ,void*,void*,PetscInt*);
352014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGetAdicMFArray(DM,PetscBool ,void*,void*,PetscInt*);
353014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGetAdicMFArray4(DM,PetscBool ,void*,void*,PetscInt*);
354014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGetAdicMFArray9(DM,PetscBool ,void*,void*,PetscInt*);
355014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGetAdicMFArrayb(DM,PetscBool ,void*,void*,PetscInt*);
356014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDARestoreAdicMFArray(DM,PetscBool ,void*,void*,PetscInt*);
357014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGetArray(DM,PetscBool ,void*);
358014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDARestoreArray(DM,PetscBool ,void*);
359014dd563SJed Brown PETSC_EXTERN PetscErrorCode ad_DAGetArray(DM,PetscBool ,void*);
360014dd563SJed Brown PETSC_EXTERN PetscErrorCode ad_DARestoreArray(DM,PetscBool ,void*);
361014dd563SJed Brown PETSC_EXTERN PetscErrorCode admf_DAGetArray(DM,PetscBool ,void*);
362014dd563SJed Brown PETSC_EXTERN PetscErrorCode admf_DARestoreArray(DM,PetscBool ,void*);
3633c48a1e8SJed Brown 
364014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDACreatePF(DM,PF*);
3653c48a1e8SJed Brown 
36657459e9aSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMDAGetNumCells(DM, PetscInt *);
36757459e9aSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMDAGetNumVertices(DM, PetscInt *, PetscInt *, PetscInt *, PetscInt *);
36857459e9aSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMDAGetNumFaces(DM, PetscInt *, PetscInt *, PetscInt *, PetscInt *, PetscInt *, PetscInt *);
36957459e9aSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMDAGetHeightStratum(DM, PetscInt, PetscInt *, PetscInt *);
370014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDACreateSection(DM, PetscInt[], PetscInt[], PetscInt[], PetscInt[]);
37157459e9aSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMDAComputeCellGeometry(DM, PetscInt, PetscQuadrature *, PetscReal [], PetscReal [], PetscReal [], PetscReal []);
37257459e9aSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMDAVecGetClosure(DM, PetscSection, Vec, PetscInt, const PetscScalar **);
37357459e9aSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMDAVecSetClosure(DM, PetscSection, Vec, PetscInt, const PetscScalar *, InsertMode);
374aa1993deSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMDAGetClosure(DM,PetscSection,PetscInt,PetscInt*,const PetscInt**);
375aa1993deSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMDARestoreClosure(DM,PetscSection,PetscInt,PetscInt*,const PetscInt**);
376aa1993deSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMDAGetClosureScalar(DM,PetscSection,PetscInt,PetscScalar*,const PetscScalar**);
377aa1993deSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMDARestoreClosureScalar(DM,PetscSection,PetscInt,PetscScalar*,const PetscScalar**);
378aa1993deSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMDASetClosureScalar(DM,PetscSection,PetscInt,PetscScalar*,const PetscScalar*,InsertMode);
379f0e3914cSSatish Balay PETSC_EXTERN PetscErrorCode DMDAConvertToCell(DM, MatStencil, PetscInt *);
38080800b1aSMatthew G Knepley 
381b859378eSBarry Smith #define DMDA_FILE_CLASSID 1211220
3823c48a1e8SJed Brown #endif
383