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