xref: /petsc/include/petscdmda.h (revision 785e854f82a3c614b452fca2cf5ad4f2afe8bdde)
1 #if !defined(__PETSCDMDA_H)
2 #define __PETSCDMDA_H
3 
4 #include <petscdm.h>
5 #include <petscdt.h>
6 #include <petscdmdatypes.h>
7 #include <petscpf.h>
8 #include <petscao.h>
9 
10 /*MC
11      DMDA_STENCIL_STAR - "Star"-type stencil. In logical grid coordinates, only (i,j,k), (i+s,j,k), (i,j+s,k),
12                        (i,j,k+s) are in the stencil  NOT, for example, (i+s,j+s,k)
13 
14      Level: beginner
15 
16 .seealso: DMDA_STENCIL_BOX, DMDAStencilType, DMDASetStencilType()
17 M*/
18 
19 /*MC
20      DMDA_STENCIL_BOX - "Box"-type stencil. In logical grid coordinates, any of (i,j,k), (i+s,j+r,k+t) may
21                       be in the stencil.
22 
23      Level: beginner
24 
25 .seealso: DMDA_STENCIL_STAR, DMDAStencilType, DMDASetStencilType()
26 M*/
27 
28 PETSC_EXTERN const char *const DMDABoundaryTypes[];
29 
30 PETSC_EXTERN PetscErrorCode DMDASetInterpolationType(DM,DMDAInterpolationType);
31 PETSC_EXTERN PetscErrorCode DMDAGetInterpolationType(DM,DMDAInterpolationType*);
32 
33 PETSC_EXTERN PetscErrorCode DMDASetElementType(DM,DMDAElementType);
34 PETSC_EXTERN PetscErrorCode DMDAGetElementType(DM,DMDAElementType*);
35 PETSC_EXTERN PetscErrorCode DMDAGetElements(DM,PetscInt *,PetscInt *,const PetscInt*[]);
36 PETSC_EXTERN PetscErrorCode DMDARestoreElements(DM,PetscInt *,PetscInt *,const PetscInt*[]);
37 
38 typedef enum { DMDA_X,DMDA_Y,DMDA_Z } DMDADirection;
39 
40 #define MATSEQUSFFT        "sequsfft"
41 
42 PETSC_EXTERN PetscErrorCode DMDACreate(MPI_Comm,DM*);
43 PETSC_EXTERN PetscErrorCode DMDASetDim(DM,PetscInt);
44 PETSC_EXTERN PetscErrorCode DMDASetSizes(DM,PetscInt,PetscInt,PetscInt);
45 PETSC_EXTERN PetscErrorCode DMDACreate1d(MPI_Comm,DMDABoundaryType,PetscInt,PetscInt,PetscInt,const PetscInt[],DM *);
46 PETSC_EXTERN PetscErrorCode DMDACreate2d(MPI_Comm,DMDABoundaryType,DMDABoundaryType,DMDAStencilType,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],DM*);
47 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*);
48 
49 PETSC_EXTERN PetscErrorCode DMDAGlobalToNaturalBegin(DM,Vec,InsertMode,Vec);
50 PETSC_EXTERN PetscErrorCode DMDAGlobalToNaturalEnd(DM,Vec,InsertMode,Vec);
51 PETSC_EXTERN PetscErrorCode DMDANaturalToGlobalBegin(DM,Vec,InsertMode,Vec);
52 PETSC_EXTERN PetscErrorCode DMDANaturalToGlobalEnd(DM,Vec,InsertMode,Vec);
53 PETSC_DEPRECATED("Use DMLocalToLocalBegin()") PETSC_STATIC_INLINE PetscErrorCode DMDALocalToLocalBegin(DM dm,Vec g,InsertMode mode,Vec l) {return DMLocalToLocalBegin(dm,g,mode,l);}
54 PETSC_DEPRECATED("Use DMLocalToLocalEnd()") PETSC_STATIC_INLINE PetscErrorCode DMDALocalToLocalEnd(DM dm,Vec g,InsertMode mode,Vec l) {return DMLocalToLocalEnd(dm,g,mode,l);}
55 PETSC_EXTERN PetscErrorCode DMDACreateNaturalVector(DM,Vec *);
56 
57 PETSC_EXTERN PetscErrorCode DMDAGetCorners(DM,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
58 PETSC_EXTERN PetscErrorCode DMDAGetGhostCorners(DM,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
59 PETSC_EXTERN PetscErrorCode DMDAGetInfo(DM,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,DMDABoundaryType*,DMDABoundaryType*,DMDABoundaryType*,DMDAStencilType*);
60 PETSC_EXTERN PetscErrorCode DMDAGetProcessorSubset(DM,DMDADirection,PetscInt,MPI_Comm*);
61 PETSC_EXTERN PetscErrorCode DMDAGetProcessorSubsets(DM,DMDADirection,MPI_Comm*);
62 PETSC_EXTERN PetscErrorCode DMDAGetRay(DM,DMDADirection,PetscInt,Vec*,VecScatter*);
63 
64 PETSC_EXTERN PetscErrorCode DMDAGlobalToNaturalAllCreate(DM,VecScatter*);
65 PETSC_EXTERN PetscErrorCode DMDANaturalAllToGlobalCreate(DM,VecScatter*);
66 
67 PETSC_EXTERN PetscErrorCode DMDAGetGlobalIndices(DM,PetscInt*,PetscInt**);
68 
69 PETSC_EXTERN PetscErrorCode DMDAGetScatter(DM,VecScatter*,VecScatter*,VecScatter*);
70 PETSC_EXTERN PetscErrorCode DMDAGetNeighbors(DM,const PetscMPIInt**);
71 
72 PETSC_EXTERN PetscErrorCode DMDAGetAO(DM,AO*);
73 PETSC_EXTERN PetscErrorCode DMDASetUniformCoordinates(DM,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal);
74 PETSC_EXTERN PetscErrorCode DMDAGetBoundingBox(DM,PetscReal[],PetscReal[]);
75 PETSC_EXTERN PetscErrorCode DMDAGetLocalBoundingBox(DM,PetscReal[],PetscReal[]);
76 PETSC_EXTERN PetscErrorCode DMDAGetLogicalCoordinate(DM,PetscScalar,PetscScalar,PetscScalar,PetscInt*,PetscInt*,PetscInt*,PetscScalar*,PetscScalar*,PetscScalar*);
77 /* function to wrap coordinates around boundary */
78 PETSC_EXTERN PetscErrorCode DMDAMapCoordsToPeriodicDomain(DM,PetscScalar*,PetscScalar*);
79 
80 PETSC_EXTERN PetscErrorCode DMDAGetReducedDMDA(DM,PetscInt,DM*);
81 
82 PETSC_EXTERN PetscErrorCode DMDASetFieldName(DM,PetscInt,const char[]);
83 PETSC_EXTERN PetscErrorCode DMDAGetFieldName(DM,PetscInt,const char**);
84 PETSC_EXTERN PetscErrorCode DMDASetCoordinateName(DM,PetscInt,const char[]);
85 PETSC_EXTERN PetscErrorCode DMDAGetCoordinateName(DM,PetscInt,const char**);
86 
87 PETSC_EXTERN PetscErrorCode DMDASetBoundaryType(DM,DMDABoundaryType,DMDABoundaryType,DMDABoundaryType);
88 PETSC_EXTERN PetscErrorCode DMDASetDof(DM, PetscInt);
89 PETSC_EXTERN PetscErrorCode DMDASetOverlap(DM,PetscInt,PetscInt,PetscInt);
90 PETSC_EXTERN PetscErrorCode DMDAGetOverlap(DM,PetscInt*,PetscInt*,PetscInt*);
91 PETSC_EXTERN PetscErrorCode DMDASetNumLocalSubDomains(DM,PetscInt);
92 PETSC_EXTERN PetscErrorCode DMDAGetNumLocalSubDomains(DM,PetscInt*);
93 PETSC_EXTERN PetscErrorCode DMDAGetOffset(DM,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
94 PETSC_EXTERN PetscErrorCode DMDASetOffset(DM,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt);
95 PETSC_EXTERN PetscErrorCode DMDAGetNonOverlappingRegion(DM,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
96 PETSC_EXTERN PetscErrorCode DMDASetNonOverlappingRegion(DM,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt);
97 PETSC_EXTERN PetscErrorCode DMDASetStencilWidth(DM, PetscInt);
98 PETSC_EXTERN PetscErrorCode DMDASetOwnershipRanges(DM,const PetscInt[],const PetscInt[],const PetscInt[]);
99 PETSC_EXTERN PetscErrorCode DMDAGetOwnershipRanges(DM,const PetscInt**,const PetscInt**,const PetscInt**);
100 PETSC_EXTERN PetscErrorCode DMDASetNumProcs(DM, PetscInt, PetscInt, PetscInt);
101 PETSC_EXTERN PetscErrorCode DMDASetStencilType(DM, DMDAStencilType);
102 
103 PETSC_EXTERN PetscErrorCode DMDAVecGetArray(DM,Vec,void *);
104 PETSC_EXTERN PetscErrorCode DMDAVecRestoreArray(DM,Vec,void *);
105 
106 PETSC_EXTERN PetscErrorCode DMDAVecGetArrayDOF(DM,Vec,void *);
107 PETSC_EXTERN PetscErrorCode DMDAVecRestoreArrayDOF(DM,Vec,void *);
108 
109 PETSC_EXTERN PetscErrorCode DMDASplitComm2d(MPI_Comm,PetscInt,PetscInt,PetscInt,MPI_Comm*);
110 
111 PETSC_EXTERN PetscErrorCode DMDACreatePatchIS(DM,MatStencil*,MatStencil*,IS*);
112 
113 
114 /*MC
115       DMDACoor2d - Structure for holding 2d (x and y) coordinates.
116 
117     Level: intermediate
118 
119     Sample Usage:
120       DMDACoor2d **coors;
121       Vec      vcoors;
122       DM       cda;
123 
124       DMGetCoordinates(da,&vcoors);
125       DMDAGetCoordinateDA(da,&cda);
126       DMDAVecGetArray(cda,vcoors,&coors);
127       DMDAGetCorners(cda,&mstart,&nstart,0,&m,&n,0)
128       for (i=mstart; i<mstart+m; i++) {
129         for (j=nstart; j<nstart+n; j++) {
130           x = coors[j][i].x;
131           y = coors[j][i].y;
132           ......
133         }
134       }
135       DMDAVecRestoreArray(dac,vcoors,&coors);
136 
137 .seealso: DMDACoor3d, DMDAForEachPointBegin(), DMDAGetCoordinateDA(), DMGetCoordinates(), DMDAGetGhostCoordinates()
138 M*/
139 typedef struct {PetscScalar x,y;} DMDACoor2d;
140 
141 /*MC
142       DMDACoor3d - Structure for holding 3d (x, y and z) coordinates.
143 
144     Level: intermediate
145 
146     Sample Usage:
147       DMDACoor3d ***coors;
148       Vec      vcoors;
149       DM       cda;
150 
151       DMGetCoordinates(da,&vcoors);
152       DMDAGetCoordinateDA(da,&cda);
153       DMDAVecGetArray(cda,vcoors,&coors);
154       DMDAGetCorners(cda,&mstart,&nstart,&pstart,&m,&n,&p)
155       for (i=mstart; i<mstart+m; i++) {
156         for (j=nstart; j<nstart+n; j++) {
157           for (k=pstart; k<pstart+p; k++) {
158             x = coors[k][j][i].x;
159             y = coors[k][j][i].y;
160             z = coors[k][j][i].z;
161           ......
162         }
163       }
164       DMDAVecRestoreArray(dac,vcoors,&coors);
165 
166 .seealso: DMDACoor2d, DMDAForEachPointBegin(), DMDAGetCoordinateDA(), DMGetCoordinates(), DMDAGetGhostCoordinates()
167 M*/
168 typedef struct {PetscScalar x,y,z;} DMDACoor3d;
169 
170 PETSC_EXTERN PetscErrorCode DMDAGetLocalInfo(DM,DMDALocalInfo*);
171 
172 PETSC_EXTERN PetscErrorCode MatRegisterDAAD(void);
173 PETSC_EXTERN PetscErrorCode MatCreateDAAD(DM,Mat*);
174 PETSC_EXTERN PetscErrorCode MatCreateSeqUSFFT(Vec,DM,Mat*);
175 
176 PETSC_EXTERN PetscErrorCode DMDASetGetMatrix(DM,PetscErrorCode (*)(DM, Mat *));
177 PETSC_EXTERN PetscErrorCode DMDASetBlockFills(DM,const PetscInt*,const PetscInt*);
178 PETSC_EXTERN PetscErrorCode DMDASetRefinementFactor(DM,PetscInt,PetscInt,PetscInt);
179 PETSC_EXTERN PetscErrorCode DMDAGetRefinementFactor(DM,PetscInt*,PetscInt*,PetscInt*);
180 
181 PETSC_EXTERN PetscErrorCode DMDAGetArray(DM,PetscBool ,void*);
182 PETSC_EXTERN PetscErrorCode DMDARestoreArray(DM,PetscBool ,void*);
183 
184 PETSC_EXTERN PetscErrorCode DMDACreatePF(DM,PF*);
185 
186 PETSC_EXTERN PetscErrorCode DMDAGetNumCells(DM, PetscInt *);
187 PETSC_EXTERN PetscErrorCode DMDAGetNumVertices(DM, PetscInt *, PetscInt *, PetscInt *, PetscInt *);
188 PETSC_EXTERN PetscErrorCode DMDAGetNumFaces(DM, PetscInt *, PetscInt *, PetscInt *, PetscInt *, PetscInt *, PetscInt *);
189 PETSC_EXTERN PetscErrorCode DMDAGetHeightStratum(DM, PetscInt, PetscInt *, PetscInt *);
190 PETSC_EXTERN PetscErrorCode DMDACreateSection(DM, PetscInt[], PetscInt[], PetscInt[], PetscInt[]);
191 PETSC_EXTERN PetscErrorCode DMDAComputeCellGeometry(DM, PetscInt, PetscQuadrature *, PetscReal [], PetscReal [], PetscReal [], PetscReal []);
192 PETSC_EXTERN PetscErrorCode DMDAVecGetClosure(DM, PetscSection, Vec, PetscInt, const PetscScalar **);
193 PETSC_EXTERN PetscErrorCode DMDAVecSetClosure(DM, PetscSection, Vec, PetscInt, const PetscScalar *, InsertMode);
194 PETSC_EXTERN PetscErrorCode DMDAGetClosure(DM,PetscSection,PetscInt,PetscInt*,const PetscInt**);
195 PETSC_EXTERN PetscErrorCode DMDARestoreClosure(DM,PetscSection,PetscInt,PetscInt*,const PetscInt**);
196 PETSC_EXTERN PetscErrorCode DMDAGetClosureScalar(DM,PetscSection,PetscInt,PetscScalar*,const PetscScalar**);
197 PETSC_EXTERN PetscErrorCode DMDARestoreClosureScalar(DM,PetscSection,PetscInt,PetscScalar*,const PetscScalar**);
198 PETSC_EXTERN PetscErrorCode DMDASetClosureScalar(DM,PetscSection,PetscInt,PetscScalar*,const PetscScalar*,InsertMode);
199 PETSC_EXTERN PetscErrorCode DMDAConvertToCell(DM, MatStencil, PetscInt *);
200 
201 #endif
202