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() 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 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 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), 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 109 PETSC_EXTERN PetscErrorCode DMDAGlobalToNaturalAllCreate(DM,VecScatter*); 110 PETSC_EXTERN PetscErrorCode DMDANaturalAllToGlobalCreate(DM,VecScatter*); 111 112 PETSC_EXTERN PetscErrorCode DMDAGetGlobalIndices(DM,PetscInt*,PetscInt**); 113 114 PETSC_EXTERN PetscErrorCode DMDAGetScatter(DM,VecScatter*,VecScatter*,VecScatter*); 115 PETSC_EXTERN PetscErrorCode DMDAGetNeighbors(DM,const PetscMPIInt**); 116 117 PETSC_EXTERN PetscErrorCode DMDAGetAO(DM,AO*); 118 PETSC_EXTERN PetscErrorCode DMDASetCoordinates(DM,Vec); 119 PETSC_EXTERN PetscErrorCode DMDASetGhostedCoordinates(DM,Vec); 120 PETSC_EXTERN PetscErrorCode DMDAGetCoordinates(DM,Vec *); 121 PETSC_EXTERN PetscErrorCode DMDAGetGhostedCoordinates(DM,Vec *); 122 PETSC_EXTERN PetscErrorCode DMDAGetCoordinateDA(DM,DM *); 123 PETSC_EXTERN PetscErrorCode DMDASetUniformCoordinates(DM,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal); 124 PETSC_EXTERN PetscErrorCode DMDAGetBoundingBox(DM,PetscReal[],PetscReal[]); 125 PETSC_EXTERN PetscErrorCode DMDAGetLocalBoundingBox(DM,PetscReal[],PetscReal[]); 126 /* function to wrap coordinates around boundary */ 127 PETSC_EXTERN PetscErrorCode DMDAMapCoordsToPeriodicDomain(DM,PetscScalar*,PetscScalar*); 128 129 PETSC_EXTERN PetscErrorCode DMDAGetReducedDA(DM,PetscInt,DM*); 130 131 PETSC_EXTERN PetscErrorCode DMDASetFieldName(DM,PetscInt,const char[]); 132 PETSC_EXTERN PetscErrorCode DMDAGetFieldName(DM,PetscInt,const char**); 133 134 PETSC_EXTERN PetscErrorCode DMDASetBoundaryType(DM,DMDABoundaryType,DMDABoundaryType,DMDABoundaryType); 135 PETSC_EXTERN PetscErrorCode DMDASetDof(DM, PetscInt); 136 PETSC_EXTERN PetscErrorCode DMDASetOverlap(DM, 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 /*S 152 DMDALocalInfo - C struct that contains information about a structured grid and a processors logical 153 location in it. 154 155 Level: beginner 156 157 Concepts: distributed array 158 159 Developer note: Then entries in this struct are int instead of PetscInt so that the elements may 160 be extracted in Fortran as if from an integer array 161 162 .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DM, DMDAGetLocalInfo(), DMDAGetInfo() 163 S*/ 164 typedef struct { 165 PetscInt dim,dof,sw; 166 PetscInt mx,my,mz; /* global number of grid points in each direction */ 167 PetscInt xs,ys,zs; /* starting point of this processor, excluding ghosts */ 168 PetscInt xm,ym,zm; /* number of grid points on this processor, excluding ghosts */ 169 PetscInt gxs,gys,gzs; /* starting point of this processor including ghosts */ 170 PetscInt gxm,gym,gzm; /* number of grid points on this processor including ghosts */ 171 DMDABoundaryType bx,by,bz; /* type of ghost nodes at boundary */ 172 DMDAStencilType st; 173 DM da; 174 } DMDALocalInfo; 175 176 /*MC 177 DMDAForEachPointBegin2d - Starts a loop over the local part of a two dimensional DMDA 178 179 Synopsis: 180 void DMDAForEachPointBegin2d(DALocalInfo *info,PetscInt i,PetscInt j); 181 182 Not Collective 183 184 Level: intermediate 185 186 .seealso: DMDAForEachPointEnd2d(), DMDAVecGetArray() 187 M*/ 188 #define DMDAForEachPointBegin2d(info,i,j) {\ 189 PetscInt _xints = info->xs,_xinte = info->xs+info->xm,_yints = info->ys,_yinte = info->ys+info->ym;\ 190 for (j=_yints; j<_yinte; j++) {\ 191 for (i=_xints; i<_xinte; i++) {\ 192 193 /*MC 194 DMDAForEachPointEnd2d - Ends a loop over the local part of a two dimensional DMDA 195 196 Synopsis: 197 void DMDAForEachPointEnd2d; 198 199 Not Collective 200 201 Level: intermediate 202 203 .seealso: DMDAForEachPointBegin2d(), DMDAVecGetArray() 204 M*/ 205 #define DMDAForEachPointEnd2d }}} 206 207 /*MC 208 DMDACoor2d - Structure for holding 2d (x and y) coordinates. 209 210 Level: intermediate 211 212 Sample Usage: 213 DMDACoor2d **coors; 214 Vec vcoors; 215 DM cda; 216 217 DMDAGetCoordinates(da,&vcoors); 218 DMDAGetCoordinateDA(da,&cda); 219 DMDAVecGetArray(cda,vcoors,&coors); 220 DMDAGetCorners(cda,&mstart,&nstart,0,&m,&n,0) 221 for (i=mstart; i<mstart+m; i++) { 222 for (j=nstart; j<nstart+n; j++) { 223 x = coors[j][i].x; 224 y = coors[j][i].y; 225 ...... 226 } 227 } 228 DMDAVecRestoreArray(dac,vcoors,&coors); 229 230 .seealso: DMDACoor3d, DMDAForEachPointBegin(), DMDAGetCoordinateDA(), DMDAGetCoordinates(), DMDAGetGhostCoordinates() 231 M*/ 232 typedef struct {PetscScalar x,y;} DMDACoor2d; 233 234 /*MC 235 DMDACoor3d - Structure for holding 3d (x, y and z) coordinates. 236 237 Level: intermediate 238 239 Sample Usage: 240 DMDACoor3d ***coors; 241 Vec vcoors; 242 DM cda; 243 244 DMDAGetCoordinates(da,&vcoors); 245 DMDAGetCoordinateDA(da,&cda); 246 DMDAVecGetArray(cda,vcoors,&coors); 247 DMDAGetCorners(cda,&mstart,&nstart,&pstart,&m,&n,&p) 248 for (i=mstart; i<mstart+m; i++) { 249 for (j=nstart; j<nstart+n; j++) { 250 for (k=pstart; k<pstart+p; k++) { 251 x = coors[k][j][i].x; 252 y = coors[k][j][i].y; 253 z = coors[k][j][i].z; 254 ...... 255 } 256 } 257 DMDAVecRestoreArray(dac,vcoors,&coors); 258 259 .seealso: DMDACoor2d, DMDAForEachPointBegin(), DMDAGetCoordinateDA(), DMDAGetCoordinates(), DMDAGetGhostCoordinates() 260 M*/ 261 typedef struct {PetscScalar x,y,z;} DMDACoor3d; 262 263 PETSC_EXTERN PetscErrorCode DMDAGetLocalInfo(DM,DMDALocalInfo*); 264 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDALocalFunction1)(DMDALocalInfo*,void*,void*,void*); 265 PETSC_EXTERN PetscErrorCode DMDAComputeFunctionLocal(DM, DMDALocalFunction1, Vec, Vec, void *); 266 PETSC_EXTERN PetscErrorCode DMDAComputeFunctionLocalGhost(DM, DMDALocalFunction1, Vec, Vec, void *); 267 PETSC_EXTERN PetscErrorCode DMDAFormJacobianLocal(DM, DMDALocalFunction1, Vec, Mat, void *); 268 PETSC_EXTERN PetscErrorCode DMDAComputeFunction1(DM,Vec,Vec,void*); 269 PETSC_EXTERN PetscErrorCode DMDAComputeFunction(DM,PetscErrorCode (*)(void),Vec,Vec,void*); 270 PETSC_EXTERN PetscErrorCode DMDAComputeFunctioni1(DM,PetscInt,Vec,PetscScalar*,void*); 271 PETSC_EXTERN PetscErrorCode DMDAComputeFunctionib1(DM,PetscInt,Vec,PetscScalar*,void*); 272 PETSC_EXTERN PetscErrorCode DMDAComputeJacobian1WithAdic(DM,Vec,Mat,void*); 273 PETSC_EXTERN PetscErrorCode DMDAComputeJacobian1WithAdifor(DM,Vec,Mat,void*); 274 PETSC_EXTERN PetscErrorCode DMDAMultiplyByJacobian1WithAdic(DM,Vec,Vec,Vec,void*); 275 PETSC_EXTERN PetscErrorCode DMDAMultiplyByJacobian1WithAdifor(DM,Vec,Vec,Vec,void*); 276 PETSC_EXTERN PetscErrorCode DMDAMultiplyByJacobian1WithAD(DM,Vec,Vec,Vec,void*); 277 PETSC_EXTERN PetscErrorCode DMDAComputeJacobian1(DM,Vec,Mat,void*); 278 PETSC_EXTERN PetscErrorCode DMDAGetLocalFunction(DM,DMDALocalFunction1*); 279 PETSC_EXTERN PetscErrorCode DMDASetLocalFunction(DM,DMDALocalFunction1); 280 PETSC_EXTERN PetscErrorCode DMDASetLocalFunctioni(DM,PetscErrorCode (*)(DMDALocalInfo*,MatStencil*,void*,PetscScalar*,void*)); 281 PETSC_EXTERN PetscErrorCode DMDASetLocalFunctionib(DM,PetscErrorCode (*)(DMDALocalInfo*,MatStencil*,void*,PetscScalar*,void*)); 282 PETSC_EXTERN PetscErrorCode DMDAGetLocalJacobian(DM,DMDALocalFunction1*); 283 PETSC_EXTERN PetscErrorCode DMDASetLocalJacobian(DM,DMDALocalFunction1); 284 PETSC_EXTERN PetscErrorCode DMDASetLocalAdicFunction_Private(DM,DMDALocalFunction1); 285 286 PETSC_EXTERN PetscErrorCode MatSetDM(Mat,DM); 287 PETSC_EXTERN PetscErrorCode MatRegisterDAAD(void); 288 PETSC_EXTERN PetscErrorCode MatCreateDAAD(DM,Mat*); 289 PETSC_EXTERN PetscErrorCode MatCreateSeqUSFFT(Vec,DM,Mat*); 290 291 /*MC 292 DMDASetLocalAdicFunction - Caches in a DM a local function computed by ADIC/ADIFOR 293 294 Synopsis: 295 PetscErrorCode DMDASetLocalAdicFunction(DM da,DMDALocalFunction1 ad_lf) 296 297 Logically Collective on DM 298 299 Input Parameter: 300 + da - initial distributed array 301 - ad_lf - the local function as computed by ADIC/ADIFOR 302 303 Level: intermediate 304 305 .keywords: distributed array, refine 306 307 .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDAGetLocalFunction(), DMDASetLocalFunction(), 308 DMDASetLocalJacobian() 309 M*/ 310 #if defined(PETSC_HAVE_ADIC) 311 # define DMDASetLocalAdicFunction(a,d) DMDASetLocalAdicFunction_Private(a,(DMDALocalFunction1)d) 312 #else 313 # define DMDASetLocalAdicFunction(a,d) DMDASetLocalAdicFunction_Private(a,0) 314 #endif 315 316 PETSC_EXTERN PetscErrorCode DMDASetLocalAdicMFFunction_Private(DM,DMDALocalFunction1); 317 #if defined(PETSC_HAVE_ADIC) 318 # define DMDASetLocalAdicMFFunction(a,d) DMDASetLocalAdicMFFunction_Private(a,(DMDALocalFunction1)d) 319 #else 320 # define DMDASetLocalAdicMFFunction(a,d) DMDASetLocalAdicMFFunction_Private(a,0) 321 #endif 322 PETSC_EXTERN PetscErrorCode DMDASetLocalAdicFunctioni_Private(DM,PetscErrorCode (*)(DMDALocalInfo*,MatStencil*,void*,void*,void*)); 323 #if defined(PETSC_HAVE_ADIC) 324 # define DMDASetLocalAdicFunctioni(a,d) DMDASetLocalAdicFunctioni_Private(a,(PetscErrorCode (*)(DMDALocalInfo*,MatStencil*,void*,void*,void*))d) 325 #else 326 # define DMDASetLocalAdicFunctioni(a,d) DMDASetLocalAdicFunctioni_Private(a,0) 327 #endif 328 PETSC_EXTERN PetscErrorCode DMDASetLocalAdicMFFunctioni_Private(DM,PetscErrorCode (*)(DMDALocalInfo*,MatStencil*,void*,void*,void*)); 329 #if defined(PETSC_HAVE_ADIC) 330 # define DMDASetLocalAdicMFFunctioni(a,d) DMDASetLocalAdicMFFunctioni_Private(a,(PetscErrorCode (*)(DMDALocalInfo*,MatStencil*,void*,void*,void*))d) 331 #else 332 # define DMDASetLocalAdicMFFunctioni(a,d) DMDASetLocalAdicMFFunctioni_Private(a,0) 333 #endif 334 335 PETSC_EXTERN PetscErrorCode DMDASetLocalAdicFunctionib_Private(DM,PetscErrorCode (*)(DMDALocalInfo*,MatStencil*,void*,void*,void*)); 336 #if defined(PETSC_HAVE_ADIC) 337 # define DMDASetLocalAdicFunctionib(a,d) DMDASetLocalAdicFunctionib_Private(a,(PetscErrorCode (*)(DMDALocalInfo*,MatStencil*,void*,void*,void*))d) 338 #else 339 # define DMDASetLocalAdicFunctionib(a,d) DMDASetLocalAdicFunctionib_Private(a,0) 340 #endif 341 PETSC_EXTERN PetscErrorCode DMDASetLocalAdicMFFunctionib_Private(DM,PetscErrorCode (*)(DMDALocalInfo*,MatStencil*,void*,void*,void*)); 342 #if defined(PETSC_HAVE_ADIC) 343 # define DMDASetLocalAdicMFFunctionib(a,d) DMDASetLocalAdicMFFunctionib_Private(a,(PetscErrorCode (*)(DMDALocalInfo*,MatStencil*,void*,void*,void*))d) 344 #else 345 # define DMDASetLocalAdicMFFunctionib(a,d) DMDASetLocalAdicMFFunctionib_Private(a,0) 346 #endif 347 348 PETSC_EXTERN PetscErrorCode DMDAComputeFunctioniTest1(DM,void*); 349 PETSC_EXTERN PetscErrorCode DMDASetGetMatrix(DM,PetscErrorCode (*)(DM, const MatType,Mat *)); 350 PETSC_EXTERN PetscErrorCode DMDASetBlockFills(DM,PetscInt*,PetscInt*); 351 PETSC_EXTERN PetscErrorCode DMDASetRefinementFactor(DM,PetscInt,PetscInt,PetscInt); 352 PETSC_EXTERN PetscErrorCode DMDAGetRefinementFactor(DM,PetscInt*,PetscInt*,PetscInt*); 353 354 PETSC_EXTERN PetscErrorCode DMDAGetAdicArray(DM,PetscBool ,void*,void*,PetscInt*); 355 PETSC_EXTERN PetscErrorCode DMDARestoreAdicArray(DM,PetscBool ,void*,void*,PetscInt*); 356 PETSC_EXTERN PetscErrorCode DMDAGetAdicMFArray(DM,PetscBool ,void*,void*,PetscInt*); 357 PETSC_EXTERN PetscErrorCode DMDAGetAdicMFArray4(DM,PetscBool ,void*,void*,PetscInt*); 358 PETSC_EXTERN PetscErrorCode DMDAGetAdicMFArray9(DM,PetscBool ,void*,void*,PetscInt*); 359 PETSC_EXTERN PetscErrorCode DMDAGetAdicMFArrayb(DM,PetscBool ,void*,void*,PetscInt*); 360 PETSC_EXTERN PetscErrorCode DMDARestoreAdicMFArray(DM,PetscBool ,void*,void*,PetscInt*); 361 PETSC_EXTERN PetscErrorCode DMDAGetArray(DM,PetscBool ,void*); 362 PETSC_EXTERN PetscErrorCode DMDARestoreArray(DM,PetscBool ,void*); 363 PETSC_EXTERN PetscErrorCode ad_DAGetArray(DM,PetscBool ,void*); 364 PETSC_EXTERN PetscErrorCode ad_DARestoreArray(DM,PetscBool ,void*); 365 PETSC_EXTERN PetscErrorCode admf_DAGetArray(DM,PetscBool ,void*); 366 PETSC_EXTERN PetscErrorCode admf_DARestoreArray(DM,PetscBool ,void*); 367 368 PETSC_EXTERN PetscErrorCode DMDACreatePF(DM,PF*); 369 370 PETSC_EXTERN PetscErrorCode DMDAGetNumCells(DM, PetscInt *); 371 PETSC_EXTERN PetscErrorCode DMDAGetNumVertices(DM, PetscInt *, PetscInt *, PetscInt *, PetscInt *); 372 PETSC_EXTERN PetscErrorCode DMDAGetNumFaces(DM, PetscInt *, PetscInt *, PetscInt *, PetscInt *, PetscInt *, PetscInt *); 373 PETSC_EXTERN PetscErrorCode DMDAGetHeightStratum(DM, PetscInt, PetscInt *, PetscInt *); 374 PETSC_EXTERN PetscErrorCode DMDACreateSection(DM, PetscInt[], PetscInt[], PetscInt[], PetscInt[]); 375 PETSC_EXTERN PetscErrorCode DMDAComputeCellGeometry(DM, PetscInt, PetscQuadrature *, PetscReal [], PetscReal [], PetscReal [], PetscReal []); 376 PETSC_EXTERN PetscErrorCode DMDAVecGetClosure(DM, PetscSection, Vec, PetscInt, const PetscScalar **); 377 PETSC_EXTERN PetscErrorCode DMDAVecSetClosure(DM, PetscSection, Vec, PetscInt, const PetscScalar *, InsertMode); 378 PETSC_EXTERN PetscErrorCode DMDAGetClosure(DM,PetscSection,PetscInt,PetscInt*,const PetscInt**); 379 PETSC_EXTERN PetscErrorCode DMDARestoreClosure(DM,PetscSection,PetscInt,PetscInt*,const PetscInt**); 380 PETSC_EXTERN PetscErrorCode DMDAGetClosureScalar(DM,PetscSection,PetscInt,PetscScalar*,const PetscScalar**); 381 PETSC_EXTERN PetscErrorCode DMDARestoreClosureScalar(DM,PetscSection,PetscInt,PetscScalar*,const PetscScalar**); 382 PETSC_EXTERN PetscErrorCode DMDASetClosureScalar(DM,PetscSection,PetscInt,PetscScalar*,const PetscScalar*,InsertMode); 383 PETSC_EXTERN PetscErrorCode DMDAConvertToCell(DM, MatStencil, PetscInt *); 384 385 #define DMDA_FILE_CLASSID 1211220 386 #endif 387