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