xref: /petsc/include/petscdmda.h (revision 9fbee5477fd88ea4536ebb185f3c80da15fb55c0)
126bd1501SBarry Smith #if !defined(PETSCDMDA_H)
226bd1501SBarry Smith #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*);
3397779f9aSLisandro Dalcin PETSC_EXTERN PetscErrorCode DMDACreateAggregates(DM,DM,Mat*);
343c48a1e8SJed Brown 
35d4a6ed37SStefano Zampini /* FEM */
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*[]);
40d4a6ed37SStefano Zampini PETSC_EXTERN PetscErrorCode DMDAGetElementsSizes(DM,PetscInt*,PetscInt*,PetscInt*);
41d4a6ed37SStefano Zampini PETSC_EXTERN PetscErrorCode DMDAGetElementsCorners(DM,PetscInt*,PetscInt*,PetscInt*);
42d4a6ed37SStefano Zampini PETSC_EXTERN PetscErrorCode DMDAGetSubdomainCornersIS(DM,IS*);
43d4a6ed37SStefano Zampini PETSC_EXTERN PetscErrorCode DMDARestoreSubdomainCornersIS(DM,IS*);
443c48a1e8SJed Brown 
453c48a1e8SJed Brown #define MATSEQUSFFT        "sequsfft"
463c48a1e8SJed Brown 
47014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDACreate(MPI_Comm,DM*);
48014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDASetSizes(DM,PetscInt,PetscInt,PetscInt);
49bff4a2f0SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMDACreate1d(MPI_Comm,DMBoundaryType,PetscInt,PetscInt,PetscInt,const PetscInt[],DM *);
50bff4a2f0SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMDACreate2d(MPI_Comm,DMBoundaryType,DMBoundaryType,DMDAStencilType,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],DM*);
51bff4a2f0SMatthew 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*);
523c48a1e8SJed Brown 
53014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGlobalToNaturalBegin(DM,Vec,InsertMode,Vec);
54014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGlobalToNaturalEnd(DM,Vec,InsertMode,Vec);
55014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDANaturalToGlobalBegin(DM,Vec,InsertMode,Vec);
56014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDANaturalToGlobalEnd(DM,Vec,InsertMode,Vec);
57*9fbee547SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use DMLocalToLocalBegin() (since version 3.5)") static inline PetscErrorCode DMDALocalToLocalBegin(DM dm,Vec g,InsertMode mode,Vec l) {return DMLocalToLocalBegin(dm,g,mode,l);}
58*9fbee547SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use DMLocalToLocalEnd() (since version 3.5)") static inline PetscErrorCode DMDALocalToLocalEnd(DM dm,Vec g,InsertMode mode,Vec l) {return DMLocalToLocalEnd(dm,g,mode,l);}
59014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDACreateNaturalVector(DM,Vec *);
603c48a1e8SJed Brown 
61014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGetCorners(DM,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
62014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGetGhostCorners(DM,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
63bff4a2f0SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMDAGetInfo(DM,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,DMBoundaryType*,DMBoundaryType*,DMBoundaryType*,DMDAStencilType*);
643ee9839eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMDAGetProcessorSubset(DM,DMDirection,PetscInt,MPI_Comm*);
653ee9839eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMDAGetProcessorSubsets(DM,DMDirection,MPI_Comm*);
663ee9839eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMDAGetRay(DM,DMDirection,PetscInt,Vec*,VecScatter*);
673c48a1e8SJed Brown 
68014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGlobalToNaturalAllCreate(DM,VecScatter*);
69014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDANaturalAllToGlobalCreate(DM,VecScatter*);
703c48a1e8SJed Brown 
71bd1fc5aeSBarry Smith PETSC_EXTERN PetscErrorCode DMDAGetScatter(DM,VecScatter*,VecScatter*);
72014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGetNeighbors(DM,const PetscMPIInt**);
733c48a1e8SJed Brown 
749db3d8bcSStefano Zampini PETSC_EXTERN PetscErrorCode DMDASetAOType(DM,AOType);
75014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGetAO(DM,AO*);
76014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDASetUniformCoordinates(DM,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal);
778272889dSSatish Balay PETSC_EXTERN PetscErrorCode DMDASetGLLCoordinates(DM,PetscInt,PetscReal*);
78c593f006SBarry Smith PETSC_EXTERN PetscErrorCode DMDAGetCoordinateArray(DM,void*);
79c593f006SBarry Smith PETSC_EXTERN PetscErrorCode DMDARestoreCoordinateArray(DM,void*);
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 
84211d3afbSPatrick Sanan PETSC_EXTERN PetscErrorCode DMDACreateCompatibleDMDA(DM,PetscInt,DM*);
8525ef9dfeSBarry Smith PETSC_EXTERN PETSC_DEPRECATED_FUNCTION("Use DMDACreateCompatibleDMDA()  (since version 3.10)") PetscErrorCode DMDAGetReducedDMDA(DM,PetscInt,DM*);
863c48a1e8SJed Brown 
87014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDASetFieldName(DM,PetscInt,const char[]);
88014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGetFieldName(DM,PetscInt,const char**);
89c629b14aSBarry Smith PETSC_EXTERN PetscErrorCode DMDASetFieldNames(DM,const char * const *);
90c629b14aSBarry Smith PETSC_EXTERN PetscErrorCode DMDAGetFieldNames(DM,const char * const **);
91109c9344SBarry Smith PETSC_EXTERN PetscErrorCode DMDASetCoordinateName(DM,PetscInt,const char[]);
92109c9344SBarry Smith PETSC_EXTERN PetscErrorCode DMDAGetCoordinateName(DM,PetscInt,const char**);
933c48a1e8SJed Brown 
94bff4a2f0SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMDASetBoundaryType(DM,DMBoundaryType,DMBoundaryType,DMBoundaryType);
95014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDASetDof(DM, PetscInt);
96fb6725baSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMDAGetDof(DM, PetscInt*);
977ddda789SPeter Brune PETSC_EXTERN PetscErrorCode DMDASetOverlap(DM,PetscInt,PetscInt,PetscInt);
987ddda789SPeter Brune PETSC_EXTERN PetscErrorCode DMDAGetOverlap(DM,PetscInt*,PetscInt*,PetscInt*);
993e7870d2SPeter Brune PETSC_EXTERN PetscErrorCode DMDASetNumLocalSubDomains(DM,PetscInt);
1003e7870d2SPeter Brune PETSC_EXTERN PetscErrorCode DMDAGetNumLocalSubDomains(DM,PetscInt*);
10195c13181SPeter Brune PETSC_EXTERN PetscErrorCode DMDAGetOffset(DM,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
10295c13181SPeter Brune PETSC_EXTERN PetscErrorCode DMDASetOffset(DM,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt);
10340234c92SPeter Brune PETSC_EXTERN PetscErrorCode DMDAGetNonOverlappingRegion(DM,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
10440234c92SPeter Brune PETSC_EXTERN PetscErrorCode DMDASetNonOverlappingRegion(DM,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt);
105014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDASetStencilWidth(DM, PetscInt);
106fb6725baSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMDAGetStencilWidth(DM, PetscInt*);
107014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDASetOwnershipRanges(DM,const PetscInt[],const PetscInt[],const PetscInt[]);
108014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGetOwnershipRanges(DM,const PetscInt**,const PetscInt**,const PetscInt**);
109014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDASetNumProcs(DM, PetscInt, PetscInt, PetscInt);
110014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDASetStencilType(DM, DMDAStencilType);
111fb6725baSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMDAGetStencilType(DM, DMDAStencilType*);
1123c48a1e8SJed Brown 
113014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAVecGetArray(DM,Vec,void *);
114014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAVecRestoreArray(DM,Vec,void *);
115fdc842d1SBarry Smith PETSC_EXTERN PetscErrorCode DMDAVecGetArrayWrite(DM,Vec,void *);
116fdc842d1SBarry Smith PETSC_EXTERN PetscErrorCode DMDAVecRestoreArrayWrite(DM,Vec,void *);
1173c48a1e8SJed Brown 
118014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAVecGetArrayDOF(DM,Vec,void *);
119014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAVecRestoreArrayDOF(DM,Vec,void *);
1203c48a1e8SJed Brown 
1215edff71fSBarry Smith PETSC_EXTERN PetscErrorCode DMDAVecGetArrayRead(DM,Vec,void *);
1225edff71fSBarry Smith PETSC_EXTERN PetscErrorCode DMDAVecRestoreArrayRead(DM,Vec,void *);
1235edff71fSBarry Smith 
1245edff71fSBarry Smith PETSC_EXTERN PetscErrorCode DMDAVecGetArrayDOFRead(DM,Vec,void *);
1255edff71fSBarry Smith PETSC_EXTERN PetscErrorCode DMDAVecRestoreArrayDOFRead(DM,Vec,void *);
1265edff71fSBarry Smith 
1271e5d2365SBarry Smith PETSC_EXTERN PetscErrorCode DMDAVecGetArrayDOFWrite(DM,Vec,void *);
1281e5d2365SBarry Smith PETSC_EXTERN PetscErrorCode DMDAVecRestoreArrayDOFWrite(DM,Vec,void *);
1291e5d2365SBarry Smith 
1303b5e53cdSSajid Ali PETSC_EXTERN PetscErrorCode DMDACreatePatchIS(DM,MatStencil*,MatStencil*,IS*,PetscBool);
131ff9846d9SPeter Brune 
1323c48a1e8SJed Brown /*MC
1333c48a1e8SJed Brown       DMDACoor2d - Structure for holding 2d (x and y) coordinates.
1343c48a1e8SJed Brown 
1353c48a1e8SJed Brown     Level: intermediate
1363c48a1e8SJed Brown 
1378e897e54SSatish Balay     Synopsis:
1383c48a1e8SJed Brown       DMDACoor2d **coors;
1393c48a1e8SJed Brown       Vec      vcoors;
1403c48a1e8SJed Brown       DM       cda;
1412150357eSBarry Smith       DMGetCoordinates(da,&vcoors);
1421957e957SBarry Smith       DMGetCoordinateDM(da,&cda);
1433c48a1e8SJed Brown       DMDAVecGetArray(cda,vcoors,&coors);
1443c48a1e8SJed Brown       DMDAGetCorners(cda,&mstart,&nstart,0,&m,&n,0)
1453c48a1e8SJed Brown       for (i=mstart; i<mstart+m; i++) {
1463c48a1e8SJed Brown         for (j=nstart; j<nstart+n; j++) {
1473c48a1e8SJed Brown           x = coors[j][i].x;
1483c48a1e8SJed Brown           y = coors[j][i].y;
1493c48a1e8SJed Brown           ......
1503c48a1e8SJed Brown         }
1513c48a1e8SJed Brown       }
1523c48a1e8SJed Brown       DMDAVecRestoreArray(dac,vcoors,&coors);
1533c48a1e8SJed Brown 
1548e897e54SSatish Balay .seealso: DMDACoor3d, DMGetCoordinateDM(), DMGetCoordinates()
1553c48a1e8SJed Brown M*/
1563c48a1e8SJed Brown typedef struct {PetscScalar x,y;} DMDACoor2d;
1573c48a1e8SJed Brown 
1583c48a1e8SJed Brown /*MC
1593c48a1e8SJed Brown       DMDACoor3d - Structure for holding 3d (x, y and z) coordinates.
1603c48a1e8SJed Brown 
1613c48a1e8SJed Brown     Level: intermediate
1623c48a1e8SJed Brown 
1638e897e54SSatish Balay     Synopsis:
1643c48a1e8SJed Brown       DMDACoor3d ***coors;
1653c48a1e8SJed Brown       Vec      vcoors;
1663c48a1e8SJed Brown       DM       cda;
1672150357eSBarry Smith       DMGetCoordinates(da,&vcoors);
1681957e957SBarry Smith       DMGetCoordinateDM(da,&cda);
1693c48a1e8SJed Brown       DMDAVecGetArray(cda,vcoors,&coors);
1703c48a1e8SJed Brown       DMDAGetCorners(cda,&mstart,&nstart,&pstart,&m,&n,&p)
1713c48a1e8SJed Brown       for (i=mstart; i<mstart+m; i++) {
1723c48a1e8SJed Brown         for (j=nstart; j<nstart+n; j++) {
1733c48a1e8SJed Brown           for (k=pstart; k<pstart+p; k++) {
1743c48a1e8SJed Brown             x = coors[k][j][i].x;
1753c48a1e8SJed Brown             y = coors[k][j][i].y;
1763c48a1e8SJed Brown             z = coors[k][j][i].z;
1773c48a1e8SJed Brown           ......
1783c48a1e8SJed Brown         }
1793c48a1e8SJed Brown       }
1803c48a1e8SJed Brown       DMDAVecRestoreArray(dac,vcoors,&coors);
1813c48a1e8SJed Brown 
1828e897e54SSatish Balay .seealso: DMDACoor2d, DMGetCoordinateDM(), DMGetCoordinates()
1833c48a1e8SJed Brown M*/
1843c48a1e8SJed Brown typedef struct {PetscScalar x,y,z;} DMDACoor3d;
1853c48a1e8SJed Brown 
186014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDAGetLocalInfo(DM,DMDALocalInfo*);
1873c48a1e8SJed Brown 
188014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRegisterDAAD(void);
189014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateDAAD(DM,Mat*);
190014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqUSFFT(Vec,DM,Mat*);
1913c48a1e8SJed Brown 
192b412c318SBarry Smith PETSC_EXTERN PetscErrorCode DMDASetGetMatrix(DM,PetscErrorCode (*)(DM, Mat *));
193ce308e1dSBarry Smith PETSC_EXTERN PetscErrorCode DMDASetBlockFills(DM,const PetscInt*,const PetscInt*);
19409e28618SBarry Smith PETSC_EXTERN PetscErrorCode DMDASetBlockFillsSparse(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 *);
2108e0841e0SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMDAComputeCellGeometryFEM(DM, PetscInt, PetscQuadrature, PetscReal [], PetscReal [], PetscReal [], PetscReal []);
2116693a731SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMDAGetTransitiveClosure(DM, PetscInt, PetscBool, PetscInt *, PetscInt **);
2126693a731SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMDARestoreTransitiveClosure(DM, PetscInt, PetscBool, PetscInt *, PetscInt **);
213f0e3914cSSatish Balay PETSC_EXTERN PetscErrorCode DMDAConvertToCell(DM, MatStencil, PetscInt *);
2143385ff7eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMDASetVertexCoordinates(DM,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal);
2158ca2f13cSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMDASetPreallocationCenterDimension(DM, PetscInt);
2168ca2f13cSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMDAGetPreallocationCenterDimension(DM, PetscInt*);
21780800b1aSMatthew G Knepley 
2183c48a1e8SJed Brown #endif
219