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