1 #if !defined(__PETSCDMDA_H) 2 #define __PETSCDMDA_H 3 4 #include <petscdm.h> 5 #include <petscdmdatypes.h> 6 #include <petscpf.h> 7 #include <petscao.h> 8 #include <petscfe.h> 9 10 /*MC 11 DMDA_STENCIL_STAR - "Star"-type stencil. In logical grid coordinates, only (i,j,k), (i+s,j,k), (i,j+s,k), 12 (i,j,k+s) are in the stencil NOT, for example, (i+s,j+s,k) 13 14 Level: beginner 15 16 Determines what ghost point values are brought over to each process; in this case the "corner" values are not 17 brought over and hence should not be accessed locally 18 19 .seealso: DMDA_STENCIL_BOX, DMDAStencilType, DMDASetStencilType() 20 M*/ 21 22 /*MC 23 DMDA_STENCIL_BOX - "Box"-type stencil. In logical grid coordinates, any of (i,j,k), (i+s,j+r,k+t) may 24 be in the stencil. 25 26 Level: beginner 27 28 .seealso: DMDA_STENCIL_STAR, DMDAStencilType, DMDASetStencilType() 29 M*/ 30 31 PETSC_EXTERN PetscErrorCode DMDASetInterpolationType(DM,DMDAInterpolationType); 32 PETSC_EXTERN PetscErrorCode DMDAGetInterpolationType(DM,DMDAInterpolationType*); 33 34 PETSC_EXTERN PetscErrorCode DMDASetElementType(DM,DMDAElementType); 35 PETSC_EXTERN PetscErrorCode DMDAGetElementType(DM,DMDAElementType*); 36 PETSC_EXTERN PetscErrorCode DMDAGetElements(DM,PetscInt *,PetscInt *,const PetscInt*[]); 37 PETSC_EXTERN PetscErrorCode DMDARestoreElements(DM,PetscInt *,PetscInt *,const PetscInt*[]); 38 39 typedef enum { DMDA_X,DMDA_Y,DMDA_Z } DMDADirection; 40 41 #define MATSEQUSFFT "sequsfft" 42 43 PETSC_EXTERN PetscErrorCode DMDACreate(MPI_Comm,DM*); 44 PETSC_EXTERN PetscErrorCode DMDASetDim(DM,PetscInt); 45 PETSC_EXTERN PetscErrorCode DMDASetSizes(DM,PetscInt,PetscInt,PetscInt); 46 PETSC_EXTERN PetscErrorCode DMDACreate1d(MPI_Comm,DMBoundaryType,PetscInt,PetscInt,PetscInt,const PetscInt[],DM *); 47 PETSC_EXTERN PetscErrorCode DMDACreate2d(MPI_Comm,DMBoundaryType,DMBoundaryType,DMDAStencilType,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],DM*); 48 PETSC_EXTERN PetscErrorCode DMDACreate3d(MPI_Comm,DMBoundaryType,DMBoundaryType,DMBoundaryType,DMDAStencilType,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscInt[],DM*); 49 50 PETSC_EXTERN PetscErrorCode DMDAGlobalToNaturalBegin(DM,Vec,InsertMode,Vec); 51 PETSC_EXTERN PetscErrorCode DMDAGlobalToNaturalEnd(DM,Vec,InsertMode,Vec); 52 PETSC_EXTERN PetscErrorCode DMDANaturalToGlobalBegin(DM,Vec,InsertMode,Vec); 53 PETSC_EXTERN PetscErrorCode DMDANaturalToGlobalEnd(DM,Vec,InsertMode,Vec); 54 PETSC_DEPRECATED("Use DMLocalToLocalBegin()") PETSC_STATIC_INLINE PetscErrorCode DMDALocalToLocalBegin(DM dm,Vec g,InsertMode mode,Vec l) {return DMLocalToLocalBegin(dm,g,mode,l);} 55 PETSC_DEPRECATED("Use DMLocalToLocalEnd()") PETSC_STATIC_INLINE PetscErrorCode DMDALocalToLocalEnd(DM dm,Vec g,InsertMode mode,Vec l) {return DMLocalToLocalEnd(dm,g,mode,l);} 56 PETSC_EXTERN PetscErrorCode DMDACreateNaturalVector(DM,Vec *); 57 58 PETSC_EXTERN PetscErrorCode DMDAGetCorners(DM,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*); 59 PETSC_EXTERN PetscErrorCode DMDAGetGhostCorners(DM,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*); 60 PETSC_EXTERN PetscErrorCode DMDAGetInfo(DM,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,DMBoundaryType*,DMBoundaryType*,DMBoundaryType*,DMDAStencilType*); 61 PETSC_EXTERN PetscErrorCode DMDAGetProcessorSubset(DM,DMDADirection,PetscInt,MPI_Comm*); 62 PETSC_EXTERN PetscErrorCode DMDAGetProcessorSubsets(DM,DMDADirection,MPI_Comm*); 63 PETSC_EXTERN PetscErrorCode DMDAGetRay(DM,DMDADirection,PetscInt,Vec*,VecScatter*); 64 65 PETSC_EXTERN PetscErrorCode DMDAGlobalToNaturalAllCreate(DM,VecScatter*); 66 PETSC_EXTERN PetscErrorCode DMDANaturalAllToGlobalCreate(DM,VecScatter*); 67 68 PETSC_EXTERN PetscErrorCode DMDAGetGlobalIndices(DM,PetscInt*,const PetscInt*[]); 69 PETSC_EXTERN PetscErrorCode DMDARestoreGlobalIndices(DM,PetscInt*,const PetscInt*[]); 70 71 PETSC_EXTERN PetscErrorCode DMDAGetScatter(DM,VecScatter*,VecScatter*,VecScatter*); 72 PETSC_EXTERN PetscErrorCode DMDAGetNeighbors(DM,const PetscMPIInt**); 73 74 PETSC_EXTERN PetscErrorCode DMDAGetAO(DM,AO*); 75 PETSC_EXTERN PetscErrorCode DMDASetUniformCoordinates(DM,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal); 76 PETSC_EXTERN PetscErrorCode DMDAGetBoundingBox(DM,PetscReal[],PetscReal[]); 77 PETSC_EXTERN PetscErrorCode DMDAGetLocalBoundingBox(DM,PetscReal[],PetscReal[]); 78 PETSC_EXTERN PetscErrorCode DMDAGetLogicalCoordinate(DM,PetscScalar,PetscScalar,PetscScalar,PetscInt*,PetscInt*,PetscInt*,PetscScalar*,PetscScalar*,PetscScalar*); 79 /* function to wrap coordinates around boundary */ 80 PETSC_EXTERN PetscErrorCode DMDAMapCoordsToPeriodicDomain(DM,PetscScalar*,PetscScalar*); 81 82 PETSC_EXTERN PetscErrorCode DMDAGetReducedDMDA(DM,PetscInt,DM*); 83 84 PETSC_EXTERN PetscErrorCode DMDASetFieldName(DM,PetscInt,const char[]); 85 PETSC_EXTERN PetscErrorCode DMDAGetFieldName(DM,PetscInt,const char**); 86 PETSC_EXTERN PetscErrorCode DMDASetCoordinateName(DM,PetscInt,const char[]); 87 PETSC_EXTERN PetscErrorCode DMDAGetCoordinateName(DM,PetscInt,const char**); 88 89 PETSC_EXTERN PetscErrorCode DMDASetBoundaryType(DM,DMBoundaryType,DMBoundaryType,DMBoundaryType); 90 PETSC_EXTERN PetscErrorCode DMDASetDof(DM, PetscInt); 91 PETSC_EXTERN PetscErrorCode DMDASetOverlap(DM,PetscInt,PetscInt,PetscInt); 92 PETSC_EXTERN PetscErrorCode DMDAGetOverlap(DM,PetscInt*,PetscInt*,PetscInt*); 93 PETSC_EXTERN PetscErrorCode DMDASetNumLocalSubDomains(DM,PetscInt); 94 PETSC_EXTERN PetscErrorCode DMDAGetNumLocalSubDomains(DM,PetscInt*); 95 PETSC_EXTERN PetscErrorCode DMDAGetOffset(DM,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*); 96 PETSC_EXTERN PetscErrorCode DMDASetOffset(DM,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt); 97 PETSC_EXTERN PetscErrorCode DMDAGetNonOverlappingRegion(DM,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*); 98 PETSC_EXTERN PetscErrorCode DMDASetNonOverlappingRegion(DM,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt); 99 PETSC_EXTERN PetscErrorCode DMDASetStencilWidth(DM, PetscInt); 100 PETSC_EXTERN PetscErrorCode DMDASetOwnershipRanges(DM,const PetscInt[],const PetscInt[],const PetscInt[]); 101 PETSC_EXTERN PetscErrorCode DMDAGetOwnershipRanges(DM,const PetscInt**,const PetscInt**,const PetscInt**); 102 PETSC_EXTERN PetscErrorCode DMDASetNumProcs(DM, PetscInt, PetscInt, PetscInt); 103 PETSC_EXTERN PetscErrorCode DMDASetStencilType(DM, DMDAStencilType); 104 105 PETSC_EXTERN PetscErrorCode DMDAVecGetArray(DM,Vec,void *); 106 PETSC_EXTERN PetscErrorCode DMDAVecRestoreArray(DM,Vec,void *); 107 108 PETSC_EXTERN PetscErrorCode DMDAVecGetArrayDOF(DM,Vec,void *); 109 PETSC_EXTERN PetscErrorCode DMDAVecRestoreArrayDOF(DM,Vec,void *); 110 111 PETSC_EXTERN PetscErrorCode DMDASplitComm2d(MPI_Comm,PetscInt,PetscInt,PetscInt,MPI_Comm*); 112 113 PETSC_EXTERN PetscErrorCode DMDACreatePatchIS(DM,MatStencil*,MatStencil*,IS*); 114 115 116 /*MC 117 DMDACoor2d - Structure for holding 2d (x and y) coordinates. 118 119 Level: intermediate 120 121 Sample Usage: 122 DMDACoor2d **coors; 123 Vec vcoors; 124 DM cda; 125 126 DMGetCoordinates(da,&vcoors); 127 DMGetCoordinateDM(da,&cda); 128 DMDAVecGetArray(cda,vcoors,&coors); 129 DMDAGetCorners(cda,&mstart,&nstart,0,&m,&n,0) 130 for (i=mstart; i<mstart+m; i++) { 131 for (j=nstart; j<nstart+n; j++) { 132 x = coors[j][i].x; 133 y = coors[j][i].y; 134 ...... 135 } 136 } 137 DMDAVecRestoreArray(dac,vcoors,&coors); 138 139 .seealso: DMDACoor3d, DMGetCoordinateDM(), DMGetCoordinates(), DMDAGetGhostCoordinates() 140 M*/ 141 typedef struct {PetscScalar x,y;} DMDACoor2d; 142 143 /*MC 144 DMDACoor3d - Structure for holding 3d (x, y and z) coordinates. 145 146 Level: intermediate 147 148 Sample Usage: 149 DMDACoor3d ***coors; 150 Vec vcoors; 151 DM cda; 152 153 DMGetCoordinates(da,&vcoors); 154 DMGetCoordinateDM(da,&cda); 155 DMDAVecGetArray(cda,vcoors,&coors); 156 DMDAGetCorners(cda,&mstart,&nstart,&pstart,&m,&n,&p) 157 for (i=mstart; i<mstart+m; i++) { 158 for (j=nstart; j<nstart+n; j++) { 159 for (k=pstart; k<pstart+p; k++) { 160 x = coors[k][j][i].x; 161 y = coors[k][j][i].y; 162 z = coors[k][j][i].z; 163 ...... 164 } 165 } 166 DMDAVecRestoreArray(dac,vcoors,&coors); 167 168 .seealso: DMDACoor2d, DMGetCoordinateDM(), DMGetCoordinates(), DMDAGetGhostCoordinates() 169 M*/ 170 typedef struct {PetscScalar x,y,z;} DMDACoor3d; 171 172 PETSC_EXTERN PetscErrorCode DMDAGetLocalInfo(DM,DMDALocalInfo*); 173 174 PETSC_EXTERN PetscErrorCode MatRegisterDAAD(void); 175 PETSC_EXTERN PetscErrorCode MatCreateDAAD(DM,Mat*); 176 PETSC_EXTERN PetscErrorCode MatCreateSeqUSFFT(Vec,DM,Mat*); 177 178 PETSC_EXTERN PetscErrorCode DMDASetGetMatrix(DM,PetscErrorCode (*)(DM, Mat *)); 179 PETSC_EXTERN PetscErrorCode DMDASetBlockFills(DM,const PetscInt*,const PetscInt*); 180 PETSC_EXTERN PetscErrorCode DMDASetRefinementFactor(DM,PetscInt,PetscInt,PetscInt); 181 PETSC_EXTERN PetscErrorCode DMDAGetRefinementFactor(DM,PetscInt*,PetscInt*,PetscInt*); 182 183 PETSC_EXTERN PetscErrorCode DMDAGetArray(DM,PetscBool ,void*); 184 PETSC_EXTERN PetscErrorCode DMDARestoreArray(DM,PetscBool ,void*); 185 186 PETSC_EXTERN PetscErrorCode DMDACreatePF(DM,PF*); 187 188 /* Emulation of DMPlex */ 189 PETSC_EXTERN PetscErrorCode DMDAGetNumCells(DM, PetscInt *, PetscInt *, PetscInt *, PetscInt *); 190 PETSC_EXTERN PetscErrorCode DMDAGetCellPoint(DM, PetscInt, PetscInt, PetscInt, PetscInt *); 191 PETSC_EXTERN PetscErrorCode DMDAGetNumVertices(DM, PetscInt *, PetscInt *, PetscInt *, PetscInt *); 192 PETSC_EXTERN PetscErrorCode DMDAGetNumFaces(DM, PetscInt *, PetscInt *, PetscInt *, PetscInt *, PetscInt *, PetscInt *); 193 PETSC_EXTERN PetscErrorCode DMDAGetHeightStratum(DM, PetscInt, PetscInt *, PetscInt *); 194 PETSC_EXTERN PetscErrorCode DMDACreateSection(DM, const PetscInt[], const PetscInt[], const PetscInt[], PetscSection *); 195 PETSC_EXTERN PetscErrorCode DMDAComputeCellGeometry(DM, PetscInt, PetscQuadrature, PetscReal [], PetscReal [], PetscReal [], PetscReal []); 196 PETSC_EXTERN PetscErrorCode DMDAGetTransitiveClosure(DM, PetscInt, PetscBool, PetscInt *, PetscInt **); 197 PETSC_EXTERN PetscErrorCode DMDARestoreTransitiveClosure(DM, PetscInt, PetscBool, PetscInt *, PetscInt **); 198 PETSC_EXTERN PetscErrorCode DMDAVecGetClosure(DM, PetscSection, Vec, PetscInt, PetscInt *, PetscScalar **); 199 PETSC_EXTERN PetscErrorCode DMDAVecRestoreClosure(DM, PetscSection, Vec, PetscInt, PetscInt *, PetscScalar **); 200 PETSC_EXTERN PetscErrorCode DMDAVecSetClosure(DM, PetscSection, Vec, PetscInt, const PetscScalar *, InsertMode); 201 PETSC_EXTERN PetscErrorCode DMDAGetClosure(DM, PetscSection, PetscInt, PetscInt*, const PetscInt**); 202 PETSC_EXTERN PetscErrorCode DMDARestoreClosure(DM, PetscSection, PetscInt, PetscInt*, const PetscInt**); 203 PETSC_EXTERN PetscErrorCode DMDAGetClosureScalar(DM, PetscSection, PetscInt, PetscScalar*, PetscInt*, PetscScalar**); 204 PETSC_EXTERN PetscErrorCode DMDARestoreClosureScalar(DM, PetscSection, PetscInt, PetscScalar*, PetscInt*, PetscScalar**); 205 PETSC_EXTERN PetscErrorCode DMDASetClosureScalar(DM,PetscSection,PetscInt,PetscScalar*,const PetscScalar*,InsertMode); 206 PETSC_EXTERN PetscErrorCode DMDAConvertToCell(DM, MatStencil, PetscInt *); 207 PETSC_EXTERN PetscErrorCode DMDASetVertexCoordinates(DM,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal); 208 PETSC_EXTERN PetscErrorCode DMDASetPreallocationCenterDimension(DM, PetscInt); 209 PETSC_EXTERN PetscErrorCode DMDAGetPreallocationCenterDimension(DM, PetscInt*); 210 211 PETSC_EXTERN PetscErrorCode DMDAProjectFunction(DM, PetscFE[], void (**)(const PetscReal [], PetscScalar *, void *), void **, InsertMode, Vec); 212 PETSC_EXTERN PetscErrorCode DMDAProjectFunctionLocal(DM, PetscFE[], void (**)(const PetscReal [], PetscScalar *, void *), void **, InsertMode, Vec); 213 PETSC_EXTERN PetscErrorCode DMDAComputeL2Diff(DM, PetscFE[], void (**)(const PetscReal [], PetscScalar *, void *), void **, Vec, PetscReal *); 214 PETSC_EXTERN PetscErrorCode DMDAComputeL2GradientDiff(DM, PetscFE[], void (**)(const PetscReal [], const PetscReal [], PetscScalar *, void *), void **, Vec, const PetscReal[], PetscReal *); 215 216 #endif 217