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