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 ierr = DMSetLabelValue(dm, "Face Sets", cone[1], markerRight);CHKERRQ(ierr); 229 } 230 } else if (vx == 0) { 231 ierr = DMSetLabelValue(dm, "marker", edge, markerLeft);CHKERRQ(ierr); 232 ierr = DMSetLabelValue(dm, "marker", cone[0], markerLeft);CHKERRQ(ierr); 233 if (ey == edges[1]-1) { 234 ierr = DMSetLabelValue(dm, "marker", cone[1], markerLeft);CHKERRQ(ierr); 235 ierr = DMSetLabelValue(dm, "Face Sets", cone[1], markerLeft);CHKERRQ(ierr); 236 } 237 } 238 } 239 } 240 for (vy = 0; vy <= edges[1]; vy++) { 241 for (ex = 0; ex < edges[0]; ex++) { 242 PetscInt edge = vy*edges[0] + ex; 243 PetscInt vertex = vy*(edges[0]+1) + ex + numEdges; 244 PetscInt cone[2]; 245 246 cone[0] = vertex; cone[1] = vertex+1; 247 ierr = DMPlexSetCone(dm, edge, cone);CHKERRQ(ierr); 248 if (vy == edges[1]) { 249 ierr = DMSetLabelValue(dm, "marker", edge, markerTop);CHKERRQ(ierr); 250 ierr = DMSetLabelValue(dm, "marker", cone[0], markerTop);CHKERRQ(ierr); 251 if (ex == edges[0]-1) { 252 ierr = DMSetLabelValue(dm, "marker", cone[1], markerTop);CHKERRQ(ierr); 253 ierr = DMSetLabelValue(dm, "Face Sets", cone[1], markerTop);CHKERRQ(ierr); 254 } 255 } else if (vy == 0) { 256 ierr = DMSetLabelValue(dm, "marker", edge, markerBottom);CHKERRQ(ierr); 257 ierr = DMSetLabelValue(dm, "marker", cone[0], markerBottom);CHKERRQ(ierr); 258 if (ex == edges[0]-1) { 259 ierr = DMSetLabelValue(dm, "marker", cone[1], markerBottom);CHKERRQ(ierr); 260 ierr = DMSetLabelValue(dm, "Face Sets", cone[1], markerBottom);CHKERRQ(ierr); 261 } 262 } 263 } 264 } 265 } 266 ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr); 267 ierr = DMPlexStratify(dm);CHKERRQ(ierr); 268 /* Build coordinates */ 269 ierr = DMSetCoordinateDim(dm, 2);CHKERRQ(ierr); 270 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 271 ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr); 272 ierr = PetscSectionSetChart(coordSection, numEdges, numEdges + numVertices);CHKERRQ(ierr); 273 ierr = PetscSectionSetFieldComponents(coordSection, 0, 2);CHKERRQ(ierr); 274 for (v = numEdges; v < numEdges+numVertices; ++v) { 275 ierr = PetscSectionSetDof(coordSection, v, 2);CHKERRQ(ierr); 276 ierr = PetscSectionSetFieldDof(coordSection, v, 0, 2);CHKERRQ(ierr); 277 } 278 ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 279 ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 280 ierr = VecCreate(PETSC_COMM_SELF, &coordinates);CHKERRQ(ierr); 281 ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr); 282 ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 283 ierr = VecSetBlockSize(coordinates, 2);CHKERRQ(ierr); 284 ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr); 285 ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 286 for (vy = 0; vy <= edges[1]; ++vy) { 287 for (vx = 0; vx <= edges[0]; ++vx) { 288 coords[(vy*(edges[0]+1)+vx)*2+0] = lower[0] + ((upper[0] - lower[0])/edges[0])*vx; 289 coords[(vy*(edges[0]+1)+vx)*2+1] = lower[1] + ((upper[1] - lower[1])/edges[1])*vy; 290 } 291 } 292 ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 293 ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr); 294 ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 295 PetscFunctionReturn(0); 296 } 297 298 /*@ 299 DMPlexCreateCubeBoundary - Creates a 2D mesh the is the boundary of a cubic lattice. 300 301 Collective on MPI_Comm 302 303 Input Parameters: 304 + comm - The communicator for the DM object 305 . lower - The lower left front corner coordinates 306 . upper - The upper right back corner coordinates 307 - edges - The number of cells in each direction 308 309 Output Parameter: 310 . dm - The DM object 311 312 Level: beginner 313 314 .keywords: DM, create 315 .seealso: DMPlexCreateBoxMesh(), DMPlexCreateSquareBoundary(), DMSetType(), DMCreate() 316 @*/ 317 PetscErrorCode DMPlexCreateCubeBoundary(DM dm, const PetscReal lower[], const PetscReal upper[], const PetscInt faces[]) 318 { 319 PetscInt vertices[3], numVertices; 320 PetscInt numFaces = 2*faces[0]*faces[1] + 2*faces[1]*faces[2] + 2*faces[0]*faces[2]; 321 Vec coordinates; 322 PetscSection coordSection; 323 PetscScalar *coords; 324 PetscInt coordSize; 325 PetscMPIInt rank; 326 PetscInt v, vx, vy, vz; 327 PetscInt voffset, iface=0, cone[4]; 328 PetscErrorCode ierr; 329 330 PetscFunctionBegin; 331 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"); 332 ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);CHKERRQ(ierr); 333 vertices[0] = faces[0]+1; vertices[1] = faces[1]+1; vertices[2] = faces[2]+1; 334 numVertices = vertices[0]*vertices[1]*vertices[2]; 335 if (!rank) { 336 PetscInt f; 337 338 ierr = DMPlexSetChart(dm, 0, numFaces+numVertices);CHKERRQ(ierr); 339 for (f = 0; f < numFaces; ++f) { 340 ierr = DMPlexSetConeSize(dm, f, 4);CHKERRQ(ierr); 341 } 342 ierr = DMSetUp(dm);CHKERRQ(ierr); /* Allocate space for cones */ 343 for (v = 0; v < numFaces+numVertices; ++v) { 344 ierr = DMSetLabelValue(dm, "marker", v, 1);CHKERRQ(ierr); 345 } 346 347 /* Side 0 (Top) */ 348 for (vy = 0; vy < faces[1]; vy++) { 349 for (vx = 0; vx < faces[0]; vx++) { 350 voffset = numFaces + vertices[0]*vertices[1]*(vertices[2]-1) + vy*vertices[0] + vx; 351 cone[0] = voffset; cone[1] = voffset+1; cone[2] = voffset+vertices[0]+1; cone[3] = voffset+vertices[0]; 352 ierr = DMPlexSetCone(dm, iface, cone);CHKERRQ(ierr); 353 iface++; 354 } 355 } 356 357 /* Side 1 (Bottom) */ 358 for (vy = 0; vy < faces[1]; vy++) { 359 for (vx = 0; vx < faces[0]; vx++) { 360 voffset = numFaces + vy*(faces[0]+1) + vx; 361 cone[0] = voffset+1; cone[1] = voffset; cone[2] = voffset+vertices[0]; cone[3] = voffset+vertices[0]+1; 362 ierr = DMPlexSetCone(dm, iface, cone);CHKERRQ(ierr); 363 iface++; 364 } 365 } 366 367 /* Side 2 (Front) */ 368 for (vz = 0; vz < faces[2]; vz++) { 369 for (vx = 0; vx < faces[0]; vx++) { 370 voffset = numFaces + vz*vertices[0]*vertices[1] + vx; 371 cone[0] = voffset; cone[1] = voffset+1; cone[2] = voffset+vertices[0]*vertices[1]+1; cone[3] = voffset+vertices[0]*vertices[1]; 372 ierr = DMPlexSetCone(dm, iface, cone);CHKERRQ(ierr); 373 iface++; 374 } 375 } 376 377 /* Side 3 (Back) */ 378 for (vz = 0; vz < faces[2]; vz++) { 379 for (vx = 0; vx < faces[0]; vx++) { 380 voffset = numFaces + vz*vertices[0]*vertices[1] + vertices[0]*(vertices[1]-1) + vx; 381 cone[0] = voffset+vertices[0]*vertices[1]; cone[1] = voffset+vertices[0]*vertices[1]+1; 382 cone[2] = voffset+1; cone[3] = voffset; 383 ierr = DMPlexSetCone(dm, iface, cone);CHKERRQ(ierr); 384 iface++; 385 } 386 } 387 388 /* Side 4 (Left) */ 389 for (vz = 0; vz < faces[2]; vz++) { 390 for (vy = 0; vy < faces[1]; vy++) { 391 voffset = numFaces + vz*vertices[0]*vertices[1] + vy*vertices[0]; 392 cone[0] = voffset; cone[1] = voffset+vertices[0]*vertices[1]; 393 cone[2] = voffset+vertices[0]*vertices[1]+vertices[0]; cone[3] = voffset+vertices[0]; 394 ierr = DMPlexSetCone(dm, iface, cone);CHKERRQ(ierr); 395 iface++; 396 } 397 } 398 399 /* Side 5 (Right) */ 400 for (vz = 0; vz < faces[2]; vz++) { 401 for (vy = 0; vy < faces[1]; vy++) { 402 voffset = numFaces + vz*vertices[0]*vertices[1] + vy*vertices[0] + faces[0]; 403 cone[0] = voffset+vertices[0]*vertices[1]; cone[1] = voffset; 404 cone[2] = voffset+vertices[0]; cone[3] = voffset+vertices[0]*vertices[1]+vertices[0]; 405 ierr = DMPlexSetCone(dm, iface, cone);CHKERRQ(ierr); 406 iface++; 407 } 408 } 409 } 410 ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr); 411 ierr = DMPlexStratify(dm);CHKERRQ(ierr); 412 /* Build coordinates */ 413 ierr = DMSetCoordinateDim(dm, 3);CHKERRQ(ierr); 414 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 415 ierr = PetscSectionSetChart(coordSection, numFaces, numFaces + numVertices);CHKERRQ(ierr); 416 for (v = numFaces; v < numFaces+numVertices; ++v) { 417 ierr = PetscSectionSetDof(coordSection, v, 3);CHKERRQ(ierr); 418 } 419 ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 420 ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 421 ierr = VecCreate(PETSC_COMM_SELF, &coordinates);CHKERRQ(ierr); 422 ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr); 423 ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 424 ierr = VecSetBlockSize(coordinates, 3);CHKERRQ(ierr); 425 ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr); 426 ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 427 for (vz = 0; vz <= faces[2]; ++vz) { 428 for (vy = 0; vy <= faces[1]; ++vy) { 429 for (vx = 0; vx <= faces[0]; ++vx) { 430 coords[((vz*(faces[1]+1)+vy)*(faces[0]+1)+vx)*3+0] = lower[0] + ((upper[0] - lower[0])/faces[0])*vx; 431 coords[((vz*(faces[1]+1)+vy)*(faces[0]+1)+vx)*3+1] = lower[1] + ((upper[1] - lower[1])/faces[1])*vy; 432 coords[((vz*(faces[1]+1)+vy)*(faces[0]+1)+vx)*3+2] = lower[2] + ((upper[2] - lower[2])/faces[2])*vz; 433 } 434 } 435 } 436 ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 437 ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr); 438 ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 439 PetscFunctionReturn(0); 440 } 441 442 static PetscErrorCode DMPlexCreateBoxMesh_Simplex_Internal(MPI_Comm comm, PetscInt dim, const PetscInt faces[], const PetscReal lower[], const PetscReal upper[], const DMBoundaryType periodicity[], PetscBool interpolate, DM *dm) 443 { 444 DM boundary; 445 PetscInt i; 446 PetscErrorCode ierr; 447 448 PetscFunctionBegin; 449 PetscValidPointer(dm, 4); 450 for (i = 0; i < dim; ++i) if (periodicity[i] != DM_BOUNDARY_NONE) SETERRQ(comm, PETSC_ERR_SUP, "Periodicity is not supported for simplex meshes"); 451 ierr = DMCreate(comm, &boundary);CHKERRQ(ierr); 452 PetscValidLogicalCollectiveInt(boundary,dim,2); 453 ierr = DMSetType(boundary, DMPLEX);CHKERRQ(ierr); 454 ierr = DMSetDimension(boundary, dim-1);CHKERRQ(ierr); 455 ierr = DMSetCoordinateDim(boundary, dim);CHKERRQ(ierr); 456 switch (dim) { 457 case 2: ierr = DMPlexCreateSquareBoundary(boundary, lower, upper, faces);CHKERRQ(ierr);break; 458 case 3: ierr = DMPlexCreateCubeBoundary(boundary, lower, upper, faces);CHKERRQ(ierr);break; 459 default: SETERRQ1(comm, PETSC_ERR_SUP, "Dimension not supported: %d", dim); 460 } 461 ierr = DMPlexGenerate(boundary, NULL, interpolate, dm);CHKERRQ(ierr); 462 ierr = DMDestroy(&boundary);CHKERRQ(ierr); 463 PetscFunctionReturn(0); 464 } 465 466 static PetscErrorCode DMPlexCreateCubeMesh_Internal(DM dm, const PetscReal lower[], const PetscReal upper[], const PetscInt edges[], DMBoundaryType bdX, DMBoundaryType bdY, DMBoundaryType bdZ) 467 { 468 DMLabel cutLabel = NULL; 469 PetscInt markerTop = 1, faceMarkerTop = 1; 470 PetscInt markerBottom = 1, faceMarkerBottom = 1; 471 PetscInt markerFront = 1, faceMarkerFront = 1; 472 PetscInt markerBack = 1, faceMarkerBack = 1; 473 PetscInt markerRight = 1, faceMarkerRight = 1; 474 PetscInt markerLeft = 1, faceMarkerLeft = 1; 475 PetscInt dim; 476 PetscBool markerSeparate = PETSC_FALSE, cutMarker = PETSC_FALSE; 477 PetscMPIInt rank; 478 PetscErrorCode ierr; 479 480 PetscFunctionBegin; 481 ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr); 482 ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);CHKERRQ(ierr); 483 ierr = DMCreateLabel(dm,"marker");CHKERRQ(ierr); 484 ierr = DMCreateLabel(dm,"Face Sets");CHKERRQ(ierr); 485 if (bdX == DM_BOUNDARY_PERIODIC || bdX == DM_BOUNDARY_TWIST || 486 bdY == DM_BOUNDARY_PERIODIC || bdY == DM_BOUNDARY_TWIST || 487 bdZ == DM_BOUNDARY_PERIODIC || bdZ == DM_BOUNDARY_TWIST) { 488 ierr = PetscOptionsGetBool(((PetscObject) dm)->options,((PetscObject) dm)->prefix, "-dm_plex_periodic_cut", &cutMarker, NULL);CHKERRQ(ierr); 489 if (cutMarker) {ierr = DMCreateLabel(dm, "periodic_cut");CHKERRQ(ierr); ierr = DMGetLabel(dm, "periodic_cut", &cutLabel);CHKERRQ(ierr);} 490 } 491 switch (dim) { 492 case 2: 493 faceMarkerTop = 3; 494 faceMarkerBottom = 1; 495 faceMarkerRight = 2; 496 faceMarkerLeft = 4; 497 break; 498 case 3: 499 faceMarkerBottom = 1; 500 faceMarkerTop = 2; 501 faceMarkerFront = 3; 502 faceMarkerBack = 4; 503 faceMarkerRight = 5; 504 faceMarkerLeft = 6; 505 break; 506 default: 507 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Dimension %d not supported",dim); 508 break; 509 } 510 ierr = PetscOptionsGetBool(((PetscObject) dm)->options,((PetscObject) dm)->prefix, "-dm_plex_separate_marker", &markerSeparate, NULL);CHKERRQ(ierr); 511 if (markerSeparate) { 512 markerBottom = faceMarkerBottom; 513 markerTop = faceMarkerTop; 514 markerFront = faceMarkerFront; 515 markerBack = faceMarkerBack; 516 markerRight = faceMarkerRight; 517 markerLeft = faceMarkerLeft; 518 } 519 { 520 const PetscInt numXEdges = !rank ? edges[0] : 0; 521 const PetscInt numYEdges = !rank ? edges[1] : 0; 522 const PetscInt numZEdges = !rank ? edges[2] : 0; 523 const PetscInt numXVertices = !rank ? (bdX == DM_BOUNDARY_PERIODIC || bdX == DM_BOUNDARY_TWIST ? edges[0] : edges[0]+1) : 0; 524 const PetscInt numYVertices = !rank ? (bdY == DM_BOUNDARY_PERIODIC || bdY == DM_BOUNDARY_TWIST ? edges[1] : edges[1]+1) : 0; 525 const PetscInt numZVertices = !rank ? (bdZ == DM_BOUNDARY_PERIODIC || bdZ == DM_BOUNDARY_TWIST ? edges[2] : edges[2]+1) : 0; 526 const PetscInt numCells = numXEdges*numYEdges*numZEdges; 527 const PetscInt numXFaces = numYEdges*numZEdges; 528 const PetscInt numYFaces = numXEdges*numZEdges; 529 const PetscInt numZFaces = numXEdges*numYEdges; 530 const PetscInt numTotXFaces = numXVertices*numXFaces; 531 const PetscInt numTotYFaces = numYVertices*numYFaces; 532 const PetscInt numTotZFaces = numZVertices*numZFaces; 533 const PetscInt numFaces = numTotXFaces + numTotYFaces + numTotZFaces; 534 const PetscInt numTotXEdges = numXEdges*numYVertices*numZVertices; 535 const PetscInt numTotYEdges = numYEdges*numXVertices*numZVertices; 536 const PetscInt numTotZEdges = numZEdges*numXVertices*numYVertices; 537 const PetscInt numVertices = numXVertices*numYVertices*numZVertices; 538 const PetscInt numEdges = numTotXEdges + numTotYEdges + numTotZEdges; 539 const PetscInt firstVertex = (dim == 2) ? numFaces : numCells; 540 const PetscInt firstXFace = (dim == 2) ? 0 : numCells + numVertices; 541 const PetscInt firstYFace = firstXFace + numTotXFaces; 542 const PetscInt firstZFace = firstYFace + numTotYFaces; 543 const PetscInt firstXEdge = numCells + numFaces + numVertices; 544 const PetscInt firstYEdge = firstXEdge + numTotXEdges; 545 const PetscInt firstZEdge = firstYEdge + numTotYEdges; 546 Vec coordinates; 547 PetscSection coordSection; 548 PetscScalar *coords; 549 PetscInt coordSize; 550 PetscInt v, vx, vy, vz; 551 PetscInt c, f, fx, fy, fz, e, ex, ey, ez; 552 553 ierr = DMPlexSetChart(dm, 0, numCells+numFaces+numEdges+numVertices);CHKERRQ(ierr); 554 for (c = 0; c < numCells; c++) { 555 ierr = DMPlexSetConeSize(dm, c, 6);CHKERRQ(ierr); 556 } 557 for (f = firstXFace; f < firstXFace+numFaces; ++f) { 558 ierr = DMPlexSetConeSize(dm, f, 4);CHKERRQ(ierr); 559 } 560 for (e = firstXEdge; e < firstXEdge+numEdges; ++e) { 561 ierr = DMPlexSetConeSize(dm, e, 2);CHKERRQ(ierr); 562 } 563 ierr = DMSetUp(dm);CHKERRQ(ierr); /* Allocate space for cones */ 564 /* Build cells */ 565 for (fz = 0; fz < numZEdges; ++fz) { 566 for (fy = 0; fy < numYEdges; ++fy) { 567 for (fx = 0; fx < numXEdges; ++fx) { 568 PetscInt cell = (fz*numYEdges + fy)*numXEdges + fx; 569 PetscInt faceB = firstZFace + (fy*numXEdges+fx)*numZVertices + fz; 570 PetscInt faceT = firstZFace + (fy*numXEdges+fx)*numZVertices + ((fz+1)%numZVertices); 571 PetscInt faceF = firstYFace + (fz*numXEdges+fx)*numYVertices + fy; 572 PetscInt faceK = firstYFace + (fz*numXEdges+fx)*numYVertices + ((fy+1)%numYVertices); 573 PetscInt faceL = firstXFace + (fz*numYEdges+fy)*numXVertices + fx; 574 PetscInt faceR = firstXFace + (fz*numYEdges+fy)*numXVertices + ((fx+1)%numXVertices); 575 /* B, T, F, K, R, L */ 576 PetscInt ornt[6] = {-4, 0, 0, -1, 0, -4}; /* ??? */ 577 PetscInt cone[6]; 578 579 /* no boundary twisting in 3D */ 580 cone[0] = faceB; cone[1] = faceT; cone[2] = faceF; cone[3] = faceK; cone[4] = faceR; cone[5] = faceL; 581 ierr = DMPlexSetCone(dm, cell, cone);CHKERRQ(ierr); 582 ierr = DMPlexSetConeOrientation(dm, cell, ornt);CHKERRQ(ierr); 583 if (bdX != DM_BOUNDARY_NONE && fx == numXEdges-1 && cutLabel) {ierr = DMLabelSetValue(cutLabel, cell, 2);CHKERRQ(ierr);} 584 if (bdY != DM_BOUNDARY_NONE && fy == numYEdges-1 && cutLabel) {ierr = DMLabelSetValue(cutLabel, cell, 2);CHKERRQ(ierr);} 585 if (bdZ != DM_BOUNDARY_NONE && fz == numZEdges-1 && cutLabel) {ierr = DMLabelSetValue(cutLabel, cell, 2);CHKERRQ(ierr);} 586 } 587 } 588 } 589 /* Build x faces */ 590 for (fz = 0; fz < numZEdges; ++fz) { 591 for (fy = 0; fy < numYEdges; ++fy) { 592 for (fx = 0; fx < numXVertices; ++fx) { 593 PetscInt face = firstXFace + (fz*numYEdges+fy)*numXVertices + fx; 594 PetscInt edgeL = firstZEdge + ( fy* numXVertices+fx)*numZEdges + fz; 595 PetscInt edgeR = firstZEdge + (((fy+1)%numYVertices)*numXVertices+fx)*numZEdges + fz; 596 PetscInt edgeB = firstYEdge + ( fz* numXVertices+fx)*numYEdges + fy; 597 PetscInt edgeT = firstYEdge + (((fz+1)%numZVertices)*numXVertices+fx)*numYEdges + fy; 598 PetscInt ornt[4] = {0, 0, -2, -2}; 599 PetscInt cone[4]; 600 601 if (dim == 3) { 602 /* markers */ 603 if (bdX != DM_BOUNDARY_PERIODIC) { 604 if (fx == numXVertices-1) { 605 ierr = DMSetLabelValue(dm, "Face Sets", face, faceMarkerRight);CHKERRQ(ierr); 606 ierr = DMSetLabelValue(dm, "marker", face, markerRight);CHKERRQ(ierr); 607 } 608 else if (fx == 0) { 609 ierr = DMSetLabelValue(dm, "Face Sets", face, faceMarkerLeft);CHKERRQ(ierr); 610 ierr = DMSetLabelValue(dm, "marker", face, markerLeft);CHKERRQ(ierr); 611 } 612 } 613 } 614 cone[0] = edgeB; cone[1] = edgeR; cone[2] = edgeT; cone[3] = edgeL; 615 ierr = DMPlexSetCone(dm, face, cone);CHKERRQ(ierr); 616 ierr = DMPlexSetConeOrientation(dm, face, ornt);CHKERRQ(ierr); 617 } 618 } 619 } 620 /* Build y faces */ 621 for (fz = 0; fz < numZEdges; ++fz) { 622 for (fx = 0; fx < numXEdges; ++fx) { 623 for (fy = 0; fy < numYVertices; ++fy) { 624 PetscInt face = firstYFace + (fz*numXEdges+fx)*numYVertices + fy; 625 PetscInt edgeL = firstZEdge + (fy*numXVertices+ fx )*numZEdges + fz; 626 PetscInt edgeR = firstZEdge + (fy*numXVertices+((fx+1)%numXVertices))*numZEdges + fz; 627 PetscInt edgeB = firstXEdge + ( fz *numYVertices+fy)*numXEdges + fx; 628 PetscInt edgeT = firstXEdge + (((fz+1)%numZVertices)*numYVertices+fy)*numXEdges + fx; 629 PetscInt ornt[4] = {0, 0, -2, -2}; 630 PetscInt cone[4]; 631 632 if (dim == 3) { 633 /* markers */ 634 if (bdY != DM_BOUNDARY_PERIODIC) { 635 if (fy == numYVertices-1) { 636 ierr = DMSetLabelValue(dm, "Face Sets", face, faceMarkerBack);CHKERRQ(ierr); 637 ierr = DMSetLabelValue(dm, "marker", face, markerBack);CHKERRQ(ierr); 638 } 639 else if (fy == 0) { 640 ierr = DMSetLabelValue(dm, "Face Sets", face, faceMarkerFront);CHKERRQ(ierr); 641 ierr = DMSetLabelValue(dm, "marker", face, markerFront);CHKERRQ(ierr); 642 } 643 } 644 } 645 cone[0] = edgeB; cone[1] = edgeR; cone[2] = edgeT; cone[3] = edgeL; 646 ierr = DMPlexSetCone(dm, face, cone);CHKERRQ(ierr); 647 ierr = DMPlexSetConeOrientation(dm, face, ornt);CHKERRQ(ierr); 648 } 649 } 650 } 651 /* Build z faces */ 652 for (fy = 0; fy < numYEdges; ++fy) { 653 for (fx = 0; fx < numXEdges; ++fx) { 654 for (fz = 0; fz < numZVertices; fz++) { 655 PetscInt face = firstZFace + (fy*numXEdges+fx)*numZVertices + fz; 656 PetscInt edgeL = firstYEdge + (fz*numXVertices+ fx )*numYEdges + fy; 657 PetscInt edgeR = firstYEdge + (fz*numXVertices+((fx+1)%numXVertices))*numYEdges + fy; 658 PetscInt edgeB = firstXEdge + (fz*numYVertices+ fy )*numXEdges + fx; 659 PetscInt edgeT = firstXEdge + (fz*numYVertices+((fy+1)%numYVertices))*numXEdges + fx; 660 PetscInt ornt[4] = {0, 0, -2, -2}; 661 PetscInt cone[4]; 662 663 if (dim == 2) { 664 if (bdX == DM_BOUNDARY_TWIST && fx == numXEdges-1) {edgeR += numYEdges-1-2*fy; ornt[1] = -2;} 665 if (bdY == DM_BOUNDARY_TWIST && fy == numYEdges-1) {edgeT += numXEdges-1-2*fx; ornt[2] = 0;} 666 if (bdX != DM_BOUNDARY_NONE && fx == numXEdges-1 && cutLabel) {ierr = DMLabelSetValue(cutLabel, face, 2);CHKERRQ(ierr);} 667 if (bdY != DM_BOUNDARY_NONE && fy == numYEdges-1 && cutLabel) {ierr = DMLabelSetValue(cutLabel, face, 2);CHKERRQ(ierr);} 668 } else { 669 /* markers */ 670 if (bdZ != DM_BOUNDARY_PERIODIC) { 671 if (fz == numZVertices-1) { 672 ierr = DMSetLabelValue(dm, "Face Sets", face, faceMarkerTop);CHKERRQ(ierr); 673 ierr = DMSetLabelValue(dm, "marker", face, markerTop);CHKERRQ(ierr); 674 } 675 else if (fz == 0) { 676 ierr = DMSetLabelValue(dm, "Face Sets", face, faceMarkerBottom);CHKERRQ(ierr); 677 ierr = DMSetLabelValue(dm, "marker", face, markerBottom);CHKERRQ(ierr); 678 } 679 } 680 } 681 cone[0] = edgeB; cone[1] = edgeR; cone[2] = edgeT; cone[3] = edgeL; 682 ierr = DMPlexSetCone(dm, face, cone);CHKERRQ(ierr); 683 ierr = DMPlexSetConeOrientation(dm, face, ornt);CHKERRQ(ierr); 684 } 685 } 686 } 687 /* Build Z edges*/ 688 for (vy = 0; vy < numYVertices; vy++) { 689 for (vx = 0; vx < numXVertices; vx++) { 690 for (ez = 0; ez < numZEdges; ez++) { 691 const PetscInt edge = firstZEdge + (vy*numXVertices+vx)*numZEdges + ez; 692 const PetscInt vertexB = firstVertex + ( ez *numYVertices+vy)*numXVertices + vx; 693 const PetscInt vertexT = firstVertex + (((ez+1)%numZVertices)*numYVertices+vy)*numXVertices + vx; 694 PetscInt cone[2]; 695 696 if (dim == 3) { 697 if (bdX != DM_BOUNDARY_PERIODIC) { 698 if (vx == numXVertices-1) { 699 ierr = DMSetLabelValue(dm, "marker", edge, markerRight);CHKERRQ(ierr); 700 } 701 else if (vx == 0) { 702 ierr = DMSetLabelValue(dm, "marker", edge, markerLeft);CHKERRQ(ierr); 703 } 704 } 705 if (bdY != DM_BOUNDARY_PERIODIC) { 706 if (vy == numYVertices-1) { 707 ierr = DMSetLabelValue(dm, "marker", edge, markerBack);CHKERRQ(ierr); 708 } 709 else if (vy == 0) { 710 ierr = DMSetLabelValue(dm, "marker", edge, markerFront);CHKERRQ(ierr); 711 } 712 } 713 } 714 cone[0] = vertexB; cone[1] = vertexT; 715 ierr = DMPlexSetCone(dm, edge, cone);CHKERRQ(ierr); 716 } 717 } 718 } 719 /* Build Y edges*/ 720 for (vz = 0; vz < numZVertices; vz++) { 721 for (vx = 0; vx < numXVertices; vx++) { 722 for (ey = 0; ey < numYEdges; ey++) { 723 const PetscInt nextv = (dim == 2 && bdY == DM_BOUNDARY_TWIST && ey == numYEdges-1) ? (numXVertices-vx-1) : (vz*numYVertices+((ey+1)%numYVertices))*numXVertices + vx; 724 const PetscInt edge = firstYEdge + (vz*numXVertices+vx)*numYEdges + ey; 725 const PetscInt vertexF = firstVertex + (vz*numYVertices+ey)*numXVertices + vx; 726 const PetscInt vertexK = firstVertex + nextv; 727 PetscInt cone[2]; 728 729 cone[0] = vertexF; cone[1] = vertexK; 730 ierr = DMPlexSetCone(dm, edge, cone);CHKERRQ(ierr); 731 if (dim == 2) { 732 if ((bdX != DM_BOUNDARY_PERIODIC) && (bdX != DM_BOUNDARY_TWIST)) { 733 if (vx == numXVertices-1) { 734 ierr = DMSetLabelValue(dm, "Face Sets", edge, faceMarkerRight);CHKERRQ(ierr); 735 ierr = DMSetLabelValue(dm, "marker", edge, markerRight);CHKERRQ(ierr); 736 ierr = DMSetLabelValue(dm, "marker", cone[0], markerRight);CHKERRQ(ierr); 737 if (ey == numYEdges-1) { 738 ierr = DMSetLabelValue(dm, "marker", cone[1], markerRight);CHKERRQ(ierr); 739 } 740 } else if (vx == 0) { 741 ierr = DMSetLabelValue(dm, "Face Sets", edge, faceMarkerLeft);CHKERRQ(ierr); 742 ierr = DMSetLabelValue(dm, "marker", edge, markerLeft);CHKERRQ(ierr); 743 ierr = DMSetLabelValue(dm, "marker", cone[0], markerLeft);CHKERRQ(ierr); 744 if (ey == numYEdges-1) { 745 ierr = DMSetLabelValue(dm, "marker", cone[1], markerLeft);CHKERRQ(ierr); 746 } 747 } 748 } else { 749 if (vx == 0 && cutMarker) { 750 ierr = DMLabelSetValue(cutLabel, edge, 1);CHKERRQ(ierr); 751 ierr = DMLabelSetValue(cutLabel, cone[0], 1);CHKERRQ(ierr); 752 if (ey == numYEdges-1) { 753 ierr = DMLabelSetValue(cutLabel, cone[1], 1);CHKERRQ(ierr); 754 } 755 } 756 } 757 } else { 758 if (bdX != DM_BOUNDARY_PERIODIC) { 759 if (vx == numXVertices-1) { 760 ierr = DMSetLabelValue(dm, "marker", edge, markerRight);CHKERRQ(ierr); 761 } else if (vx == 0) { 762 ierr = DMSetLabelValue(dm, "marker", edge, markerLeft);CHKERRQ(ierr); 763 } 764 } 765 if (bdZ != DM_BOUNDARY_PERIODIC) { 766 if (vz == numZVertices-1) { 767 ierr = DMSetLabelValue(dm, "marker", edge, markerTop);CHKERRQ(ierr); 768 } else if (vz == 0) { 769 ierr = DMSetLabelValue(dm, "marker", edge, markerBottom);CHKERRQ(ierr); 770 } 771 } 772 } 773 } 774 } 775 } 776 /* Build X edges*/ 777 for (vz = 0; vz < numZVertices; vz++) { 778 for (vy = 0; vy < numYVertices; vy++) { 779 for (ex = 0; ex < numXEdges; ex++) { 780 const PetscInt nextv = (dim == 2 && bdX == DM_BOUNDARY_TWIST && ex == numXEdges-1) ? (numYVertices-vy-1)*numXVertices : (vz*numYVertices+vy)*numXVertices + (ex+1)%numXVertices; 781 const PetscInt edge = firstXEdge + (vz*numYVertices+vy)*numXEdges + ex; 782 const PetscInt vertexL = firstVertex + (vz*numYVertices+vy)*numXVertices + ex; 783 const PetscInt vertexR = firstVertex + nextv; 784 PetscInt cone[2]; 785 786 cone[0] = vertexL; cone[1] = vertexR; 787 ierr = DMPlexSetCone(dm, edge, cone);CHKERRQ(ierr); 788 if (dim == 2) { 789 if ((bdY != DM_BOUNDARY_PERIODIC) && (bdY != DM_BOUNDARY_TWIST)) { 790 if (vy == numYVertices-1) { 791 ierr = DMSetLabelValue(dm, "Face Sets", edge, faceMarkerTop);CHKERRQ(ierr); 792 ierr = DMSetLabelValue(dm, "marker", edge, markerTop);CHKERRQ(ierr); 793 ierr = DMSetLabelValue(dm, "marker", cone[0], markerTop);CHKERRQ(ierr); 794 if (ex == numXEdges-1) { 795 ierr = DMSetLabelValue(dm, "marker", cone[1], markerTop);CHKERRQ(ierr); 796 } 797 } else if (vy == 0) { 798 ierr = DMSetLabelValue(dm, "Face Sets", edge, faceMarkerBottom);CHKERRQ(ierr); 799 ierr = DMSetLabelValue(dm, "marker", edge, markerBottom);CHKERRQ(ierr); 800 ierr = DMSetLabelValue(dm, "marker", cone[0], markerBottom);CHKERRQ(ierr); 801 if (ex == numXEdges-1) { 802 ierr = DMSetLabelValue(dm, "marker", cone[1], markerBottom);CHKERRQ(ierr); 803 } 804 } 805 } else { 806 if (vy == 0 && cutMarker) { 807 ierr = DMLabelSetValue(cutLabel, edge, 1);CHKERRQ(ierr); 808 ierr = DMLabelSetValue(cutLabel, cone[0], 1);CHKERRQ(ierr); 809 if (ex == numXEdges-1) { 810 ierr = DMLabelSetValue(cutLabel, cone[1], 1);CHKERRQ(ierr); 811 } 812 } 813 } 814 } else { 815 if (bdY != DM_BOUNDARY_PERIODIC) { 816 if (vy == numYVertices-1) { 817 ierr = DMSetLabelValue(dm, "marker", edge, markerBack);CHKERRQ(ierr); 818 } 819 else if (vy == 0) { 820 ierr = DMSetLabelValue(dm, "marker", edge, markerFront);CHKERRQ(ierr); 821 } 822 } 823 if (bdZ != DM_BOUNDARY_PERIODIC) { 824 if (vz == numZVertices-1) { 825 ierr = DMSetLabelValue(dm, "marker", edge, markerTop);CHKERRQ(ierr); 826 } 827 else if (vz == 0) { 828 ierr = DMSetLabelValue(dm, "marker", edge, markerBottom);CHKERRQ(ierr); 829 } 830 } 831 } 832 } 833 } 834 } 835 ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr); 836 ierr = DMPlexStratify(dm);CHKERRQ(ierr); 837 /* Build coordinates */ 838 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 839 ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr); 840 ierr = PetscSectionSetFieldComponents(coordSection, 0, dim);CHKERRQ(ierr); 841 ierr = PetscSectionSetChart(coordSection, firstVertex, firstVertex+numVertices);CHKERRQ(ierr); 842 for (v = firstVertex; v < firstVertex+numVertices; ++v) { 843 ierr = PetscSectionSetDof(coordSection, v, dim);CHKERRQ(ierr); 844 ierr = PetscSectionSetFieldDof(coordSection, v, 0, dim);CHKERRQ(ierr); 845 } 846 ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 847 ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 848 ierr = VecCreate(PETSC_COMM_SELF, &coordinates);CHKERRQ(ierr); 849 ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr); 850 ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 851 ierr = VecSetBlockSize(coordinates, dim);CHKERRQ(ierr); 852 ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr); 853 ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 854 for (vz = 0; vz < numZVertices; ++vz) { 855 for (vy = 0; vy < numYVertices; ++vy) { 856 for (vx = 0; vx < numXVertices; ++vx) { 857 coords[((vz*numYVertices+vy)*numXVertices+vx)*dim+0] = lower[0] + ((upper[0] - lower[0])/numXEdges)*vx; 858 coords[((vz*numYVertices+vy)*numXVertices+vx)*dim+1] = lower[1] + ((upper[1] - lower[1])/numYEdges)*vy; 859 if (dim == 3) { 860 coords[((vz*numYVertices+vy)*numXVertices+vx)*dim+2] = lower[2] + ((upper[2] - lower[2])/numZEdges)*vz; 861 } 862 } 863 } 864 } 865 ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 866 ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr); 867 ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 868 } 869 PetscFunctionReturn(0); 870 } 871 872 static PetscErrorCode DMPlexCreateBoxMesh_Tensor_Internal(MPI_Comm comm, PetscInt dim, const PetscInt faces[], const PetscReal lower[], const PetscReal upper[], const DMBoundaryType periodicity[], PetscBool interpolate, DM *dm) 873 { 874 PetscInt i; 875 PetscErrorCode ierr; 876 877 PetscFunctionBegin; 878 PetscValidPointer(dm, 7); 879 ierr = DMCreate(comm, dm);CHKERRQ(ierr); 880 PetscValidLogicalCollectiveInt(*dm,dim,2); 881 ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 882 ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr); 883 switch (dim) { 884 case 2: {ierr = DMPlexCreateCubeMesh_Internal(*dm, lower, upper, faces, periodicity[0], periodicity[1], DM_BOUNDARY_NONE);CHKERRQ(ierr);break;} 885 case 3: {ierr = DMPlexCreateCubeMesh_Internal(*dm, lower, upper, faces, periodicity[0], periodicity[1], periodicity[2]);CHKERRQ(ierr);break;} 886 default: SETERRQ1(comm, PETSC_ERR_SUP, "Dimension not supported: %d", dim); 887 } 888 if (periodicity[0] == DM_BOUNDARY_PERIODIC || periodicity[0] == DM_BOUNDARY_TWIST || 889 periodicity[1] == DM_BOUNDARY_PERIODIC || periodicity[1] == DM_BOUNDARY_TWIST || 890 (dim > 2 && (periodicity[2] == DM_BOUNDARY_PERIODIC || periodicity[2] == DM_BOUNDARY_TWIST))) { 891 PetscReal L[3]; 892 PetscReal maxCell[3]; 893 894 for (i = 0; i < dim; i++) { 895 L[i] = upper[i] - lower[i]; 896 maxCell[i] = 1.1 * (L[i] / PetscMax(1,faces[i])); 897 } 898 ierr = DMSetPeriodicity(*dm,PETSC_TRUE,maxCell,L,periodicity);CHKERRQ(ierr); 899 } 900 if (!interpolate) { 901 DM udm; 902 903 ierr = DMPlexUninterpolate(*dm, &udm);CHKERRQ(ierr); 904 ierr = DMDestroy(dm);CHKERRQ(ierr); 905 *dm = udm; 906 } 907 PetscFunctionReturn(0); 908 } 909 910 /*@C 911 DMPlexCreateBoxMesh - Creates a mesh on the tensor product of unit intervals (box) using simplices or tensor cells (hexahedra). 912 913 Collective on MPI_Comm 914 915 Input Parameters: 916 + comm - The communicator for the DM object 917 . dim - The spatial dimension 918 . simplex - PETSC_TRUE for simplices, PETSC_FALSE for tensor cells 919 . faces - Number of faces per dimension, or NULL for (2, 2) in 2D and (1, 1, 1) in 3D 920 . lower - The lower left corner, or NULL for (0, 0, 0) 921 . upper - The upper right corner, or NULL for (1, 1, 1) 922 . periodicity - The boundary type for the X,Y,Z direction, or NULL for DM_BOUDNARY_NONE 923 - interpolate - Flag to create intermediate mesh pieces (edges, faces) 924 925 Output Parameter: 926 . dm - The DM object 927 928 Note: Here is the numbering returned for 2 faces in each direction for tensor cells: 929 $ 10---17---11---18----12 930 $ | | | 931 $ | | | 932 $ 20 2 22 3 24 933 $ | | | 934 $ | | | 935 $ 7---15----8---16----9 936 $ | | | 937 $ | | | 938 $ 19 0 21 1 23 939 $ | | | 940 $ | | | 941 $ 4---13----5---14----6 942 943 and for simplicial cells 944 945 $ 14----8---15----9----16 946 $ |\ 5 |\ 7 | 947 $ | \ | \ | 948 $ 13 2 14 3 15 949 $ | 4 \ | 6 \ | 950 $ | \ | \ | 951 $ 11----6---12----7----13 952 $ |\ |\ | 953 $ | \ 1 | \ 3 | 954 $ 10 0 11 1 12 955 $ | 0 \ | 2 \ | 956 $ | \ | \ | 957 $ 8----4----9----5----10 958 959 Level: beginner 960 961 .keywords: DM, create 962 .seealso: DMPlexCreateFromFile(), DMPlexCreateHexCylinderMesh(), DMSetType(), DMCreate() 963 @*/ 964 PetscErrorCode DMPlexCreateBoxMesh(MPI_Comm comm, PetscInt dim, PetscBool simplex, const PetscInt faces[], const PetscReal lower[], const PetscReal upper[], const DMBoundaryType periodicity[], PetscBool interpolate, DM *dm) 965 { 966 PetscReal low[3]; 967 PetscReal upp[3]; 968 PetscInt fac[3]; 969 DMBoundaryType bdt[3]; 970 PetscErrorCode ierr; 971 972 PetscFunctionBegin; 973 low[0] = lower ? lower[0] : 0.0; low[1] = lower ? lower[1] : 0.0; low[2] = lower ? lower[2] : 0.0; 974 upp[0] = upper ? upper[0] : 1.0; upp[1] = upper ? upper[1] : 1.0; upp[2] = upper ? upper[2] : 1.0; 975 fac[0] = faces ? faces[0] : 4-dim; fac[1] = faces ? faces[1] : 4-dim; fac[2] = dim > 2 ? (faces ? faces[2] : 4-dim) : 0; 976 bdt[0] = periodicity ? periodicity[0] : DM_BOUNDARY_NONE; bdt[1] = periodicity ? periodicity[1] : DM_BOUNDARY_NONE; bdt[2] = periodicity ? periodicity[2] : DM_BOUNDARY_NONE; 977 if (simplex) {ierr = DMPlexCreateBoxMesh_Simplex_Internal(comm, dim, fac, low, upp, bdt, interpolate, dm);CHKERRQ(ierr);} 978 else {ierr = DMPlexCreateBoxMesh_Tensor_Internal(comm, dim, fac, low, upp, bdt, interpolate, dm);CHKERRQ(ierr);} 979 PetscFunctionReturn(0); 980 } 981 982 /*@C 983 DMPlexSetOptionsPrefix - Sets the prefix used for searching for all DM options in the database. 984 985 Logically Collective on DM 986 987 Input Parameters: 988 + dm - the DM context 989 - prefix - the prefix to prepend to all option names 990 991 Notes: 992 A hyphen (-) must NOT be given at the beginning of the prefix name. 993 The first character of all runtime options is AUTOMATICALLY the hyphen. 994 995 Level: advanced 996 997 .seealso: SNESSetFromOptions() 998 @*/ 999 PetscErrorCode DMPlexSetOptionsPrefix(DM dm, const char prefix[]) 1000 { 1001 DM_Plex *mesh = (DM_Plex *) dm->data; 1002 PetscErrorCode ierr; 1003 1004 PetscFunctionBegin; 1005 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1006 ierr = PetscObjectSetOptionsPrefix((PetscObject) dm, prefix);CHKERRQ(ierr); 1007 ierr = PetscObjectSetOptionsPrefix((PetscObject) mesh->partitioner, prefix);CHKERRQ(ierr); 1008 PetscFunctionReturn(0); 1009 } 1010 1011 /*@ 1012 DMPlexCreateHexCylinderMesh - Creates a mesh on the tensor product of the unit interval with the circle (cylinder) using hexahedra. 1013 1014 Collective on MPI_Comm 1015 1016 Input Parameters: 1017 + comm - The communicator for the DM object 1018 . numRefine - The number of regular refinements to the basic 5 cell structure 1019 - periodicZ - The boundary type for the Z direction 1020 1021 Output Parameter: 1022 . dm - The DM object 1023 1024 Note: Here is the output numbering looking from the bottom of the cylinder: 1025 $ 17-----14 1026 $ | | 1027 $ | 2 | 1028 $ | | 1029 $ 17-----8-----7-----14 1030 $ | | | | 1031 $ | 3 | 0 | 1 | 1032 $ | | | | 1033 $ 19-----5-----6-----13 1034 $ | | 1035 $ | 4 | 1036 $ | | 1037 $ 19-----13 1038 $ 1039 $ and up through the top 1040 $ 1041 $ 18-----16 1042 $ | | 1043 $ | 2 | 1044 $ | | 1045 $ 18----10----11-----16 1046 $ | | | | 1047 $ | 3 | 0 | 1 | 1048 $ | | | | 1049 $ 20-----9----12-----15 1050 $ | | 1051 $ | 4 | 1052 $ | | 1053 $ 20-----15 1054 1055 Level: beginner 1056 1057 .keywords: DM, create 1058 .seealso: DMPlexCreateBoxMesh(), DMSetType(), DMCreate() 1059 @*/ 1060 PetscErrorCode DMPlexCreateHexCylinderMesh(MPI_Comm comm, PetscInt numRefine, DMBoundaryType periodicZ, DM *dm) 1061 { 1062 const PetscInt dim = 3; 1063 PetscInt numCells, numVertices, r; 1064 PetscErrorCode ierr; 1065 1066 PetscFunctionBegin; 1067 PetscValidPointer(dm, 4); 1068 if (numRefine < 0) SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Number of refinements %D cannot be negative", numRefine); 1069 ierr = DMCreate(comm, dm);CHKERRQ(ierr); 1070 ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 1071 ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr); 1072 /* Create topology */ 1073 { 1074 PetscInt cone[8], c; 1075 1076 numCells = 5; 1077 numVertices = 16; 1078 if (periodicZ == DM_BOUNDARY_PERIODIC) { 1079 numCells *= 3; 1080 numVertices = 24; 1081 } 1082 ierr = DMPlexSetChart(*dm, 0, numCells+numVertices);CHKERRQ(ierr); 1083 for (c = 0; c < numCells; c++) {ierr = DMPlexSetConeSize(*dm, c, 8);CHKERRQ(ierr);} 1084 ierr = DMSetUp(*dm);CHKERRQ(ierr); 1085 if (periodicZ == DM_BOUNDARY_PERIODIC) { 1086 cone[0] = 15; cone[1] = 18; cone[2] = 17; cone[3] = 16; 1087 cone[4] = 31; cone[5] = 32; cone[6] = 33; cone[7] = 34; 1088 ierr = DMPlexSetCone(*dm, 0, cone);CHKERRQ(ierr); 1089 cone[0] = 16; cone[1] = 17; cone[2] = 24; cone[3] = 23; 1090 cone[4] = 32; cone[5] = 36; cone[6] = 37; cone[7] = 33; /* 22 25 26 21 */ 1091 ierr = DMPlexSetCone(*dm, 1, cone);CHKERRQ(ierr); 1092 cone[0] = 18; cone[1] = 27; cone[2] = 24; cone[3] = 17; 1093 cone[4] = 34; cone[5] = 33; cone[6] = 37; cone[7] = 38; 1094 ierr = DMPlexSetCone(*dm, 2, cone);CHKERRQ(ierr); 1095 cone[0] = 29; cone[1] = 27; cone[2] = 18; cone[3] = 15; 1096 cone[4] = 35; cone[5] = 31; cone[6] = 34; cone[7] = 38; 1097 ierr = DMPlexSetCone(*dm, 3, cone);CHKERRQ(ierr); 1098 cone[0] = 29; cone[1] = 15; cone[2] = 16; cone[3] = 23; 1099 cone[4] = 35; cone[5] = 36; cone[6] = 32; cone[7] = 31; 1100 ierr = DMPlexSetCone(*dm, 4, cone);CHKERRQ(ierr); 1101 1102 cone[0] = 31; cone[1] = 34; cone[2] = 33; cone[3] = 32; 1103 cone[4] = 19; cone[5] = 22; cone[6] = 21; cone[7] = 20; 1104 ierr = DMPlexSetCone(*dm, 5, cone);CHKERRQ(ierr); 1105 cone[0] = 32; cone[1] = 33; cone[2] = 37; cone[3] = 36; 1106 cone[4] = 22; cone[5] = 25; cone[6] = 26; cone[7] = 21; 1107 ierr = DMPlexSetCone(*dm, 6, cone);CHKERRQ(ierr); 1108 cone[0] = 34; cone[1] = 38; cone[2] = 37; cone[3] = 33; 1109 cone[4] = 20; cone[5] = 21; cone[6] = 26; cone[7] = 28; 1110 ierr = DMPlexSetCone(*dm, 7, cone);CHKERRQ(ierr); 1111 cone[0] = 35; cone[1] = 38; cone[2] = 34; cone[3] = 31; 1112 cone[4] = 30; cone[5] = 19; cone[6] = 20; cone[7] = 28; 1113 ierr = DMPlexSetCone(*dm, 8, cone);CHKERRQ(ierr); 1114 cone[0] = 35; cone[1] = 31; cone[2] = 32; cone[3] = 36; 1115 cone[4] = 30; cone[5] = 25; cone[6] = 22; cone[7] = 19; 1116 ierr = DMPlexSetCone(*dm, 9, cone);CHKERRQ(ierr); 1117 1118 cone[0] = 19; cone[1] = 20; cone[2] = 21; cone[3] = 22; 1119 cone[4] = 15; cone[5] = 16; cone[6] = 17; cone[7] = 18; 1120 ierr = DMPlexSetCone(*dm, 10, cone);CHKERRQ(ierr); 1121 cone[0] = 22; cone[1] = 21; cone[2] = 26; cone[3] = 25; 1122 cone[4] = 16; cone[5] = 23; cone[6] = 24; cone[7] = 17; 1123 ierr = DMPlexSetCone(*dm, 11, cone);CHKERRQ(ierr); 1124 cone[0] = 20; cone[1] = 28; cone[2] = 26; cone[3] = 21; 1125 cone[4] = 18; cone[5] = 17; cone[6] = 24; cone[7] = 27; 1126 ierr = DMPlexSetCone(*dm, 12, cone);CHKERRQ(ierr); 1127 cone[0] = 30; cone[1] = 28; cone[2] = 20; cone[3] = 19; 1128 cone[4] = 29; cone[5] = 15; cone[6] = 18; cone[7] = 27; 1129 ierr = DMPlexSetCone(*dm, 13, cone);CHKERRQ(ierr); 1130 cone[0] = 30; cone[1] = 19; cone[2] = 22; cone[3] = 25; 1131 cone[4] = 29; cone[5] = 23; cone[6] = 16; cone[7] = 15; 1132 ierr = DMPlexSetCone(*dm, 14, cone);CHKERRQ(ierr); 1133 } else { 1134 cone[0] = 5; cone[1] = 8; cone[2] = 7; cone[3] = 6; 1135 cone[4] = 9; cone[5] = 12; cone[6] = 11; cone[7] = 10; 1136 ierr = DMPlexSetCone(*dm, 0, cone);CHKERRQ(ierr); 1137 cone[0] = 6; cone[1] = 7; cone[2] = 14; cone[3] = 13; 1138 cone[4] = 12; cone[5] = 15; cone[6] = 16; cone[7] = 11; 1139 ierr = DMPlexSetCone(*dm, 1, cone);CHKERRQ(ierr); 1140 cone[0] = 8; cone[1] = 17; cone[2] = 14; cone[3] = 7; 1141 cone[4] = 10; cone[5] = 11; cone[6] = 16; cone[7] = 18; 1142 ierr = DMPlexSetCone(*dm, 2, cone);CHKERRQ(ierr); 1143 cone[0] = 19; cone[1] = 17; cone[2] = 8; cone[3] = 5; 1144 cone[4] = 20; cone[5] = 9; cone[6] = 10; cone[7] = 18; 1145 ierr = DMPlexSetCone(*dm, 3, cone);CHKERRQ(ierr); 1146 cone[0] = 19; cone[1] = 5; cone[2] = 6; cone[3] = 13; 1147 cone[4] = 20; cone[5] = 15; cone[6] = 12; cone[7] = 9; 1148 ierr = DMPlexSetCone(*dm, 4, cone);CHKERRQ(ierr); 1149 } 1150 ierr = DMPlexSymmetrize(*dm);CHKERRQ(ierr); 1151 ierr = DMPlexStratify(*dm);CHKERRQ(ierr); 1152 } 1153 /* Interpolate */ 1154 { 1155 DM idm = NULL; 1156 1157 ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr); 1158 ierr = DMDestroy(dm);CHKERRQ(ierr); 1159 *dm = idm; 1160 } 1161 /* Create cube geometry */ 1162 { 1163 Vec coordinates; 1164 PetscSection coordSection; 1165 PetscScalar *coords; 1166 PetscInt coordSize, v; 1167 const PetscReal dis = 1.0/PetscSqrtReal(2.0); 1168 const PetscReal ds2 = dis/2.0; 1169 1170 /* Build coordinates */ 1171 ierr = DMGetCoordinateSection(*dm, &coordSection);CHKERRQ(ierr); 1172 ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr); 1173 ierr = PetscSectionSetFieldComponents(coordSection, 0, dim);CHKERRQ(ierr); 1174 ierr = PetscSectionSetChart(coordSection, numCells, numCells+numVertices);CHKERRQ(ierr); 1175 for (v = numCells; v < numCells+numVertices; ++v) { 1176 ierr = PetscSectionSetDof(coordSection, v, dim);CHKERRQ(ierr); 1177 ierr = PetscSectionSetFieldDof(coordSection, v, 0, dim);CHKERRQ(ierr); 1178 } 1179 ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 1180 ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 1181 ierr = VecCreate(PETSC_COMM_SELF, &coordinates);CHKERRQ(ierr); 1182 ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr); 1183 ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 1184 ierr = VecSetBlockSize(coordinates, dim);CHKERRQ(ierr); 1185 ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr); 1186 ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 1187 coords[0*dim+0] = -ds2; coords[0*dim+1] = -ds2; coords[0*dim+2] = 0.0; 1188 coords[1*dim+0] = ds2; coords[1*dim+1] = -ds2; coords[1*dim+2] = 0.0; 1189 coords[2*dim+0] = ds2; coords[2*dim+1] = ds2; coords[2*dim+2] = 0.0; 1190 coords[3*dim+0] = -ds2; coords[3*dim+1] = ds2; coords[3*dim+2] = 0.0; 1191 coords[4*dim+0] = -ds2; coords[4*dim+1] = -ds2; coords[4*dim+2] = 1.0; 1192 coords[5*dim+0] = -ds2; coords[5*dim+1] = ds2; coords[5*dim+2] = 1.0; 1193 coords[6*dim+0] = ds2; coords[6*dim+1] = ds2; coords[6*dim+2] = 1.0; 1194 coords[7*dim+0] = ds2; coords[7*dim+1] = -ds2; coords[7*dim+2] = 1.0; 1195 coords[ 8*dim+0] = dis; coords[ 8*dim+1] = -dis; coords[ 8*dim+2] = 0.0; 1196 coords[ 9*dim+0] = dis; coords[ 9*dim+1] = dis; coords[ 9*dim+2] = 0.0; 1197 coords[10*dim+0] = dis; coords[10*dim+1] = -dis; coords[10*dim+2] = 1.0; 1198 coords[11*dim+0] = dis; coords[11*dim+1] = dis; coords[11*dim+2] = 1.0; 1199 coords[12*dim+0] = -dis; coords[12*dim+1] = dis; coords[12*dim+2] = 0.0; 1200 coords[13*dim+0] = -dis; coords[13*dim+1] = dis; coords[13*dim+2] = 1.0; 1201 coords[14*dim+0] = -dis; coords[14*dim+1] = -dis; coords[14*dim+2] = 0.0; 1202 coords[15*dim+0] = -dis; coords[15*dim+1] = -dis; coords[15*dim+2] = 1.0; 1203 if (periodicZ == DM_BOUNDARY_PERIODIC) { 1204 /* 15 31 19 */ coords[16*dim+0] = -ds2; coords[16*dim+1] = -ds2; coords[16*dim+2] = 0.5; 1205 /* 16 32 22 */ coords[17*dim+0] = ds2; coords[17*dim+1] = -ds2; coords[17*dim+2] = 0.5; 1206 /* 17 33 21 */ coords[18*dim+0] = ds2; coords[18*dim+1] = ds2; coords[18*dim+2] = 0.5; 1207 /* 18 34 20 */ coords[19*dim+0] = -ds2; coords[19*dim+1] = ds2; coords[19*dim+2] = 0.5; 1208 /* 29 35 30 */ coords[20*dim+0] = -dis; coords[20*dim+1] = -dis; coords[20*dim+2] = 0.5; 1209 /* 23 36 25 */ coords[21*dim+0] = dis; coords[21*dim+1] = -dis; coords[21*dim+2] = 0.5; 1210 /* 24 37 26 */ coords[22*dim+0] = dis; coords[22*dim+1] = dis; coords[22*dim+2] = 0.5; 1211 /* 27 38 28 */ coords[23*dim+0] = -dis; coords[23*dim+1] = dis; coords[23*dim+2] = 0.5; 1212 } 1213 ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 1214 ierr = DMSetCoordinatesLocal(*dm, coordinates);CHKERRQ(ierr); 1215 ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 1216 } 1217 /* Create periodicity */ 1218 if (periodicZ == DM_BOUNDARY_PERIODIC || periodicZ == DM_BOUNDARY_TWIST) { 1219 PetscReal L[3]; 1220 PetscReal maxCell[3]; 1221 DMBoundaryType bdType[3]; 1222 PetscReal lower[3] = {0.0, 0.0, 0.0}; 1223 PetscReal upper[3] = {1.0, 1.0, 1.5}; 1224 PetscInt i, numZCells = 3; 1225 1226 bdType[0] = DM_BOUNDARY_NONE; 1227 bdType[1] = DM_BOUNDARY_NONE; 1228 bdType[2] = periodicZ; 1229 for (i = 0; i < dim; i++) { 1230 L[i] = upper[i] - lower[i]; 1231 maxCell[i] = 1.1 * (L[i] / numZCells); 1232 } 1233 ierr = DMSetPeriodicity(*dm, PETSC_TRUE, maxCell, L, bdType);CHKERRQ(ierr); 1234 } 1235 /* Refine topology */ 1236 for (r = 0; r < numRefine; ++r) { 1237 DM rdm = NULL; 1238 1239 ierr = DMRefine(*dm, comm, &rdm);CHKERRQ(ierr); 1240 ierr = DMDestroy(dm);CHKERRQ(ierr); 1241 *dm = rdm; 1242 } 1243 /* Remap geometry to cylinder 1244 Interior square: Linear interpolation is correct 1245 The other cells all have vertices on rays from the origin. We want to uniformly expand the spacing 1246 such that the last vertex is on the unit circle. So the closest and farthest vertices are at distance 1247 1248 phi = arctan(y/x) 1249 d_close = sqrt(1/8 + 1/4 sin^2(phi)) 1250 d_far = sqrt(1/2 + sin^2(phi)) 1251 1252 so we remap them using 1253 1254 x_new = x_close + (x - x_close) (1 - d_close) / (d_far - d_close) 1255 y_new = y_close + (y - y_close) (1 - d_close) / (d_far - d_close) 1256 1257 If pi/4 < phi < 3pi/4 or -3pi/4 < phi < -pi/4, then we switch x and y. 1258 */ 1259 { 1260 Vec coordinates; 1261 PetscSection coordSection; 1262 PetscScalar *coords; 1263 PetscInt vStart, vEnd, v; 1264 const PetscReal dis = 1.0/PetscSqrtReal(2.0); 1265 const PetscReal ds2 = 0.5*dis; 1266 1267 ierr = DMPlexGetDepthStratum(*dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 1268 ierr = DMGetCoordinateSection(*dm, &coordSection);CHKERRQ(ierr); 1269 ierr = DMGetCoordinatesLocal(*dm, &coordinates);CHKERRQ(ierr); 1270 ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 1271 for (v = vStart; v < vEnd; ++v) { 1272 PetscReal phi, sinp, cosp, dc, df, x, y, xc, yc; 1273 PetscInt off; 1274 1275 ierr = PetscSectionGetOffset(coordSection, v, &off);CHKERRQ(ierr); 1276 if ((PetscAbsScalar(coords[off+0]) <= ds2) && (PetscAbsScalar(coords[off+1]) <= ds2)) continue; 1277 x = PetscRealPart(coords[off]); 1278 y = PetscRealPart(coords[off+1]); 1279 phi = PetscAtan2Real(y, x); 1280 sinp = PetscSinReal(phi); 1281 cosp = PetscCosReal(phi); 1282 if ((PetscAbsReal(phi) > PETSC_PI/4.0) && (PetscAbsReal(phi) < 3.0*PETSC_PI/4.0)) { 1283 dc = PetscAbsReal(ds2/sinp); 1284 df = PetscAbsReal(dis/sinp); 1285 xc = ds2*x/PetscAbsReal(y); 1286 yc = ds2*PetscSignReal(y); 1287 } else { 1288 dc = PetscAbsReal(ds2/cosp); 1289 df = PetscAbsReal(dis/cosp); 1290 xc = ds2*PetscSignReal(x); 1291 yc = ds2*y/PetscAbsReal(x); 1292 } 1293 coords[off+0] = xc + (coords[off+0] - xc)*(1.0 - dc)/(df - dc); 1294 coords[off+1] = yc + (coords[off+1] - yc)*(1.0 - dc)/(df - dc); 1295 } 1296 ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 1297 if (periodicZ == DM_BOUNDARY_PERIODIC || periodicZ == DM_BOUNDARY_TWIST) { 1298 ierr = DMLocalizeCoordinates(*dm);CHKERRQ(ierr); 1299 } 1300 } 1301 PetscFunctionReturn(0); 1302 } 1303 1304 /*@ 1305 DMPlexCreateWedgeCylinderMesh - Creates a mesh on the tensor product of the unit interval with the circle (cylinder) using wedges. 1306 1307 Collective on MPI_Comm 1308 1309 Input Parameters: 1310 + comm - The communicator for the DM object 1311 . n - The number of wedges around the origin 1312 - interpolate - Create edges and faces 1313 1314 Output Parameter: 1315 . dm - The DM object 1316 1317 Level: beginner 1318 1319 .keywords: DM, create 1320 .seealso: DMPlexCreateHexCylinderMesh(), DMPlexCreateBoxMesh(), DMSetType(), DMCreate() 1321 @*/ 1322 PetscErrorCode DMPlexCreateWedgeCylinderMesh(MPI_Comm comm, PetscInt n, PetscBool interpolate, DM *dm) 1323 { 1324 const PetscInt dim = 3; 1325 PetscInt numCells, numVertices; 1326 PetscErrorCode ierr; 1327 1328 PetscFunctionBegin; 1329 PetscValidPointer(dm, 3); 1330 if (n < 0) SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Number of wedges %D cannot be negative", n); 1331 ierr = DMCreate(comm, dm);CHKERRQ(ierr); 1332 ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 1333 ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr); 1334 /* Create topology */ 1335 { 1336 PetscInt cone[6], c; 1337 1338 numCells = n; 1339 numVertices = 2*(n+1); 1340 ierr = DMPlexSetChart(*dm, 0, numCells+numVertices);CHKERRQ(ierr); 1341 for (c = 0; c < numCells; c++) {ierr = DMPlexSetConeSize(*dm, c, 6);CHKERRQ(ierr);} 1342 ierr = DMSetUp(*dm);CHKERRQ(ierr); 1343 for (c = 0; c < numCells; c++) { 1344 cone[0] = c+n*1; cone[1] = (c+1)%n+n*1; cone[2] = 0+3*n; 1345 cone[3] = c+n*2; cone[4] = (c+1)%n+n*2; cone[5] = 1+3*n; 1346 ierr = DMPlexSetCone(*dm, c, cone);CHKERRQ(ierr); 1347 } 1348 ierr = DMPlexSymmetrize(*dm);CHKERRQ(ierr); 1349 ierr = DMPlexStratify(*dm);CHKERRQ(ierr); 1350 } 1351 /* Interpolate */ 1352 if (interpolate) { 1353 DM idm = NULL; 1354 1355 ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr); 1356 ierr = DMDestroy(dm);CHKERRQ(ierr); 1357 *dm = idm; 1358 } 1359 /* Create cylinder geometry */ 1360 { 1361 Vec coordinates; 1362 PetscSection coordSection; 1363 PetscScalar *coords; 1364 PetscInt coordSize, v, c; 1365 1366 /* Build coordinates */ 1367 ierr = DMGetCoordinateSection(*dm, &coordSection);CHKERRQ(ierr); 1368 ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr); 1369 ierr = PetscSectionSetFieldComponents(coordSection, 0, dim);CHKERRQ(ierr); 1370 ierr = PetscSectionSetChart(coordSection, numCells, numCells+numVertices);CHKERRQ(ierr); 1371 for (v = numCells; v < numCells+numVertices; ++v) { 1372 ierr = PetscSectionSetDof(coordSection, v, dim);CHKERRQ(ierr); 1373 ierr = PetscSectionSetFieldDof(coordSection, v, 0, dim);CHKERRQ(ierr); 1374 } 1375 ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 1376 ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 1377 ierr = VecCreate(PETSC_COMM_SELF, &coordinates);CHKERRQ(ierr); 1378 ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr); 1379 ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 1380 ierr = VecSetBlockSize(coordinates, dim);CHKERRQ(ierr); 1381 ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr); 1382 ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 1383 for (c = 0; c < numCells; c++) { 1384 coords[(c+0*n)*dim+0] = PetscCosReal(2.0*c*PETSC_PI/n); coords[(c+0*n)*dim+1] = PetscSinReal(2.0*c*PETSC_PI/n); coords[(c+0*n)*dim+2] = 1.0; 1385 coords[(c+1*n)*dim+0] = PetscCosReal(2.0*c*PETSC_PI/n); coords[(c+1*n)*dim+1] = PetscSinReal(2.0*c*PETSC_PI/n); coords[(c+1*n)*dim+2] = 0.0; 1386 } 1387 coords[(2*n+0)*dim+0] = 0.0; coords[(2*n+0)*dim+1] = 0.0; coords[(2*n+0)*dim+2] = 1.0; 1388 coords[(2*n+1)*dim+0] = 0.0; coords[(2*n+1)*dim+1] = 0.0; coords[(2*n+1)*dim+2] = 0.0; 1389 ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 1390 ierr = DMSetCoordinatesLocal(*dm, coordinates);CHKERRQ(ierr); 1391 ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 1392 } 1393 PetscFunctionReturn(0); 1394 } 1395 1396 PETSC_STATIC_INLINE PetscReal DiffNormReal(PetscInt dim, const PetscReal x[], const PetscReal y[]) 1397 { 1398 PetscReal prod = 0.0; 1399 PetscInt i; 1400 for (i = 0; i < dim; ++i) prod += PetscSqr(x[i] - y[i]); 1401 return PetscSqrtReal(prod); 1402 } 1403 PETSC_STATIC_INLINE PetscReal DotReal(PetscInt dim, const PetscReal x[], const PetscReal y[]) 1404 { 1405 PetscReal prod = 0.0; 1406 PetscInt i; 1407 for (i = 0; i < dim; ++i) prod += x[i]*y[i]; 1408 return prod; 1409 } 1410 1411 /*@ 1412 DMPlexCreateSphereMesh - Creates a mesh on the d-dimensional sphere, S^d. 1413 1414 Collective on MPI_Comm 1415 1416 Input Parameters: 1417 . comm - The communicator for the DM object 1418 . dim - The dimension 1419 - simplex - Use simplices, or tensor product cells 1420 1421 Output Parameter: 1422 . dm - The DM object 1423 1424 Level: beginner 1425 1426 .keywords: DM, create 1427 .seealso: DMPlexCreateBoxMesh(), DMSetType(), DMCreate() 1428 @*/ 1429 PetscErrorCode DMPlexCreateSphereMesh(MPI_Comm comm, PetscInt dim, PetscBool simplex, DM *dm) 1430 { 1431 const PetscInt embedDim = dim+1; 1432 PetscSection coordSection; 1433 Vec coordinates; 1434 PetscScalar *coords; 1435 PetscReal *coordsIn; 1436 PetscInt numCells, numEdges, numVerts, firstVertex, v, firstEdge, coordSize, d, c, e; 1437 PetscMPIInt rank; 1438 PetscErrorCode ierr; 1439 1440 PetscFunctionBegin; 1441 PetscValidPointer(dm, 4); 1442 ierr = DMCreate(comm, dm);CHKERRQ(ierr); 1443 ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 1444 ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr); 1445 ierr = DMSetCoordinateDim(*dm, dim+1);CHKERRQ(ierr); 1446 ierr = MPI_Comm_rank(PetscObjectComm((PetscObject) *dm), &rank);CHKERRQ(ierr); 1447 switch (dim) { 1448 case 2: 1449 if (simplex) { 1450 DM idm = NULL; 1451 const PetscReal edgeLen = 2.0/(1.0 + PETSC_PHI); 1452 const PetscReal vertex[3] = {0.0, 1.0/(1.0 + PETSC_PHI), PETSC_PHI/(1.0 + PETSC_PHI)}; 1453 const PetscInt degree = 5; 1454 PetscInt s[3] = {1, 1, 1}; 1455 PetscInt cone[3]; 1456 PetscInt *graph, p, i, j, k; 1457 1458 numCells = !rank ? 20 : 0; 1459 numVerts = !rank ? 12 : 0; 1460 firstVertex = numCells; 1461 /* Use icosahedron, which for a unit sphere has coordinates which are all cyclic permutations of 1462 1463 (0, \pm 1/\phi+1, \pm \phi/\phi+1) 1464 1465 where \phi^2 - \phi - 1 = 0, meaning \phi is the golden ratio \frac{1 + \sqrt{5}}{2}. The edge 1466 length is then given by 2/\phi = 2 * 2.73606 = 5.47214. 1467 */ 1468 /* Construct vertices */ 1469 ierr = PetscCalloc1(numVerts * embedDim, &coordsIn);CHKERRQ(ierr); 1470 for (p = 0, i = 0; p < embedDim; ++p) { 1471 for (s[1] = -1; s[1] < 2; s[1] += 2) { 1472 for (s[2] = -1; s[2] < 2; s[2] += 2) { 1473 for (d = 0; d < embedDim; ++d) coordsIn[i*embedDim+d] = s[(d+p)%embedDim]*vertex[(d+p)%embedDim]; 1474 ++i; 1475 } 1476 } 1477 } 1478 /* Construct graph */ 1479 ierr = PetscCalloc1(numVerts * numVerts, &graph);CHKERRQ(ierr); 1480 for (i = 0; i < numVerts; ++i) { 1481 for (j = 0, k = 0; j < numVerts; ++j) { 1482 if (PetscAbsReal(DiffNormReal(embedDim, &coordsIn[i*embedDim], &coordsIn[j*embedDim]) - edgeLen) < PETSC_SMALL) {graph[i*numVerts+j] = 1; ++k;} 1483 } 1484 if (k != degree) SETERRQ3(comm, PETSC_ERR_PLIB, "Invalid icosahedron, vertex %D degree %D != %D", i, k, degree); 1485 } 1486 /* Build Topology */ 1487 ierr = DMPlexSetChart(*dm, 0, numCells+numVerts);CHKERRQ(ierr); 1488 for (c = 0; c < numCells; c++) { 1489 ierr = DMPlexSetConeSize(*dm, c, embedDim);CHKERRQ(ierr); 1490 } 1491 ierr = DMSetUp(*dm);CHKERRQ(ierr); /* Allocate space for cones */ 1492 /* Cells */ 1493 for (i = 0, c = 0; i < numVerts; ++i) { 1494 for (j = 0; j < i; ++j) { 1495 for (k = 0; k < j; ++k) { 1496 if (graph[i*numVerts+j] && graph[j*numVerts+k] && graph[k*numVerts+i]) { 1497 cone[0] = firstVertex+i; cone[1] = firstVertex+j; cone[2] = firstVertex+k; 1498 /* Check orientation */ 1499 { 1500 const PetscInt epsilon[3][3][3] = {{{0, 0, 0}, {0, 0, 1}, {0, -1, 0}}, {{0, 0, -1}, {0, 0, 0}, {1, 0, 0}}, {{0, 1, 0}, {-1, 0, 0}, {0, 0, 0}}}; 1501 PetscReal normal[3]; 1502 PetscInt e, f; 1503 1504 for (d = 0; d < embedDim; ++d) { 1505 normal[d] = 0.0; 1506 for (e = 0; e < embedDim; ++e) { 1507 for (f = 0; f < embedDim; ++f) { 1508 normal[d] += epsilon[d][e][f]*(coordsIn[j*embedDim+e] - coordsIn[i*embedDim+e])*(coordsIn[k*embedDim+f] - coordsIn[i*embedDim+f]); 1509 } 1510 } 1511 } 1512 if (DotReal(embedDim, normal, &coordsIn[i*embedDim]) < 0) {PetscInt tmp = cone[1]; cone[1] = cone[2]; cone[2] = tmp;} 1513 } 1514 ierr = DMPlexSetCone(*dm, c++, cone);CHKERRQ(ierr); 1515 } 1516 } 1517 } 1518 } 1519 ierr = DMPlexSymmetrize(*dm);CHKERRQ(ierr); 1520 ierr = DMPlexStratify(*dm);CHKERRQ(ierr); 1521 ierr = PetscFree(graph);CHKERRQ(ierr); 1522 /* Interpolate mesh */ 1523 ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr); 1524 ierr = DMDestroy(dm);CHKERRQ(ierr); 1525 *dm = idm; 1526 } else { 1527 /* 1528 12-21--13 1529 | | 1530 25 4 24 1531 | | 1532 12-25--9-16--8-24--13 1533 | | | | 1534 23 5 17 0 15 3 22 1535 | | | | 1536 10-20--6-14--7-19--11 1537 | | 1538 20 1 19 1539 | | 1540 10-18--11 1541 | | 1542 23 2 22 1543 | | 1544 12-21--13 1545 */ 1546 const PetscReal dist = 1.0/PetscSqrtReal(3.0); 1547 PetscInt cone[4], ornt[4]; 1548 1549 numCells = !rank ? 6 : 0; 1550 numEdges = !rank ? 12 : 0; 1551 numVerts = !rank ? 8 : 0; 1552 firstVertex = numCells; 1553 firstEdge = numCells + numVerts; 1554 /* Build Topology */ 1555 ierr = DMPlexSetChart(*dm, 0, numCells+numEdges+numVerts);CHKERRQ(ierr); 1556 for (c = 0; c < numCells; c++) { 1557 ierr = DMPlexSetConeSize(*dm, c, 4);CHKERRQ(ierr); 1558 } 1559 for (e = firstEdge; e < firstEdge+numEdges; ++e) { 1560 ierr = DMPlexSetConeSize(*dm, e, 2);CHKERRQ(ierr); 1561 } 1562 ierr = DMSetUp(*dm);CHKERRQ(ierr); /* Allocate space for cones */ 1563 /* Cell 0 */ 1564 cone[0] = 14; cone[1] = 15; cone[2] = 16; cone[3] = 17; 1565 ierr = DMPlexSetCone(*dm, 0, cone);CHKERRQ(ierr); 1566 ornt[0] = 0; ornt[1] = 0; ornt[2] = 0; ornt[3] = 0; 1567 ierr = DMPlexSetConeOrientation(*dm, 0, ornt);CHKERRQ(ierr); 1568 /* Cell 1 */ 1569 cone[0] = 18; cone[1] = 19; cone[2] = 14; cone[3] = 20; 1570 ierr = DMPlexSetCone(*dm, 1, cone);CHKERRQ(ierr); 1571 ornt[0] = 0; ornt[1] = 0; ornt[2] = -2; ornt[3] = 0; 1572 ierr = DMPlexSetConeOrientation(*dm, 1, ornt);CHKERRQ(ierr); 1573 /* Cell 2 */ 1574 cone[0] = 21; cone[1] = 22; cone[2] = 18; cone[3] = 23; 1575 ierr = DMPlexSetCone(*dm, 2, cone);CHKERRQ(ierr); 1576 ornt[0] = 0; ornt[1] = 0; ornt[2] = -2; ornt[3] = 0; 1577 ierr = DMPlexSetConeOrientation(*dm, 2, ornt);CHKERRQ(ierr); 1578 /* Cell 3 */ 1579 cone[0] = 19; cone[1] = 22; cone[2] = 24; cone[3] = 15; 1580 ierr = DMPlexSetCone(*dm, 3, cone);CHKERRQ(ierr); 1581 ornt[0] = -2; ornt[1] = -2; ornt[2] = 0; ornt[3] = -2; 1582 ierr = DMPlexSetConeOrientation(*dm, 3, ornt);CHKERRQ(ierr); 1583 /* Cell 4 */ 1584 cone[0] = 16; cone[1] = 24; cone[2] = 21; cone[3] = 25; 1585 ierr = DMPlexSetCone(*dm, 4, cone);CHKERRQ(ierr); 1586 ornt[0] = -2; ornt[1] = -2; ornt[2] = -2; ornt[3] = 0; 1587 ierr = DMPlexSetConeOrientation(*dm, 4, ornt);CHKERRQ(ierr); 1588 /* Cell 5 */ 1589 cone[0] = 20; cone[1] = 17; cone[2] = 25; cone[3] = 23; 1590 ierr = DMPlexSetCone(*dm, 5, cone);CHKERRQ(ierr); 1591 ornt[0] = -2; ornt[1] = -2; ornt[2] = -2; ornt[3] = -2; 1592 ierr = DMPlexSetConeOrientation(*dm, 5, ornt);CHKERRQ(ierr); 1593 /* Edges */ 1594 cone[0] = 6; cone[1] = 7; 1595 ierr = DMPlexSetCone(*dm, 14, cone);CHKERRQ(ierr); 1596 cone[0] = 7; cone[1] = 8; 1597 ierr = DMPlexSetCone(*dm, 15, cone);CHKERRQ(ierr); 1598 cone[0] = 8; cone[1] = 9; 1599 ierr = DMPlexSetCone(*dm, 16, cone);CHKERRQ(ierr); 1600 cone[0] = 9; cone[1] = 6; 1601 ierr = DMPlexSetCone(*dm, 17, cone);CHKERRQ(ierr); 1602 cone[0] = 10; cone[1] = 11; 1603 ierr = DMPlexSetCone(*dm, 18, cone);CHKERRQ(ierr); 1604 cone[0] = 11; cone[1] = 7; 1605 ierr = DMPlexSetCone(*dm, 19, cone);CHKERRQ(ierr); 1606 cone[0] = 6; cone[1] = 10; 1607 ierr = DMPlexSetCone(*dm, 20, cone);CHKERRQ(ierr); 1608 cone[0] = 12; cone[1] = 13; 1609 ierr = DMPlexSetCone(*dm, 21, cone);CHKERRQ(ierr); 1610 cone[0] = 13; cone[1] = 11; 1611 ierr = DMPlexSetCone(*dm, 22, cone);CHKERRQ(ierr); 1612 cone[0] = 10; cone[1] = 12; 1613 ierr = DMPlexSetCone(*dm, 23, cone);CHKERRQ(ierr); 1614 cone[0] = 13; cone[1] = 8; 1615 ierr = DMPlexSetCone(*dm, 24, cone);CHKERRQ(ierr); 1616 cone[0] = 12; cone[1] = 9; 1617 ierr = DMPlexSetCone(*dm, 25, cone);CHKERRQ(ierr); 1618 ierr = DMPlexSymmetrize(*dm);CHKERRQ(ierr); 1619 ierr = DMPlexStratify(*dm);CHKERRQ(ierr); 1620 /* Build coordinates */ 1621 ierr = PetscCalloc1(numVerts * embedDim, &coordsIn);CHKERRQ(ierr); 1622 coordsIn[0*embedDim+0] = -dist; coordsIn[0*embedDim+1] = dist; coordsIn[0*embedDim+2] = -dist; 1623 coordsIn[1*embedDim+0] = dist; coordsIn[1*embedDim+1] = dist; coordsIn[1*embedDim+2] = -dist; 1624 coordsIn[2*embedDim+0] = dist; coordsIn[2*embedDim+1] = -dist; coordsIn[2*embedDim+2] = -dist; 1625 coordsIn[3*embedDim+0] = -dist; coordsIn[3*embedDim+1] = -dist; coordsIn[3*embedDim+2] = -dist; 1626 coordsIn[4*embedDim+0] = -dist; coordsIn[4*embedDim+1] = dist; coordsIn[4*embedDim+2] = dist; 1627 coordsIn[5*embedDim+0] = dist; coordsIn[5*embedDim+1] = dist; coordsIn[5*embedDim+2] = dist; 1628 coordsIn[6*embedDim+0] = -dist; coordsIn[6*embedDim+1] = -dist; coordsIn[6*embedDim+2] = dist; 1629 coordsIn[7*embedDim+0] = dist; coordsIn[7*embedDim+1] = -dist; coordsIn[7*embedDim+2] = dist; 1630 } 1631 break; 1632 case 3: 1633 if (simplex) { 1634 DM idm = NULL; 1635 const PetscReal edgeLen = 1.0/PETSC_PHI; 1636 const PetscReal vertexA[4] = {0.5, 0.5, 0.5, 0.5}; 1637 const PetscReal vertexB[4] = {1.0, 0.0, 0.0, 0.0}; 1638 const PetscReal vertexC[4] = {0.5, 0.5*PETSC_PHI, 0.5/PETSC_PHI, 0.0}; 1639 const PetscInt degree = 12; 1640 PetscInt s[4] = {1, 1, 1}; 1641 PetscInt evenPerm[12][4] = {{0, 1, 2, 3}, {0, 2, 3, 1}, {0, 3, 1, 2}, {1, 0, 3, 2}, {1, 2, 0, 3}, {1, 3, 2, 0}, 1642 {2, 0, 1, 3}, {2, 1, 3, 0}, {2, 3, 0, 1}, {3, 0, 2, 1}, {3, 1, 0, 2}, {3, 2, 1, 0}}; 1643 PetscInt cone[4]; 1644 PetscInt *graph, p, i, j, k, l; 1645 1646 numCells = !rank ? 600 : 0; 1647 numVerts = !rank ? 120 : 0; 1648 firstVertex = numCells; 1649 /* Use the 600-cell, which for a unit sphere has coordinates which are 1650 1651 1/2 (\pm 1, \pm 1, \pm 1, \pm 1) 16 1652 (\pm 1, 0, 0, 0) all cyclic permutations 8 1653 1/2 (\pm 1, \pm phi, \pm 1/phi, 0) all even permutations 96 1654 1655 where \phi^2 - \phi - 1 = 0, meaning \phi is the golden ratio \frac{1 + \sqrt{5}}{2}. The edge 1656 length is then given by 1/\phi = 2.73606. 1657 1658 http://buzzard.pugetsound.edu/sage-practice/ch03s03.html 1659 http://mathworld.wolfram.com/600-Cell.html 1660 */ 1661 /* Construct vertices */ 1662 ierr = PetscCalloc1(numVerts * embedDim, &coordsIn);CHKERRQ(ierr); 1663 i = 0; 1664 for (s[0] = -1; s[0] < 2; s[0] += 2) { 1665 for (s[1] = -1; s[1] < 2; s[1] += 2) { 1666 for (s[2] = -1; s[2] < 2; s[2] += 2) { 1667 for (s[3] = -1; s[3] < 2; s[3] += 2) { 1668 for (d = 0; d < embedDim; ++d) coordsIn[i*embedDim+d] = s[d]*vertexA[d]; 1669 ++i; 1670 } 1671 } 1672 } 1673 } 1674 for (p = 0; p < embedDim; ++p) { 1675 s[1] = s[2] = s[3] = 1; 1676 for (s[0] = -1; s[0] < 2; s[0] += 2) { 1677 for (d = 0; d < embedDim; ++d) coordsIn[i*embedDim+d] = s[(d+p)%embedDim]*vertexB[(d+p)%embedDim]; 1678 ++i; 1679 } 1680 } 1681 for (p = 0; p < 12; ++p) { 1682 s[3] = 1; 1683 for (s[0] = -1; s[0] < 2; s[0] += 2) { 1684 for (s[1] = -1; s[1] < 2; s[1] += 2) { 1685 for (s[2] = -1; s[2] < 2; s[2] += 2) { 1686 for (d = 0; d < embedDim; ++d) coordsIn[i*embedDim+d] = s[evenPerm[p][d]]*vertexC[evenPerm[p][d]]; 1687 ++i; 1688 } 1689 } 1690 } 1691 } 1692 if (i != numVerts) SETERRQ2(comm, PETSC_ERR_PLIB, "Invalid 600-cell, vertices %D != %D", i, numVerts); 1693 /* Construct graph */ 1694 ierr = PetscCalloc1(numVerts * numVerts, &graph);CHKERRQ(ierr); 1695 for (i = 0; i < numVerts; ++i) { 1696 for (j = 0, k = 0; j < numVerts; ++j) { 1697 if (PetscAbsReal(DiffNormReal(embedDim, &coordsIn[i*embedDim], &coordsIn[j*embedDim]) - edgeLen) < PETSC_SMALL) {graph[i*numVerts+j] = 1; ++k;} 1698 } 1699 if (k != degree) SETERRQ3(comm, PETSC_ERR_PLIB, "Invalid 600-cell, vertex %D degree %D != %D", i, k, degree); 1700 } 1701 /* Build Topology */ 1702 ierr = DMPlexSetChart(*dm, 0, numCells+numVerts);CHKERRQ(ierr); 1703 for (c = 0; c < numCells; c++) { 1704 ierr = DMPlexSetConeSize(*dm, c, embedDim);CHKERRQ(ierr); 1705 } 1706 ierr = DMSetUp(*dm);CHKERRQ(ierr); /* Allocate space for cones */ 1707 /* Cells */ 1708 for (i = 0, c = 0; i < numVerts; ++i) { 1709 for (j = 0; j < i; ++j) { 1710 for (k = 0; k < j; ++k) { 1711 for (l = 0; l < k; ++l) { 1712 if (graph[i*numVerts+j] && graph[j*numVerts+k] && graph[k*numVerts+i] && 1713 graph[l*numVerts+i] && graph[l*numVerts+j] && graph[l*numVerts+k]) { 1714 cone[0] = firstVertex+i; cone[1] = firstVertex+j; cone[2] = firstVertex+k; cone[3] = firstVertex+l; 1715 /* Check orientation: https://ef.gy/linear-algebra:normal-vectors-in-higher-dimensional-spaces */ 1716 { 1717 const PetscInt epsilon[4][4][4][4] = {{{{0, 0, 0, 0}, { 0, 0, 0, 0}, { 0, 0, 0, 0}, { 0, 0, 0, 0}}, 1718 {{0, 0, 0, 0}, { 0, 0, 0, 0}, { 0, 0, 0, 1}, { 0, 0, -1, 0}}, 1719 {{0, 0, 0, 0}, { 0, 0, 0, -1}, { 0, 0, 0, 0}, { 0, 1, 0, 0}}, 1720 {{0, 0, 0, 0}, { 0, 0, 1, 0}, { 0, -1, 0, 0}, { 0, 0, 0, 0}}}, 1721 1722 {{{0, 0, 0, 0}, { 0, 0, 0, 0}, { 0, 0, 0, -1}, { 0, 0, 1, 0}}, 1723 {{0, 0, 0, 0}, { 0, 0, 0, 0}, { 0, 0, 0, 0}, { 0, 0, 0, 0}}, 1724 {{0, 0, 0, 1}, { 0, 0, 0, 0}, { 0, 0, 0, 0}, {-1, 0, 0, 0}}, 1725 {{0, 0, -1, 0}, { 0, 0, 0, 0}, { 1, 0, 0, 0}, { 0, 0, 0, 0}}}, 1726 1727 {{{0, 0, 0, 0}, { 0, 0, 0, 1}, { 0, 0, 0, 0}, { 0, -1, 0, 0}}, 1728 {{0, 0, 0, -1}, { 0, 0, 0, 0}, { 0, 0, 0, 0}, { 1, 0, 0, 0}}, 1729 {{0, 0, 0, 0}, { 0, 0, 0, 0}, { 0, 0, 0, 0}, { 0, 0, 0, 0}}, 1730 {{0, 1, 0, 0}, {-1, 0, 0, 0}, { 0, 0, 0, 0}, { 0, 0, 0, 0}}}, 1731 1732 {{{0, 0, 0, 0}, { 0, 0, -1, 0}, { 0, 1, 0, 0}, { 0, 0, 0, 0}}, 1733 {{0, 0, 1, 0}, { 0, 0, 0, 0}, {-1, 0, 0, 0}, { 0, 0, 0, 0}}, 1734 {{0, -1, 0, 0}, { 1, 0, 0, 0}, { 0, 0, 0, 0}, { 0, 0, 0, 0}}, 1735 {{0, 0, 0, 0}, { 0, 0, 0, 0}, { 0, 0, 0, 0}, { 0, 0, 0, 0}}}}; 1736 PetscReal normal[4]; 1737 PetscInt e, f, g; 1738 1739 for (d = 0; d < embedDim; ++d) { 1740 normal[d] = 0.0; 1741 for (e = 0; e < embedDim; ++e) { 1742 for (f = 0; f < embedDim; ++f) { 1743 for (g = 0; g < embedDim; ++g) { 1744 normal[d] += epsilon[d][e][f][g]*(coordsIn[j*embedDim+e] - coordsIn[i*embedDim+e])*(coordsIn[k*embedDim+f] - coordsIn[i*embedDim+f])*(coordsIn[l*embedDim+f] - coordsIn[i*embedDim+f]); 1745 } 1746 } 1747 } 1748 } 1749 if (DotReal(embedDim, normal, &coordsIn[i*embedDim]) < 0) {PetscInt tmp = cone[1]; cone[1] = cone[2]; cone[2] = tmp;} 1750 } 1751 ierr = DMPlexSetCone(*dm, c++, cone);CHKERRQ(ierr); 1752 } 1753 } 1754 } 1755 } 1756 } 1757 ierr = DMPlexSymmetrize(*dm);CHKERRQ(ierr); 1758 ierr = DMPlexStratify(*dm);CHKERRQ(ierr); 1759 ierr = PetscFree(graph);CHKERRQ(ierr); 1760 /* Interpolate mesh */ 1761 ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr); 1762 ierr = DMDestroy(dm);CHKERRQ(ierr); 1763 *dm = idm; 1764 break; 1765 } 1766 default: SETERRQ1(comm, PETSC_ERR_SUP, "Unsupported dimension for sphere: %D", dim); 1767 } 1768 /* Create coordinates */ 1769 ierr = DMGetCoordinateSection(*dm, &coordSection);CHKERRQ(ierr); 1770 ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr); 1771 ierr = PetscSectionSetFieldComponents(coordSection, 0, embedDim);CHKERRQ(ierr); 1772 ierr = PetscSectionSetChart(coordSection, firstVertex, firstVertex+numVerts);CHKERRQ(ierr); 1773 for (v = firstVertex; v < firstVertex+numVerts; ++v) { 1774 ierr = PetscSectionSetDof(coordSection, v, embedDim);CHKERRQ(ierr); 1775 ierr = PetscSectionSetFieldDof(coordSection, v, 0, embedDim);CHKERRQ(ierr); 1776 } 1777 ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 1778 ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 1779 ierr = VecCreate(PETSC_COMM_SELF, &coordinates);CHKERRQ(ierr); 1780 ierr = VecSetBlockSize(coordinates, embedDim);CHKERRQ(ierr); 1781 ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr); 1782 ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 1783 ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr); 1784 ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 1785 for (v = 0; v < numVerts; ++v) for (d = 0; d < embedDim; ++d) {coords[v*embedDim+d] = coordsIn[v*embedDim+d];} 1786 ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 1787 ierr = DMSetCoordinatesLocal(*dm, coordinates);CHKERRQ(ierr); 1788 ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 1789 ierr = PetscFree(coordsIn);CHKERRQ(ierr); 1790 PetscFunctionReturn(0); 1791 } 1792 1793 /* External function declarations here */ 1794 extern PetscErrorCode DMCreateInterpolation_Plex(DM dmCoarse, DM dmFine, Mat *interpolation, Vec *scaling); 1795 extern PetscErrorCode DMCreateInjection_Plex(DM dmCoarse, DM dmFine, Mat *mat); 1796 extern PetscErrorCode DMCreateDefaultSection_Plex(DM dm); 1797 extern PetscErrorCode DMCreateDefaultConstraints_Plex(DM dm); 1798 extern PetscErrorCode DMCreateMatrix_Plex(DM dm, Mat *J); 1799 extern PetscErrorCode DMCreateCoordinateDM_Plex(DM dm, DM *cdm); 1800 PETSC_INTERN PetscErrorCode DMClone_Plex(DM dm, DM *newdm); 1801 extern PetscErrorCode DMSetUp_Plex(DM dm); 1802 extern PetscErrorCode DMDestroy_Plex(DM dm); 1803 extern PetscErrorCode DMView_Plex(DM dm, PetscViewer viewer); 1804 extern PetscErrorCode DMLoad_Plex(DM dm, PetscViewer viewer); 1805 extern PetscErrorCode DMCreateSubDM_Plex(DM dm, PetscInt numFields, PetscInt fields[], IS *is, DM *subdm); 1806 static PetscErrorCode DMInitialize_Plex(DM dm); 1807 1808 /* Replace dm with the contents of dmNew 1809 - Share the DM_Plex structure 1810 - Share the coordinates 1811 - Share the SF 1812 */ 1813 static PetscErrorCode DMPlexReplace_Static(DM dm, DM dmNew) 1814 { 1815 PetscSF sf; 1816 DM coordDM, coarseDM; 1817 Vec coords; 1818 PetscBool isper; 1819 const PetscReal *maxCell, *L; 1820 const DMBoundaryType *bd; 1821 PetscErrorCode ierr; 1822 1823 PetscFunctionBegin; 1824 ierr = DMGetPointSF(dmNew, &sf);CHKERRQ(ierr); 1825 ierr = DMSetPointSF(dm, sf);CHKERRQ(ierr); 1826 ierr = DMGetCoordinateDM(dmNew, &coordDM);CHKERRQ(ierr); 1827 ierr = DMGetCoordinatesLocal(dmNew, &coords);CHKERRQ(ierr); 1828 ierr = DMSetCoordinateDM(dm, coordDM);CHKERRQ(ierr); 1829 ierr = DMSetCoordinatesLocal(dm, coords);CHKERRQ(ierr); 1830 ierr = DMGetPeriodicity(dm, &isper, &maxCell, &L, &bd);CHKERRQ(ierr); 1831 ierr = DMSetPeriodicity(dmNew, isper, maxCell, L, bd);CHKERRQ(ierr); 1832 ierr = DMDestroy_Plex(dm);CHKERRQ(ierr); 1833 ierr = DMInitialize_Plex(dm);CHKERRQ(ierr); 1834 dm->data = dmNew->data; 1835 ((DM_Plex *) dmNew->data)->refct++; 1836 dmNew->labels->refct++; 1837 if (!--(dm->labels->refct)) { 1838 DMLabelLink next = dm->labels->next; 1839 1840 /* destroy the labels */ 1841 while (next) { 1842 DMLabelLink tmp = next->next; 1843 1844 ierr = DMLabelDestroy(&next->label);CHKERRQ(ierr); 1845 ierr = PetscFree(next);CHKERRQ(ierr); 1846 next = tmp; 1847 } 1848 ierr = PetscFree(dm->labels);CHKERRQ(ierr); 1849 } 1850 dm->labels = dmNew->labels; 1851 dm->depthLabel = dmNew->depthLabel; 1852 ierr = DMGetCoarseDM(dmNew,&coarseDM);CHKERRQ(ierr); 1853 ierr = DMSetCoarseDM(dm,coarseDM);CHKERRQ(ierr); 1854 PetscFunctionReturn(0); 1855 } 1856 1857 /* Swap dm with the contents of dmNew 1858 - Swap the DM_Plex structure 1859 - Swap the coordinates 1860 - Swap the point PetscSF 1861 */ 1862 static PetscErrorCode DMPlexSwap_Static(DM dmA, DM dmB) 1863 { 1864 DM coordDMA, coordDMB; 1865 Vec coordsA, coordsB; 1866 PetscSF sfA, sfB; 1867 void *tmp; 1868 DMLabelLinkList listTmp; 1869 DMLabel depthTmp; 1870 PetscErrorCode ierr; 1871 1872 PetscFunctionBegin; 1873 ierr = DMGetPointSF(dmA, &sfA);CHKERRQ(ierr); 1874 ierr = DMGetPointSF(dmB, &sfB);CHKERRQ(ierr); 1875 ierr = PetscObjectReference((PetscObject) sfA);CHKERRQ(ierr); 1876 ierr = DMSetPointSF(dmA, sfB);CHKERRQ(ierr); 1877 ierr = DMSetPointSF(dmB, sfA);CHKERRQ(ierr); 1878 ierr = PetscObjectDereference((PetscObject) sfA);CHKERRQ(ierr); 1879 1880 ierr = DMGetCoordinateDM(dmA, &coordDMA);CHKERRQ(ierr); 1881 ierr = DMGetCoordinateDM(dmB, &coordDMB);CHKERRQ(ierr); 1882 ierr = PetscObjectReference((PetscObject) coordDMA);CHKERRQ(ierr); 1883 ierr = DMSetCoordinateDM(dmA, coordDMB);CHKERRQ(ierr); 1884 ierr = DMSetCoordinateDM(dmB, coordDMA);CHKERRQ(ierr); 1885 ierr = PetscObjectDereference((PetscObject) coordDMA);CHKERRQ(ierr); 1886 1887 ierr = DMGetCoordinatesLocal(dmA, &coordsA);CHKERRQ(ierr); 1888 ierr = DMGetCoordinatesLocal(dmB, &coordsB);CHKERRQ(ierr); 1889 ierr = PetscObjectReference((PetscObject) coordsA);CHKERRQ(ierr); 1890 ierr = DMSetCoordinatesLocal(dmA, coordsB);CHKERRQ(ierr); 1891 ierr = DMSetCoordinatesLocal(dmB, coordsA);CHKERRQ(ierr); 1892 ierr = PetscObjectDereference((PetscObject) coordsA);CHKERRQ(ierr); 1893 1894 tmp = dmA->data; 1895 dmA->data = dmB->data; 1896 dmB->data = tmp; 1897 listTmp = dmA->labels; 1898 dmA->labels = dmB->labels; 1899 dmB->labels = listTmp; 1900 depthTmp = dmA->depthLabel; 1901 dmA->depthLabel = dmB->depthLabel; 1902 dmB->depthLabel = depthTmp; 1903 PetscFunctionReturn(0); 1904 } 1905 1906 PetscErrorCode DMSetFromOptions_NonRefinement_Plex(PetscOptionItems *PetscOptionsObject,DM dm) 1907 { 1908 DM_Plex *mesh = (DM_Plex*) dm->data; 1909 PetscErrorCode ierr; 1910 1911 PetscFunctionBegin; 1912 /* Handle viewing */ 1913 ierr = PetscOptionsBool("-dm_plex_print_set_values", "Output all set values info", "DMView", PETSC_FALSE, &mesh->printSetValues, NULL);CHKERRQ(ierr); 1914 ierr = PetscOptionsInt("-dm_plex_print_fem", "Debug output level all fem computations", "DMView", 0, &mesh->printFEM, NULL);CHKERRQ(ierr); 1915 ierr = PetscOptionsReal("-dm_plex_print_tol", "Tolerance for FEM output", "DMView", mesh->printTol, &mesh->printTol, NULL);CHKERRQ(ierr); 1916 /* Point Location */ 1917 ierr = PetscOptionsBool("-dm_plex_hash_location", "Use grid hashing for point location", "DMView", PETSC_FALSE, &mesh->useHashLocation, NULL);CHKERRQ(ierr); 1918 /* Generation and remeshing */ 1919 ierr = PetscOptionsBool("-dm_plex_remesh_bd", "Allow changes to the boundary on remeshing", "DMView", PETSC_FALSE, &mesh->remeshBd, NULL);CHKERRQ(ierr); 1920 /* Projection behavior */ 1921 ierr = PetscOptionsInt("-dm_plex_max_projection_height", "Maxmimum mesh point height used to project locally", "DMPlexSetMaxProjectionHeight", 0, &mesh->maxProjectionHeight, NULL);CHKERRQ(ierr); 1922 ierr = PetscOptionsBool("-dm_plex_regular_refinement", "Use special nested projection algorithm for regular refinement", "DMPlexSetRegularRefinement", mesh->regularRefinement, &mesh->regularRefinement, NULL);CHKERRQ(ierr); 1923 1924 ierr = PetscPartitionerSetFromOptions(mesh->partitioner);CHKERRQ(ierr); 1925 PetscFunctionReturn(0); 1926 } 1927 1928 static PetscErrorCode DMSetFromOptions_Plex(PetscOptionItems *PetscOptionsObject,DM dm) 1929 { 1930 PetscInt refine = 0, coarsen = 0, r; 1931 PetscBool isHierarchy; 1932 PetscErrorCode ierr; 1933 1934 PetscFunctionBegin; 1935 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1936 ierr = PetscOptionsHead(PetscOptionsObject,"DMPlex Options");CHKERRQ(ierr); 1937 /* Handle DMPlex refinement */ 1938 ierr = PetscOptionsInt("-dm_refine", "The number of uniform refinements", "DMCreate", refine, &refine, NULL);CHKERRQ(ierr); 1939 ierr = PetscOptionsInt("-dm_refine_hierarchy", "The number of uniform refinements", "DMCreate", refine, &refine, &isHierarchy);CHKERRQ(ierr); 1940 if (refine) {ierr = DMPlexSetRefinementUniform(dm, PETSC_TRUE);CHKERRQ(ierr);} 1941 if (refine && isHierarchy) { 1942 DM *dms, coarseDM; 1943 1944 ierr = DMGetCoarseDM(dm, &coarseDM);CHKERRQ(ierr); 1945 ierr = PetscObjectReference((PetscObject)coarseDM);CHKERRQ(ierr); 1946 ierr = PetscMalloc1(refine,&dms);CHKERRQ(ierr); 1947 ierr = DMRefineHierarchy(dm, refine, dms);CHKERRQ(ierr); 1948 /* Total hack since we do not pass in a pointer */ 1949 ierr = DMPlexSwap_Static(dm, dms[refine-1]);CHKERRQ(ierr); 1950 if (refine == 1) { 1951 ierr = DMSetCoarseDM(dm, dms[0]);CHKERRQ(ierr); 1952 ierr = DMPlexSetRegularRefinement(dm, PETSC_TRUE);CHKERRQ(ierr); 1953 } else { 1954 ierr = DMSetCoarseDM(dm, dms[refine-2]);CHKERRQ(ierr); 1955 ierr = DMPlexSetRegularRefinement(dm, PETSC_TRUE);CHKERRQ(ierr); 1956 ierr = DMSetCoarseDM(dms[0], dms[refine-1]);CHKERRQ(ierr); 1957 ierr = DMPlexSetRegularRefinement(dms[0], PETSC_TRUE);CHKERRQ(ierr); 1958 } 1959 ierr = DMSetCoarseDM(dms[refine-1], coarseDM);CHKERRQ(ierr); 1960 ierr = PetscObjectDereference((PetscObject)coarseDM);CHKERRQ(ierr); 1961 /* Free DMs */ 1962 for (r = 0; r < refine; ++r) { 1963 ierr = DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dms[r]);CHKERRQ(ierr); 1964 ierr = DMDestroy(&dms[r]);CHKERRQ(ierr); 1965 } 1966 ierr = PetscFree(dms);CHKERRQ(ierr); 1967 } else { 1968 for (r = 0; r < refine; ++r) { 1969 DM refinedMesh; 1970 1971 ierr = DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dm);CHKERRQ(ierr); 1972 ierr = DMRefine(dm, PetscObjectComm((PetscObject) dm), &refinedMesh);CHKERRQ(ierr); 1973 /* Total hack since we do not pass in a pointer */ 1974 ierr = DMPlexReplace_Static(dm, refinedMesh);CHKERRQ(ierr); 1975 ierr = DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dm);CHKERRQ(ierr); 1976 ierr = DMDestroy(&refinedMesh);CHKERRQ(ierr); 1977 } 1978 } 1979 /* Handle DMPlex coarsening */ 1980 ierr = PetscOptionsInt("-dm_coarsen", "Coarsen the mesh", "DMCreate", coarsen, &coarsen, NULL);CHKERRQ(ierr); 1981 ierr = PetscOptionsInt("-dm_coarsen_hierarchy", "The number of coarsenings", "DMCreate", coarsen, &coarsen, &isHierarchy);CHKERRQ(ierr); 1982 if (coarsen && isHierarchy) { 1983 DM *dms; 1984 1985 ierr = PetscMalloc1(coarsen, &dms);CHKERRQ(ierr); 1986 ierr = DMCoarsenHierarchy(dm, coarsen, dms);CHKERRQ(ierr); 1987 /* Free DMs */ 1988 for (r = 0; r < coarsen; ++r) { 1989 ierr = DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dms[r]);CHKERRQ(ierr); 1990 ierr = DMDestroy(&dms[r]);CHKERRQ(ierr); 1991 } 1992 ierr = PetscFree(dms);CHKERRQ(ierr); 1993 } else { 1994 for (r = 0; r < coarsen; ++r) { 1995 DM coarseMesh; 1996 1997 ierr = DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dm);CHKERRQ(ierr); 1998 ierr = DMCoarsen(dm, PetscObjectComm((PetscObject) dm), &coarseMesh);CHKERRQ(ierr); 1999 /* Total hack since we do not pass in a pointer */ 2000 ierr = DMPlexReplace_Static(dm, coarseMesh);CHKERRQ(ierr); 2001 ierr = DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dm);CHKERRQ(ierr); 2002 ierr = DMDestroy(&coarseMesh);CHKERRQ(ierr); 2003 } 2004 } 2005 /* Handle */ 2006 ierr = DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dm);CHKERRQ(ierr); 2007 ierr = PetscOptionsTail();CHKERRQ(ierr); 2008 PetscFunctionReturn(0); 2009 } 2010 2011 static PetscErrorCode DMCreateGlobalVector_Plex(DM dm,Vec *vec) 2012 { 2013 PetscErrorCode ierr; 2014 2015 PetscFunctionBegin; 2016 ierr = DMCreateGlobalVector_Section_Private(dm,vec);CHKERRQ(ierr); 2017 /* ierr = VecSetOperation(*vec, VECOP_DUPLICATE, (void(*)(void)) VecDuplicate_MPI_DM);CHKERRQ(ierr); */ 2018 ierr = VecSetOperation(*vec, VECOP_VIEW, (void (*)(void)) VecView_Plex);CHKERRQ(ierr); 2019 ierr = VecSetOperation(*vec, VECOP_VIEWNATIVE, (void (*)(void)) VecView_Plex_Native);CHKERRQ(ierr); 2020 ierr = VecSetOperation(*vec, VECOP_LOAD, (void (*)(void)) VecLoad_Plex);CHKERRQ(ierr); 2021 ierr = VecSetOperation(*vec, VECOP_LOADNATIVE, (void (*)(void)) VecLoad_Plex_Native);CHKERRQ(ierr); 2022 PetscFunctionReturn(0); 2023 } 2024 2025 static PetscErrorCode DMCreateLocalVector_Plex(DM dm,Vec *vec) 2026 { 2027 PetscErrorCode ierr; 2028 2029 PetscFunctionBegin; 2030 ierr = DMCreateLocalVector_Section_Private(dm,vec);CHKERRQ(ierr); 2031 ierr = VecSetOperation(*vec, VECOP_VIEW, (void (*)(void)) VecView_Plex_Local);CHKERRQ(ierr); 2032 ierr = VecSetOperation(*vec, VECOP_LOAD, (void (*)(void)) VecLoad_Plex_Local);CHKERRQ(ierr); 2033 PetscFunctionReturn(0); 2034 } 2035 2036 static PetscErrorCode DMGetDimPoints_Plex(DM dm, PetscInt dim, PetscInt *pStart, PetscInt *pEnd) 2037 { 2038 PetscInt depth, d; 2039 PetscErrorCode ierr; 2040 2041 PetscFunctionBegin; 2042 ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 2043 if (depth == 1) { 2044 ierr = DMGetDimension(dm, &d);CHKERRQ(ierr); 2045 if (dim == 0) {ierr = DMPlexGetDepthStratum(dm, dim, pStart, pEnd);CHKERRQ(ierr);} 2046 else if (dim == d) {ierr = DMPlexGetDepthStratum(dm, 1, pStart, pEnd);CHKERRQ(ierr);} 2047 else {*pStart = 0; *pEnd = 0;} 2048 } else { 2049 ierr = DMPlexGetDepthStratum(dm, dim, pStart, pEnd);CHKERRQ(ierr); 2050 } 2051 PetscFunctionReturn(0); 2052 } 2053 2054 static PetscErrorCode DMGetNeighors_Plex(DM dm, PetscInt *nranks, const PetscMPIInt *ranks[]) 2055 { 2056 PetscSF sf; 2057 PetscErrorCode ierr; 2058 2059 PetscFunctionBegin; 2060 ierr = DMGetPointSF(dm, &sf);CHKERRQ(ierr); 2061 ierr = PetscSFGetRanks(sf, nranks, ranks, NULL, NULL, NULL);CHKERRQ(ierr); 2062 PetscFunctionReturn(0); 2063 } 2064 2065 static PetscErrorCode DMInitialize_Plex(DM dm) 2066 { 2067 PetscErrorCode ierr; 2068 2069 PetscFunctionBegin; 2070 dm->ops->view = DMView_Plex; 2071 dm->ops->load = DMLoad_Plex; 2072 dm->ops->setfromoptions = DMSetFromOptions_Plex; 2073 dm->ops->clone = DMClone_Plex; 2074 dm->ops->setup = DMSetUp_Plex; 2075 dm->ops->createdefaultsection = DMCreateDefaultSection_Plex; 2076 dm->ops->createdefaultconstraints = DMCreateDefaultConstraints_Plex; 2077 dm->ops->createglobalvector = DMCreateGlobalVector_Plex; 2078 dm->ops->createlocalvector = DMCreateLocalVector_Plex; 2079 dm->ops->getlocaltoglobalmapping = NULL; 2080 dm->ops->createfieldis = NULL; 2081 dm->ops->createcoordinatedm = DMCreateCoordinateDM_Plex; 2082 dm->ops->getcoloring = NULL; 2083 dm->ops->creatematrix = DMCreateMatrix_Plex; 2084 dm->ops->createinterpolation = DMCreateInterpolation_Plex; 2085 dm->ops->getaggregates = NULL; 2086 dm->ops->getinjection = DMCreateInjection_Plex; 2087 dm->ops->refine = DMRefine_Plex; 2088 dm->ops->coarsen = DMCoarsen_Plex; 2089 dm->ops->refinehierarchy = DMRefineHierarchy_Plex; 2090 dm->ops->coarsenhierarchy = DMCoarsenHierarchy_Plex; 2091 dm->ops->adaptlabel = DMAdaptLabel_Plex; 2092 dm->ops->adaptmetric = DMAdaptMetric_Plex; 2093 dm->ops->globaltolocalbegin = NULL; 2094 dm->ops->globaltolocalend = NULL; 2095 dm->ops->localtoglobalbegin = NULL; 2096 dm->ops->localtoglobalend = NULL; 2097 dm->ops->destroy = DMDestroy_Plex; 2098 dm->ops->createsubdm = DMCreateSubDM_Plex; 2099 dm->ops->getdimpoints = DMGetDimPoints_Plex; 2100 dm->ops->locatepoints = DMLocatePoints_Plex; 2101 dm->ops->projectfunctionlocal = DMProjectFunctionLocal_Plex; 2102 dm->ops->projectfunctionlabellocal = DMProjectFunctionLabelLocal_Plex; 2103 dm->ops->projectfieldlocal = DMProjectFieldLocal_Plex; 2104 dm->ops->projectfieldlabellocal = DMProjectFieldLabelLocal_Plex; 2105 dm->ops->computel2diff = DMComputeL2Diff_Plex; 2106 dm->ops->computel2gradientdiff = DMComputeL2GradientDiff_Plex; 2107 dm->ops->computel2fielddiff = DMComputeL2FieldDiff_Plex; 2108 dm->ops->getneighbors = DMGetNeighors_Plex; 2109 ierr = PetscObjectComposeFunction((PetscObject)dm,"DMPlexInsertBoundaryValues_C",DMPlexInsertBoundaryValues_Plex);CHKERRQ(ierr); 2110 ierr = PetscObjectComposeFunction((PetscObject)dm,"DMSetUpGLVisViewer_C",DMSetUpGLVisViewer_Plex);CHKERRQ(ierr); 2111 PetscFunctionReturn(0); 2112 } 2113 2114 PETSC_INTERN PetscErrorCode DMClone_Plex(DM dm, DM *newdm) 2115 { 2116 DM_Plex *mesh = (DM_Plex *) dm->data; 2117 PetscErrorCode ierr; 2118 2119 PetscFunctionBegin; 2120 mesh->refct++; 2121 (*newdm)->data = mesh; 2122 ierr = PetscObjectChangeTypeName((PetscObject) *newdm, DMPLEX);CHKERRQ(ierr); 2123 ierr = DMInitialize_Plex(*newdm);CHKERRQ(ierr); 2124 PetscFunctionReturn(0); 2125 } 2126 2127 /*MC 2128 DMPLEX = "plex" - A DM object that encapsulates an unstructured mesh, or CW Complex, which can be expressed using a Hasse Diagram. 2129 In the local representation, Vecs contain all unknowns in the interior and shared boundary. This is 2130 specified by a PetscSection object. Ownership in the global representation is determined by 2131 ownership of the underlying DMPlex points. This is specified by another PetscSection object. 2132 2133 Level: intermediate 2134 2135 .seealso: DMType, DMPlexCreate(), DMCreate(), DMSetType() 2136 M*/ 2137 2138 PETSC_EXTERN PetscErrorCode DMCreate_Plex(DM dm) 2139 { 2140 DM_Plex *mesh; 2141 PetscInt unit, d; 2142 PetscErrorCode ierr; 2143 2144 PetscFunctionBegin; 2145 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 2146 ierr = PetscNewLog(dm,&mesh);CHKERRQ(ierr); 2147 dm->dim = 0; 2148 dm->data = mesh; 2149 2150 mesh->refct = 1; 2151 ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm), &mesh->coneSection);CHKERRQ(ierr); 2152 mesh->maxConeSize = 0; 2153 mesh->cones = NULL; 2154 mesh->coneOrientations = NULL; 2155 ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm), &mesh->supportSection);CHKERRQ(ierr); 2156 mesh->maxSupportSize = 0; 2157 mesh->supports = NULL; 2158 mesh->refinementUniform = PETSC_TRUE; 2159 mesh->refinementLimit = -1.0; 2160 2161 mesh->facesTmp = NULL; 2162 2163 mesh->tetgenOpts = NULL; 2164 mesh->triangleOpts = NULL; 2165 ierr = PetscPartitionerCreate(PetscObjectComm((PetscObject)dm), &mesh->partitioner);CHKERRQ(ierr); 2166 mesh->remeshBd = PETSC_FALSE; 2167 2168 mesh->subpointMap = NULL; 2169 2170 for (unit = 0; unit < NUM_PETSC_UNITS; ++unit) mesh->scale[unit] = 1.0; 2171 2172 mesh->regularRefinement = PETSC_FALSE; 2173 mesh->depthState = -1; 2174 mesh->globalVertexNumbers = NULL; 2175 mesh->globalCellNumbers = NULL; 2176 mesh->anchorSection = NULL; 2177 mesh->anchorIS = NULL; 2178 mesh->createanchors = NULL; 2179 mesh->computeanchormatrix = NULL; 2180 mesh->parentSection = NULL; 2181 mesh->parents = NULL; 2182 mesh->childIDs = NULL; 2183 mesh->childSection = NULL; 2184 mesh->children = NULL; 2185 mesh->referenceTree = NULL; 2186 mesh->getchildsymmetry = NULL; 2187 for (d = 0; d < 8; ++d) mesh->hybridPointMax[d] = PETSC_DETERMINE; 2188 mesh->vtkCellHeight = 0; 2189 mesh->useCone = PETSC_FALSE; 2190 mesh->useClosure = PETSC_TRUE; 2191 mesh->useAnchors = PETSC_FALSE; 2192 2193 mesh->maxProjectionHeight = 0; 2194 2195 mesh->printSetValues = PETSC_FALSE; 2196 mesh->printFEM = 0; 2197 mesh->printTol = 1.0e-10; 2198 2199 ierr = DMInitialize_Plex(dm);CHKERRQ(ierr); 2200 PetscFunctionReturn(0); 2201 } 2202 2203 /*@ 2204 DMPlexCreate - Creates a DMPlex object, which encapsulates an unstructured mesh, or CW complex, which can be expressed using a Hasse Diagram. 2205 2206 Collective on MPI_Comm 2207 2208 Input Parameter: 2209 . comm - The communicator for the DMPlex object 2210 2211 Output Parameter: 2212 . mesh - The DMPlex object 2213 2214 Level: beginner 2215 2216 .keywords: DMPlex, create 2217 @*/ 2218 PetscErrorCode DMPlexCreate(MPI_Comm comm, DM *mesh) 2219 { 2220 PetscErrorCode ierr; 2221 2222 PetscFunctionBegin; 2223 PetscValidPointer(mesh,2); 2224 ierr = DMCreate(comm, mesh);CHKERRQ(ierr); 2225 ierr = DMSetType(*mesh, DMPLEX);CHKERRQ(ierr); 2226 PetscFunctionReturn(0); 2227 } 2228 2229 /* 2230 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 2231 */ 2232 static PetscErrorCode DMPlexBuildFromCellList_Parallel_Private(DM dm, PetscInt numCells, PetscInt numVertices, PetscInt numCorners, const int cells[], PetscSF *sfVert) 2233 { 2234 PetscSF sfPoint; 2235 PetscLayout vLayout; 2236 PetscHashI vhash; 2237 PetscSFNode *remoteVerticesAdj, *vertexLocal, *vertexOwner, *remoteVertex; 2238 const PetscInt *vrange; 2239 PetscInt numVerticesAdj, off, *verticesAdj, numVerticesGhost = 0, *localVertex, *cone, c, p, v, g; 2240 PETSC_UNUSED PetscHashIIter ret, iter; 2241 PetscMPIInt rank, size; 2242 PetscErrorCode ierr; 2243 2244 PetscFunctionBegin; 2245 ierr = MPI_Comm_rank(PetscObjectComm((PetscObject) dm), &rank);CHKERRQ(ierr); 2246 ierr = MPI_Comm_size(PetscObjectComm((PetscObject) dm), &size);CHKERRQ(ierr); 2247 /* Partition vertices */ 2248 ierr = PetscLayoutCreate(PetscObjectComm((PetscObject) dm), &vLayout);CHKERRQ(ierr); 2249 ierr = PetscLayoutSetLocalSize(vLayout, numVertices);CHKERRQ(ierr); 2250 ierr = PetscLayoutSetBlockSize(vLayout, 1);CHKERRQ(ierr); 2251 ierr = PetscLayoutSetUp(vLayout);CHKERRQ(ierr); 2252 ierr = PetscLayoutGetRanges(vLayout, &vrange);CHKERRQ(ierr); 2253 /* Count vertices and map them to procs */ 2254 PetscHashICreate(vhash); 2255 for (c = 0; c < numCells; ++c) { 2256 for (p = 0; p < numCorners; ++p) { 2257 PetscHashIPut(vhash, cells[c*numCorners+p], ret, iter); 2258 } 2259 } 2260 PetscHashISize(vhash, numVerticesAdj); 2261 ierr = PetscMalloc1(numVerticesAdj, &verticesAdj);CHKERRQ(ierr); 2262 off = 0; ierr = PetscHashIGetKeys(vhash, &off, verticesAdj);CHKERRQ(ierr); 2263 if (off != numVerticesAdj) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid number of local vertices %D should be %D", off, numVerticesAdj); 2264 ierr = PetscSortInt(numVerticesAdj, verticesAdj);CHKERRQ(ierr); 2265 ierr = PetscMalloc1(numVerticesAdj, &remoteVerticesAdj);CHKERRQ(ierr); 2266 for (v = 0; v < numVerticesAdj; ++v) { 2267 const PetscInt gv = verticesAdj[v]; 2268 PetscInt vrank; 2269 2270 ierr = PetscFindInt(gv, size+1, vrange, &vrank);CHKERRQ(ierr); 2271 vrank = vrank < 0 ? -(vrank+2) : vrank; 2272 remoteVerticesAdj[v].index = gv - vrange[vrank]; 2273 remoteVerticesAdj[v].rank = vrank; 2274 } 2275 /* Create cones */ 2276 ierr = DMPlexSetChart(dm, 0, numCells+numVerticesAdj);CHKERRQ(ierr); 2277 for (c = 0; c < numCells; ++c) {ierr = DMPlexSetConeSize(dm, c, numCorners);CHKERRQ(ierr);} 2278 ierr = DMSetUp(dm);CHKERRQ(ierr); 2279 ierr = DMGetWorkArray(dm, numCorners, PETSC_INT, &cone);CHKERRQ(ierr); 2280 for (c = 0; c < numCells; ++c) { 2281 for (p = 0; p < numCorners; ++p) { 2282 const PetscInt gv = cells[c*numCorners+p]; 2283 PetscInt lv; 2284 2285 ierr = PetscFindInt(gv, numVerticesAdj, verticesAdj, &lv);CHKERRQ(ierr); 2286 if (lv < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Could not find global vertex %D in local connectivity", gv); 2287 cone[p] = lv+numCells; 2288 } 2289 ierr = DMPlexSetCone(dm, c, cone);CHKERRQ(ierr); 2290 } 2291 ierr = DMRestoreWorkArray(dm, numCorners, PETSC_INT, &cone);CHKERRQ(ierr); 2292 /* Create SF for vertices */ 2293 ierr = PetscSFCreate(PetscObjectComm((PetscObject)dm), sfVert);CHKERRQ(ierr); 2294 ierr = PetscObjectSetName((PetscObject) *sfVert, "Vertex Ownership SF");CHKERRQ(ierr); 2295 ierr = PetscSFSetFromOptions(*sfVert);CHKERRQ(ierr); 2296 ierr = PetscSFSetGraph(*sfVert, numVertices, numVerticesAdj, NULL, PETSC_OWN_POINTER, remoteVerticesAdj, PETSC_OWN_POINTER);CHKERRQ(ierr); 2297 ierr = PetscFree(verticesAdj);CHKERRQ(ierr); 2298 /* Build pointSF */ 2299 ierr = PetscMalloc2(numVerticesAdj, &vertexLocal, numVertices, &vertexOwner);CHKERRQ(ierr); 2300 for (v = 0; v < numVerticesAdj; ++v) {vertexLocal[v].index = v+numCells; vertexLocal[v].rank = rank;} 2301 for (v = 0; v < numVertices; ++v) {vertexOwner[v].index = -1; vertexOwner[v].rank = -1;} 2302 ierr = PetscSFReduceBegin(*sfVert, MPIU_2INT, vertexLocal, vertexOwner, MPI_MAXLOC);CHKERRQ(ierr); 2303 ierr = PetscSFReduceEnd(*sfVert, MPIU_2INT, vertexLocal, vertexOwner, MPI_MAXLOC);CHKERRQ(ierr); 2304 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); 2305 ierr = PetscSFBcastBegin(*sfVert, MPIU_2INT, vertexOwner, vertexLocal);CHKERRQ(ierr); 2306 ierr = PetscSFBcastEnd(*sfVert, MPIU_2INT, vertexOwner, vertexLocal);CHKERRQ(ierr); 2307 for (v = 0; v < numVerticesAdj; ++v) if (vertexLocal[v].rank != rank) ++numVerticesGhost; 2308 ierr = PetscMalloc1(numVerticesGhost, &localVertex);CHKERRQ(ierr); 2309 ierr = PetscMalloc1(numVerticesGhost, &remoteVertex);CHKERRQ(ierr); 2310 for (v = 0, g = 0; v < numVerticesAdj; ++v) { 2311 if (vertexLocal[v].rank != rank) { 2312 localVertex[g] = v+numCells; 2313 remoteVertex[g].index = vertexLocal[v].index; 2314 remoteVertex[g].rank = vertexLocal[v].rank; 2315 ++g; 2316 } 2317 } 2318 ierr = PetscFree2(vertexLocal, vertexOwner);CHKERRQ(ierr); 2319 if (g != numVerticesGhost) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid number of vertex ghosts %D should be %D", g, numVerticesGhost); 2320 ierr = DMGetPointSF(dm, &sfPoint);CHKERRQ(ierr); 2321 ierr = PetscObjectSetName((PetscObject) sfPoint, "point SF");CHKERRQ(ierr); 2322 ierr = PetscSFSetGraph(sfPoint, numCells+numVerticesAdj, numVerticesGhost, localVertex, PETSC_OWN_POINTER, remoteVertex, PETSC_OWN_POINTER);CHKERRQ(ierr); 2323 ierr = PetscLayoutDestroy(&vLayout);CHKERRQ(ierr); 2324 PetscHashIDestroy(vhash); 2325 /* Fill in the rest of the topology structure */ 2326 ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr); 2327 ierr = DMPlexStratify(dm);CHKERRQ(ierr); 2328 PetscFunctionReturn(0); 2329 } 2330 2331 /* 2332 This takes as input the coordinates for each owned vertex 2333 */ 2334 static PetscErrorCode DMPlexBuildCoordinates_Parallel_Private(DM dm, PetscInt spaceDim, PetscInt numCells, PetscInt numV, PetscSF sfVert, const PetscReal vertexCoords[]) 2335 { 2336 PetscSection coordSection; 2337 Vec coordinates; 2338 PetscScalar *coords; 2339 PetscInt numVertices, numVerticesAdj, coordSize, v; 2340 PetscErrorCode ierr; 2341 2342 PetscFunctionBegin; 2343 ierr = DMSetCoordinateDim(dm, spaceDim);CHKERRQ(ierr); 2344 ierr = PetscSFGetGraph(sfVert, &numVertices, &numVerticesAdj, NULL, NULL);CHKERRQ(ierr); 2345 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 2346 ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr); 2347 ierr = PetscSectionSetFieldComponents(coordSection, 0, spaceDim);CHKERRQ(ierr); 2348 ierr = PetscSectionSetChart(coordSection, numCells, numCells + numVerticesAdj);CHKERRQ(ierr); 2349 for (v = numCells; v < numCells+numVerticesAdj; ++v) { 2350 ierr = PetscSectionSetDof(coordSection, v, spaceDim);CHKERRQ(ierr); 2351 ierr = PetscSectionSetFieldDof(coordSection, v, 0, spaceDim);CHKERRQ(ierr); 2352 } 2353 ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 2354 ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 2355 ierr = VecCreate(PetscObjectComm((PetscObject)dm), &coordinates);CHKERRQ(ierr); 2356 ierr = VecSetBlockSize(coordinates, spaceDim);CHKERRQ(ierr); 2357 ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr); 2358 ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 2359 ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr); 2360 ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 2361 { 2362 MPI_Datatype coordtype; 2363 2364 /* Need a temp buffer for coords if we have complex/single */ 2365 ierr = MPI_Type_contiguous(spaceDim, MPIU_SCALAR, &coordtype);CHKERRQ(ierr); 2366 ierr = MPI_Type_commit(&coordtype);CHKERRQ(ierr); 2367 #if defined(PETSC_USE_COMPLEX) 2368 { 2369 PetscScalar *svertexCoords; 2370 PetscInt i; 2371 ierr = PetscMalloc1(numV*spaceDim,&svertexCoords);CHKERRQ(ierr); 2372 for (i=0; i<numV*spaceDim; i++) svertexCoords[i] = vertexCoords[i]; 2373 ierr = PetscSFBcastBegin(sfVert, coordtype, svertexCoords, coords);CHKERRQ(ierr); 2374 ierr = PetscSFBcastEnd(sfVert, coordtype, svertexCoords, coords);CHKERRQ(ierr); 2375 ierr = PetscFree(svertexCoords);CHKERRQ(ierr); 2376 } 2377 #else 2378 ierr = PetscSFBcastBegin(sfVert, coordtype, vertexCoords, coords);CHKERRQ(ierr); 2379 ierr = PetscSFBcastEnd(sfVert, coordtype, vertexCoords, coords);CHKERRQ(ierr); 2380 #endif 2381 ierr = MPI_Type_free(&coordtype);CHKERRQ(ierr); 2382 } 2383 ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 2384 ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr); 2385 ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 2386 PetscFunctionReturn(0); 2387 } 2388 2389 /*@C 2390 DMPlexCreateFromCellListParallel - This takes as input common mesh generator output, a list of the vertices for each cell, and produces a DM 2391 2392 Input Parameters: 2393 + comm - The communicator 2394 . dim - The topological dimension of the mesh 2395 . numCells - The number of cells owned by this process 2396 . numVertices - The number of vertices owned by this process 2397 . numCorners - The number of vertices for each cell 2398 . interpolate - Flag indicating that intermediate mesh entities (faces, edges) should be created automatically 2399 . cells - An array of numCells*numCorners numbers, the global vertex numbers for each cell 2400 . spaceDim - The spatial dimension used for coordinates 2401 - vertexCoords - An array of numVertices*spaceDim numbers, the coordinates of each vertex 2402 2403 Output Parameter: 2404 + dm - The DM 2405 - vertexSF - Optional, SF describing complete vertex ownership 2406 2407 Note: Two triangles sharing a face 2408 $ 2409 $ 2 2410 $ / | \ 2411 $ / | \ 2412 $ / | \ 2413 $ 0 0 | 1 3 2414 $ \ | / 2415 $ \ | / 2416 $ \ | / 2417 $ 1 2418 would have input 2419 $ numCells = 2, numVertices = 4 2420 $ cells = [0 1 2 1 3 2] 2421 $ 2422 which would result in the DMPlex 2423 $ 2424 $ 4 2425 $ / | \ 2426 $ / | \ 2427 $ / | \ 2428 $ 2 0 | 1 5 2429 $ \ | / 2430 $ \ | / 2431 $ \ | / 2432 $ 3 2433 2434 Level: beginner 2435 2436 .seealso: DMPlexCreateFromCellList(), DMPlexCreateFromDAG(), DMPlexCreate() 2437 @*/ 2438 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) 2439 { 2440 PetscSF sfVert; 2441 PetscErrorCode ierr; 2442 2443 PetscFunctionBegin; 2444 ierr = DMCreate(comm, dm);CHKERRQ(ierr); 2445 ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 2446 PetscValidLogicalCollectiveInt(*dm, dim, 2); 2447 PetscValidLogicalCollectiveInt(*dm, spaceDim, 8); 2448 ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr); 2449 ierr = DMPlexBuildFromCellList_Parallel_Private(*dm, numCells, numVertices, numCorners, cells, &sfVert);CHKERRQ(ierr); 2450 if (interpolate) { 2451 DM idm = NULL; 2452 2453 ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr); 2454 ierr = DMDestroy(dm);CHKERRQ(ierr); 2455 *dm = idm; 2456 } 2457 ierr = DMPlexBuildCoordinates_Parallel_Private(*dm, spaceDim, numCells, numVertices,sfVert, vertexCoords);CHKERRQ(ierr); 2458 if (vertexSF) *vertexSF = sfVert; 2459 else {ierr = PetscSFDestroy(&sfVert);CHKERRQ(ierr);} 2460 PetscFunctionReturn(0); 2461 } 2462 2463 /* 2464 This takes as input the common mesh generator output, a list of the vertices for each cell 2465 */ 2466 static PetscErrorCode DMPlexBuildFromCellList_Private(DM dm, PetscInt numCells, PetscInt numVertices, PetscInt numCorners, const int cells[]) 2467 { 2468 PetscInt *cone, c, p; 2469 PetscErrorCode ierr; 2470 2471 PetscFunctionBegin; 2472 ierr = DMPlexSetChart(dm, 0, numCells+numVertices);CHKERRQ(ierr); 2473 for (c = 0; c < numCells; ++c) { 2474 ierr = DMPlexSetConeSize(dm, c, numCorners);CHKERRQ(ierr); 2475 } 2476 ierr = DMSetUp(dm);CHKERRQ(ierr); 2477 ierr = DMGetWorkArray(dm, numCorners, PETSC_INT, &cone);CHKERRQ(ierr); 2478 for (c = 0; c < numCells; ++c) { 2479 for (p = 0; p < numCorners; ++p) { 2480 cone[p] = cells[c*numCorners+p]+numCells; 2481 } 2482 ierr = DMPlexSetCone(dm, c, cone);CHKERRQ(ierr); 2483 } 2484 ierr = DMRestoreWorkArray(dm, numCorners, PETSC_INT, &cone);CHKERRQ(ierr); 2485 ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr); 2486 ierr = DMPlexStratify(dm);CHKERRQ(ierr); 2487 PetscFunctionReturn(0); 2488 } 2489 2490 /* 2491 This takes as input the coordinates for each vertex 2492 */ 2493 static PetscErrorCode DMPlexBuildCoordinates_Private(DM dm, PetscInt spaceDim, PetscInt numCells, PetscInt numVertices, const double vertexCoords[]) 2494 { 2495 PetscSection coordSection; 2496 Vec coordinates; 2497 DM cdm; 2498 PetscScalar *coords; 2499 PetscInt v, d; 2500 PetscErrorCode ierr; 2501 2502 PetscFunctionBegin; 2503 ierr = DMSetCoordinateDim(dm, spaceDim);CHKERRQ(ierr); 2504 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 2505 ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr); 2506 ierr = PetscSectionSetFieldComponents(coordSection, 0, spaceDim);CHKERRQ(ierr); 2507 ierr = PetscSectionSetChart(coordSection, numCells, numCells + numVertices);CHKERRQ(ierr); 2508 for (v = numCells; v < numCells+numVertices; ++v) { 2509 ierr = PetscSectionSetDof(coordSection, v, spaceDim);CHKERRQ(ierr); 2510 ierr = PetscSectionSetFieldDof(coordSection, v, 0, spaceDim);CHKERRQ(ierr); 2511 } 2512 ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 2513 2514 ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr); 2515 ierr = DMCreateLocalVector(cdm, &coordinates);CHKERRQ(ierr); 2516 ierr = VecSetBlockSize(coordinates, spaceDim);CHKERRQ(ierr); 2517 ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr); 2518 ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 2519 for (v = 0; v < numVertices; ++v) { 2520 for (d = 0; d < spaceDim; ++d) { 2521 coords[v*spaceDim+d] = vertexCoords[v*spaceDim+d]; 2522 } 2523 } 2524 ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 2525 ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr); 2526 ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 2527 PetscFunctionReturn(0); 2528 } 2529 2530 /*@C 2531 DMPlexCreateFromCellList - This takes as input common mesh generator output, a list of the vertices for each cell, and produces a DM 2532 2533 Input Parameters: 2534 + comm - The communicator 2535 . dim - The topological dimension of the mesh 2536 . numCells - The number of cells 2537 . numVertices - The number of vertices 2538 . numCorners - The number of vertices for each cell 2539 . interpolate - Flag indicating that intermediate mesh entities (faces, edges) should be created automatically 2540 . cells - An array of numCells*numCorners numbers, the vertices for each cell 2541 . spaceDim - The spatial dimension used for coordinates 2542 - vertexCoords - An array of numVertices*spaceDim numbers, the coordinates of each vertex 2543 2544 Output Parameter: 2545 . dm - The DM 2546 2547 Note: Two triangles sharing a face 2548 $ 2549 $ 2 2550 $ / | \ 2551 $ / | \ 2552 $ / | \ 2553 $ 0 0 | 1 3 2554 $ \ | / 2555 $ \ | / 2556 $ \ | / 2557 $ 1 2558 would have input 2559 $ numCells = 2, numVertices = 4 2560 $ cells = [0 1 2 1 3 2] 2561 $ 2562 which would result in the DMPlex 2563 $ 2564 $ 4 2565 $ / | \ 2566 $ / | \ 2567 $ / | \ 2568 $ 2 0 | 1 5 2569 $ \ | / 2570 $ \ | / 2571 $ \ | / 2572 $ 3 2573 2574 Level: beginner 2575 2576 .seealso: DMPlexCreateFromDAG(), DMPlexCreate() 2577 @*/ 2578 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) 2579 { 2580 PetscErrorCode ierr; 2581 2582 PetscFunctionBegin; 2583 ierr = DMCreate(comm, dm);CHKERRQ(ierr); 2584 ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 2585 ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr); 2586 ierr = DMPlexBuildFromCellList_Private(*dm, numCells, numVertices, numCorners, cells);CHKERRQ(ierr); 2587 if (interpolate) { 2588 DM idm = NULL; 2589 2590 ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr); 2591 ierr = DMDestroy(dm);CHKERRQ(ierr); 2592 *dm = idm; 2593 } 2594 ierr = DMPlexBuildCoordinates_Private(*dm, spaceDim, numCells, numVertices, vertexCoords);CHKERRQ(ierr); 2595 PetscFunctionReturn(0); 2596 } 2597 2598 /*@ 2599 DMPlexCreateFromDAG - This takes as input the adjacency-list representation of the Directed Acyclic Graph (Hasse Diagram) encoding a mesh, and produces a DM 2600 2601 Input Parameters: 2602 + dm - The empty DM object, usually from DMCreate() and DMSetDimension() 2603 . depth - The depth of the DAG 2604 . numPoints - The number of points at each depth 2605 . coneSize - The cone size of each point 2606 . cones - The concatenation of the cone points for each point, the cone list must be oriented correctly for each point 2607 . coneOrientations - The orientation of each cone point 2608 - vertexCoords - An array of numVertices*dim numbers, the coordinates of each vertex 2609 2610 Output Parameter: 2611 . dm - The DM 2612 2613 Note: Two triangles sharing a face would have input 2614 $ depth = 1, numPoints = [4 2], coneSize = [3 3 0 0 0 0] 2615 $ cones = [2 3 4 3 5 4], coneOrientations = [0 0 0 0 0 0] 2616 $ vertexCoords = [-1.0 0.0 0.0 -1.0 0.0 1.0 1.0 0.0] 2617 $ 2618 which would result in the DMPlex 2619 $ 2620 $ 4 2621 $ / | \ 2622 $ / | \ 2623 $ / | \ 2624 $ 2 0 | 1 5 2625 $ \ | / 2626 $ \ | / 2627 $ \ | / 2628 $ 3 2629 $ 2630 $ Notice that all points are numbered consecutively, unlikely DMPlexCreateFromCellList() 2631 2632 Level: advanced 2633 2634 .seealso: DMPlexCreateFromCellList(), DMPlexCreate() 2635 @*/ 2636 PetscErrorCode DMPlexCreateFromDAG(DM dm, PetscInt depth, const PetscInt numPoints[], const PetscInt coneSize[], const PetscInt cones[], const PetscInt coneOrientations[], const PetscScalar vertexCoords[]) 2637 { 2638 Vec coordinates; 2639 PetscSection coordSection; 2640 PetscScalar *coords; 2641 PetscInt coordSize, firstVertex = -1, pStart = 0, pEnd = 0, p, v, dim, dimEmbed, d, off; 2642 PetscErrorCode ierr; 2643 2644 PetscFunctionBegin; 2645 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 2646 ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr); 2647 if (dimEmbed < dim) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Embedding dimension %d cannot be less than intrinsic dimension %d",dimEmbed,dim); 2648 for (d = 0; d <= depth; ++d) pEnd += numPoints[d]; 2649 ierr = DMPlexSetChart(dm, pStart, pEnd);CHKERRQ(ierr); 2650 for (p = pStart; p < pEnd; ++p) { 2651 ierr = DMPlexSetConeSize(dm, p, coneSize[p-pStart]);CHKERRQ(ierr); 2652 if (firstVertex < 0 && !coneSize[p - pStart]) { 2653 firstVertex = p - pStart; 2654 } 2655 } 2656 if (firstVertex < 0 && numPoints[0]) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Expected %d vertices but could not find any", numPoints[0]); 2657 ierr = DMSetUp(dm);CHKERRQ(ierr); /* Allocate space for cones */ 2658 for (p = pStart, off = 0; p < pEnd; off += coneSize[p-pStart], ++p) { 2659 ierr = DMPlexSetCone(dm, p, &cones[off]);CHKERRQ(ierr); 2660 ierr = DMPlexSetConeOrientation(dm, p, &coneOrientations[off]);CHKERRQ(ierr); 2661 } 2662 ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr); 2663 ierr = DMPlexStratify(dm);CHKERRQ(ierr); 2664 /* Build coordinates */ 2665 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 2666 ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr); 2667 ierr = PetscSectionSetFieldComponents(coordSection, 0, dimEmbed);CHKERRQ(ierr); 2668 ierr = PetscSectionSetChart(coordSection, firstVertex, firstVertex+numPoints[0]);CHKERRQ(ierr); 2669 for (v = firstVertex; v < firstVertex+numPoints[0]; ++v) { 2670 ierr = PetscSectionSetDof(coordSection, v, dimEmbed);CHKERRQ(ierr); 2671 ierr = PetscSectionSetFieldDof(coordSection, v, 0, dimEmbed);CHKERRQ(ierr); 2672 } 2673 ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 2674 ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 2675 ierr = VecCreate(PETSC_COMM_SELF, &coordinates);CHKERRQ(ierr); 2676 ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr); 2677 ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 2678 ierr = VecSetBlockSize(coordinates, dimEmbed);CHKERRQ(ierr); 2679 ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr); 2680 ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 2681 for (v = 0; v < numPoints[0]; ++v) { 2682 PetscInt off; 2683 2684 ierr = PetscSectionGetOffset(coordSection, v+firstVertex, &off);CHKERRQ(ierr); 2685 for (d = 0; d < dimEmbed; ++d) { 2686 coords[off+d] = vertexCoords[v*dimEmbed+d]; 2687 } 2688 } 2689 ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 2690 ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr); 2691 ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 2692 PetscFunctionReturn(0); 2693 } 2694 2695 /*@C 2696 DMPlexCreateFromFile - This takes a filename and produces a DM 2697 2698 Input Parameters: 2699 + comm - The communicator 2700 . filename - A file name 2701 - interpolate - Flag to create intermediate mesh pieces (edges, faces) 2702 2703 Output Parameter: 2704 . dm - The DM 2705 2706 Level: beginner 2707 2708 .seealso: DMPlexCreateFromDAG(), DMPlexCreateFromCellList(), DMPlexCreate() 2709 @*/ 2710 PetscErrorCode DMPlexCreateFromFile(MPI_Comm comm, const char filename[], PetscBool interpolate, DM *dm) 2711 { 2712 const char *extGmsh = ".msh"; 2713 const char *extCGNS = ".cgns"; 2714 const char *extExodus = ".exo"; 2715 const char *extGenesis = ".gen"; 2716 const char *extFluent = ".cas"; 2717 const char *extHDF5 = ".h5"; 2718 const char *extMed = ".med"; 2719 const char *extPLY = ".ply"; 2720 size_t len; 2721 PetscBool isGmsh, isCGNS, isExodus, isGenesis, isFluent, isHDF5, isMed, isPLY; 2722 PetscMPIInt rank; 2723 PetscErrorCode ierr; 2724 2725 PetscFunctionBegin; 2726 PetscValidPointer(filename, 2); 2727 PetscValidPointer(dm, 4); 2728 ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 2729 ierr = PetscStrlen(filename, &len);CHKERRQ(ierr); 2730 if (!len) SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Filename must be a valid path"); 2731 ierr = PetscStrncmp(&filename[PetscMax(0,len-4)], extGmsh, 4, &isGmsh);CHKERRQ(ierr); 2732 ierr = PetscStrncmp(&filename[PetscMax(0,len-5)], extCGNS, 5, &isCGNS);CHKERRQ(ierr); 2733 ierr = PetscStrncmp(&filename[PetscMax(0,len-4)], extExodus, 4, &isExodus);CHKERRQ(ierr); 2734 ierr = PetscStrncmp(&filename[PetscMax(0,len-4)], extGenesis, 4, &isGenesis);CHKERRQ(ierr); 2735 ierr = PetscStrncmp(&filename[PetscMax(0,len-4)], extFluent, 4, &isFluent);CHKERRQ(ierr); 2736 ierr = PetscStrncmp(&filename[PetscMax(0,len-3)], extHDF5, 3, &isHDF5);CHKERRQ(ierr); 2737 ierr = PetscStrncmp(&filename[PetscMax(0,len-4)], extMed, 4, &isMed);CHKERRQ(ierr); 2738 ierr = PetscStrncmp(&filename[PetscMax(0,len-4)], extPLY, 4, &isPLY);CHKERRQ(ierr); 2739 if (isGmsh) { 2740 ierr = DMPlexCreateGmshFromFile(comm, filename, interpolate, dm);CHKERRQ(ierr); 2741 } else if (isCGNS) { 2742 ierr = DMPlexCreateCGNSFromFile(comm, filename, interpolate, dm);CHKERRQ(ierr); 2743 } else if (isExodus || isGenesis) { 2744 ierr = DMPlexCreateExodusFromFile(comm, filename, interpolate, dm);CHKERRQ(ierr); 2745 } else if (isFluent) { 2746 ierr = DMPlexCreateFluentFromFile(comm, filename, interpolate, dm);CHKERRQ(ierr); 2747 } else if (isHDF5) { 2748 PetscViewer viewer; 2749 2750 ierr = PetscViewerCreate(comm, &viewer);CHKERRQ(ierr); 2751 ierr = PetscViewerSetType(viewer, PETSCVIEWERHDF5);CHKERRQ(ierr); 2752 ierr = PetscViewerFileSetMode(viewer, FILE_MODE_READ);CHKERRQ(ierr); 2753 ierr = PetscViewerFileSetName(viewer, filename);CHKERRQ(ierr); 2754 ierr = DMCreate(comm, dm);CHKERRQ(ierr); 2755 ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 2756 ierr = DMLoad(*dm, viewer);CHKERRQ(ierr); 2757 ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 2758 } else if (isMed) { 2759 ierr = DMPlexCreateMedFromFile(comm, filename, interpolate, dm);CHKERRQ(ierr); 2760 } else if (isPLY) { 2761 ierr = DMPlexCreatePLYFromFile(comm, filename, interpolate, dm);CHKERRQ(ierr); 2762 } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cannot load file %s: unrecognized extension", filename); 2763 PetscFunctionReturn(0); 2764 } 2765 2766 /*@ 2767 DMPlexCreateReferenceCell - Create a DMPLEX with the appropriate FEM reference cell 2768 2769 Collective on comm 2770 2771 Input Parameters: 2772 + comm - The communicator 2773 . dim - The spatial dimension 2774 - simplex - Flag for simplex, otherwise use a tensor-product cell 2775 2776 Output Parameter: 2777 . refdm - The reference cell 2778 2779 Level: intermediate 2780 2781 .keywords: reference cell 2782 .seealso: 2783 @*/ 2784 PetscErrorCode DMPlexCreateReferenceCell(MPI_Comm comm, PetscInt dim, PetscBool simplex, DM *refdm) 2785 { 2786 DM rdm; 2787 Vec coords; 2788 PetscErrorCode ierr; 2789 2790 PetscFunctionBeginUser; 2791 ierr = DMCreate(comm, &rdm);CHKERRQ(ierr); 2792 ierr = DMSetType(rdm, DMPLEX);CHKERRQ(ierr); 2793 ierr = DMSetDimension(rdm, dim);CHKERRQ(ierr); 2794 switch (dim) { 2795 case 0: 2796 { 2797 PetscInt numPoints[1] = {1}; 2798 PetscInt coneSize[1] = {0}; 2799 PetscInt cones[1] = {0}; 2800 PetscInt coneOrientations[1] = {0}; 2801 PetscScalar vertexCoords[1] = {0.0}; 2802 2803 ierr = DMPlexCreateFromDAG(rdm, 0, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr); 2804 } 2805 break; 2806 case 1: 2807 { 2808 PetscInt numPoints[2] = {2, 1}; 2809 PetscInt coneSize[3] = {2, 0, 0}; 2810 PetscInt cones[2] = {1, 2}; 2811 PetscInt coneOrientations[2] = {0, 0}; 2812 PetscScalar vertexCoords[2] = {-1.0, 1.0}; 2813 2814 ierr = DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr); 2815 } 2816 break; 2817 case 2: 2818 if (simplex) { 2819 PetscInt numPoints[2] = {3, 1}; 2820 PetscInt coneSize[4] = {3, 0, 0, 0}; 2821 PetscInt cones[3] = {1, 2, 3}; 2822 PetscInt coneOrientations[3] = {0, 0, 0}; 2823 PetscScalar vertexCoords[6] = {-1.0, -1.0, 1.0, -1.0, -1.0, 1.0}; 2824 2825 ierr = DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr); 2826 } else { 2827 PetscInt numPoints[2] = {4, 1}; 2828 PetscInt coneSize[5] = {4, 0, 0, 0, 0}; 2829 PetscInt cones[4] = {1, 2, 3, 4}; 2830 PetscInt coneOrientations[4] = {0, 0, 0, 0}; 2831 PetscScalar vertexCoords[8] = {-1.0, -1.0, 1.0, -1.0, 1.0, 1.0, -1.0, 1.0}; 2832 2833 ierr = DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr); 2834 } 2835 break; 2836 case 3: 2837 if (simplex) { 2838 PetscInt numPoints[2] = {4, 1}; 2839 PetscInt coneSize[5] = {4, 0, 0, 0, 0}; 2840 PetscInt cones[4] = {1, 3, 2, 4}; 2841 PetscInt coneOrientations[4] = {0, 0, 0, 0}; 2842 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}; 2843 2844 ierr = DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr); 2845 } else { 2846 PetscInt numPoints[2] = {8, 1}; 2847 PetscInt coneSize[9] = {8, 0, 0, 0, 0, 0, 0, 0, 0}; 2848 PetscInt cones[8] = {1, 4, 3, 2, 5, 6, 7, 8}; 2849 PetscInt coneOrientations[8] = {0, 0, 0, 0, 0, 0, 0, 0}; 2850 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, 2851 -1.0, -1.0, 1.0, 1.0, -1.0, 1.0, 1.0, 1.0, 1.0, -1.0, 1.0, 1.0}; 2852 2853 ierr = DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr); 2854 } 2855 break; 2856 default: 2857 SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Cannot create reference cell for dimension %d", dim); 2858 } 2859 *refdm = NULL; 2860 ierr = DMPlexInterpolate(rdm, refdm);CHKERRQ(ierr); 2861 if (rdm->coordinateDM) { 2862 DM ncdm; 2863 PetscSection cs; 2864 PetscInt pEnd = -1; 2865 2866 ierr = DMGetDefaultSection(rdm->coordinateDM, &cs);CHKERRQ(ierr); 2867 if (cs) {ierr = PetscSectionGetChart(cs, NULL, &pEnd);CHKERRQ(ierr);} 2868 if (pEnd >= 0) { 2869 ierr = DMClone(rdm->coordinateDM, &ncdm);CHKERRQ(ierr); 2870 ierr = DMSetDefaultSection(ncdm, cs);CHKERRQ(ierr); 2871 ierr = DMSetCoordinateDM(*refdm, ncdm);CHKERRQ(ierr); 2872 ierr = DMDestroy(&ncdm);CHKERRQ(ierr); 2873 } 2874 } 2875 ierr = DMGetCoordinatesLocal(rdm, &coords);CHKERRQ(ierr); 2876 if (coords) { 2877 ierr = DMSetCoordinatesLocal(*refdm, coords);CHKERRQ(ierr); 2878 } else { 2879 ierr = DMGetCoordinates(rdm, &coords);CHKERRQ(ierr); 2880 if (coords) {ierr = DMSetCoordinates(*refdm, coords);CHKERRQ(ierr);} 2881 } 2882 ierr = DMDestroy(&rdm);CHKERRQ(ierr); 2883 PetscFunctionReturn(0); 2884 } 2885