1 #if !defined(__PETSCDMDA_H) 2 #define __PETSCDMDA_H 3 4 #include "petscdm.h" 5 #include "petscpf.h" 6 #include "petscao.h" 7 PETSC_EXTERN_CXX_BEGIN 8 9 /*E 10 DMDAStencilType - Determines if the stencil extends only along the coordinate directions, or also 11 to the northeast, northwest etc 12 13 Level: beginner 14 15 .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDACreate() 16 E*/ 17 typedef enum { DMDA_STENCIL_STAR,DMDA_STENCIL_BOX } DMDAStencilType; 18 19 /*MC 20 DMDA_STENCIL_STAR - "Star"-type stencil. In logical grid coordinates, only (i,j,k), (i+s,j,k), (i,j+s,k), 21 (i,j,k+s) are in the stencil NOT, for example, (i+s,j+s,k) 22 23 Level: beginner 24 25 .seealso: DMDA_STENCIL_BOX, DMDAStencilType 26 M*/ 27 28 /*MC 29 DMDA_STENCIL_BOX - "Box"-type stencil. In logical grid coordinates, any of (i,j,k), (i+s,j+r,k+t) may 30 be in the stencil. 31 32 Level: beginner 33 34 .seealso: DMDA_STENCIL_STAR, DMDAStencilType 35 M*/ 36 37 /*E 38 DMDABoundaryType - Describes the choice for fill of ghost cells on domain boundaries. 39 40 Level: beginner 41 42 A boundary may be of type DMDA_BOUNDARY_NONE (no ghost nodes), DMDA_BOUNDARY_GHOST (ghost nodes 43 exist but aren't filled), DMDA_BOUNDARY_MIRROR (not yet implemented), or DMDA_BOUNDARY_PERIODIC 44 (ghost nodes filled by the opposite edge of the domain). 45 46 .seealso: DMDASetBoundaryType(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDACreate() 47 E*/ 48 typedef enum { DMDA_BOUNDARY_NONE, DMDA_BOUNDARY_GHOSTED, DMDA_BOUNDARY_MIRROR, DMDA_BOUNDARY_PERIODIC } DMDABoundaryType; 49 50 extern const char *DMDABoundaryTypes[]; 51 52 /*E 53 DMDAInterpolationType - Defines the type of interpolation that will be returned by 54 DMGetInterpolation. 55 56 Level: beginner 57 58 .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMGetInterpolation(), DMDASetInterpolationType(), DMDACreate() 59 E*/ 60 typedef enum { DMDA_Q0, DMDA_Q1 } DMDAInterpolationType; 61 62 extern PetscErrorCode DMDASetInterpolationType(DM,DMDAInterpolationType); 63 extern PetscErrorCode DMDAGetInterpolationType(DM,DMDAInterpolationType*); 64 65 /*E 66 DMDAElementType - Defines the type of elements that will be returned by 67 DMDAGetElements() 68 69 Level: beginner 70 71 .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMGetInterpolation(), DMDASetInterpolationType(), 72 DMDASetElementType(), DMDAGetElements(), DMDARestoreElements(), DMDACreate() 73 E*/ 74 typedef enum { DMDA_ELEMENT_P1, DMDA_ELEMENT_Q1 } DMDAElementType; 75 76 extern PetscErrorCode DMDASetElementType(DM,DMDAElementType); 77 extern PetscErrorCode DMDAGetElementType(DM,DMDAElementType*); 78 extern PetscErrorCode DMDAGetElements(DM,PetscInt *,PetscInt *,const PetscInt*[]); 79 extern PetscErrorCode DMDARestoreElements(DM,PetscInt *,PetscInt *,const PetscInt*[]); 80 81 typedef enum { DMDA_X,DMDA_Y,DMDA_Z } DMDADirection; 82 83 #define MATSEQUSFFT "sequsfft" 84 85 extern PetscErrorCode DMDACreate(MPI_Comm,DM*); 86 extern PetscErrorCode DMDASetDim(DM,PetscInt); 87 extern PetscErrorCode DMDASetSizes(DM,PetscInt,PetscInt,PetscInt); 88 extern PetscErrorCode DMDACreate1d(MPI_Comm,DMDABoundaryType,PetscInt,PetscInt,PetscInt,const PetscInt[],DM *); 89 extern PetscErrorCode DMDACreate2d(MPI_Comm,DMDABoundaryType,DMDABoundaryType,DMDAStencilType,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],DM*); 90 extern PetscErrorCode DMDACreate3d(MPI_Comm,DMDABoundaryType,DMDABoundaryType,DMDABoundaryType,DMDAStencilType,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscInt[],DM*); 91 92 extern PetscErrorCode DMDAGlobalToNaturalBegin(DM,Vec,InsertMode,Vec); 93 extern PetscErrorCode DMDAGlobalToNaturalEnd(DM,Vec,InsertMode,Vec); 94 extern PetscErrorCode DMDANaturalToGlobalBegin(DM,Vec,InsertMode,Vec); 95 extern PetscErrorCode DMDANaturalToGlobalEnd(DM,Vec,InsertMode,Vec); 96 extern PetscErrorCode DMDALocalToLocalBegin(DM,Vec,InsertMode,Vec); 97 extern PetscErrorCode DMDALocalToLocalEnd(DM,Vec,InsertMode,Vec); 98 extern PetscErrorCode DMDACreateNaturalVector(DM,Vec *); 99 100 extern PetscErrorCode DMDALoad(PetscViewer,PetscInt,PetscInt,PetscInt,DM *); 101 extern PetscErrorCode DMDAGetCorners(DM,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*); 102 extern PetscErrorCode DMDAGetGhostCorners(DM,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*); 103 extern PetscErrorCode DMDAGetInfo(DM,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,DMDABoundaryType*,DMDABoundaryType*,DMDABoundaryType*,DMDAStencilType*); 104 extern PetscErrorCode DMDAGetProcessorSubset(DM,DMDADirection,PetscInt,MPI_Comm*); 105 extern PetscErrorCode DMDAGetProcessorSubsets(DM,DMDADirection,MPI_Comm*); 106 107 extern PetscErrorCode DMDAGlobalToNaturalAllCreate(DM,VecScatter*); 108 extern PetscErrorCode DMDANaturalAllToGlobalCreate(DM,VecScatter*); 109 110 extern PetscErrorCode DMDAGetGlobalIndices(DM,PetscInt*,PetscInt**); 111 112 extern PetscErrorCode DMDAGetScatter(DM,VecScatter*,VecScatter*,VecScatter*); 113 extern PetscErrorCode DMDAGetNeighbors(DM,const PetscMPIInt**); 114 115 extern PetscErrorCode DMDAGetAO(DM,AO*); 116 extern PetscErrorCode DMDASetCoordinates(DM,Vec); 117 extern PetscErrorCode DMDASetGhostedCoordinates(DM,Vec); 118 extern PetscErrorCode DMDAGetCoordinates(DM,Vec *); 119 extern PetscErrorCode DMDAGetGhostedCoordinates(DM,Vec *); 120 extern PetscErrorCode DMDAGetCoordinateDA(DM,DM *); 121 extern PetscErrorCode DMDASetUniformCoordinates(DM,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal); 122 extern PetscErrorCode DMDAGetBoundingBox(DM,PetscReal[],PetscReal[]); 123 extern PetscErrorCode DMDAGetLocalBoundingBox(DM,PetscReal[],PetscReal[]); 124 /* function to wrap coordinates around boundary */ 125 extern PetscErrorCode DMDAMapCoordsToPeriodicDomain(DM,PetscScalar*,PetscScalar*); 126 127 extern PetscErrorCode DMDAGetReducedDA(DM,PetscInt,DM*); 128 129 extern PetscErrorCode DMDASetFieldName(DM,PetscInt,const char[]); 130 extern PetscErrorCode DMDAGetFieldName(DM,PetscInt,const char**); 131 132 extern PetscErrorCode DMDASetBoundaryType(DM,DMDABoundaryType,DMDABoundaryType,DMDABoundaryType); 133 extern PetscErrorCode DMDASetDof(DM, PetscInt); 134 extern PetscErrorCode DMDASetStencilWidth(DM, PetscInt); 135 extern PetscErrorCode DMDASetOwnershipRanges(DM,const PetscInt[],const PetscInt[],const PetscInt[]); 136 extern PetscErrorCode DMDAGetOwnershipRanges(DM,const PetscInt**,const PetscInt**,const PetscInt**); 137 extern PetscErrorCode DMDASetNumProcs(DM, PetscInt, PetscInt, PetscInt); 138 extern PetscErrorCode DMDASetStencilType(DM, DMDAStencilType); 139 140 extern PetscErrorCode DMDAVecGetArray(DM,Vec,void *); 141 extern PetscErrorCode DMDAVecRestoreArray(DM,Vec,void *); 142 143 extern PetscErrorCode DMDAVecGetArrayDOF(DM,Vec,void *); 144 extern PetscErrorCode DMDAVecRestoreArrayDOF(DM,Vec,void *); 145 146 extern PetscErrorCode DMDASplitComm2d(MPI_Comm,PetscInt,PetscInt,PetscInt,MPI_Comm*); 147 148 /*S 149 DMDALocalInfo - C struct that contains information about a structured grid and a processors logical 150 location in it. 151 152 Level: beginner 153 154 Concepts: distributed array 155 156 Developer note: Then entries in this struct are int instead of PetscInt so that the elements may 157 be extracted in Fortran as if from an integer array 158 159 .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DM, DMDAGetLocalInfo(), DMDAGetInfo() 160 S*/ 161 typedef struct { 162 PetscInt dim,dof,sw; 163 PetscInt mx,my,mz; /* global number of grid points in each direction */ 164 PetscInt xs,ys,zs; /* starting point of this processor, excluding ghosts */ 165 PetscInt xm,ym,zm; /* number of grid points on this processor, excluding ghosts */ 166 PetscInt gxs,gys,gzs; /* starting point of this processor including ghosts */ 167 PetscInt gxm,gym,gzm; /* number of grid points on this processor including ghosts */ 168 DMDABoundaryType bx,by,bz; /* type of ghost nodes at boundary */ 169 DMDAStencilType st; 170 DM da; 171 } DMDALocalInfo; 172 173 /*MC 174 DMDAForEachPointBegin2d - Starts a loop over the local part of a two dimensional DMDA 175 176 Synopsis: 177 void DMDAForEachPointBegin2d(DALocalInfo *info,PetscInt i,PetscInt j); 178 179 Not Collective 180 181 Level: intermediate 182 183 .seealso: DMDAForEachPointEnd2d(), DMDAVecGetArray() 184 M*/ 185 #define DMDAForEachPointBegin2d(info,i,j) {\ 186 PetscInt _xints = info->xs,_xinte = info->xs+info->xm,_yints = info->ys,_yinte = info->ys+info->ym;\ 187 for (j=_yints; j<_yinte; j++) {\ 188 for (i=_xints; i<_xinte; i++) {\ 189 190 /*MC 191 DMDAForEachPointEnd2d - Ends a loop over the local part of a two dimensional DMDA 192 193 Synopsis: 194 void DMDAForEachPointEnd2d; 195 196 Not Collective 197 198 Level: intermediate 199 200 .seealso: DMDAForEachPointBegin2d(), DMDAVecGetArray() 201 M*/ 202 #define DMDAForEachPointEnd2d }}} 203 204 /*MC 205 DMDACoor2d - Structure for holding 2d (x and y) coordinates. 206 207 Level: intermediate 208 209 Sample Usage: 210 DMDACoor2d **coors; 211 Vec vcoors; 212 DM cda; 213 214 DMDAGetCoordinates(da,&vcoors); 215 DMDAGetCoordinateDA(da,&cda); 216 DMDAVecGetArray(cda,vcoors,&coors); 217 DMDAGetCorners(cda,&mstart,&nstart,0,&m,&n,0) 218 for (i=mstart; i<mstart+m; i++) { 219 for (j=nstart; j<nstart+n; j++) { 220 x = coors[j][i].x; 221 y = coors[j][i].y; 222 ...... 223 } 224 } 225 DMDAVecRestoreArray(dac,vcoors,&coors); 226 227 .seealso: DMDACoor3d, DMDAForEachPointBegin(), DMDAGetCoordinateDA(), DMDAGetCoordinates(), DMDAGetGhostCoordinates() 228 M*/ 229 typedef struct {PetscScalar x,y;} DMDACoor2d; 230 231 /*MC 232 DMDACoor3d - Structure for holding 3d (x, y and z) coordinates. 233 234 Level: intermediate 235 236 Sample Usage: 237 DMDACoor3d ***coors; 238 Vec vcoors; 239 DM cda; 240 241 DMDAGetCoordinates(da,&vcoors); 242 DMDAGetCoordinateDA(da,&cda); 243 DMDAVecGetArray(cda,vcoors,&coors); 244 DMDAGetCorners(cda,&mstart,&nstart,&pstart,&m,&n,&p) 245 for (i=mstart; i<mstart+m; i++) { 246 for (j=nstart; j<nstart+n; j++) { 247 for (k=pstart; k<pstart+p; k++) { 248 x = coors[k][j][i].x; 249 y = coors[k][j][i].y; 250 z = coors[k][j][i].z; 251 ...... 252 } 253 } 254 DMDAVecRestoreArray(dac,vcoors,&coors); 255 256 .seealso: DMDACoor2d, DMDAForEachPointBegin(), DMDAGetCoordinateDA(), DMDAGetCoordinates(), DMDAGetGhostCoordinates() 257 M*/ 258 typedef struct {PetscScalar x,y,z;} DMDACoor3d; 259 260 extern PetscErrorCode DMDAGetLocalInfo(DM,DMDALocalInfo*); 261 typedef PetscErrorCode (*DMDALocalFunction1)(DMDALocalInfo*,void*,void*,void*); 262 extern PetscErrorCode DMDAFormFunctionLocal(DM, DMDALocalFunction1, Vec, Vec, void *); 263 extern PetscErrorCode DMDAFormFunctionLocalGhost(DM, DMDALocalFunction1, Vec, Vec, void *); 264 extern PetscErrorCode DMDAFormJacobianLocal(DM, DMDALocalFunction1, Vec, Mat, void *); 265 extern PetscErrorCode DMDAFormFunction1(DM,Vec,Vec,void*); 266 extern PetscErrorCode DMDAFormFunction(DM,PetscErrorCode (*)(void),Vec,Vec,void*); 267 extern PetscErrorCode DMDAFormFunctioni1(DM,PetscInt,Vec,PetscScalar*,void*); 268 extern PetscErrorCode DMDAFormFunctionib1(DM,PetscInt,Vec,PetscScalar*,void*); 269 extern PetscErrorCode DMDAComputeJacobian1WithAdic(DM,Vec,Mat,void*); 270 extern PetscErrorCode DMDAComputeJacobian1WithAdifor(DM,Vec,Mat,void*); 271 extern PetscErrorCode DMDAMultiplyByJacobian1WithAdic(DM,Vec,Vec,Vec,void*); 272 extern PetscErrorCode DMDAMultiplyByJacobian1WithAdifor(DM,Vec,Vec,Vec,void*); 273 extern PetscErrorCode DMDAMultiplyByJacobian1WithAD(DM,Vec,Vec,Vec,void*); 274 extern PetscErrorCode DMDAComputeJacobian1(DM,Vec,Mat,void*); 275 extern PetscErrorCode DMDAGetLocalFunction(DM,DMDALocalFunction1*); 276 extern PetscErrorCode DMDASetLocalFunction(DM,DMDALocalFunction1); 277 extern PetscErrorCode DMDASetLocalFunctioni(DM,PetscErrorCode (*)(DMDALocalInfo*,MatStencil*,void*,PetscScalar*,void*)); 278 extern PetscErrorCode DMDASetLocalFunctionib(DM,PetscErrorCode (*)(DMDALocalInfo*,MatStencil*,void*,PetscScalar*,void*)); 279 extern PetscErrorCode DMDAGetLocalJacobian(DM,DMDALocalFunction1*); 280 extern PetscErrorCode DMDASetLocalJacobian(DM,DMDALocalFunction1); 281 extern PetscErrorCode DMDASetLocalAdicFunction_Private(DM,DMDALocalFunction1); 282 283 extern PetscErrorCode MatSetDM(Mat,DM); 284 extern PetscErrorCode MatRegisterDAAD(void); 285 extern PetscErrorCode MatCreateDAAD(DM,Mat*); 286 extern PetscErrorCode MatCreateSeqUSFFT(Vec,DM,Mat*); 287 288 /*MC 289 DMDASetLocalAdicFunction - Caches in a DM a local function computed by ADIC/ADIFOR 290 291 Synopsis: 292 PetscErrorCode DMDASetLocalAdicFunction(DM da,DMDALocalFunction1 ad_lf) 293 294 Logically Collective on DM 295 296 Input Parameter: 297 + da - initial distributed array 298 - ad_lf - the local function as computed by ADIC/ADIFOR 299 300 Level: intermediate 301 302 .keywords: distributed array, refine 303 304 .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDAGetLocalFunction(), DMDASetLocalFunction(), 305 DMDASetLocalJacobian() 306 M*/ 307 #if defined(PETSC_HAVE_ADIC) 308 # define DMDASetLocalAdicFunction(a,d) DMDASetLocalAdicFunction_Private(a,(DMDALocalFunction1)d) 309 #else 310 # define DMDASetLocalAdicFunction(a,d) DMDASetLocalAdicFunction_Private(a,0) 311 #endif 312 313 extern PetscErrorCode DMDASetLocalAdicMFFunction_Private(DM,DMDALocalFunction1); 314 #if defined(PETSC_HAVE_ADIC) 315 # define DMDASetLocalAdicMFFunction(a,d) DMDASetLocalAdicMFFunction_Private(a,(DMDALocalFunction1)d) 316 #else 317 # define DMDASetLocalAdicMFFunction(a,d) DMDASetLocalAdicMFFunction_Private(a,0) 318 #endif 319 extern PetscErrorCode DMDASetLocalAdicFunctioni_Private(DM,PetscErrorCode (*)(DMDALocalInfo*,MatStencil*,void*,void*,void*)); 320 #if defined(PETSC_HAVE_ADIC) 321 # define DMDASetLocalAdicFunctioni(a,d) DMDASetLocalAdicFunctioni_Private(a,(PetscErrorCode (*)(DMDALocalInfo*,MatStencil*,void*,void*,void*))d) 322 #else 323 # define DMDASetLocalAdicFunctioni(a,d) DMDASetLocalAdicFunctioni_Private(a,0) 324 #endif 325 extern PetscErrorCode DMDASetLocalAdicMFFunctioni_Private(DM,PetscErrorCode (*)(DMDALocalInfo*,MatStencil*,void*,void*,void*)); 326 #if defined(PETSC_HAVE_ADIC) 327 # define DMDASetLocalAdicMFFunctioni(a,d) DMDASetLocalAdicMFFunctioni_Private(a,(PetscErrorCode (*)(DMDALocalInfo*,MatStencil*,void*,void*,void*))d) 328 #else 329 # define DMDASetLocalAdicMFFunctioni(a,d) DMDASetLocalAdicMFFunctioni_Private(a,0) 330 #endif 331 332 extern PetscErrorCode DMDASetLocalAdicFunctionib_Private(DM,PetscErrorCode (*)(DMDALocalInfo*,MatStencil*,void*,void*,void*)); 333 #if defined(PETSC_HAVE_ADIC) 334 # define DMDASetLocalAdicFunctionib(a,d) DMDASetLocalAdicFunctionib_Private(a,(PetscErrorCode (*)(DMDALocalInfo*,MatStencil*,void*,void*,void*))d) 335 #else 336 # define DMDASetLocalAdicFunctionib(a,d) DMDASetLocalAdicFunctionib_Private(a,0) 337 #endif 338 extern PetscErrorCode DMDASetLocalAdicMFFunctionib_Private(DM,PetscErrorCode (*)(DMDALocalInfo*,MatStencil*,void*,void*,void*)); 339 #if defined(PETSC_HAVE_ADIC) 340 # define DMDASetLocalAdicMFFunctionib(a,d) DMDASetLocalAdicMFFunctionib_Private(a,(PetscErrorCode (*)(DMDALocalInfo*,MatStencil*,void*,void*,void*))d) 341 #else 342 # define DMDASetLocalAdicMFFunctionib(a,d) DMDASetLocalAdicMFFunctionib_Private(a,0) 343 #endif 344 345 extern PetscErrorCode DMDAFormFunctioniTest1(DM,void*); 346 extern PetscErrorCode DMDASetGetMatrix(DM,PetscErrorCode (*)(DM, const MatType,Mat *)); 347 extern PetscErrorCode DMDASetBlockFills(DM,PetscInt*,PetscInt*); 348 extern PetscErrorCode DMDASetRefinementFactor(DM,PetscInt,PetscInt,PetscInt); 349 extern PetscErrorCode DMDAGetRefinementFactor(DM,PetscInt*,PetscInt*,PetscInt*); 350 351 extern PetscErrorCode DMDAGetAdicArray(DM,PetscBool ,void*,void*,PetscInt*); 352 extern PetscErrorCode DMDARestoreAdicArray(DM,PetscBool ,void*,void*,PetscInt*); 353 extern PetscErrorCode DMDAGetAdicMFArray(DM,PetscBool ,void*,void*,PetscInt*); 354 extern PetscErrorCode DMDAGetAdicMFArray4(DM,PetscBool ,void*,void*,PetscInt*); 355 extern PetscErrorCode DMDAGetAdicMFArray9(DM,PetscBool ,void*,void*,PetscInt*); 356 extern PetscErrorCode DMDAGetAdicMFArrayb(DM,PetscBool ,void*,void*,PetscInt*); 357 extern PetscErrorCode DMDARestoreAdicMFArray(DM,PetscBool ,void*,void*,PetscInt*); 358 extern PetscErrorCode DMDAGetArray(DM,PetscBool ,void*); 359 extern PetscErrorCode DMDARestoreArray(DM,PetscBool ,void*); 360 extern PetscErrorCode ad_DAGetArray(DM,PetscBool ,void*); 361 extern PetscErrorCode ad_DARestoreArray(DM,PetscBool ,void*); 362 extern PetscErrorCode admf_DAGetArray(DM,PetscBool ,void*); 363 extern PetscErrorCode admf_DARestoreArray(DM,PetscBool ,void*); 364 365 extern PetscErrorCode DMDACreatePF(DM,PF*); 366 367 #define DMDA_FILE_CLASSID 1211220 368 PETSC_EXTERN_CXX_END 369 #endif 370