1 #define PETSCDM_DLL 2 #include <petsc-private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 3 #include <petscdmda.h> 4 #include <petscsf.h> 5 6 #undef __FUNCT__ 7 #define __FUNCT__ "DMPlexCreateDoublet" 8 /*@ 9 DMPlexCreateDoublet - Creates a mesh of two cells of the specified type, optionally with later refinement. 10 11 Collective on MPI_Comm 12 13 Input Parameters: 14 + comm - The communicator for the DM object 15 . dim - The spatial dimension 16 . simplex - Flag for simplicial cells, otherwise they are tensor product cells 17 . interpolate - Flag to create intermediate mesh pieces (edges, faces) 18 . refinementUniform - Flag for uniform parallel refinement 19 - refinementLimit - A nonzero number indicates the largest admissible volume for a refined cell 20 21 Output Parameter: 22 . dm - The DM object 23 24 Level: beginner 25 26 .keywords: DM, create 27 .seealso: DMSetType(), DMCreate() 28 @*/ 29 PetscErrorCode DMPlexCreateDoublet(MPI_Comm comm, PetscInt dim, PetscBool simplex, PetscBool interpolate, PetscBool refinementUniform, PetscReal refinementLimit, DM *newdm) 30 { 31 DM dm; 32 PetscInt p; 33 PetscMPIInt rank; 34 PetscErrorCode ierr; 35 36 PetscFunctionBegin; 37 ierr = DMCreate(comm, &dm);CHKERRQ(ierr); 38 ierr = DMSetType(dm, DMPLEX);CHKERRQ(ierr); 39 ierr = DMPlexSetDimension(dm, dim);CHKERRQ(ierr); 40 ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 41 switch (dim) { 42 case 2: 43 if (simplex) {ierr = PetscObjectSetName((PetscObject) dm, "triangular");CHKERRQ(ierr);} 44 else {ierr = PetscObjectSetName((PetscObject) dm, "quadrilateral");CHKERRQ(ierr);} 45 break; 46 case 3: 47 if (simplex) {ierr = PetscObjectSetName((PetscObject) dm, "tetrahedral");CHKERRQ(ierr);} 48 else {ierr = PetscObjectSetName((PetscObject) dm, "hexahedral");CHKERRQ(ierr);} 49 break; 50 default: 51 SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Cannot make meshes for dimension %d", dim); 52 } 53 if (rank) { 54 PetscInt numPoints[2] = {0, 0}; 55 ierr = DMPlexCreateFromDAG(dm, 1, numPoints, NULL, NULL, NULL, NULL);CHKERRQ(ierr); 56 } else { 57 switch (dim) { 58 case 2: 59 if (simplex) { 60 PetscInt numPoints[2] = {4, 2}; 61 PetscInt coneSize[6] = {3, 3, 0, 0, 0, 0}; 62 PetscInt cones[6] = {2, 3, 4, 5, 4, 3}; 63 PetscInt coneOrientations[6] = {0, 0, 0, 0, 0, 0}; 64 PetscScalar vertexCoords[8] = {-0.5, 0.5, 0.0, 0.0, 0.0, 1.0, 0.5, 0.5}; 65 PetscInt markerPoints[8] = {2, 1, 3, 1, 4, 1, 5, 1}; 66 67 ierr = DMPlexCreateFromDAG(dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr); 68 for (p = 0; p < 4; ++p) {ierr = DMPlexSetLabelValue(dm, "marker", markerPoints[p*2], markerPoints[p*2+1]);CHKERRQ(ierr);} 69 } else { 70 PetscInt numPoints[2] = {6, 2}; 71 PetscInt coneSize[8] = {4, 4, 0, 0, 0, 0, 0, 0}; 72 PetscInt cones[8] = {2, 3, 4, 5, 3, 6, 7, 4}; 73 PetscInt coneOrientations[8] = {0, 0, 0, 0, 0, 0, 0, 0}; 74 PetscScalar vertexCoords[12] = {-1.0, -0.5, 0.0, -0.5, 0.0, 0.5, -1.0, 0.5, 1.0, -0.5, 1.0, 0.5}; 75 76 ierr = DMPlexCreateFromDAG(dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr); 77 } 78 break; 79 case 3: 80 if (simplex) { 81 PetscInt numPoints[2] = {5, 2}; 82 PetscInt coneSize[7] = {4, 4, 0, 0, 0, 0, 0}; 83 PetscInt cones[8] = {4, 3, 5, 2, 5, 3, 4, 6}; 84 PetscInt coneOrientations[8] = {0, 0, 0, 0, 0, 0, 0, 0}; 85 PetscScalar vertexCoords[15] = {-1.0, 0.0, 0.0, 0.0, -1.0, 0.0, 0.0, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0, 0.0, 0.0}; 86 PetscInt markerPoints[10] = {2, 1, 3, 1, 4, 1, 5, 1, 6, 1}; 87 88 ierr = DMPlexCreateFromDAG(dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr); 89 for (p = 0; p < 5; ++p) {ierr = DMPlexSetLabelValue(dm, "marker", markerPoints[p*2], markerPoints[p*2+1]);CHKERRQ(ierr);} 90 } else { 91 PetscInt numPoints[2] = {12, 2}; 92 PetscInt coneSize[14] = {8, 8, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0}; 93 PetscInt cones[16] = {2, 3, 4, 5, 6, 7, 8, 9, 5, 4, 10, 11, 7, 12, 13, 8}; 94 PetscInt coneOrientations[16] = {0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0}; 95 PetscScalar vertexCoords[36] = {-1.0, -0.5, -0.5, -1.0, 0.5, -0.5, 0.0, 0.5, -0.5, 0.0, -0.5, -0.5, 96 -1.0, -0.5, 0.5, 0.0, -0.5, 0.5, 0.0, 0.5, 0.5, -1.0, 0.5, 0.5, 97 1.0, 0.5, -0.5, 1.0, -0.5, -0.5, 1.0, -0.5, 0.5, 1.0, 0.5, 0.5}; 98 99 ierr = DMPlexCreateFromDAG(dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr); 100 } 101 break; 102 default: 103 SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Cannot make meshes for dimension %d", dim); 104 } 105 } 106 *newdm = dm; 107 if (refinementLimit > 0.0) { 108 DM rdm; 109 const char *name; 110 111 ierr = DMPlexSetRefinementUniform(*newdm, PETSC_FALSE);CHKERRQ(ierr); 112 ierr = DMPlexSetRefinementLimit(*newdm, refinementLimit);CHKERRQ(ierr); 113 ierr = DMRefine(*newdm, comm, &rdm);CHKERRQ(ierr); 114 ierr = PetscObjectGetName((PetscObject) *newdm, &name);CHKERRQ(ierr); 115 ierr = PetscObjectSetName((PetscObject) rdm, name);CHKERRQ(ierr); 116 ierr = DMDestroy(newdm);CHKERRQ(ierr); 117 *newdm = rdm; 118 } 119 if (interpolate) { 120 DM idm; 121 const char *name; 122 123 ierr = DMPlexInterpolate(*newdm, &idm);CHKERRQ(ierr); 124 ierr = PetscObjectGetName((PetscObject) *newdm, &name);CHKERRQ(ierr); 125 ierr = PetscObjectSetName((PetscObject) idm, name);CHKERRQ(ierr); 126 ierr = DMPlexCopyCoordinates(*newdm, idm);CHKERRQ(ierr); 127 ierr = DMPlexCopyLabels(*newdm, idm);CHKERRQ(ierr); 128 ierr = DMDestroy(newdm);CHKERRQ(ierr); 129 *newdm = idm; 130 } 131 { 132 DM refinedMesh = NULL; 133 DM distributedMesh = NULL; 134 135 /* Distribute mesh over processes */ 136 ierr = DMPlexDistribute(*newdm, NULL, 0, NULL, &distributedMesh);CHKERRQ(ierr); 137 if (distributedMesh) { 138 ierr = DMDestroy(newdm);CHKERRQ(ierr); 139 *newdm = distributedMesh; 140 } 141 if (refinementUniform) { 142 ierr = DMPlexSetRefinementUniform(*newdm, refinementUniform);CHKERRQ(ierr); 143 ierr = DMRefine(*newdm, comm, &refinedMesh);CHKERRQ(ierr); 144 if (refinedMesh) { 145 ierr = DMDestroy(newdm);CHKERRQ(ierr); 146 *newdm = refinedMesh; 147 } 148 } 149 } 150 PetscFunctionReturn(0); 151 } 152 153 #undef __FUNCT__ 154 #define __FUNCT__ "DMPlexCreateSquareBoundary" 155 /*@ 156 DMPlexCreateSquareBoundary - Creates a 1D mesh the is the boundary of a square lattice. 157 158 Collective on MPI_Comm 159 160 Input Parameters: 161 + comm - The communicator for the DM object 162 . lower - The lower left corner coordinates 163 . upper - The upper right corner coordinates 164 - edges - The number of cells in each direction 165 166 Output Parameter: 167 . dm - The DM object 168 169 Note: Here is the numbering returned for 2 cells in each direction: 170 $ 18--5-17--4--16 171 $ | | | 172 $ 6 10 3 173 $ | | | 174 $ 19-11-20--9--15 175 $ | | | 176 $ 7 8 2 177 $ | | | 178 $ 12--0-13--1--14 179 180 Level: beginner 181 182 .keywords: DM, create 183 .seealso: DMPlexCreateBoxMesh(), DMPlexCreateCubeBoundary(), DMSetType(), DMCreate() 184 @*/ 185 PetscErrorCode DMPlexCreateSquareBoundary(DM dm, const PetscReal lower[], const PetscReal upper[], const PetscInt edges[]) 186 { 187 PetscInt numVertices = (edges[0]+1)*(edges[1]+1); 188 PetscInt numEdges = edges[0]*(edges[1]+1) + (edges[0]+1)*edges[1]; 189 PetscInt markerTop = 1; 190 PetscInt markerBottom = 1; 191 PetscInt markerRight = 1; 192 PetscInt markerLeft = 1; 193 PetscBool markerSeparate = PETSC_FALSE; 194 Vec coordinates; 195 PetscSection coordSection; 196 PetscScalar *coords; 197 PetscInt coordSize; 198 PetscMPIInt rank; 199 PetscInt v, vx, vy; 200 PetscErrorCode ierr; 201 202 PetscFunctionBegin; 203 ierr = PetscOptionsGetBool(((PetscObject) dm)->prefix, "-dm_plex_separate_marker", &markerSeparate, NULL);CHKERRQ(ierr); 204 if (markerSeparate) { 205 markerTop = 1; 206 markerBottom = 0; 207 markerRight = 0; 208 markerLeft = 0; 209 } 210 ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);CHKERRQ(ierr); 211 if (!rank) { 212 PetscInt e, ex, ey; 213 214 ierr = DMPlexSetChart(dm, 0, numEdges+numVertices);CHKERRQ(ierr); 215 for (e = 0; e < numEdges; ++e) { 216 ierr = DMPlexSetConeSize(dm, e, 2);CHKERRQ(ierr); 217 } 218 ierr = DMSetUp(dm);CHKERRQ(ierr); /* Allocate space for cones */ 219 for (vx = 0; vx <= edges[0]; vx++) { 220 for (ey = 0; ey < edges[1]; ey++) { 221 PetscInt edge = vx*edges[1] + ey + edges[0]*(edges[1]+1); 222 PetscInt vertex = ey*(edges[0]+1) + vx + numEdges; 223 PetscInt cone[2]; 224 225 cone[0] = vertex; cone[1] = vertex+edges[0]+1; 226 ierr = DMPlexSetCone(dm, edge, cone);CHKERRQ(ierr); 227 if (vx == edges[0]) { 228 ierr = DMPlexSetLabelValue(dm, "marker", edge, markerRight);CHKERRQ(ierr); 229 ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerRight);CHKERRQ(ierr); 230 if (ey == edges[1]-1) { 231 ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerRight);CHKERRQ(ierr); 232 } 233 } else if (vx == 0) { 234 ierr = DMPlexSetLabelValue(dm, "marker", edge, markerLeft);CHKERRQ(ierr); 235 ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerLeft);CHKERRQ(ierr); 236 if (ey == edges[1]-1) { 237 ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerLeft);CHKERRQ(ierr); 238 } 239 } 240 } 241 } 242 for (vy = 0; vy <= edges[1]; vy++) { 243 for (ex = 0; ex < edges[0]; ex++) { 244 PetscInt edge = vy*edges[0] + ex; 245 PetscInt vertex = vy*(edges[0]+1) + ex + numEdges; 246 PetscInt cone[2]; 247 248 cone[0] = vertex; cone[1] = vertex+1; 249 ierr = DMPlexSetCone(dm, edge, cone);CHKERRQ(ierr); 250 if (vy == edges[1]) { 251 ierr = DMPlexSetLabelValue(dm, "marker", edge, markerTop);CHKERRQ(ierr); 252 ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerTop);CHKERRQ(ierr); 253 if (ex == edges[0]-1) { 254 ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerTop);CHKERRQ(ierr); 255 } 256 } else if (vy == 0) { 257 ierr = DMPlexSetLabelValue(dm, "marker", edge, markerBottom);CHKERRQ(ierr); 258 ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerBottom);CHKERRQ(ierr); 259 if (ex == edges[0]-1) { 260 ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerBottom);CHKERRQ(ierr); 261 } 262 } 263 } 264 } 265 } 266 ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr); 267 ierr = DMPlexStratify(dm);CHKERRQ(ierr); 268 /* Build coordinates */ 269 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 270 ierr = PetscSectionSetChart(coordSection, numEdges, numEdges + numVertices);CHKERRQ(ierr); 271 for (v = numEdges; v < numEdges+numVertices; ++v) { 272 ierr = PetscSectionSetDof(coordSection, v, 2);CHKERRQ(ierr); 273 } 274 ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 275 ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 276 ierr = VecCreate(PetscObjectComm((PetscObject)dm), &coordinates);CHKERRQ(ierr); 277 ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 278 ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr); 279 ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 280 for (vy = 0; vy <= edges[1]; ++vy) { 281 for (vx = 0; vx <= edges[0]; ++vx) { 282 coords[(vy*(edges[0]+1)+vx)*2+0] = lower[0] + ((upper[0] - lower[0])/edges[0])*vx; 283 coords[(vy*(edges[0]+1)+vx)*2+1] = lower[1] + ((upper[1] - lower[1])/edges[1])*vy; 284 } 285 } 286 ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 287 ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr); 288 ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 289 PetscFunctionReturn(0); 290 } 291 292 #undef __FUNCT__ 293 #define __FUNCT__ "DMPlexCreateCubeBoundary" 294 /*@ 295 DMPlexCreateCubeBoundary - Creates a 2D mesh the is the boundary of a cubic lattice. 296 297 Collective on MPI_Comm 298 299 Input Parameters: 300 + comm - The communicator for the DM object 301 . lower - The lower left front corner coordinates 302 . upper - The upper right back corner coordinates 303 - edges - The number of cells in each direction 304 305 Output Parameter: 306 . dm - The DM object 307 308 Level: beginner 309 310 .keywords: DM, create 311 .seealso: DMPlexCreateBoxMesh(), DMPlexCreateSquareBoundary(), DMSetType(), DMCreate() 312 @*/ 313 PetscErrorCode DMPlexCreateCubeBoundary(DM dm, const PetscReal lower[], const PetscReal upper[], const PetscInt faces[]) 314 { 315 PetscInt numVertices = (faces[0]+1)*(faces[1]+1)*(faces[2]+1); 316 PetscInt numFaces = 6; 317 Vec coordinates; 318 PetscSection coordSection; 319 PetscScalar *coords; 320 PetscInt coordSize; 321 PetscMPIInt rank; 322 PetscInt v, vx, vy, vz; 323 PetscErrorCode ierr; 324 325 PetscFunctionBegin; 326 if ((faces[0] < 1) || (faces[1] < 1) || (faces[2] < 1)) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Must have at least 1 face per side"); 327 if ((faces[0] > 1) || (faces[1] > 1) || (faces[2] > 1)) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Currently can't handle more than 1 face per side"); 328 ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);CHKERRQ(ierr); 329 if (!rank) { 330 PetscInt f; 331 332 ierr = DMPlexSetChart(dm, 0, numFaces+numVertices);CHKERRQ(ierr); 333 for (f = 0; f < numFaces; ++f) { 334 ierr = DMPlexSetConeSize(dm, f, 4);CHKERRQ(ierr); 335 } 336 ierr = DMSetUp(dm);CHKERRQ(ierr); /* Allocate space for cones */ 337 for (v = 0; v < numFaces+numVertices; ++v) { 338 ierr = DMPlexSetLabelValue(dm, "marker", v, 1);CHKERRQ(ierr); 339 } 340 { /* Side 0 (Front) */ 341 PetscInt cone[4]; 342 cone[0] = numFaces+4; cone[1] = numFaces+5; cone[2] = numFaces+7; cone[3] = numFaces+6; 343 ierr = DMPlexSetCone(dm, 0, cone);CHKERRQ(ierr); 344 } 345 { /* Side 1 (Back) */ 346 PetscInt cone[4]; 347 cone[0] = numFaces+1; cone[1] = numFaces+0; cone[2] = numFaces+2; cone[3] = numFaces+3; 348 ierr = DMPlexSetCone(dm, 1, cone);CHKERRQ(ierr); 349 } 350 { /* Side 2 (Bottom) */ 351 PetscInt cone[4]; 352 cone[0] = numFaces+0; cone[1] = numFaces+1; cone[2] = numFaces+5; cone[3] = numFaces+4; 353 ierr = DMPlexSetCone(dm, 2, cone);CHKERRQ(ierr); 354 } 355 { /* Side 3 (Top) */ 356 PetscInt cone[4]; 357 cone[0] = numFaces+6; cone[1] = numFaces+7; cone[2] = numFaces+3; cone[3] = numFaces+2; 358 ierr = DMPlexSetCone(dm, 3, cone);CHKERRQ(ierr); 359 } 360 { /* Side 4 (Left) */ 361 PetscInt cone[4]; 362 cone[0] = numFaces+0; cone[1] = numFaces+4; cone[2] = numFaces+6; cone[3] = numFaces+2; 363 ierr = DMPlexSetCone(dm, 4, cone);CHKERRQ(ierr); 364 } 365 { /* Side 5 (Right) */ 366 PetscInt cone[4]; 367 cone[0] = numFaces+5; cone[1] = numFaces+1; cone[2] = numFaces+3; cone[3] = numFaces+7; 368 ierr = DMPlexSetCone(dm, 5, cone);CHKERRQ(ierr); 369 } 370 } 371 ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr); 372 ierr = DMPlexStratify(dm);CHKERRQ(ierr); 373 /* Build coordinates */ 374 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 375 ierr = PetscSectionSetChart(coordSection, numFaces, numFaces + numVertices);CHKERRQ(ierr); 376 for (v = numFaces; v < numFaces+numVertices; ++v) { 377 ierr = PetscSectionSetDof(coordSection, v, 3);CHKERRQ(ierr); 378 } 379 ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 380 ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 381 ierr = VecCreate(PetscObjectComm((PetscObject)dm), &coordinates);CHKERRQ(ierr); 382 ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr); 383 ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 384 ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr); 385 ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 386 for (vz = 0; vz <= faces[2]; ++vz) { 387 for (vy = 0; vy <= faces[1]; ++vy) { 388 for (vx = 0; vx <= faces[0]; ++vx) { 389 coords[((vz*(faces[1]+1)+vy)*(faces[0]+1)+vx)*3+0] = lower[0] + ((upper[0] - lower[0])/faces[0])*vx; 390 coords[((vz*(faces[1]+1)+vy)*(faces[0]+1)+vx)*3+1] = lower[1] + ((upper[1] - lower[1])/faces[1])*vy; 391 coords[((vz*(faces[1]+1)+vy)*(faces[0]+1)+vx)*3+2] = lower[2] + ((upper[2] - lower[2])/faces[2])*vz; 392 } 393 } 394 } 395 ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 396 ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr); 397 ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 398 PetscFunctionReturn(0); 399 } 400 401 #undef __FUNCT__ 402 #define __FUNCT__ "DMPlexCreateSquareMesh" 403 /*@ 404 DMPlexCreateSquareMesh - Creates a 2D mesh for a square lattice. 405 406 Collective on MPI_Comm 407 408 Input Parameters: 409 + comm - The communicator for the DM object 410 . lower - The lower left corner coordinates 411 . upper - The upper right corner coordinates 412 . edges - The number of cells in each direction 413 . bdX - The boundary type for the X direction 414 - bdY - The boundary type for the Y direction 415 416 Output Parameter: 417 . dm - The DM object 418 419 Note: Here is the numbering returned for 2 cells in each direction: 420 $ 22--8-23--9--24 421 $ | | | 422 $ 13 2 14 3 15 423 $ | | | 424 $ 19--6-20--7--21 425 $ | | | 426 $ 10 0 11 1 12 427 $ | | | 428 $ 16--4-17--5--18 429 430 Level: beginner 431 432 .keywords: DM, create 433 .seealso: DMPlexCreateBoxMesh(), DMPlexCreateSquareBoundary(), DMPlexCreateCubeBoundary(), DMSetType(), DMCreate() 434 @*/ 435 PetscErrorCode DMPlexCreateSquareMesh(DM dm, const PetscReal lower[], const PetscReal upper[], const PetscInt edges[], DMBoundaryType bdX, DMBoundaryType bdY) 436 { 437 PetscInt markerTop = 1; 438 PetscInt markerBottom = 1; 439 PetscInt markerRight = 1; 440 PetscInt markerLeft = 1; 441 PetscBool markerSeparate = PETSC_FALSE; 442 PetscMPIInt rank; 443 PetscErrorCode ierr; 444 445 PetscFunctionBegin; 446 ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);CHKERRQ(ierr); 447 ierr = PetscOptionsGetBool(((PetscObject) dm)->prefix, "-dm_plex_separate_marker", &markerSeparate, NULL);CHKERRQ(ierr); 448 if (markerSeparate) { 449 markerTop = 3; 450 markerBottom = 1; 451 markerRight = 2; 452 markerLeft = 4; 453 } 454 { 455 const PetscInt numXEdges = !rank ? edges[0] : 0; 456 const PetscInt numYEdges = !rank ? edges[1] : 0; 457 const PetscInt numXVertices = !rank ? (bdX == DM_BOUNDARY_PERIODIC || bdX == DM_BOUNDARY_TWIST ? edges[0] : edges[0]+1) : 0; 458 const PetscInt numYVertices = !rank ? (bdY == DM_BOUNDARY_PERIODIC || bdY == DM_BOUNDARY_TWIST ? edges[1] : edges[1]+1) : 0; 459 const PetscInt numTotXEdges = numXEdges*numYVertices; 460 const PetscInt numTotYEdges = numYEdges*numXVertices; 461 const PetscInt numVertices = numXVertices*numYVertices; 462 const PetscInt numEdges = numTotXEdges + numTotYEdges; 463 const PetscInt numFaces = numXEdges*numYEdges; 464 const PetscInt firstVertex = numFaces; 465 const PetscInt firstXEdge = numFaces + numVertices; 466 const PetscInt firstYEdge = numFaces + numVertices + numTotXEdges; 467 Vec coordinates; 468 PetscSection coordSection; 469 PetscScalar *coords; 470 PetscInt coordSize; 471 PetscInt v, vx, vy; 472 PetscInt f, fx, fy, e, ex, ey; 473 474 ierr = DMPlexSetChart(dm, 0, numFaces+numEdges+numVertices);CHKERRQ(ierr); 475 for (f = 0; f < numFaces; ++f) { 476 ierr = DMPlexSetConeSize(dm, f, 4);CHKERRQ(ierr); 477 } 478 for (e = firstXEdge; e < firstXEdge+numEdges; ++e) { 479 ierr = DMPlexSetConeSize(dm, e, 2);CHKERRQ(ierr); 480 } 481 ierr = DMSetUp(dm);CHKERRQ(ierr); /* Allocate space for cones */ 482 /* Build faces */ 483 for (fy = 0; fy < numYEdges; fy++) { 484 for (fx = 0; fx < numXEdges; fx++) { 485 const PetscInt face = fy*numXEdges + fx; 486 const PetscInt edgeL = firstYEdge + fx *numYEdges + fy; 487 const PetscInt edgeR = firstYEdge + ((fx+1)%numXVertices)*numYEdges + fy; 488 const PetscInt edgeB = firstXEdge + fy *numXEdges + fx; 489 const PetscInt edgeT = firstXEdge + ((fy+1)%numYVertices)*numXEdges + fx; 490 const PetscInt ornt[4] = {0, 0, -2, -2}; 491 PetscInt cone[4]; 492 493 cone[0] = edgeB; cone[1] = edgeR; cone[2] = edgeT; cone[3] = edgeL; 494 ierr = DMPlexSetCone(dm, face, cone);CHKERRQ(ierr); 495 ierr = DMPlexSetConeOrientation(dm, face, ornt);CHKERRQ(ierr); 496 } 497 } 498 /* Build Y edges*/ 499 for (vx = 0; vx < numXVertices; vx++) { 500 for (ey = 0; ey < numYEdges; ey++) { 501 const PetscInt nextv = (bdY == DM_BOUNDARY_TWIST && ey == numYEdges-1) ? (numXVertices-vx-1) : ((ey+1)%numYVertices)*numXVertices + vx; 502 const PetscInt edge = firstYEdge + vx*numYEdges + ey; 503 const PetscInt vertexB = firstVertex + ey*numXVertices + vx; 504 const PetscInt vertexT = firstVertex + nextv; 505 PetscInt cone[2]; 506 507 cone[0] = vertexB; cone[1] = vertexT; 508 ierr = DMPlexSetCone(dm, edge, cone);CHKERRQ(ierr); 509 if ((bdX != DM_BOUNDARY_PERIODIC) && (bdX != DM_BOUNDARY_TWIST)) { 510 if (vx == numXVertices-1) { 511 ierr = DMPlexSetLabelValue(dm, "marker", edge, markerRight);CHKERRQ(ierr); 512 ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerRight);CHKERRQ(ierr); 513 if (ey == numYEdges-1) { 514 ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerRight);CHKERRQ(ierr); 515 } 516 } else if (vx == 0) { 517 ierr = DMPlexSetLabelValue(dm, "marker", edge, markerLeft);CHKERRQ(ierr); 518 ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerLeft);CHKERRQ(ierr); 519 if (ey == numYEdges-1) { 520 ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerLeft);CHKERRQ(ierr); 521 } 522 } 523 } 524 } 525 } 526 /* Build X edges*/ 527 for (vy = 0; vy < numYVertices; vy++) { 528 for (ex = 0; ex < numXEdges; ex++) { 529 const PetscInt nextv = (bdX == DM_BOUNDARY_TWIST && ex == numXEdges-1) ? (numYVertices-vy-1)*numXVertices : vy*numXVertices + (ex+1)%numXVertices; 530 const PetscInt edge = firstXEdge + vy*numXEdges + ex; 531 const PetscInt vertexL = firstVertex + vy*numXVertices + ex; 532 const PetscInt vertexR = firstVertex + nextv; 533 PetscInt cone[2]; 534 535 cone[0] = vertexL; cone[1] = vertexR; 536 ierr = DMPlexSetCone(dm, edge, cone);CHKERRQ(ierr); 537 if ((bdY != DM_BOUNDARY_PERIODIC) && (bdY != DM_BOUNDARY_TWIST)) { 538 if (vy == numYVertices-1) { 539 ierr = DMPlexSetLabelValue(dm, "marker", edge, markerTop);CHKERRQ(ierr); 540 ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerTop);CHKERRQ(ierr); 541 if (ex == numXEdges-1) { 542 ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerTop);CHKERRQ(ierr); 543 } 544 } else if (vy == 0) { 545 ierr = DMPlexSetLabelValue(dm, "marker", edge, markerBottom);CHKERRQ(ierr); 546 ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerBottom);CHKERRQ(ierr); 547 if (ex == numXEdges-1) { 548 ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerBottom);CHKERRQ(ierr); 549 } 550 } 551 } 552 } 553 } 554 ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr); 555 ierr = DMPlexStratify(dm);CHKERRQ(ierr); 556 /* Build coordinates */ 557 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 558 ierr = PetscSectionSetChart(coordSection, firstVertex, firstVertex+numVertices);CHKERRQ(ierr); 559 for (v = firstVertex; v < firstVertex+numVertices; ++v) { 560 ierr = PetscSectionSetDof(coordSection, v, 2);CHKERRQ(ierr); 561 } 562 ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 563 ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 564 ierr = VecCreate(PetscObjectComm((PetscObject)dm), &coordinates);CHKERRQ(ierr); 565 ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 566 ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr); 567 ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 568 for (vy = 0; vy < numYVertices; ++vy) { 569 for (vx = 0; vx < numXVertices; ++vx) { 570 coords[(vy*numXVertices+vx)*2+0] = lower[0] + ((upper[0] - lower[0])/numXEdges)*vx; 571 coords[(vy*numXVertices+vx)*2+1] = lower[1] + ((upper[1] - lower[1])/numYEdges)*vy; 572 } 573 } 574 ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 575 ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr); 576 ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 577 } 578 PetscFunctionReturn(0); 579 } 580 581 #undef __FUNCT__ 582 #define __FUNCT__ "DMPlexCreateBoxMesh" 583 /*@ 584 DMPlexCreateBoxMesh - Creates a mesh on the tensor product of unit intervals (box) using simplices. 585 586 Collective on MPI_Comm 587 588 Input Parameters: 589 + comm - The communicator for the DM object 590 . dim - The spatial dimension 591 - interpolate - Flag to create intermediate mesh pieces (edges, faces) 592 593 Output Parameter: 594 . dm - The DM object 595 596 Level: beginner 597 598 .keywords: DM, create 599 .seealso: DMPlexCreateHexBoxMesh(), DMSetType(), DMCreate() 600 @*/ 601 PetscErrorCode DMPlexCreateBoxMesh(MPI_Comm comm, PetscInt dim, PetscBool interpolate, DM *dm) 602 { 603 DM boundary; 604 PetscErrorCode ierr; 605 606 PetscFunctionBegin; 607 PetscValidPointer(dm, 4); 608 ierr = DMCreate(comm, &boundary);CHKERRQ(ierr); 609 PetscValidLogicalCollectiveInt(boundary,dim,2); 610 ierr = DMSetType(boundary, DMPLEX);CHKERRQ(ierr); 611 ierr = DMPlexSetDimension(boundary, dim-1);CHKERRQ(ierr); 612 switch (dim) { 613 case 2: 614 { 615 PetscReal lower[2] = {0.0, 0.0}; 616 PetscReal upper[2] = {1.0, 1.0}; 617 PetscInt edges[2] = {2, 2}; 618 619 ierr = DMPlexCreateSquareBoundary(boundary, lower, upper, edges);CHKERRQ(ierr); 620 break; 621 } 622 case 3: 623 { 624 PetscReal lower[3] = {0.0, 0.0, 0.0}; 625 PetscReal upper[3] = {1.0, 1.0, 1.0}; 626 PetscInt faces[3] = {1, 1, 1}; 627 628 ierr = DMPlexCreateCubeBoundary(boundary, lower, upper, faces);CHKERRQ(ierr); 629 break; 630 } 631 default: 632 SETERRQ1(comm, PETSC_ERR_SUP, "Dimension not supported: %d", dim); 633 } 634 ierr = DMPlexGenerate(boundary, NULL, interpolate, dm);CHKERRQ(ierr); 635 ierr = DMDestroy(&boundary);CHKERRQ(ierr); 636 PetscFunctionReturn(0); 637 } 638 639 #undef __FUNCT__ 640 #define __FUNCT__ "DMPlexCreateHexBoxMesh" 641 /*@ 642 DMPlexCreateHexBoxMesh - Creates a mesh on the tensor product of unit intervals (box) using hexahedra. 643 644 Collective on MPI_Comm 645 646 Input Parameters: 647 + comm - The communicator for the DM object 648 . dim - The spatial dimension 649 . periodicX - The boundary type for the X direction 650 . periodicY - The boundary type for the Y direction 651 . periodicZ - The boundary type for the Z direction 652 - cells - The number of cells in each direction 653 654 Output Parameter: 655 . dm - The DM object 656 657 Level: beginner 658 659 .keywords: DM, create 660 .seealso: DMPlexCreateBoxMesh(), DMSetType(), DMCreate() 661 @*/ 662 PetscErrorCode DMPlexCreateHexBoxMesh(MPI_Comm comm, PetscInt dim, const PetscInt cells[], DMBoundaryType periodicX, DMBoundaryType periodicY, DMBoundaryType periodicZ, DM *dm) 663 { 664 PetscErrorCode ierr; 665 666 PetscFunctionBegin; 667 PetscValidPointer(dm, 4); 668 ierr = DMCreate(comm, dm);CHKERRQ(ierr); 669 PetscValidLogicalCollectiveInt(*dm,dim,2); 670 ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 671 ierr = DMPlexSetDimension(*dm, dim);CHKERRQ(ierr); 672 switch (dim) { 673 case 2: 674 { 675 PetscReal lower[2] = {0.0, 0.0}; 676 PetscReal upper[2] = {1.0, 1.0}; 677 678 ierr = DMPlexCreateSquareMesh(*dm, lower, upper, cells, periodicX, periodicY);CHKERRQ(ierr); 679 break; 680 } 681 #if 0 682 case 3: 683 { 684 PetscReal lower[3] = {0.0, 0.0, 0.0}; 685 PetscReal upper[3] = {1.0, 1.0, 1.0}; 686 687 ierr = DMPlexCreateCubeMesh(boundary, lower, upper, cells, periodicX, periodicY, periodicZ);CHKERRQ(ierr); 688 break; 689 } 690 #endif 691 default: 692 SETERRQ1(comm, PETSC_ERR_SUP, "Dimension not supported: %d", dim); 693 } 694 PetscFunctionReturn(0); 695 } 696 697 /* External function declarations here */ 698 extern PetscErrorCode DMCreateInterpolation_Plex(DM dmCoarse, DM dmFine, Mat *interpolation, Vec *scaling); 699 extern PetscErrorCode DMCreateInjection_Plex(DM dmCoarse, DM dmFine, VecScatter *ctx); 700 extern PetscErrorCode DMCreateMatrix_Plex(DM dm, Mat *J); 701 extern PetscErrorCode DMCreateCoordinateDM_Plex(DM dm, DM *cdm); 702 extern PetscErrorCode DMRefine_Plex(DM dm, MPI_Comm comm, DM *dmRefined); 703 extern PetscErrorCode DMCoarsen_Plex(DM dm, MPI_Comm comm, DM *dmCoarsened); 704 extern PetscErrorCode DMRefineHierarchy_Plex(DM dm, PetscInt nlevels, DM dmRefined[]); 705 extern PetscErrorCode DMClone_Plex(DM dm, DM *newdm); 706 extern PetscErrorCode DMSetUp_Plex(DM dm); 707 extern PetscErrorCode DMDestroy_Plex(DM dm); 708 extern PetscErrorCode DMView_Plex(DM dm, PetscViewer viewer); 709 extern PetscErrorCode DMCreateSubDM_Plex(DM dm, PetscInt numFields, PetscInt fields[], IS *is, DM *subdm); 710 extern PetscErrorCode DMLocatePoints_Plex(DM dm, Vec v, IS *cellIS); 711 712 #undef __FUNCT__ 713 #define __FUNCT__ "DMPlexReplace_Static" 714 /* Replace dm with the contents of dmNew 715 - Share the DM_Plex structure 716 - Share the coordinates 717 */ 718 static PetscErrorCode DMPlexReplace_Static(DM dm, DM dmNew) 719 { 720 PetscSection coordSection; 721 Vec coords; 722 PetscErrorCode ierr; 723 724 PetscFunctionBegin; 725 ierr = DMGetCoordinateSection(dmNew, &coordSection);CHKERRQ(ierr); 726 ierr = DMGetCoordinatesLocal(dmNew, &coords);CHKERRQ(ierr); 727 ierr = DMSetCoordinateSection(dm, coordSection);CHKERRQ(ierr); 728 ierr = DMSetCoordinatesLocal(dm, coords);CHKERRQ(ierr); 729 ierr = DMDestroy_Plex(dm);CHKERRQ(ierr); 730 dm->data = dmNew->data; 731 ((DM_Plex *) dmNew->data)->refct++; 732 PetscFunctionReturn(0); 733 } 734 735 #undef __FUNCT__ 736 #define __FUNCT__ "DMPlexSwap_Static" 737 /* Swap dm with the contents of dmNew 738 - Swap the DM_Plex structure 739 - Swap the coordinates 740 */ 741 static PetscErrorCode DMPlexSwap_Static(DM dmA, DM dmB) 742 { 743 DM coordDMA, coordDMB; 744 Vec coordsA, coordsB; 745 void *tmp; 746 PetscErrorCode ierr; 747 748 PetscFunctionBegin; 749 ierr = DMGetCoordinateDM(dmA, &coordDMA);CHKERRQ(ierr); 750 ierr = DMGetCoordinateDM(dmB, &coordDMB);CHKERRQ(ierr); 751 ierr = PetscObjectReference((PetscObject) coordDMA);CHKERRQ(ierr); 752 ierr = DMSetCoordinateDM(dmA, coordDMB);CHKERRQ(ierr); 753 ierr = DMSetCoordinateDM(dmB, coordDMA);CHKERRQ(ierr); 754 ierr = PetscObjectDereference((PetscObject) coordDMA);CHKERRQ(ierr); 755 756 ierr = DMGetCoordinatesLocal(dmA, &coordsA);CHKERRQ(ierr); 757 ierr = DMGetCoordinatesLocal(dmB, &coordsB);CHKERRQ(ierr); 758 ierr = PetscObjectReference((PetscObject) coordsA);CHKERRQ(ierr); 759 ierr = DMSetCoordinatesLocal(dmA, coordsB);CHKERRQ(ierr); 760 ierr = DMSetCoordinatesLocal(dmB, coordsA);CHKERRQ(ierr); 761 ierr = PetscObjectDereference((PetscObject) coordsA);CHKERRQ(ierr); 762 tmp = dmA->data; 763 dmA->data = dmB->data; 764 dmB->data = tmp; 765 PetscFunctionReturn(0); 766 } 767 768 #undef __FUNCT__ 769 #define __FUNCT__ "DMSetFromOptions_Plex" 770 PetscErrorCode DMSetFromOptions_Plex(DM dm) 771 { 772 DM_Plex *mesh = (DM_Plex*) dm->data; 773 PetscInt refine = 0, r; 774 PetscBool isHierarchy; 775 PetscErrorCode ierr; 776 777 PetscFunctionBegin; 778 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 779 ierr = PetscOptionsHead("DMPlex Options");CHKERRQ(ierr); 780 /* Handle DMPlex refinement */ 781 ierr = PetscOptionsInt("-dm_refine", "The number of uniform refinements", "DMCreate", refine, &refine, NULL);CHKERRQ(ierr); 782 ierr = PetscOptionsInt("-dm_refine_hierarchy", "The number of uniform refinements", "DMCreate", refine, &refine, &isHierarchy);CHKERRQ(ierr); 783 ierr = DMPlexSetRefinementUniform(dm, refine ? PETSC_TRUE : PETSC_FALSE);CHKERRQ(ierr); 784 if (refine && isHierarchy) { 785 DM *dms; 786 787 ierr = PetscMalloc1(refine,&dms);CHKERRQ(ierr); 788 ierr = DMRefineHierarchy(dm, refine, dms);CHKERRQ(ierr); 789 /* Total hack since we do not pass in a pointer */ 790 ierr = DMPlexSwap_Static(dm, dms[refine-1]);CHKERRQ(ierr); 791 if (refine == 1) { 792 ierr = DMPlexSetCoarseDM(dm, dms[0]);CHKERRQ(ierr); 793 } else { 794 ierr = DMPlexSetCoarseDM(dm, dms[refine-2]);CHKERRQ(ierr); 795 ierr = DMPlexSetCoarseDM(dms[0], dms[refine-1]);CHKERRQ(ierr); 796 } 797 /* Free DMs */ 798 for (r = 0; r < refine; ++r) {ierr = DMDestroy(&dms[r]);CHKERRQ(ierr);} 799 ierr = PetscFree(dms);CHKERRQ(ierr); 800 } else { 801 for (r = 0; r < refine; ++r) { 802 DM refinedMesh; 803 804 ierr = DMRefine(dm, PetscObjectComm((PetscObject) dm), &refinedMesh);CHKERRQ(ierr); 805 /* Total hack since we do not pass in a pointer */ 806 ierr = DMPlexReplace_Static(dm, refinedMesh);CHKERRQ(ierr); 807 ierr = DMDestroy(&refinedMesh);CHKERRQ(ierr); 808 } 809 } 810 /* Handle associated vectors */ 811 /* Handle viewing */ 812 ierr = PetscOptionsBool("-dm_plex_print_set_values", "Output all set values info", "DMView", PETSC_FALSE, &mesh->printSetValues, NULL);CHKERRQ(ierr); 813 ierr = PetscOptionsInt("-dm_plex_print_fem", "Debug output level all fem computations", "DMView", 0, &mesh->printFEM, NULL);CHKERRQ(ierr); 814 ierr = PetscOptionsReal("-dm_plex_print_tol", "Tolerance for FEM output", "DMView", PETSC_FALSE, &mesh->printTol, NULL);CHKERRQ(ierr); 815 ierr = PetscOptionsTail();CHKERRQ(ierr); 816 PetscFunctionReturn(0); 817 } 818 819 #undef __FUNCT__ 820 #define __FUNCT__ "DMCreateGlobalVector_Plex" 821 static PetscErrorCode DMCreateGlobalVector_Plex(DM dm,Vec *vec) 822 { 823 PetscErrorCode ierr; 824 825 PetscFunctionBegin; 826 ierr = DMCreateGlobalVector_Section_Private(dm,vec);CHKERRQ(ierr); 827 /* ierr = VecSetOperation(*vec, VECOP_DUPLICATE, (void(*)(void)) VecDuplicate_MPI_DM);CHKERRQ(ierr); */ 828 ierr = VecSetOperation(*vec, VECOP_VIEW, (void (*)(void))VecView_Plex);CHKERRQ(ierr); 829 PetscFunctionReturn(0); 830 } 831 832 #undef __FUNCT__ 833 #define __FUNCT__ "DMCreateLocalVector_Plex" 834 static PetscErrorCode DMCreateLocalVector_Plex(DM dm,Vec *vec) 835 { 836 PetscErrorCode ierr; 837 838 PetscFunctionBegin; 839 ierr = DMCreateLocalVector_Section_Private(dm,vec);CHKERRQ(ierr); 840 ierr = VecSetOperation(*vec, VECOP_VIEW, (void(*)(void)) VecView_Plex_Local);CHKERRQ(ierr); 841 PetscFunctionReturn(0); 842 } 843 844 #undef __FUNCT__ 845 #define __FUNCT__ "DMInitialize_Plex" 846 PetscErrorCode DMInitialize_Plex(DM dm) 847 { 848 PetscFunctionBegin; 849 dm->ops->view = DMView_Plex; 850 dm->ops->setfromoptions = DMSetFromOptions_Plex; 851 dm->ops->clone = DMClone_Plex; 852 dm->ops->setup = DMSetUp_Plex; 853 dm->ops->createglobalvector = DMCreateGlobalVector_Plex; 854 dm->ops->createlocalvector = DMCreateLocalVector_Plex; 855 dm->ops->getlocaltoglobalmapping = NULL; 856 dm->ops->getlocaltoglobalmappingblock = NULL; 857 dm->ops->createfieldis = NULL; 858 dm->ops->createcoordinatedm = DMCreateCoordinateDM_Plex; 859 dm->ops->getcoloring = NULL; 860 dm->ops->creatematrix = DMCreateMatrix_Plex; 861 dm->ops->createinterpolation = DMCreateInterpolation_Plex; 862 dm->ops->getaggregates = NULL; 863 dm->ops->getinjection = DMCreateInjection_Plex; 864 dm->ops->refine = DMRefine_Plex; 865 dm->ops->coarsen = DMCoarsen_Plex; 866 dm->ops->refinehierarchy = DMRefineHierarchy_Plex; 867 dm->ops->coarsenhierarchy = NULL; 868 dm->ops->globaltolocalbegin = NULL; 869 dm->ops->globaltolocalend = NULL; 870 dm->ops->localtoglobalbegin = NULL; 871 dm->ops->localtoglobalend = NULL; 872 dm->ops->destroy = DMDestroy_Plex; 873 dm->ops->createsubdm = DMCreateSubDM_Plex; 874 dm->ops->locatepoints = DMLocatePoints_Plex; 875 PetscFunctionReturn(0); 876 } 877 878 #undef __FUNCT__ 879 #define __FUNCT__ "DMClone_Plex" 880 PetscErrorCode DMClone_Plex(DM dm, DM *newdm) 881 { 882 DM_Plex *mesh = (DM_Plex *) dm->data; 883 PetscErrorCode ierr; 884 885 PetscFunctionBegin; 886 mesh->refct++; 887 (*newdm)->data = mesh; 888 ierr = PetscObjectChangeTypeName((PetscObject) *newdm, DMPLEX);CHKERRQ(ierr); 889 ierr = DMInitialize_Plex(*newdm);CHKERRQ(ierr); 890 PetscFunctionReturn(0); 891 } 892 893 /*MC 894 DMPLEX = "plex" - A DM object that encapsulates an unstructured mesh, or CW Complex, which can be expressed using a Hasse Diagram. 895 In the local representation, Vecs contain all unknowns in the interior and shared boundary. This is 896 specified by a PetscSection object. Ownership in the global representation is determined by 897 ownership of the underlying DMPlex points. This is specified by another PetscSection object. 898 899 Level: intermediate 900 901 .seealso: DMType, DMPlexCreate(), DMCreate(), DMSetType() 902 M*/ 903 904 #undef __FUNCT__ 905 #define __FUNCT__ "DMCreate_Plex" 906 PETSC_EXTERN PetscErrorCode DMCreate_Plex(DM dm) 907 { 908 DM_Plex *mesh; 909 PetscInt unit, d; 910 PetscErrorCode ierr; 911 912 PetscFunctionBegin; 913 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 914 ierr = PetscNewLog(dm,&mesh);CHKERRQ(ierr); 915 dm->data = mesh; 916 917 mesh->refct = 1; 918 mesh->dim = 0; 919 ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm), &mesh->coneSection);CHKERRQ(ierr); 920 mesh->maxConeSize = 0; 921 mesh->cones = NULL; 922 mesh->coneOrientations = NULL; 923 ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm), &mesh->supportSection);CHKERRQ(ierr); 924 mesh->maxSupportSize = 0; 925 mesh->supports = NULL; 926 mesh->refinementUniform = PETSC_TRUE; 927 mesh->refinementLimit = -1.0; 928 929 mesh->facesTmp = NULL; 930 931 mesh->subpointMap = NULL; 932 933 for (unit = 0; unit < NUM_PETSC_UNITS; ++unit) mesh->scale[unit] = 1.0; 934 935 mesh->labels = NULL; 936 mesh->depthLabel = NULL; 937 mesh->globalVertexNumbers = NULL; 938 mesh->globalCellNumbers = NULL; 939 for (d = 0; d < 8; ++d) mesh->hybridPointMax[d] = PETSC_DETERMINE; 940 mesh->vtkCellHeight = 0; 941 mesh->preallocCenterDim = -1; 942 943 mesh->printSetValues = PETSC_FALSE; 944 mesh->printFEM = 0; 945 mesh->printTol = 1.0e-10; 946 947 ierr = DMInitialize_Plex(dm);CHKERRQ(ierr); 948 PetscFunctionReturn(0); 949 } 950 951 #undef __FUNCT__ 952 #define __FUNCT__ "DMPlexCreate" 953 /*@ 954 DMPlexCreate - Creates a DMPlex object, which encapsulates an unstructured mesh, or CW complex, which can be expressed using a Hasse Diagram. 955 956 Collective on MPI_Comm 957 958 Input Parameter: 959 . comm - The communicator for the DMPlex object 960 961 Output Parameter: 962 . mesh - The DMPlex object 963 964 Level: beginner 965 966 .keywords: DMPlex, create 967 @*/ 968 PetscErrorCode DMPlexCreate(MPI_Comm comm, DM *mesh) 969 { 970 PetscErrorCode ierr; 971 972 PetscFunctionBegin; 973 PetscValidPointer(mesh,2); 974 ierr = DMCreate(comm, mesh);CHKERRQ(ierr); 975 ierr = DMSetType(*mesh, DMPLEX);CHKERRQ(ierr); 976 PetscFunctionReturn(0); 977 } 978 979 #undef __FUNCT__ 980 #define __FUNCT__ "DMPlexBuildFromCellList_Private" 981 /* 982 This takes as input the common mesh generator output, a list of the vertices for each cell 983 */ 984 PetscErrorCode DMPlexBuildFromCellList_Private(DM dm, PetscInt numCells, PetscInt numVertices, PetscInt numCorners, const int cells[]) 985 { 986 PetscInt *cone, c, p; 987 PetscErrorCode ierr; 988 989 PetscFunctionBegin; 990 ierr = DMPlexSetChart(dm, 0, numCells+numVertices);CHKERRQ(ierr); 991 for (c = 0; c < numCells; ++c) { 992 ierr = DMPlexSetConeSize(dm, c, numCorners);CHKERRQ(ierr); 993 } 994 ierr = DMSetUp(dm);CHKERRQ(ierr); 995 ierr = DMGetWorkArray(dm, numCorners, PETSC_INT, &cone);CHKERRQ(ierr); 996 for (c = 0; c < numCells; ++c) { 997 for (p = 0; p < numCorners; ++p) { 998 cone[p] = cells[c*numCorners+p]+numCells; 999 } 1000 ierr = DMPlexSetCone(dm, c, cone);CHKERRQ(ierr); 1001 } 1002 ierr = DMRestoreWorkArray(dm, numCorners, PETSC_INT, &cone);CHKERRQ(ierr); 1003 ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr); 1004 ierr = DMPlexStratify(dm);CHKERRQ(ierr); 1005 PetscFunctionReturn(0); 1006 } 1007 1008 #undef __FUNCT__ 1009 #define __FUNCT__ "DMPlexBuildCoordinates_Private" 1010 /* 1011 This takes as input the coordinates for each vertex 1012 */ 1013 PetscErrorCode DMPlexBuildCoordinates_Private(DM dm, PetscInt spaceDim, PetscInt numCells, PetscInt numVertices, const double vertexCoords[]) 1014 { 1015 PetscSection coordSection; 1016 Vec coordinates; 1017 PetscScalar *coords; 1018 PetscInt coordSize, v, d; 1019 PetscErrorCode ierr; 1020 1021 PetscFunctionBegin; 1022 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 1023 ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr); 1024 ierr = PetscSectionSetFieldComponents(coordSection, 0, spaceDim);CHKERRQ(ierr); 1025 ierr = PetscSectionSetChart(coordSection, numCells, numCells + numVertices);CHKERRQ(ierr); 1026 for (v = numCells; v < numCells+numVertices; ++v) { 1027 ierr = PetscSectionSetDof(coordSection, v, spaceDim);CHKERRQ(ierr); 1028 ierr = PetscSectionSetFieldDof(coordSection, v, 0, spaceDim);CHKERRQ(ierr); 1029 } 1030 ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 1031 ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 1032 ierr = VecCreate(PetscObjectComm((PetscObject)dm), &coordinates);CHKERRQ(ierr); 1033 ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr); 1034 ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 1035 ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr); 1036 ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 1037 for (v = 0; v < numVertices; ++v) { 1038 for (d = 0; d < spaceDim; ++d) { 1039 coords[v*spaceDim+d] = vertexCoords[v*spaceDim+d]; 1040 } 1041 } 1042 ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 1043 ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr); 1044 ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 1045 PetscFunctionReturn(0); 1046 } 1047 1048 #undef __FUNCT__ 1049 #define __FUNCT__ "DMPlexCreateFromCellList" 1050 /*@C 1051 DMPlexCreateFromCellList - This takes as input common mesh generator output, a list of the vertices for each cell, and produces a DM 1052 1053 Input Parameters: 1054 + comm - The communicator 1055 . dim - The topological dimension of the mesh 1056 . numCells - The number of cells 1057 . numVertices - The number of vertices 1058 . numCorners - The number of vertices for each cell 1059 . interpolate - Flag indicating that intermediate mesh entities (faces, edges) should be created automatically 1060 . cells - An array of numCells*numCorners numbers, the vertices for each cell 1061 . spaceDim - The spatial dimension used for coordinates 1062 - vertexCoords - An array of numVertices*spaceDim numbers, the coordinates of each vertex 1063 1064 Output Parameter: 1065 . dm - The DM 1066 1067 Note: Two triangles sharing a face 1068 $ 1069 $ 2 1070 $ / | \ 1071 $ / | \ 1072 $ / | \ 1073 $ 0 0 | 1 3 1074 $ \ | / 1075 $ \ | / 1076 $ \ | / 1077 $ 1 1078 would have input 1079 $ numCells = 2, numVertices = 4 1080 $ cells = [0 1 2 1 3 2] 1081 $ 1082 which would result in the DMPlex 1083 $ 1084 $ 4 1085 $ / | \ 1086 $ / | \ 1087 $ / | \ 1088 $ 2 0 | 1 5 1089 $ \ | / 1090 $ \ | / 1091 $ \ | / 1092 $ 3 1093 1094 Level: beginner 1095 1096 .seealso: DMPlexCreateFromDAG(), DMPlexCreate() 1097 @*/ 1098 PetscErrorCode DMPlexCreateFromCellList(MPI_Comm comm, PetscInt dim, PetscInt numCells, PetscInt numVertices, PetscInt numCorners, PetscBool interpolate, const int cells[], PetscInt spaceDim, const double vertexCoords[], DM *dm) 1099 { 1100 PetscErrorCode ierr; 1101 1102 PetscFunctionBegin; 1103 ierr = DMCreate(comm, dm);CHKERRQ(ierr); 1104 ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 1105 ierr = DMPlexSetDimension(*dm, dim);CHKERRQ(ierr); 1106 ierr = DMPlexBuildFromCellList_Private(*dm, numCells, numVertices, numCorners, cells);CHKERRQ(ierr); 1107 if (interpolate) { 1108 DM idm; 1109 1110 ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr); 1111 ierr = DMDestroy(dm);CHKERRQ(ierr); 1112 *dm = idm; 1113 } 1114 ierr = DMPlexBuildCoordinates_Private(*dm, spaceDim, numCells, numVertices, vertexCoords);CHKERRQ(ierr); 1115 PetscFunctionReturn(0); 1116 } 1117 1118 #undef __FUNCT__ 1119 #define __FUNCT__ "DMPlexCreateFromDAG" 1120 /*@ 1121 DMPlexCreateFromDAG - This takes as input the adjacency-list representation of the Directed Acyclic Graph (Hasse Diagram) encoding a mesh, and produces a DM 1122 1123 Input Parameters: 1124 + dm - The empty DM object, usually from DMCreate() and DMPlexSetDimension() 1125 . depth - The depth of the DAG 1126 . numPoints - The number of points at each depth 1127 . coneSize - The cone size of each point 1128 . cones - The concatenation of the cone points for each point, the cone list must be oriented correctly for each point 1129 . coneOrientations - The orientation of each cone point 1130 - vertexCoords - An array of numVertices*dim numbers, the coordinates of each vertex 1131 1132 Output Parameter: 1133 . dm - The DM 1134 1135 Note: Two triangles sharing a face would have input 1136 $ depth = 1, numPoints = [4 2], coneSize = [3 3 0 0 0 0] 1137 $ cones = [2 3 4 3 5 4], coneOrientations = [0 0 0 0 0 0] 1138 $ vertexCoords = [-1.0 0.0 0.0 -1.0 0.0 1.0 1.0 0.0] 1139 $ 1140 which would result in the DMPlex 1141 $ 1142 $ 4 1143 $ / | \ 1144 $ / | \ 1145 $ / | \ 1146 $ 2 0 | 1 5 1147 $ \ | / 1148 $ \ | / 1149 $ \ | / 1150 $ 3 1151 $ 1152 $ Notice that all points are numbered consecutively, unlikely DMPlexCreateFromCellList() 1153 1154 Level: advanced 1155 1156 .seealso: DMPlexCreateFromCellList(), DMPlexCreate() 1157 @*/ 1158 PetscErrorCode DMPlexCreateFromDAG(DM dm, PetscInt depth, const PetscInt numPoints[], const PetscInt coneSize[], const PetscInt cones[], const PetscInt coneOrientations[], const PetscScalar vertexCoords[]) 1159 { 1160 Vec coordinates; 1161 PetscSection coordSection; 1162 PetscScalar *coords; 1163 PetscInt coordSize, firstVertex = numPoints[depth], pStart = 0, pEnd = 0, p, v, dim, d, off; 1164 PetscErrorCode ierr; 1165 1166 PetscFunctionBegin; 1167 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 1168 for (d = 0; d <= depth; ++d) pEnd += numPoints[d]; 1169 ierr = DMPlexSetChart(dm, pStart, pEnd);CHKERRQ(ierr); 1170 for (p = pStart; p < pEnd; ++p) { 1171 ierr = DMPlexSetConeSize(dm, p, coneSize[p-pStart]);CHKERRQ(ierr); 1172 } 1173 ierr = DMSetUp(dm);CHKERRQ(ierr); /* Allocate space for cones */ 1174 for (p = pStart, off = 0; p < pEnd; off += coneSize[p-pStart], ++p) { 1175 ierr = DMPlexSetCone(dm, p, &cones[off]);CHKERRQ(ierr); 1176 ierr = DMPlexSetConeOrientation(dm, p, &coneOrientations[off]);CHKERRQ(ierr); 1177 } 1178 ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr); 1179 ierr = DMPlexStratify(dm);CHKERRQ(ierr); 1180 /* Build coordinates */ 1181 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 1182 ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr); 1183 ierr = PetscSectionSetFieldComponents(coordSection, 0, dim);CHKERRQ(ierr); 1184 ierr = PetscSectionSetChart(coordSection, firstVertex, firstVertex+numPoints[0]);CHKERRQ(ierr); 1185 for (v = firstVertex; v < firstVertex+numPoints[0]; ++v) { 1186 ierr = PetscSectionSetDof(coordSection, v, dim);CHKERRQ(ierr); 1187 ierr = PetscSectionSetFieldDof(coordSection, v, 0, dim);CHKERRQ(ierr); 1188 } 1189 ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 1190 ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 1191 ierr = VecCreate(PetscObjectComm((PetscObject)dm), &coordinates);CHKERRQ(ierr); 1192 ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr); 1193 ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 1194 ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr); 1195 ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 1196 for (v = 0; v < numPoints[0]; ++v) { 1197 PetscInt off; 1198 1199 ierr = PetscSectionGetOffset(coordSection, v+firstVertex, &off);CHKERRQ(ierr); 1200 for (d = 0; d < dim; ++d) { 1201 coords[off+d] = vertexCoords[v*dim+d]; 1202 } 1203 } 1204 ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 1205 ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr); 1206 ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 1207 PetscFunctionReturn(0); 1208 } 1209