1 #define PETSCDM_DLL 2 #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 3 #include <petscdmda.h> 4 #include <petscsf.h> 5 6 /*@ 7 DMPlexCreateDoublet - Creates a mesh of two cells of the specified type, optionally with later refinement. 8 9 Collective on MPI_Comm 10 11 Input Parameters: 12 + comm - The communicator for the DM object 13 . dim - The spatial dimension 14 . simplex - Flag for simplicial cells, otherwise they are tensor product cells 15 . interpolate - Flag to create intermediate mesh pieces (edges, faces) 16 . refinementUniform - Flag for uniform parallel refinement 17 - refinementLimit - A nonzero number indicates the largest admissible volume for a refined cell 18 19 Output Parameter: 20 . dm - The DM object 21 22 Level: beginner 23 24 .keywords: DM, create 25 .seealso: DMSetType(), DMCreate() 26 @*/ 27 PetscErrorCode DMPlexCreateDoublet(MPI_Comm comm, PetscInt dim, PetscBool simplex, PetscBool interpolate, PetscBool refinementUniform, PetscReal refinementLimit, DM *newdm) 28 { 29 DM dm; 30 PetscInt p; 31 PetscMPIInt rank; 32 PetscErrorCode ierr; 33 34 PetscFunctionBegin; 35 ierr = DMCreate(comm, &dm);CHKERRQ(ierr); 36 ierr = DMSetType(dm, DMPLEX);CHKERRQ(ierr); 37 ierr = DMSetDimension(dm, dim);CHKERRQ(ierr); 38 ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 39 switch (dim) { 40 case 2: 41 if (simplex) {ierr = PetscObjectSetName((PetscObject) dm, "triangular");CHKERRQ(ierr);} 42 else {ierr = PetscObjectSetName((PetscObject) dm, "quadrilateral");CHKERRQ(ierr);} 43 break; 44 case 3: 45 if (simplex) {ierr = PetscObjectSetName((PetscObject) dm, "tetrahedral");CHKERRQ(ierr);} 46 else {ierr = PetscObjectSetName((PetscObject) dm, "hexahedral");CHKERRQ(ierr);} 47 break; 48 default: 49 SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Cannot make meshes for dimension %d", dim); 50 } 51 if (rank) { 52 PetscInt numPoints[2] = {0, 0}; 53 ierr = DMPlexCreateFromDAG(dm, 1, numPoints, NULL, NULL, NULL, NULL);CHKERRQ(ierr); 54 } else { 55 switch (dim) { 56 case 2: 57 if (simplex) { 58 PetscInt numPoints[2] = {4, 2}; 59 PetscInt coneSize[6] = {3, 3, 0, 0, 0, 0}; 60 PetscInt cones[6] = {2, 3, 4, 5, 4, 3}; 61 PetscInt coneOrientations[6] = {0, 0, 0, 0, 0, 0}; 62 PetscScalar vertexCoords[8] = {-0.5, 0.5, 0.0, 0.0, 0.0, 1.0, 0.5, 0.5}; 63 PetscInt markerPoints[8] = {2, 1, 3, 1, 4, 1, 5, 1}; 64 65 ierr = DMPlexCreateFromDAG(dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr); 66 for (p = 0; p < 4; ++p) {ierr = DMSetLabelValue(dm, "marker", markerPoints[p*2], markerPoints[p*2+1]);CHKERRQ(ierr);} 67 } else { 68 PetscInt numPoints[2] = {6, 2}; 69 PetscInt coneSize[8] = {4, 4, 0, 0, 0, 0, 0, 0}; 70 PetscInt cones[8] = {2, 3, 4, 5, 3, 6, 7, 4}; 71 PetscInt coneOrientations[8] = {0, 0, 0, 0, 0, 0, 0, 0}; 72 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}; 73 74 ierr = DMPlexCreateFromDAG(dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr); 75 } 76 break; 77 case 3: 78 if (simplex) { 79 PetscInt numPoints[2] = {5, 2}; 80 PetscInt coneSize[7] = {4, 4, 0, 0, 0, 0, 0}; 81 PetscInt cones[8] = {4, 3, 5, 2, 5, 3, 4, 6}; 82 PetscInt coneOrientations[8] = {0, 0, 0, 0, 0, 0, 0, 0}; 83 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}; 84 PetscInt markerPoints[10] = {2, 1, 3, 1, 4, 1, 5, 1, 6, 1}; 85 86 ierr = DMPlexCreateFromDAG(dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr); 87 for (p = 0; p < 5; ++p) {ierr = DMSetLabelValue(dm, "marker", markerPoints[p*2], markerPoints[p*2+1]);CHKERRQ(ierr);} 88 } else { 89 PetscInt numPoints[2] = {12, 2}; 90 PetscInt coneSize[14] = {8, 8, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0}; 91 PetscInt cones[16] = {2, 3, 4, 5, 6, 7, 8, 9, 5, 4, 10, 11, 7, 12, 13, 8}; 92 PetscInt coneOrientations[16] = {0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0}; 93 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, 94 -1.0, -0.5, 0.5, 0.0, -0.5, 0.5, 0.0, 0.5, 0.5, -1.0, 0.5, 0.5, 95 1.0, 0.5, -0.5, 1.0, -0.5, -0.5, 1.0, -0.5, 0.5, 1.0, 0.5, 0.5}; 96 97 ierr = DMPlexCreateFromDAG(dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr); 98 } 99 break; 100 default: 101 SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Cannot make meshes for dimension %d", dim); 102 } 103 } 104 *newdm = dm; 105 if (refinementLimit > 0.0) { 106 DM rdm; 107 const char *name; 108 109 ierr = DMPlexSetRefinementUniform(*newdm, PETSC_FALSE);CHKERRQ(ierr); 110 ierr = DMPlexSetRefinementLimit(*newdm, refinementLimit);CHKERRQ(ierr); 111 ierr = DMRefine(*newdm, comm, &rdm);CHKERRQ(ierr); 112 ierr = PetscObjectGetName((PetscObject) *newdm, &name);CHKERRQ(ierr); 113 ierr = PetscObjectSetName((PetscObject) rdm, name);CHKERRQ(ierr); 114 ierr = DMDestroy(newdm);CHKERRQ(ierr); 115 *newdm = rdm; 116 } 117 if (interpolate) { 118 DM idm = NULL; 119 const char *name; 120 121 ierr = DMPlexInterpolate(*newdm, &idm);CHKERRQ(ierr); 122 ierr = PetscObjectGetName((PetscObject) *newdm, &name);CHKERRQ(ierr); 123 ierr = PetscObjectSetName((PetscObject) idm, name);CHKERRQ(ierr); 124 ierr = DMPlexCopyCoordinates(*newdm, idm);CHKERRQ(ierr); 125 ierr = DMCopyLabels(*newdm, idm);CHKERRQ(ierr); 126 ierr = DMDestroy(newdm);CHKERRQ(ierr); 127 *newdm = idm; 128 } 129 { 130 DM refinedMesh = NULL; 131 DM distributedMesh = NULL; 132 133 /* Distribute mesh over processes */ 134 ierr = DMPlexDistribute(*newdm, 0, NULL, &distributedMesh);CHKERRQ(ierr); 135 if (distributedMesh) { 136 ierr = DMDestroy(newdm);CHKERRQ(ierr); 137 *newdm = distributedMesh; 138 } 139 if (refinementUniform) { 140 ierr = DMPlexSetRefinementUniform(*newdm, refinementUniform);CHKERRQ(ierr); 141 ierr = DMRefine(*newdm, comm, &refinedMesh);CHKERRQ(ierr); 142 if (refinedMesh) { 143 ierr = DMDestroy(newdm);CHKERRQ(ierr); 144 *newdm = refinedMesh; 145 } 146 } 147 } 148 PetscFunctionReturn(0); 149 } 150 151 /*@ 152 DMPlexCreateSquareBoundary - Creates a 1D mesh the is the boundary of a square lattice. 153 154 Collective on MPI_Comm 155 156 Input Parameters: 157 + comm - The communicator for the DM object 158 . lower - The lower left corner coordinates 159 . upper - The upper right corner coordinates 160 - edges - The number of cells in each direction 161 162 Output Parameter: 163 . dm - The DM object 164 165 Note: Here is the numbering returned for 2 cells in each direction: 166 $ 18--5-17--4--16 167 $ | | | 168 $ 6 10 3 169 $ | | | 170 $ 19-11-20--9--15 171 $ | | | 172 $ 7 8 2 173 $ | | | 174 $ 12--0-13--1--14 175 176 Level: beginner 177 178 .keywords: DM, create 179 .seealso: DMPlexCreateBoxMesh(), DMPlexCreateCubeBoundary(), DMSetType(), DMCreate() 180 @*/ 181 PetscErrorCode DMPlexCreateSquareBoundary(DM dm, const PetscReal lower[], const PetscReal upper[], const PetscInt edges[]) 182 { 183 const PetscInt numVertices = (edges[0]+1)*(edges[1]+1); 184 const PetscInt numEdges = edges[0]*(edges[1]+1) + (edges[0]+1)*edges[1]; 185 PetscInt markerTop = 1; 186 PetscInt markerBottom = 1; 187 PetscInt markerRight = 1; 188 PetscInt markerLeft = 1; 189 PetscBool markerSeparate = PETSC_FALSE; 190 Vec coordinates; 191 PetscSection coordSection; 192 PetscScalar *coords; 193 PetscInt coordSize; 194 PetscMPIInt rank; 195 PetscInt v, vx, vy; 196 PetscErrorCode ierr; 197 198 PetscFunctionBegin; 199 ierr = PetscOptionsGetBool(((PetscObject) dm)->options,((PetscObject) dm)->prefix, "-dm_plex_separate_marker", &markerSeparate, NULL);CHKERRQ(ierr); 200 if (markerSeparate) { 201 markerTop = 3; 202 markerBottom = 1; 203 markerRight = 2; 204 markerLeft = 4; 205 } 206 ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);CHKERRQ(ierr); 207 if (!rank) { 208 PetscInt e, ex, ey; 209 210 ierr = DMPlexSetChart(dm, 0, numEdges+numVertices);CHKERRQ(ierr); 211 for (e = 0; e < numEdges; ++e) { 212 ierr = DMPlexSetConeSize(dm, e, 2);CHKERRQ(ierr); 213 } 214 ierr = DMSetUp(dm);CHKERRQ(ierr); /* Allocate space for cones */ 215 for (vx = 0; vx <= edges[0]; vx++) { 216 for (ey = 0; ey < edges[1]; ey++) { 217 PetscInt edge = vx*edges[1] + ey + edges[0]*(edges[1]+1); 218 PetscInt vertex = ey*(edges[0]+1) + vx + numEdges; 219 PetscInt cone[2]; 220 221 cone[0] = vertex; cone[1] = vertex+edges[0]+1; 222 ierr = DMPlexSetCone(dm, edge, cone);CHKERRQ(ierr); 223 if (vx == edges[0]) { 224 ierr = DMSetLabelValue(dm, "marker", edge, markerRight);CHKERRQ(ierr); 225 ierr = DMSetLabelValue(dm, "marker", cone[0], markerRight);CHKERRQ(ierr); 226 if (ey == edges[1]-1) { 227 ierr = DMSetLabelValue(dm, "marker", cone[1], markerRight);CHKERRQ(ierr); 228 } 229 } else if (vx == 0) { 230 ierr = DMSetLabelValue(dm, "marker", edge, markerLeft);CHKERRQ(ierr); 231 ierr = DMSetLabelValue(dm, "marker", cone[0], markerLeft);CHKERRQ(ierr); 232 if (ey == edges[1]-1) { 233 ierr = DMSetLabelValue(dm, "marker", cone[1], markerLeft);CHKERRQ(ierr); 234 } 235 } 236 } 237 } 238 for (vy = 0; vy <= edges[1]; vy++) { 239 for (ex = 0; ex < edges[0]; ex++) { 240 PetscInt edge = vy*edges[0] + ex; 241 PetscInt vertex = vy*(edges[0]+1) + ex + numEdges; 242 PetscInt cone[2]; 243 244 cone[0] = vertex; cone[1] = vertex+1; 245 ierr = DMPlexSetCone(dm, edge, cone);CHKERRQ(ierr); 246 if (vy == edges[1]) { 247 ierr = DMSetLabelValue(dm, "marker", edge, markerTop);CHKERRQ(ierr); 248 ierr = DMSetLabelValue(dm, "marker", cone[0], markerTop);CHKERRQ(ierr); 249 if (ex == edges[0]-1) { 250 ierr = DMSetLabelValue(dm, "marker", cone[1], markerTop);CHKERRQ(ierr); 251 } 252 } else if (vy == 0) { 253 ierr = DMSetLabelValue(dm, "marker", edge, markerBottom);CHKERRQ(ierr); 254 ierr = DMSetLabelValue(dm, "marker", cone[0], markerBottom);CHKERRQ(ierr); 255 if (ex == edges[0]-1) { 256 ierr = DMSetLabelValue(dm, "marker", cone[1], markerBottom);CHKERRQ(ierr); 257 } 258 } 259 } 260 } 261 } 262 ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr); 263 ierr = DMPlexStratify(dm);CHKERRQ(ierr); 264 /* Build coordinates */ 265 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 266 ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr); 267 ierr = PetscSectionSetChart(coordSection, numEdges, numEdges + numVertices);CHKERRQ(ierr); 268 ierr = PetscSectionSetFieldComponents(coordSection, 0, 2);CHKERRQ(ierr); 269 for (v = numEdges; v < numEdges+numVertices; ++v) { 270 ierr = PetscSectionSetDof(coordSection, v, 2);CHKERRQ(ierr); 271 ierr = PetscSectionSetFieldDof(coordSection, v, 0, 2);CHKERRQ(ierr); 272 } 273 ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 274 ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 275 ierr = VecCreate(PETSC_COMM_SELF, &coordinates);CHKERRQ(ierr); 276 ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr); 277 ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 278 ierr = VecSetBlockSize(coordinates, 2);CHKERRQ(ierr); 279 ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr); 280 ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 281 for (vy = 0; vy <= edges[1]; ++vy) { 282 for (vx = 0; vx <= edges[0]; ++vx) { 283 coords[(vy*(edges[0]+1)+vx)*2+0] = lower[0] + ((upper[0] - lower[0])/edges[0])*vx; 284 coords[(vy*(edges[0]+1)+vx)*2+1] = lower[1] + ((upper[1] - lower[1])/edges[1])*vy; 285 } 286 } 287 ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 288 ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr); 289 ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 290 PetscFunctionReturn(0); 291 } 292 293 /*@ 294 DMPlexCreateCubeBoundary - Creates a 2D mesh the is the boundary of a cubic lattice. 295 296 Collective on MPI_Comm 297 298 Input Parameters: 299 + comm - The communicator for the DM object 300 . lower - The lower left front corner coordinates 301 . upper - The upper right back corner coordinates 302 - edges - The number of cells in each direction 303 304 Output Parameter: 305 . dm - The DM object 306 307 Level: beginner 308 309 .keywords: DM, create 310 .seealso: DMPlexCreateBoxMesh(), DMPlexCreateSquareBoundary(), DMSetType(), DMCreate() 311 @*/ 312 PetscErrorCode DMPlexCreateCubeBoundary(DM dm, const PetscReal lower[], const PetscReal upper[], const PetscInt faces[]) 313 { 314 PetscInt vertices[3], numVertices; 315 PetscInt numFaces = 2*faces[0]*faces[1] + 2*faces[1]*faces[2] + 2*faces[0]*faces[2]; 316 Vec coordinates; 317 PetscSection coordSection; 318 PetscScalar *coords; 319 PetscInt coordSize; 320 PetscMPIInt rank; 321 PetscInt v, vx, vy, vz; 322 PetscInt voffset, iface=0, cone[4]; 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 ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);CHKERRQ(ierr); 328 vertices[0] = faces[0]+1; vertices[1] = faces[1]+1; vertices[2] = faces[2]+1; 329 numVertices = vertices[0]*vertices[1]*vertices[2]; 330 if (!rank) { 331 PetscInt f; 332 333 ierr = DMPlexSetChart(dm, 0, numFaces+numVertices);CHKERRQ(ierr); 334 for (f = 0; f < numFaces; ++f) { 335 ierr = DMPlexSetConeSize(dm, f, 4);CHKERRQ(ierr); 336 } 337 ierr = DMSetUp(dm);CHKERRQ(ierr); /* Allocate space for cones */ 338 for (v = 0; v < numFaces+numVertices; ++v) { 339 ierr = DMSetLabelValue(dm, "marker", v, 1);CHKERRQ(ierr); 340 } 341 342 /* Side 0 (Top) */ 343 for (vy = 0; vy < faces[1]; vy++) { 344 for (vx = 0; vx < faces[0]; vx++) { 345 voffset = numFaces + vertices[0]*vertices[1]*(vertices[2]-1) + vy*vertices[0] + vx; 346 cone[0] = voffset; cone[1] = voffset+1; cone[2] = voffset+vertices[0]+1; cone[3] = voffset+vertices[0]; 347 ierr = DMPlexSetCone(dm, iface, cone);CHKERRQ(ierr); 348 iface++; 349 } 350 } 351 352 /* Side 1 (Bottom) */ 353 for (vy = 0; vy < faces[1]; vy++) { 354 for (vx = 0; vx < faces[0]; vx++) { 355 voffset = numFaces + vy*(faces[0]+1) + vx; 356 cone[0] = voffset+1; cone[1] = voffset; cone[2] = voffset+vertices[0]; cone[3] = voffset+vertices[0]+1; 357 ierr = DMPlexSetCone(dm, iface, cone);CHKERRQ(ierr); 358 iface++; 359 } 360 } 361 362 /* Side 2 (Front) */ 363 for (vz = 0; vz < faces[2]; vz++) { 364 for (vx = 0; vx < faces[0]; vx++) { 365 voffset = numFaces + vz*vertices[0]*vertices[1] + vx; 366 cone[0] = voffset; cone[1] = voffset+1; cone[2] = voffset+vertices[0]*vertices[1]+1; cone[3] = voffset+vertices[0]*vertices[1]; 367 ierr = DMPlexSetCone(dm, iface, cone);CHKERRQ(ierr); 368 iface++; 369 } 370 } 371 372 /* Side 3 (Back) */ 373 for (vz = 0; vz < faces[2]; vz++) { 374 for (vx = 0; vx < faces[0]; vx++) { 375 voffset = numFaces + vz*vertices[0]*vertices[1] + vertices[0]*(vertices[1]-1) + vx; 376 cone[0] = voffset+vertices[0]*vertices[1]; cone[1] = voffset+vertices[0]*vertices[1]+1; 377 cone[2] = voffset+1; cone[3] = voffset; 378 ierr = DMPlexSetCone(dm, iface, cone);CHKERRQ(ierr); 379 iface++; 380 } 381 } 382 383 /* Side 4 (Left) */ 384 for (vz = 0; vz < faces[2]; vz++) { 385 for (vy = 0; vy < faces[1]; vy++) { 386 voffset = numFaces + vz*vertices[0]*vertices[1] + vy*vertices[0]; 387 cone[0] = voffset; cone[1] = voffset+vertices[0]*vertices[1]; 388 cone[2] = voffset+vertices[0]*vertices[1]+vertices[0]; cone[3] = voffset+vertices[0]; 389 ierr = DMPlexSetCone(dm, iface, cone);CHKERRQ(ierr); 390 iface++; 391 } 392 } 393 394 /* Side 5 (Right) */ 395 for (vz = 0; vz < faces[2]; vz++) { 396 for (vy = 0; vy < faces[1]; vy++) { 397 voffset = numFaces + vz*vertices[0]*vertices[1] + vy*vertices[0] + vx; 398 cone[0] = voffset+vertices[0]*vertices[1]; cone[1] = voffset; 399 cone[2] = voffset+vertices[0]; cone[3] = voffset+vertices[0]*vertices[1]+vertices[0]; 400 ierr = DMPlexSetCone(dm, iface, cone);CHKERRQ(ierr); 401 iface++; 402 } 403 } 404 } 405 ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr); 406 ierr = DMPlexStratify(dm);CHKERRQ(ierr); 407 /* Build coordinates */ 408 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 409 ierr = PetscSectionSetChart(coordSection, numFaces, numFaces + numVertices);CHKERRQ(ierr); 410 for (v = numFaces; v < numFaces+numVertices; ++v) { 411 ierr = PetscSectionSetDof(coordSection, v, 3);CHKERRQ(ierr); 412 } 413 ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 414 ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 415 ierr = VecCreate(PETSC_COMM_SELF, &coordinates);CHKERRQ(ierr); 416 ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr); 417 ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 418 ierr = VecSetBlockSize(coordinates, 3);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 /*@ 437 DMPlexCreateBoxMesh - Creates a mesh on the tensor product of unit intervals (box) using simplices. 438 439 Collective on MPI_Comm 440 441 Input Parameters: 442 + comm - The communicator for the DM object 443 . dim - The spatial dimension 444 . numFaces - Number of faces per dimension 445 - interpolate - Flag to create intermediate mesh pieces (edges, faces) 446 447 Output Parameter: 448 . dm - The DM object 449 450 Level: beginner 451 452 .keywords: DM, create 453 .seealso: DMPlexCreateHexBoxMesh(), DMSetType(), DMCreate() 454 @*/ 455 PetscErrorCode DMPlexCreateBoxMesh(MPI_Comm comm, PetscInt dim, PetscInt numFaces, PetscBool interpolate, DM *dm) 456 { 457 DM boundary; 458 PetscErrorCode ierr; 459 460 PetscFunctionBegin; 461 PetscValidPointer(dm, 4); 462 ierr = DMCreate(comm, &boundary);CHKERRQ(ierr); 463 PetscValidLogicalCollectiveInt(boundary,dim,2); 464 ierr = DMSetType(boundary, DMPLEX);CHKERRQ(ierr); 465 ierr = DMSetDimension(boundary, dim-1);CHKERRQ(ierr); 466 ierr = DMSetCoordinateDim(boundary, dim);CHKERRQ(ierr); 467 switch (dim) { 468 case 2: 469 { 470 PetscReal lower[2] = {0.0, 0.0}; 471 PetscReal upper[2] = {1.0, 1.0}; 472 PetscInt edges[2]; 473 474 edges[0] = numFaces; edges[1] = numFaces; 475 ierr = DMPlexCreateSquareBoundary(boundary, lower, upper, edges);CHKERRQ(ierr); 476 break; 477 } 478 case 3: 479 { 480 PetscReal lower[3] = {0.0, 0.0, 0.0}; 481 PetscReal upper[3] = {1.0, 1.0, 1.0}; 482 PetscInt faces[3]; 483 484 faces[0] = numFaces; faces[1] = numFaces; faces[2] = numFaces; 485 ierr = DMPlexCreateCubeBoundary(boundary, lower, upper, faces);CHKERRQ(ierr); 486 break; 487 } 488 default: 489 SETERRQ1(comm, PETSC_ERR_SUP, "Dimension not supported: %d", dim); 490 } 491 ierr = DMPlexGenerate(boundary, NULL, interpolate, dm);CHKERRQ(ierr); 492 ierr = DMDestroy(&boundary);CHKERRQ(ierr); 493 PetscFunctionReturn(0); 494 } 495 496 static PetscErrorCode DMPlexCreateCubeMesh_Internal(DM dm, const PetscReal lower[], const PetscReal upper[], const PetscInt edges[], DMBoundaryType bdX, DMBoundaryType bdY, DMBoundaryType bdZ) 497 { 498 DMLabel cutLabel = NULL; 499 PetscInt markerTop = 1, faceMarkerTop = 1; 500 PetscInt markerBottom = 1, faceMarkerBottom = 1; 501 PetscInt markerFront = 1, faceMarkerFront = 1; 502 PetscInt markerBack = 1, faceMarkerBack = 1; 503 PetscInt markerRight = 1, faceMarkerRight = 1; 504 PetscInt markerLeft = 1, faceMarkerLeft = 1; 505 PetscInt dim; 506 PetscBool markerSeparate = PETSC_FALSE, cutMarker = PETSC_FALSE; 507 PetscMPIInt rank; 508 PetscErrorCode ierr; 509 510 PetscFunctionBegin; 511 ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr); 512 ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);CHKERRQ(ierr); 513 ierr = DMCreateLabel(dm,"marker");CHKERRQ(ierr); 514 ierr = DMCreateLabel(dm,"Face Sets");CHKERRQ(ierr); 515 if (bdX == DM_BOUNDARY_PERIODIC || bdX == DM_BOUNDARY_TWIST || 516 bdY == DM_BOUNDARY_PERIODIC || bdY == DM_BOUNDARY_TWIST || 517 bdZ == DM_BOUNDARY_PERIODIC || bdZ == DM_BOUNDARY_TWIST) { 518 ierr = PetscOptionsGetBool(((PetscObject) dm)->options,((PetscObject) dm)->prefix, "-dm_plex_periodic_cut", &cutMarker, NULL);CHKERRQ(ierr); 519 if (cutMarker) {ierr = DMCreateLabel(dm, "periodic_cut");CHKERRQ(ierr); ierr = DMGetLabel(dm, "periodic_cut", &cutLabel);CHKERRQ(ierr);} 520 } 521 switch (dim) { 522 case 2: 523 faceMarkerTop = 3; 524 faceMarkerBottom = 1; 525 faceMarkerRight = 2; 526 faceMarkerLeft = 4; 527 break; 528 case 3: 529 faceMarkerBottom = 1; 530 faceMarkerTop = 2; 531 faceMarkerFront = 3; 532 faceMarkerBack = 4; 533 faceMarkerRight = 5; 534 faceMarkerLeft = 6; 535 break; 536 default: 537 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Dimension %d not supported",dim); 538 break; 539 } 540 ierr = PetscOptionsGetBool(((PetscObject) dm)->options,((PetscObject) dm)->prefix, "-dm_plex_separate_marker", &markerSeparate, NULL);CHKERRQ(ierr); 541 if (markerSeparate) { 542 markerBottom = faceMarkerBottom; 543 markerTop = faceMarkerTop; 544 markerFront = faceMarkerFront; 545 markerBack = faceMarkerBack; 546 markerRight = faceMarkerRight; 547 markerLeft = faceMarkerLeft; 548 } 549 { 550 const PetscInt numXEdges = !rank ? edges[0] : 0; 551 const PetscInt numYEdges = !rank ? edges[1] : 0; 552 const PetscInt numZEdges = !rank ? edges[2] : 0; 553 const PetscInt numXVertices = !rank ? (bdX == DM_BOUNDARY_PERIODIC || bdX == DM_BOUNDARY_TWIST ? edges[0] : edges[0]+1) : 0; 554 const PetscInt numYVertices = !rank ? (bdY == DM_BOUNDARY_PERIODIC || bdY == DM_BOUNDARY_TWIST ? edges[1] : edges[1]+1) : 0; 555 const PetscInt numZVertices = !rank ? (bdZ == DM_BOUNDARY_PERIODIC || bdZ == DM_BOUNDARY_TWIST ? edges[2] : edges[2]+1) : 0; 556 const PetscInt numCells = numXEdges*numYEdges*numZEdges; 557 const PetscInt numXFaces = numYEdges*numZEdges; 558 const PetscInt numYFaces = numXEdges*numZEdges; 559 const PetscInt numZFaces = numXEdges*numYEdges; 560 const PetscInt numTotXFaces = numXVertices*numXFaces; 561 const PetscInt numTotYFaces = numYVertices*numYFaces; 562 const PetscInt numTotZFaces = numZVertices*numZFaces; 563 const PetscInt numFaces = numTotXFaces + numTotYFaces + numTotZFaces; 564 const PetscInt numTotXEdges = numXEdges*numYVertices*numZVertices; 565 const PetscInt numTotYEdges = numYEdges*numXVertices*numZVertices; 566 const PetscInt numTotZEdges = numZEdges*numXVertices*numYVertices; 567 const PetscInt numVertices = numXVertices*numYVertices*numZVertices; 568 const PetscInt numEdges = numTotXEdges + numTotYEdges + numTotZEdges; 569 const PetscInt firstVertex = (dim == 2) ? numFaces : numCells; 570 const PetscInt firstXFace = (dim == 2) ? 0 : numCells + numVertices; 571 const PetscInt firstYFace = firstXFace + numTotXFaces; 572 const PetscInt firstZFace = firstYFace + numTotYFaces; 573 const PetscInt firstXEdge = numCells + numFaces + numVertices; 574 const PetscInt firstYEdge = firstXEdge + numTotXEdges; 575 const PetscInt firstZEdge = firstYEdge + numTotYEdges; 576 Vec coordinates; 577 PetscSection coordSection; 578 PetscScalar *coords; 579 PetscInt coordSize; 580 PetscInt v, vx, vy, vz; 581 PetscInt c, f, fx, fy, fz, e, ex, ey, ez; 582 583 ierr = DMPlexSetChart(dm, 0, numCells+numFaces+numEdges+numVertices);CHKERRQ(ierr); 584 for (c = 0; c < numCells; c++) { 585 ierr = DMPlexSetConeSize(dm, c, 6);CHKERRQ(ierr); 586 } 587 for (f = firstXFace; f < firstXFace+numFaces; ++f) { 588 ierr = DMPlexSetConeSize(dm, f, 4);CHKERRQ(ierr); 589 } 590 for (e = firstXEdge; e < firstXEdge+numEdges; ++e) { 591 ierr = DMPlexSetConeSize(dm, e, 2);CHKERRQ(ierr); 592 } 593 ierr = DMSetUp(dm);CHKERRQ(ierr); /* Allocate space for cones */ 594 /* Build cells */ 595 for (fz = 0; fz < numZEdges; ++fz) { 596 for (fy = 0; fy < numYEdges; ++fy) { 597 for (fx = 0; fx < numXEdges; ++fx) { 598 PetscInt cell = (fz*numYEdges + fy)*numXEdges + fx; 599 PetscInt faceB = firstZFace + (fy*numXEdges+fx)*numZVertices + fz; 600 PetscInt faceT = firstZFace + (fy*numXEdges+fx)*numZVertices + ((fz+1)%numZVertices); 601 PetscInt faceF = firstYFace + (fz*numXEdges+fx)*numYVertices + fy; 602 PetscInt faceK = firstYFace + (fz*numXEdges+fx)*numYVertices + ((fy+1)%numYVertices); 603 PetscInt faceL = firstXFace + (fz*numYEdges+fy)*numXVertices + fx; 604 PetscInt faceR = firstXFace + (fz*numYEdges+fy)*numXVertices + ((fx+1)%numXVertices); 605 /* B, T, F, K, R, L */ 606 PetscInt ornt[6] = {-4, 0, 0, -1, 0, -4}; /* ??? */ 607 PetscInt cone[6]; 608 609 /* no boundary twisting in 3D */ 610 cone[0] = faceB; cone[1] = faceT; cone[2] = faceF; cone[3] = faceK; cone[4] = faceR; cone[5] = faceL; 611 ierr = DMPlexSetCone(dm, cell, cone);CHKERRQ(ierr); 612 ierr = DMPlexSetConeOrientation(dm, cell, ornt);CHKERRQ(ierr); 613 if (bdX != DM_BOUNDARY_NONE && fx == numXEdges-1 && cutLabel) {ierr = DMLabelSetValue(cutLabel, cell, 2);CHKERRQ(ierr);} 614 if (bdY != DM_BOUNDARY_NONE && fy == numYEdges-1 && cutLabel) {ierr = DMLabelSetValue(cutLabel, cell, 2);CHKERRQ(ierr);} 615 if (bdZ != DM_BOUNDARY_NONE && fz == numZEdges-1 && cutLabel) {ierr = DMLabelSetValue(cutLabel, cell, 2);CHKERRQ(ierr);} 616 } 617 } 618 } 619 /* Build x faces */ 620 for (fz = 0; fz < numZEdges; ++fz) { 621 for (fy = 0; fy < numYEdges; ++fy) { 622 for (fx = 0; fx < numXVertices; ++fx) { 623 PetscInt face = firstXFace + (fz*numYEdges+fy)*numXVertices + fx; 624 PetscInt edgeL = firstZEdge + ( fy* numXVertices+fx)*numZEdges + fz; 625 PetscInt edgeR = firstZEdge + (((fy+1)%numYVertices)*numXVertices+fx)*numZEdges + fz; 626 PetscInt edgeB = firstYEdge + ( fz* numXVertices+fx)*numYEdges + fy; 627 PetscInt edgeT = firstYEdge + (((fz+1)%numZVertices)*numXVertices+fx)*numYEdges + fy; 628 PetscInt ornt[4] = {0, 0, -2, -2}; 629 PetscInt cone[4]; 630 631 if (dim == 3) { 632 /* markers */ 633 if (bdX != DM_BOUNDARY_PERIODIC) { 634 if (fx == numXVertices-1) { 635 ierr = DMSetLabelValue(dm, "Face Sets", face, faceMarkerRight);CHKERRQ(ierr); 636 ierr = DMSetLabelValue(dm, "marker", face, markerRight);CHKERRQ(ierr); 637 } 638 else if (fx == 0) { 639 ierr = DMSetLabelValue(dm, "Face Sets", face, faceMarkerLeft);CHKERRQ(ierr); 640 ierr = DMSetLabelValue(dm, "marker", face, markerLeft);CHKERRQ(ierr); 641 } 642 } 643 } 644 cone[0] = edgeB; cone[1] = edgeR; cone[2] = edgeT; cone[3] = edgeL; 645 ierr = DMPlexSetCone(dm, face, cone);CHKERRQ(ierr); 646 ierr = DMPlexSetConeOrientation(dm, face, ornt);CHKERRQ(ierr); 647 } 648 } 649 } 650 /* Build y faces */ 651 for (fz = 0; fz < numZEdges; ++fz) { 652 for (fx = 0; fx < numXEdges; ++fx) { 653 for (fy = 0; fy < numYVertices; ++fy) { 654 PetscInt face = firstYFace + (fz*numXEdges+fx)*numYVertices + fy; 655 PetscInt edgeL = firstZEdge + (fy*numXVertices+ fx )*numZEdges + fz; 656 PetscInt edgeR = firstZEdge + (fy*numXVertices+((fx+1)%numXVertices))*numZEdges + fz; 657 PetscInt edgeB = firstXEdge + ( fz *numYVertices+fy)*numXEdges + fx; 658 PetscInt edgeT = firstXEdge + (((fz+1)%numZVertices)*numYVertices+fy)*numXEdges + fx; 659 PetscInt ornt[4] = {0, 0, -2, -2}; 660 PetscInt cone[4]; 661 662 if (dim == 3) { 663 /* markers */ 664 if (bdY != DM_BOUNDARY_PERIODIC) { 665 if (fy == numYVertices-1) { 666 ierr = DMSetLabelValue(dm, "Face Sets", face, faceMarkerBack);CHKERRQ(ierr); 667 ierr = DMSetLabelValue(dm, "marker", face, markerBack);CHKERRQ(ierr); 668 } 669 else if (fy == 0) { 670 ierr = DMSetLabelValue(dm, "Face Sets", face, faceMarkerFront);CHKERRQ(ierr); 671 ierr = DMSetLabelValue(dm, "marker", face, markerFront);CHKERRQ(ierr); 672 } 673 } 674 } 675 cone[0] = edgeB; cone[1] = edgeR; cone[2] = edgeT; cone[3] = edgeL; 676 ierr = DMPlexSetCone(dm, face, cone);CHKERRQ(ierr); 677 ierr = DMPlexSetConeOrientation(dm, face, ornt);CHKERRQ(ierr); 678 } 679 } 680 } 681 /* Build z faces */ 682 for (fy = 0; fy < numYEdges; ++fy) { 683 for (fx = 0; fx < numXEdges; ++fx) { 684 for (fz = 0; fz < numZVertices; fz++) { 685 PetscInt face = firstZFace + (fy*numXEdges+fx)*numZVertices + fz; 686 PetscInt edgeL = firstYEdge + (fz*numXVertices+ fx )*numYEdges + fy; 687 PetscInt edgeR = firstYEdge + (fz*numXVertices+((fx+1)%numXVertices))*numYEdges + fy; 688 PetscInt edgeB = firstXEdge + (fz*numYVertices+ fy )*numXEdges + fx; 689 PetscInt edgeT = firstXEdge + (fz*numYVertices+((fy+1)%numYVertices))*numXEdges + fx; 690 PetscInt ornt[4] = {0, 0, -2, -2}; 691 PetscInt cone[4]; 692 693 if (dim == 2) { 694 if (bdX == DM_BOUNDARY_TWIST && fx == numXEdges-1) {edgeR += numYEdges-1-2*fy; ornt[1] = -2;} 695 if (bdY == DM_BOUNDARY_TWIST && fy == numYEdges-1) {edgeT += numXEdges-1-2*fx; ornt[2] = 0;} 696 if (bdX != DM_BOUNDARY_NONE && fx == numXEdges-1 && cutLabel) {ierr = DMLabelSetValue(cutLabel, face, 2);CHKERRQ(ierr);} 697 if (bdY != DM_BOUNDARY_NONE && fy == numYEdges-1 && cutLabel) {ierr = DMLabelSetValue(cutLabel, face, 2);CHKERRQ(ierr);} 698 } else { 699 /* markers */ 700 if (bdZ != DM_BOUNDARY_PERIODIC) { 701 if (fz == numZVertices-1) { 702 ierr = DMSetLabelValue(dm, "Face Sets", face, faceMarkerTop);CHKERRQ(ierr); 703 ierr = DMSetLabelValue(dm, "marker", face, markerTop);CHKERRQ(ierr); 704 } 705 else if (fz == 0) { 706 ierr = DMSetLabelValue(dm, "Face Sets", face, faceMarkerBottom);CHKERRQ(ierr); 707 ierr = DMSetLabelValue(dm, "marker", face, markerBottom);CHKERRQ(ierr); 708 } 709 } 710 } 711 cone[0] = edgeB; cone[1] = edgeR; cone[2] = edgeT; cone[3] = edgeL; 712 ierr = DMPlexSetCone(dm, face, cone);CHKERRQ(ierr); 713 ierr = DMPlexSetConeOrientation(dm, face, ornt);CHKERRQ(ierr); 714 } 715 } 716 } 717 /* Build Z edges*/ 718 for (vy = 0; vy < numYVertices; vy++) { 719 for (vx = 0; vx < numXVertices; vx++) { 720 for (ez = 0; ez < numZEdges; ez++) { 721 const PetscInt edge = firstZEdge + (vy*numXVertices+vx)*numZEdges + ez; 722 const PetscInt vertexB = firstVertex + ( ez *numYVertices+vy)*numXVertices + vx; 723 const PetscInt vertexT = firstVertex + (((ez+1)%numZVertices)*numYVertices+vy)*numXVertices + vx; 724 PetscInt cone[2]; 725 726 if (dim == 3) { 727 if (bdX != DM_BOUNDARY_PERIODIC) { 728 if (vx == numXVertices-1) { 729 ierr = DMSetLabelValue(dm, "marker", edge, markerRight);CHKERRQ(ierr); 730 } 731 else if (vx == 0) { 732 ierr = DMSetLabelValue(dm, "marker", edge, markerLeft);CHKERRQ(ierr); 733 } 734 } 735 if (bdY != DM_BOUNDARY_PERIODIC) { 736 if (vy == numYVertices-1) { 737 ierr = DMSetLabelValue(dm, "marker", edge, markerBack);CHKERRQ(ierr); 738 } 739 else if (vy == 0) { 740 ierr = DMSetLabelValue(dm, "marker", edge, markerFront);CHKERRQ(ierr); 741 } 742 } 743 } 744 cone[0] = vertexB; cone[1] = vertexT; 745 ierr = DMPlexSetCone(dm, edge, cone);CHKERRQ(ierr); 746 } 747 } 748 } 749 /* Build Y edges*/ 750 for (vz = 0; vz < numZVertices; vz++) { 751 for (vx = 0; vx < numXVertices; vx++) { 752 for (ey = 0; ey < numYEdges; ey++) { 753 const PetscInt nextv = (dim == 2 && bdY == DM_BOUNDARY_TWIST && ey == numYEdges-1) ? (numXVertices-vx-1) : (vz*numYVertices+((ey+1)%numYVertices))*numXVertices + vx; 754 const PetscInt edge = firstYEdge + (vz*numXVertices+vx)*numYEdges + ey; 755 const PetscInt vertexF = firstVertex + (vz*numYVertices+ey)*numXVertices + vx; 756 const PetscInt vertexK = firstVertex + nextv; 757 PetscInt cone[2]; 758 759 cone[0] = vertexF; cone[1] = vertexK; 760 ierr = DMPlexSetCone(dm, edge, cone);CHKERRQ(ierr); 761 if (dim == 2) { 762 if ((bdX != DM_BOUNDARY_PERIODIC) && (bdX != DM_BOUNDARY_TWIST)) { 763 if (vx == numXVertices-1) { 764 ierr = DMSetLabelValue(dm, "Face Sets", edge, faceMarkerRight);CHKERRQ(ierr); 765 ierr = DMSetLabelValue(dm, "marker", edge, markerRight);CHKERRQ(ierr); 766 ierr = DMSetLabelValue(dm, "marker", cone[0], markerRight);CHKERRQ(ierr); 767 if (ey == numYEdges-1) { 768 ierr = DMSetLabelValue(dm, "marker", cone[1], markerRight);CHKERRQ(ierr); 769 } 770 } else if (vx == 0) { 771 ierr = DMSetLabelValue(dm, "Face Sets", edge, faceMarkerLeft);CHKERRQ(ierr); 772 ierr = DMSetLabelValue(dm, "marker", edge, markerLeft);CHKERRQ(ierr); 773 ierr = DMSetLabelValue(dm, "marker", cone[0], markerLeft);CHKERRQ(ierr); 774 if (ey == numYEdges-1) { 775 ierr = DMSetLabelValue(dm, "marker", cone[1], markerLeft);CHKERRQ(ierr); 776 } 777 } 778 } else { 779 if (vx == 0 && cutMarker) { 780 ierr = DMLabelSetValue(cutLabel, edge, 1);CHKERRQ(ierr); 781 ierr = DMLabelSetValue(cutLabel, cone[0], 1);CHKERRQ(ierr); 782 if (ey == numYEdges-1) { 783 ierr = DMLabelSetValue(cutLabel, cone[1], 1);CHKERRQ(ierr); 784 } 785 } 786 } 787 } else { 788 if (bdX != DM_BOUNDARY_PERIODIC) { 789 if (vx == numXVertices-1) { 790 ierr = DMSetLabelValue(dm, "marker", edge, markerRight);CHKERRQ(ierr); 791 } else if (vx == 0) { 792 ierr = DMSetLabelValue(dm, "marker", edge, markerLeft);CHKERRQ(ierr); 793 } 794 } 795 if (bdZ != DM_BOUNDARY_PERIODIC) { 796 if (vz == numZVertices-1) { 797 ierr = DMSetLabelValue(dm, "marker", edge, markerTop);CHKERRQ(ierr); 798 } else if (vz == 0) { 799 ierr = DMSetLabelValue(dm, "marker", edge, markerBottom);CHKERRQ(ierr); 800 } 801 } 802 } 803 } 804 } 805 } 806 /* Build X edges*/ 807 for (vz = 0; vz < numZVertices; vz++) { 808 for (vy = 0; vy < numYVertices; vy++) { 809 for (ex = 0; ex < numXEdges; ex++) { 810 const PetscInt nextv = (dim == 2 && bdX == DM_BOUNDARY_TWIST && ex == numXEdges-1) ? (numYVertices-vy-1)*numXVertices : (vz*numYVertices+vy)*numXVertices + (ex+1)%numXVertices; 811 const PetscInt edge = firstXEdge + (vz*numYVertices+vy)*numXEdges + ex; 812 const PetscInt vertexL = firstVertex + (vz*numYVertices+vy)*numXVertices + ex; 813 const PetscInt vertexR = firstVertex + nextv; 814 PetscInt cone[2]; 815 816 cone[0] = vertexL; cone[1] = vertexR; 817 ierr = DMPlexSetCone(dm, edge, cone);CHKERRQ(ierr); 818 if (dim == 2) { 819 if ((bdY != DM_BOUNDARY_PERIODIC) && (bdY != DM_BOUNDARY_TWIST)) { 820 if (vy == numYVertices-1) { 821 ierr = DMSetLabelValue(dm, "Face Sets", edge, faceMarkerTop);CHKERRQ(ierr); 822 ierr = DMSetLabelValue(dm, "marker", edge, markerTop);CHKERRQ(ierr); 823 ierr = DMSetLabelValue(dm, "marker", cone[0], markerTop);CHKERRQ(ierr); 824 if (ex == numXEdges-1) { 825 ierr = DMSetLabelValue(dm, "marker", cone[1], markerTop);CHKERRQ(ierr); 826 } 827 } else if (vy == 0) { 828 ierr = DMSetLabelValue(dm, "Face Sets", edge, faceMarkerBottom);CHKERRQ(ierr); 829 ierr = DMSetLabelValue(dm, "marker", edge, markerBottom);CHKERRQ(ierr); 830 ierr = DMSetLabelValue(dm, "marker", cone[0], markerBottom);CHKERRQ(ierr); 831 if (ex == numXEdges-1) { 832 ierr = DMSetLabelValue(dm, "marker", cone[1], markerBottom);CHKERRQ(ierr); 833 } 834 } 835 } else { 836 if (vy == 0 && cutMarker) { 837 ierr = DMLabelSetValue(cutLabel, edge, 1);CHKERRQ(ierr); 838 ierr = DMLabelSetValue(cutLabel, cone[0], 1);CHKERRQ(ierr); 839 if (ex == numXEdges-1) { 840 ierr = DMLabelSetValue(cutLabel, cone[1], 1);CHKERRQ(ierr); 841 } 842 } 843 } 844 } else { 845 if (bdY != DM_BOUNDARY_PERIODIC) { 846 if (vy == numYVertices-1) { 847 ierr = DMSetLabelValue(dm, "marker", edge, markerBack);CHKERRQ(ierr); 848 } 849 else if (vy == 0) { 850 ierr = DMSetLabelValue(dm, "marker", edge, markerFront);CHKERRQ(ierr); 851 } 852 } 853 if (bdZ != DM_BOUNDARY_PERIODIC) { 854 if (vz == numZVertices-1) { 855 ierr = DMSetLabelValue(dm, "marker", edge, markerTop);CHKERRQ(ierr); 856 } 857 else if (vz == 0) { 858 ierr = DMSetLabelValue(dm, "marker", edge, markerBottom);CHKERRQ(ierr); 859 } 860 } 861 } 862 } 863 } 864 } 865 ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr); 866 ierr = DMPlexStratify(dm);CHKERRQ(ierr); 867 /* Build coordinates */ 868 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 869 ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr); 870 ierr = PetscSectionSetFieldComponents(coordSection, 0, dim);CHKERRQ(ierr); 871 ierr = PetscSectionSetChart(coordSection, firstVertex, firstVertex+numVertices);CHKERRQ(ierr); 872 for (v = firstVertex; v < firstVertex+numVertices; ++v) { 873 ierr = PetscSectionSetDof(coordSection, v, dim);CHKERRQ(ierr); 874 ierr = PetscSectionSetFieldDof(coordSection, v, 0, dim);CHKERRQ(ierr); 875 } 876 ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 877 ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 878 ierr = VecCreate(PETSC_COMM_SELF, &coordinates);CHKERRQ(ierr); 879 ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr); 880 ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 881 ierr = VecSetBlockSize(coordinates, dim);CHKERRQ(ierr); 882 ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr); 883 ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 884 for (vz = 0; vz < numZVertices; ++vz) { 885 for (vy = 0; vy < numYVertices; ++vy) { 886 for (vx = 0; vx < numXVertices; ++vx) { 887 coords[((vz*numYVertices+vy)*numXVertices+vx)*dim+0] = lower[0] + ((upper[0] - lower[0])/numXEdges)*vx; 888 coords[((vz*numYVertices+vy)*numXVertices+vx)*dim+1] = lower[1] + ((upper[1] - lower[1])/numYEdges)*vy; 889 if (dim == 3) { 890 coords[((vz*numYVertices+vy)*numXVertices+vx)*dim+2] = lower[2] + ((upper[2] - lower[2])/numZEdges)*vz; 891 } 892 } 893 } 894 } 895 ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 896 ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr); 897 ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 898 } 899 PetscFunctionReturn(0); 900 } 901 902 /*@ 903 DMPlexCreateHexBoxMesh - Creates a mesh on the tensor product of unit intervals (box) using hexahedra. 904 905 Collective on MPI_Comm 906 907 Input Parameters: 908 + comm - The communicator for the DM object 909 . dim - The spatial dimension 910 . periodicX - The boundary type for the X direction 911 . periodicY - The boundary type for the Y direction 912 . periodicZ - The boundary type for the Z direction 913 - cells - The number of cells in each direction 914 915 Output Parameter: 916 . dm - The DM object 917 918 Note: Here is the numbering returned for 2 cells in each direction: 919 $ 22--8-23--9--24 920 $ | | | 921 $ 13 2 14 3 15 922 $ | | | 923 $ 19--6-20--7--21 924 $ | | | 925 $ 10 0 11 1 12 926 $ | | | 927 $ 16--4-17--5--18 928 929 Level: beginner 930 931 .keywords: DM, create 932 .seealso: DMPlexCreateBoxMesh(), DMSetType(), DMCreate() 933 @*/ 934 PetscErrorCode DMPlexCreateHexBoxMesh(MPI_Comm comm, PetscInt dim, const PetscInt cells[], DMBoundaryType periodicX, DMBoundaryType periodicY, DMBoundaryType periodicZ, DM *dm) 935 { 936 PetscInt i; 937 PetscErrorCode ierr; 938 939 PetscFunctionBegin; 940 PetscValidPointer(dm, 7); 941 ierr = DMCreate(comm, dm);CHKERRQ(ierr); 942 PetscValidLogicalCollectiveInt(*dm,dim,2); 943 ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 944 ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr); 945 switch (dim) { 946 case 2: 947 { 948 PetscReal lower[3] = {0.0, 0.0, 0.0}; 949 PetscReal upper[3] = {1.0, 1.0, 0.0}; 950 PetscInt edges[3]; 951 952 edges[0] = cells[0]; 953 edges[1] = cells[1]; 954 edges[2] = 0; 955 956 ierr = DMPlexCreateCubeMesh_Internal(*dm, lower, upper, edges, periodicX, periodicY, DM_BOUNDARY_NONE);CHKERRQ(ierr); 957 if (periodicX == DM_BOUNDARY_PERIODIC || periodicX == DM_BOUNDARY_TWIST || 958 periodicY == DM_BOUNDARY_PERIODIC || periodicY == DM_BOUNDARY_TWIST) { 959 PetscReal L[2]; 960 PetscReal maxCell[2]; 961 DMBoundaryType bdType[2]; 962 963 bdType[0] = periodicX; 964 bdType[1] = periodicY; 965 for (i = 0; i < dim; i++) { 966 L[i] = upper[i] - lower[i]; 967 maxCell[i] = 1.1 * (L[i] / PetscMax(1,cells[i])); 968 } 969 970 ierr = DMSetPeriodicity(*dm,maxCell,L,bdType);CHKERRQ(ierr); 971 } 972 break; 973 } 974 case 3: 975 { 976 PetscReal lower[3] = {0.0, 0.0, 0.0}; 977 PetscReal upper[3] = {1.0, 1.0, 1.0}; 978 979 ierr = DMPlexCreateCubeMesh_Internal(*dm, lower, upper, cells, periodicX, periodicY, periodicZ);CHKERRQ(ierr); 980 if (periodicX == DM_BOUNDARY_PERIODIC || periodicX == DM_BOUNDARY_TWIST || 981 periodicY == DM_BOUNDARY_PERIODIC || periodicY == DM_BOUNDARY_TWIST || 982 periodicZ == DM_BOUNDARY_PERIODIC || periodicZ == DM_BOUNDARY_TWIST) { 983 PetscReal L[3]; 984 PetscReal maxCell[3]; 985 DMBoundaryType bdType[3]; 986 987 bdType[0] = periodicX; 988 bdType[1] = periodicY; 989 bdType[2] = periodicZ; 990 for (i = 0; i < dim; i++) { 991 L[i] = upper[i] - lower[i]; 992 maxCell[i] = 1.1 * (L[i] / cells[i]); 993 } 994 995 ierr = DMSetPeriodicity(*dm,maxCell,L,bdType);CHKERRQ(ierr); 996 } 997 break; 998 } 999 default: 1000 SETERRQ1(comm, PETSC_ERR_SUP, "Dimension not supported: %d", dim); 1001 } 1002 PetscFunctionReturn(0); 1003 } 1004 1005 /*@C 1006 DMPlexSetOptionsPrefix - Sets the prefix used for searching for all DM options in the database. 1007 1008 Logically Collective on DM 1009 1010 Input Parameters: 1011 + dm - the DM context 1012 - prefix - the prefix to prepend to all option names 1013 1014 Notes: 1015 A hyphen (-) must NOT be given at the beginning of the prefix name. 1016 The first character of all runtime options is AUTOMATICALLY the hyphen. 1017 1018 Level: advanced 1019 1020 .seealso: SNESSetFromOptions() 1021 @*/ 1022 PetscErrorCode DMPlexSetOptionsPrefix(DM dm, const char prefix[]) 1023 { 1024 DM_Plex *mesh = (DM_Plex *) dm->data; 1025 PetscErrorCode ierr; 1026 1027 PetscFunctionBegin; 1028 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1029 ierr = PetscObjectSetOptionsPrefix((PetscObject) dm, prefix);CHKERRQ(ierr); 1030 ierr = PetscObjectSetOptionsPrefix((PetscObject) mesh->partitioner, prefix);CHKERRQ(ierr); 1031 PetscFunctionReturn(0); 1032 } 1033 1034 /*@ 1035 DMPlexCreateHexCylinderMesh - Creates a mesh on the tensor product of the unit interval with the circle (cylinder) using hexahedra. 1036 1037 Collective on MPI_Comm 1038 1039 Input Parameters: 1040 + comm - The communicator for the DM object 1041 . numRefine - The number of regular refinements to the basic 5 cell structure 1042 - periodicZ - The boundary type for the Z direction 1043 1044 Output Parameter: 1045 . dm - The DM object 1046 1047 Note: Here is the output numbering looking from the bottom of the cylinder: 1048 $ 17-----14 1049 $ | | 1050 $ | 2 | 1051 $ | | 1052 $ 17-----8-----7-----14 1053 $ | | | | 1054 $ | 3 | 0 | 1 | 1055 $ | | | | 1056 $ 19-----5-----6-----13 1057 $ | | 1058 $ | 4 | 1059 $ | | 1060 $ 19-----13 1061 $ 1062 $ and up through the top 1063 $ 1064 $ 18-----16 1065 $ | | 1066 $ | 2 | 1067 $ | | 1068 $ 18----10----11-----16 1069 $ | | | | 1070 $ | 3 | 0 | 1 | 1071 $ | | | | 1072 $ 20-----9----12-----15 1073 $ | | 1074 $ | 4 | 1075 $ | | 1076 $ 20-----15 1077 1078 Level: beginner 1079 1080 .keywords: DM, create 1081 .seealso: DMPlexCreateHexBoxMesh(), DMPlexCreateBoxMesh(), DMSetType(), DMCreate() 1082 @*/ 1083 PetscErrorCode DMPlexCreateHexCylinderMesh(MPI_Comm comm, PetscInt numRefine, DMBoundaryType periodicZ, DM *dm) 1084 { 1085 const PetscInt dim = 3; 1086 PetscInt numCells, numVertices, r; 1087 PetscErrorCode ierr; 1088 1089 PetscFunctionBegin; 1090 PetscValidPointer(dm, 4); 1091 if (numRefine < 0) SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Number of refinements %D cannot be negative", numRefine); 1092 ierr = DMCreate(comm, dm);CHKERRQ(ierr); 1093 ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 1094 ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr); 1095 /* Create topology */ 1096 { 1097 PetscInt cone[8], c; 1098 1099 numCells = 5; 1100 numVertices = 16; 1101 if (periodicZ == DM_BOUNDARY_PERIODIC) { 1102 numCells *= 3; 1103 numVertices = 24; 1104 } 1105 ierr = DMPlexSetChart(*dm, 0, numCells+numVertices);CHKERRQ(ierr); 1106 for (c = 0; c < numCells; c++) {ierr = DMPlexSetConeSize(*dm, c, 8);CHKERRQ(ierr);} 1107 ierr = DMSetUp(*dm);CHKERRQ(ierr); 1108 if (periodicZ == DM_BOUNDARY_PERIODIC) { 1109 cone[0] = 15; cone[1] = 18; cone[2] = 17; cone[3] = 16; 1110 cone[4] = 31; cone[5] = 32; cone[6] = 33; cone[7] = 34; 1111 ierr = DMPlexSetCone(*dm, 0, cone);CHKERRQ(ierr); 1112 cone[0] = 16; cone[1] = 17; cone[2] = 24; cone[3] = 23; 1113 cone[4] = 32; cone[5] = 36; cone[6] = 37; cone[7] = 33; /* 22 25 26 21 */ 1114 ierr = DMPlexSetCone(*dm, 1, cone);CHKERRQ(ierr); 1115 cone[0] = 18; cone[1] = 27; cone[2] = 24; cone[3] = 17; 1116 cone[4] = 34; cone[5] = 33; cone[6] = 37; cone[7] = 38; 1117 ierr = DMPlexSetCone(*dm, 2, cone);CHKERRQ(ierr); 1118 cone[0] = 29; cone[1] = 27; cone[2] = 18; cone[3] = 15; 1119 cone[4] = 35; cone[5] = 31; cone[6] = 34; cone[7] = 38; 1120 ierr = DMPlexSetCone(*dm, 3, cone);CHKERRQ(ierr); 1121 cone[0] = 29; cone[1] = 15; cone[2] = 16; cone[3] = 23; 1122 cone[4] = 35; cone[5] = 36; cone[6] = 32; cone[7] = 31; 1123 ierr = DMPlexSetCone(*dm, 4, cone);CHKERRQ(ierr); 1124 1125 cone[0] = 31; cone[1] = 34; cone[2] = 33; cone[3] = 32; 1126 cone[4] = 19; cone[5] = 22; cone[6] = 21; cone[7] = 20; 1127 ierr = DMPlexSetCone(*dm, 5, cone);CHKERRQ(ierr); 1128 cone[0] = 32; cone[1] = 33; cone[2] = 37; cone[3] = 36; 1129 cone[4] = 22; cone[5] = 25; cone[6] = 26; cone[7] = 21; 1130 ierr = DMPlexSetCone(*dm, 6, cone);CHKERRQ(ierr); 1131 cone[0] = 34; cone[1] = 38; cone[2] = 37; cone[3] = 33; 1132 cone[4] = 20; cone[5] = 21; cone[6] = 26; cone[7] = 28; 1133 ierr = DMPlexSetCone(*dm, 7, cone);CHKERRQ(ierr); 1134 cone[0] = 35; cone[1] = 38; cone[2] = 34; cone[3] = 31; 1135 cone[4] = 30; cone[5] = 19; cone[6] = 20; cone[7] = 28; 1136 ierr = DMPlexSetCone(*dm, 8, cone);CHKERRQ(ierr); 1137 cone[0] = 35; cone[1] = 31; cone[2] = 32; cone[3] = 36; 1138 cone[4] = 30; cone[5] = 25; cone[6] = 22; cone[7] = 19; 1139 ierr = DMPlexSetCone(*dm, 9, cone);CHKERRQ(ierr); 1140 1141 cone[0] = 19; cone[1] = 20; cone[2] = 21; cone[3] = 22; 1142 cone[4] = 15; cone[5] = 16; cone[6] = 17; cone[7] = 18; 1143 ierr = DMPlexSetCone(*dm, 10, cone);CHKERRQ(ierr); 1144 cone[0] = 22; cone[1] = 21; cone[2] = 26; cone[3] = 25; 1145 cone[4] = 16; cone[5] = 23; cone[6] = 24; cone[7] = 17; 1146 ierr = DMPlexSetCone(*dm, 11, cone);CHKERRQ(ierr); 1147 cone[0] = 20; cone[1] = 28; cone[2] = 26; cone[3] = 21; 1148 cone[4] = 18; cone[5] = 17; cone[6] = 24; cone[7] = 27; 1149 ierr = DMPlexSetCone(*dm, 12, cone);CHKERRQ(ierr); 1150 cone[0] = 30; cone[1] = 28; cone[2] = 20; cone[3] = 19; 1151 cone[4] = 29; cone[5] = 15; cone[6] = 18; cone[7] = 27; 1152 ierr = DMPlexSetCone(*dm, 13, cone);CHKERRQ(ierr); 1153 cone[0] = 30; cone[1] = 19; cone[2] = 22; cone[3] = 25; 1154 cone[4] = 29; cone[5] = 23; cone[6] = 16; cone[7] = 15; 1155 ierr = DMPlexSetCone(*dm, 14, cone);CHKERRQ(ierr); 1156 } else { 1157 cone[0] = 5; cone[1] = 8; cone[2] = 7; cone[3] = 6; 1158 cone[4] = 9; cone[5] = 12; cone[6] = 11; cone[7] = 10; 1159 ierr = DMPlexSetCone(*dm, 0, cone);CHKERRQ(ierr); 1160 cone[0] = 6; cone[1] = 7; cone[2] = 14; cone[3] = 13; 1161 cone[4] = 12; cone[5] = 15; cone[6] = 16; cone[7] = 11; 1162 ierr = DMPlexSetCone(*dm, 1, cone);CHKERRQ(ierr); 1163 cone[0] = 8; cone[1] = 17; cone[2] = 14; cone[3] = 7; 1164 cone[4] = 10; cone[5] = 11; cone[6] = 16; cone[7] = 18; 1165 ierr = DMPlexSetCone(*dm, 2, cone);CHKERRQ(ierr); 1166 cone[0] = 19; cone[1] = 17; cone[2] = 8; cone[3] = 5; 1167 cone[4] = 20; cone[5] = 9; cone[6] = 10; cone[7] = 18; 1168 ierr = DMPlexSetCone(*dm, 3, cone);CHKERRQ(ierr); 1169 cone[0] = 19; cone[1] = 5; cone[2] = 6; cone[3] = 13; 1170 cone[4] = 20; cone[5] = 15; cone[6] = 12; cone[7] = 9; 1171 ierr = DMPlexSetCone(*dm, 4, cone);CHKERRQ(ierr); 1172 } 1173 ierr = DMPlexSymmetrize(*dm);CHKERRQ(ierr); 1174 ierr = DMPlexStratify(*dm);CHKERRQ(ierr); 1175 } 1176 /* Interpolate */ 1177 { 1178 DM idm = NULL; 1179 1180 ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr); 1181 ierr = DMDestroy(dm);CHKERRQ(ierr); 1182 *dm = idm; 1183 } 1184 /* Create cube geometry */ 1185 { 1186 Vec coordinates; 1187 PetscSection coordSection; 1188 PetscScalar *coords; 1189 PetscInt coordSize, v; 1190 const PetscReal dis = 1.0/PetscSqrtReal(2.0); 1191 const PetscReal ds2 = dis/2.0; 1192 1193 /* Build coordinates */ 1194 ierr = DMGetCoordinateSection(*dm, &coordSection);CHKERRQ(ierr); 1195 ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr); 1196 ierr = PetscSectionSetFieldComponents(coordSection, 0, dim);CHKERRQ(ierr); 1197 ierr = PetscSectionSetChart(coordSection, numCells, numCells+numVertices);CHKERRQ(ierr); 1198 for (v = numCells; v < numCells+numVertices; ++v) { 1199 ierr = PetscSectionSetDof(coordSection, v, dim);CHKERRQ(ierr); 1200 ierr = PetscSectionSetFieldDof(coordSection, v, 0, dim);CHKERRQ(ierr); 1201 } 1202 ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 1203 ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 1204 ierr = VecCreate(PETSC_COMM_SELF, &coordinates);CHKERRQ(ierr); 1205 ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr); 1206 ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 1207 ierr = VecSetBlockSize(coordinates, dim);CHKERRQ(ierr); 1208 ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr); 1209 ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 1210 coords[0*dim+0] = -ds2; coords[0*dim+1] = -ds2; coords[0*dim+2] = 0.0; 1211 coords[1*dim+0] = ds2; coords[1*dim+1] = -ds2; coords[1*dim+2] = 0.0; 1212 coords[2*dim+0] = ds2; coords[2*dim+1] = ds2; coords[2*dim+2] = 0.0; 1213 coords[3*dim+0] = -ds2; coords[3*dim+1] = ds2; coords[3*dim+2] = 0.0; 1214 coords[4*dim+0] = -ds2; coords[4*dim+1] = -ds2; coords[4*dim+2] = 1.0; 1215 coords[5*dim+0] = -ds2; coords[5*dim+1] = ds2; coords[5*dim+2] = 1.0; 1216 coords[6*dim+0] = ds2; coords[6*dim+1] = ds2; coords[6*dim+2] = 1.0; 1217 coords[7*dim+0] = ds2; coords[7*dim+1] = -ds2; coords[7*dim+2] = 1.0; 1218 coords[ 8*dim+0] = dis; coords[ 8*dim+1] = -dis; coords[ 8*dim+2] = 0.0; 1219 coords[ 9*dim+0] = dis; coords[ 9*dim+1] = dis; coords[ 9*dim+2] = 0.0; 1220 coords[10*dim+0] = dis; coords[10*dim+1] = -dis; coords[10*dim+2] = 1.0; 1221 coords[11*dim+0] = dis; coords[11*dim+1] = dis; coords[11*dim+2] = 1.0; 1222 coords[12*dim+0] = -dis; coords[12*dim+1] = dis; coords[12*dim+2] = 0.0; 1223 coords[13*dim+0] = -dis; coords[13*dim+1] = dis; coords[13*dim+2] = 1.0; 1224 coords[14*dim+0] = -dis; coords[14*dim+1] = -dis; coords[14*dim+2] = 0.0; 1225 coords[15*dim+0] = -dis; coords[15*dim+1] = -dis; coords[15*dim+2] = 1.0; 1226 if (periodicZ == DM_BOUNDARY_PERIODIC) { 1227 /* 15 31 19 */ coords[16*dim+0] = -ds2; coords[16*dim+1] = -ds2; coords[16*dim+2] = 0.5; 1228 /* 16 32 22 */ coords[17*dim+0] = ds2; coords[17*dim+1] = -ds2; coords[17*dim+2] = 0.5; 1229 /* 17 33 21 */ coords[18*dim+0] = ds2; coords[18*dim+1] = ds2; coords[18*dim+2] = 0.5; 1230 /* 18 34 20 */ coords[19*dim+0] = -ds2; coords[19*dim+1] = ds2; coords[19*dim+2] = 0.5; 1231 /* 29 35 30 */ coords[20*dim+0] = -dis; coords[20*dim+1] = -dis; coords[20*dim+2] = 0.5; 1232 /* 23 36 25 */ coords[21*dim+0] = dis; coords[21*dim+1] = -dis; coords[21*dim+2] = 0.5; 1233 /* 24 37 26 */ coords[22*dim+0] = dis; coords[22*dim+1] = dis; coords[22*dim+2] = 0.5; 1234 /* 27 38 28 */ coords[23*dim+0] = -dis; coords[23*dim+1] = dis; coords[23*dim+2] = 0.5; 1235 } 1236 ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 1237 ierr = DMSetCoordinatesLocal(*dm, coordinates);CHKERRQ(ierr); 1238 ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 1239 } 1240 /* Create periodicity */ 1241 if (periodicZ == DM_BOUNDARY_PERIODIC || periodicZ == DM_BOUNDARY_TWIST) { 1242 PetscReal L[3]; 1243 PetscReal maxCell[3]; 1244 DMBoundaryType bdType[3]; 1245 PetscReal lower[3] = {0.0, 0.0, 0.0}; 1246 PetscReal upper[3] = {1.0, 1.0, 1.5}; 1247 PetscInt i, numZCells = 3; 1248 1249 bdType[0] = DM_BOUNDARY_NONE; 1250 bdType[1] = DM_BOUNDARY_NONE; 1251 bdType[2] = periodicZ; 1252 for (i = 0; i < dim; i++) { 1253 L[i] = upper[i] - lower[i]; 1254 maxCell[i] = 1.1 * (L[i] / numZCells); 1255 } 1256 ierr = DMSetPeriodicity(*dm, maxCell, L, bdType);CHKERRQ(ierr); 1257 } 1258 /* Refine topology */ 1259 for (r = 0; r < numRefine; ++r) { 1260 DM rdm = NULL; 1261 1262 ierr = DMRefine(*dm, comm, &rdm);CHKERRQ(ierr); 1263 ierr = DMDestroy(dm);CHKERRQ(ierr); 1264 *dm = rdm; 1265 } 1266 /* Remap geometry to cylinder 1267 Interior square: Linear interpolation is correct 1268 The other cells all have vertices on rays from the origin. We want to uniformly expand the spacing 1269 such that the last vertex is on the unit circle. So the closest and farthest vertices are at distance 1270 1271 phi = arctan(y/x) 1272 d_close = sqrt(1/8 + 1/4 sin^2(phi)) 1273 d_far = sqrt(1/2 + sin^2(phi)) 1274 1275 so we remap them using 1276 1277 x_new = x_close + (x - x_close) (1 - d_close) / (d_far - d_close) 1278 y_new = y_close + (y - y_close) (1 - d_close) / (d_far - d_close) 1279 1280 If pi/4 < phi < 3pi/4 or -3pi/4 < phi < -pi/4, then we switch x and y. 1281 */ 1282 { 1283 Vec coordinates; 1284 PetscSection coordSection; 1285 PetscScalar *coords; 1286 PetscInt vStart, vEnd, v; 1287 const PetscReal dis = 1.0/PetscSqrtReal(2.0); 1288 const PetscReal ds2 = 0.5*dis; 1289 1290 ierr = DMPlexGetDepthStratum(*dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 1291 ierr = DMGetCoordinateSection(*dm, &coordSection);CHKERRQ(ierr); 1292 ierr = DMGetCoordinatesLocal(*dm, &coordinates);CHKERRQ(ierr); 1293 ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 1294 for (v = vStart; v < vEnd; ++v) { 1295 PetscReal phi, sinp, cosp, dc, df, x, y, xc, yc; 1296 PetscInt off; 1297 1298 ierr = PetscSectionGetOffset(coordSection, v, &off);CHKERRQ(ierr); 1299 if ((PetscAbsScalar(coords[off+0]) <= ds2) && (PetscAbsScalar(coords[off+1]) <= ds2)) continue; 1300 x = PetscRealPart(coords[off]); 1301 y = PetscRealPart(coords[off+1]); 1302 phi = PetscAtan2Real(y, x); 1303 sinp = PetscSinReal(phi); 1304 cosp = PetscCosReal(phi); 1305 if ((PetscAbsReal(phi) > PETSC_PI/4.0) && (PetscAbsReal(phi) < 3.0*PETSC_PI/4.0)) { 1306 dc = PetscAbsReal(ds2/sinp); 1307 df = PetscAbsReal(dis/sinp); 1308 xc = ds2*x/PetscAbsScalar(y); 1309 yc = ds2*PetscSignReal(y); 1310 } else { 1311 dc = PetscAbsReal(ds2/cosp); 1312 df = PetscAbsReal(dis/cosp); 1313 xc = ds2*PetscSignReal(x); 1314 yc = ds2*y/PetscAbsScalar(x); 1315 } 1316 coords[off+0] = xc + (coords[off+0] - xc)*(1.0 - dc)/(df - dc); 1317 coords[off+1] = yc + (coords[off+1] - yc)*(1.0 - dc)/(df - dc); 1318 } 1319 ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 1320 if (periodicZ == DM_BOUNDARY_PERIODIC || periodicZ == DM_BOUNDARY_TWIST) { 1321 ierr = DMLocalizeCoordinates(*dm);CHKERRQ(ierr); 1322 } 1323 } 1324 PetscFunctionReturn(0); 1325 } 1326 1327 /* External function declarations here */ 1328 extern PetscErrorCode DMCreateInterpolation_Plex(DM dmCoarse, DM dmFine, Mat *interpolation, Vec *scaling); 1329 extern PetscErrorCode DMCreateInjection_Plex(DM dmCoarse, DM dmFine, Mat *mat); 1330 extern PetscErrorCode DMCreateDefaultSection_Plex(DM dm); 1331 extern PetscErrorCode DMCreateDefaultConstraints_Plex(DM dm); 1332 extern PetscErrorCode DMCreateMatrix_Plex(DM dm, Mat *J); 1333 extern PetscErrorCode DMCreateCoordinateDM_Plex(DM dm, DM *cdm); 1334 PETSC_INTERN PetscErrorCode DMClone_Plex(DM dm, DM *newdm); 1335 extern PetscErrorCode DMSetUp_Plex(DM dm); 1336 extern PetscErrorCode DMDestroy_Plex(DM dm); 1337 extern PetscErrorCode DMView_Plex(DM dm, PetscViewer viewer); 1338 extern PetscErrorCode DMLoad_Plex(DM dm, PetscViewer viewer); 1339 extern PetscErrorCode DMCreateSubDM_Plex(DM dm, PetscInt numFields, PetscInt fields[], IS *is, DM *subdm); 1340 1341 /* Replace dm with the contents of dmNew 1342 - Share the DM_Plex structure 1343 - Share the coordinates 1344 - Share the SF 1345 */ 1346 static PetscErrorCode DMPlexReplace_Static(DM dm, DM dmNew) 1347 { 1348 PetscSF sf; 1349 DM coordDM, coarseDM; 1350 Vec coords; 1351 const PetscReal *maxCell, *L; 1352 const DMBoundaryType *bd; 1353 PetscErrorCode ierr; 1354 1355 PetscFunctionBegin; 1356 ierr = DMGetPointSF(dmNew, &sf);CHKERRQ(ierr); 1357 ierr = DMSetPointSF(dm, sf);CHKERRQ(ierr); 1358 ierr = DMGetCoordinateDM(dmNew, &coordDM);CHKERRQ(ierr); 1359 ierr = DMGetCoordinatesLocal(dmNew, &coords);CHKERRQ(ierr); 1360 ierr = DMSetCoordinateDM(dm, coordDM);CHKERRQ(ierr); 1361 ierr = DMSetCoordinatesLocal(dm, coords);CHKERRQ(ierr); 1362 ierr = DMGetPeriodicity(dm, &maxCell, &L, &bd);CHKERRQ(ierr); 1363 if (L) {ierr = DMSetPeriodicity(dmNew, maxCell, L, bd);CHKERRQ(ierr);} 1364 ierr = DMDestroy_Plex(dm);CHKERRQ(ierr); 1365 dm->data = dmNew->data; 1366 ((DM_Plex *) dmNew->data)->refct++; 1367 dmNew->labels->refct++; 1368 if (!--(dm->labels->refct)) { 1369 DMLabelLink next = dm->labels->next; 1370 1371 /* destroy the labels */ 1372 while (next) { 1373 DMLabelLink tmp = next->next; 1374 1375 ierr = DMLabelDestroy(&next->label);CHKERRQ(ierr); 1376 ierr = PetscFree(next);CHKERRQ(ierr); 1377 next = tmp; 1378 } 1379 ierr = PetscFree(dm->labels);CHKERRQ(ierr); 1380 } 1381 dm->labels = dmNew->labels; 1382 dm->depthLabel = dmNew->depthLabel; 1383 ierr = DMGetCoarseDM(dmNew,&coarseDM);CHKERRQ(ierr); 1384 ierr = DMSetCoarseDM(dm,coarseDM);CHKERRQ(ierr); 1385 PetscFunctionReturn(0); 1386 } 1387 1388 /* Swap dm with the contents of dmNew 1389 - Swap the DM_Plex structure 1390 - Swap the coordinates 1391 - Swap the point PetscSF 1392 */ 1393 static PetscErrorCode DMPlexSwap_Static(DM dmA, DM dmB) 1394 { 1395 DM coordDMA, coordDMB; 1396 Vec coordsA, coordsB; 1397 PetscSF sfA, sfB; 1398 void *tmp; 1399 DMLabelLinkList listTmp; 1400 DMLabel depthTmp; 1401 PetscErrorCode ierr; 1402 1403 PetscFunctionBegin; 1404 ierr = DMGetPointSF(dmA, &sfA);CHKERRQ(ierr); 1405 ierr = DMGetPointSF(dmB, &sfB);CHKERRQ(ierr); 1406 ierr = PetscObjectReference((PetscObject) sfA);CHKERRQ(ierr); 1407 ierr = DMSetPointSF(dmA, sfB);CHKERRQ(ierr); 1408 ierr = DMSetPointSF(dmB, sfA);CHKERRQ(ierr); 1409 ierr = PetscObjectDereference((PetscObject) sfA);CHKERRQ(ierr); 1410 1411 ierr = DMGetCoordinateDM(dmA, &coordDMA);CHKERRQ(ierr); 1412 ierr = DMGetCoordinateDM(dmB, &coordDMB);CHKERRQ(ierr); 1413 ierr = PetscObjectReference((PetscObject) coordDMA);CHKERRQ(ierr); 1414 ierr = DMSetCoordinateDM(dmA, coordDMB);CHKERRQ(ierr); 1415 ierr = DMSetCoordinateDM(dmB, coordDMA);CHKERRQ(ierr); 1416 ierr = PetscObjectDereference((PetscObject) coordDMA);CHKERRQ(ierr); 1417 1418 ierr = DMGetCoordinatesLocal(dmA, &coordsA);CHKERRQ(ierr); 1419 ierr = DMGetCoordinatesLocal(dmB, &coordsB);CHKERRQ(ierr); 1420 ierr = PetscObjectReference((PetscObject) coordsA);CHKERRQ(ierr); 1421 ierr = DMSetCoordinatesLocal(dmA, coordsB);CHKERRQ(ierr); 1422 ierr = DMSetCoordinatesLocal(dmB, coordsA);CHKERRQ(ierr); 1423 ierr = PetscObjectDereference((PetscObject) coordsA);CHKERRQ(ierr); 1424 1425 tmp = dmA->data; 1426 dmA->data = dmB->data; 1427 dmB->data = tmp; 1428 listTmp = dmA->labels; 1429 dmA->labels = dmB->labels; 1430 dmB->labels = listTmp; 1431 depthTmp = dmA->depthLabel; 1432 dmA->depthLabel = dmB->depthLabel; 1433 dmB->depthLabel = depthTmp; 1434 PetscFunctionReturn(0); 1435 } 1436 1437 PetscErrorCode DMSetFromOptions_NonRefinement_Plex(PetscOptionItems *PetscOptionsObject,DM dm) 1438 { 1439 DM_Plex *mesh = (DM_Plex*) dm->data; 1440 PetscErrorCode ierr; 1441 1442 PetscFunctionBegin; 1443 /* Handle viewing */ 1444 ierr = PetscOptionsBool("-dm_plex_print_set_values", "Output all set values info", "DMView", PETSC_FALSE, &mesh->printSetValues, NULL);CHKERRQ(ierr); 1445 ierr = PetscOptionsInt("-dm_plex_print_fem", "Debug output level all fem computations", "DMView", 0, &mesh->printFEM, NULL);CHKERRQ(ierr); 1446 ierr = PetscOptionsReal("-dm_plex_print_tol", "Tolerance for FEM output", "DMView", mesh->printTol, &mesh->printTol, NULL);CHKERRQ(ierr); 1447 /* Point Location */ 1448 ierr = PetscOptionsBool("-dm_plex_hash_location", "Use grid hashing for point location", "DMView", PETSC_FALSE, &mesh->useHashLocation, NULL);CHKERRQ(ierr); 1449 /* Generation and remeshing */ 1450 ierr = PetscOptionsBool("-dm_plex_remesh_bd", "Allow changes to the boundary on remeshing", "DMView", PETSC_FALSE, &mesh->remeshBd, NULL);CHKERRQ(ierr); 1451 /* Projection behavior */ 1452 ierr = PetscOptionsInt("-dm_plex_max_projection_height", "Maxmimum mesh point height used to project locally", "DMPlexSetMaxProjectionHeight", 0, &mesh->maxProjectionHeight, NULL);CHKERRQ(ierr); 1453 ierr = PetscOptionsBool("-dm_plex_regular_refinement", "Use special nested projection algorithm for regular refinement", "DMPlexSetRegularRefinement", mesh->regularRefinement, &mesh->regularRefinement, NULL);CHKERRQ(ierr); 1454 1455 ierr = PetscPartitionerSetFromOptions(mesh->partitioner);CHKERRQ(ierr); 1456 PetscFunctionReturn(0); 1457 } 1458 1459 static PetscErrorCode DMSetFromOptions_Plex(PetscOptionItems *PetscOptionsObject,DM dm) 1460 { 1461 PetscInt refine = 0, coarsen = 0, r; 1462 PetscBool isHierarchy; 1463 PetscErrorCode ierr; 1464 1465 PetscFunctionBegin; 1466 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1467 ierr = PetscOptionsHead(PetscOptionsObject,"DMPlex Options");CHKERRQ(ierr); 1468 /* Handle DMPlex refinement */ 1469 ierr = PetscOptionsInt("-dm_refine", "The number of uniform refinements", "DMCreate", refine, &refine, NULL);CHKERRQ(ierr); 1470 ierr = PetscOptionsInt("-dm_refine_hierarchy", "The number of uniform refinements", "DMCreate", refine, &refine, &isHierarchy);CHKERRQ(ierr); 1471 if (refine) {ierr = DMPlexSetRefinementUniform(dm, PETSC_TRUE);CHKERRQ(ierr);} 1472 if (refine && isHierarchy) { 1473 DM *dms, coarseDM; 1474 1475 ierr = DMGetCoarseDM(dm, &coarseDM);CHKERRQ(ierr); 1476 ierr = PetscObjectReference((PetscObject)coarseDM);CHKERRQ(ierr); 1477 ierr = PetscMalloc1(refine,&dms);CHKERRQ(ierr); 1478 ierr = DMRefineHierarchy(dm, refine, dms);CHKERRQ(ierr); 1479 /* Total hack since we do not pass in a pointer */ 1480 ierr = DMPlexSwap_Static(dm, dms[refine-1]);CHKERRQ(ierr); 1481 if (refine == 1) { 1482 ierr = DMSetCoarseDM(dm, dms[0]);CHKERRQ(ierr); 1483 ierr = DMPlexSetRegularRefinement(dm, PETSC_TRUE);CHKERRQ(ierr); 1484 } else { 1485 ierr = DMSetCoarseDM(dm, dms[refine-2]);CHKERRQ(ierr); 1486 ierr = DMPlexSetRegularRefinement(dm, PETSC_TRUE);CHKERRQ(ierr); 1487 ierr = DMSetCoarseDM(dms[0], dms[refine-1]);CHKERRQ(ierr); 1488 ierr = DMPlexSetRegularRefinement(dms[0], PETSC_TRUE);CHKERRQ(ierr); 1489 } 1490 ierr = DMSetCoarseDM(dms[refine-1], coarseDM);CHKERRQ(ierr); 1491 ierr = PetscObjectDereference((PetscObject)coarseDM);CHKERRQ(ierr); 1492 /* Free DMs */ 1493 for (r = 0; r < refine; ++r) { 1494 ierr = DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dms[r]);CHKERRQ(ierr); 1495 ierr = DMDestroy(&dms[r]);CHKERRQ(ierr); 1496 } 1497 ierr = PetscFree(dms);CHKERRQ(ierr); 1498 } else { 1499 for (r = 0; r < refine; ++r) { 1500 DM refinedMesh; 1501 1502 ierr = DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dm);CHKERRQ(ierr); 1503 ierr = DMRefine(dm, PetscObjectComm((PetscObject) dm), &refinedMesh);CHKERRQ(ierr); 1504 /* Total hack since we do not pass in a pointer */ 1505 ierr = DMPlexReplace_Static(dm, refinedMesh);CHKERRQ(ierr); 1506 ierr = DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dm);CHKERRQ(ierr); 1507 ierr = DMDestroy(&refinedMesh);CHKERRQ(ierr); 1508 } 1509 } 1510 /* Handle DMPlex coarsening */ 1511 ierr = PetscOptionsInt("-dm_coarsen", "Coarsen the mesh", "DMCreate", coarsen, &coarsen, NULL);CHKERRQ(ierr); 1512 ierr = PetscOptionsInt("-dm_coarsen_hierarchy", "The number of coarsenings", "DMCreate", coarsen, &coarsen, &isHierarchy);CHKERRQ(ierr); 1513 if (coarsen && isHierarchy) { 1514 DM *dms; 1515 1516 ierr = PetscMalloc1(coarsen, &dms);CHKERRQ(ierr); 1517 ierr = DMCoarsenHierarchy(dm, coarsen, dms);CHKERRQ(ierr); 1518 /* Free DMs */ 1519 for (r = 0; r < coarsen; ++r) { 1520 ierr = DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dms[r]);CHKERRQ(ierr); 1521 ierr = DMDestroy(&dms[r]);CHKERRQ(ierr); 1522 } 1523 ierr = PetscFree(dms);CHKERRQ(ierr); 1524 } else { 1525 for (r = 0; r < coarsen; ++r) { 1526 DM coarseMesh; 1527 1528 ierr = DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dm);CHKERRQ(ierr); 1529 ierr = DMCoarsen(dm, PetscObjectComm((PetscObject) dm), &coarseMesh);CHKERRQ(ierr); 1530 /* Total hack since we do not pass in a pointer */ 1531 ierr = DMPlexReplace_Static(dm, coarseMesh);CHKERRQ(ierr); 1532 ierr = DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dm);CHKERRQ(ierr); 1533 ierr = DMDestroy(&coarseMesh);CHKERRQ(ierr); 1534 } 1535 } 1536 /* Handle */ 1537 ierr = DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dm);CHKERRQ(ierr); 1538 ierr = PetscOptionsTail();CHKERRQ(ierr); 1539 PetscFunctionReturn(0); 1540 } 1541 1542 static PetscErrorCode DMCreateGlobalVector_Plex(DM dm,Vec *vec) 1543 { 1544 PetscErrorCode ierr; 1545 1546 PetscFunctionBegin; 1547 ierr = DMCreateGlobalVector_Section_Private(dm,vec);CHKERRQ(ierr); 1548 /* ierr = VecSetOperation(*vec, VECOP_DUPLICATE, (void(*)(void)) VecDuplicate_MPI_DM);CHKERRQ(ierr); */ 1549 ierr = VecSetOperation(*vec, VECOP_VIEW, (void (*)(void)) VecView_Plex);CHKERRQ(ierr); 1550 ierr = VecSetOperation(*vec, VECOP_VIEWNATIVE, (void (*)(void)) VecView_Plex_Native);CHKERRQ(ierr); 1551 ierr = VecSetOperation(*vec, VECOP_LOAD, (void (*)(void)) VecLoad_Plex);CHKERRQ(ierr); 1552 ierr = VecSetOperation(*vec, VECOP_LOADNATIVE, (void (*)(void)) VecLoad_Plex_Native);CHKERRQ(ierr); 1553 PetscFunctionReturn(0); 1554 } 1555 1556 static PetscErrorCode DMCreateLocalVector_Plex(DM dm,Vec *vec) 1557 { 1558 PetscErrorCode ierr; 1559 1560 PetscFunctionBegin; 1561 ierr = DMCreateLocalVector_Section_Private(dm,vec);CHKERRQ(ierr); 1562 ierr = VecSetOperation(*vec, VECOP_VIEW, (void (*)(void)) VecView_Plex_Local);CHKERRQ(ierr); 1563 ierr = VecSetOperation(*vec, VECOP_LOAD, (void (*)(void)) VecLoad_Plex_Local);CHKERRQ(ierr); 1564 PetscFunctionReturn(0); 1565 } 1566 1567 static PetscErrorCode DMGetDimPoints_Plex(DM dm, PetscInt dim, PetscInt *pStart, PetscInt *pEnd) 1568 { 1569 PetscInt depth, d; 1570 PetscErrorCode ierr; 1571 1572 PetscFunctionBegin; 1573 ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 1574 if (depth == 1) { 1575 ierr = DMGetDimension(dm, &d);CHKERRQ(ierr); 1576 if (dim == 0) {ierr = DMPlexGetDepthStratum(dm, dim, pStart, pEnd);CHKERRQ(ierr);} 1577 else if (dim == d) {ierr = DMPlexGetDepthStratum(dm, 1, pStart, pEnd);CHKERRQ(ierr);} 1578 else {*pStart = 0; *pEnd = 0;} 1579 } else { 1580 ierr = DMPlexGetDepthStratum(dm, dim, pStart, pEnd);CHKERRQ(ierr); 1581 } 1582 PetscFunctionReturn(0); 1583 } 1584 1585 static PetscErrorCode DMGetNeighors_Plex(DM dm, PetscInt *nranks, const PetscMPIInt *ranks[]) 1586 { 1587 PetscSF sf; 1588 PetscErrorCode ierr; 1589 1590 PetscFunctionBegin; 1591 ierr = DMGetPointSF(dm, &sf);CHKERRQ(ierr); 1592 ierr = PetscSFGetRanks(sf, nranks, ranks, NULL, NULL, NULL);CHKERRQ(ierr); 1593 PetscFunctionReturn(0); 1594 } 1595 1596 static PetscErrorCode DMInitialize_Plex(DM dm) 1597 { 1598 PetscErrorCode ierr; 1599 1600 PetscFunctionBegin; 1601 dm->ops->view = DMView_Plex; 1602 dm->ops->load = DMLoad_Plex; 1603 dm->ops->setfromoptions = DMSetFromOptions_Plex; 1604 dm->ops->clone = DMClone_Plex; 1605 dm->ops->setup = DMSetUp_Plex; 1606 dm->ops->createdefaultsection = DMCreateDefaultSection_Plex; 1607 dm->ops->createdefaultconstraints = DMCreateDefaultConstraints_Plex; 1608 dm->ops->createglobalvector = DMCreateGlobalVector_Plex; 1609 dm->ops->createlocalvector = DMCreateLocalVector_Plex; 1610 dm->ops->getlocaltoglobalmapping = NULL; 1611 dm->ops->createfieldis = NULL; 1612 dm->ops->createcoordinatedm = DMCreateCoordinateDM_Plex; 1613 dm->ops->getcoloring = NULL; 1614 dm->ops->creatematrix = DMCreateMatrix_Plex; 1615 dm->ops->createinterpolation = DMCreateInterpolation_Plex; 1616 dm->ops->getaggregates = NULL; 1617 dm->ops->getinjection = DMCreateInjection_Plex; 1618 dm->ops->refine = DMRefine_Plex; 1619 dm->ops->coarsen = DMCoarsen_Plex; 1620 dm->ops->refinehierarchy = DMRefineHierarchy_Plex; 1621 dm->ops->coarsenhierarchy = DMCoarsenHierarchy_Plex; 1622 dm->ops->globaltolocalbegin = NULL; 1623 dm->ops->globaltolocalend = NULL; 1624 dm->ops->localtoglobalbegin = NULL; 1625 dm->ops->localtoglobalend = NULL; 1626 dm->ops->destroy = DMDestroy_Plex; 1627 dm->ops->createsubdm = DMCreateSubDM_Plex; 1628 dm->ops->getdimpoints = DMGetDimPoints_Plex; 1629 dm->ops->locatepoints = DMLocatePoints_Plex; 1630 dm->ops->projectfunctionlocal = DMProjectFunctionLocal_Plex; 1631 dm->ops->projectfunctionlabellocal = DMProjectFunctionLabelLocal_Plex; 1632 dm->ops->projectfieldlocal = DMProjectFieldLocal_Plex; 1633 dm->ops->projectfieldlabellocal = DMProjectFieldLabelLocal_Plex; 1634 dm->ops->computel2diff = DMComputeL2Diff_Plex; 1635 dm->ops->computel2gradientdiff = DMComputeL2GradientDiff_Plex; 1636 dm->ops->computel2fielddiff = DMComputeL2FieldDiff_Plex; 1637 dm->ops->getneighbors = DMGetNeighors_Plex; 1638 ierr = PetscObjectComposeFunction((PetscObject)dm,"DMAdaptLabel_C",DMAdaptLabel_Plex);CHKERRQ(ierr); 1639 PetscFunctionReturn(0); 1640 } 1641 1642 PETSC_INTERN PetscErrorCode DMClone_Plex(DM dm, DM *newdm) 1643 { 1644 DM_Plex *mesh = (DM_Plex *) dm->data; 1645 PetscErrorCode ierr; 1646 1647 PetscFunctionBegin; 1648 mesh->refct++; 1649 (*newdm)->data = mesh; 1650 ierr = PetscObjectChangeTypeName((PetscObject) *newdm, DMPLEX);CHKERRQ(ierr); 1651 ierr = DMInitialize_Plex(*newdm);CHKERRQ(ierr); 1652 PetscFunctionReturn(0); 1653 } 1654 1655 /*MC 1656 DMPLEX = "plex" - A DM object that encapsulates an unstructured mesh, or CW Complex, which can be expressed using a Hasse Diagram. 1657 In the local representation, Vecs contain all unknowns in the interior and shared boundary. This is 1658 specified by a PetscSection object. Ownership in the global representation is determined by 1659 ownership of the underlying DMPlex points. This is specified by another PetscSection object. 1660 1661 Level: intermediate 1662 1663 .seealso: DMType, DMPlexCreate(), DMCreate(), DMSetType() 1664 M*/ 1665 1666 PETSC_EXTERN PetscErrorCode DMCreate_Plex(DM dm) 1667 { 1668 DM_Plex *mesh; 1669 PetscInt unit, d; 1670 PetscErrorCode ierr; 1671 1672 PetscFunctionBegin; 1673 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1674 ierr = PetscNewLog(dm,&mesh);CHKERRQ(ierr); 1675 dm->dim = 0; 1676 dm->data = mesh; 1677 1678 mesh->refct = 1; 1679 ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm), &mesh->coneSection);CHKERRQ(ierr); 1680 mesh->maxConeSize = 0; 1681 mesh->cones = NULL; 1682 mesh->coneOrientations = NULL; 1683 ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm), &mesh->supportSection);CHKERRQ(ierr); 1684 mesh->maxSupportSize = 0; 1685 mesh->supports = NULL; 1686 mesh->refinementUniform = PETSC_TRUE; 1687 mesh->refinementLimit = -1.0; 1688 1689 mesh->facesTmp = NULL; 1690 1691 mesh->tetgenOpts = NULL; 1692 mesh->triangleOpts = NULL; 1693 ierr = PetscPartitionerCreate(PetscObjectComm((PetscObject)dm), &mesh->partitioner);CHKERRQ(ierr); 1694 mesh->remeshBd = PETSC_FALSE; 1695 1696 mesh->subpointMap = NULL; 1697 1698 for (unit = 0; unit < NUM_PETSC_UNITS; ++unit) mesh->scale[unit] = 1.0; 1699 1700 mesh->regularRefinement = PETSC_FALSE; 1701 mesh->depthState = -1; 1702 mesh->globalVertexNumbers = NULL; 1703 mesh->globalCellNumbers = NULL; 1704 mesh->anchorSection = NULL; 1705 mesh->anchorIS = NULL; 1706 mesh->createanchors = NULL; 1707 mesh->computeanchormatrix = NULL; 1708 mesh->parentSection = NULL; 1709 mesh->parents = NULL; 1710 mesh->childIDs = NULL; 1711 mesh->childSection = NULL; 1712 mesh->children = NULL; 1713 mesh->referenceTree = NULL; 1714 mesh->getchildsymmetry = NULL; 1715 for (d = 0; d < 8; ++d) mesh->hybridPointMax[d] = PETSC_DETERMINE; 1716 mesh->vtkCellHeight = 0; 1717 mesh->useCone = PETSC_FALSE; 1718 mesh->useClosure = PETSC_TRUE; 1719 mesh->useAnchors = PETSC_FALSE; 1720 1721 mesh->maxProjectionHeight = 0; 1722 1723 mesh->printSetValues = PETSC_FALSE; 1724 mesh->printFEM = 0; 1725 mesh->printTol = 1.0e-10; 1726 1727 ierr = DMInitialize_Plex(dm);CHKERRQ(ierr); 1728 PetscFunctionReturn(0); 1729 } 1730 1731 /*@ 1732 DMPlexCreate - Creates a DMPlex object, which encapsulates an unstructured mesh, or CW complex, which can be expressed using a Hasse Diagram. 1733 1734 Collective on MPI_Comm 1735 1736 Input Parameter: 1737 . comm - The communicator for the DMPlex object 1738 1739 Output Parameter: 1740 . mesh - The DMPlex object 1741 1742 Level: beginner 1743 1744 .keywords: DMPlex, create 1745 @*/ 1746 PetscErrorCode DMPlexCreate(MPI_Comm comm, DM *mesh) 1747 { 1748 PetscErrorCode ierr; 1749 1750 PetscFunctionBegin; 1751 PetscValidPointer(mesh,2); 1752 ierr = DMCreate(comm, mesh);CHKERRQ(ierr); 1753 ierr = DMSetType(*mesh, DMPLEX);CHKERRQ(ierr); 1754 PetscFunctionReturn(0); 1755 } 1756 1757 /* 1758 This takes as input the common mesh generator output, a list of the vertices for each cell, but vertex numbers are global and an SF is built for them 1759 */ 1760 static PetscErrorCode DMPlexBuildFromCellList_Parallel_Private(DM dm, PetscInt numCells, PetscInt numVertices, PetscInt numCorners, const int cells[], PetscSF *sfVert) 1761 { 1762 PetscSF sfPoint; 1763 PetscLayout vLayout; 1764 PetscHashI vhash; 1765 PetscSFNode *remoteVerticesAdj, *vertexLocal, *vertexOwner, *remoteVertex; 1766 const PetscInt *vrange; 1767 PetscInt numVerticesAdj, off, *verticesAdj, numVerticesGhost = 0, *localVertex, *cone, c, p, v, g; 1768 PETSC_UNUSED PetscHashIIter ret, iter; 1769 PetscMPIInt rank, size; 1770 PetscErrorCode ierr; 1771 1772 PetscFunctionBegin; 1773 ierr = MPI_Comm_rank(PetscObjectComm((PetscObject) dm), &rank);CHKERRQ(ierr); 1774 ierr = MPI_Comm_size(PetscObjectComm((PetscObject) dm), &size);CHKERRQ(ierr); 1775 /* Partition vertices */ 1776 ierr = PetscLayoutCreate(PetscObjectComm((PetscObject) dm), &vLayout);CHKERRQ(ierr); 1777 ierr = PetscLayoutSetLocalSize(vLayout, numVertices);CHKERRQ(ierr); 1778 ierr = PetscLayoutSetBlockSize(vLayout, 1);CHKERRQ(ierr); 1779 ierr = PetscLayoutSetUp(vLayout);CHKERRQ(ierr); 1780 ierr = PetscLayoutGetRanges(vLayout, &vrange);CHKERRQ(ierr); 1781 /* Count vertices and map them to procs */ 1782 PetscHashICreate(vhash); 1783 for (c = 0; c < numCells; ++c) { 1784 for (p = 0; p < numCorners; ++p) { 1785 PetscHashIPut(vhash, cells[c*numCorners+p], ret, iter); 1786 } 1787 } 1788 PetscHashISize(vhash, numVerticesAdj); 1789 ierr = PetscMalloc1(numVerticesAdj, &verticesAdj);CHKERRQ(ierr); 1790 off = 0; ierr = PetscHashIGetKeys(vhash, &off, verticesAdj);CHKERRQ(ierr); 1791 if (off != numVerticesAdj) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid number of local vertices %D should be %D", off, numVerticesAdj); 1792 ierr = PetscSortInt(numVerticesAdj, verticesAdj);CHKERRQ(ierr); 1793 ierr = PetscMalloc1(numVerticesAdj, &remoteVerticesAdj);CHKERRQ(ierr); 1794 for (v = 0; v < numVerticesAdj; ++v) { 1795 const PetscInt gv = verticesAdj[v]; 1796 PetscInt vrank; 1797 1798 ierr = PetscFindInt(gv, size+1, vrange, &vrank);CHKERRQ(ierr); 1799 vrank = vrank < 0 ? -(vrank+2) : vrank; 1800 remoteVerticesAdj[v].index = gv - vrange[vrank]; 1801 remoteVerticesAdj[v].rank = vrank; 1802 } 1803 /* Create cones */ 1804 ierr = DMPlexSetChart(dm, 0, numCells+numVerticesAdj);CHKERRQ(ierr); 1805 for (c = 0; c < numCells; ++c) {ierr = DMPlexSetConeSize(dm, c, numCorners);CHKERRQ(ierr);} 1806 ierr = DMSetUp(dm);CHKERRQ(ierr); 1807 ierr = DMGetWorkArray(dm, numCorners, PETSC_INT, &cone);CHKERRQ(ierr); 1808 for (c = 0; c < numCells; ++c) { 1809 for (p = 0; p < numCorners; ++p) { 1810 const PetscInt gv = cells[c*numCorners+p]; 1811 PetscInt lv; 1812 1813 ierr = PetscFindInt(gv, numVerticesAdj, verticesAdj, &lv);CHKERRQ(ierr); 1814 if (lv < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Could not find global vertex %D in local connectivity", gv); 1815 cone[p] = lv+numCells; 1816 } 1817 ierr = DMPlexSetCone(dm, c, cone);CHKERRQ(ierr); 1818 } 1819 ierr = DMRestoreWorkArray(dm, numCorners, PETSC_INT, &cone);CHKERRQ(ierr); 1820 /* Create SF for vertices */ 1821 ierr = PetscSFCreate(PetscObjectComm((PetscObject)dm), sfVert);CHKERRQ(ierr); 1822 ierr = PetscObjectSetName((PetscObject) *sfVert, "Vertex Ownership SF");CHKERRQ(ierr); 1823 ierr = PetscSFSetFromOptions(*sfVert);CHKERRQ(ierr); 1824 ierr = PetscSFSetGraph(*sfVert, numVertices, numVerticesAdj, NULL, PETSC_OWN_POINTER, remoteVerticesAdj, PETSC_OWN_POINTER);CHKERRQ(ierr); 1825 ierr = PetscFree(verticesAdj);CHKERRQ(ierr); 1826 /* Build pointSF */ 1827 ierr = PetscMalloc2(numVerticesAdj, &vertexLocal, numVertices, &vertexOwner);CHKERRQ(ierr); 1828 for (v = 0; v < numVerticesAdj; ++v) {vertexLocal[v].index = v+numCells; vertexLocal[v].rank = rank;} 1829 for (v = 0; v < numVertices; ++v) {vertexOwner[v].index = -1; vertexOwner[v].rank = -1;} 1830 ierr = PetscSFReduceBegin(*sfVert, MPIU_2INT, vertexLocal, vertexOwner, MPI_MAXLOC);CHKERRQ(ierr); 1831 ierr = PetscSFReduceEnd(*sfVert, MPIU_2INT, vertexLocal, vertexOwner, MPI_MAXLOC);CHKERRQ(ierr); 1832 for (v = 0; v < numVertices; ++v) if (vertexOwner[v].rank < 0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Global vertex %d on rank %d was unclaimed", v + vrange[rank], rank); 1833 ierr = PetscSFBcastBegin(*sfVert, MPIU_2INT, vertexOwner, vertexLocal);CHKERRQ(ierr); 1834 ierr = PetscSFBcastEnd(*sfVert, MPIU_2INT, vertexOwner, vertexLocal);CHKERRQ(ierr); 1835 for (v = 0; v < numVerticesAdj; ++v) if (vertexLocal[v].rank != rank) ++numVerticesGhost; 1836 ierr = PetscMalloc1(numVerticesGhost, &localVertex);CHKERRQ(ierr); 1837 ierr = PetscMalloc1(numVerticesGhost, &remoteVertex);CHKERRQ(ierr); 1838 for (v = 0, g = 0; v < numVerticesAdj; ++v) { 1839 if (vertexLocal[v].rank != rank) { 1840 localVertex[g] = v+numCells; 1841 remoteVertex[g].index = vertexLocal[v].index; 1842 remoteVertex[g].rank = vertexLocal[v].rank; 1843 ++g; 1844 } 1845 } 1846 ierr = PetscFree2(vertexLocal, vertexOwner);CHKERRQ(ierr); 1847 if (g != numVerticesGhost) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid number of vertex ghosts %D should be %D", g, numVerticesGhost); 1848 ierr = DMGetPointSF(dm, &sfPoint);CHKERRQ(ierr); 1849 ierr = PetscObjectSetName((PetscObject) sfPoint, "point SF");CHKERRQ(ierr); 1850 ierr = PetscSFSetGraph(sfPoint, numCells+numVerticesAdj, numVerticesGhost, localVertex, PETSC_OWN_POINTER, remoteVertex, PETSC_OWN_POINTER);CHKERRQ(ierr); 1851 ierr = PetscLayoutDestroy(&vLayout);CHKERRQ(ierr); 1852 PetscHashIDestroy(vhash); 1853 /* Fill in the rest of the topology structure */ 1854 ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr); 1855 ierr = DMPlexStratify(dm);CHKERRQ(ierr); 1856 PetscFunctionReturn(0); 1857 } 1858 1859 /* 1860 This takes as input the coordinates for each owned vertex 1861 */ 1862 static PetscErrorCode DMPlexBuildCoordinates_Parallel_Private(DM dm, PetscInt spaceDim, PetscInt numCells, PetscInt numV, PetscSF sfVert, const PetscReal vertexCoords[]) 1863 { 1864 PetscSection coordSection; 1865 Vec coordinates; 1866 PetscScalar *coords; 1867 PetscInt numVertices, numVerticesAdj, coordSize, v; 1868 PetscErrorCode ierr; 1869 1870 PetscFunctionBegin; 1871 ierr = PetscSFGetGraph(sfVert, &numVertices, &numVerticesAdj, NULL, NULL);CHKERRQ(ierr); 1872 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 1873 ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr); 1874 ierr = PetscSectionSetFieldComponents(coordSection, 0, spaceDim);CHKERRQ(ierr); 1875 ierr = PetscSectionSetChart(coordSection, numCells, numCells + numVerticesAdj);CHKERRQ(ierr); 1876 for (v = numCells; v < numCells+numVerticesAdj; ++v) { 1877 ierr = PetscSectionSetDof(coordSection, v, spaceDim);CHKERRQ(ierr); 1878 ierr = PetscSectionSetFieldDof(coordSection, v, 0, spaceDim);CHKERRQ(ierr); 1879 } 1880 ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 1881 ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 1882 ierr = VecCreate(PetscObjectComm((PetscObject)dm), &coordinates);CHKERRQ(ierr); 1883 ierr = VecSetBlockSize(coordinates, spaceDim);CHKERRQ(ierr); 1884 ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr); 1885 ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 1886 ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr); 1887 ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 1888 { 1889 MPI_Datatype coordtype; 1890 1891 /* Need a temp buffer for coords if we have complex/single */ 1892 ierr = MPI_Type_contiguous(spaceDim, MPIU_SCALAR, &coordtype);CHKERRQ(ierr); 1893 ierr = MPI_Type_commit(&coordtype);CHKERRQ(ierr); 1894 #if defined(PETSC_USE_COMPLEX) 1895 { 1896 PetscScalar *svertexCoords; 1897 PetscInt i; 1898 ierr = PetscMalloc1(numV*spaceDim,&svertexCoords);CHKERRQ(ierr); 1899 for (i=0; i<numV*spaceDim; i++) svertexCoords[i] = vertexCoords[i]; 1900 ierr = PetscSFBcastBegin(sfVert, coordtype, svertexCoords, coords);CHKERRQ(ierr); 1901 ierr = PetscSFBcastEnd(sfVert, coordtype, svertexCoords, coords);CHKERRQ(ierr); 1902 ierr = PetscFree(svertexCoords);CHKERRQ(ierr); 1903 } 1904 #else 1905 ierr = PetscSFBcastBegin(sfVert, coordtype, vertexCoords, coords);CHKERRQ(ierr); 1906 ierr = PetscSFBcastEnd(sfVert, coordtype, vertexCoords, coords);CHKERRQ(ierr); 1907 #endif 1908 ierr = MPI_Type_free(&coordtype);CHKERRQ(ierr); 1909 } 1910 ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 1911 ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr); 1912 ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 1913 PetscFunctionReturn(0); 1914 } 1915 1916 /*@C 1917 DMPlexCreateFromCellListParallel - This takes as input common mesh generator output, a list of the vertices for each cell, and produces a DM 1918 1919 Input Parameters: 1920 + comm - The communicator 1921 . dim - The topological dimension of the mesh 1922 . numCells - The number of cells owned by this process 1923 . numVertices - The number of vertices owned by this process 1924 . numCorners - The number of vertices for each cell 1925 . interpolate - Flag indicating that intermediate mesh entities (faces, edges) should be created automatically 1926 . cells - An array of numCells*numCorners numbers, the global vertex numbers for each cell 1927 . spaceDim - The spatial dimension used for coordinates 1928 - vertexCoords - An array of numVertices*spaceDim numbers, the coordinates of each vertex 1929 1930 Output Parameter: 1931 + dm - The DM 1932 - vertexSF - Optional, SF describing complete vertex ownership 1933 1934 Note: Two triangles sharing a face 1935 $ 1936 $ 2 1937 $ / | \ 1938 $ / | \ 1939 $ / | \ 1940 $ 0 0 | 1 3 1941 $ \ | / 1942 $ \ | / 1943 $ \ | / 1944 $ 1 1945 would have input 1946 $ numCells = 2, numVertices = 4 1947 $ cells = [0 1 2 1 3 2] 1948 $ 1949 which would result in the DMPlex 1950 $ 1951 $ 4 1952 $ / | \ 1953 $ / | \ 1954 $ / | \ 1955 $ 2 0 | 1 5 1956 $ \ | / 1957 $ \ | / 1958 $ \ | / 1959 $ 3 1960 1961 Level: beginner 1962 1963 .seealso: DMPlexCreateFromCellList(), DMPlexCreateFromDAG(), DMPlexCreate() 1964 @*/ 1965 PetscErrorCode DMPlexCreateFromCellListParallel(MPI_Comm comm, PetscInt dim, PetscInt numCells, PetscInt numVertices, PetscInt numCorners, PetscBool interpolate, const int cells[], PetscInt spaceDim, const PetscReal vertexCoords[], PetscSF *vertexSF, DM *dm) 1966 { 1967 PetscSF sfVert; 1968 PetscErrorCode ierr; 1969 1970 PetscFunctionBegin; 1971 ierr = DMCreate(comm, dm);CHKERRQ(ierr); 1972 ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 1973 PetscValidLogicalCollectiveInt(*dm, dim, 2); 1974 PetscValidLogicalCollectiveInt(*dm, spaceDim, 8); 1975 ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr); 1976 ierr = DMPlexBuildFromCellList_Parallel_Private(*dm, numCells, numVertices, numCorners, cells, &sfVert);CHKERRQ(ierr); 1977 if (interpolate) { 1978 DM idm = NULL; 1979 1980 ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr); 1981 ierr = DMDestroy(dm);CHKERRQ(ierr); 1982 *dm = idm; 1983 } 1984 ierr = DMPlexBuildCoordinates_Parallel_Private(*dm, spaceDim, numCells, numVertices,sfVert, vertexCoords);CHKERRQ(ierr); 1985 if (vertexSF) *vertexSF = sfVert; 1986 else {ierr = PetscSFDestroy(&sfVert);CHKERRQ(ierr);} 1987 PetscFunctionReturn(0); 1988 } 1989 1990 /* 1991 This takes as input the common mesh generator output, a list of the vertices for each cell 1992 */ 1993 static PetscErrorCode DMPlexBuildFromCellList_Private(DM dm, PetscInt numCells, PetscInt numVertices, PetscInt numCorners, const int cells[]) 1994 { 1995 PetscInt *cone, c, p; 1996 PetscErrorCode ierr; 1997 1998 PetscFunctionBegin; 1999 ierr = DMPlexSetChart(dm, 0, numCells+numVertices);CHKERRQ(ierr); 2000 for (c = 0; c < numCells; ++c) { 2001 ierr = DMPlexSetConeSize(dm, c, numCorners);CHKERRQ(ierr); 2002 } 2003 ierr = DMSetUp(dm);CHKERRQ(ierr); 2004 ierr = DMGetWorkArray(dm, numCorners, PETSC_INT, &cone);CHKERRQ(ierr); 2005 for (c = 0; c < numCells; ++c) { 2006 for (p = 0; p < numCorners; ++p) { 2007 cone[p] = cells[c*numCorners+p]+numCells; 2008 } 2009 ierr = DMPlexSetCone(dm, c, cone);CHKERRQ(ierr); 2010 } 2011 ierr = DMRestoreWorkArray(dm, numCorners, PETSC_INT, &cone);CHKERRQ(ierr); 2012 ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr); 2013 ierr = DMPlexStratify(dm);CHKERRQ(ierr); 2014 PetscFunctionReturn(0); 2015 } 2016 2017 /* 2018 This takes as input the coordinates for each vertex 2019 */ 2020 static PetscErrorCode DMPlexBuildCoordinates_Private(DM dm, PetscInt spaceDim, PetscInt numCells, PetscInt numVertices, const double vertexCoords[]) 2021 { 2022 PetscSection coordSection; 2023 Vec coordinates; 2024 DM cdm; 2025 PetscScalar *coords; 2026 PetscInt v, d; 2027 PetscErrorCode ierr; 2028 2029 PetscFunctionBegin; 2030 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 2031 ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr); 2032 ierr = PetscSectionSetFieldComponents(coordSection, 0, spaceDim);CHKERRQ(ierr); 2033 ierr = PetscSectionSetChart(coordSection, numCells, numCells + numVertices);CHKERRQ(ierr); 2034 for (v = numCells; v < numCells+numVertices; ++v) { 2035 ierr = PetscSectionSetDof(coordSection, v, spaceDim);CHKERRQ(ierr); 2036 ierr = PetscSectionSetFieldDof(coordSection, v, 0, spaceDim);CHKERRQ(ierr); 2037 } 2038 ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 2039 2040 ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr); 2041 ierr = DMCreateLocalVector(cdm, &coordinates);CHKERRQ(ierr); 2042 ierr = VecSetBlockSize(coordinates, spaceDim);CHKERRQ(ierr); 2043 ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr); 2044 ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 2045 for (v = 0; v < numVertices; ++v) { 2046 for (d = 0; d < spaceDim; ++d) { 2047 coords[v*spaceDim+d] = vertexCoords[v*spaceDim+d]; 2048 } 2049 } 2050 ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 2051 ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr); 2052 ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 2053 PetscFunctionReturn(0); 2054 } 2055 2056 /*@C 2057 DMPlexCreateFromCellList - This takes as input common mesh generator output, a list of the vertices for each cell, and produces a DM 2058 2059 Input Parameters: 2060 + comm - The communicator 2061 . dim - The topological dimension of the mesh 2062 . numCells - The number of cells 2063 . numVertices - The number of vertices 2064 . numCorners - The number of vertices for each cell 2065 . interpolate - Flag indicating that intermediate mesh entities (faces, edges) should be created automatically 2066 . cells - An array of numCells*numCorners numbers, the vertices for each cell 2067 . spaceDim - The spatial dimension used for coordinates 2068 - vertexCoords - An array of numVertices*spaceDim numbers, the coordinates of each vertex 2069 2070 Output Parameter: 2071 . dm - The DM 2072 2073 Note: Two triangles sharing a face 2074 $ 2075 $ 2 2076 $ / | \ 2077 $ / | \ 2078 $ / | \ 2079 $ 0 0 | 1 3 2080 $ \ | / 2081 $ \ | / 2082 $ \ | / 2083 $ 1 2084 would have input 2085 $ numCells = 2, numVertices = 4 2086 $ cells = [0 1 2 1 3 2] 2087 $ 2088 which would result in the DMPlex 2089 $ 2090 $ 4 2091 $ / | \ 2092 $ / | \ 2093 $ / | \ 2094 $ 2 0 | 1 5 2095 $ \ | / 2096 $ \ | / 2097 $ \ | / 2098 $ 3 2099 2100 Level: beginner 2101 2102 .seealso: DMPlexCreateFromDAG(), DMPlexCreate() 2103 @*/ 2104 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) 2105 { 2106 PetscErrorCode ierr; 2107 2108 PetscFunctionBegin; 2109 ierr = DMCreate(comm, dm);CHKERRQ(ierr); 2110 ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 2111 ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr); 2112 ierr = DMPlexBuildFromCellList_Private(*dm, numCells, numVertices, numCorners, cells);CHKERRQ(ierr); 2113 if (interpolate) { 2114 DM idm = NULL; 2115 2116 ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr); 2117 ierr = DMDestroy(dm);CHKERRQ(ierr); 2118 *dm = idm; 2119 } 2120 ierr = DMPlexBuildCoordinates_Private(*dm, spaceDim, numCells, numVertices, vertexCoords);CHKERRQ(ierr); 2121 PetscFunctionReturn(0); 2122 } 2123 2124 /*@ 2125 DMPlexCreateFromDAG - This takes as input the adjacency-list representation of the Directed Acyclic Graph (Hasse Diagram) encoding a mesh, and produces a DM 2126 2127 Input Parameters: 2128 + dm - The empty DM object, usually from DMCreate() and DMSetDimension() 2129 . depth - The depth of the DAG 2130 . numPoints - The number of points at each depth 2131 . coneSize - The cone size of each point 2132 . cones - The concatenation of the cone points for each point, the cone list must be oriented correctly for each point 2133 . coneOrientations - The orientation of each cone point 2134 - vertexCoords - An array of numVertices*dim numbers, the coordinates of each vertex 2135 2136 Output Parameter: 2137 . dm - The DM 2138 2139 Note: Two triangles sharing a face would have input 2140 $ depth = 1, numPoints = [4 2], coneSize = [3 3 0 0 0 0] 2141 $ cones = [2 3 4 3 5 4], coneOrientations = [0 0 0 0 0 0] 2142 $ vertexCoords = [-1.0 0.0 0.0 -1.0 0.0 1.0 1.0 0.0] 2143 $ 2144 which would result in the DMPlex 2145 $ 2146 $ 4 2147 $ / | \ 2148 $ / | \ 2149 $ / | \ 2150 $ 2 0 | 1 5 2151 $ \ | / 2152 $ \ | / 2153 $ \ | / 2154 $ 3 2155 $ 2156 $ Notice that all points are numbered consecutively, unlikely DMPlexCreateFromCellList() 2157 2158 Level: advanced 2159 2160 .seealso: DMPlexCreateFromCellList(), DMPlexCreate() 2161 @*/ 2162 PetscErrorCode DMPlexCreateFromDAG(DM dm, PetscInt depth, const PetscInt numPoints[], const PetscInt coneSize[], const PetscInt cones[], const PetscInt coneOrientations[], const PetscScalar vertexCoords[]) 2163 { 2164 Vec coordinates; 2165 PetscSection coordSection; 2166 PetscScalar *coords; 2167 PetscInt coordSize, firstVertex = -1, pStart = 0, pEnd = 0, p, v, dim, dimEmbed, d, off; 2168 PetscErrorCode ierr; 2169 2170 PetscFunctionBegin; 2171 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 2172 ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr); 2173 if (dimEmbed < dim) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Embedding dimension %d cannot be less than intrinsic dimension %d",dimEmbed,dim); 2174 for (d = 0; d <= depth; ++d) pEnd += numPoints[d]; 2175 ierr = DMPlexSetChart(dm, pStart, pEnd);CHKERRQ(ierr); 2176 for (p = pStart; p < pEnd; ++p) { 2177 ierr = DMPlexSetConeSize(dm, p, coneSize[p-pStart]);CHKERRQ(ierr); 2178 if (firstVertex < 0 && !coneSize[p - pStart]) { 2179 firstVertex = p - pStart; 2180 } 2181 } 2182 if (firstVertex < 0 && numPoints[0]) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Expected %d vertices but could not find any", numPoints[0]); 2183 ierr = DMSetUp(dm);CHKERRQ(ierr); /* Allocate space for cones */ 2184 for (p = pStart, off = 0; p < pEnd; off += coneSize[p-pStart], ++p) { 2185 ierr = DMPlexSetCone(dm, p, &cones[off]);CHKERRQ(ierr); 2186 ierr = DMPlexSetConeOrientation(dm, p, &coneOrientations[off]);CHKERRQ(ierr); 2187 } 2188 ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr); 2189 ierr = DMPlexStratify(dm);CHKERRQ(ierr); 2190 /* Build coordinates */ 2191 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 2192 ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr); 2193 ierr = PetscSectionSetFieldComponents(coordSection, 0, dimEmbed);CHKERRQ(ierr); 2194 ierr = PetscSectionSetChart(coordSection, firstVertex, firstVertex+numPoints[0]);CHKERRQ(ierr); 2195 for (v = firstVertex; v < firstVertex+numPoints[0]; ++v) { 2196 ierr = PetscSectionSetDof(coordSection, v, dimEmbed);CHKERRQ(ierr); 2197 ierr = PetscSectionSetFieldDof(coordSection, v, 0, dimEmbed);CHKERRQ(ierr); 2198 } 2199 ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 2200 ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 2201 ierr = VecCreate(PETSC_COMM_SELF, &coordinates);CHKERRQ(ierr); 2202 ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr); 2203 ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 2204 ierr = VecSetBlockSize(coordinates, dimEmbed);CHKERRQ(ierr); 2205 ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr); 2206 ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 2207 for (v = 0; v < numPoints[0]; ++v) { 2208 PetscInt off; 2209 2210 ierr = PetscSectionGetOffset(coordSection, v+firstVertex, &off);CHKERRQ(ierr); 2211 for (d = 0; d < dimEmbed; ++d) { 2212 coords[off+d] = vertexCoords[v*dimEmbed+d]; 2213 } 2214 } 2215 ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 2216 ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr); 2217 ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 2218 PetscFunctionReturn(0); 2219 } 2220 2221 /*@C 2222 DMPlexCreateFromFile - This takes a filename and produces a DM 2223 2224 Input Parameters: 2225 + comm - The communicator 2226 . filename - A file name 2227 - interpolate - Flag to create intermediate mesh pieces (edges, faces) 2228 2229 Output Parameter: 2230 . dm - The DM 2231 2232 Level: beginner 2233 2234 .seealso: DMPlexCreateFromDAG(), DMPlexCreateFromCellList(), DMPlexCreate() 2235 @*/ 2236 PetscErrorCode DMPlexCreateFromFile(MPI_Comm comm, const char filename[], PetscBool interpolate, DM *dm) 2237 { 2238 const char *extGmsh = ".msh"; 2239 const char *extCGNS = ".cgns"; 2240 const char *extExodus = ".exo"; 2241 const char *extFluent = ".cas"; 2242 const char *extHDF5 = ".h5"; 2243 const char *extMed = ".med"; 2244 const char *extPLY = ".ply"; 2245 size_t len; 2246 PetscBool isGmsh, isCGNS, isExodus, isFluent, isHDF5, isMed, isPLY; 2247 PetscMPIInt rank; 2248 PetscErrorCode ierr; 2249 2250 PetscFunctionBegin; 2251 PetscValidPointer(filename, 2); 2252 PetscValidPointer(dm, 4); 2253 ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 2254 ierr = PetscStrlen(filename, &len);CHKERRQ(ierr); 2255 if (!len) SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Filename must be a valid path"); 2256 ierr = PetscStrncmp(&filename[PetscMax(0,len-4)], extGmsh, 4, &isGmsh);CHKERRQ(ierr); 2257 ierr = PetscStrncmp(&filename[PetscMax(0,len-5)], extCGNS, 5, &isCGNS);CHKERRQ(ierr); 2258 ierr = PetscStrncmp(&filename[PetscMax(0,len-4)], extExodus, 4, &isExodus);CHKERRQ(ierr); 2259 ierr = PetscStrncmp(&filename[PetscMax(0,len-4)], extFluent, 4, &isFluent);CHKERRQ(ierr); 2260 ierr = PetscStrncmp(&filename[PetscMax(0,len-3)], extHDF5, 3, &isHDF5);CHKERRQ(ierr); 2261 ierr = PetscStrncmp(&filename[PetscMax(0,len-4)], extMed, 4, &isMed);CHKERRQ(ierr); 2262 ierr = PetscStrncmp(&filename[PetscMax(0,len-4)], extPLY, 4, &isPLY);CHKERRQ(ierr); 2263 if (isGmsh) { 2264 ierr = DMPlexCreateGmshFromFile(comm, filename, interpolate, dm);CHKERRQ(ierr); 2265 } else if (isCGNS) { 2266 ierr = DMPlexCreateCGNSFromFile(comm, filename, interpolate, dm);CHKERRQ(ierr); 2267 } else if (isExodus) { 2268 ierr = DMPlexCreateExodusFromFile(comm, filename, interpolate, dm);CHKERRQ(ierr); 2269 } else if (isFluent) { 2270 ierr = DMPlexCreateFluentFromFile(comm, filename, interpolate, dm);CHKERRQ(ierr); 2271 } else if (isHDF5) { 2272 PetscViewer viewer; 2273 2274 ierr = PetscViewerCreate(comm, &viewer);CHKERRQ(ierr); 2275 ierr = PetscViewerSetType(viewer, PETSCVIEWERHDF5);CHKERRQ(ierr); 2276 ierr = PetscViewerFileSetMode(viewer, FILE_MODE_READ);CHKERRQ(ierr); 2277 ierr = PetscViewerFileSetName(viewer, filename);CHKERRQ(ierr); 2278 ierr = DMCreate(comm, dm);CHKERRQ(ierr); 2279 ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 2280 ierr = DMLoad(*dm, viewer);CHKERRQ(ierr); 2281 ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 2282 } else if (isMed) { 2283 ierr = DMPlexCreateMedFromFile(comm, filename, interpolate, dm);CHKERRQ(ierr); 2284 } else if (isPLY) { 2285 ierr = DMPlexCreatePLYFromFile(comm, filename, interpolate, dm);CHKERRQ(ierr); 2286 } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cannot load file %s: unrecognized extension", filename); 2287 PetscFunctionReturn(0); 2288 } 2289 2290 /*@ 2291 DMPlexCreateReferenceCell - Create a DMPLEX with the appropriate FEM reference cell 2292 2293 Collective on comm 2294 2295 Input Parameters: 2296 + comm - The communicator 2297 . dim - The spatial dimension 2298 - simplex - Flag for simplex, otherwise use a tensor-product cell 2299 2300 Output Parameter: 2301 . refdm - The reference cell 2302 2303 Level: intermediate 2304 2305 .keywords: reference cell 2306 .seealso: 2307 @*/ 2308 PetscErrorCode DMPlexCreateReferenceCell(MPI_Comm comm, PetscInt dim, PetscBool simplex, DM *refdm) 2309 { 2310 DM rdm; 2311 Vec coords; 2312 PetscErrorCode ierr; 2313 2314 PetscFunctionBeginUser; 2315 ierr = DMCreate(comm, &rdm);CHKERRQ(ierr); 2316 ierr = DMSetType(rdm, DMPLEX);CHKERRQ(ierr); 2317 ierr = DMSetDimension(rdm, dim);CHKERRQ(ierr); 2318 switch (dim) { 2319 case 0: 2320 { 2321 PetscInt numPoints[1] = {1}; 2322 PetscInt coneSize[1] = {0}; 2323 PetscInt cones[1] = {0}; 2324 PetscInt coneOrientations[1] = {0}; 2325 PetscScalar vertexCoords[1] = {0.0}; 2326 2327 ierr = DMPlexCreateFromDAG(rdm, 0, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr); 2328 } 2329 break; 2330 case 1: 2331 { 2332 PetscInt numPoints[2] = {2, 1}; 2333 PetscInt coneSize[3] = {2, 0, 0}; 2334 PetscInt cones[2] = {1, 2}; 2335 PetscInt coneOrientations[2] = {0, 0}; 2336 PetscScalar vertexCoords[2] = {-1.0, 1.0}; 2337 2338 ierr = DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr); 2339 } 2340 break; 2341 case 2: 2342 if (simplex) { 2343 PetscInt numPoints[2] = {3, 1}; 2344 PetscInt coneSize[4] = {3, 0, 0, 0}; 2345 PetscInt cones[3] = {1, 2, 3}; 2346 PetscInt coneOrientations[3] = {0, 0, 0}; 2347 PetscScalar vertexCoords[6] = {-1.0, -1.0, 1.0, -1.0, -1.0, 1.0}; 2348 2349 ierr = DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr); 2350 } else { 2351 PetscInt numPoints[2] = {4, 1}; 2352 PetscInt coneSize[5] = {4, 0, 0, 0, 0}; 2353 PetscInt cones[4] = {1, 2, 3, 4}; 2354 PetscInt coneOrientations[4] = {0, 0, 0, 0}; 2355 PetscScalar vertexCoords[8] = {-1.0, -1.0, 1.0, -1.0, 1.0, 1.0, -1.0, 1.0}; 2356 2357 ierr = DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr); 2358 } 2359 break; 2360 case 3: 2361 if (simplex) { 2362 PetscInt numPoints[2] = {4, 1}; 2363 PetscInt coneSize[5] = {4, 0, 0, 0, 0}; 2364 PetscInt cones[4] = {1, 3, 2, 4}; 2365 PetscInt coneOrientations[4] = {0, 0, 0, 0}; 2366 PetscScalar vertexCoords[12] = {-1.0, -1.0, -1.0, 1.0, -1.0, -1.0, -1.0, 1.0, -1.0, -1.0, -1.0, 1.0}; 2367 2368 ierr = DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr); 2369 } else { 2370 PetscInt numPoints[2] = {8, 1}; 2371 PetscInt coneSize[9] = {8, 0, 0, 0, 0, 0, 0, 0, 0}; 2372 PetscInt cones[8] = {1, 4, 3, 2, 5, 6, 7, 8}; 2373 PetscInt coneOrientations[8] = {0, 0, 0, 0, 0, 0, 0, 0}; 2374 PetscScalar vertexCoords[24] = {-1.0, -1.0, -1.0, 1.0, -1.0, -1.0, 1.0, 1.0, -1.0, -1.0, 1.0, -1.0, 2375 -1.0, -1.0, 1.0, 1.0, -1.0, 1.0, 1.0, 1.0, 1.0, -1.0, 1.0, 1.0}; 2376 2377 ierr = DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr); 2378 } 2379 break; 2380 default: 2381 SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Cannot create reference cell for dimension %d", dim); 2382 } 2383 *refdm = NULL; 2384 ierr = DMPlexInterpolate(rdm, refdm);CHKERRQ(ierr); 2385 if (rdm->coordinateDM) { 2386 DM ncdm; 2387 PetscSection cs; 2388 PetscInt pEnd = -1; 2389 2390 ierr = DMGetDefaultSection(rdm->coordinateDM, &cs);CHKERRQ(ierr); 2391 if (cs) {ierr = PetscSectionGetChart(cs, NULL, &pEnd);CHKERRQ(ierr);} 2392 if (pEnd >= 0) { 2393 ierr = DMClone(rdm->coordinateDM, &ncdm);CHKERRQ(ierr); 2394 ierr = DMSetDefaultSection(ncdm, cs);CHKERRQ(ierr); 2395 ierr = DMSetCoordinateDM(*refdm, ncdm);CHKERRQ(ierr); 2396 ierr = DMDestroy(&ncdm);CHKERRQ(ierr); 2397 } 2398 } 2399 ierr = DMGetCoordinatesLocal(rdm, &coords);CHKERRQ(ierr); 2400 if (coords) { 2401 ierr = DMSetCoordinatesLocal(*refdm, coords);CHKERRQ(ierr); 2402 } else { 2403 ierr = DMGetCoordinates(rdm, &coords);CHKERRQ(ierr); 2404 if (coords) {ierr = DMSetCoordinates(*refdm, coords);CHKERRQ(ierr);} 2405 } 2406 ierr = DMDestroy(&rdm);CHKERRQ(ierr); 2407 PetscFunctionReturn(0); 2408 } 2409