1 #ifndef 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 /* SUBMANSEC = DMDA */ 11 12 /*MC 13 DMDA_STENCIL_STAR - "Star"-type stencil. In logical grid coordinates, only (i,j,k), (i+s,j,k), (i,j+s,k), 14 (i,j,k+s) are in the stencil NOT, for example, (i+s,j+s,k) 15 16 Level: beginner 17 18 Note: 19 Determines what ghost point values are brought over to each process; in this case the "corner" values are not 20 brought over and hence should not be accessed locally 21 22 .seealso: `DMDA_STENCIL_BOX`, `DMDAStencilType`, `DMDASetStencilType()` 23 M*/ 24 25 /*MC 26 DMDA_STENCIL_BOX - "Box"-type stencil. In logical grid coordinates, any of (i,j,k), (i+s,j+r,k+t) may 27 be in the stencil. 28 29 Level: beginner 30 31 .seealso: `DMDA_STENCIL_STAR`, `DMDAStencilType`, `DMDASetStencilType()` 32 M*/ 33 34 PETSC_EXTERN PetscErrorCode DMDASetInterpolationType(DM, DMDAInterpolationType); 35 PETSC_EXTERN PetscErrorCode DMDAGetInterpolationType(DM, DMDAInterpolationType *); 36 PETSC_EXTERN PetscErrorCode DMDACreateAggregates(DM, DM, Mat *); 37 38 /* FEM */ 39 PETSC_EXTERN PetscErrorCode DMDASetElementType(DM, DMDAElementType); 40 PETSC_EXTERN PetscErrorCode DMDAGetElementType(DM, DMDAElementType *); 41 PETSC_EXTERN PetscErrorCode DMDAGetElements(DM, PetscInt *, PetscInt *, const PetscInt *[]); 42 PETSC_EXTERN PetscErrorCode DMDARestoreElements(DM, PetscInt *, PetscInt *, const PetscInt *[]); 43 PETSC_EXTERN PetscErrorCode DMDAGetElementsSizes(DM, PetscInt *, PetscInt *, PetscInt *); 44 PETSC_EXTERN PetscErrorCode DMDAGetElementsCorners(DM, PetscInt *, PetscInt *, PetscInt *); 45 PETSC_EXTERN PetscErrorCode DMDAGetSubdomainCornersIS(DM, IS *); 46 PETSC_EXTERN PetscErrorCode DMDARestoreSubdomainCornersIS(DM, IS *); 47 48 #define MATSEQUSFFT "sequsfft" 49 50 PETSC_EXTERN PetscErrorCode DMDACreate(MPI_Comm, DM *); 51 PETSC_EXTERN PetscErrorCode DMDASetSizes(DM, PetscInt, PetscInt, PetscInt); 52 PETSC_EXTERN PetscErrorCode DMDACreate1d(MPI_Comm, DMBoundaryType, PetscInt, PetscInt, PetscInt, const PetscInt[], DM *); 53 PETSC_EXTERN PetscErrorCode DMDACreate2d(MPI_Comm, DMBoundaryType, DMBoundaryType, DMDAStencilType, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], DM *); 54 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 *); 55 56 PETSC_EXTERN PetscErrorCode DMDAGlobalToNaturalBegin(DM, Vec, InsertMode, Vec); 57 PETSC_EXTERN PetscErrorCode DMDAGlobalToNaturalEnd(DM, Vec, InsertMode, Vec); 58 PETSC_EXTERN PetscErrorCode DMDANaturalToGlobalBegin(DM, Vec, InsertMode, Vec); 59 PETSC_EXTERN PetscErrorCode DMDANaturalToGlobalEnd(DM, Vec, InsertMode, Vec); 60 PETSC_DEPRECATED_FUNCTION(3, 5, 0, "DMLocalToLocalBegin()", ) static inline PetscErrorCode DMDALocalToLocalBegin(DM dm, Vec g, InsertMode mode, Vec l) 61 { 62 return DMLocalToLocalBegin(dm, g, mode, l); 63 } 64 PETSC_DEPRECATED_FUNCTION(3, 5, 0, "DMLocalToLocalEnd()", ) static inline PetscErrorCode DMDALocalToLocalEnd(DM dm, Vec g, InsertMode mode, Vec l) 65 { 66 return DMLocalToLocalEnd(dm, g, mode, l); 67 } 68 PETSC_EXTERN PetscErrorCode DMDACreateNaturalVector(DM, Vec *); 69 70 PETSC_EXTERN PetscErrorCode DMDAGetCorners(DM, PetscInt *, PetscInt *, PetscInt *, PetscInt *, PetscInt *, PetscInt *); 71 PETSC_EXTERN PetscErrorCode DMDAGetGhostCorners(DM, PetscInt *, PetscInt *, PetscInt *, PetscInt *, PetscInt *, PetscInt *); 72 PETSC_EXTERN PetscErrorCode DMDAGetInfo(DM, PetscInt *, PetscInt *, PetscInt *, PetscInt *, PetscInt *, PetscInt *, PetscInt *, PetscInt *, PetscInt *, DMBoundaryType *, DMBoundaryType *, DMBoundaryType *, DMDAStencilType *); 73 PETSC_EXTERN PetscErrorCode DMDAGetProcessorSubset(DM, DMDirection, PetscInt, MPI_Comm *); 74 PETSC_EXTERN PetscErrorCode DMDAGetProcessorSubsets(DM, DMDirection, MPI_Comm *); 75 PETSC_EXTERN PetscErrorCode DMDAGetRay(DM, DMDirection, PetscInt, Vec *, VecScatter *); 76 77 PETSC_EXTERN PetscErrorCode DMDAGlobalToNaturalAllCreate(DM, VecScatter *); 78 PETSC_EXTERN PetscErrorCode DMDANaturalAllToGlobalCreate(DM, VecScatter *); 79 80 PETSC_EXTERN PetscErrorCode DMDAGetScatter(DM, VecScatter *, VecScatter *); 81 PETSC_EXTERN PetscErrorCode DMDAGetNeighbors(DM, const PetscMPIInt **); 82 83 PETSC_EXTERN PetscErrorCode DMDASetAOType(DM, AOType); 84 PETSC_EXTERN PetscErrorCode DMDAGetAO(DM, AO *); 85 PETSC_EXTERN PetscErrorCode DMDASetUniformCoordinates(DM, PetscReal, PetscReal, PetscReal, PetscReal, PetscReal, PetscReal); 86 PETSC_EXTERN PetscErrorCode DMDASetGLLCoordinates(DM, PetscInt, PetscReal *); 87 PETSC_EXTERN PetscErrorCode DMDAGetCoordinateArray(DM, void *); 88 PETSC_EXTERN PetscErrorCode DMDARestoreCoordinateArray(DM, void *); 89 PETSC_EXTERN PetscErrorCode DMDAGetLogicalCoordinate(DM, PetscScalar, PetscScalar, PetscScalar, PetscInt *, PetscInt *, PetscInt *, PetscScalar *, PetscScalar *, PetscScalar *); 90 /* function to wrap coordinates around boundary */ 91 PETSC_EXTERN PetscErrorCode DMDAMapCoordsToPeriodicDomain(DM, PetscScalar *, PetscScalar *); 92 93 PETSC_EXTERN PetscErrorCode DMDACreateCompatibleDMDA(DM, PetscInt, DM *); 94 PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 10, 0, "DMDACreateCompatibleDMDA()", ) PetscErrorCode DMDAGetReducedDMDA(DM, PetscInt, DM *); 95 96 PETSC_EXTERN PetscErrorCode DMDASetFieldName(DM, PetscInt, const char[]); 97 PETSC_EXTERN PetscErrorCode DMDAGetFieldName(DM, PetscInt, const char **); 98 PETSC_EXTERN PetscErrorCode DMDASetFieldNames(DM, const char *const *); 99 PETSC_EXTERN PetscErrorCode DMDAGetFieldNames(DM, const char *const **); 100 PETSC_EXTERN PetscErrorCode DMDASetCoordinateName(DM, PetscInt, const char[]); 101 PETSC_EXTERN PetscErrorCode DMDAGetCoordinateName(DM, PetscInt, const char **); 102 103 PETSC_EXTERN PetscErrorCode DMDASetBoundaryType(DM, DMBoundaryType, DMBoundaryType, DMBoundaryType); 104 PETSC_EXTERN PetscErrorCode DMDASetDof(DM, PetscInt); 105 PETSC_EXTERN PetscErrorCode DMDAGetDof(DM, PetscInt *); 106 PETSC_EXTERN PetscErrorCode DMDASetOverlap(DM, PetscInt, PetscInt, PetscInt); 107 PETSC_EXTERN PetscErrorCode DMDAGetOverlap(DM, PetscInt *, PetscInt *, PetscInt *); 108 PETSC_EXTERN PetscErrorCode DMDASetNumLocalSubDomains(DM, PetscInt); 109 PETSC_EXTERN PetscErrorCode DMDAGetNumLocalSubDomains(DM, PetscInt *); 110 PETSC_EXTERN PetscErrorCode DMDAGetOffset(DM, PetscInt *, PetscInt *, PetscInt *, PetscInt *, PetscInt *, PetscInt *); 111 PETSC_EXTERN PetscErrorCode DMDASetOffset(DM, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt); 112 PETSC_EXTERN PetscErrorCode DMDAGetNonOverlappingRegion(DM, PetscInt *, PetscInt *, PetscInt *, PetscInt *, PetscInt *, PetscInt *); 113 PETSC_EXTERN PetscErrorCode DMDASetNonOverlappingRegion(DM, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt); 114 PETSC_EXTERN PetscErrorCode DMDASetStencilWidth(DM, PetscInt); 115 PETSC_EXTERN PetscErrorCode DMDAGetStencilWidth(DM, PetscInt *); 116 PETSC_EXTERN PetscErrorCode DMDAMapMatStencilToGlobal(DM, PetscInt, const MatStencil[], PetscInt[]); 117 PETSC_EXTERN PetscErrorCode DMDASetOwnershipRanges(DM, const PetscInt[], const PetscInt[], const PetscInt[]); 118 PETSC_EXTERN PetscErrorCode DMDAGetOwnershipRanges(DM, const PetscInt **, const PetscInt **, const PetscInt **); 119 PETSC_EXTERN PetscErrorCode DMDASetNumProcs(DM, PetscInt, PetscInt, PetscInt); 120 PETSC_EXTERN PetscErrorCode DMDASetStencilType(DM, DMDAStencilType); 121 PETSC_EXTERN PetscErrorCode DMDAGetStencilType(DM, DMDAStencilType *); 122 123 PETSC_EXTERN PetscErrorCode DMDAVecGetArray(DM, Vec, void *); 124 PETSC_EXTERN PetscErrorCode DMDAVecRestoreArray(DM, Vec, void *); 125 PETSC_EXTERN PetscErrorCode DMDAVecGetArrayWrite(DM, Vec, void *); 126 PETSC_EXTERN PetscErrorCode DMDAVecRestoreArrayWrite(DM, Vec, void *); 127 128 PETSC_EXTERN PetscErrorCode DMDAVecGetArrayDOF(DM, Vec, void *); 129 PETSC_EXTERN PetscErrorCode DMDAVecRestoreArrayDOF(DM, Vec, void *); 130 131 PETSC_EXTERN PetscErrorCode DMDAVecGetArrayRead(DM, Vec, void *); 132 PETSC_EXTERN PetscErrorCode DMDAVecRestoreArrayRead(DM, Vec, void *); 133 134 PETSC_EXTERN PetscErrorCode DMDAVecGetArrayDOFRead(DM, Vec, void *); 135 PETSC_EXTERN PetscErrorCode DMDAVecRestoreArrayDOFRead(DM, Vec, void *); 136 137 PETSC_EXTERN PetscErrorCode DMDAVecGetArrayDOFWrite(DM, Vec, void *); 138 PETSC_EXTERN PetscErrorCode DMDAVecRestoreArrayDOFWrite(DM, Vec, void *); 139 140 PETSC_EXTERN PetscErrorCode DMDACreatePatchIS(DM, MatStencil *, MatStencil *, IS *, PetscBool); 141 142 /*MC 143 DMDACoor2d - Structure for holding 2d (x and y) coordinates when working with `DMDA` 144 145 Synopsis: 146 .vb 147 DMDACoor2d **coors; 148 Vec vcoors; 149 DM cda; 150 DMGetCoordinates(da,&vcoors); 151 DMGetCoordinateDM(da,&cda); 152 DMDAVecGetArray(cda,vcoors,&coors); 153 DMDAGetCorners(cda,&mstart,&nstart,0,&m,&n,0) 154 for (i=mstart; i<mstart+m; i++) { 155 for (j=nstart; j<nstart+n; j++) { 156 x = coors[j][i].x; 157 y = coors[j][i].y; 158 ...... 159 } 160 } 161 DMDAVecRestoreArray(dac,vcoors,&coors); 162 .ve 163 164 Level: intermediate 165 166 .seealso: `DMDACoor3d`, `DMDAVecRestoreArray()`, `DMDAVecGetArray()`, `DMGetCoordinateDM()`, `DMGetCoordinates()` 167 M*/ 168 typedef struct { 169 PetscScalar x, y; 170 } DMDACoor2d; 171 172 /*MC 173 DMDACoor3d - Structure for holding 3d (x, y and z) coordinates coordinates when working with `DMDA` 174 175 Synopsis: 176 .vb 177 DMDACoor3d ***coors; 178 Vec vcoors; 179 DM cda; 180 DMGetCoordinates(da,&vcoors); 181 DMGetCoordinateDM(da,&cda); 182 DMDAVecGetArray(cda,vcoors,&coors); 183 DMDAGetCorners(cda,&mstart,&nstart,&pstart,&m,&n,&p) 184 for (i=mstart; i<mstart+m; i++) { 185 for (j=nstart; j<nstart+n; j++) { 186 for (k=pstart; k<pstart+p; k++) { 187 x = coors[k][j][i].x; 188 y = coors[k][j][i].y; 189 z = coors[k][j][i].z; 190 ...... 191 } 192 } 193 DMDAVecRestoreArray(dac,vcoors,&coors); 194 .ve 195 196 Level: intermediate 197 198 .seealso: `DMDACoor2d`, `DMDAVecRestoreArray()`, `DMDAVecGetArray()`, `DMGetCoordinateDM()`, `DMGetCoordinates()` 199 M*/ 200 typedef struct { 201 PetscScalar x, y, z; 202 } DMDACoor3d; 203 204 PETSC_EXTERN PetscErrorCode DMDAGetLocalInfo(DM, DMDALocalInfo *); 205 206 PETSC_EXTERN PetscErrorCode MatRegisterDAAD(void); 207 PETSC_EXTERN PetscErrorCode MatCreateDAAD(DM, Mat *); 208 PETSC_EXTERN PetscErrorCode MatCreateSeqUSFFT(Vec, DM, Mat *); 209 210 PETSC_EXTERN PetscErrorCode DMDASetGetMatrix(DM, PetscErrorCode (*)(DM, Mat *)); 211 PETSC_EXTERN PetscErrorCode DMDASetBlockFills(DM, const PetscInt *, const PetscInt *); 212 PETSC_EXTERN PetscErrorCode DMDASetBlockFillsSparse(DM, const PetscInt *, const PetscInt *); 213 PETSC_EXTERN PetscErrorCode DMDASetRefinementFactor(DM, PetscInt, PetscInt, PetscInt); 214 PETSC_EXTERN PetscErrorCode DMDAGetRefinementFactor(DM, PetscInt *, PetscInt *, PetscInt *); 215 216 PETSC_EXTERN PetscErrorCode DMDAGetArray(DM, PetscBool, void *); 217 PETSC_EXTERN PetscErrorCode DMDARestoreArray(DM, PetscBool, void *); 218 219 PETSC_EXTERN PetscErrorCode DMDACreatePF(DM, PF *); 220 221 /* Emulation of DMPlex */ 222 PETSC_EXTERN PetscErrorCode DMDAGetNumCells(DM, PetscInt *, PetscInt *, PetscInt *, PetscInt *); 223 PETSC_EXTERN PetscErrorCode DMDAGetCellPoint(DM, PetscInt, PetscInt, PetscInt, PetscInt *); 224 PETSC_EXTERN PetscErrorCode DMDAGetNumVertices(DM, PetscInt *, PetscInt *, PetscInt *, PetscInt *); 225 PETSC_EXTERN PetscErrorCode DMDAGetNumFaces(DM, PetscInt *, PetscInt *, PetscInt *, PetscInt *, PetscInt *, PetscInt *); 226 PETSC_EXTERN PetscErrorCode DMDAGetHeightStratum(DM, PetscInt, PetscInt *, PetscInt *); 227 PETSC_EXTERN PetscErrorCode DMDAGetDepthStratum(DM, PetscInt, PetscInt *, PetscInt *); 228 PETSC_EXTERN PetscErrorCode DMDAComputeCellGeometryFEM(DM, PetscInt, PetscQuadrature, PetscReal[], PetscReal[], PetscReal[], PetscReal[]); 229 PETSC_EXTERN PetscErrorCode DMDAGetTransitiveClosure(DM, PetscInt, PetscBool, PetscInt *, PetscInt **); 230 PETSC_EXTERN PetscErrorCode DMDARestoreTransitiveClosure(DM, PetscInt, PetscBool, PetscInt *, PetscInt **); 231 PETSC_EXTERN PetscErrorCode DMDAConvertToCell(DM, MatStencil, PetscInt *); 232 PETSC_EXTERN PetscErrorCode DMDASetVertexCoordinates(DM, PetscReal, PetscReal, PetscReal, PetscReal, PetscReal, PetscReal); 233 PETSC_EXTERN PetscErrorCode DMDASetPreallocationCenterDimension(DM, PetscInt); 234 PETSC_EXTERN PetscErrorCode DMDAGetPreallocationCenterDimension(DM, PetscInt *); 235 236 #endif 237