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] = {lower ? lower[0] : 0.0, lower ? lower[1] : 0.0, lower ? lower[2] : 0.0}; 967 PetscReal upp[3] = {upper ? upper[0] : 1.0, upper ? upper[1] : 1.0, upper ? upper[2] : 1.0}; 968 PetscInt fac[3] = {faces ? faces[0] : 4-dim, faces ? faces[1] : 4-dim, faces ? faces[2] : (dim > 2 ? 4-dim : 0)}; 969 DMBoundaryType bdt[3] = {periodicity ? periodicity[0] : DM_BOUNDARY_NONE, periodicity ? periodicity[1] : DM_BOUNDARY_NONE, periodicity ? periodicity[2] : DM_BOUNDARY_NONE}; 970 PetscErrorCode ierr; 971 972 PetscFunctionBegin; 973 if (simplex) {ierr = DMPlexCreateBoxMesh_Simplex_Internal(comm, dim, fac, low, upp, bdt, interpolate, dm);CHKERRQ(ierr);} 974 else {ierr = DMPlexCreateBoxMesh_Tensor_Internal(comm, dim, fac, low, upp, bdt, interpolate, dm);CHKERRQ(ierr);} 975 PetscFunctionReturn(0); 976 } 977 978 /*@C 979 DMPlexSetOptionsPrefix - Sets the prefix used for searching for all DM options in the database. 980 981 Logically Collective on DM 982 983 Input Parameters: 984 + dm - the DM context 985 - prefix - the prefix to prepend to all option names 986 987 Notes: 988 A hyphen (-) must NOT be given at the beginning of the prefix name. 989 The first character of all runtime options is AUTOMATICALLY the hyphen. 990 991 Level: advanced 992 993 .seealso: SNESSetFromOptions() 994 @*/ 995 PetscErrorCode DMPlexSetOptionsPrefix(DM dm, const char prefix[]) 996 { 997 DM_Plex *mesh = (DM_Plex *) dm->data; 998 PetscErrorCode ierr; 999 1000 PetscFunctionBegin; 1001 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1002 ierr = PetscObjectSetOptionsPrefix((PetscObject) dm, prefix);CHKERRQ(ierr); 1003 ierr = PetscObjectSetOptionsPrefix((PetscObject) mesh->partitioner, prefix);CHKERRQ(ierr); 1004 PetscFunctionReturn(0); 1005 } 1006 1007 /*@ 1008 DMPlexCreateHexCylinderMesh - Creates a mesh on the tensor product of the unit interval with the circle (cylinder) using hexahedra. 1009 1010 Collective on MPI_Comm 1011 1012 Input Parameters: 1013 + comm - The communicator for the DM object 1014 . numRefine - The number of regular refinements to the basic 5 cell structure 1015 - periodicZ - The boundary type for the Z direction 1016 1017 Output Parameter: 1018 . dm - The DM object 1019 1020 Note: Here is the output numbering looking from the bottom of the cylinder: 1021 $ 17-----14 1022 $ | | 1023 $ | 2 | 1024 $ | | 1025 $ 17-----8-----7-----14 1026 $ | | | | 1027 $ | 3 | 0 | 1 | 1028 $ | | | | 1029 $ 19-----5-----6-----13 1030 $ | | 1031 $ | 4 | 1032 $ | | 1033 $ 19-----13 1034 $ 1035 $ and up through the top 1036 $ 1037 $ 18-----16 1038 $ | | 1039 $ | 2 | 1040 $ | | 1041 $ 18----10----11-----16 1042 $ | | | | 1043 $ | 3 | 0 | 1 | 1044 $ | | | | 1045 $ 20-----9----12-----15 1046 $ | | 1047 $ | 4 | 1048 $ | | 1049 $ 20-----15 1050 1051 Level: beginner 1052 1053 .keywords: DM, create 1054 .seealso: DMPlexCreateBoxMesh(), DMSetType(), DMCreate() 1055 @*/ 1056 PetscErrorCode DMPlexCreateHexCylinderMesh(MPI_Comm comm, PetscInt numRefine, DMBoundaryType periodicZ, DM *dm) 1057 { 1058 const PetscInt dim = 3; 1059 PetscInt numCells, numVertices, r; 1060 PetscErrorCode ierr; 1061 1062 PetscFunctionBegin; 1063 PetscValidPointer(dm, 4); 1064 if (numRefine < 0) SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Number of refinements %D cannot be negative", numRefine); 1065 ierr = DMCreate(comm, dm);CHKERRQ(ierr); 1066 ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 1067 ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr); 1068 /* Create topology */ 1069 { 1070 PetscInt cone[8], c; 1071 1072 numCells = 5; 1073 numVertices = 16; 1074 if (periodicZ == DM_BOUNDARY_PERIODIC) { 1075 numCells *= 3; 1076 numVertices = 24; 1077 } 1078 ierr = DMPlexSetChart(*dm, 0, numCells+numVertices);CHKERRQ(ierr); 1079 for (c = 0; c < numCells; c++) {ierr = DMPlexSetConeSize(*dm, c, 8);CHKERRQ(ierr);} 1080 ierr = DMSetUp(*dm);CHKERRQ(ierr); 1081 if (periodicZ == DM_BOUNDARY_PERIODIC) { 1082 cone[0] = 15; cone[1] = 18; cone[2] = 17; cone[3] = 16; 1083 cone[4] = 31; cone[5] = 32; cone[6] = 33; cone[7] = 34; 1084 ierr = DMPlexSetCone(*dm, 0, cone);CHKERRQ(ierr); 1085 cone[0] = 16; cone[1] = 17; cone[2] = 24; cone[3] = 23; 1086 cone[4] = 32; cone[5] = 36; cone[6] = 37; cone[7] = 33; /* 22 25 26 21 */ 1087 ierr = DMPlexSetCone(*dm, 1, cone);CHKERRQ(ierr); 1088 cone[0] = 18; cone[1] = 27; cone[2] = 24; cone[3] = 17; 1089 cone[4] = 34; cone[5] = 33; cone[6] = 37; cone[7] = 38; 1090 ierr = DMPlexSetCone(*dm, 2, cone);CHKERRQ(ierr); 1091 cone[0] = 29; cone[1] = 27; cone[2] = 18; cone[3] = 15; 1092 cone[4] = 35; cone[5] = 31; cone[6] = 34; cone[7] = 38; 1093 ierr = DMPlexSetCone(*dm, 3, cone);CHKERRQ(ierr); 1094 cone[0] = 29; cone[1] = 15; cone[2] = 16; cone[3] = 23; 1095 cone[4] = 35; cone[5] = 36; cone[6] = 32; cone[7] = 31; 1096 ierr = DMPlexSetCone(*dm, 4, cone);CHKERRQ(ierr); 1097 1098 cone[0] = 31; cone[1] = 34; cone[2] = 33; cone[3] = 32; 1099 cone[4] = 19; cone[5] = 22; cone[6] = 21; cone[7] = 20; 1100 ierr = DMPlexSetCone(*dm, 5, cone);CHKERRQ(ierr); 1101 cone[0] = 32; cone[1] = 33; cone[2] = 37; cone[3] = 36; 1102 cone[4] = 22; cone[5] = 25; cone[6] = 26; cone[7] = 21; 1103 ierr = DMPlexSetCone(*dm, 6, cone);CHKERRQ(ierr); 1104 cone[0] = 34; cone[1] = 38; cone[2] = 37; cone[3] = 33; 1105 cone[4] = 20; cone[5] = 21; cone[6] = 26; cone[7] = 28; 1106 ierr = DMPlexSetCone(*dm, 7, cone);CHKERRQ(ierr); 1107 cone[0] = 35; cone[1] = 38; cone[2] = 34; cone[3] = 31; 1108 cone[4] = 30; cone[5] = 19; cone[6] = 20; cone[7] = 28; 1109 ierr = DMPlexSetCone(*dm, 8, cone);CHKERRQ(ierr); 1110 cone[0] = 35; cone[1] = 31; cone[2] = 32; cone[3] = 36; 1111 cone[4] = 30; cone[5] = 25; cone[6] = 22; cone[7] = 19; 1112 ierr = DMPlexSetCone(*dm, 9, cone);CHKERRQ(ierr); 1113 1114 cone[0] = 19; cone[1] = 20; cone[2] = 21; cone[3] = 22; 1115 cone[4] = 15; cone[5] = 16; cone[6] = 17; cone[7] = 18; 1116 ierr = DMPlexSetCone(*dm, 10, cone);CHKERRQ(ierr); 1117 cone[0] = 22; cone[1] = 21; cone[2] = 26; cone[3] = 25; 1118 cone[4] = 16; cone[5] = 23; cone[6] = 24; cone[7] = 17; 1119 ierr = DMPlexSetCone(*dm, 11, cone);CHKERRQ(ierr); 1120 cone[0] = 20; cone[1] = 28; cone[2] = 26; cone[3] = 21; 1121 cone[4] = 18; cone[5] = 17; cone[6] = 24; cone[7] = 27; 1122 ierr = DMPlexSetCone(*dm, 12, cone);CHKERRQ(ierr); 1123 cone[0] = 30; cone[1] = 28; cone[2] = 20; cone[3] = 19; 1124 cone[4] = 29; cone[5] = 15; cone[6] = 18; cone[7] = 27; 1125 ierr = DMPlexSetCone(*dm, 13, cone);CHKERRQ(ierr); 1126 cone[0] = 30; cone[1] = 19; cone[2] = 22; cone[3] = 25; 1127 cone[4] = 29; cone[5] = 23; cone[6] = 16; cone[7] = 15; 1128 ierr = DMPlexSetCone(*dm, 14, cone);CHKERRQ(ierr); 1129 } else { 1130 cone[0] = 5; cone[1] = 8; cone[2] = 7; cone[3] = 6; 1131 cone[4] = 9; cone[5] = 12; cone[6] = 11; cone[7] = 10; 1132 ierr = DMPlexSetCone(*dm, 0, cone);CHKERRQ(ierr); 1133 cone[0] = 6; cone[1] = 7; cone[2] = 14; cone[3] = 13; 1134 cone[4] = 12; cone[5] = 15; cone[6] = 16; cone[7] = 11; 1135 ierr = DMPlexSetCone(*dm, 1, cone);CHKERRQ(ierr); 1136 cone[0] = 8; cone[1] = 17; cone[2] = 14; cone[3] = 7; 1137 cone[4] = 10; cone[5] = 11; cone[6] = 16; cone[7] = 18; 1138 ierr = DMPlexSetCone(*dm, 2, cone);CHKERRQ(ierr); 1139 cone[0] = 19; cone[1] = 17; cone[2] = 8; cone[3] = 5; 1140 cone[4] = 20; cone[5] = 9; cone[6] = 10; cone[7] = 18; 1141 ierr = DMPlexSetCone(*dm, 3, cone);CHKERRQ(ierr); 1142 cone[0] = 19; cone[1] = 5; cone[2] = 6; cone[3] = 13; 1143 cone[4] = 20; cone[5] = 15; cone[6] = 12; cone[7] = 9; 1144 ierr = DMPlexSetCone(*dm, 4, cone);CHKERRQ(ierr); 1145 } 1146 ierr = DMPlexSymmetrize(*dm);CHKERRQ(ierr); 1147 ierr = DMPlexStratify(*dm);CHKERRQ(ierr); 1148 } 1149 /* Interpolate */ 1150 { 1151 DM idm = NULL; 1152 1153 ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr); 1154 ierr = DMDestroy(dm);CHKERRQ(ierr); 1155 *dm = idm; 1156 } 1157 /* Create cube geometry */ 1158 { 1159 Vec coordinates; 1160 PetscSection coordSection; 1161 PetscScalar *coords; 1162 PetscInt coordSize, v; 1163 const PetscReal dis = 1.0/PetscSqrtReal(2.0); 1164 const PetscReal ds2 = dis/2.0; 1165 1166 /* Build coordinates */ 1167 ierr = DMGetCoordinateSection(*dm, &coordSection);CHKERRQ(ierr); 1168 ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr); 1169 ierr = PetscSectionSetFieldComponents(coordSection, 0, dim);CHKERRQ(ierr); 1170 ierr = PetscSectionSetChart(coordSection, numCells, numCells+numVertices);CHKERRQ(ierr); 1171 for (v = numCells; v < numCells+numVertices; ++v) { 1172 ierr = PetscSectionSetDof(coordSection, v, dim);CHKERRQ(ierr); 1173 ierr = PetscSectionSetFieldDof(coordSection, v, 0, dim);CHKERRQ(ierr); 1174 } 1175 ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 1176 ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 1177 ierr = VecCreate(PETSC_COMM_SELF, &coordinates);CHKERRQ(ierr); 1178 ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr); 1179 ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 1180 ierr = VecSetBlockSize(coordinates, dim);CHKERRQ(ierr); 1181 ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr); 1182 ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 1183 coords[0*dim+0] = -ds2; coords[0*dim+1] = -ds2; coords[0*dim+2] = 0.0; 1184 coords[1*dim+0] = ds2; coords[1*dim+1] = -ds2; coords[1*dim+2] = 0.0; 1185 coords[2*dim+0] = ds2; coords[2*dim+1] = ds2; coords[2*dim+2] = 0.0; 1186 coords[3*dim+0] = -ds2; coords[3*dim+1] = ds2; coords[3*dim+2] = 0.0; 1187 coords[4*dim+0] = -ds2; coords[4*dim+1] = -ds2; coords[4*dim+2] = 1.0; 1188 coords[5*dim+0] = -ds2; coords[5*dim+1] = ds2; coords[5*dim+2] = 1.0; 1189 coords[6*dim+0] = ds2; coords[6*dim+1] = ds2; coords[6*dim+2] = 1.0; 1190 coords[7*dim+0] = ds2; coords[7*dim+1] = -ds2; coords[7*dim+2] = 1.0; 1191 coords[ 8*dim+0] = dis; coords[ 8*dim+1] = -dis; coords[ 8*dim+2] = 0.0; 1192 coords[ 9*dim+0] = dis; coords[ 9*dim+1] = dis; coords[ 9*dim+2] = 0.0; 1193 coords[10*dim+0] = dis; coords[10*dim+1] = -dis; coords[10*dim+2] = 1.0; 1194 coords[11*dim+0] = dis; coords[11*dim+1] = dis; coords[11*dim+2] = 1.0; 1195 coords[12*dim+0] = -dis; coords[12*dim+1] = dis; coords[12*dim+2] = 0.0; 1196 coords[13*dim+0] = -dis; coords[13*dim+1] = dis; coords[13*dim+2] = 1.0; 1197 coords[14*dim+0] = -dis; coords[14*dim+1] = -dis; coords[14*dim+2] = 0.0; 1198 coords[15*dim+0] = -dis; coords[15*dim+1] = -dis; coords[15*dim+2] = 1.0; 1199 if (periodicZ == DM_BOUNDARY_PERIODIC) { 1200 /* 15 31 19 */ coords[16*dim+0] = -ds2; coords[16*dim+1] = -ds2; coords[16*dim+2] = 0.5; 1201 /* 16 32 22 */ coords[17*dim+0] = ds2; coords[17*dim+1] = -ds2; coords[17*dim+2] = 0.5; 1202 /* 17 33 21 */ coords[18*dim+0] = ds2; coords[18*dim+1] = ds2; coords[18*dim+2] = 0.5; 1203 /* 18 34 20 */ coords[19*dim+0] = -ds2; coords[19*dim+1] = ds2; coords[19*dim+2] = 0.5; 1204 /* 29 35 30 */ coords[20*dim+0] = -dis; coords[20*dim+1] = -dis; coords[20*dim+2] = 0.5; 1205 /* 23 36 25 */ coords[21*dim+0] = dis; coords[21*dim+1] = -dis; coords[21*dim+2] = 0.5; 1206 /* 24 37 26 */ coords[22*dim+0] = dis; coords[22*dim+1] = dis; coords[22*dim+2] = 0.5; 1207 /* 27 38 28 */ coords[23*dim+0] = -dis; coords[23*dim+1] = dis; coords[23*dim+2] = 0.5; 1208 } 1209 ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 1210 ierr = DMSetCoordinatesLocal(*dm, coordinates);CHKERRQ(ierr); 1211 ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 1212 } 1213 /* Create periodicity */ 1214 if (periodicZ == DM_BOUNDARY_PERIODIC || periodicZ == DM_BOUNDARY_TWIST) { 1215 PetscReal L[3]; 1216 PetscReal maxCell[3]; 1217 DMBoundaryType bdType[3]; 1218 PetscReal lower[3] = {0.0, 0.0, 0.0}; 1219 PetscReal upper[3] = {1.0, 1.0, 1.5}; 1220 PetscInt i, numZCells = 3; 1221 1222 bdType[0] = DM_BOUNDARY_NONE; 1223 bdType[1] = DM_BOUNDARY_NONE; 1224 bdType[2] = periodicZ; 1225 for (i = 0; i < dim; i++) { 1226 L[i] = upper[i] - lower[i]; 1227 maxCell[i] = 1.1 * (L[i] / numZCells); 1228 } 1229 ierr = DMSetPeriodicity(*dm, PETSC_TRUE, maxCell, L, bdType);CHKERRQ(ierr); 1230 } 1231 /* Refine topology */ 1232 for (r = 0; r < numRefine; ++r) { 1233 DM rdm = NULL; 1234 1235 ierr = DMRefine(*dm, comm, &rdm);CHKERRQ(ierr); 1236 ierr = DMDestroy(dm);CHKERRQ(ierr); 1237 *dm = rdm; 1238 } 1239 /* Remap geometry to cylinder 1240 Interior square: Linear interpolation is correct 1241 The other cells all have vertices on rays from the origin. We want to uniformly expand the spacing 1242 such that the last vertex is on the unit circle. So the closest and farthest vertices are at distance 1243 1244 phi = arctan(y/x) 1245 d_close = sqrt(1/8 + 1/4 sin^2(phi)) 1246 d_far = sqrt(1/2 + sin^2(phi)) 1247 1248 so we remap them using 1249 1250 x_new = x_close + (x - x_close) (1 - d_close) / (d_far - d_close) 1251 y_new = y_close + (y - y_close) (1 - d_close) / (d_far - d_close) 1252 1253 If pi/4 < phi < 3pi/4 or -3pi/4 < phi < -pi/4, then we switch x and y. 1254 */ 1255 { 1256 Vec coordinates; 1257 PetscSection coordSection; 1258 PetscScalar *coords; 1259 PetscInt vStart, vEnd, v; 1260 const PetscReal dis = 1.0/PetscSqrtReal(2.0); 1261 const PetscReal ds2 = 0.5*dis; 1262 1263 ierr = DMPlexGetDepthStratum(*dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 1264 ierr = DMGetCoordinateSection(*dm, &coordSection);CHKERRQ(ierr); 1265 ierr = DMGetCoordinatesLocal(*dm, &coordinates);CHKERRQ(ierr); 1266 ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 1267 for (v = vStart; v < vEnd; ++v) { 1268 PetscReal phi, sinp, cosp, dc, df, x, y, xc, yc; 1269 PetscInt off; 1270 1271 ierr = PetscSectionGetOffset(coordSection, v, &off);CHKERRQ(ierr); 1272 if ((PetscAbsScalar(coords[off+0]) <= ds2) && (PetscAbsScalar(coords[off+1]) <= ds2)) continue; 1273 x = PetscRealPart(coords[off]); 1274 y = PetscRealPart(coords[off+1]); 1275 phi = PetscAtan2Real(y, x); 1276 sinp = PetscSinReal(phi); 1277 cosp = PetscCosReal(phi); 1278 if ((PetscAbsReal(phi) > PETSC_PI/4.0) && (PetscAbsReal(phi) < 3.0*PETSC_PI/4.0)) { 1279 dc = PetscAbsReal(ds2/sinp); 1280 df = PetscAbsReal(dis/sinp); 1281 xc = ds2*x/PetscAbsReal(y); 1282 yc = ds2*PetscSignReal(y); 1283 } else { 1284 dc = PetscAbsReal(ds2/cosp); 1285 df = PetscAbsReal(dis/cosp); 1286 xc = ds2*PetscSignReal(x); 1287 yc = ds2*y/PetscAbsReal(x); 1288 } 1289 coords[off+0] = xc + (coords[off+0] - xc)*(1.0 - dc)/(df - dc); 1290 coords[off+1] = yc + (coords[off+1] - yc)*(1.0 - dc)/(df - dc); 1291 } 1292 ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 1293 if (periodicZ == DM_BOUNDARY_PERIODIC || periodicZ == DM_BOUNDARY_TWIST) { 1294 ierr = DMLocalizeCoordinates(*dm);CHKERRQ(ierr); 1295 } 1296 } 1297 PetscFunctionReturn(0); 1298 } 1299 1300 /*@ 1301 DMPlexCreateWedgeCylinderMesh - Creates a mesh on the tensor product of the unit interval with the circle (cylinder) using wedges. 1302 1303 Collective on MPI_Comm 1304 1305 Input Parameters: 1306 + comm - The communicator for the DM object 1307 . n - The number of wedges around the origin 1308 - interpolate - Create edges and faces 1309 1310 Output Parameter: 1311 . dm - The DM object 1312 1313 Level: beginner 1314 1315 .keywords: DM, create 1316 .seealso: DMPlexCreateHexCylinderMesh(), DMPlexCreateBoxMesh(), DMSetType(), DMCreate() 1317 @*/ 1318 PetscErrorCode DMPlexCreateWedgeCylinderMesh(MPI_Comm comm, PetscInt n, PetscBool interpolate, DM *dm) 1319 { 1320 const PetscInt dim = 3; 1321 PetscInt numCells, numVertices; 1322 PetscErrorCode ierr; 1323 1324 PetscFunctionBegin; 1325 PetscValidPointer(dm, 3); 1326 if (n < 0) SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Number of wedges %D cannot be negative", n); 1327 ierr = DMCreate(comm, dm);CHKERRQ(ierr); 1328 ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 1329 ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr); 1330 /* Create topology */ 1331 { 1332 PetscInt cone[6], c; 1333 1334 numCells = n; 1335 numVertices = 2*(n+1); 1336 ierr = DMPlexSetChart(*dm, 0, numCells+numVertices);CHKERRQ(ierr); 1337 for (c = 0; c < numCells; c++) {ierr = DMPlexSetConeSize(*dm, c, 6);CHKERRQ(ierr);} 1338 ierr = DMSetUp(*dm);CHKERRQ(ierr); 1339 for (c = 0; c < numCells; c++) { 1340 cone[0] = c+n*1; cone[1] = (c+1)%n+n*1; cone[2] = 0+3*n; 1341 cone[3] = c+n*2; cone[4] = (c+1)%n+n*2; cone[5] = 1+3*n; 1342 ierr = DMPlexSetCone(*dm, c, cone);CHKERRQ(ierr); 1343 } 1344 ierr = DMPlexSymmetrize(*dm);CHKERRQ(ierr); 1345 ierr = DMPlexStratify(*dm);CHKERRQ(ierr); 1346 } 1347 /* Interpolate */ 1348 if (interpolate) { 1349 DM idm = NULL; 1350 1351 ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr); 1352 ierr = DMDestroy(dm);CHKERRQ(ierr); 1353 *dm = idm; 1354 } 1355 /* Create cylinder geometry */ 1356 { 1357 Vec coordinates; 1358 PetscSection coordSection; 1359 PetscScalar *coords; 1360 PetscInt coordSize, v, c; 1361 1362 /* Build coordinates */ 1363 ierr = DMGetCoordinateSection(*dm, &coordSection);CHKERRQ(ierr); 1364 ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr); 1365 ierr = PetscSectionSetFieldComponents(coordSection, 0, dim);CHKERRQ(ierr); 1366 ierr = PetscSectionSetChart(coordSection, numCells, numCells+numVertices);CHKERRQ(ierr); 1367 for (v = numCells; v < numCells+numVertices; ++v) { 1368 ierr = PetscSectionSetDof(coordSection, v, dim);CHKERRQ(ierr); 1369 ierr = PetscSectionSetFieldDof(coordSection, v, 0, dim);CHKERRQ(ierr); 1370 } 1371 ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 1372 ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 1373 ierr = VecCreate(PETSC_COMM_SELF, &coordinates);CHKERRQ(ierr); 1374 ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr); 1375 ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 1376 ierr = VecSetBlockSize(coordinates, dim);CHKERRQ(ierr); 1377 ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr); 1378 ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 1379 for (c = 0; c < numCells; c++) { 1380 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; 1381 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; 1382 } 1383 coords[(2*n+0)*dim+0] = 0.0; coords[(2*n+0)*dim+1] = 0.0; coords[(2*n+0)*dim+2] = 1.0; 1384 coords[(2*n+1)*dim+0] = 0.0; coords[(2*n+1)*dim+1] = 0.0; coords[(2*n+1)*dim+2] = 0.0; 1385 ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 1386 ierr = DMSetCoordinatesLocal(*dm, coordinates);CHKERRQ(ierr); 1387 ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 1388 } 1389 PetscFunctionReturn(0); 1390 } 1391 1392 PETSC_STATIC_INLINE PetscReal DiffNormReal(PetscInt dim, const PetscReal x[], const PetscReal y[]) 1393 { 1394 PetscReal prod = 0.0; 1395 PetscInt i; 1396 for (i = 0; i < dim; ++i) prod += PetscSqr(x[i] - y[i]); 1397 return PetscSqrtReal(prod); 1398 } 1399 PETSC_STATIC_INLINE PetscReal DotReal(PetscInt dim, const PetscReal x[], const PetscReal y[]) 1400 { 1401 PetscReal prod = 0.0; 1402 PetscInt i; 1403 for (i = 0; i < dim; ++i) prod += x[i]*y[i]; 1404 return prod; 1405 } 1406 1407 /*@ 1408 DMPlexCreateSphereMesh - Creates a mesh on the d-dimensional sphere, S^d. 1409 1410 Collective on MPI_Comm 1411 1412 Input Parameters: 1413 . comm - The communicator for the DM object 1414 . dim - The dimension 1415 - simplex - Use simplices, or tensor product cells 1416 1417 Output Parameter: 1418 . dm - The DM object 1419 1420 Level: beginner 1421 1422 .keywords: DM, create 1423 .seealso: DMPlexCreateBoxMesh(), DMSetType(), DMCreate() 1424 @*/ 1425 PetscErrorCode DMPlexCreateSphereMesh(MPI_Comm comm, PetscInt dim, PetscBool simplex, DM *dm) 1426 { 1427 const PetscInt embedDim = dim+1; 1428 PetscSection coordSection; 1429 Vec coordinates; 1430 PetscScalar *coords; 1431 PetscReal *coordsIn; 1432 PetscInt numCells, numEdges, numVerts, firstVertex, v, firstEdge, coordSize, d, c, e; 1433 PetscMPIInt rank; 1434 PetscErrorCode ierr; 1435 1436 PetscFunctionBegin; 1437 PetscValidPointer(dm, 4); 1438 ierr = DMCreate(comm, dm);CHKERRQ(ierr); 1439 ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 1440 ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr); 1441 ierr = DMSetCoordinateDim(*dm, dim+1);CHKERRQ(ierr); 1442 ierr = MPI_Comm_rank(PetscObjectComm((PetscObject) *dm), &rank);CHKERRQ(ierr); 1443 switch (dim) { 1444 case 2: 1445 if (simplex) { 1446 DM idm = NULL; 1447 const PetscReal edgeLen = 2.0/(1.0 + PETSC_PHI); 1448 const PetscReal vertex[3] = {0.0, 1.0/(1.0 + PETSC_PHI), PETSC_PHI/(1.0 + PETSC_PHI)}; 1449 const PetscInt degree = 5; 1450 PetscInt s[3] = {1, 1, 1}; 1451 PetscInt cone[3]; 1452 PetscInt *graph, p, i, j, k; 1453 1454 numCells = !rank ? 20 : 0; 1455 numVerts = !rank ? 12 : 0; 1456 firstVertex = numCells; 1457 /* Use icosahedron, which for a unit sphere has coordinates which are all cyclic permutations of 1458 1459 (0, \pm 1/\phi+1, \pm \phi/\phi+1) 1460 1461 where \phi^2 - \phi - 1 = 0, meaning \phi is the golden ratio \frac{1 + \sqrt{5}}{2}. The edge 1462 length is then given by 2/\phi = 2 * 2.73606 = 5.47214. 1463 */ 1464 /* Construct vertices */ 1465 ierr = PetscCalloc1(numVerts * embedDim, &coordsIn);CHKERRQ(ierr); 1466 for (p = 0, i = 0; p < embedDim; ++p) { 1467 for (s[1] = -1; s[1] < 2; s[1] += 2) { 1468 for (s[2] = -1; s[2] < 2; s[2] += 2) { 1469 for (d = 0; d < embedDim; ++d) coordsIn[i*embedDim+d] = s[(d+p)%embedDim]*vertex[(d+p)%embedDim]; 1470 ++i; 1471 } 1472 } 1473 } 1474 /* Construct graph */ 1475 ierr = PetscCalloc1(numVerts * numVerts, &graph);CHKERRQ(ierr); 1476 for (i = 0; i < numVerts; ++i) { 1477 for (j = 0, k = 0; j < numVerts; ++j) { 1478 if (PetscAbsReal(DiffNormReal(embedDim, &coordsIn[i*embedDim], &coordsIn[j*embedDim]) - edgeLen) < PETSC_SMALL) {graph[i*numVerts+j] = 1; ++k;} 1479 } 1480 if (k != degree) SETERRQ3(comm, PETSC_ERR_PLIB, "Invalid icosahedron, vertex %D degree %D != %D", i, k, degree); 1481 } 1482 /* Build Topology */ 1483 ierr = DMPlexSetChart(*dm, 0, numCells+numVerts);CHKERRQ(ierr); 1484 for (c = 0; c < numCells; c++) { 1485 ierr = DMPlexSetConeSize(*dm, c, embedDim);CHKERRQ(ierr); 1486 } 1487 ierr = DMSetUp(*dm);CHKERRQ(ierr); /* Allocate space for cones */ 1488 /* Cells */ 1489 for (i = 0, c = 0; i < numVerts; ++i) { 1490 for (j = 0; j < i; ++j) { 1491 for (k = 0; k < j; ++k) { 1492 if (graph[i*numVerts+j] && graph[j*numVerts+k] && graph[k*numVerts+i]) { 1493 cone[0] = firstVertex+i; cone[1] = firstVertex+j; cone[2] = firstVertex+k; 1494 /* Check orientation */ 1495 { 1496 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}}}; 1497 PetscReal normal[3]; 1498 PetscInt e, f; 1499 1500 for (d = 0; d < embedDim; ++d) { 1501 normal[d] = 0.0; 1502 for (e = 0; e < embedDim; ++e) { 1503 for (f = 0; f < embedDim; ++f) { 1504 normal[d] += epsilon[d][e][f]*(coordsIn[j*embedDim+e] - coordsIn[i*embedDim+e])*(coordsIn[k*embedDim+f] - coordsIn[i*embedDim+f]); 1505 } 1506 } 1507 } 1508 if (DotReal(embedDim, normal, &coordsIn[i*embedDim]) < 0) {PetscInt tmp = cone[1]; cone[1] = cone[2]; cone[2] = tmp;} 1509 } 1510 ierr = DMPlexSetCone(*dm, c++, cone);CHKERRQ(ierr); 1511 } 1512 } 1513 } 1514 } 1515 ierr = DMPlexSymmetrize(*dm);CHKERRQ(ierr); 1516 ierr = DMPlexStratify(*dm);CHKERRQ(ierr); 1517 ierr = PetscFree(graph);CHKERRQ(ierr); 1518 /* Interpolate mesh */ 1519 ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr); 1520 ierr = DMDestroy(dm);CHKERRQ(ierr); 1521 *dm = idm; 1522 } else { 1523 /* 1524 12-21--13 1525 | | 1526 25 4 24 1527 | | 1528 12-25--9-16--8-24--13 1529 | | | | 1530 23 5 17 0 15 3 22 1531 | | | | 1532 10-20--6-14--7-19--11 1533 | | 1534 20 1 19 1535 | | 1536 10-18--11 1537 | | 1538 23 2 22 1539 | | 1540 12-21--13 1541 */ 1542 const PetscReal dist = 1.0/PetscSqrtReal(3.0); 1543 PetscInt cone[4], ornt[4]; 1544 1545 numCells = !rank ? 6 : 0; 1546 numEdges = !rank ? 12 : 0; 1547 numVerts = !rank ? 8 : 0; 1548 firstVertex = numCells; 1549 firstEdge = numCells + numVerts; 1550 /* Build Topology */ 1551 ierr = DMPlexSetChart(*dm, 0, numCells+numEdges+numVerts);CHKERRQ(ierr); 1552 for (c = 0; c < numCells; c++) { 1553 ierr = DMPlexSetConeSize(*dm, c, 4);CHKERRQ(ierr); 1554 } 1555 for (e = firstEdge; e < firstEdge+numEdges; ++e) { 1556 ierr = DMPlexSetConeSize(*dm, e, 2);CHKERRQ(ierr); 1557 } 1558 ierr = DMSetUp(*dm);CHKERRQ(ierr); /* Allocate space for cones */ 1559 /* Cell 0 */ 1560 cone[0] = 14; cone[1] = 15; cone[2] = 16; cone[3] = 17; 1561 ierr = DMPlexSetCone(*dm, 0, cone);CHKERRQ(ierr); 1562 ornt[0] = 0; ornt[1] = 0; ornt[2] = 0; ornt[3] = 0; 1563 ierr = DMPlexSetConeOrientation(*dm, 0, ornt);CHKERRQ(ierr); 1564 /* Cell 1 */ 1565 cone[0] = 18; cone[1] = 19; cone[2] = 14; cone[3] = 20; 1566 ierr = DMPlexSetCone(*dm, 1, cone);CHKERRQ(ierr); 1567 ornt[0] = 0; ornt[1] = 0; ornt[2] = -2; ornt[3] = 0; 1568 ierr = DMPlexSetConeOrientation(*dm, 1, ornt);CHKERRQ(ierr); 1569 /* Cell 2 */ 1570 cone[0] = 21; cone[1] = 22; cone[2] = 18; cone[3] = 23; 1571 ierr = DMPlexSetCone(*dm, 2, cone);CHKERRQ(ierr); 1572 ornt[0] = 0; ornt[1] = 0; ornt[2] = -2; ornt[3] = 0; 1573 ierr = DMPlexSetConeOrientation(*dm, 2, ornt);CHKERRQ(ierr); 1574 /* Cell 3 */ 1575 cone[0] = 19; cone[1] = 22; cone[2] = 24; cone[3] = 15; 1576 ierr = DMPlexSetCone(*dm, 3, cone);CHKERRQ(ierr); 1577 ornt[0] = -2; ornt[1] = -2; ornt[2] = 0; ornt[3] = -2; 1578 ierr = DMPlexSetConeOrientation(*dm, 3, ornt);CHKERRQ(ierr); 1579 /* Cell 4 */ 1580 cone[0] = 16; cone[1] = 24; cone[2] = 21; cone[3] = 25; 1581 ierr = DMPlexSetCone(*dm, 4, cone);CHKERRQ(ierr); 1582 ornt[0] = -2; ornt[1] = -2; ornt[2] = -2; ornt[3] = 0; 1583 ierr = DMPlexSetConeOrientation(*dm, 4, ornt);CHKERRQ(ierr); 1584 /* Cell 5 */ 1585 cone[0] = 20; cone[1] = 17; cone[2] = 25; cone[3] = 23; 1586 ierr = DMPlexSetCone(*dm, 5, cone);CHKERRQ(ierr); 1587 ornt[0] = -2; ornt[1] = -2; ornt[2] = -2; ornt[3] = -2; 1588 ierr = DMPlexSetConeOrientation(*dm, 5, ornt);CHKERRQ(ierr); 1589 /* Edges */ 1590 cone[0] = 6; cone[1] = 7; 1591 ierr = DMPlexSetCone(*dm, 14, cone);CHKERRQ(ierr); 1592 cone[0] = 7; cone[1] = 8; 1593 ierr = DMPlexSetCone(*dm, 15, cone);CHKERRQ(ierr); 1594 cone[0] = 8; cone[1] = 9; 1595 ierr = DMPlexSetCone(*dm, 16, cone);CHKERRQ(ierr); 1596 cone[0] = 9; cone[1] = 6; 1597 ierr = DMPlexSetCone(*dm, 17, cone);CHKERRQ(ierr); 1598 cone[0] = 10; cone[1] = 11; 1599 ierr = DMPlexSetCone(*dm, 18, cone);CHKERRQ(ierr); 1600 cone[0] = 11; cone[1] = 7; 1601 ierr = DMPlexSetCone(*dm, 19, cone);CHKERRQ(ierr); 1602 cone[0] = 6; cone[1] = 10; 1603 ierr = DMPlexSetCone(*dm, 20, cone);CHKERRQ(ierr); 1604 cone[0] = 12; cone[1] = 13; 1605 ierr = DMPlexSetCone(*dm, 21, cone);CHKERRQ(ierr); 1606 cone[0] = 13; cone[1] = 11; 1607 ierr = DMPlexSetCone(*dm, 22, cone);CHKERRQ(ierr); 1608 cone[0] = 10; cone[1] = 12; 1609 ierr = DMPlexSetCone(*dm, 23, cone);CHKERRQ(ierr); 1610 cone[0] = 13; cone[1] = 8; 1611 ierr = DMPlexSetCone(*dm, 24, cone);CHKERRQ(ierr); 1612 cone[0] = 12; cone[1] = 9; 1613 ierr = DMPlexSetCone(*dm, 25, cone);CHKERRQ(ierr); 1614 ierr = DMPlexSymmetrize(*dm);CHKERRQ(ierr); 1615 ierr = DMPlexStratify(*dm);CHKERRQ(ierr); 1616 /* Build coordinates */ 1617 ierr = PetscCalloc1(numVerts * embedDim, &coordsIn);CHKERRQ(ierr); 1618 coordsIn[0*embedDim+0] = -dist; coordsIn[0*embedDim+1] = dist; coordsIn[0*embedDim+2] = -dist; 1619 coordsIn[1*embedDim+0] = dist; coordsIn[1*embedDim+1] = dist; coordsIn[1*embedDim+2] = -dist; 1620 coordsIn[2*embedDim+0] = dist; coordsIn[2*embedDim+1] = -dist; coordsIn[2*embedDim+2] = -dist; 1621 coordsIn[3*embedDim+0] = -dist; coordsIn[3*embedDim+1] = -dist; coordsIn[3*embedDim+2] = -dist; 1622 coordsIn[4*embedDim+0] = -dist; coordsIn[4*embedDim+1] = dist; coordsIn[4*embedDim+2] = dist; 1623 coordsIn[5*embedDim+0] = dist; coordsIn[5*embedDim+1] = dist; coordsIn[5*embedDim+2] = dist; 1624 coordsIn[6*embedDim+0] = -dist; coordsIn[6*embedDim+1] = -dist; coordsIn[6*embedDim+2] = dist; 1625 coordsIn[7*embedDim+0] = dist; coordsIn[7*embedDim+1] = -dist; coordsIn[7*embedDim+2] = dist; 1626 } 1627 break; 1628 case 3: 1629 if (simplex) { 1630 DM idm = NULL; 1631 const PetscReal edgeLen = 1.0/PETSC_PHI; 1632 const PetscReal vertexA[4] = {0.5, 0.5, 0.5, 0.5}; 1633 const PetscReal vertexB[4] = {1.0, 0.0, 0.0, 0.0}; 1634 const PetscReal vertexC[4] = {0.5, 0.5*PETSC_PHI, 0.5/PETSC_PHI, 0.0}; 1635 const PetscInt degree = 12; 1636 PetscInt s[4] = {1, 1, 1}; 1637 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}, 1638 {2, 0, 1, 3}, {2, 1, 3, 0}, {2, 3, 0, 1}, {3, 0, 2, 1}, {3, 1, 0, 2}, {3, 2, 1, 0}}; 1639 PetscInt cone[4]; 1640 PetscInt *graph, p, i, j, k, l; 1641 1642 numCells = !rank ? 600 : 0; 1643 numVerts = !rank ? 120 : 0; 1644 firstVertex = numCells; 1645 /* Use the 600-cell, which for a unit sphere has coordinates which are 1646 1647 1/2 (\pm 1, \pm 1, \pm 1, \pm 1) 16 1648 (\pm 1, 0, 0, 0) all cyclic permutations 8 1649 1/2 (\pm 1, \pm phi, \pm 1/phi, 0) all even permutations 96 1650 1651 where \phi^2 - \phi - 1 = 0, meaning \phi is the golden ratio \frac{1 + \sqrt{5}}{2}. The edge 1652 length is then given by 1/\phi = 2.73606. 1653 1654 http://buzzard.pugetsound.edu/sage-practice/ch03s03.html 1655 http://mathworld.wolfram.com/600-Cell.html 1656 */ 1657 /* Construct vertices */ 1658 ierr = PetscCalloc1(numVerts * embedDim, &coordsIn);CHKERRQ(ierr); 1659 i = 0; 1660 for (s[0] = -1; s[0] < 2; s[0] += 2) { 1661 for (s[1] = -1; s[1] < 2; s[1] += 2) { 1662 for (s[2] = -1; s[2] < 2; s[2] += 2) { 1663 for (s[3] = -1; s[3] < 2; s[3] += 2) { 1664 for (d = 0; d < embedDim; ++d) coordsIn[i*embedDim+d] = s[d]*vertexA[d]; 1665 ++i; 1666 } 1667 } 1668 } 1669 } 1670 for (p = 0; p < embedDim; ++p) { 1671 s[1] = s[2] = s[3] = 1; 1672 for (s[0] = -1; s[0] < 2; s[0] += 2) { 1673 for (d = 0; d < embedDim; ++d) coordsIn[i*embedDim+d] = s[(d+p)%embedDim]*vertexB[(d+p)%embedDim]; 1674 ++i; 1675 } 1676 } 1677 for (p = 0; p < 12; ++p) { 1678 s[3] = 1; 1679 for (s[0] = -1; s[0] < 2; s[0] += 2) { 1680 for (s[1] = -1; s[1] < 2; s[1] += 2) { 1681 for (s[2] = -1; s[2] < 2; s[2] += 2) { 1682 for (d = 0; d < embedDim; ++d) coordsIn[i*embedDim+d] = s[evenPerm[p][d]]*vertexC[evenPerm[p][d]]; 1683 ++i; 1684 } 1685 } 1686 } 1687 } 1688 if (i != numVerts) SETERRQ2(comm, PETSC_ERR_PLIB, "Invalid 600-cell, vertices %D != %D", i, numVerts); 1689 /* Construct graph */ 1690 ierr = PetscCalloc1(numVerts * numVerts, &graph);CHKERRQ(ierr); 1691 for (i = 0; i < numVerts; ++i) { 1692 for (j = 0, k = 0; j < numVerts; ++j) { 1693 if (PetscAbsReal(DiffNormReal(embedDim, &coordsIn[i*embedDim], &coordsIn[j*embedDim]) - edgeLen) < PETSC_SMALL) {graph[i*numVerts+j] = 1; ++k;} 1694 } 1695 if (k != degree) SETERRQ3(comm, PETSC_ERR_PLIB, "Invalid 600-cell, vertex %D degree %D != %D", i, k, degree); 1696 } 1697 /* Build Topology */ 1698 ierr = DMPlexSetChart(*dm, 0, numCells+numVerts);CHKERRQ(ierr); 1699 for (c = 0; c < numCells; c++) { 1700 ierr = DMPlexSetConeSize(*dm, c, embedDim);CHKERRQ(ierr); 1701 } 1702 ierr = DMSetUp(*dm);CHKERRQ(ierr); /* Allocate space for cones */ 1703 /* Cells */ 1704 for (i = 0, c = 0; i < numVerts; ++i) { 1705 for (j = 0; j < i; ++j) { 1706 for (k = 0; k < j; ++k) { 1707 for (l = 0; l < k; ++l) { 1708 if (graph[i*numVerts+j] && graph[j*numVerts+k] && graph[k*numVerts+i] && 1709 graph[l*numVerts+i] && graph[l*numVerts+j] && graph[l*numVerts+k]) { 1710 cone[0] = firstVertex+i; cone[1] = firstVertex+j; cone[2] = firstVertex+k; cone[3] = firstVertex+l; 1711 /* Check orientation: https://ef.gy/linear-algebra:normal-vectors-in-higher-dimensional-spaces */ 1712 { 1713 const PetscInt epsilon[4][4][4][4] = {{{{0, 0, 0, 0}, { 0, 0, 0, 0}, { 0, 0, 0, 0}, { 0, 0, 0, 0}}, 1714 {{0, 0, 0, 0}, { 0, 0, 0, 0}, { 0, 0, 0, 1}, { 0, 0, -1, 0}}, 1715 {{0, 0, 0, 0}, { 0, 0, 0, -1}, { 0, 0, 0, 0}, { 0, 1, 0, 0}}, 1716 {{0, 0, 0, 0}, { 0, 0, 1, 0}, { 0, -1, 0, 0}, { 0, 0, 0, 0}}}, 1717 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, 0}, { 0, 0, 0, 0}, { 0, 0, 0, 0}}, 1720 {{0, 0, 0, 1}, { 0, 0, 0, 0}, { 0, 0, 0, 0}, {-1, 0, 0, 0}}, 1721 {{0, 0, -1, 0}, { 0, 0, 0, 0}, { 1, 0, 0, 0}, { 0, 0, 0, 0}}}, 1722 1723 {{{0, 0, 0, 0}, { 0, 0, 0, 1}, { 0, 0, 0, 0}, { 0, -1, 0, 0}}, 1724 {{0, 0, 0, -1}, { 0, 0, 0, 0}, { 0, 0, 0, 0}, { 1, 0, 0, 0}}, 1725 {{0, 0, 0, 0}, { 0, 0, 0, 0}, { 0, 0, 0, 0}, { 0, 0, 0, 0}}, 1726 {{0, 1, 0, 0}, {-1, 0, 0, 0}, { 0, 0, 0, 0}, { 0, 0, 0, 0}}}, 1727 1728 {{{0, 0, 0, 0}, { 0, 0, -1, 0}, { 0, 1, 0, 0}, { 0, 0, 0, 0}}, 1729 {{0, 0, 1, 0}, { 0, 0, 0, 0}, {-1, 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 {{0, 0, 0, 0}, { 0, 0, 0, 0}, { 0, 0, 0, 0}, { 0, 0, 0, 0}}}}; 1732 PetscReal normal[4]; 1733 PetscInt e, f, g; 1734 1735 for (d = 0; d < embedDim; ++d) { 1736 normal[d] = 0.0; 1737 for (e = 0; e < embedDim; ++e) { 1738 for (f = 0; f < embedDim; ++f) { 1739 for (g = 0; g < embedDim; ++g) { 1740 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]); 1741 } 1742 } 1743 } 1744 } 1745 if (DotReal(embedDim, normal, &coordsIn[i*embedDim]) < 0) {PetscInt tmp = cone[1]; cone[1] = cone[2]; cone[2] = tmp;} 1746 } 1747 ierr = DMPlexSetCone(*dm, c++, cone);CHKERRQ(ierr); 1748 } 1749 } 1750 } 1751 } 1752 } 1753 ierr = DMPlexSymmetrize(*dm);CHKERRQ(ierr); 1754 ierr = DMPlexStratify(*dm);CHKERRQ(ierr); 1755 ierr = PetscFree(graph);CHKERRQ(ierr); 1756 /* Interpolate mesh */ 1757 ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr); 1758 ierr = DMDestroy(dm);CHKERRQ(ierr); 1759 *dm = idm; 1760 break; 1761 } 1762 default: SETERRQ1(comm, PETSC_ERR_SUP, "Unsupported dimension for sphere: %D", dim); 1763 } 1764 /* Create coordinates */ 1765 ierr = DMGetCoordinateSection(*dm, &coordSection);CHKERRQ(ierr); 1766 ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr); 1767 ierr = PetscSectionSetFieldComponents(coordSection, 0, embedDim);CHKERRQ(ierr); 1768 ierr = PetscSectionSetChart(coordSection, firstVertex, firstVertex+numVerts);CHKERRQ(ierr); 1769 for (v = firstVertex; v < firstVertex+numVerts; ++v) { 1770 ierr = PetscSectionSetDof(coordSection, v, embedDim);CHKERRQ(ierr); 1771 ierr = PetscSectionSetFieldDof(coordSection, v, 0, embedDim);CHKERRQ(ierr); 1772 } 1773 ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 1774 ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 1775 ierr = VecCreate(PETSC_COMM_SELF, &coordinates);CHKERRQ(ierr); 1776 ierr = VecSetBlockSize(coordinates, embedDim);CHKERRQ(ierr); 1777 ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr); 1778 ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 1779 ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr); 1780 ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 1781 for (v = 0; v < numVerts; ++v) for (d = 0; d < embedDim; ++d) {coords[v*embedDim+d] = coordsIn[v*embedDim+d];} 1782 ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 1783 ierr = DMSetCoordinatesLocal(*dm, coordinates);CHKERRQ(ierr); 1784 ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 1785 ierr = PetscFree(coordsIn);CHKERRQ(ierr); 1786 PetscFunctionReturn(0); 1787 } 1788 1789 /* External function declarations here */ 1790 extern PetscErrorCode DMCreateInterpolation_Plex(DM dmCoarse, DM dmFine, Mat *interpolation, Vec *scaling); 1791 extern PetscErrorCode DMCreateInjection_Plex(DM dmCoarse, DM dmFine, Mat *mat); 1792 extern PetscErrorCode DMCreateDefaultSection_Plex(DM dm); 1793 extern PetscErrorCode DMCreateDefaultConstraints_Plex(DM dm); 1794 extern PetscErrorCode DMCreateMatrix_Plex(DM dm, Mat *J); 1795 extern PetscErrorCode DMCreateCoordinateDM_Plex(DM dm, DM *cdm); 1796 PETSC_INTERN PetscErrorCode DMClone_Plex(DM dm, DM *newdm); 1797 extern PetscErrorCode DMSetUp_Plex(DM dm); 1798 extern PetscErrorCode DMDestroy_Plex(DM dm); 1799 extern PetscErrorCode DMView_Plex(DM dm, PetscViewer viewer); 1800 extern PetscErrorCode DMLoad_Plex(DM dm, PetscViewer viewer); 1801 extern PetscErrorCode DMCreateSubDM_Plex(DM dm, PetscInt numFields, PetscInt fields[], IS *is, DM *subdm); 1802 static PetscErrorCode DMInitialize_Plex(DM dm); 1803 1804 /* Replace dm with the contents of dmNew 1805 - Share the DM_Plex structure 1806 - Share the coordinates 1807 - Share the SF 1808 */ 1809 static PetscErrorCode DMPlexReplace_Static(DM dm, DM dmNew) 1810 { 1811 PetscSF sf; 1812 DM coordDM, coarseDM; 1813 Vec coords; 1814 PetscBool isper; 1815 const PetscReal *maxCell, *L; 1816 const DMBoundaryType *bd; 1817 PetscErrorCode ierr; 1818 1819 PetscFunctionBegin; 1820 ierr = DMGetPointSF(dmNew, &sf);CHKERRQ(ierr); 1821 ierr = DMSetPointSF(dm, sf);CHKERRQ(ierr); 1822 ierr = DMGetCoordinateDM(dmNew, &coordDM);CHKERRQ(ierr); 1823 ierr = DMGetCoordinatesLocal(dmNew, &coords);CHKERRQ(ierr); 1824 ierr = DMSetCoordinateDM(dm, coordDM);CHKERRQ(ierr); 1825 ierr = DMSetCoordinatesLocal(dm, coords);CHKERRQ(ierr); 1826 ierr = DMGetPeriodicity(dm, &isper, &maxCell, &L, &bd);CHKERRQ(ierr); 1827 ierr = DMSetPeriodicity(dmNew, isper, maxCell, L, bd);CHKERRQ(ierr); 1828 ierr = DMDestroy_Plex(dm);CHKERRQ(ierr); 1829 ierr = DMInitialize_Plex(dm);CHKERRQ(ierr); 1830 dm->data = dmNew->data; 1831 ((DM_Plex *) dmNew->data)->refct++; 1832 dmNew->labels->refct++; 1833 if (!--(dm->labels->refct)) { 1834 DMLabelLink next = dm->labels->next; 1835 1836 /* destroy the labels */ 1837 while (next) { 1838 DMLabelLink tmp = next->next; 1839 1840 ierr = DMLabelDestroy(&next->label);CHKERRQ(ierr); 1841 ierr = PetscFree(next);CHKERRQ(ierr); 1842 next = tmp; 1843 } 1844 ierr = PetscFree(dm->labels);CHKERRQ(ierr); 1845 } 1846 dm->labels = dmNew->labels; 1847 dm->depthLabel = dmNew->depthLabel; 1848 ierr = DMGetCoarseDM(dmNew,&coarseDM);CHKERRQ(ierr); 1849 ierr = DMSetCoarseDM(dm,coarseDM);CHKERRQ(ierr); 1850 PetscFunctionReturn(0); 1851 } 1852 1853 /* Swap dm with the contents of dmNew 1854 - Swap the DM_Plex structure 1855 - Swap the coordinates 1856 - Swap the point PetscSF 1857 */ 1858 static PetscErrorCode DMPlexSwap_Static(DM dmA, DM dmB) 1859 { 1860 DM coordDMA, coordDMB; 1861 Vec coordsA, coordsB; 1862 PetscSF sfA, sfB; 1863 void *tmp; 1864 DMLabelLinkList listTmp; 1865 DMLabel depthTmp; 1866 PetscErrorCode ierr; 1867 1868 PetscFunctionBegin; 1869 ierr = DMGetPointSF(dmA, &sfA);CHKERRQ(ierr); 1870 ierr = DMGetPointSF(dmB, &sfB);CHKERRQ(ierr); 1871 ierr = PetscObjectReference((PetscObject) sfA);CHKERRQ(ierr); 1872 ierr = DMSetPointSF(dmA, sfB);CHKERRQ(ierr); 1873 ierr = DMSetPointSF(dmB, sfA);CHKERRQ(ierr); 1874 ierr = PetscObjectDereference((PetscObject) sfA);CHKERRQ(ierr); 1875 1876 ierr = DMGetCoordinateDM(dmA, &coordDMA);CHKERRQ(ierr); 1877 ierr = DMGetCoordinateDM(dmB, &coordDMB);CHKERRQ(ierr); 1878 ierr = PetscObjectReference((PetscObject) coordDMA);CHKERRQ(ierr); 1879 ierr = DMSetCoordinateDM(dmA, coordDMB);CHKERRQ(ierr); 1880 ierr = DMSetCoordinateDM(dmB, coordDMA);CHKERRQ(ierr); 1881 ierr = PetscObjectDereference((PetscObject) coordDMA);CHKERRQ(ierr); 1882 1883 ierr = DMGetCoordinatesLocal(dmA, &coordsA);CHKERRQ(ierr); 1884 ierr = DMGetCoordinatesLocal(dmB, &coordsB);CHKERRQ(ierr); 1885 ierr = PetscObjectReference((PetscObject) coordsA);CHKERRQ(ierr); 1886 ierr = DMSetCoordinatesLocal(dmA, coordsB);CHKERRQ(ierr); 1887 ierr = DMSetCoordinatesLocal(dmB, coordsA);CHKERRQ(ierr); 1888 ierr = PetscObjectDereference((PetscObject) coordsA);CHKERRQ(ierr); 1889 1890 tmp = dmA->data; 1891 dmA->data = dmB->data; 1892 dmB->data = tmp; 1893 listTmp = dmA->labels; 1894 dmA->labels = dmB->labels; 1895 dmB->labels = listTmp; 1896 depthTmp = dmA->depthLabel; 1897 dmA->depthLabel = dmB->depthLabel; 1898 dmB->depthLabel = depthTmp; 1899 PetscFunctionReturn(0); 1900 } 1901 1902 PetscErrorCode DMSetFromOptions_NonRefinement_Plex(PetscOptionItems *PetscOptionsObject,DM dm) 1903 { 1904 DM_Plex *mesh = (DM_Plex*) dm->data; 1905 PetscErrorCode ierr; 1906 1907 PetscFunctionBegin; 1908 /* Handle viewing */ 1909 ierr = PetscOptionsBool("-dm_plex_print_set_values", "Output all set values info", "DMView", PETSC_FALSE, &mesh->printSetValues, NULL);CHKERRQ(ierr); 1910 ierr = PetscOptionsInt("-dm_plex_print_fem", "Debug output level all fem computations", "DMView", 0, &mesh->printFEM, NULL);CHKERRQ(ierr); 1911 ierr = PetscOptionsReal("-dm_plex_print_tol", "Tolerance for FEM output", "DMView", mesh->printTol, &mesh->printTol, NULL);CHKERRQ(ierr); 1912 /* Point Location */ 1913 ierr = PetscOptionsBool("-dm_plex_hash_location", "Use grid hashing for point location", "DMView", PETSC_FALSE, &mesh->useHashLocation, NULL);CHKERRQ(ierr); 1914 /* Generation and remeshing */ 1915 ierr = PetscOptionsBool("-dm_plex_remesh_bd", "Allow changes to the boundary on remeshing", "DMView", PETSC_FALSE, &mesh->remeshBd, NULL);CHKERRQ(ierr); 1916 /* Projection behavior */ 1917 ierr = PetscOptionsInt("-dm_plex_max_projection_height", "Maxmimum mesh point height used to project locally", "DMPlexSetMaxProjectionHeight", 0, &mesh->maxProjectionHeight, NULL);CHKERRQ(ierr); 1918 ierr = PetscOptionsBool("-dm_plex_regular_refinement", "Use special nested projection algorithm for regular refinement", "DMPlexSetRegularRefinement", mesh->regularRefinement, &mesh->regularRefinement, NULL);CHKERRQ(ierr); 1919 1920 ierr = PetscPartitionerSetFromOptions(mesh->partitioner);CHKERRQ(ierr); 1921 PetscFunctionReturn(0); 1922 } 1923 1924 static PetscErrorCode DMSetFromOptions_Plex(PetscOptionItems *PetscOptionsObject,DM dm) 1925 { 1926 PetscInt refine = 0, coarsen = 0, r; 1927 PetscBool isHierarchy; 1928 PetscErrorCode ierr; 1929 1930 PetscFunctionBegin; 1931 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1932 ierr = PetscOptionsHead(PetscOptionsObject,"DMPlex Options");CHKERRQ(ierr); 1933 /* Handle DMPlex refinement */ 1934 ierr = PetscOptionsInt("-dm_refine", "The number of uniform refinements", "DMCreate", refine, &refine, NULL);CHKERRQ(ierr); 1935 ierr = PetscOptionsInt("-dm_refine_hierarchy", "The number of uniform refinements", "DMCreate", refine, &refine, &isHierarchy);CHKERRQ(ierr); 1936 if (refine) {ierr = DMPlexSetRefinementUniform(dm, PETSC_TRUE);CHKERRQ(ierr);} 1937 if (refine && isHierarchy) { 1938 DM *dms, coarseDM; 1939 1940 ierr = DMGetCoarseDM(dm, &coarseDM);CHKERRQ(ierr); 1941 ierr = PetscObjectReference((PetscObject)coarseDM);CHKERRQ(ierr); 1942 ierr = PetscMalloc1(refine,&dms);CHKERRQ(ierr); 1943 ierr = DMRefineHierarchy(dm, refine, dms);CHKERRQ(ierr); 1944 /* Total hack since we do not pass in a pointer */ 1945 ierr = DMPlexSwap_Static(dm, dms[refine-1]);CHKERRQ(ierr); 1946 if (refine == 1) { 1947 ierr = DMSetCoarseDM(dm, dms[0]);CHKERRQ(ierr); 1948 ierr = DMPlexSetRegularRefinement(dm, PETSC_TRUE);CHKERRQ(ierr); 1949 } else { 1950 ierr = DMSetCoarseDM(dm, dms[refine-2]);CHKERRQ(ierr); 1951 ierr = DMPlexSetRegularRefinement(dm, PETSC_TRUE);CHKERRQ(ierr); 1952 ierr = DMSetCoarseDM(dms[0], dms[refine-1]);CHKERRQ(ierr); 1953 ierr = DMPlexSetRegularRefinement(dms[0], PETSC_TRUE);CHKERRQ(ierr); 1954 } 1955 ierr = DMSetCoarseDM(dms[refine-1], coarseDM);CHKERRQ(ierr); 1956 ierr = PetscObjectDereference((PetscObject)coarseDM);CHKERRQ(ierr); 1957 /* Free DMs */ 1958 for (r = 0; r < refine; ++r) { 1959 ierr = DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dms[r]);CHKERRQ(ierr); 1960 ierr = DMDestroy(&dms[r]);CHKERRQ(ierr); 1961 } 1962 ierr = PetscFree(dms);CHKERRQ(ierr); 1963 } else { 1964 for (r = 0; r < refine; ++r) { 1965 DM refinedMesh; 1966 1967 ierr = DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dm);CHKERRQ(ierr); 1968 ierr = DMRefine(dm, PetscObjectComm((PetscObject) dm), &refinedMesh);CHKERRQ(ierr); 1969 /* Total hack since we do not pass in a pointer */ 1970 ierr = DMPlexReplace_Static(dm, refinedMesh);CHKERRQ(ierr); 1971 ierr = DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dm);CHKERRQ(ierr); 1972 ierr = DMDestroy(&refinedMesh);CHKERRQ(ierr); 1973 } 1974 } 1975 /* Handle DMPlex coarsening */ 1976 ierr = PetscOptionsInt("-dm_coarsen", "Coarsen the mesh", "DMCreate", coarsen, &coarsen, NULL);CHKERRQ(ierr); 1977 ierr = PetscOptionsInt("-dm_coarsen_hierarchy", "The number of coarsenings", "DMCreate", coarsen, &coarsen, &isHierarchy);CHKERRQ(ierr); 1978 if (coarsen && isHierarchy) { 1979 DM *dms; 1980 1981 ierr = PetscMalloc1(coarsen, &dms);CHKERRQ(ierr); 1982 ierr = DMCoarsenHierarchy(dm, coarsen, dms);CHKERRQ(ierr); 1983 /* Free DMs */ 1984 for (r = 0; r < coarsen; ++r) { 1985 ierr = DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dms[r]);CHKERRQ(ierr); 1986 ierr = DMDestroy(&dms[r]);CHKERRQ(ierr); 1987 } 1988 ierr = PetscFree(dms);CHKERRQ(ierr); 1989 } else { 1990 for (r = 0; r < coarsen; ++r) { 1991 DM coarseMesh; 1992 1993 ierr = DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dm);CHKERRQ(ierr); 1994 ierr = DMCoarsen(dm, PetscObjectComm((PetscObject) dm), &coarseMesh);CHKERRQ(ierr); 1995 /* Total hack since we do not pass in a pointer */ 1996 ierr = DMPlexReplace_Static(dm, coarseMesh);CHKERRQ(ierr); 1997 ierr = DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dm);CHKERRQ(ierr); 1998 ierr = DMDestroy(&coarseMesh);CHKERRQ(ierr); 1999 } 2000 } 2001 /* Handle */ 2002 ierr = DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dm);CHKERRQ(ierr); 2003 ierr = PetscOptionsTail();CHKERRQ(ierr); 2004 PetscFunctionReturn(0); 2005 } 2006 2007 static PetscErrorCode DMCreateGlobalVector_Plex(DM dm,Vec *vec) 2008 { 2009 PetscErrorCode ierr; 2010 2011 PetscFunctionBegin; 2012 ierr = DMCreateGlobalVector_Section_Private(dm,vec);CHKERRQ(ierr); 2013 /* ierr = VecSetOperation(*vec, VECOP_DUPLICATE, (void(*)(void)) VecDuplicate_MPI_DM);CHKERRQ(ierr); */ 2014 ierr = VecSetOperation(*vec, VECOP_VIEW, (void (*)(void)) VecView_Plex);CHKERRQ(ierr); 2015 ierr = VecSetOperation(*vec, VECOP_VIEWNATIVE, (void (*)(void)) VecView_Plex_Native);CHKERRQ(ierr); 2016 ierr = VecSetOperation(*vec, VECOP_LOAD, (void (*)(void)) VecLoad_Plex);CHKERRQ(ierr); 2017 ierr = VecSetOperation(*vec, VECOP_LOADNATIVE, (void (*)(void)) VecLoad_Plex_Native);CHKERRQ(ierr); 2018 PetscFunctionReturn(0); 2019 } 2020 2021 static PetscErrorCode DMCreateLocalVector_Plex(DM dm,Vec *vec) 2022 { 2023 PetscErrorCode ierr; 2024 2025 PetscFunctionBegin; 2026 ierr = DMCreateLocalVector_Section_Private(dm,vec);CHKERRQ(ierr); 2027 ierr = VecSetOperation(*vec, VECOP_VIEW, (void (*)(void)) VecView_Plex_Local);CHKERRQ(ierr); 2028 ierr = VecSetOperation(*vec, VECOP_LOAD, (void (*)(void)) VecLoad_Plex_Local);CHKERRQ(ierr); 2029 PetscFunctionReturn(0); 2030 } 2031 2032 static PetscErrorCode DMGetDimPoints_Plex(DM dm, PetscInt dim, PetscInt *pStart, PetscInt *pEnd) 2033 { 2034 PetscInt depth, d; 2035 PetscErrorCode ierr; 2036 2037 PetscFunctionBegin; 2038 ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 2039 if (depth == 1) { 2040 ierr = DMGetDimension(dm, &d);CHKERRQ(ierr); 2041 if (dim == 0) {ierr = DMPlexGetDepthStratum(dm, dim, pStart, pEnd);CHKERRQ(ierr);} 2042 else if (dim == d) {ierr = DMPlexGetDepthStratum(dm, 1, pStart, pEnd);CHKERRQ(ierr);} 2043 else {*pStart = 0; *pEnd = 0;} 2044 } else { 2045 ierr = DMPlexGetDepthStratum(dm, dim, pStart, pEnd);CHKERRQ(ierr); 2046 } 2047 PetscFunctionReturn(0); 2048 } 2049 2050 static PetscErrorCode DMGetNeighors_Plex(DM dm, PetscInt *nranks, const PetscMPIInt *ranks[]) 2051 { 2052 PetscSF sf; 2053 PetscErrorCode ierr; 2054 2055 PetscFunctionBegin; 2056 ierr = DMGetPointSF(dm, &sf);CHKERRQ(ierr); 2057 ierr = PetscSFGetRanks(sf, nranks, ranks, NULL, NULL, NULL);CHKERRQ(ierr); 2058 PetscFunctionReturn(0); 2059 } 2060 2061 static PetscErrorCode DMInitialize_Plex(DM dm) 2062 { 2063 PetscErrorCode ierr; 2064 2065 PetscFunctionBegin; 2066 dm->ops->view = DMView_Plex; 2067 dm->ops->load = DMLoad_Plex; 2068 dm->ops->setfromoptions = DMSetFromOptions_Plex; 2069 dm->ops->clone = DMClone_Plex; 2070 dm->ops->setup = DMSetUp_Plex; 2071 dm->ops->createdefaultsection = DMCreateDefaultSection_Plex; 2072 dm->ops->createdefaultconstraints = DMCreateDefaultConstraints_Plex; 2073 dm->ops->createglobalvector = DMCreateGlobalVector_Plex; 2074 dm->ops->createlocalvector = DMCreateLocalVector_Plex; 2075 dm->ops->getlocaltoglobalmapping = NULL; 2076 dm->ops->createfieldis = NULL; 2077 dm->ops->createcoordinatedm = DMCreateCoordinateDM_Plex; 2078 dm->ops->getcoloring = NULL; 2079 dm->ops->creatematrix = DMCreateMatrix_Plex; 2080 dm->ops->createinterpolation = DMCreateInterpolation_Plex; 2081 dm->ops->getaggregates = NULL; 2082 dm->ops->getinjection = DMCreateInjection_Plex; 2083 dm->ops->refine = DMRefine_Plex; 2084 dm->ops->coarsen = DMCoarsen_Plex; 2085 dm->ops->refinehierarchy = DMRefineHierarchy_Plex; 2086 dm->ops->coarsenhierarchy = DMCoarsenHierarchy_Plex; 2087 dm->ops->adaptlabel = DMAdaptLabel_Plex; 2088 dm->ops->adaptmetric = DMAdaptMetric_Plex; 2089 dm->ops->globaltolocalbegin = NULL; 2090 dm->ops->globaltolocalend = NULL; 2091 dm->ops->localtoglobalbegin = NULL; 2092 dm->ops->localtoglobalend = NULL; 2093 dm->ops->destroy = DMDestroy_Plex; 2094 dm->ops->createsubdm = DMCreateSubDM_Plex; 2095 dm->ops->getdimpoints = DMGetDimPoints_Plex; 2096 dm->ops->locatepoints = DMLocatePoints_Plex; 2097 dm->ops->projectfunctionlocal = DMProjectFunctionLocal_Plex; 2098 dm->ops->projectfunctionlabellocal = DMProjectFunctionLabelLocal_Plex; 2099 dm->ops->projectfieldlocal = DMProjectFieldLocal_Plex; 2100 dm->ops->projectfieldlabellocal = DMProjectFieldLabelLocal_Plex; 2101 dm->ops->computel2diff = DMComputeL2Diff_Plex; 2102 dm->ops->computel2gradientdiff = DMComputeL2GradientDiff_Plex; 2103 dm->ops->computel2fielddiff = DMComputeL2FieldDiff_Plex; 2104 dm->ops->getneighbors = DMGetNeighors_Plex; 2105 ierr = PetscObjectComposeFunction((PetscObject)dm,"DMPlexInsertBoundaryValues_C",DMPlexInsertBoundaryValues_Plex);CHKERRQ(ierr); 2106 ierr = PetscObjectComposeFunction((PetscObject)dm,"DMSetUpGLVisViewer_C",DMSetUpGLVisViewer_Plex);CHKERRQ(ierr); 2107 PetscFunctionReturn(0); 2108 } 2109 2110 PETSC_INTERN PetscErrorCode DMClone_Plex(DM dm, DM *newdm) 2111 { 2112 DM_Plex *mesh = (DM_Plex *) dm->data; 2113 PetscErrorCode ierr; 2114 2115 PetscFunctionBegin; 2116 mesh->refct++; 2117 (*newdm)->data = mesh; 2118 ierr = PetscObjectChangeTypeName((PetscObject) *newdm, DMPLEX);CHKERRQ(ierr); 2119 ierr = DMInitialize_Plex(*newdm);CHKERRQ(ierr); 2120 PetscFunctionReturn(0); 2121 } 2122 2123 /*MC 2124 DMPLEX = "plex" - A DM object that encapsulates an unstructured mesh, or CW Complex, which can be expressed using a Hasse Diagram. 2125 In the local representation, Vecs contain all unknowns in the interior and shared boundary. This is 2126 specified by a PetscSection object. Ownership in the global representation is determined by 2127 ownership of the underlying DMPlex points. This is specified by another PetscSection object. 2128 2129 Level: intermediate 2130 2131 .seealso: DMType, DMPlexCreate(), DMCreate(), DMSetType() 2132 M*/ 2133 2134 PETSC_EXTERN PetscErrorCode DMCreate_Plex(DM dm) 2135 { 2136 DM_Plex *mesh; 2137 PetscInt unit, d; 2138 PetscErrorCode ierr; 2139 2140 PetscFunctionBegin; 2141 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 2142 ierr = PetscNewLog(dm,&mesh);CHKERRQ(ierr); 2143 dm->dim = 0; 2144 dm->data = mesh; 2145 2146 mesh->refct = 1; 2147 ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm), &mesh->coneSection);CHKERRQ(ierr); 2148 mesh->maxConeSize = 0; 2149 mesh->cones = NULL; 2150 mesh->coneOrientations = NULL; 2151 ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm), &mesh->supportSection);CHKERRQ(ierr); 2152 mesh->maxSupportSize = 0; 2153 mesh->supports = NULL; 2154 mesh->refinementUniform = PETSC_TRUE; 2155 mesh->refinementLimit = -1.0; 2156 2157 mesh->facesTmp = NULL; 2158 2159 mesh->tetgenOpts = NULL; 2160 mesh->triangleOpts = NULL; 2161 ierr = PetscPartitionerCreate(PetscObjectComm((PetscObject)dm), &mesh->partitioner);CHKERRQ(ierr); 2162 mesh->remeshBd = PETSC_FALSE; 2163 2164 mesh->subpointMap = NULL; 2165 2166 for (unit = 0; unit < NUM_PETSC_UNITS; ++unit) mesh->scale[unit] = 1.0; 2167 2168 mesh->regularRefinement = PETSC_FALSE; 2169 mesh->depthState = -1; 2170 mesh->globalVertexNumbers = NULL; 2171 mesh->globalCellNumbers = NULL; 2172 mesh->anchorSection = NULL; 2173 mesh->anchorIS = NULL; 2174 mesh->createanchors = NULL; 2175 mesh->computeanchormatrix = NULL; 2176 mesh->parentSection = NULL; 2177 mesh->parents = NULL; 2178 mesh->childIDs = NULL; 2179 mesh->childSection = NULL; 2180 mesh->children = NULL; 2181 mesh->referenceTree = NULL; 2182 mesh->getchildsymmetry = NULL; 2183 for (d = 0; d < 8; ++d) mesh->hybridPointMax[d] = PETSC_DETERMINE; 2184 mesh->vtkCellHeight = 0; 2185 mesh->useCone = PETSC_FALSE; 2186 mesh->useClosure = PETSC_TRUE; 2187 mesh->useAnchors = PETSC_FALSE; 2188 2189 mesh->maxProjectionHeight = 0; 2190 2191 mesh->printSetValues = PETSC_FALSE; 2192 mesh->printFEM = 0; 2193 mesh->printTol = 1.0e-10; 2194 2195 ierr = DMInitialize_Plex(dm);CHKERRQ(ierr); 2196 PetscFunctionReturn(0); 2197 } 2198 2199 /*@ 2200 DMPlexCreate - Creates a DMPlex object, which encapsulates an unstructured mesh, or CW complex, which can be expressed using a Hasse Diagram. 2201 2202 Collective on MPI_Comm 2203 2204 Input Parameter: 2205 . comm - The communicator for the DMPlex object 2206 2207 Output Parameter: 2208 . mesh - The DMPlex object 2209 2210 Level: beginner 2211 2212 .keywords: DMPlex, create 2213 @*/ 2214 PetscErrorCode DMPlexCreate(MPI_Comm comm, DM *mesh) 2215 { 2216 PetscErrorCode ierr; 2217 2218 PetscFunctionBegin; 2219 PetscValidPointer(mesh,2); 2220 ierr = DMCreate(comm, mesh);CHKERRQ(ierr); 2221 ierr = DMSetType(*mesh, DMPLEX);CHKERRQ(ierr); 2222 PetscFunctionReturn(0); 2223 } 2224 2225 /* 2226 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 2227 */ 2228 static PetscErrorCode DMPlexBuildFromCellList_Parallel_Private(DM dm, PetscInt numCells, PetscInt numVertices, PetscInt numCorners, const int cells[], PetscSF *sfVert) 2229 { 2230 PetscSF sfPoint; 2231 PetscLayout vLayout; 2232 PetscHashI vhash; 2233 PetscSFNode *remoteVerticesAdj, *vertexLocal, *vertexOwner, *remoteVertex; 2234 const PetscInt *vrange; 2235 PetscInt numVerticesAdj, off, *verticesAdj, numVerticesGhost = 0, *localVertex, *cone, c, p, v, g; 2236 PETSC_UNUSED PetscHashIIter ret, iter; 2237 PetscMPIInt rank, size; 2238 PetscErrorCode ierr; 2239 2240 PetscFunctionBegin; 2241 ierr = MPI_Comm_rank(PetscObjectComm((PetscObject) dm), &rank);CHKERRQ(ierr); 2242 ierr = MPI_Comm_size(PetscObjectComm((PetscObject) dm), &size);CHKERRQ(ierr); 2243 /* Partition vertices */ 2244 ierr = PetscLayoutCreate(PetscObjectComm((PetscObject) dm), &vLayout);CHKERRQ(ierr); 2245 ierr = PetscLayoutSetLocalSize(vLayout, numVertices);CHKERRQ(ierr); 2246 ierr = PetscLayoutSetBlockSize(vLayout, 1);CHKERRQ(ierr); 2247 ierr = PetscLayoutSetUp(vLayout);CHKERRQ(ierr); 2248 ierr = PetscLayoutGetRanges(vLayout, &vrange);CHKERRQ(ierr); 2249 /* Count vertices and map them to procs */ 2250 PetscHashICreate(vhash); 2251 for (c = 0; c < numCells; ++c) { 2252 for (p = 0; p < numCorners; ++p) { 2253 PetscHashIPut(vhash, cells[c*numCorners+p], ret, iter); 2254 } 2255 } 2256 PetscHashISize(vhash, numVerticesAdj); 2257 ierr = PetscMalloc1(numVerticesAdj, &verticesAdj);CHKERRQ(ierr); 2258 off = 0; ierr = PetscHashIGetKeys(vhash, &off, verticesAdj);CHKERRQ(ierr); 2259 if (off != numVerticesAdj) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid number of local vertices %D should be %D", off, numVerticesAdj); 2260 ierr = PetscSortInt(numVerticesAdj, verticesAdj);CHKERRQ(ierr); 2261 ierr = PetscMalloc1(numVerticesAdj, &remoteVerticesAdj);CHKERRQ(ierr); 2262 for (v = 0; v < numVerticesAdj; ++v) { 2263 const PetscInt gv = verticesAdj[v]; 2264 PetscInt vrank; 2265 2266 ierr = PetscFindInt(gv, size+1, vrange, &vrank);CHKERRQ(ierr); 2267 vrank = vrank < 0 ? -(vrank+2) : vrank; 2268 remoteVerticesAdj[v].index = gv - vrange[vrank]; 2269 remoteVerticesAdj[v].rank = vrank; 2270 } 2271 /* Create cones */ 2272 ierr = DMPlexSetChart(dm, 0, numCells+numVerticesAdj);CHKERRQ(ierr); 2273 for (c = 0; c < numCells; ++c) {ierr = DMPlexSetConeSize(dm, c, numCorners);CHKERRQ(ierr);} 2274 ierr = DMSetUp(dm);CHKERRQ(ierr); 2275 ierr = DMGetWorkArray(dm, numCorners, PETSC_INT, &cone);CHKERRQ(ierr); 2276 for (c = 0; c < numCells; ++c) { 2277 for (p = 0; p < numCorners; ++p) { 2278 const PetscInt gv = cells[c*numCorners+p]; 2279 PetscInt lv; 2280 2281 ierr = PetscFindInt(gv, numVerticesAdj, verticesAdj, &lv);CHKERRQ(ierr); 2282 if (lv < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Could not find global vertex %D in local connectivity", gv); 2283 cone[p] = lv+numCells; 2284 } 2285 ierr = DMPlexSetCone(dm, c, cone);CHKERRQ(ierr); 2286 } 2287 ierr = DMRestoreWorkArray(dm, numCorners, PETSC_INT, &cone);CHKERRQ(ierr); 2288 /* Create SF for vertices */ 2289 ierr = PetscSFCreate(PetscObjectComm((PetscObject)dm), sfVert);CHKERRQ(ierr); 2290 ierr = PetscObjectSetName((PetscObject) *sfVert, "Vertex Ownership SF");CHKERRQ(ierr); 2291 ierr = PetscSFSetFromOptions(*sfVert);CHKERRQ(ierr); 2292 ierr = PetscSFSetGraph(*sfVert, numVertices, numVerticesAdj, NULL, PETSC_OWN_POINTER, remoteVerticesAdj, PETSC_OWN_POINTER);CHKERRQ(ierr); 2293 ierr = PetscFree(verticesAdj);CHKERRQ(ierr); 2294 /* Build pointSF */ 2295 ierr = PetscMalloc2(numVerticesAdj, &vertexLocal, numVertices, &vertexOwner);CHKERRQ(ierr); 2296 for (v = 0; v < numVerticesAdj; ++v) {vertexLocal[v].index = v+numCells; vertexLocal[v].rank = rank;} 2297 for (v = 0; v < numVertices; ++v) {vertexOwner[v].index = -1; vertexOwner[v].rank = -1;} 2298 ierr = PetscSFReduceBegin(*sfVert, MPIU_2INT, vertexLocal, vertexOwner, MPI_MAXLOC);CHKERRQ(ierr); 2299 ierr = PetscSFReduceEnd(*sfVert, MPIU_2INT, vertexLocal, vertexOwner, MPI_MAXLOC);CHKERRQ(ierr); 2300 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); 2301 ierr = PetscSFBcastBegin(*sfVert, MPIU_2INT, vertexOwner, vertexLocal);CHKERRQ(ierr); 2302 ierr = PetscSFBcastEnd(*sfVert, MPIU_2INT, vertexOwner, vertexLocal);CHKERRQ(ierr); 2303 for (v = 0; v < numVerticesAdj; ++v) if (vertexLocal[v].rank != rank) ++numVerticesGhost; 2304 ierr = PetscMalloc1(numVerticesGhost, &localVertex);CHKERRQ(ierr); 2305 ierr = PetscMalloc1(numVerticesGhost, &remoteVertex);CHKERRQ(ierr); 2306 for (v = 0, g = 0; v < numVerticesAdj; ++v) { 2307 if (vertexLocal[v].rank != rank) { 2308 localVertex[g] = v+numCells; 2309 remoteVertex[g].index = vertexLocal[v].index; 2310 remoteVertex[g].rank = vertexLocal[v].rank; 2311 ++g; 2312 } 2313 } 2314 ierr = PetscFree2(vertexLocal, vertexOwner);CHKERRQ(ierr); 2315 if (g != numVerticesGhost) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid number of vertex ghosts %D should be %D", g, numVerticesGhost); 2316 ierr = DMGetPointSF(dm, &sfPoint);CHKERRQ(ierr); 2317 ierr = PetscObjectSetName((PetscObject) sfPoint, "point SF");CHKERRQ(ierr); 2318 ierr = PetscSFSetGraph(sfPoint, numCells+numVerticesAdj, numVerticesGhost, localVertex, PETSC_OWN_POINTER, remoteVertex, PETSC_OWN_POINTER);CHKERRQ(ierr); 2319 ierr = PetscLayoutDestroy(&vLayout);CHKERRQ(ierr); 2320 PetscHashIDestroy(vhash); 2321 /* Fill in the rest of the topology structure */ 2322 ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr); 2323 ierr = DMPlexStratify(dm);CHKERRQ(ierr); 2324 PetscFunctionReturn(0); 2325 } 2326 2327 /* 2328 This takes as input the coordinates for each owned vertex 2329 */ 2330 static PetscErrorCode DMPlexBuildCoordinates_Parallel_Private(DM dm, PetscInt spaceDim, PetscInt numCells, PetscInt numV, PetscSF sfVert, const PetscReal vertexCoords[]) 2331 { 2332 PetscSection coordSection; 2333 Vec coordinates; 2334 PetscScalar *coords; 2335 PetscInt numVertices, numVerticesAdj, coordSize, v; 2336 PetscErrorCode ierr; 2337 2338 PetscFunctionBegin; 2339 ierr = DMSetCoordinateDim(dm, spaceDim);CHKERRQ(ierr); 2340 ierr = PetscSFGetGraph(sfVert, &numVertices, &numVerticesAdj, NULL, NULL);CHKERRQ(ierr); 2341 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 2342 ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr); 2343 ierr = PetscSectionSetFieldComponents(coordSection, 0, spaceDim);CHKERRQ(ierr); 2344 ierr = PetscSectionSetChart(coordSection, numCells, numCells + numVerticesAdj);CHKERRQ(ierr); 2345 for (v = numCells; v < numCells+numVerticesAdj; ++v) { 2346 ierr = PetscSectionSetDof(coordSection, v, spaceDim);CHKERRQ(ierr); 2347 ierr = PetscSectionSetFieldDof(coordSection, v, 0, spaceDim);CHKERRQ(ierr); 2348 } 2349 ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 2350 ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 2351 ierr = VecCreate(PetscObjectComm((PetscObject)dm), &coordinates);CHKERRQ(ierr); 2352 ierr = VecSetBlockSize(coordinates, spaceDim);CHKERRQ(ierr); 2353 ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr); 2354 ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 2355 ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr); 2356 ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 2357 { 2358 MPI_Datatype coordtype; 2359 2360 /* Need a temp buffer for coords if we have complex/single */ 2361 ierr = MPI_Type_contiguous(spaceDim, MPIU_SCALAR, &coordtype);CHKERRQ(ierr); 2362 ierr = MPI_Type_commit(&coordtype);CHKERRQ(ierr); 2363 #if defined(PETSC_USE_COMPLEX) 2364 { 2365 PetscScalar *svertexCoords; 2366 PetscInt i; 2367 ierr = PetscMalloc1(numV*spaceDim,&svertexCoords);CHKERRQ(ierr); 2368 for (i=0; i<numV*spaceDim; i++) svertexCoords[i] = vertexCoords[i]; 2369 ierr = PetscSFBcastBegin(sfVert, coordtype, svertexCoords, coords);CHKERRQ(ierr); 2370 ierr = PetscSFBcastEnd(sfVert, coordtype, svertexCoords, coords);CHKERRQ(ierr); 2371 ierr = PetscFree(svertexCoords);CHKERRQ(ierr); 2372 } 2373 #else 2374 ierr = PetscSFBcastBegin(sfVert, coordtype, vertexCoords, coords);CHKERRQ(ierr); 2375 ierr = PetscSFBcastEnd(sfVert, coordtype, vertexCoords, coords);CHKERRQ(ierr); 2376 #endif 2377 ierr = MPI_Type_free(&coordtype);CHKERRQ(ierr); 2378 } 2379 ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 2380 ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr); 2381 ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 2382 PetscFunctionReturn(0); 2383 } 2384 2385 /*@C 2386 DMPlexCreateFromCellListParallel - This takes as input common mesh generator output, a list of the vertices for each cell, and produces a DM 2387 2388 Input Parameters: 2389 + comm - The communicator 2390 . dim - The topological dimension of the mesh 2391 . numCells - The number of cells owned by this process 2392 . numVertices - The number of vertices owned by this process 2393 . numCorners - The number of vertices for each cell 2394 . interpolate - Flag indicating that intermediate mesh entities (faces, edges) should be created automatically 2395 . cells - An array of numCells*numCorners numbers, the global vertex numbers for each cell 2396 . spaceDim - The spatial dimension used for coordinates 2397 - vertexCoords - An array of numVertices*spaceDim numbers, the coordinates of each vertex 2398 2399 Output Parameter: 2400 + dm - The DM 2401 - vertexSF - Optional, SF describing complete vertex ownership 2402 2403 Note: Two triangles sharing a face 2404 $ 2405 $ 2 2406 $ / | \ 2407 $ / | \ 2408 $ / | \ 2409 $ 0 0 | 1 3 2410 $ \ | / 2411 $ \ | / 2412 $ \ | / 2413 $ 1 2414 would have input 2415 $ numCells = 2, numVertices = 4 2416 $ cells = [0 1 2 1 3 2] 2417 $ 2418 which would result in the DMPlex 2419 $ 2420 $ 4 2421 $ / | \ 2422 $ / | \ 2423 $ / | \ 2424 $ 2 0 | 1 5 2425 $ \ | / 2426 $ \ | / 2427 $ \ | / 2428 $ 3 2429 2430 Level: beginner 2431 2432 .seealso: DMPlexCreateFromCellList(), DMPlexCreateFromDAG(), DMPlexCreate() 2433 @*/ 2434 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) 2435 { 2436 PetscSF sfVert; 2437 PetscErrorCode ierr; 2438 2439 PetscFunctionBegin; 2440 ierr = DMCreate(comm, dm);CHKERRQ(ierr); 2441 ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 2442 PetscValidLogicalCollectiveInt(*dm, dim, 2); 2443 PetscValidLogicalCollectiveInt(*dm, spaceDim, 8); 2444 ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr); 2445 ierr = DMPlexBuildFromCellList_Parallel_Private(*dm, numCells, numVertices, numCorners, cells, &sfVert);CHKERRQ(ierr); 2446 if (interpolate) { 2447 DM idm = NULL; 2448 2449 ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr); 2450 ierr = DMDestroy(dm);CHKERRQ(ierr); 2451 *dm = idm; 2452 } 2453 ierr = DMPlexBuildCoordinates_Parallel_Private(*dm, spaceDim, numCells, numVertices,sfVert, vertexCoords);CHKERRQ(ierr); 2454 if (vertexSF) *vertexSF = sfVert; 2455 else {ierr = PetscSFDestroy(&sfVert);CHKERRQ(ierr);} 2456 PetscFunctionReturn(0); 2457 } 2458 2459 /* 2460 This takes as input the common mesh generator output, a list of the vertices for each cell 2461 */ 2462 static PetscErrorCode DMPlexBuildFromCellList_Private(DM dm, PetscInt numCells, PetscInt numVertices, PetscInt numCorners, const int cells[]) 2463 { 2464 PetscInt *cone, c, p; 2465 PetscErrorCode ierr; 2466 2467 PetscFunctionBegin; 2468 ierr = DMPlexSetChart(dm, 0, numCells+numVertices);CHKERRQ(ierr); 2469 for (c = 0; c < numCells; ++c) { 2470 ierr = DMPlexSetConeSize(dm, c, numCorners);CHKERRQ(ierr); 2471 } 2472 ierr = DMSetUp(dm);CHKERRQ(ierr); 2473 ierr = DMGetWorkArray(dm, numCorners, PETSC_INT, &cone);CHKERRQ(ierr); 2474 for (c = 0; c < numCells; ++c) { 2475 for (p = 0; p < numCorners; ++p) { 2476 cone[p] = cells[c*numCorners+p]+numCells; 2477 } 2478 ierr = DMPlexSetCone(dm, c, cone);CHKERRQ(ierr); 2479 } 2480 ierr = DMRestoreWorkArray(dm, numCorners, PETSC_INT, &cone);CHKERRQ(ierr); 2481 ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr); 2482 ierr = DMPlexStratify(dm);CHKERRQ(ierr); 2483 PetscFunctionReturn(0); 2484 } 2485 2486 /* 2487 This takes as input the coordinates for each vertex 2488 */ 2489 static PetscErrorCode DMPlexBuildCoordinates_Private(DM dm, PetscInt spaceDim, PetscInt numCells, PetscInt numVertices, const double vertexCoords[]) 2490 { 2491 PetscSection coordSection; 2492 Vec coordinates; 2493 DM cdm; 2494 PetscScalar *coords; 2495 PetscInt v, d; 2496 PetscErrorCode ierr; 2497 2498 PetscFunctionBegin; 2499 ierr = DMSetCoordinateDim(dm, spaceDim);CHKERRQ(ierr); 2500 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 2501 ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr); 2502 ierr = PetscSectionSetFieldComponents(coordSection, 0, spaceDim);CHKERRQ(ierr); 2503 ierr = PetscSectionSetChart(coordSection, numCells, numCells + numVertices);CHKERRQ(ierr); 2504 for (v = numCells; v < numCells+numVertices; ++v) { 2505 ierr = PetscSectionSetDof(coordSection, v, spaceDim);CHKERRQ(ierr); 2506 ierr = PetscSectionSetFieldDof(coordSection, v, 0, spaceDim);CHKERRQ(ierr); 2507 } 2508 ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 2509 2510 ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr); 2511 ierr = DMCreateLocalVector(cdm, &coordinates);CHKERRQ(ierr); 2512 ierr = VecSetBlockSize(coordinates, spaceDim);CHKERRQ(ierr); 2513 ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr); 2514 ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 2515 for (v = 0; v < numVertices; ++v) { 2516 for (d = 0; d < spaceDim; ++d) { 2517 coords[v*spaceDim+d] = vertexCoords[v*spaceDim+d]; 2518 } 2519 } 2520 ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 2521 ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr); 2522 ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 2523 PetscFunctionReturn(0); 2524 } 2525 2526 /*@C 2527 DMPlexCreateFromCellList - This takes as input common mesh generator output, a list of the vertices for each cell, and produces a DM 2528 2529 Input Parameters: 2530 + comm - The communicator 2531 . dim - The topological dimension of the mesh 2532 . numCells - The number of cells 2533 . numVertices - The number of vertices 2534 . numCorners - The number of vertices for each cell 2535 . interpolate - Flag indicating that intermediate mesh entities (faces, edges) should be created automatically 2536 . cells - An array of numCells*numCorners numbers, the vertices for each cell 2537 . spaceDim - The spatial dimension used for coordinates 2538 - vertexCoords - An array of numVertices*spaceDim numbers, the coordinates of each vertex 2539 2540 Output Parameter: 2541 . dm - The DM 2542 2543 Note: Two triangles sharing a face 2544 $ 2545 $ 2 2546 $ / | \ 2547 $ / | \ 2548 $ / | \ 2549 $ 0 0 | 1 3 2550 $ \ | / 2551 $ \ | / 2552 $ \ | / 2553 $ 1 2554 would have input 2555 $ numCells = 2, numVertices = 4 2556 $ cells = [0 1 2 1 3 2] 2557 $ 2558 which would result in the DMPlex 2559 $ 2560 $ 4 2561 $ / | \ 2562 $ / | \ 2563 $ / | \ 2564 $ 2 0 | 1 5 2565 $ \ | / 2566 $ \ | / 2567 $ \ | / 2568 $ 3 2569 2570 Level: beginner 2571 2572 .seealso: DMPlexCreateFromDAG(), DMPlexCreate() 2573 @*/ 2574 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) 2575 { 2576 PetscErrorCode ierr; 2577 2578 PetscFunctionBegin; 2579 ierr = DMCreate(comm, dm);CHKERRQ(ierr); 2580 ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 2581 ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr); 2582 ierr = DMPlexBuildFromCellList_Private(*dm, numCells, numVertices, numCorners, cells);CHKERRQ(ierr); 2583 if (interpolate) { 2584 DM idm = NULL; 2585 2586 ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr); 2587 ierr = DMDestroy(dm);CHKERRQ(ierr); 2588 *dm = idm; 2589 } 2590 ierr = DMPlexBuildCoordinates_Private(*dm, spaceDim, numCells, numVertices, vertexCoords);CHKERRQ(ierr); 2591 PetscFunctionReturn(0); 2592 } 2593 2594 /*@ 2595 DMPlexCreateFromDAG - This takes as input the adjacency-list representation of the Directed Acyclic Graph (Hasse Diagram) encoding a mesh, and produces a DM 2596 2597 Input Parameters: 2598 + dm - The empty DM object, usually from DMCreate() and DMSetDimension() 2599 . depth - The depth of the DAG 2600 . numPoints - The number of points at each depth 2601 . coneSize - The cone size of each point 2602 . cones - The concatenation of the cone points for each point, the cone list must be oriented correctly for each point 2603 . coneOrientations - The orientation of each cone point 2604 - vertexCoords - An array of numVertices*dim numbers, the coordinates of each vertex 2605 2606 Output Parameter: 2607 . dm - The DM 2608 2609 Note: Two triangles sharing a face would have input 2610 $ depth = 1, numPoints = [4 2], coneSize = [3 3 0 0 0 0] 2611 $ cones = [2 3 4 3 5 4], coneOrientations = [0 0 0 0 0 0] 2612 $ vertexCoords = [-1.0 0.0 0.0 -1.0 0.0 1.0 1.0 0.0] 2613 $ 2614 which would result in the DMPlex 2615 $ 2616 $ 4 2617 $ / | \ 2618 $ / | \ 2619 $ / | \ 2620 $ 2 0 | 1 5 2621 $ \ | / 2622 $ \ | / 2623 $ \ | / 2624 $ 3 2625 $ 2626 $ Notice that all points are numbered consecutively, unlikely DMPlexCreateFromCellList() 2627 2628 Level: advanced 2629 2630 .seealso: DMPlexCreateFromCellList(), DMPlexCreate() 2631 @*/ 2632 PetscErrorCode DMPlexCreateFromDAG(DM dm, PetscInt depth, const PetscInt numPoints[], const PetscInt coneSize[], const PetscInt cones[], const PetscInt coneOrientations[], const PetscScalar vertexCoords[]) 2633 { 2634 Vec coordinates; 2635 PetscSection coordSection; 2636 PetscScalar *coords; 2637 PetscInt coordSize, firstVertex = -1, pStart = 0, pEnd = 0, p, v, dim, dimEmbed, d, off; 2638 PetscErrorCode ierr; 2639 2640 PetscFunctionBegin; 2641 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 2642 ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr); 2643 if (dimEmbed < dim) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Embedding dimension %d cannot be less than intrinsic dimension %d",dimEmbed,dim); 2644 for (d = 0; d <= depth; ++d) pEnd += numPoints[d]; 2645 ierr = DMPlexSetChart(dm, pStart, pEnd);CHKERRQ(ierr); 2646 for (p = pStart; p < pEnd; ++p) { 2647 ierr = DMPlexSetConeSize(dm, p, coneSize[p-pStart]);CHKERRQ(ierr); 2648 if (firstVertex < 0 && !coneSize[p - pStart]) { 2649 firstVertex = p - pStart; 2650 } 2651 } 2652 if (firstVertex < 0 && numPoints[0]) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Expected %d vertices but could not find any", numPoints[0]); 2653 ierr = DMSetUp(dm);CHKERRQ(ierr); /* Allocate space for cones */ 2654 for (p = pStart, off = 0; p < pEnd; off += coneSize[p-pStart], ++p) { 2655 ierr = DMPlexSetCone(dm, p, &cones[off]);CHKERRQ(ierr); 2656 ierr = DMPlexSetConeOrientation(dm, p, &coneOrientations[off]);CHKERRQ(ierr); 2657 } 2658 ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr); 2659 ierr = DMPlexStratify(dm);CHKERRQ(ierr); 2660 /* Build coordinates */ 2661 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 2662 ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr); 2663 ierr = PetscSectionSetFieldComponents(coordSection, 0, dimEmbed);CHKERRQ(ierr); 2664 ierr = PetscSectionSetChart(coordSection, firstVertex, firstVertex+numPoints[0]);CHKERRQ(ierr); 2665 for (v = firstVertex; v < firstVertex+numPoints[0]; ++v) { 2666 ierr = PetscSectionSetDof(coordSection, v, dimEmbed);CHKERRQ(ierr); 2667 ierr = PetscSectionSetFieldDof(coordSection, v, 0, dimEmbed);CHKERRQ(ierr); 2668 } 2669 ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 2670 ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 2671 ierr = VecCreate(PETSC_COMM_SELF, &coordinates);CHKERRQ(ierr); 2672 ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr); 2673 ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 2674 ierr = VecSetBlockSize(coordinates, dimEmbed);CHKERRQ(ierr); 2675 ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr); 2676 ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 2677 for (v = 0; v < numPoints[0]; ++v) { 2678 PetscInt off; 2679 2680 ierr = PetscSectionGetOffset(coordSection, v+firstVertex, &off);CHKERRQ(ierr); 2681 for (d = 0; d < dimEmbed; ++d) { 2682 coords[off+d] = vertexCoords[v*dimEmbed+d]; 2683 } 2684 } 2685 ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 2686 ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr); 2687 ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 2688 PetscFunctionReturn(0); 2689 } 2690 2691 /*@C 2692 DMPlexCreateFromFile - This takes a filename and produces a DM 2693 2694 Input Parameters: 2695 + comm - The communicator 2696 . filename - A file name 2697 - interpolate - Flag to create intermediate mesh pieces (edges, faces) 2698 2699 Output Parameter: 2700 . dm - The DM 2701 2702 Level: beginner 2703 2704 .seealso: DMPlexCreateFromDAG(), DMPlexCreateFromCellList(), DMPlexCreate() 2705 @*/ 2706 PetscErrorCode DMPlexCreateFromFile(MPI_Comm comm, const char filename[], PetscBool interpolate, DM *dm) 2707 { 2708 const char *extGmsh = ".msh"; 2709 const char *extCGNS = ".cgns"; 2710 const char *extExodus = ".exo"; 2711 const char *extGenesis = ".gen"; 2712 const char *extFluent = ".cas"; 2713 const char *extHDF5 = ".h5"; 2714 const char *extMed = ".med"; 2715 const char *extPLY = ".ply"; 2716 size_t len; 2717 PetscBool isGmsh, isCGNS, isExodus, isGenesis, isFluent, isHDF5, isMed, isPLY; 2718 PetscMPIInt rank; 2719 PetscErrorCode ierr; 2720 2721 PetscFunctionBegin; 2722 PetscValidPointer(filename, 2); 2723 PetscValidPointer(dm, 4); 2724 ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 2725 ierr = PetscStrlen(filename, &len);CHKERRQ(ierr); 2726 if (!len) SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Filename must be a valid path"); 2727 ierr = PetscStrncmp(&filename[PetscMax(0,len-4)], extGmsh, 4, &isGmsh);CHKERRQ(ierr); 2728 ierr = PetscStrncmp(&filename[PetscMax(0,len-5)], extCGNS, 5, &isCGNS);CHKERRQ(ierr); 2729 ierr = PetscStrncmp(&filename[PetscMax(0,len-4)], extExodus, 4, &isExodus);CHKERRQ(ierr); 2730 ierr = PetscStrncmp(&filename[PetscMax(0,len-4)], extGenesis, 4, &isGenesis);CHKERRQ(ierr); 2731 ierr = PetscStrncmp(&filename[PetscMax(0,len-4)], extFluent, 4, &isFluent);CHKERRQ(ierr); 2732 ierr = PetscStrncmp(&filename[PetscMax(0,len-3)], extHDF5, 3, &isHDF5);CHKERRQ(ierr); 2733 ierr = PetscStrncmp(&filename[PetscMax(0,len-4)], extMed, 4, &isMed);CHKERRQ(ierr); 2734 ierr = PetscStrncmp(&filename[PetscMax(0,len-4)], extPLY, 4, &isPLY);CHKERRQ(ierr); 2735 if (isGmsh) { 2736 ierr = DMPlexCreateGmshFromFile(comm, filename, interpolate, dm);CHKERRQ(ierr); 2737 } else if (isCGNS) { 2738 ierr = DMPlexCreateCGNSFromFile(comm, filename, interpolate, dm);CHKERRQ(ierr); 2739 } else if (isExodus || isGenesis) { 2740 ierr = DMPlexCreateExodusFromFile(comm, filename, interpolate, dm);CHKERRQ(ierr); 2741 } else if (isFluent) { 2742 ierr = DMPlexCreateFluentFromFile(comm, filename, interpolate, dm);CHKERRQ(ierr); 2743 } else if (isHDF5) { 2744 PetscViewer viewer; 2745 2746 ierr = PetscViewerCreate(comm, &viewer);CHKERRQ(ierr); 2747 ierr = PetscViewerSetType(viewer, PETSCVIEWERHDF5);CHKERRQ(ierr); 2748 ierr = PetscViewerFileSetMode(viewer, FILE_MODE_READ);CHKERRQ(ierr); 2749 ierr = PetscViewerFileSetName(viewer, filename);CHKERRQ(ierr); 2750 ierr = DMCreate(comm, dm);CHKERRQ(ierr); 2751 ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 2752 ierr = DMLoad(*dm, viewer);CHKERRQ(ierr); 2753 ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 2754 } else if (isMed) { 2755 ierr = DMPlexCreateMedFromFile(comm, filename, interpolate, dm);CHKERRQ(ierr); 2756 } else if (isPLY) { 2757 ierr = DMPlexCreatePLYFromFile(comm, filename, interpolate, dm);CHKERRQ(ierr); 2758 } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cannot load file %s: unrecognized extension", filename); 2759 PetscFunctionReturn(0); 2760 } 2761 2762 /*@ 2763 DMPlexCreateReferenceCell - Create a DMPLEX with the appropriate FEM reference cell 2764 2765 Collective on comm 2766 2767 Input Parameters: 2768 + comm - The communicator 2769 . dim - The spatial dimension 2770 - simplex - Flag for simplex, otherwise use a tensor-product cell 2771 2772 Output Parameter: 2773 . refdm - The reference cell 2774 2775 Level: intermediate 2776 2777 .keywords: reference cell 2778 .seealso: 2779 @*/ 2780 PetscErrorCode DMPlexCreateReferenceCell(MPI_Comm comm, PetscInt dim, PetscBool simplex, DM *refdm) 2781 { 2782 DM rdm; 2783 Vec coords; 2784 PetscErrorCode ierr; 2785 2786 PetscFunctionBeginUser; 2787 ierr = DMCreate(comm, &rdm);CHKERRQ(ierr); 2788 ierr = DMSetType(rdm, DMPLEX);CHKERRQ(ierr); 2789 ierr = DMSetDimension(rdm, dim);CHKERRQ(ierr); 2790 switch (dim) { 2791 case 0: 2792 { 2793 PetscInt numPoints[1] = {1}; 2794 PetscInt coneSize[1] = {0}; 2795 PetscInt cones[1] = {0}; 2796 PetscInt coneOrientations[1] = {0}; 2797 PetscScalar vertexCoords[1] = {0.0}; 2798 2799 ierr = DMPlexCreateFromDAG(rdm, 0, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr); 2800 } 2801 break; 2802 case 1: 2803 { 2804 PetscInt numPoints[2] = {2, 1}; 2805 PetscInt coneSize[3] = {2, 0, 0}; 2806 PetscInt cones[2] = {1, 2}; 2807 PetscInt coneOrientations[2] = {0, 0}; 2808 PetscScalar vertexCoords[2] = {-1.0, 1.0}; 2809 2810 ierr = DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr); 2811 } 2812 break; 2813 case 2: 2814 if (simplex) { 2815 PetscInt numPoints[2] = {3, 1}; 2816 PetscInt coneSize[4] = {3, 0, 0, 0}; 2817 PetscInt cones[3] = {1, 2, 3}; 2818 PetscInt coneOrientations[3] = {0, 0, 0}; 2819 PetscScalar vertexCoords[6] = {-1.0, -1.0, 1.0, -1.0, -1.0, 1.0}; 2820 2821 ierr = DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr); 2822 } else { 2823 PetscInt numPoints[2] = {4, 1}; 2824 PetscInt coneSize[5] = {4, 0, 0, 0, 0}; 2825 PetscInt cones[4] = {1, 2, 3, 4}; 2826 PetscInt coneOrientations[4] = {0, 0, 0, 0}; 2827 PetscScalar vertexCoords[8] = {-1.0, -1.0, 1.0, -1.0, 1.0, 1.0, -1.0, 1.0}; 2828 2829 ierr = DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr); 2830 } 2831 break; 2832 case 3: 2833 if (simplex) { 2834 PetscInt numPoints[2] = {4, 1}; 2835 PetscInt coneSize[5] = {4, 0, 0, 0, 0}; 2836 PetscInt cones[4] = {1, 3, 2, 4}; 2837 PetscInt coneOrientations[4] = {0, 0, 0, 0}; 2838 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}; 2839 2840 ierr = DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr); 2841 } else { 2842 PetscInt numPoints[2] = {8, 1}; 2843 PetscInt coneSize[9] = {8, 0, 0, 0, 0, 0, 0, 0, 0}; 2844 PetscInt cones[8] = {1, 4, 3, 2, 5, 6, 7, 8}; 2845 PetscInt coneOrientations[8] = {0, 0, 0, 0, 0, 0, 0, 0}; 2846 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, 2847 -1.0, -1.0, 1.0, 1.0, -1.0, 1.0, 1.0, 1.0, 1.0, -1.0, 1.0, 1.0}; 2848 2849 ierr = DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr); 2850 } 2851 break; 2852 default: 2853 SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Cannot create reference cell for dimension %d", dim); 2854 } 2855 *refdm = NULL; 2856 ierr = DMPlexInterpolate(rdm, refdm);CHKERRQ(ierr); 2857 if (rdm->coordinateDM) { 2858 DM ncdm; 2859 PetscSection cs; 2860 PetscInt pEnd = -1; 2861 2862 ierr = DMGetDefaultSection(rdm->coordinateDM, &cs);CHKERRQ(ierr); 2863 if (cs) {ierr = PetscSectionGetChart(cs, NULL, &pEnd);CHKERRQ(ierr);} 2864 if (pEnd >= 0) { 2865 ierr = DMClone(rdm->coordinateDM, &ncdm);CHKERRQ(ierr); 2866 ierr = DMSetDefaultSection(ncdm, cs);CHKERRQ(ierr); 2867 ierr = DMSetCoordinateDM(*refdm, ncdm);CHKERRQ(ierr); 2868 ierr = DMDestroy(&ncdm);CHKERRQ(ierr); 2869 } 2870 } 2871 ierr = DMGetCoordinatesLocal(rdm, &coords);CHKERRQ(ierr); 2872 if (coords) { 2873 ierr = DMSetCoordinatesLocal(*refdm, coords);CHKERRQ(ierr); 2874 } else { 2875 ierr = DMGetCoordinates(rdm, &coords);CHKERRQ(ierr); 2876 if (coords) {ierr = DMSetCoordinates(*refdm, coords);CHKERRQ(ierr);} 2877 } 2878 ierr = DMDestroy(&rdm);CHKERRQ(ierr); 2879 PetscFunctionReturn(0); 2880 } 2881