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