xref: /petsc/include/petscdmda.h (revision 6ffe77eaecce1557e50d00ca5347a7f48e598865)
1 #if !defined(PETSCDMDA_H)
2 #define PETSCDMDA_H
3 
4 #include <petscdm.h>
5 #include <petscdmdatypes.h>
6 #include <petscpf.h>
7 #include <petscao.h>
8 #include <petscfe.h>
9 
10 /* SUBMANSEC = DMDA */
11 
12 /*MC
13      DMDA_STENCIL_STAR - "Star"-type stencil. In logical grid coordinates, only (i,j,k), (i+s,j,k), (i,j+s,k),
14                        (i,j,k+s) are in the stencil  NOT, for example, (i+s,j+s,k)
15 
16      Level: beginner
17 
18      Determines what ghost point values are brought over to each process; in this case the "corner" values are not
19      brought over and hence should not be accessed locally
20 
21 .seealso: `DMDA_STENCIL_BOX`, `DMDAStencilType`, `DMDASetStencilType()`
22 M*/
23 
24 /*MC
25      DMDA_STENCIL_BOX - "Box"-type stencil. In logical grid coordinates, any of (i,j,k), (i+s,j+r,k+t) may
26                       be in the stencil.
27 
28      Level: beginner
29 
30 .seealso: `DMDA_STENCIL_STAR`, `DMDAStencilType`, `DMDASetStencilType()`
31 M*/
32 
33 PETSC_EXTERN PetscErrorCode DMDASetInterpolationType(DM,DMDAInterpolationType);
34 PETSC_EXTERN PetscErrorCode DMDAGetInterpolationType(DM,DMDAInterpolationType*);
35 PETSC_EXTERN PetscErrorCode DMDACreateAggregates(DM,DM,Mat*);
36 
37 /* FEM */
38 PETSC_EXTERN PetscErrorCode DMDASetElementType(DM,DMDAElementType);
39 PETSC_EXTERN PetscErrorCode DMDAGetElementType(DM,DMDAElementType*);
40 PETSC_EXTERN PetscErrorCode DMDAGetElements(DM,PetscInt*,PetscInt*,const PetscInt*[]);
41 PETSC_EXTERN PetscErrorCode DMDARestoreElements(DM,PetscInt*,PetscInt*,const PetscInt*[]);
42 PETSC_EXTERN PetscErrorCode DMDAGetElementsSizes(DM,PetscInt*,PetscInt*,PetscInt*);
43 PETSC_EXTERN PetscErrorCode DMDAGetElementsCorners(DM,PetscInt*,PetscInt*,PetscInt*);
44 PETSC_EXTERN PetscErrorCode DMDAGetSubdomainCornersIS(DM,IS*);
45 PETSC_EXTERN PetscErrorCode DMDARestoreSubdomainCornersIS(DM,IS*);
46 
47 #define MATSEQUSFFT        "sequsfft"
48 
49 PETSC_EXTERN PetscErrorCode DMDACreate(MPI_Comm,DM*);
50 PETSC_EXTERN PetscErrorCode DMDASetSizes(DM,PetscInt,PetscInt,PetscInt);
51 PETSC_EXTERN PetscErrorCode DMDACreate1d(MPI_Comm,DMBoundaryType,PetscInt,PetscInt,PetscInt,const PetscInt[],DM *);
52 PETSC_EXTERN PetscErrorCode DMDACreate2d(MPI_Comm,DMBoundaryType,DMBoundaryType,DMDAStencilType,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],DM*);
53 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*);
54 
55 PETSC_EXTERN PetscErrorCode DMDAGlobalToNaturalBegin(DM,Vec,InsertMode,Vec);
56 PETSC_EXTERN PetscErrorCode DMDAGlobalToNaturalEnd(DM,Vec,InsertMode,Vec);
57 PETSC_EXTERN PetscErrorCode DMDANaturalToGlobalBegin(DM,Vec,InsertMode,Vec);
58 PETSC_EXTERN PetscErrorCode DMDANaturalToGlobalEnd(DM,Vec,InsertMode,Vec);
59 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);}
60 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);}
61 PETSC_EXTERN PetscErrorCode DMDACreateNaturalVector(DM,Vec *);
62 
63 PETSC_EXTERN PetscErrorCode DMDAGetCorners(DM,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
64 PETSC_EXTERN PetscErrorCode DMDAGetGhostCorners(DM,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
65 PETSC_EXTERN PetscErrorCode DMDAGetInfo(DM,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,DMBoundaryType*,DMBoundaryType*,DMBoundaryType*,DMDAStencilType*);
66 PETSC_EXTERN PetscErrorCode DMDAGetProcessorSubset(DM,DMDirection,PetscInt,MPI_Comm*);
67 PETSC_EXTERN PetscErrorCode DMDAGetProcessorSubsets(DM,DMDirection,MPI_Comm*);
68 PETSC_EXTERN PetscErrorCode DMDAGetRay(DM,DMDirection,PetscInt,Vec*,VecScatter*);
69 
70 PETSC_EXTERN PetscErrorCode DMDAGlobalToNaturalAllCreate(DM,VecScatter*);
71 PETSC_EXTERN PetscErrorCode DMDANaturalAllToGlobalCreate(DM,VecScatter*);
72 
73 PETSC_EXTERN PetscErrorCode DMDAGetScatter(DM,VecScatter*,VecScatter*);
74 PETSC_EXTERN PetscErrorCode DMDAGetNeighbors(DM,const PetscMPIInt**);
75 
76 PETSC_EXTERN PetscErrorCode DMDASetAOType(DM,AOType);
77 PETSC_EXTERN PetscErrorCode DMDAGetAO(DM,AO*);
78 PETSC_EXTERN PetscErrorCode DMDASetUniformCoordinates(DM,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal);
79 PETSC_EXTERN PetscErrorCode DMDASetGLLCoordinates(DM,PetscInt,PetscReal*);
80 PETSC_EXTERN PetscErrorCode DMDAGetCoordinateArray(DM,void*);
81 PETSC_EXTERN PetscErrorCode DMDARestoreCoordinateArray(DM,void*);
82 PETSC_EXTERN PetscErrorCode DMDAGetLogicalCoordinate(DM,PetscScalar,PetscScalar,PetscScalar,PetscInt*,PetscInt*,PetscInt*,PetscScalar*,PetscScalar*,PetscScalar*);
83 /* function to wrap coordinates around boundary */
84 PETSC_EXTERN PetscErrorCode DMDAMapCoordsToPeriodicDomain(DM,PetscScalar*,PetscScalar*);
85 
86 PETSC_EXTERN PetscErrorCode DMDACreateCompatibleDMDA(DM,PetscInt,DM*);
87 PETSC_EXTERN PETSC_DEPRECATED_FUNCTION("Use DMDACreateCompatibleDMDA()  (since version 3.10)") PetscErrorCode DMDAGetReducedDMDA(DM,PetscInt,DM*);
88 
89 PETSC_EXTERN PetscErrorCode DMDASetFieldName(DM,PetscInt,const char[]);
90 PETSC_EXTERN PetscErrorCode DMDAGetFieldName(DM,PetscInt,const char**);
91 PETSC_EXTERN PetscErrorCode DMDASetFieldNames(DM,const char * const *);
92 PETSC_EXTERN PetscErrorCode DMDAGetFieldNames(DM,const char * const **);
93 PETSC_EXTERN PetscErrorCode DMDASetCoordinateName(DM,PetscInt,const char[]);
94 PETSC_EXTERN PetscErrorCode DMDAGetCoordinateName(DM,PetscInt,const char**);
95 
96 PETSC_EXTERN PetscErrorCode DMDASetBoundaryType(DM,DMBoundaryType,DMBoundaryType,DMBoundaryType);
97 PETSC_EXTERN PetscErrorCode DMDASetDof(DM, PetscInt);
98 PETSC_EXTERN PetscErrorCode DMDAGetDof(DM, PetscInt*);
99 PETSC_EXTERN PetscErrorCode DMDASetOverlap(DM,PetscInt,PetscInt,PetscInt);
100 PETSC_EXTERN PetscErrorCode DMDAGetOverlap(DM,PetscInt*,PetscInt*,PetscInt*);
101 PETSC_EXTERN PetscErrorCode DMDASetNumLocalSubDomains(DM,PetscInt);
102 PETSC_EXTERN PetscErrorCode DMDAGetNumLocalSubDomains(DM,PetscInt*);
103 PETSC_EXTERN PetscErrorCode DMDAGetOffset(DM,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
104 PETSC_EXTERN PetscErrorCode DMDASetOffset(DM,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt);
105 PETSC_EXTERN PetscErrorCode DMDAGetNonOverlappingRegion(DM,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
106 PETSC_EXTERN PetscErrorCode DMDASetNonOverlappingRegion(DM,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt);
107 PETSC_EXTERN PetscErrorCode DMDASetStencilWidth(DM, PetscInt);
108 PETSC_EXTERN PetscErrorCode DMDAGetStencilWidth(DM, PetscInt*);
109 PETSC_EXTERN PetscErrorCode DMDASetOwnershipRanges(DM,const PetscInt[],const PetscInt[],const PetscInt[]);
110 PETSC_EXTERN PetscErrorCode DMDAGetOwnershipRanges(DM,const PetscInt**,const PetscInt**,const PetscInt**);
111 PETSC_EXTERN PetscErrorCode DMDASetNumProcs(DM, PetscInt, PetscInt, PetscInt);
112 PETSC_EXTERN PetscErrorCode DMDASetStencilType(DM, DMDAStencilType);
113 PETSC_EXTERN PetscErrorCode DMDAGetStencilType(DM, DMDAStencilType*);
114 
115 PETSC_EXTERN PetscErrorCode DMDAVecGetArray(DM,Vec,void *);
116 PETSC_EXTERN PetscErrorCode DMDAVecRestoreArray(DM,Vec,void *);
117 PETSC_EXTERN PetscErrorCode DMDAVecGetArrayWrite(DM,Vec,void *);
118 PETSC_EXTERN PetscErrorCode DMDAVecRestoreArrayWrite(DM,Vec,void *);
119 
120 PETSC_EXTERN PetscErrorCode DMDAVecGetArrayDOF(DM,Vec,void *);
121 PETSC_EXTERN PetscErrorCode DMDAVecRestoreArrayDOF(DM,Vec,void *);
122 
123 PETSC_EXTERN PetscErrorCode DMDAVecGetArrayRead(DM,Vec,void *);
124 PETSC_EXTERN PetscErrorCode DMDAVecRestoreArrayRead(DM,Vec,void *);
125 
126 PETSC_EXTERN PetscErrorCode DMDAVecGetArrayDOFRead(DM,Vec,void *);
127 PETSC_EXTERN PetscErrorCode DMDAVecRestoreArrayDOFRead(DM,Vec,void *);
128 
129 PETSC_EXTERN PetscErrorCode DMDAVecGetArrayDOFWrite(DM,Vec,void *);
130 PETSC_EXTERN PetscErrorCode DMDAVecRestoreArrayDOFWrite(DM,Vec,void *);
131 
132 PETSC_EXTERN PetscErrorCode DMDACreatePatchIS(DM,MatStencil*,MatStencil*,IS*,PetscBool);
133 
134 /*MC
135       DMDACoor2d - Structure for holding 2d (x and y) coordinates.
136 
137     Level: intermediate
138 
139     Synopsis:
140       DMDACoor2d **coors;
141       Vec      vcoors;
142       DM       cda;
143       DMGetCoordinates(da,&vcoors);
144       DMGetCoordinateDM(da,&cda);
145       DMDAVecGetArray(cda,vcoors,&coors);
146       DMDAGetCorners(cda,&mstart,&nstart,0,&m,&n,0)
147       for (i=mstart; i<mstart+m; i++) {
148         for (j=nstart; j<nstart+n; j++) {
149           x = coors[j][i].x;
150           y = coors[j][i].y;
151           ......
152         }
153       }
154       DMDAVecRestoreArray(dac,vcoors,&coors);
155 
156 .seealso: `DMDACoor3d`, `DMGetCoordinateDM()`, `DMGetCoordinates()`
157 M*/
158 typedef struct {PetscScalar x,y;} DMDACoor2d;
159 
160 /*MC
161       DMDACoor3d - Structure for holding 3d (x, y and z) coordinates.
162 
163     Level: intermediate
164 
165     Synopsis:
166       DMDACoor3d ***coors;
167       Vec      vcoors;
168       DM       cda;
169       DMGetCoordinates(da,&vcoors);
170       DMGetCoordinateDM(da,&cda);
171       DMDAVecGetArray(cda,vcoors,&coors);
172       DMDAGetCorners(cda,&mstart,&nstart,&pstart,&m,&n,&p)
173       for (i=mstart; i<mstart+m; i++) {
174         for (j=nstart; j<nstart+n; j++) {
175           for (k=pstart; k<pstart+p; k++) {
176             x = coors[k][j][i].x;
177             y = coors[k][j][i].y;
178             z = coors[k][j][i].z;
179           ......
180         }
181       }
182       DMDAVecRestoreArray(dac,vcoors,&coors);
183 
184 .seealso: `DMDACoor2d`, `DMGetCoordinateDM()`, `DMGetCoordinates()`
185 M*/
186 typedef struct {PetscScalar x,y,z;} DMDACoor3d;
187 
188 PETSC_EXTERN PetscErrorCode DMDAGetLocalInfo(DM,DMDALocalInfo*);
189 
190 PETSC_EXTERN PetscErrorCode MatRegisterDAAD(void);
191 PETSC_EXTERN PetscErrorCode MatCreateDAAD(DM,Mat*);
192 PETSC_EXTERN PetscErrorCode MatCreateSeqUSFFT(Vec,DM,Mat*);
193 
194 PETSC_EXTERN PetscErrorCode DMDASetGetMatrix(DM,PetscErrorCode (*)(DM, Mat *));
195 PETSC_EXTERN PetscErrorCode DMDASetBlockFills(DM,const PetscInt*,const PetscInt*);
196 PETSC_EXTERN PetscErrorCode DMDASetBlockFillsSparse(DM,const PetscInt*,const PetscInt*);
197 PETSC_EXTERN PetscErrorCode DMDASetRefinementFactor(DM,PetscInt,PetscInt,PetscInt);
198 PETSC_EXTERN PetscErrorCode DMDAGetRefinementFactor(DM,PetscInt*,PetscInt*,PetscInt*);
199 
200 PETSC_EXTERN PetscErrorCode DMDAGetArray(DM,PetscBool ,void*);
201 PETSC_EXTERN PetscErrorCode DMDARestoreArray(DM,PetscBool ,void*);
202 
203 PETSC_EXTERN PetscErrorCode DMDACreatePF(DM,PF*);
204 
205 /* Emulation of DMPlex */
206 PETSC_EXTERN PetscErrorCode DMDAGetNumCells(DM, PetscInt *, PetscInt *, PetscInt *, PetscInt *);
207 PETSC_EXTERN PetscErrorCode DMDAGetCellPoint(DM, PetscInt, PetscInt, PetscInt, PetscInt *);
208 PETSC_EXTERN PetscErrorCode DMDAGetNumVertices(DM, PetscInt *, PetscInt *, PetscInt *, PetscInt *);
209 PETSC_EXTERN PetscErrorCode DMDAGetNumFaces(DM, PetscInt *, PetscInt *, PetscInt *, PetscInt *, PetscInt *, PetscInt *);
210 PETSC_EXTERN PetscErrorCode DMDAGetHeightStratum(DM, PetscInt, PetscInt *, PetscInt *);
211 PETSC_EXTERN PetscErrorCode DMDAGetDepthStratum(DM, PetscInt, PetscInt *, PetscInt *);
212 PETSC_EXTERN PetscErrorCode DMDAComputeCellGeometryFEM(DM, PetscInt, PetscQuadrature, PetscReal [], PetscReal [], PetscReal [], PetscReal []);
213 PETSC_EXTERN PetscErrorCode DMDAGetTransitiveClosure(DM, PetscInt, PetscBool, PetscInt *, PetscInt **);
214 PETSC_EXTERN PetscErrorCode DMDARestoreTransitiveClosure(DM, PetscInt, PetscBool, PetscInt *, PetscInt **);
215 PETSC_EXTERN PetscErrorCode DMDAConvertToCell(DM, MatStencil, PetscInt *);
216 PETSC_EXTERN PetscErrorCode DMDASetVertexCoordinates(DM,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal);
217 PETSC_EXTERN PetscErrorCode DMDASetPreallocationCenterDimension(DM, PetscInt);
218 PETSC_EXTERN PetscErrorCode DMDAGetPreallocationCenterDimension(DM, PetscInt*);
219 
220 #endif
221