xref: /petsc/include/petscdmda.h (revision 2b8d69ca7ea5fe9190df62c1dce3bbd66fce84dd)
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 /*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      Determines what ghost point values are brought over to each process; in this case the "corner" values are not
17      brought over and hence should not be accessed locally
18 
19 .seealso: DMDA_STENCIL_BOX, DMDAStencilType, DMDASetStencilType()
20 M*/
21 
22 /*MC
23      DMDA_STENCIL_BOX - "Box"-type stencil. In logical grid coordinates, any of (i,j,k), (i+s,j+r,k+t) may
24                       be in the stencil.
25 
26      Level: beginner
27 
28 .seealso: DMDA_STENCIL_STAR, DMDAStencilType, DMDASetStencilType()
29 M*/
30 
31 PETSC_EXTERN PetscErrorCode DMDASetInterpolationType(DM,DMDAInterpolationType);
32 PETSC_EXTERN PetscErrorCode DMDAGetInterpolationType(DM,DMDAInterpolationType*);
33 
34 PETSC_EXTERN PetscErrorCode DMDASetElementType(DM,DMDAElementType);
35 PETSC_EXTERN PetscErrorCode DMDAGetElementType(DM,DMDAElementType*);
36 PETSC_EXTERN PetscErrorCode DMDAGetElements(DM,PetscInt *,PetscInt *,const PetscInt*[]);
37 PETSC_EXTERN PetscErrorCode DMDARestoreElements(DM,PetscInt *,PetscInt *,const PetscInt*[]);
38 
39 typedef enum { DMDA_X,DMDA_Y,DMDA_Z } DMDADirection;
40 
41 #define MATSEQUSFFT        "sequsfft"
42 
43 PETSC_EXTERN PetscErrorCode DMDACreate(MPI_Comm,DM*);
44 PETSC_EXTERN PetscErrorCode DMDASetSizes(DM,PetscInt,PetscInt,PetscInt);
45 PETSC_EXTERN PetscErrorCode DMDACreate1d(MPI_Comm,DMBoundaryType,PetscInt,PetscInt,PetscInt,const PetscInt[],DM *);
46 PETSC_EXTERN PetscErrorCode DMDACreate2d(MPI_Comm,DMBoundaryType,DMBoundaryType,DMDAStencilType,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],DM*);
47 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*);
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*,DMBoundaryType*,DMBoundaryType*,DMBoundaryType*,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 DMDAGetScatter(DM,VecScatter*,VecScatter*);
68 PETSC_EXTERN PetscErrorCode DMDAGetNeighbors(DM,const PetscMPIInt**);
69 
70 PETSC_EXTERN PetscErrorCode DMDASetAOType(DM,AOType);
71 PETSC_EXTERN PetscErrorCode DMDAGetAO(DM,AO*);
72 PETSC_EXTERN PetscErrorCode DMDASetUniformCoordinates(DM,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal);
73 PETSC_EXTERN PetscErrorCode DMDAGetCoordinateArray(DM,void*);
74 PETSC_EXTERN PetscErrorCode DMDARestoreCoordinateArray(DM,void*);
75 PETSC_EXTERN PetscErrorCode DMDAGetBoundingBox(DM,PetscReal[],PetscReal[]);
76 PETSC_EXTERN PetscErrorCode DMDAGetLocalBoundingBox(DM,PetscReal[],PetscReal[]);
77 PETSC_EXTERN PetscErrorCode DMDAGetLogicalCoordinate(DM,PetscScalar,PetscScalar,PetscScalar,PetscInt*,PetscInt*,PetscInt*,PetscScalar*,PetscScalar*,PetscScalar*);
78 /* function to wrap coordinates around boundary */
79 PETSC_EXTERN PetscErrorCode DMDAMapCoordsToPeriodicDomain(DM,PetscScalar*,PetscScalar*);
80 
81 PETSC_EXTERN PetscErrorCode DMDAGetReducedDMDA(DM,PetscInt,DM*);
82 
83 PETSC_EXTERN PetscErrorCode DMDASetFieldName(DM,PetscInt,const char[]);
84 PETSC_EXTERN PetscErrorCode DMDAGetFieldName(DM,PetscInt,const char**);
85 PETSC_EXTERN PetscErrorCode DMDASetFieldNames(DM,const char * const *);
86 PETSC_EXTERN PetscErrorCode DMDAGetFieldNames(DM,const char * const **);
87 PETSC_EXTERN PetscErrorCode DMDASetCoordinateName(DM,PetscInt,const char[]);
88 PETSC_EXTERN PetscErrorCode DMDAGetCoordinateName(DM,PetscInt,const char**);
89 
90 PETSC_EXTERN PetscErrorCode DMDASetBoundaryType(DM,DMBoundaryType,DMBoundaryType,DMBoundaryType);
91 PETSC_EXTERN PetscErrorCode DMDASetDof(DM, PetscInt);
92 PETSC_EXTERN PetscErrorCode DMDAGetDof(DM, PetscInt*);
93 PETSC_EXTERN PetscErrorCode DMDASetOverlap(DM,PetscInt,PetscInt,PetscInt);
94 PETSC_EXTERN PetscErrorCode DMDAGetOverlap(DM,PetscInt*,PetscInt*,PetscInt*);
95 PETSC_EXTERN PetscErrorCode DMDASetNumLocalSubDomains(DM,PetscInt);
96 PETSC_EXTERN PetscErrorCode DMDAGetNumLocalSubDomains(DM,PetscInt*);
97 PETSC_EXTERN PetscErrorCode DMDAGetOffset(DM,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
98 PETSC_EXTERN PetscErrorCode DMDASetOffset(DM,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt);
99 PETSC_EXTERN PetscErrorCode DMDAGetNonOverlappingRegion(DM,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
100 PETSC_EXTERN PetscErrorCode DMDASetNonOverlappingRegion(DM,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt);
101 PETSC_EXTERN PetscErrorCode DMDASetStencilWidth(DM, PetscInt);
102 PETSC_EXTERN PetscErrorCode DMDAGetStencilWidth(DM, PetscInt*);
103 PETSC_EXTERN PetscErrorCode DMDASetOwnershipRanges(DM,const PetscInt[],const PetscInt[],const PetscInt[]);
104 PETSC_EXTERN PetscErrorCode DMDAGetOwnershipRanges(DM,const PetscInt**,const PetscInt**,const PetscInt**);
105 PETSC_EXTERN PetscErrorCode DMDASetNumProcs(DM, PetscInt, PetscInt, PetscInt);
106 PETSC_EXTERN PetscErrorCode DMDASetStencilType(DM, DMDAStencilType);
107 PETSC_EXTERN PetscErrorCode DMDAGetStencilType(DM, DMDAStencilType*);
108 
109 PETSC_EXTERN PetscErrorCode DMDAVecGetArray(DM,Vec,void *);
110 PETSC_EXTERN PetscErrorCode DMDAVecRestoreArray(DM,Vec,void *);
111 
112 PETSC_EXTERN PetscErrorCode DMDAVecGetArrayDOF(DM,Vec,void *);
113 PETSC_EXTERN PetscErrorCode DMDAVecRestoreArrayDOF(DM,Vec,void *);
114 
115 PETSC_EXTERN PetscErrorCode DMDAVecGetArrayRead(DM,Vec,void *);
116 PETSC_EXTERN PetscErrorCode DMDAVecRestoreArrayRead(DM,Vec,void *);
117 
118 PETSC_EXTERN PetscErrorCode DMDAVecGetArrayDOFRead(DM,Vec,void *);
119 PETSC_EXTERN PetscErrorCode DMDAVecRestoreArrayDOFRead(DM,Vec,void *);
120 
121 PETSC_EXTERN PetscErrorCode DMDASplitComm2d(MPI_Comm,PetscInt,PetscInt,PetscInt,MPI_Comm*);
122 
123 PETSC_EXTERN PetscErrorCode DMDACreatePatchIS(DM,MatStencil*,MatStencil*,IS*);
124 
125 
126 /*MC
127       DMDACoor2d - Structure for holding 2d (x and y) coordinates.
128 
129     Level: intermediate
130 
131     Synopsis:
132       DMDACoor2d **coors;
133       Vec      vcoors;
134       DM       cda;
135       DMGetCoordinates(da,&vcoors);
136       DMGetCoordinateDM(da,&cda);
137       DMDAVecGetArray(cda,vcoors,&coors);
138       DMDAGetCorners(cda,&mstart,&nstart,0,&m,&n,0)
139       for (i=mstart; i<mstart+m; i++) {
140         for (j=nstart; j<nstart+n; j++) {
141           x = coors[j][i].x;
142           y = coors[j][i].y;
143           ......
144         }
145       }
146       DMDAVecRestoreArray(dac,vcoors,&coors);
147 
148 .seealso: DMDACoor3d, DMGetCoordinateDM(), DMGetCoordinates()
149 M*/
150 typedef struct {PetscScalar x,y;} DMDACoor2d;
151 
152 /*MC
153       DMDACoor3d - Structure for holding 3d (x, y and z) coordinates.
154 
155     Level: intermediate
156 
157     Synopsis:
158       DMDACoor3d ***coors;
159       Vec      vcoors;
160       DM       cda;
161       DMGetCoordinates(da,&vcoors);
162       DMGetCoordinateDM(da,&cda);
163       DMDAVecGetArray(cda,vcoors,&coors);
164       DMDAGetCorners(cda,&mstart,&nstart,&pstart,&m,&n,&p)
165       for (i=mstart; i<mstart+m; i++) {
166         for (j=nstart; j<nstart+n; j++) {
167           for (k=pstart; k<pstart+p; k++) {
168             x = coors[k][j][i].x;
169             y = coors[k][j][i].y;
170             z = coors[k][j][i].z;
171           ......
172         }
173       }
174       DMDAVecRestoreArray(dac,vcoors,&coors);
175 
176 .seealso: DMDACoor2d, DMGetCoordinateDM(), DMGetCoordinates()
177 M*/
178 typedef struct {PetscScalar x,y,z;} DMDACoor3d;
179 
180 PETSC_EXTERN PetscErrorCode DMDAGetLocalInfo(DM,DMDALocalInfo*);
181 
182 PETSC_EXTERN PetscErrorCode MatRegisterDAAD(void);
183 PETSC_EXTERN PetscErrorCode MatCreateDAAD(DM,Mat*);
184 PETSC_EXTERN PetscErrorCode MatCreateSeqUSFFT(Vec,DM,Mat*);
185 
186 PETSC_EXTERN PetscErrorCode DMDASetGetMatrix(DM,PetscErrorCode (*)(DM, Mat *));
187 PETSC_EXTERN PetscErrorCode DMDASetBlockFills(DM,const PetscInt*,const PetscInt*);
188 PETSC_EXTERN PetscErrorCode DMDASetRefinementFactor(DM,PetscInt,PetscInt,PetscInt);
189 PETSC_EXTERN PetscErrorCode DMDAGetRefinementFactor(DM,PetscInt*,PetscInt*,PetscInt*);
190 
191 PETSC_EXTERN PetscErrorCode DMDAGetArray(DM,PetscBool ,void*);
192 PETSC_EXTERN PetscErrorCode DMDARestoreArray(DM,PetscBool ,void*);
193 
194 PETSC_EXTERN PetscErrorCode DMDACreatePF(DM,PF*);
195 
196 /* Emulation of DMPlex */
197 PETSC_EXTERN PetscErrorCode DMDAGetNumCells(DM, PetscInt *, PetscInt *, PetscInt *, PetscInt *);
198 PETSC_EXTERN PetscErrorCode DMDAGetCellPoint(DM, PetscInt, PetscInt, PetscInt, PetscInt *);
199 PETSC_EXTERN PetscErrorCode DMDAGetNumVertices(DM, PetscInt *, PetscInt *, PetscInt *, PetscInt *);
200 PETSC_EXTERN PetscErrorCode DMDAGetNumFaces(DM, PetscInt *, PetscInt *, PetscInt *, PetscInt *, PetscInt *, PetscInt *);
201 PETSC_EXTERN PetscErrorCode DMDAGetHeightStratum(DM, PetscInt, PetscInt *, PetscInt *);
202 PETSC_EXTERN PetscErrorCode DMDAGetDepthStratum(DM, PetscInt, PetscInt *, PetscInt *);
203 PETSC_EXTERN PetscErrorCode DMDACreateSection(DM, const PetscInt[], const PetscInt[], const PetscInt[], PetscSection *);
204 PETSC_EXTERN PetscErrorCode DMDAComputeCellGeometryFEM(DM, PetscInt, PetscQuadrature, PetscReal [], PetscReal [], PetscReal [], PetscReal []);
205 PETSC_EXTERN PetscErrorCode DMDAGetTransitiveClosure(DM, PetscInt, PetscBool, PetscInt *, PetscInt **);
206 PETSC_EXTERN PetscErrorCode DMDARestoreTransitiveClosure(DM, PetscInt, PetscBool, PetscInt *, PetscInt **);
207 PETSC_EXTERN PetscErrorCode DMDAVecGetClosure(DM, PetscSection, Vec, PetscInt, PetscInt *, PetscScalar **);
208 PETSC_EXTERN PetscErrorCode DMDAVecRestoreClosure(DM, PetscSection, Vec, PetscInt, PetscInt *, PetscScalar **);
209 PETSC_EXTERN PetscErrorCode DMDAVecSetClosure(DM, PetscSection, Vec, PetscInt, const PetscScalar *, InsertMode);
210 PETSC_EXTERN PetscErrorCode DMDAGetClosure(DM, PetscSection, PetscInt, PetscInt*, const PetscInt**);
211 PETSC_EXTERN PetscErrorCode DMDARestoreClosure(DM, PetscSection, PetscInt, PetscInt*, const PetscInt**);
212 PETSC_EXTERN PetscErrorCode DMDAGetClosureScalar(DM, PetscSection, PetscInt, PetscScalar*, PetscInt*, PetscScalar**);
213 PETSC_EXTERN PetscErrorCode DMDARestoreClosureScalar(DM, PetscSection, PetscInt, PetscScalar*, PetscInt*, PetscScalar**);
214 PETSC_EXTERN PetscErrorCode DMDASetClosureScalar(DM,PetscSection,PetscInt,PetscScalar*,const PetscScalar*,InsertMode);
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