xref: /petsc/include/petscdmda.h (revision 1957e957c9edf41e4162784af8cf55e6b64ba41e)
13c48a1e8SJed Brown #if !defined(__PETSCDMDA_H)
23c48a1e8SJed Brown #define __PETSCDMDA_H
33c48a1e8SJed Brown 
42c8e378dSBarry Smith #include <petscdm.h>
5a0845e3aSMatthew G. Knepley #include <petscdt.h>
6df66330eSJed Brown #include <petscfetypes.h>
71e25c274SJed Brown #include <petscdmdatypes.h>
82c8e378dSBarry Smith #include <petscpf.h>
92c8e378dSBarry Smith #include <petscao.h>
10bb619f44SMatthew G. Knepley #include <petscfe.h>
113c48a1e8SJed Brown 
123c48a1e8SJed Brown /*MC
133c48a1e8SJed Brown      DMDA_STENCIL_STAR - "Star"-type stencil. In logical grid coordinates, only (i,j,k), (i+s,j,k), (i,j+s,k),
143c48a1e8SJed Brown                        (i,j,k+s) are in the stencil  NOT, for example, (i+s,j+s,k)
153c48a1e8SJed Brown 
163c48a1e8SJed Brown      Level: beginner
173c48a1e8SJed Brown 
18*1957e957SBarry Smith      Determines what ghost point values are brought over to each process; in this case the "corner" values are not
19*1957e957SBarry Smith      brought over and hence should not be accessed locally
20*1957e957SBarry Smith 
21db05f41bSBarry Smith .seealso: 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 
30db05f41bSBarry Smith .seealso: DMDA_STENCIL_STAR, DMDAStencilType, DMDASetStencilType()
313c48a1e8SJed Brown M*/
323c48a1e8SJed Brown 
33014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDASetInterpolationType(DM,DMDAInterpolationType);
34014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGetInterpolationType(DM,DMDAInterpolationType*);
353c48a1e8SJed Brown 
36014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDASetElementType(DM,DMDAElementType);
37014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGetElementType(DM,DMDAElementType*);
38014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGetElements(DM,PetscInt *,PetscInt *,const PetscInt*[]);
39014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDARestoreElements(DM,PetscInt *,PetscInt *,const PetscInt*[]);
403c48a1e8SJed Brown 
413c48a1e8SJed Brown typedef enum { DMDA_X,DMDA_Y,DMDA_Z } DMDADirection;
423c48a1e8SJed Brown 
433c48a1e8SJed Brown #define MATSEQUSFFT        "sequsfft"
443c48a1e8SJed Brown 
45014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDACreate(MPI_Comm,DM*);
46014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDASetDim(DM,PetscInt);
47014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDASetSizes(DM,PetscInt,PetscInt,PetscInt);
48bff4a2f0SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMDACreate1d(MPI_Comm,DMBoundaryType,PetscInt,PetscInt,PetscInt,const PetscInt[],DM *);
49bff4a2f0SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMDACreate2d(MPI_Comm,DMBoundaryType,DMBoundaryType,DMDAStencilType,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],DM*);
50bff4a2f0SMatthew 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*);
513c48a1e8SJed Brown 
52014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGlobalToNaturalBegin(DM,Vec,InsertMode,Vec);
53014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGlobalToNaturalEnd(DM,Vec,InsertMode,Vec);
54014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDANaturalToGlobalBegin(DM,Vec,InsertMode,Vec);
55014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDANaturalToGlobalEnd(DM,Vec,InsertMode,Vec);
56d78e899eSRichard 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);}
57d78e899eSRichard 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);}
58014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDACreateNaturalVector(DM,Vec *);
593c48a1e8SJed Brown 
60014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGetCorners(DM,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
61014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGetGhostCorners(DM,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
62bff4a2f0SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMDAGetInfo(DM,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,DMBoundaryType*,DMBoundaryType*,DMBoundaryType*,DMDAStencilType*);
63014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGetProcessorSubset(DM,DMDADirection,PetscInt,MPI_Comm*);
64014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGetProcessorSubsets(DM,DMDADirection,MPI_Comm*);
65db05f41bSBarry Smith PETSC_EXTERN PetscErrorCode DMDAGetRay(DM,DMDADirection,PetscInt,Vec*,VecScatter*);
663c48a1e8SJed Brown 
67014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGlobalToNaturalAllCreate(DM,VecScatter*);
68014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDANaturalAllToGlobalCreate(DM,VecScatter*);
693c48a1e8SJed Brown 
708ea3bf28SBarry Smith PETSC_EXTERN PetscErrorCode DMDAGetGlobalIndices(DM,PetscInt*,const PetscInt*[]);
718ea3bf28SBarry Smith PETSC_EXTERN PetscErrorCode DMDARestoreGlobalIndices(DM,PetscInt*,const PetscInt*[]);
723c48a1e8SJed Brown 
73014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGetScatter(DM,VecScatter*,VecScatter*,VecScatter*);
74014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGetNeighbors(DM,const PetscMPIInt**);
753c48a1e8SJed Brown 
76014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGetAO(DM,AO*);
77014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDASetUniformCoordinates(DM,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal);
78014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGetBoundingBox(DM,PetscReal[],PetscReal[]);
79014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGetLocalBoundingBox(DM,PetscReal[],PetscReal[]);
80a74ba6f7SBarry Smith PETSC_EXTERN PetscErrorCode DMDAGetLogicalCoordinate(DM,PetscScalar,PetscScalar,PetscScalar,PetscInt*,PetscInt*,PetscInt*,PetscScalar*,PetscScalar*,PetscScalar*);
81f0aa4865SJed Brown /* function to wrap coordinates around boundary */
82014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAMapCoordsToPeriodicDomain(DM,PetscScalar*,PetscScalar*);
83f0aa4865SJed Brown 
84db05f41bSBarry Smith PETSC_EXTERN PetscErrorCode DMDAGetReducedDMDA(DM,PetscInt,DM*);
853c48a1e8SJed Brown 
86014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDASetFieldName(DM,PetscInt,const char[]);
87014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGetFieldName(DM,PetscInt,const char**);
88109c9344SBarry Smith PETSC_EXTERN PetscErrorCode DMDASetCoordinateName(DM,PetscInt,const char[]);
89109c9344SBarry Smith PETSC_EXTERN PetscErrorCode DMDAGetCoordinateName(DM,PetscInt,const char**);
903c48a1e8SJed Brown 
91bff4a2f0SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMDASetBoundaryType(DM,DMBoundaryType,DMBoundaryType,DMBoundaryType);
92014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDASetDof(DM, PetscInt);
937ddda789SPeter Brune PETSC_EXTERN PetscErrorCode DMDASetOverlap(DM,PetscInt,PetscInt,PetscInt);
947ddda789SPeter Brune PETSC_EXTERN PetscErrorCode DMDAGetOverlap(DM,PetscInt*,PetscInt*,PetscInt*);
953e7870d2SPeter Brune PETSC_EXTERN PetscErrorCode DMDASetNumLocalSubDomains(DM,PetscInt);
963e7870d2SPeter Brune PETSC_EXTERN PetscErrorCode DMDAGetNumLocalSubDomains(DM,PetscInt*);
9795c13181SPeter Brune PETSC_EXTERN PetscErrorCode DMDAGetOffset(DM,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
9895c13181SPeter Brune PETSC_EXTERN PetscErrorCode DMDASetOffset(DM,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt);
9940234c92SPeter Brune PETSC_EXTERN PetscErrorCode DMDAGetNonOverlappingRegion(DM,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
10040234c92SPeter Brune PETSC_EXTERN PetscErrorCode DMDASetNonOverlappingRegion(DM,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt);
101014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDASetStencilWidth(DM, PetscInt);
102014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDASetOwnershipRanges(DM,const PetscInt[],const PetscInt[],const PetscInt[]);
103014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGetOwnershipRanges(DM,const PetscInt**,const PetscInt**,const PetscInt**);
104014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDASetNumProcs(DM, PetscInt, PetscInt, PetscInt);
105014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDASetStencilType(DM, DMDAStencilType);
1063c48a1e8SJed Brown 
107014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAVecGetArray(DM,Vec,void *);
108014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAVecRestoreArray(DM,Vec,void *);
1093c48a1e8SJed Brown 
110014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAVecGetArrayDOF(DM,Vec,void *);
111014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAVecRestoreArrayDOF(DM,Vec,void *);
1123c48a1e8SJed Brown 
113014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDASplitComm2d(MPI_Comm,PetscInt,PetscInt,PetscInt,MPI_Comm*);
1143c48a1e8SJed Brown 
115ff9846d9SPeter Brune PETSC_EXTERN PetscErrorCode DMDACreatePatchIS(DM,MatStencil*,MatStencil*,IS*);
116ff9846d9SPeter Brune 
1173c48a1e8SJed Brown 
1183c48a1e8SJed Brown /*MC
1193c48a1e8SJed Brown       DMDACoor2d - Structure for holding 2d (x and y) coordinates.
1203c48a1e8SJed Brown 
1213c48a1e8SJed Brown     Level: intermediate
1223c48a1e8SJed Brown 
1233c48a1e8SJed Brown     Sample Usage:
1243c48a1e8SJed Brown       DMDACoor2d **coors;
1253c48a1e8SJed Brown       Vec      vcoors;
1263c48a1e8SJed Brown       DM       cda;
1273c48a1e8SJed Brown 
1282150357eSBarry Smith       DMGetCoordinates(da,&vcoors);
129*1957e957SBarry Smith       DMGetCoordinateDM(da,&cda);
1303c48a1e8SJed Brown       DMDAVecGetArray(cda,vcoors,&coors);
1313c48a1e8SJed Brown       DMDAGetCorners(cda,&mstart,&nstart,0,&m,&n,0)
1323c48a1e8SJed Brown       for (i=mstart; i<mstart+m; i++) {
1333c48a1e8SJed Brown         for (j=nstart; j<nstart+n; j++) {
1343c48a1e8SJed Brown           x = coors[j][i].x;
1353c48a1e8SJed Brown           y = coors[j][i].y;
1363c48a1e8SJed Brown           ......
1373c48a1e8SJed Brown         }
1383c48a1e8SJed Brown       }
1393c48a1e8SJed Brown       DMDAVecRestoreArray(dac,vcoors,&coors);
1403c48a1e8SJed Brown 
141*1957e957SBarry Smith .seealso: DMDACoor3d, DMGetCoordinateDM(), DMGetCoordinates(), DMDAGetGhostCoordinates()
1423c48a1e8SJed Brown M*/
1433c48a1e8SJed Brown typedef struct {PetscScalar x,y;} DMDACoor2d;
1443c48a1e8SJed Brown 
1453c48a1e8SJed Brown /*MC
1463c48a1e8SJed Brown       DMDACoor3d - Structure for holding 3d (x, y and z) coordinates.
1473c48a1e8SJed Brown 
1483c48a1e8SJed Brown     Level: intermediate
1493c48a1e8SJed Brown 
1503c48a1e8SJed Brown     Sample Usage:
1513c48a1e8SJed Brown       DMDACoor3d ***coors;
1523c48a1e8SJed Brown       Vec      vcoors;
1533c48a1e8SJed Brown       DM       cda;
1543c48a1e8SJed Brown 
1552150357eSBarry Smith       DMGetCoordinates(da,&vcoors);
156*1957e957SBarry Smith       DMGetCoordinateDM(da,&cda);
1573c48a1e8SJed Brown       DMDAVecGetArray(cda,vcoors,&coors);
1583c48a1e8SJed Brown       DMDAGetCorners(cda,&mstart,&nstart,&pstart,&m,&n,&p)
1593c48a1e8SJed Brown       for (i=mstart; i<mstart+m; i++) {
1603c48a1e8SJed Brown         for (j=nstart; j<nstart+n; j++) {
1613c48a1e8SJed Brown           for (k=pstart; k<pstart+p; k++) {
1623c48a1e8SJed Brown             x = coors[k][j][i].x;
1633c48a1e8SJed Brown             y = coors[k][j][i].y;
1643c48a1e8SJed Brown             z = coors[k][j][i].z;
1653c48a1e8SJed Brown           ......
1663c48a1e8SJed Brown         }
1673c48a1e8SJed Brown       }
1683c48a1e8SJed Brown       DMDAVecRestoreArray(dac,vcoors,&coors);
1693c48a1e8SJed Brown 
170*1957e957SBarry Smith .seealso: DMDACoor2d, DMGetCoordinateDM(), DMGetCoordinates(), DMDAGetGhostCoordinates()
1713c48a1e8SJed Brown M*/
1723c48a1e8SJed Brown typedef struct {PetscScalar x,y,z;} DMDACoor3d;
1733c48a1e8SJed Brown 
174014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGetLocalInfo(DM,DMDALocalInfo*);
1753c48a1e8SJed Brown 
176014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRegisterDAAD(void);
177014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateDAAD(DM,Mat*);
178014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqUSFFT(Vec,DM,Mat*);
1793c48a1e8SJed Brown 
180b412c318SBarry Smith PETSC_EXTERN PetscErrorCode DMDASetGetMatrix(DM,PetscErrorCode (*)(DM, Mat *));
181ce308e1dSBarry Smith PETSC_EXTERN PetscErrorCode DMDASetBlockFills(DM,const PetscInt*,const PetscInt*);
182014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDASetRefinementFactor(DM,PetscInt,PetscInt,PetscInt);
183014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGetRefinementFactor(DM,PetscInt*,PetscInt*,PetscInt*);
1843c48a1e8SJed Brown 
185014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGetArray(DM,PetscBool ,void*);
186014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDARestoreArray(DM,PetscBool ,void*);
1873c48a1e8SJed Brown 
188014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDACreatePF(DM,PF*);
1893c48a1e8SJed Brown 
1903582350cSMatthew G. Knepley /* Emulation of DMPlex */
1913582350cSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMDAGetNumCells(DM, PetscInt *, PetscInt *, PetscInt *, PetscInt *);
19248710be7SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMDAGetCellPoint(DM, PetscInt, PetscInt, PetscInt, PetscInt *);
19357459e9aSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMDAGetNumVertices(DM, PetscInt *, PetscInt *, PetscInt *, PetscInt *);
19457459e9aSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMDAGetNumFaces(DM, PetscInt *, PetscInt *, PetscInt *, PetscInt *, PetscInt *, PetscInt *);
19557459e9aSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMDAGetHeightStratum(DM, PetscInt, PetscInt *, PetscInt *);
196affa55c7SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMDACreateSection(DM, PetscInt[], PetscInt[], PetscInt[], PetscSection *);
19757459e9aSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMDAComputeCellGeometry(DM, PetscInt, PetscQuadrature *, PetscReal [], PetscReal [], PetscReal [], PetscReal []);
1986693a731SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMDAGetTransitiveClosure(DM, PetscInt, PetscBool, PetscInt *, PetscInt **);
1996693a731SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMDARestoreTransitiveClosure(DM, PetscInt, PetscBool, PetscInt *, PetscInt **);
200d7a12f93SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMDAVecGetClosure(DM, PetscSection, Vec, PetscInt, PetscInt *, PetscScalar **);
201d7a12f93SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMDAVecRestoreClosure(DM, PetscSection, Vec, PetscInt, PetscInt *, PetscScalar **);
20257459e9aSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMDAVecSetClosure(DM, PetscSection, Vec, PetscInt, const PetscScalar *, InsertMode);
203aa1993deSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMDAGetClosure(DM, PetscSection, PetscInt, PetscInt*, const PetscInt**);
204aa1993deSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMDARestoreClosure(DM, PetscSection, PetscInt, PetscInt*, const PetscInt**);
205d7a12f93SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMDAGetClosureScalar(DM, PetscSection, PetscInt, PetscScalar*, PetscInt*, PetscScalar**);
206d7a12f93SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMDARestoreClosureScalar(DM, PetscSection, PetscInt, PetscScalar*, PetscInt*, PetscScalar**);
207aa1993deSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMDASetClosureScalar(DM,PetscSection,PetscInt,PetscScalar*,const PetscScalar*,InsertMode);
208f0e3914cSSatish Balay PETSC_EXTERN PetscErrorCode DMDAConvertToCell(DM, MatStencil, PetscInt *);
2093385ff7eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMDASetVertexCoordinates(DM,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal);
2108ca2f13cSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMDASetPreallocationCenterDimension(DM, PetscInt);
2118ca2f13cSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMDAGetPreallocationCenterDimension(DM, PetscInt*);
21280800b1aSMatthew G Knepley 
213c110b1eeSGeoffrey Irving PETSC_EXTERN PetscErrorCode DMDAProjectFunction(DM, PetscFE[], void (**)(const PetscReal [], PetscScalar *, void *), void **, InsertMode, Vec);
214c110b1eeSGeoffrey Irving PETSC_EXTERN PetscErrorCode DMDAProjectFunctionLocal(DM, PetscFE[], void (**)(const PetscReal [], PetscScalar *, void *), void **, InsertMode, Vec);
215c110b1eeSGeoffrey Irving PETSC_EXTERN PetscErrorCode DMDAComputeL2Diff(DM, PetscFE[], void (**)(const PetscReal [], PetscScalar *, void *), void **, Vec, PetscReal *);
2169a800dd8SMatthew G. Knepley 
2173c48a1e8SJed Brown #endif
218