1 #if !defined(__PETSCDMDA_H) 2 #define __PETSCDMDA_H 3 4 #include <petscdm.h> 5 #include <petscpf.h> 6 #include <petscao.h> 7 8 /*E 9 DMDAStencilType - Determines if the stencil extends only along the coordinate directions, or also 10 to the northeast, northwest etc 11 12 Level: beginner 13 14 .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDACreate(), DMDASetStencilType() 15 E*/ 16 typedef enum { DMDA_STENCIL_STAR,DMDA_STENCIL_BOX } DMDAStencilType; 17 18 /*MC 19 DMDA_STENCIL_STAR - "Star"-type stencil. In logical grid coordinates, only (i,j,k), (i+s,j,k), (i,j+s,k), 20 (i,j,k+s) are in the stencil NOT, for example, (i+s,j+s,k) 21 22 Level: beginner 23 24 .seealso: DMDA_STENCIL_BOX, DMDAStencilType, DMDASetStencilType() 25 M*/ 26 27 /*MC 28 DMDA_STENCIL_BOX - "Box"-type stencil. In logical grid coordinates, any of (i,j,k), (i+s,j+r,k+t) may 29 be in the stencil. 30 31 Level: beginner 32 33 .seealso: DMDA_STENCIL_STAR, DMDAStencilType, DMDASetStencilType() 34 M*/ 35 36 /*E 37 DMDABoundaryType - Describes the choice for fill of ghost cells on physical domain boundaries. 38 39 Level: beginner 40 41 A boundary may be of type DMDA_BOUNDARY_NONE (no ghost nodes), DMDA_BOUNDARY_GHOST (ghost nodes 42 exist but aren't filled, you can put values into them and then apply a stencil that uses those ghost locations), 43 DMDA_BOUNDARY_MIRROR (not yet implemented for 3d), or DMDA_BOUNDARY_PERIODIC 44 (ghost nodes filled by the opposite edge of the domain). 45 46 Note: This is information for the boundary of the __PHYSICAL__ domain. It has nothing to do with boundaries between 47 processes, that width is always determined by the stencil width, see DMDASetStencilWidth(). 48 49 .seealso: DMDASetBoundaryType(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDACreate() 50 E*/ 51 typedef enum { DMDA_BOUNDARY_NONE, DMDA_BOUNDARY_GHOSTED, DMDA_BOUNDARY_MIRROR, DMDA_BOUNDARY_PERIODIC } DMDABoundaryType; 52 53 PETSC_EXTERN const char *const DMDABoundaryTypes[]; 54 55 /*E 56 DMDAInterpolationType - Defines the type of interpolation that will be returned by 57 DMCreateInterpolation. 58 59 Level: beginner 60 61 .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMCreateInterpolation(), DMDASetInterpolationType(), DMDACreate() 62 E*/ 63 typedef enum { DMDA_Q0, DMDA_Q1 } DMDAInterpolationType; 64 65 PETSC_EXTERN PetscErrorCode DMDASetInterpolationType(DM,DMDAInterpolationType); 66 PETSC_EXTERN PetscErrorCode DMDAGetInterpolationType(DM,DMDAInterpolationType*); 67 68 /*E 69 DMDAElementType - Defines the type of elements that will be returned by 70 DMDAGetElements() 71 72 Level: beginner 73 74 .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMCreateInterpolation(), DMDASetInterpolationType(), 75 DMDASetElementType(), DMDAGetElements(), DMDARestoreElements(), DMDACreate() 76 E*/ 77 typedef enum { DMDA_ELEMENT_P1, DMDA_ELEMENT_Q1 } DMDAElementType; 78 79 PETSC_EXTERN PetscErrorCode DMDASetElementType(DM,DMDAElementType); 80 PETSC_EXTERN PetscErrorCode DMDAGetElementType(DM,DMDAElementType*); 81 PETSC_EXTERN PetscErrorCode DMDAGetElements(DM,PetscInt *,PetscInt *,const PetscInt*[]); 82 PETSC_EXTERN PetscErrorCode DMDARestoreElements(DM,PetscInt *,PetscInt *,const PetscInt*[]); 83 84 typedef enum { DMDA_X,DMDA_Y,DMDA_Z } DMDADirection; 85 86 #define MATSEQUSFFT "sequsfft" 87 88 PETSC_EXTERN PetscErrorCode DMDACreate(MPI_Comm,DM*); 89 PETSC_EXTERN PetscErrorCode DMDASetDim(DM,PetscInt); 90 PETSC_EXTERN PetscErrorCode DMDASetSizes(DM,PetscInt,PetscInt,PetscInt); 91 PETSC_EXTERN PetscErrorCode DMDACreate1d(MPI_Comm,DMDABoundaryType,PetscInt,PetscInt,PetscInt,const PetscInt[],DM *); 92 PETSC_EXTERN PetscErrorCode DMDACreate2d(MPI_Comm,DMDABoundaryType,DMDABoundaryType,DMDAStencilType,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],DM*); 93 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*); 94 95 PETSC_EXTERN PetscErrorCode DMDAGlobalToNaturalBegin(DM,Vec,InsertMode,Vec); 96 PETSC_EXTERN PetscErrorCode DMDAGlobalToNaturalEnd(DM,Vec,InsertMode,Vec); 97 PETSC_EXTERN PetscErrorCode DMDANaturalToGlobalBegin(DM,Vec,InsertMode,Vec); 98 PETSC_EXTERN PetscErrorCode DMDANaturalToGlobalEnd(DM,Vec,InsertMode,Vec); 99 PETSC_EXTERN PetscErrorCode DMDALocalToLocalBegin(DM,Vec,InsertMode,Vec); 100 PETSC_EXTERN PetscErrorCode DMDALocalToLocalEnd(DM,Vec,InsertMode,Vec); 101 PETSC_EXTERN PetscErrorCode DMDACreateNaturalVector(DM,Vec *); 102 103 PETSC_EXTERN PetscErrorCode DMDAGetCorners(DM,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*); 104 PETSC_EXTERN PetscErrorCode DMDAGetGhostCorners(DM,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*); 105 PETSC_EXTERN PetscErrorCode DMDAGetInfo(DM,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,DMDABoundaryType*,DMDABoundaryType*,DMDABoundaryType*,DMDAStencilType*); 106 PETSC_EXTERN PetscErrorCode DMDAGetProcessorSubset(DM,DMDADirection,PetscInt,MPI_Comm*); 107 PETSC_EXTERN PetscErrorCode DMDAGetProcessorSubsets(DM,DMDADirection,MPI_Comm*); 108 PETSC_EXTERN PetscErrorCode DMDAGetRay(DM,DMDADirection,PetscInt,Vec*,VecScatter*); 109 110 PETSC_EXTERN PetscErrorCode DMDAGlobalToNaturalAllCreate(DM,VecScatter*); 111 PETSC_EXTERN PetscErrorCode DMDANaturalAllToGlobalCreate(DM,VecScatter*); 112 113 PETSC_EXTERN PetscErrorCode DMDAGetGlobalIndices(DM,PetscInt*,PetscInt**); 114 115 PETSC_EXTERN PetscErrorCode DMDAGetScatter(DM,VecScatter*,VecScatter*,VecScatter*); 116 PETSC_EXTERN PetscErrorCode DMDAGetNeighbors(DM,const PetscMPIInt**); 117 118 PETSC_EXTERN PetscErrorCode DMDAGetAO(DM,AO*); 119 PETSC_EXTERN PetscErrorCode DMDASetUniformCoordinates(DM,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal); 120 PETSC_EXTERN PetscErrorCode DMDAGetBoundingBox(DM,PetscReal[],PetscReal[]); 121 PETSC_EXTERN PetscErrorCode DMDAGetLocalBoundingBox(DM,PetscReal[],PetscReal[]); 122 /* function to wrap coordinates around boundary */ 123 PETSC_EXTERN PetscErrorCode DMDAMapCoordsToPeriodicDomain(DM,PetscScalar*,PetscScalar*); 124 125 PETSC_EXTERN PetscErrorCode DMDAGetReducedDMDA(DM,PetscInt,DM*); 126 127 PETSC_EXTERN PetscErrorCode DMDASetFieldName(DM,PetscInt,const char[]); 128 PETSC_EXTERN PetscErrorCode DMDAGetFieldName(DM,PetscInt,const char**); 129 PETSC_EXTERN PetscErrorCode DMDASetCoordinateName(DM,PetscInt,const char[]); 130 PETSC_EXTERN PetscErrorCode DMDAGetCoordinateName(DM,PetscInt,const char**); 131 132 PETSC_EXTERN PetscErrorCode DMDASetBoundaryType(DM,DMDABoundaryType,DMDABoundaryType,DMDABoundaryType); 133 PETSC_EXTERN PetscErrorCode DMDASetDof(DM, PetscInt); 134 PETSC_EXTERN PetscErrorCode DMDASetOverlap(DM, PetscInt); 135 PETSC_EXTERN PetscErrorCode DMDAGetOffset(DM, PetscInt*,PetscInt*,PetscInt*); 136 PETSC_EXTERN PetscErrorCode DMDASetOffset(DM, PetscInt,PetscInt,PetscInt); 137 PETSC_EXTERN PetscErrorCode DMDASetStencilWidth(DM, PetscInt); 138 PETSC_EXTERN PetscErrorCode DMDASetOwnershipRanges(DM,const PetscInt[],const PetscInt[],const PetscInt[]); 139 PETSC_EXTERN PetscErrorCode DMDAGetOwnershipRanges(DM,const PetscInt**,const PetscInt**,const PetscInt**); 140 PETSC_EXTERN PetscErrorCode DMDASetNumProcs(DM, PetscInt, PetscInt, PetscInt); 141 PETSC_EXTERN PetscErrorCode DMDASetStencilType(DM, DMDAStencilType); 142 143 PETSC_EXTERN PetscErrorCode DMDAVecGetArray(DM,Vec,void *); 144 PETSC_EXTERN PetscErrorCode DMDAVecRestoreArray(DM,Vec,void *); 145 146 PETSC_EXTERN PetscErrorCode DMDAVecGetArrayDOF(DM,Vec,void *); 147 PETSC_EXTERN PetscErrorCode DMDAVecRestoreArrayDOF(DM,Vec,void *); 148 149 PETSC_EXTERN PetscErrorCode DMDASplitComm2d(MPI_Comm,PetscInt,PetscInt,PetscInt,MPI_Comm*); 150 151 PETSC_EXTERN PetscErrorCode DMDACreatePatchIS(DM,MatStencil*,MatStencil*,IS*); 152 153 /*S 154 DMDALocalInfo - C struct that contains information about a structured grid and a processors logical 155 location in it. 156 157 Level: beginner 158 159 Concepts: distributed array 160 161 Developer note: Then entries in this struct are int instead of PetscInt so that the elements may 162 be extracted in Fortran as if from an integer array 163 164 .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DM, DMDAGetLocalInfo(), DMDAGetInfo() 165 S*/ 166 typedef struct { 167 PetscInt dim,dof,sw; 168 PetscInt mx,my,mz; /* global number of grid points in each direction */ 169 PetscInt xs,ys,zs; /* starting point of this processor, excluding ghosts */ 170 PetscInt xm,ym,zm; /* number of grid points on this processor, excluding ghosts */ 171 PetscInt gxs,gys,gzs; /* starting point of this processor including ghosts */ 172 PetscInt gxm,gym,gzm; /* number of grid points on this processor including ghosts */ 173 DMDABoundaryType bx,by,bz; /* type of ghost nodes at boundary */ 174 DMDAStencilType st; 175 DM da; 176 } DMDALocalInfo; 177 178 /*MC 179 DMDAForEachPointBegin2d - Starts a loop over the local part of a two dimensional DMDA 180 181 Synopsis: 182 void DMDAForEachPointBegin2d(DALocalInfo *info,PetscInt i,PetscInt j); 183 184 Not Collective 185 186 Level: intermediate 187 188 .seealso: DMDAForEachPointEnd2d(), DMDAVecGetArray() 189 M*/ 190 #define DMDAForEachPointBegin2d(info,i,j) {\ 191 PetscInt _xints = info->xs,_xinte = info->xs+info->xm,_yints = info->ys,_yinte = info->ys+info->ym;\ 192 for (j=_yints; j<_yinte; j++) {\ 193 for (i=_xints; i<_xinte; i++) {\ 194 195 /*MC 196 DMDAForEachPointEnd2d - Ends a loop over the local part of a two dimensional DMDA 197 198 Synopsis: 199 void DMDAForEachPointEnd2d; 200 201 Not Collective 202 203 Level: intermediate 204 205 .seealso: DMDAForEachPointBegin2d(), DMDAVecGetArray() 206 M*/ 207 #define DMDAForEachPointEnd2d }}} 208 209 /*MC 210 DMDACoor2d - Structure for holding 2d (x and y) coordinates. 211 212 Level: intermediate 213 214 Sample Usage: 215 DMDACoor2d **coors; 216 Vec vcoors; 217 DM cda; 218 219 DMGetCoordinates(da,&vcoors); 220 DMDAGetCoordinateDA(da,&cda); 221 DMDAVecGetArray(cda,vcoors,&coors); 222 DMDAGetCorners(cda,&mstart,&nstart,0,&m,&n,0) 223 for (i=mstart; i<mstart+m; i++) { 224 for (j=nstart; j<nstart+n; j++) { 225 x = coors[j][i].x; 226 y = coors[j][i].y; 227 ...... 228 } 229 } 230 DMDAVecRestoreArray(dac,vcoors,&coors); 231 232 .seealso: DMDACoor3d, DMDAForEachPointBegin(), DMDAGetCoordinateDA(), DMGetCoordinates(), DMDAGetGhostCoordinates() 233 M*/ 234 typedef struct {PetscScalar x,y;} DMDACoor2d; 235 236 /*MC 237 DMDACoor3d - Structure for holding 3d (x, y and z) coordinates. 238 239 Level: intermediate 240 241 Sample Usage: 242 DMDACoor3d ***coors; 243 Vec vcoors; 244 DM cda; 245 246 DMGetCoordinates(da,&vcoors); 247 DMDAGetCoordinateDA(da,&cda); 248 DMDAVecGetArray(cda,vcoors,&coors); 249 DMDAGetCorners(cda,&mstart,&nstart,&pstart,&m,&n,&p) 250 for (i=mstart; i<mstart+m; i++) { 251 for (j=nstart; j<nstart+n; j++) { 252 for (k=pstart; k<pstart+p; k++) { 253 x = coors[k][j][i].x; 254 y = coors[k][j][i].y; 255 z = coors[k][j][i].z; 256 ...... 257 } 258 } 259 DMDAVecRestoreArray(dac,vcoors,&coors); 260 261 .seealso: DMDACoor2d, DMDAForEachPointBegin(), DMDAGetCoordinateDA(), DMGetCoordinates(), DMDAGetGhostCoordinates() 262 M*/ 263 typedef struct {PetscScalar x,y,z;} DMDACoor3d; 264 265 PETSC_EXTERN PetscErrorCode DMDAGetLocalInfo(DM,DMDALocalInfo*); 266 267 PETSC_EXTERN PetscErrorCode MatRegisterDAAD(void); 268 PETSC_EXTERN PetscErrorCode MatCreateDAAD(DM,Mat*); 269 PETSC_EXTERN PetscErrorCode MatCreateSeqUSFFT(Vec,DM,Mat*); 270 271 PETSC_EXTERN PetscErrorCode DMDASetGetMatrix(DM,PetscErrorCode (*)(DM, MatType,Mat *)); 272 PETSC_EXTERN PetscErrorCode DMDASetBlockFills(DM,const PetscInt*,const PetscInt*); 273 PETSC_EXTERN PetscErrorCode DMDASetRefinementFactor(DM,PetscInt,PetscInt,PetscInt); 274 PETSC_EXTERN PetscErrorCode DMDAGetRefinementFactor(DM,PetscInt*,PetscInt*,PetscInt*); 275 276 PETSC_EXTERN PetscErrorCode DMDAGetArray(DM,PetscBool ,void*); 277 PETSC_EXTERN PetscErrorCode DMDARestoreArray(DM,PetscBool ,void*); 278 279 PETSC_EXTERN PetscErrorCode DMDACreatePF(DM,PF*); 280 281 PETSC_EXTERN PetscErrorCode DMDAGetNumCells(DM, PetscInt *); 282 PETSC_EXTERN PetscErrorCode DMDAGetNumVertices(DM, PetscInt *, PetscInt *, PetscInt *, PetscInt *); 283 PETSC_EXTERN PetscErrorCode DMDAGetNumFaces(DM, PetscInt *, PetscInt *, PetscInt *, PetscInt *, PetscInt *, PetscInt *); 284 PETSC_EXTERN PetscErrorCode DMDAGetHeightStratum(DM, PetscInt, PetscInt *, PetscInt *); 285 PETSC_EXTERN PetscErrorCode DMDACreateSection(DM, PetscInt[], PetscInt[], PetscInt[], PetscInt[]); 286 PETSC_EXTERN PetscErrorCode DMDAComputeCellGeometry(DM, PetscInt, PetscQuadrature *, PetscReal [], PetscReal [], PetscReal [], PetscReal []); 287 PETSC_EXTERN PetscErrorCode DMDAVecGetClosure(DM, PetscSection, Vec, PetscInt, const PetscScalar **); 288 PETSC_EXTERN PetscErrorCode DMDAVecSetClosure(DM, PetscSection, Vec, PetscInt, const PetscScalar *, InsertMode); 289 PETSC_EXTERN PetscErrorCode DMDAGetClosure(DM,PetscSection,PetscInt,PetscInt*,const PetscInt**); 290 PETSC_EXTERN PetscErrorCode DMDARestoreClosure(DM,PetscSection,PetscInt,PetscInt*,const PetscInt**); 291 PETSC_EXTERN PetscErrorCode DMDAGetClosureScalar(DM,PetscSection,PetscInt,PetscScalar*,const PetscScalar**); 292 PETSC_EXTERN PetscErrorCode DMDARestoreClosureScalar(DM,PetscSection,PetscInt,PetscScalar*,const PetscScalar**); 293 PETSC_EXTERN PetscErrorCode DMDASetClosureScalar(DM,PetscSection,PetscInt,PetscScalar*,const PetscScalar*,InsertMode); 294 PETSC_EXTERN PetscErrorCode DMDAConvertToCell(DM, MatStencil, PetscInt *); 295 296 #endif 297