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 vertices[3], numVertices; 316 PetscInt numFaces = 2*faces[0]*faces[1] + 2*faces[1]*faces[2] + 2*faces[0]*faces[2]; 317 Vec coordinates; 318 PetscSection coordSection; 319 PetscScalar *coords; 320 PetscInt coordSize; 321 PetscMPIInt rank; 322 PetscInt v, vx, vy, vz; 323 PetscInt voffset, iface=0, cone[4]; 324 PetscErrorCode ierr; 325 326 PetscFunctionBegin; 327 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"); 328 ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);CHKERRQ(ierr); 329 vertices[0] = faces[0]+1; vertices[1] = faces[1]+1; vertices[2] = faces[2]+1; 330 numVertices = vertices[0]*vertices[1]*vertices[2]; 331 if (!rank) { 332 PetscInt f; 333 334 ierr = DMPlexSetChart(dm, 0, numFaces+numVertices);CHKERRQ(ierr); 335 for (f = 0; f < numFaces; ++f) { 336 ierr = DMPlexSetConeSize(dm, f, 4);CHKERRQ(ierr); 337 } 338 ierr = DMSetUp(dm);CHKERRQ(ierr); /* Allocate space for cones */ 339 for (v = 0; v < numFaces+numVertices; ++v) { 340 ierr = DMPlexSetLabelValue(dm, "marker", v, 1);CHKERRQ(ierr); 341 } 342 343 /* Side 0 (Top) */ 344 for (vy = 0; vy < faces[1]; vy++) { 345 for (vx = 0; vx < faces[0]; vx++) { 346 voffset = numFaces + vertices[0]*vertices[1]*(vertices[2]-1) + vy*vertices[0] + vx; 347 cone[0] = voffset; cone[1] = voffset+1; cone[2] = voffset+vertices[0]+1; cone[3] = voffset+vertices[0]; 348 ierr = DMPlexSetCone(dm, iface, cone);CHKERRQ(ierr); 349 iface++; 350 } 351 } 352 353 /* Side 1 (Bottom) */ 354 for (vy = 0; vy < faces[1]; vy++) { 355 for (vx = 0; vx < faces[0]; vx++) { 356 voffset = numFaces + vy*(faces[0]+1) + vx; 357 cone[0] = voffset+1; cone[1] = voffset; cone[2] = voffset+vertices[0]; cone[3] = voffset+vertices[0]+1; 358 ierr = DMPlexSetCone(dm, iface, cone);CHKERRQ(ierr); 359 iface++; 360 } 361 } 362 363 /* Side 2 (Front) */ 364 for (vz = 0; vz < faces[2]; vz++) { 365 for (vx = 0; vx < faces[0]; vx++) { 366 voffset = numFaces + vz*vertices[0]*vertices[1] + vx; 367 cone[0] = voffset; cone[1] = voffset+1; cone[2] = voffset+vertices[0]*vertices[1]+1; cone[3] = voffset+vertices[0]*vertices[1]; 368 ierr = DMPlexSetCone(dm, iface, cone);CHKERRQ(ierr); 369 iface++; 370 } 371 } 372 373 /* Side 3 (Back) */ 374 for (vz = 0; vz < faces[2]; vz++) { 375 for (vx = 0; vx < faces[0]; vx++) { 376 voffset = numFaces + vz*vertices[0]*vertices[1] + vertices[0]*(vertices[1]-1) + vx; 377 cone[0] = voffset+vertices[0]*vertices[1]; cone[1] = voffset+vertices[0]*vertices[1]+1; 378 cone[2] = voffset+1; cone[3] = voffset; 379 ierr = DMPlexSetCone(dm, iface, cone);CHKERRQ(ierr); 380 iface++; 381 } 382 } 383 384 /* Side 4 (Left) */ 385 for (vz = 0; vz < faces[2]; vz++) { 386 for (vy = 0; vy < faces[1]; vy++) { 387 voffset = numFaces + vz*vertices[0]*vertices[1] + vy*vertices[0]; 388 cone[0] = voffset; cone[1] = voffset+vertices[0]*vertices[1]; 389 cone[2] = voffset+vertices[0]*vertices[1]+vertices[0]; cone[3] = voffset+vertices[0]; 390 ierr = DMPlexSetCone(dm, iface, cone);CHKERRQ(ierr); 391 iface++; 392 } 393 } 394 395 /* Side 5 (Right) */ 396 for (vz = 0; vz < faces[2]; vz++) { 397 for (vy = 0; vy < faces[1]; vy++) { 398 voffset = numFaces + vz*vertices[0]*vertices[1] + vy*vertices[0] + vx; 399 cone[0] = voffset+vertices[0]*vertices[1]; cone[1] = voffset; 400 cone[2] = voffset+vertices[0]; cone[3] = voffset+vertices[0]*vertices[1]+vertices[0]; 401 ierr = DMPlexSetCone(dm, iface, cone);CHKERRQ(ierr); 402 iface++; 403 } 404 } 405 } 406 ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr); 407 ierr = DMPlexStratify(dm);CHKERRQ(ierr); 408 /* Build coordinates */ 409 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 410 ierr = PetscSectionSetChart(coordSection, numFaces, numFaces + numVertices);CHKERRQ(ierr); 411 for (v = numFaces; v < numFaces+numVertices; ++v) { 412 ierr = PetscSectionSetDof(coordSection, v, 3);CHKERRQ(ierr); 413 } 414 ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 415 ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 416 ierr = VecCreate(PetscObjectComm((PetscObject)dm), &coordinates);CHKERRQ(ierr); 417 ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr); 418 ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 419 ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr); 420 ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 421 for (vz = 0; vz <= faces[2]; ++vz) { 422 for (vy = 0; vy <= faces[1]; ++vy) { 423 for (vx = 0; vx <= faces[0]; ++vx) { 424 coords[((vz*(faces[1]+1)+vy)*(faces[0]+1)+vx)*3+0] = lower[0] + ((upper[0] - lower[0])/faces[0])*vx; 425 coords[((vz*(faces[1]+1)+vy)*(faces[0]+1)+vx)*3+1] = lower[1] + ((upper[1] - lower[1])/faces[1])*vy; 426 coords[((vz*(faces[1]+1)+vy)*(faces[0]+1)+vx)*3+2] = lower[2] + ((upper[2] - lower[2])/faces[2])*vz; 427 } 428 } 429 } 430 ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 431 ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr); 432 ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 433 PetscFunctionReturn(0); 434 } 435 436 #undef __FUNCT__ 437 #define __FUNCT__ "DMPlexCreateSquareMesh" 438 /*@ 439 DMPlexCreateSquareMesh - Creates a 2D mesh for a square lattice. 440 441 Collective on MPI_Comm 442 443 Input Parameters: 444 + comm - The communicator for the DM object 445 . lower - The lower left corner coordinates 446 . upper - The upper right corner coordinates 447 . edges - The number of cells in each direction 448 . bdX - The boundary type for the X direction 449 - bdY - The boundary type for the Y direction 450 451 Output Parameter: 452 . dm - The DM object 453 454 Note: Here is the numbering returned for 2 cells in each direction: 455 $ 22--8-23--9--24 456 $ | | | 457 $ 13 2 14 3 15 458 $ | | | 459 $ 19--6-20--7--21 460 $ | | | 461 $ 10 0 11 1 12 462 $ | | | 463 $ 16--4-17--5--18 464 465 Level: beginner 466 467 .keywords: DM, create 468 .seealso: DMPlexCreateBoxMesh(), DMPlexCreateSquareBoundary(), DMPlexCreateCubeBoundary(), DMSetType(), DMCreate() 469 @*/ 470 PetscErrorCode DMPlexCreateSquareMesh(DM dm, const PetscReal lower[], const PetscReal upper[], const PetscInt edges[], DMBoundaryType bdX, DMBoundaryType bdY) 471 { 472 PetscInt markerTop = 1; 473 PetscInt markerBottom = 1; 474 PetscInt markerRight = 1; 475 PetscInt markerLeft = 1; 476 PetscBool markerSeparate = PETSC_FALSE; 477 PetscMPIInt rank; 478 PetscErrorCode ierr; 479 480 PetscFunctionBegin; 481 ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);CHKERRQ(ierr); 482 ierr = PetscOptionsGetBool(((PetscObject) dm)->prefix, "-dm_plex_separate_marker", &markerSeparate, NULL);CHKERRQ(ierr); 483 if (markerSeparate) { 484 markerTop = 3; 485 markerBottom = 1; 486 markerRight = 2; 487 markerLeft = 4; 488 } 489 { 490 const PetscInt numXEdges = !rank ? edges[0] : 0; 491 const PetscInt numYEdges = !rank ? edges[1] : 0; 492 const PetscInt numXVertices = !rank ? (bdX == DM_BOUNDARY_PERIODIC || bdX == DM_BOUNDARY_TWIST ? edges[0] : edges[0]+1) : 0; 493 const PetscInt numYVertices = !rank ? (bdY == DM_BOUNDARY_PERIODIC || bdY == DM_BOUNDARY_TWIST ? edges[1] : edges[1]+1) : 0; 494 const PetscInt numTotXEdges = numXEdges*numYVertices; 495 const PetscInt numTotYEdges = numYEdges*numXVertices; 496 const PetscInt numVertices = numXVertices*numYVertices; 497 const PetscInt numEdges = numTotXEdges + numTotYEdges; 498 const PetscInt numFaces = numXEdges*numYEdges; 499 const PetscInt firstVertex = numFaces; 500 const PetscInt firstXEdge = numFaces + numVertices; 501 const PetscInt firstYEdge = numFaces + numVertices + numTotXEdges; 502 Vec coordinates; 503 PetscSection coordSection; 504 PetscScalar *coords; 505 PetscInt coordSize; 506 PetscInt v, vx, vy; 507 PetscInt f, fx, fy, e, ex, ey; 508 509 ierr = DMPlexSetChart(dm, 0, numFaces+numEdges+numVertices);CHKERRQ(ierr); 510 for (f = 0; f < numFaces; ++f) { 511 ierr = DMPlexSetConeSize(dm, f, 4);CHKERRQ(ierr); 512 } 513 for (e = firstXEdge; e < firstXEdge+numEdges; ++e) { 514 ierr = DMPlexSetConeSize(dm, e, 2);CHKERRQ(ierr); 515 } 516 ierr = DMSetUp(dm);CHKERRQ(ierr); /* Allocate space for cones */ 517 /* Build faces */ 518 for (fy = 0; fy < numYEdges; fy++) { 519 for (fx = 0; fx < numXEdges; fx++) { 520 const PetscInt face = fy*numXEdges + fx; 521 const PetscInt edgeL = firstYEdge + fx *numYEdges + fy; 522 const PetscInt edgeR = firstYEdge + ((fx+1)%numXVertices)*numYEdges + fy; 523 const PetscInt edgeB = firstXEdge + fy *numXEdges + fx; 524 const PetscInt edgeT = firstXEdge + ((fy+1)%numYVertices)*numXEdges + fx; 525 const PetscInt ornt[4] = {0, 0, -2, -2}; 526 PetscInt cone[4]; 527 528 cone[0] = edgeB; cone[1] = edgeR; cone[2] = edgeT; cone[3] = edgeL; 529 ierr = DMPlexSetCone(dm, face, cone);CHKERRQ(ierr); 530 ierr = DMPlexSetConeOrientation(dm, face, ornt);CHKERRQ(ierr); 531 } 532 } 533 /* Build Y edges*/ 534 for (vx = 0; vx < numXVertices; vx++) { 535 for (ey = 0; ey < numYEdges; ey++) { 536 const PetscInt nextv = (bdY == DM_BOUNDARY_TWIST && ey == numYEdges-1) ? (numXVertices-vx-1) : ((ey+1)%numYVertices)*numXVertices + vx; 537 const PetscInt edge = firstYEdge + vx*numYEdges + ey; 538 const PetscInt vertexB = firstVertex + ey*numXVertices + vx; 539 const PetscInt vertexT = firstVertex + nextv; 540 PetscInt cone[2]; 541 542 cone[0] = vertexB; cone[1] = vertexT; 543 ierr = DMPlexSetCone(dm, edge, cone);CHKERRQ(ierr); 544 if ((bdX != DM_BOUNDARY_PERIODIC) && (bdX != DM_BOUNDARY_TWIST)) { 545 if (vx == numXVertices-1) { 546 ierr = DMPlexSetLabelValue(dm, "marker", edge, markerRight);CHKERRQ(ierr); 547 ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerRight);CHKERRQ(ierr); 548 if (ey == numYEdges-1) { 549 ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerRight);CHKERRQ(ierr); 550 } 551 } else if (vx == 0) { 552 ierr = DMPlexSetLabelValue(dm, "marker", edge, markerLeft);CHKERRQ(ierr); 553 ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerLeft);CHKERRQ(ierr); 554 if (ey == numYEdges-1) { 555 ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerLeft);CHKERRQ(ierr); 556 } 557 } 558 } 559 } 560 } 561 /* Build X edges*/ 562 for (vy = 0; vy < numYVertices; vy++) { 563 for (ex = 0; ex < numXEdges; ex++) { 564 const PetscInt nextv = (bdX == DM_BOUNDARY_TWIST && ex == numXEdges-1) ? (numYVertices-vy-1)*numXVertices : vy*numXVertices + (ex+1)%numXVertices; 565 const PetscInt edge = firstXEdge + vy*numXEdges + ex; 566 const PetscInt vertexL = firstVertex + vy*numXVertices + ex; 567 const PetscInt vertexR = firstVertex + nextv; 568 PetscInt cone[2]; 569 570 cone[0] = vertexL; cone[1] = vertexR; 571 ierr = DMPlexSetCone(dm, edge, cone);CHKERRQ(ierr); 572 if ((bdY != DM_BOUNDARY_PERIODIC) && (bdY != DM_BOUNDARY_TWIST)) { 573 if (vy == numYVertices-1) { 574 ierr = DMPlexSetLabelValue(dm, "marker", edge, markerTop);CHKERRQ(ierr); 575 ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerTop);CHKERRQ(ierr); 576 if (ex == numXEdges-1) { 577 ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerTop);CHKERRQ(ierr); 578 } 579 } else if (vy == 0) { 580 ierr = DMPlexSetLabelValue(dm, "marker", edge, markerBottom);CHKERRQ(ierr); 581 ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerBottom);CHKERRQ(ierr); 582 if (ex == numXEdges-1) { 583 ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerBottom);CHKERRQ(ierr); 584 } 585 } 586 } 587 } 588 } 589 ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr); 590 ierr = DMPlexStratify(dm);CHKERRQ(ierr); 591 /* Build coordinates */ 592 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 593 ierr = PetscSectionSetChart(coordSection, firstVertex, firstVertex+numVertices);CHKERRQ(ierr); 594 for (v = firstVertex; v < firstVertex+numVertices; ++v) { 595 ierr = PetscSectionSetDof(coordSection, v, 2);CHKERRQ(ierr); 596 } 597 ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 598 ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 599 ierr = VecCreate(PetscObjectComm((PetscObject)dm), &coordinates);CHKERRQ(ierr); 600 ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 601 ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr); 602 ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 603 for (vy = 0; vy < numYVertices; ++vy) { 604 for (vx = 0; vx < numXVertices; ++vx) { 605 coords[(vy*numXVertices+vx)*2+0] = lower[0] + ((upper[0] - lower[0])/numXEdges)*vx; 606 coords[(vy*numXVertices+vx)*2+1] = lower[1] + ((upper[1] - lower[1])/numYEdges)*vy; 607 } 608 } 609 ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 610 ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr); 611 ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 612 } 613 PetscFunctionReturn(0); 614 } 615 616 #undef __FUNCT__ 617 #define __FUNCT__ "DMPlexCreateBoxMesh" 618 /*@ 619 DMPlexCreateBoxMesh - Creates a mesh on the tensor product of unit intervals (box) using simplices. 620 621 Collective on MPI_Comm 622 623 Input Parameters: 624 + comm - The communicator for the DM object 625 . dim - The spatial dimension 626 - interpolate - Flag to create intermediate mesh pieces (edges, faces) 627 628 Output Parameter: 629 . dm - The DM object 630 631 Level: beginner 632 633 .keywords: DM, create 634 .seealso: DMPlexCreateHexBoxMesh(), DMSetType(), DMCreate() 635 @*/ 636 PetscErrorCode DMPlexCreateBoxMesh(MPI_Comm comm, PetscInt dim, PetscBool interpolate, DM *dm) 637 { 638 DM boundary; 639 PetscErrorCode ierr; 640 641 PetscFunctionBegin; 642 PetscValidPointer(dm, 4); 643 ierr = DMCreate(comm, &boundary);CHKERRQ(ierr); 644 PetscValidLogicalCollectiveInt(boundary,dim,2); 645 ierr = DMSetType(boundary, DMPLEX);CHKERRQ(ierr); 646 ierr = DMPlexSetDimension(boundary, dim-1);CHKERRQ(ierr); 647 switch (dim) { 648 case 2: 649 { 650 PetscReal lower[2] = {0.0, 0.0}; 651 PetscReal upper[2] = {1.0, 1.0}; 652 PetscInt edges[2] = {2, 2}; 653 654 ierr = DMPlexCreateSquareBoundary(boundary, lower, upper, edges);CHKERRQ(ierr); 655 break; 656 } 657 case 3: 658 { 659 PetscReal lower[3] = {0.0, 0.0, 0.0}; 660 PetscReal upper[3] = {1.0, 1.0, 1.0}; 661 PetscInt faces[3] = {1, 1, 1}; 662 663 ierr = DMPlexCreateCubeBoundary(boundary, lower, upper, faces);CHKERRQ(ierr); 664 break; 665 } 666 default: 667 SETERRQ1(comm, PETSC_ERR_SUP, "Dimension not supported: %d", dim); 668 } 669 ierr = DMPlexGenerate(boundary, NULL, interpolate, dm);CHKERRQ(ierr); 670 ierr = DMDestroy(&boundary);CHKERRQ(ierr); 671 PetscFunctionReturn(0); 672 } 673 674 #undef __FUNCT__ 675 #define __FUNCT__ "DMPlexCreateHexBoxMesh" 676 /*@ 677 DMPlexCreateHexBoxMesh - Creates a mesh on the tensor product of unit intervals (box) using hexahedra. 678 679 Collective on MPI_Comm 680 681 Input Parameters: 682 + comm - The communicator for the DM object 683 . dim - The spatial dimension 684 . periodicX - The boundary type for the X direction 685 . periodicY - The boundary type for the Y direction 686 . periodicZ - The boundary type for the Z direction 687 - cells - The number of cells in each direction 688 689 Output Parameter: 690 . dm - The DM object 691 692 Level: beginner 693 694 .keywords: DM, create 695 .seealso: DMPlexCreateBoxMesh(), DMSetType(), DMCreate() 696 @*/ 697 PetscErrorCode DMPlexCreateHexBoxMesh(MPI_Comm comm, PetscInt dim, const PetscInt cells[], DMBoundaryType periodicX, DMBoundaryType periodicY, DMBoundaryType periodicZ, DM *dm) 698 { 699 PetscErrorCode ierr; 700 701 PetscFunctionBegin; 702 PetscValidPointer(dm, 4); 703 ierr = DMCreate(comm, dm);CHKERRQ(ierr); 704 PetscValidLogicalCollectiveInt(*dm,dim,2); 705 ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 706 ierr = DMPlexSetDimension(*dm, dim);CHKERRQ(ierr); 707 switch (dim) { 708 case 2: 709 { 710 PetscReal lower[2] = {0.0, 0.0}; 711 PetscReal upper[2] = {1.0, 1.0}; 712 713 ierr = DMPlexCreateSquareMesh(*dm, lower, upper, cells, periodicX, periodicY);CHKERRQ(ierr); 714 break; 715 } 716 #if 0 717 case 3: 718 { 719 PetscReal lower[3] = {0.0, 0.0, 0.0}; 720 PetscReal upper[3] = {1.0, 1.0, 1.0}; 721 722 ierr = DMPlexCreateCubeMesh(boundary, lower, upper, cells, periodicX, periodicY, periodicZ);CHKERRQ(ierr); 723 break; 724 } 725 #endif 726 default: 727 SETERRQ1(comm, PETSC_ERR_SUP, "Dimension not supported: %d", dim); 728 } 729 PetscFunctionReturn(0); 730 } 731 732 /* External function declarations here */ 733 extern PetscErrorCode DMCreateInterpolation_Plex(DM dmCoarse, DM dmFine, Mat *interpolation, Vec *scaling); 734 extern PetscErrorCode DMCreateInjection_Plex(DM dmCoarse, DM dmFine, VecScatter *ctx); 735 extern PetscErrorCode DMCreateDefaultSection_Plex(DM dm); 736 extern PetscErrorCode DMCreateMatrix_Plex(DM dm, Mat *J); 737 extern PetscErrorCode DMCreateCoordinateDM_Plex(DM dm, DM *cdm); 738 extern PetscErrorCode DMRefine_Plex(DM dm, MPI_Comm comm, DM *dmRefined); 739 extern PetscErrorCode DMCoarsen_Plex(DM dm, MPI_Comm comm, DM *dmCoarsened); 740 extern PetscErrorCode DMRefineHierarchy_Plex(DM dm, PetscInt nlevels, DM dmRefined[]); 741 extern PetscErrorCode DMClone_Plex(DM dm, DM *newdm); 742 extern PetscErrorCode DMSetUp_Plex(DM dm); 743 extern PetscErrorCode DMDestroy_Plex(DM dm); 744 extern PetscErrorCode DMView_Plex(DM dm, PetscViewer viewer); 745 extern PetscErrorCode DMCreateSubDM_Plex(DM dm, PetscInt numFields, PetscInt fields[], IS *is, DM *subdm); 746 extern PetscErrorCode DMLocatePoints_Plex(DM dm, Vec v, IS *cellIS); 747 748 #undef __FUNCT__ 749 #define __FUNCT__ "DMPlexReplace_Static" 750 /* Replace dm with the contents of dmNew 751 - Share the DM_Plex structure 752 - Share the coordinates 753 */ 754 static PetscErrorCode DMPlexReplace_Static(DM dm, DM dmNew) 755 { 756 PetscSection coordSection; 757 Vec coords; 758 PetscErrorCode ierr; 759 760 PetscFunctionBegin; 761 ierr = DMGetCoordinateSection(dmNew, &coordSection);CHKERRQ(ierr); 762 ierr = DMGetCoordinatesLocal(dmNew, &coords);CHKERRQ(ierr); 763 ierr = DMSetCoordinateSection(dm, coordSection);CHKERRQ(ierr); 764 ierr = DMSetCoordinatesLocal(dm, coords);CHKERRQ(ierr); 765 ierr = DMDestroy_Plex(dm);CHKERRQ(ierr); 766 dm->data = dmNew->data; 767 ((DM_Plex *) dmNew->data)->refct++; 768 PetscFunctionReturn(0); 769 } 770 771 #undef __FUNCT__ 772 #define __FUNCT__ "DMPlexSwap_Static" 773 /* Swap dm with the contents of dmNew 774 - Swap the DM_Plex structure 775 - Swap the coordinates 776 */ 777 static PetscErrorCode DMPlexSwap_Static(DM dmA, DM dmB) 778 { 779 DM coordDMA, coordDMB; 780 Vec coordsA, coordsB; 781 void *tmp; 782 PetscErrorCode ierr; 783 784 PetscFunctionBegin; 785 ierr = DMGetCoordinateDM(dmA, &coordDMA);CHKERRQ(ierr); 786 ierr = DMGetCoordinateDM(dmB, &coordDMB);CHKERRQ(ierr); 787 ierr = PetscObjectReference((PetscObject) coordDMA);CHKERRQ(ierr); 788 ierr = DMSetCoordinateDM(dmA, coordDMB);CHKERRQ(ierr); 789 ierr = DMSetCoordinateDM(dmB, coordDMA);CHKERRQ(ierr); 790 ierr = PetscObjectDereference((PetscObject) coordDMA);CHKERRQ(ierr); 791 792 ierr = DMGetCoordinatesLocal(dmA, &coordsA);CHKERRQ(ierr); 793 ierr = DMGetCoordinatesLocal(dmB, &coordsB);CHKERRQ(ierr); 794 ierr = PetscObjectReference((PetscObject) coordsA);CHKERRQ(ierr); 795 ierr = DMSetCoordinatesLocal(dmA, coordsB);CHKERRQ(ierr); 796 ierr = DMSetCoordinatesLocal(dmB, coordsA);CHKERRQ(ierr); 797 ierr = PetscObjectDereference((PetscObject) coordsA);CHKERRQ(ierr); 798 tmp = dmA->data; 799 dmA->data = dmB->data; 800 dmB->data = tmp; 801 PetscFunctionReturn(0); 802 } 803 804 #undef __FUNCT__ 805 #define __FUNCT__ "DMSetFromOptions_NonRefinement_Plex" 806 PetscErrorCode DMSetFromOptions_NonRefinement_Plex(DM dm) 807 { 808 DM_Plex *mesh = (DM_Plex*) dm->data; 809 DMBoundary b; 810 PetscErrorCode ierr; 811 812 PetscFunctionBegin; 813 /* Handle boundary conditions */ 814 ierr = PetscOptionsBegin(PetscObjectComm((PetscObject) dm), NULL, "Boundary condition options", "");CHKERRQ(ierr); 815 for (b = mesh->boundary; b; b = b->next) { 816 char optname[1024]; 817 PetscInt ids[1024], len = 1024, i; 818 PetscBool flg; 819 820 ierr = PetscSNPrintf(optname, sizeof(optname), "-bc_%s", b->name);CHKERRQ(ierr); 821 ierr = PetscMemzero(ids, sizeof(ids));CHKERRQ(ierr); 822 ierr = PetscOptionsIntArray(optname, "List of boundary IDs", "", ids, &len, &flg);CHKERRQ(ierr); 823 if (flg) { 824 DMLabel label; 825 826 ierr = DMPlexGetLabel(dm, b->name, &label);CHKERRQ(ierr); 827 for (i = 0; i < len; ++i) { 828 PetscBool has; 829 830 ierr = DMLabelHasValue(label, ids[i], &has);CHKERRQ(ierr); 831 if (!has) SETERRQ2(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Boundary id %D is not present in the label %s", ids[i], b->name); 832 } 833 b->numids = len; 834 ierr = PetscFree(b->ids);CHKERRQ(ierr); 835 ierr = PetscMalloc1(len, &b->ids);CHKERRQ(ierr); 836 ierr = PetscMemcpy(b->ids, ids, len*sizeof(PetscInt));CHKERRQ(ierr); 837 } 838 } 839 ierr = PetscOptionsEnd();CHKERRQ(ierr); 840 /* Handle viewing */ 841 ierr = PetscOptionsBool("-dm_plex_print_set_values", "Output all set values info", "DMView", PETSC_FALSE, &mesh->printSetValues, NULL);CHKERRQ(ierr); 842 ierr = PetscOptionsInt("-dm_plex_print_fem", "Debug output level all fem computations", "DMView", 0, &mesh->printFEM, NULL);CHKERRQ(ierr); 843 ierr = PetscOptionsReal("-dm_plex_print_tol", "Tolerance for FEM output", "DMView", PETSC_FALSE, &mesh->printTol, NULL);CHKERRQ(ierr); 844 PetscFunctionReturn(0); 845 } 846 847 #undef __FUNCT__ 848 #define __FUNCT__ "DMSetFromOptions_Plex" 849 PetscErrorCode DMSetFromOptions_Plex(DM dm) 850 { 851 PetscInt refine = 0, r; 852 PetscBool isHierarchy; 853 PetscErrorCode ierr; 854 855 PetscFunctionBegin; 856 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 857 ierr = PetscOptionsHead("DMPlex Options");CHKERRQ(ierr); 858 /* Handle DMPlex refinement */ 859 ierr = PetscOptionsInt("-dm_refine", "The number of uniform refinements", "DMCreate", refine, &refine, NULL);CHKERRQ(ierr); 860 ierr = PetscOptionsInt("-dm_refine_hierarchy", "The number of uniform refinements", "DMCreate", refine, &refine, &isHierarchy);CHKERRQ(ierr); 861 ierr = DMPlexSetRefinementUniform(dm, refine ? PETSC_TRUE : PETSC_FALSE);CHKERRQ(ierr); 862 if (refine && isHierarchy) { 863 DM *dms; 864 865 ierr = PetscMalloc1(refine,&dms);CHKERRQ(ierr); 866 ierr = DMRefineHierarchy(dm, refine, dms);CHKERRQ(ierr); 867 /* Total hack since we do not pass in a pointer */ 868 ierr = DMPlexSwap_Static(dm, dms[refine-1]);CHKERRQ(ierr); 869 if (refine == 1) { 870 ierr = DMPlexSetCoarseDM(dm, dms[0]);CHKERRQ(ierr); 871 } else { 872 ierr = DMPlexSetCoarseDM(dm, dms[refine-2]);CHKERRQ(ierr); 873 ierr = DMPlexSetCoarseDM(dms[0], dms[refine-1]);CHKERRQ(ierr); 874 } 875 /* Free DMs */ 876 for (r = 0; r < refine; ++r) { 877 ierr = DMSetFromOptions_NonRefinement_Plex(dm);CHKERRQ(ierr); 878 ierr = DMDestroy(&dms[r]);CHKERRQ(ierr); 879 } 880 ierr = PetscFree(dms);CHKERRQ(ierr); 881 } else { 882 for (r = 0; r < refine; ++r) { 883 DM refinedMesh; 884 885 ierr = DMSetFromOptions_NonRefinement_Plex(dm);CHKERRQ(ierr); 886 ierr = DMRefine(dm, PetscObjectComm((PetscObject) dm), &refinedMesh);CHKERRQ(ierr); 887 /* Total hack since we do not pass in a pointer */ 888 ierr = DMPlexReplace_Static(dm, refinedMesh);CHKERRQ(ierr); 889 ierr = DMDestroy(&refinedMesh);CHKERRQ(ierr); 890 } 891 } 892 ierr = DMSetFromOptions_NonRefinement_Plex(dm);CHKERRQ(ierr); 893 ierr = PetscOptionsTail();CHKERRQ(ierr); 894 PetscFunctionReturn(0); 895 } 896 897 #undef __FUNCT__ 898 #define __FUNCT__ "DMCreateGlobalVector_Plex" 899 static PetscErrorCode DMCreateGlobalVector_Plex(DM dm,Vec *vec) 900 { 901 PetscErrorCode ierr; 902 903 PetscFunctionBegin; 904 ierr = DMCreateGlobalVector_Section_Private(dm,vec);CHKERRQ(ierr); 905 /* ierr = VecSetOperation(*vec, VECOP_DUPLICATE, (void(*)(void)) VecDuplicate_MPI_DM);CHKERRQ(ierr); */ 906 ierr = VecSetOperation(*vec, VECOP_VIEW, (void (*)(void))VecView_Plex);CHKERRQ(ierr); 907 PetscFunctionReturn(0); 908 } 909 910 #undef __FUNCT__ 911 #define __FUNCT__ "DMCreateLocalVector_Plex" 912 static PetscErrorCode DMCreateLocalVector_Plex(DM dm,Vec *vec) 913 { 914 PetscErrorCode ierr; 915 916 PetscFunctionBegin; 917 ierr = DMCreateLocalVector_Section_Private(dm,vec);CHKERRQ(ierr); 918 ierr = VecSetOperation(*vec, VECOP_VIEW, (void(*)(void)) VecView_Plex_Local);CHKERRQ(ierr); 919 PetscFunctionReturn(0); 920 } 921 922 #undef __FUNCT__ 923 #define __FUNCT__ "DMInitialize_Plex" 924 PetscErrorCode DMInitialize_Plex(DM dm) 925 { 926 PetscFunctionBegin; 927 dm->ops->view = DMView_Plex; 928 dm->ops->setfromoptions = DMSetFromOptions_Plex; 929 dm->ops->clone = DMClone_Plex; 930 dm->ops->setup = DMSetUp_Plex; 931 dm->ops->createdefaultsection = DMCreateDefaultSection_Plex; 932 dm->ops->createglobalvector = DMCreateGlobalVector_Plex; 933 dm->ops->createlocalvector = DMCreateLocalVector_Plex; 934 dm->ops->getlocaltoglobalmapping = NULL; 935 dm->ops->getlocaltoglobalmappingblock = NULL; 936 dm->ops->createfieldis = NULL; 937 dm->ops->createcoordinatedm = DMCreateCoordinateDM_Plex; 938 dm->ops->getcoloring = NULL; 939 dm->ops->creatematrix = DMCreateMatrix_Plex; 940 dm->ops->createinterpolation = DMCreateInterpolation_Plex; 941 dm->ops->getaggregates = NULL; 942 dm->ops->getinjection = DMCreateInjection_Plex; 943 dm->ops->refine = DMRefine_Plex; 944 dm->ops->coarsen = DMCoarsen_Plex; 945 dm->ops->refinehierarchy = DMRefineHierarchy_Plex; 946 dm->ops->coarsenhierarchy = NULL; 947 dm->ops->globaltolocalbegin = NULL; 948 dm->ops->globaltolocalend = NULL; 949 dm->ops->localtoglobalbegin = NULL; 950 dm->ops->localtoglobalend = NULL; 951 dm->ops->destroy = DMDestroy_Plex; 952 dm->ops->createsubdm = DMCreateSubDM_Plex; 953 dm->ops->locatepoints = DMLocatePoints_Plex; 954 PetscFunctionReturn(0); 955 } 956 957 #undef __FUNCT__ 958 #define __FUNCT__ "DMClone_Plex" 959 PetscErrorCode DMClone_Plex(DM dm, DM *newdm) 960 { 961 DM_Plex *mesh = (DM_Plex *) dm->data; 962 PetscErrorCode ierr; 963 964 PetscFunctionBegin; 965 mesh->refct++; 966 (*newdm)->data = mesh; 967 ierr = PetscObjectChangeTypeName((PetscObject) *newdm, DMPLEX);CHKERRQ(ierr); 968 ierr = DMInitialize_Plex(*newdm);CHKERRQ(ierr); 969 PetscFunctionReturn(0); 970 } 971 972 /*MC 973 DMPLEX = "plex" - A DM object that encapsulates an unstructured mesh, or CW Complex, which can be expressed using a Hasse Diagram. 974 In the local representation, Vecs contain all unknowns in the interior and shared boundary. This is 975 specified by a PetscSection object. Ownership in the global representation is determined by 976 ownership of the underlying DMPlex points. This is specified by another PetscSection object. 977 978 Level: intermediate 979 980 .seealso: DMType, DMPlexCreate(), DMCreate(), DMSetType() 981 M*/ 982 983 #undef __FUNCT__ 984 #define __FUNCT__ "DMCreate_Plex" 985 PETSC_EXTERN PetscErrorCode DMCreate_Plex(DM dm) 986 { 987 DM_Plex *mesh; 988 PetscInt unit, d; 989 PetscErrorCode ierr; 990 991 PetscFunctionBegin; 992 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 993 ierr = PetscNewLog(dm,&mesh);CHKERRQ(ierr); 994 dm->data = mesh; 995 996 mesh->refct = 1; 997 mesh->dim = 0; 998 ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm), &mesh->coneSection);CHKERRQ(ierr); 999 mesh->maxConeSize = 0; 1000 mesh->cones = NULL; 1001 mesh->coneOrientations = NULL; 1002 ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm), &mesh->supportSection);CHKERRQ(ierr); 1003 mesh->maxSupportSize = 0; 1004 mesh->supports = NULL; 1005 mesh->refinementUniform = PETSC_TRUE; 1006 mesh->refinementLimit = -1.0; 1007 1008 mesh->facesTmp = NULL; 1009 1010 mesh->subpointMap = NULL; 1011 1012 for (unit = 0; unit < NUM_PETSC_UNITS; ++unit) mesh->scale[unit] = 1.0; 1013 1014 mesh->labels = NULL; 1015 mesh->depthLabel = NULL; 1016 mesh->globalVertexNumbers = NULL; 1017 mesh->globalCellNumbers = NULL; 1018 for (d = 0; d < 8; ++d) mesh->hybridPointMax[d] = PETSC_DETERMINE; 1019 mesh->vtkCellHeight = 0; 1020 mesh->useCone = PETSC_FALSE; 1021 mesh->useClosure = PETSC_TRUE; 1022 1023 mesh->printSetValues = PETSC_FALSE; 1024 mesh->printFEM = 0; 1025 mesh->printTol = 1.0e-10; 1026 1027 ierr = DMInitialize_Plex(dm);CHKERRQ(ierr); 1028 PetscFunctionReturn(0); 1029 } 1030 1031 #undef __FUNCT__ 1032 #define __FUNCT__ "DMPlexCreate" 1033 /*@ 1034 DMPlexCreate - Creates a DMPlex object, which encapsulates an unstructured mesh, or CW complex, which can be expressed using a Hasse Diagram. 1035 1036 Collective on MPI_Comm 1037 1038 Input Parameter: 1039 . comm - The communicator for the DMPlex object 1040 1041 Output Parameter: 1042 . mesh - The DMPlex object 1043 1044 Level: beginner 1045 1046 .keywords: DMPlex, create 1047 @*/ 1048 PetscErrorCode DMPlexCreate(MPI_Comm comm, DM *mesh) 1049 { 1050 PetscErrorCode ierr; 1051 1052 PetscFunctionBegin; 1053 PetscValidPointer(mesh,2); 1054 ierr = DMCreate(comm, mesh);CHKERRQ(ierr); 1055 ierr = DMSetType(*mesh, DMPLEX);CHKERRQ(ierr); 1056 PetscFunctionReturn(0); 1057 } 1058 1059 #undef __FUNCT__ 1060 #define __FUNCT__ "DMPlexBuildFromCellList_Private" 1061 /* 1062 This takes as input the common mesh generator output, a list of the vertices for each cell 1063 */ 1064 PetscErrorCode DMPlexBuildFromCellList_Private(DM dm, PetscInt numCells, PetscInt numVertices, PetscInt numCorners, const int cells[]) 1065 { 1066 PetscInt *cone, c, p; 1067 PetscErrorCode ierr; 1068 1069 PetscFunctionBegin; 1070 ierr = DMPlexSetChart(dm, 0, numCells+numVertices);CHKERRQ(ierr); 1071 for (c = 0; c < numCells; ++c) { 1072 ierr = DMPlexSetConeSize(dm, c, numCorners);CHKERRQ(ierr); 1073 } 1074 ierr = DMSetUp(dm);CHKERRQ(ierr); 1075 ierr = DMGetWorkArray(dm, numCorners, PETSC_INT, &cone);CHKERRQ(ierr); 1076 for (c = 0; c < numCells; ++c) { 1077 for (p = 0; p < numCorners; ++p) { 1078 cone[p] = cells[c*numCorners+p]+numCells; 1079 } 1080 ierr = DMPlexSetCone(dm, c, cone);CHKERRQ(ierr); 1081 } 1082 ierr = DMRestoreWorkArray(dm, numCorners, PETSC_INT, &cone);CHKERRQ(ierr); 1083 ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr); 1084 ierr = DMPlexStratify(dm);CHKERRQ(ierr); 1085 PetscFunctionReturn(0); 1086 } 1087 1088 #undef __FUNCT__ 1089 #define __FUNCT__ "DMPlexBuildCoordinates_Private" 1090 /* 1091 This takes as input the coordinates for each vertex 1092 */ 1093 PetscErrorCode DMPlexBuildCoordinates_Private(DM dm, PetscInt spaceDim, PetscInt numCells, PetscInt numVertices, const double vertexCoords[]) 1094 { 1095 PetscSection coordSection; 1096 Vec coordinates; 1097 PetscScalar *coords; 1098 PetscInt coordSize, v, d; 1099 PetscErrorCode ierr; 1100 1101 PetscFunctionBegin; 1102 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 1103 ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr); 1104 ierr = PetscSectionSetFieldComponents(coordSection, 0, spaceDim);CHKERRQ(ierr); 1105 ierr = PetscSectionSetChart(coordSection, numCells, numCells + numVertices);CHKERRQ(ierr); 1106 for (v = numCells; v < numCells+numVertices; ++v) { 1107 ierr = PetscSectionSetDof(coordSection, v, spaceDim);CHKERRQ(ierr); 1108 ierr = PetscSectionSetFieldDof(coordSection, v, 0, spaceDim);CHKERRQ(ierr); 1109 } 1110 ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 1111 ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 1112 ierr = VecCreate(PetscObjectComm((PetscObject)dm), &coordinates);CHKERRQ(ierr); 1113 ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr); 1114 ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 1115 ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr); 1116 ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 1117 for (v = 0; v < numVertices; ++v) { 1118 for (d = 0; d < spaceDim; ++d) { 1119 coords[v*spaceDim+d] = vertexCoords[v*spaceDim+d]; 1120 } 1121 } 1122 ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 1123 ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr); 1124 ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 1125 PetscFunctionReturn(0); 1126 } 1127 1128 #undef __FUNCT__ 1129 #define __FUNCT__ "DMPlexCreateFromCellList" 1130 /*@C 1131 DMPlexCreateFromCellList - This takes as input common mesh generator output, a list of the vertices for each cell, and produces a DM 1132 1133 Input Parameters: 1134 + comm - The communicator 1135 . dim - The topological dimension of the mesh 1136 . numCells - The number of cells 1137 . numVertices - The number of vertices 1138 . numCorners - The number of vertices for each cell 1139 . interpolate - Flag indicating that intermediate mesh entities (faces, edges) should be created automatically 1140 . cells - An array of numCells*numCorners numbers, the vertices for each cell 1141 . spaceDim - The spatial dimension used for coordinates 1142 - vertexCoords - An array of numVertices*spaceDim numbers, the coordinates of each vertex 1143 1144 Output Parameter: 1145 . dm - The DM 1146 1147 Note: Two triangles sharing a face 1148 $ 1149 $ 2 1150 $ / | \ 1151 $ / | \ 1152 $ / | \ 1153 $ 0 0 | 1 3 1154 $ \ | / 1155 $ \ | / 1156 $ \ | / 1157 $ 1 1158 would have input 1159 $ numCells = 2, numVertices = 4 1160 $ cells = [0 1 2 1 3 2] 1161 $ 1162 which would result in the DMPlex 1163 $ 1164 $ 4 1165 $ / | \ 1166 $ / | \ 1167 $ / | \ 1168 $ 2 0 | 1 5 1169 $ \ | / 1170 $ \ | / 1171 $ \ | / 1172 $ 3 1173 1174 Level: beginner 1175 1176 .seealso: DMPlexCreateFromDAG(), DMPlexCreate() 1177 @*/ 1178 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) 1179 { 1180 PetscErrorCode ierr; 1181 1182 PetscFunctionBegin; 1183 ierr = DMCreate(comm, dm);CHKERRQ(ierr); 1184 ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 1185 ierr = DMPlexSetDimension(*dm, dim);CHKERRQ(ierr); 1186 ierr = DMPlexBuildFromCellList_Private(*dm, numCells, numVertices, numCorners, cells);CHKERRQ(ierr); 1187 if (interpolate) { 1188 DM idm; 1189 1190 ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr); 1191 ierr = DMDestroy(dm);CHKERRQ(ierr); 1192 *dm = idm; 1193 } 1194 ierr = DMPlexBuildCoordinates_Private(*dm, spaceDim, numCells, numVertices, vertexCoords);CHKERRQ(ierr); 1195 PetscFunctionReturn(0); 1196 } 1197 1198 #undef __FUNCT__ 1199 #define __FUNCT__ "DMPlexCreateFromDAG" 1200 /*@ 1201 DMPlexCreateFromDAG - This takes as input the adjacency-list representation of the Directed Acyclic Graph (Hasse Diagram) encoding a mesh, and produces a DM 1202 1203 Input Parameters: 1204 + dm - The empty DM object, usually from DMCreate() and DMPlexSetDimension() 1205 . depth - The depth of the DAG 1206 . numPoints - The number of points at each depth 1207 . coneSize - The cone size of each point 1208 . cones - The concatenation of the cone points for each point, the cone list must be oriented correctly for each point 1209 . coneOrientations - The orientation of each cone point 1210 - vertexCoords - An array of numVertices*dim numbers, the coordinates of each vertex 1211 1212 Output Parameter: 1213 . dm - The DM 1214 1215 Note: Two triangles sharing a face would have input 1216 $ depth = 1, numPoints = [4 2], coneSize = [3 3 0 0 0 0] 1217 $ cones = [2 3 4 3 5 4], coneOrientations = [0 0 0 0 0 0] 1218 $ vertexCoords = [-1.0 0.0 0.0 -1.0 0.0 1.0 1.0 0.0] 1219 $ 1220 which would result in the DMPlex 1221 $ 1222 $ 4 1223 $ / | \ 1224 $ / | \ 1225 $ / | \ 1226 $ 2 0 | 1 5 1227 $ \ | / 1228 $ \ | / 1229 $ \ | / 1230 $ 3 1231 $ 1232 $ Notice that all points are numbered consecutively, unlikely DMPlexCreateFromCellList() 1233 1234 Level: advanced 1235 1236 .seealso: DMPlexCreateFromCellList(), DMPlexCreate() 1237 @*/ 1238 PetscErrorCode DMPlexCreateFromDAG(DM dm, PetscInt depth, const PetscInt numPoints[], const PetscInt coneSize[], const PetscInt cones[], const PetscInt coneOrientations[], const PetscScalar vertexCoords[]) 1239 { 1240 Vec coordinates; 1241 PetscSection coordSection; 1242 PetscScalar *coords; 1243 PetscInt coordSize, firstVertex = numPoints[depth], pStart = 0, pEnd = 0, p, v, dim, d, off; 1244 PetscErrorCode ierr; 1245 1246 PetscFunctionBegin; 1247 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 1248 for (d = 0; d <= depth; ++d) pEnd += numPoints[d]; 1249 ierr = DMPlexSetChart(dm, pStart, pEnd);CHKERRQ(ierr); 1250 for (p = pStart; p < pEnd; ++p) { 1251 ierr = DMPlexSetConeSize(dm, p, coneSize[p-pStart]);CHKERRQ(ierr); 1252 } 1253 ierr = DMSetUp(dm);CHKERRQ(ierr); /* Allocate space for cones */ 1254 for (p = pStart, off = 0; p < pEnd; off += coneSize[p-pStart], ++p) { 1255 ierr = DMPlexSetCone(dm, p, &cones[off]);CHKERRQ(ierr); 1256 ierr = DMPlexSetConeOrientation(dm, p, &coneOrientations[off]);CHKERRQ(ierr); 1257 } 1258 ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr); 1259 ierr = DMPlexStratify(dm);CHKERRQ(ierr); 1260 /* Build coordinates */ 1261 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 1262 ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr); 1263 ierr = PetscSectionSetFieldComponents(coordSection, 0, dim);CHKERRQ(ierr); 1264 ierr = PetscSectionSetChart(coordSection, firstVertex, firstVertex+numPoints[0]);CHKERRQ(ierr); 1265 for (v = firstVertex; v < firstVertex+numPoints[0]; ++v) { 1266 ierr = PetscSectionSetDof(coordSection, v, dim);CHKERRQ(ierr); 1267 ierr = PetscSectionSetFieldDof(coordSection, v, 0, dim);CHKERRQ(ierr); 1268 } 1269 ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 1270 ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 1271 ierr = VecCreate(PetscObjectComm((PetscObject)dm), &coordinates);CHKERRQ(ierr); 1272 ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr); 1273 ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 1274 ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr); 1275 ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 1276 for (v = 0; v < numPoints[0]; ++v) { 1277 PetscInt off; 1278 1279 ierr = PetscSectionGetOffset(coordSection, v+firstVertex, &off);CHKERRQ(ierr); 1280 for (d = 0; d < dim; ++d) { 1281 coords[off+d] = vertexCoords[v*dim+d]; 1282 } 1283 } 1284 ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 1285 ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr); 1286 ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 1287 PetscFunctionReturn(0); 1288 } 1289