xref: /petsc/include/petscdmda.h (revision 69bbac976351bc191e2527ca4ac8ff0093de2015)
1 #if !defined(__PETSCDMDA_H)
2 #define __PETSCDMDA_H
3 
4 #include <petscdm.h>
5 #include <petscdt.h>
6 #include <petscfetypes.h>
7 #include <petscdmdatypes.h>
8 #include <petscpf.h>
9 #include <petscao.h>
10 #include <petscfe.h>
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 .seealso: DMDA_STENCIL_BOX, DMDAStencilType, DMDASetStencilType()
19 M*/
20 
21 /*MC
22      DMDA_STENCIL_BOX - "Box"-type stencil. In logical grid coordinates, any of (i,j,k), (i+s,j+r,k+t) may
23                       be in the stencil.
24 
25      Level: beginner
26 
27 .seealso: DMDA_STENCIL_STAR, DMDAStencilType, DMDASetStencilType()
28 M*/
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,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 DMDAGetGlobalIndices(DM,PetscInt*,const PetscInt*[]);
68 PETSC_EXTERN PetscErrorCode DMDARestoreGlobalIndices(DM,PetscInt*,const PetscInt*[]);
69 
70 PETSC_EXTERN PetscErrorCode DMDAGetScatter(DM,VecScatter*,VecScatter*,VecScatter*);
71 PETSC_EXTERN PetscErrorCode DMDAGetNeighbors(DM,const PetscMPIInt**);
72 
73 PETSC_EXTERN PetscErrorCode DMDAGetAO(DM,AO*);
74 PETSC_EXTERN PetscErrorCode DMDASetUniformCoordinates(DM,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal);
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 DMDASetCoordinateName(DM,PetscInt,const char[]);
86 PETSC_EXTERN PetscErrorCode DMDAGetCoordinateName(DM,PetscInt,const char**);
87 
88 PETSC_EXTERN PetscErrorCode DMDASetBoundaryType(DM,DMBoundaryType,DMBoundaryType,DMBoundaryType);
89 PETSC_EXTERN PetscErrorCode DMDASetDof(DM, PetscInt);
90 PETSC_EXTERN PetscErrorCode DMDASetOverlap(DM,PetscInt,PetscInt,PetscInt);
91 PETSC_EXTERN PetscErrorCode DMDAGetOverlap(DM,PetscInt*,PetscInt*,PetscInt*);
92 PETSC_EXTERN PetscErrorCode DMDASetNumLocalSubDomains(DM,PetscInt);
93 PETSC_EXTERN PetscErrorCode DMDAGetNumLocalSubDomains(DM,PetscInt*);
94 PETSC_EXTERN PetscErrorCode DMDAGetOffset(DM,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
95 PETSC_EXTERN PetscErrorCode DMDASetOffset(DM,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt);
96 PETSC_EXTERN PetscErrorCode DMDAGetNonOverlappingRegion(DM,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
97 PETSC_EXTERN PetscErrorCode DMDASetNonOverlappingRegion(DM,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt);
98 PETSC_EXTERN PetscErrorCode DMDASetStencilWidth(DM, PetscInt);
99 PETSC_EXTERN PetscErrorCode DMDASetOwnershipRanges(DM,const PetscInt[],const PetscInt[],const PetscInt[]);
100 PETSC_EXTERN PetscErrorCode DMDAGetOwnershipRanges(DM,const PetscInt**,const PetscInt**,const PetscInt**);
101 PETSC_EXTERN PetscErrorCode DMDASetNumProcs(DM, PetscInt, PetscInt, PetscInt);
102 PETSC_EXTERN PetscErrorCode DMDASetStencilType(DM, DMDAStencilType);
103 
104 PETSC_EXTERN PetscErrorCode DMDAVecGetArray(DM,Vec,void *);
105 PETSC_EXTERN PetscErrorCode DMDAVecRestoreArray(DM,Vec,void *);
106 
107 PETSC_EXTERN PetscErrorCode DMDAVecGetArrayDOF(DM,Vec,void *);
108 PETSC_EXTERN PetscErrorCode DMDAVecRestoreArrayDOF(DM,Vec,void *);
109 
110 PETSC_EXTERN PetscErrorCode DMDASplitComm2d(MPI_Comm,PetscInt,PetscInt,PetscInt,MPI_Comm*);
111 
112 PETSC_EXTERN PetscErrorCode DMDACreatePatchIS(DM,MatStencil*,MatStencil*,IS*);
113 
114 
115 /*MC
116       DMDACoor2d - Structure for holding 2d (x and y) coordinates.
117 
118     Level: intermediate
119 
120     Sample Usage:
121       DMDACoor2d **coors;
122       Vec      vcoors;
123       DM       cda;
124 
125       DMGetCoordinates(da,&vcoors);
126       DMDAGetCoordinateDA(da,&cda);
127       DMDAVecGetArray(cda,vcoors,&coors);
128       DMDAGetCorners(cda,&mstart,&nstart,0,&m,&n,0)
129       for (i=mstart; i<mstart+m; i++) {
130         for (j=nstart; j<nstart+n; j++) {
131           x = coors[j][i].x;
132           y = coors[j][i].y;
133           ......
134         }
135       }
136       DMDAVecRestoreArray(dac,vcoors,&coors);
137 
138 .seealso: DMDACoor3d, DMDAForEachPointBegin(), DMDAGetCoordinateDA(), DMGetCoordinates(), DMDAGetGhostCoordinates()
139 M*/
140 typedef struct {PetscScalar x,y;} DMDACoor2d;
141 
142 /*MC
143       DMDACoor3d - Structure for holding 3d (x, y and z) coordinates.
144 
145     Level: intermediate
146 
147     Sample Usage:
148       DMDACoor3d ***coors;
149       Vec      vcoors;
150       DM       cda;
151 
152       DMGetCoordinates(da,&vcoors);
153       DMDAGetCoordinateDA(da,&cda);
154       DMDAVecGetArray(cda,vcoors,&coors);
155       DMDAGetCorners(cda,&mstart,&nstart,&pstart,&m,&n,&p)
156       for (i=mstart; i<mstart+m; i++) {
157         for (j=nstart; j<nstart+n; j++) {
158           for (k=pstart; k<pstart+p; k++) {
159             x = coors[k][j][i].x;
160             y = coors[k][j][i].y;
161             z = coors[k][j][i].z;
162           ......
163         }
164       }
165       DMDAVecRestoreArray(dac,vcoors,&coors);
166 
167 .seealso: DMDACoor2d, DMDAForEachPointBegin(), DMDAGetCoordinateDA(), DMGetCoordinates(), DMDAGetGhostCoordinates()
168 M*/
169 typedef struct {PetscScalar x,y,z;} DMDACoor3d;
170 
171 PETSC_EXTERN PetscErrorCode DMDAGetLocalInfo(DM,DMDALocalInfo*);
172 
173 PETSC_EXTERN PetscErrorCode MatRegisterDAAD(void);
174 PETSC_EXTERN PetscErrorCode MatCreateDAAD(DM,Mat*);
175 PETSC_EXTERN PetscErrorCode MatCreateSeqUSFFT(Vec,DM,Mat*);
176 
177 PETSC_EXTERN PetscErrorCode DMDASetGetMatrix(DM,PetscErrorCode (*)(DM, Mat *));
178 PETSC_EXTERN PetscErrorCode DMDASetBlockFills(DM,const PetscInt*,const PetscInt*);
179 PETSC_EXTERN PetscErrorCode DMDASetRefinementFactor(DM,PetscInt,PetscInt,PetscInt);
180 PETSC_EXTERN PetscErrorCode DMDAGetRefinementFactor(DM,PetscInt*,PetscInt*,PetscInt*);
181 
182 PETSC_EXTERN PetscErrorCode DMDAGetArray(DM,PetscBool ,void*);
183 PETSC_EXTERN PetscErrorCode DMDARestoreArray(DM,PetscBool ,void*);
184 
185 PETSC_EXTERN PetscErrorCode DMDACreatePF(DM,PF*);
186 
187 /* Emulation of DMPlex */
188 PETSC_EXTERN PetscErrorCode DMDAGetNumCells(DM, PetscInt *, PetscInt *, PetscInt *, PetscInt *);
189 PETSC_EXTERN PetscErrorCode DMDAGetCellPoint(DM, PetscInt, PetscInt, PetscInt, PetscInt *);
190 PETSC_EXTERN PetscErrorCode DMDAGetNumVertices(DM, PetscInt *, PetscInt *, PetscInt *, PetscInt *);
191 PETSC_EXTERN PetscErrorCode DMDAGetNumFaces(DM, PetscInt *, PetscInt *, PetscInt *, PetscInt *, PetscInt *, PetscInt *);
192 PETSC_EXTERN PetscErrorCode DMDAGetHeightStratum(DM, PetscInt, PetscInt *, PetscInt *);
193 PETSC_EXTERN PetscErrorCode DMDACreateSection(DM, PetscInt[], PetscInt[], PetscInt[], PetscSection *);
194 PETSC_EXTERN PetscErrorCode DMDAComputeCellGeometry(DM, PetscInt, PetscQuadrature *, PetscReal [], PetscReal [], PetscReal [], PetscReal []);
195 PETSC_EXTERN PetscErrorCode DMDAGetTransitiveClosure(DM, PetscInt, PetscBool, PetscInt *, PetscInt **);
196 PETSC_EXTERN PetscErrorCode DMDARestoreTransitiveClosure(DM, PetscInt, PetscBool, PetscInt *, PetscInt **);
197 PETSC_EXTERN PetscErrorCode DMDAVecGetClosure(DM, PetscSection, Vec, PetscInt, PetscInt *, PetscScalar **);
198 PETSC_EXTERN PetscErrorCode DMDAVecRestoreClosure(DM, PetscSection, Vec, PetscInt, PetscInt *, PetscScalar **);
199 PETSC_EXTERN PetscErrorCode DMDAVecSetClosure(DM, PetscSection, Vec, PetscInt, const PetscScalar *, InsertMode);
200 PETSC_EXTERN PetscErrorCode DMDAGetClosure(DM, PetscSection, PetscInt, PetscInt*, const PetscInt**);
201 PETSC_EXTERN PetscErrorCode DMDARestoreClosure(DM, PetscSection, PetscInt, PetscInt*, const PetscInt**);
202 PETSC_EXTERN PetscErrorCode DMDAGetClosureScalar(DM, PetscSection, PetscInt, PetscScalar*, PetscInt*, PetscScalar**);
203 PETSC_EXTERN PetscErrorCode DMDARestoreClosureScalar(DM, PetscSection, PetscInt, PetscScalar*, PetscInt*, PetscScalar**);
204 PETSC_EXTERN PetscErrorCode DMDASetClosureScalar(DM,PetscSection,PetscInt,PetscScalar*,const PetscScalar*,InsertMode);
205 PETSC_EXTERN PetscErrorCode DMDAConvertToCell(DM, MatStencil, PetscInt *);
206 PETSC_EXTERN PetscErrorCode DMDASetVertexCoordinates(DM,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal);
207 PETSC_EXTERN PetscErrorCode DMDASetPreallocationCenterDimension(DM, PetscInt);
208 PETSC_EXTERN PetscErrorCode DMDAGetPreallocationCenterDimension(DM, PetscInt*);
209 
210 PETSC_EXTERN PetscErrorCode DMDAProjectFunction(DM, PetscFE[], void (**)(const PetscReal [], PetscScalar *, void *), void **, InsertMode, Vec);
211 PETSC_EXTERN PetscErrorCode DMDAProjectFunctionLocal(DM, PetscFE[], void (**)(const PetscReal [], PetscScalar *, void *), void **, InsertMode, Vec);
212 PETSC_EXTERN PetscErrorCode DMDAComputeL2Diff(DM, PetscFE[], void (**)(const PetscReal [], PetscScalar *, void *), void **, Vec, PetscReal *);
213 
214 #endif
215