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