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 /* External function declarations here */ 980 extern PetscErrorCode DMCreateInterpolation_Plex(DM dmCoarse, DM dmFine, Mat *interpolation, Vec *scaling); 981 extern PetscErrorCode DMCreateInjection_Plex(DM dmCoarse, DM dmFine, Mat *mat); 982 extern PetscErrorCode DMCreateDefaultSection_Plex(DM dm); 983 extern PetscErrorCode DMCreateDefaultConstraints_Plex(DM dm); 984 extern PetscErrorCode DMCreateMatrix_Plex(DM dm, Mat *J); 985 extern PetscErrorCode DMCreateCoordinateDM_Plex(DM dm, DM *cdm); 986 PETSC_INTERN PetscErrorCode DMClone_Plex(DM dm, DM *newdm); 987 extern PetscErrorCode DMSetUp_Plex(DM dm); 988 extern PetscErrorCode DMDestroy_Plex(DM dm); 989 extern PetscErrorCode DMView_Plex(DM dm, PetscViewer viewer); 990 extern PetscErrorCode DMLoad_Plex(DM dm, PetscViewer viewer); 991 extern PetscErrorCode DMCreateSubDM_Plex(DM dm, PetscInt numFields, PetscInt fields[], IS *is, DM *subdm); 992 993 /* Replace dm with the contents of dmNew 994 - Share the DM_Plex structure 995 - Share the coordinates 996 - Share the SF 997 */ 998 static PetscErrorCode DMPlexReplace_Static(DM dm, DM dmNew) 999 { 1000 PetscSF sf; 1001 DM coordDM, coarseDM; 1002 Vec coords; 1003 const PetscReal *maxCell, *L; 1004 const DMBoundaryType *bd; 1005 PetscErrorCode ierr; 1006 1007 PetscFunctionBegin; 1008 ierr = DMGetPointSF(dmNew, &sf);CHKERRQ(ierr); 1009 ierr = DMSetPointSF(dm, sf);CHKERRQ(ierr); 1010 ierr = DMGetCoordinateDM(dmNew, &coordDM);CHKERRQ(ierr); 1011 ierr = DMGetCoordinatesLocal(dmNew, &coords);CHKERRQ(ierr); 1012 ierr = DMSetCoordinateDM(dm, coordDM);CHKERRQ(ierr); 1013 ierr = DMSetCoordinatesLocal(dm, coords);CHKERRQ(ierr); 1014 ierr = DMGetPeriodicity(dm, &maxCell, &L, &bd);CHKERRQ(ierr); 1015 if (L) {ierr = DMSetPeriodicity(dmNew, maxCell, L, bd);CHKERRQ(ierr);} 1016 ierr = DMDestroy_Plex(dm);CHKERRQ(ierr); 1017 dm->data = dmNew->data; 1018 ((DM_Plex *) dmNew->data)->refct++; 1019 dmNew->labels->refct++; 1020 if (!--(dm->labels->refct)) { 1021 DMLabelLink next = dm->labels->next; 1022 1023 /* destroy the labels */ 1024 while (next) { 1025 DMLabelLink tmp = next->next; 1026 1027 ierr = DMLabelDestroy(&next->label);CHKERRQ(ierr); 1028 ierr = PetscFree(next);CHKERRQ(ierr); 1029 next = tmp; 1030 } 1031 ierr = PetscFree(dm->labels);CHKERRQ(ierr); 1032 } 1033 dm->labels = dmNew->labels; 1034 dm->depthLabel = dmNew->depthLabel; 1035 ierr = DMGetCoarseDM(dmNew,&coarseDM);CHKERRQ(ierr); 1036 ierr = DMSetCoarseDM(dm,coarseDM);CHKERRQ(ierr); 1037 PetscFunctionReturn(0); 1038 } 1039 1040 /* Swap dm with the contents of dmNew 1041 - Swap the DM_Plex structure 1042 - Swap the coordinates 1043 - Swap the point PetscSF 1044 */ 1045 static PetscErrorCode DMPlexSwap_Static(DM dmA, DM dmB) 1046 { 1047 DM coordDMA, coordDMB; 1048 Vec coordsA, coordsB; 1049 PetscSF sfA, sfB; 1050 void *tmp; 1051 DMLabelLinkList listTmp; 1052 DMLabel depthTmp; 1053 PetscErrorCode ierr; 1054 1055 PetscFunctionBegin; 1056 ierr = DMGetPointSF(dmA, &sfA);CHKERRQ(ierr); 1057 ierr = DMGetPointSF(dmB, &sfB);CHKERRQ(ierr); 1058 ierr = PetscObjectReference((PetscObject) sfA);CHKERRQ(ierr); 1059 ierr = DMSetPointSF(dmA, sfB);CHKERRQ(ierr); 1060 ierr = DMSetPointSF(dmB, sfA);CHKERRQ(ierr); 1061 ierr = PetscObjectDereference((PetscObject) sfA);CHKERRQ(ierr); 1062 1063 ierr = DMGetCoordinateDM(dmA, &coordDMA);CHKERRQ(ierr); 1064 ierr = DMGetCoordinateDM(dmB, &coordDMB);CHKERRQ(ierr); 1065 ierr = PetscObjectReference((PetscObject) coordDMA);CHKERRQ(ierr); 1066 ierr = DMSetCoordinateDM(dmA, coordDMB);CHKERRQ(ierr); 1067 ierr = DMSetCoordinateDM(dmB, coordDMA);CHKERRQ(ierr); 1068 ierr = PetscObjectDereference((PetscObject) coordDMA);CHKERRQ(ierr); 1069 1070 ierr = DMGetCoordinatesLocal(dmA, &coordsA);CHKERRQ(ierr); 1071 ierr = DMGetCoordinatesLocal(dmB, &coordsB);CHKERRQ(ierr); 1072 ierr = PetscObjectReference((PetscObject) coordsA);CHKERRQ(ierr); 1073 ierr = DMSetCoordinatesLocal(dmA, coordsB);CHKERRQ(ierr); 1074 ierr = DMSetCoordinatesLocal(dmB, coordsA);CHKERRQ(ierr); 1075 ierr = PetscObjectDereference((PetscObject) coordsA);CHKERRQ(ierr); 1076 1077 tmp = dmA->data; 1078 dmA->data = dmB->data; 1079 dmB->data = tmp; 1080 listTmp = dmA->labels; 1081 dmA->labels = dmB->labels; 1082 dmB->labels = listTmp; 1083 depthTmp = dmA->depthLabel; 1084 dmA->depthLabel = dmB->depthLabel; 1085 dmB->depthLabel = depthTmp; 1086 PetscFunctionReturn(0); 1087 } 1088 1089 PetscErrorCode DMSetFromOptions_NonRefinement_Plex(PetscOptionItems *PetscOptionsObject,DM dm) 1090 { 1091 DM_Plex *mesh = (DM_Plex*) dm->data; 1092 PetscErrorCode ierr; 1093 1094 PetscFunctionBegin; 1095 /* Handle viewing */ 1096 ierr = PetscOptionsBool("-dm_plex_print_set_values", "Output all set values info", "DMView", PETSC_FALSE, &mesh->printSetValues, NULL);CHKERRQ(ierr); 1097 ierr = PetscOptionsInt("-dm_plex_print_fem", "Debug output level all fem computations", "DMView", 0, &mesh->printFEM, NULL);CHKERRQ(ierr); 1098 ierr = PetscOptionsReal("-dm_plex_print_tol", "Tolerance for FEM output", "DMView", mesh->printTol, &mesh->printTol, NULL);CHKERRQ(ierr); 1099 /* Point Location */ 1100 ierr = PetscOptionsBool("-dm_plex_hash_location", "Use grid hashing for point location", "DMView", PETSC_FALSE, &mesh->useHashLocation, NULL);CHKERRQ(ierr); 1101 /* Generation and remeshing */ 1102 ierr = PetscOptionsBool("-dm_plex_remesh_bd", "Allow changes to the boundary on remeshing", "DMView", PETSC_FALSE, &mesh->remeshBd, NULL);CHKERRQ(ierr); 1103 /* Projection behavior */ 1104 ierr = PetscOptionsInt("-dm_plex_max_projection_height", "Maxmimum mesh point height used to project locally", "DMPlexSetMaxProjectionHeight", 0, &mesh->maxProjectionHeight, NULL);CHKERRQ(ierr); 1105 ierr = PetscOptionsBool("-dm_plex_regular_refinement", "Use special nested projection algorithm for regular refinement", "DMPlexSetRegularRefinement", mesh->regularRefinement, &mesh->regularRefinement, NULL);CHKERRQ(ierr); 1106 PetscFunctionReturn(0); 1107 } 1108 1109 static PetscErrorCode DMSetFromOptions_Plex(PetscOptionItems *PetscOptionsObject,DM dm) 1110 { 1111 PetscInt refine = 0, coarsen = 0, r; 1112 PetscBool isHierarchy; 1113 PetscErrorCode ierr; 1114 1115 PetscFunctionBegin; 1116 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1117 ierr = PetscOptionsHead(PetscOptionsObject,"DMPlex Options");CHKERRQ(ierr); 1118 /* Handle DMPlex refinement */ 1119 ierr = PetscOptionsInt("-dm_refine", "The number of uniform refinements", "DMCreate", refine, &refine, NULL);CHKERRQ(ierr); 1120 ierr = PetscOptionsInt("-dm_refine_hierarchy", "The number of uniform refinements", "DMCreate", refine, &refine, &isHierarchy);CHKERRQ(ierr); 1121 if (refine) {ierr = DMPlexSetRefinementUniform(dm, PETSC_TRUE);CHKERRQ(ierr);} 1122 if (refine && isHierarchy) { 1123 DM *dms, coarseDM; 1124 1125 ierr = DMGetCoarseDM(dm, &coarseDM);CHKERRQ(ierr); 1126 ierr = PetscObjectReference((PetscObject)coarseDM);CHKERRQ(ierr); 1127 ierr = PetscMalloc1(refine,&dms);CHKERRQ(ierr); 1128 ierr = DMRefineHierarchy(dm, refine, dms);CHKERRQ(ierr); 1129 /* Total hack since we do not pass in a pointer */ 1130 ierr = DMPlexSwap_Static(dm, dms[refine-1]);CHKERRQ(ierr); 1131 if (refine == 1) { 1132 ierr = DMSetCoarseDM(dm, dms[0]);CHKERRQ(ierr); 1133 ierr = DMPlexSetRegularRefinement(dm, PETSC_TRUE);CHKERRQ(ierr); 1134 } else { 1135 ierr = DMSetCoarseDM(dm, dms[refine-2]);CHKERRQ(ierr); 1136 ierr = DMPlexSetRegularRefinement(dm, PETSC_TRUE);CHKERRQ(ierr); 1137 ierr = DMSetCoarseDM(dms[0], dms[refine-1]);CHKERRQ(ierr); 1138 ierr = DMPlexSetRegularRefinement(dms[0], PETSC_TRUE);CHKERRQ(ierr); 1139 } 1140 ierr = DMSetCoarseDM(dms[refine-1], coarseDM);CHKERRQ(ierr); 1141 ierr = PetscObjectDereference((PetscObject)coarseDM);CHKERRQ(ierr); 1142 /* Free DMs */ 1143 for (r = 0; r < refine; ++r) { 1144 ierr = DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dms[r]);CHKERRQ(ierr); 1145 ierr = DMDestroy(&dms[r]);CHKERRQ(ierr); 1146 } 1147 ierr = PetscFree(dms);CHKERRQ(ierr); 1148 } else { 1149 for (r = 0; r < refine; ++r) { 1150 DM refinedMesh; 1151 1152 ierr = DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dm);CHKERRQ(ierr); 1153 ierr = DMRefine(dm, PetscObjectComm((PetscObject) dm), &refinedMesh);CHKERRQ(ierr); 1154 /* Total hack since we do not pass in a pointer */ 1155 ierr = DMPlexReplace_Static(dm, refinedMesh);CHKERRQ(ierr); 1156 ierr = DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dm);CHKERRQ(ierr); 1157 ierr = DMDestroy(&refinedMesh);CHKERRQ(ierr); 1158 } 1159 } 1160 /* Handle DMPlex coarsening */ 1161 ierr = PetscOptionsInt("-dm_coarsen", "Coarsen the mesh", "DMCreate", coarsen, &coarsen, NULL);CHKERRQ(ierr); 1162 ierr = PetscOptionsInt("-dm_coarsen_hierarchy", "The number of coarsenings", "DMCreate", coarsen, &coarsen, &isHierarchy);CHKERRQ(ierr); 1163 if (coarsen && isHierarchy) { 1164 DM *dms; 1165 1166 ierr = PetscMalloc1(coarsen, &dms);CHKERRQ(ierr); 1167 ierr = DMCoarsenHierarchy(dm, coarsen, dms);CHKERRQ(ierr); 1168 /* Free DMs */ 1169 for (r = 0; r < coarsen; ++r) { 1170 ierr = DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dms[r]);CHKERRQ(ierr); 1171 ierr = DMDestroy(&dms[r]);CHKERRQ(ierr); 1172 } 1173 ierr = PetscFree(dms);CHKERRQ(ierr); 1174 } else { 1175 for (r = 0; r < coarsen; ++r) { 1176 DM coarseMesh; 1177 1178 ierr = DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dm);CHKERRQ(ierr); 1179 ierr = DMCoarsen(dm, PetscObjectComm((PetscObject) dm), &coarseMesh);CHKERRQ(ierr); 1180 /* Total hack since we do not pass in a pointer */ 1181 ierr = DMPlexReplace_Static(dm, coarseMesh);CHKERRQ(ierr); 1182 ierr = DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dm);CHKERRQ(ierr); 1183 ierr = DMDestroy(&coarseMesh);CHKERRQ(ierr); 1184 } 1185 } 1186 /* Handle */ 1187 ierr = DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dm);CHKERRQ(ierr); 1188 ierr = PetscOptionsTail();CHKERRQ(ierr); 1189 PetscFunctionReturn(0); 1190 } 1191 1192 static PetscErrorCode DMCreateGlobalVector_Plex(DM dm,Vec *vec) 1193 { 1194 PetscErrorCode ierr; 1195 1196 PetscFunctionBegin; 1197 ierr = DMCreateGlobalVector_Section_Private(dm,vec);CHKERRQ(ierr); 1198 /* ierr = VecSetOperation(*vec, VECOP_DUPLICATE, (void(*)(void)) VecDuplicate_MPI_DM);CHKERRQ(ierr); */ 1199 ierr = VecSetOperation(*vec, VECOP_VIEW, (void (*)(void)) VecView_Plex);CHKERRQ(ierr); 1200 ierr = VecSetOperation(*vec, VECOP_VIEWNATIVE, (void (*)(void)) VecView_Plex_Native);CHKERRQ(ierr); 1201 ierr = VecSetOperation(*vec, VECOP_LOAD, (void (*)(void)) VecLoad_Plex);CHKERRQ(ierr); 1202 ierr = VecSetOperation(*vec, VECOP_LOADNATIVE, (void (*)(void)) VecLoad_Plex_Native);CHKERRQ(ierr); 1203 PetscFunctionReturn(0); 1204 } 1205 1206 static PetscErrorCode DMCreateLocalVector_Plex(DM dm,Vec *vec) 1207 { 1208 PetscErrorCode ierr; 1209 1210 PetscFunctionBegin; 1211 ierr = DMCreateLocalVector_Section_Private(dm,vec);CHKERRQ(ierr); 1212 ierr = VecSetOperation(*vec, VECOP_VIEW, (void (*)(void)) VecView_Plex_Local);CHKERRQ(ierr); 1213 ierr = VecSetOperation(*vec, VECOP_LOAD, (void (*)(void)) VecLoad_Plex_Local);CHKERRQ(ierr); 1214 PetscFunctionReturn(0); 1215 } 1216 1217 static PetscErrorCode DMGetDimPoints_Plex(DM dm, PetscInt dim, PetscInt *pStart, PetscInt *pEnd) 1218 { 1219 PetscInt depth, d; 1220 PetscErrorCode ierr; 1221 1222 PetscFunctionBegin; 1223 ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 1224 if (depth == 1) { 1225 ierr = DMGetDimension(dm, &d);CHKERRQ(ierr); 1226 if (dim == 0) {ierr = DMPlexGetDepthStratum(dm, dim, pStart, pEnd);CHKERRQ(ierr);} 1227 else if (dim == d) {ierr = DMPlexGetDepthStratum(dm, 1, pStart, pEnd);CHKERRQ(ierr);} 1228 else {*pStart = 0; *pEnd = 0;} 1229 } else { 1230 ierr = DMPlexGetDepthStratum(dm, dim, pStart, pEnd);CHKERRQ(ierr); 1231 } 1232 PetscFunctionReturn(0); 1233 } 1234 1235 static PetscErrorCode DMGetNeighors_Plex(DM dm, PetscInt *nranks, const PetscMPIInt *ranks[]) 1236 { 1237 PetscSF sf; 1238 PetscErrorCode ierr; 1239 1240 PetscFunctionBegin; 1241 ierr = DMGetPointSF(dm, &sf);CHKERRQ(ierr); 1242 ierr = PetscSFGetRanks(sf, nranks, ranks, NULL, NULL, NULL);CHKERRQ(ierr); 1243 PetscFunctionReturn(0); 1244 } 1245 1246 static PetscErrorCode DMInitialize_Plex(DM dm) 1247 { 1248 PetscErrorCode ierr; 1249 1250 PetscFunctionBegin; 1251 dm->ops->view = DMView_Plex; 1252 dm->ops->load = DMLoad_Plex; 1253 dm->ops->setfromoptions = DMSetFromOptions_Plex; 1254 dm->ops->clone = DMClone_Plex; 1255 dm->ops->setup = DMSetUp_Plex; 1256 dm->ops->createdefaultsection = DMCreateDefaultSection_Plex; 1257 dm->ops->createdefaultconstraints = DMCreateDefaultConstraints_Plex; 1258 dm->ops->createglobalvector = DMCreateGlobalVector_Plex; 1259 dm->ops->createlocalvector = DMCreateLocalVector_Plex; 1260 dm->ops->getlocaltoglobalmapping = NULL; 1261 dm->ops->createfieldis = NULL; 1262 dm->ops->createcoordinatedm = DMCreateCoordinateDM_Plex; 1263 dm->ops->getcoloring = NULL; 1264 dm->ops->creatematrix = DMCreateMatrix_Plex; 1265 dm->ops->createinterpolation = DMCreateInterpolation_Plex; 1266 dm->ops->getaggregates = NULL; 1267 dm->ops->getinjection = DMCreateInjection_Plex; 1268 dm->ops->refine = DMRefine_Plex; 1269 dm->ops->coarsen = DMCoarsen_Plex; 1270 dm->ops->refinehierarchy = DMRefineHierarchy_Plex; 1271 dm->ops->coarsenhierarchy = DMCoarsenHierarchy_Plex; 1272 dm->ops->globaltolocalbegin = NULL; 1273 dm->ops->globaltolocalend = NULL; 1274 dm->ops->localtoglobalbegin = NULL; 1275 dm->ops->localtoglobalend = NULL; 1276 dm->ops->destroy = DMDestroy_Plex; 1277 dm->ops->createsubdm = DMCreateSubDM_Plex; 1278 dm->ops->getdimpoints = DMGetDimPoints_Plex; 1279 dm->ops->locatepoints = DMLocatePoints_Plex; 1280 dm->ops->projectfunctionlocal = DMProjectFunctionLocal_Plex; 1281 dm->ops->projectfunctionlabellocal = DMProjectFunctionLabelLocal_Plex; 1282 dm->ops->projectfieldlocal = DMProjectFieldLocal_Plex; 1283 dm->ops->projectfieldlabellocal = DMProjectFieldLabelLocal_Plex; 1284 dm->ops->computel2diff = DMComputeL2Diff_Plex; 1285 dm->ops->computel2gradientdiff = DMComputeL2GradientDiff_Plex; 1286 dm->ops->computel2fielddiff = DMComputeL2FieldDiff_Plex; 1287 dm->ops->getneighbors = DMGetNeighors_Plex; 1288 ierr = PetscObjectComposeFunction((PetscObject)dm,"DMAdaptLabel_C",DMAdaptLabel_Plex);CHKERRQ(ierr); 1289 PetscFunctionReturn(0); 1290 } 1291 1292 PETSC_INTERN PetscErrorCode DMClone_Plex(DM dm, DM *newdm) 1293 { 1294 DM_Plex *mesh = (DM_Plex *) dm->data; 1295 PetscErrorCode ierr; 1296 1297 PetscFunctionBegin; 1298 mesh->refct++; 1299 (*newdm)->data = mesh; 1300 ierr = PetscObjectChangeTypeName((PetscObject) *newdm, DMPLEX);CHKERRQ(ierr); 1301 ierr = DMInitialize_Plex(*newdm);CHKERRQ(ierr); 1302 PetscFunctionReturn(0); 1303 } 1304 1305 /*MC 1306 DMPLEX = "plex" - A DM object that encapsulates an unstructured mesh, or CW Complex, which can be expressed using a Hasse Diagram. 1307 In the local representation, Vecs contain all unknowns in the interior and shared boundary. This is 1308 specified by a PetscSection object. Ownership in the global representation is determined by 1309 ownership of the underlying DMPlex points. This is specified by another PetscSection object. 1310 1311 Level: intermediate 1312 1313 .seealso: DMType, DMPlexCreate(), DMCreate(), DMSetType() 1314 M*/ 1315 1316 PETSC_EXTERN PetscErrorCode DMCreate_Plex(DM dm) 1317 { 1318 DM_Plex *mesh; 1319 PetscInt unit, d; 1320 PetscErrorCode ierr; 1321 1322 PetscFunctionBegin; 1323 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1324 ierr = PetscNewLog(dm,&mesh);CHKERRQ(ierr); 1325 dm->dim = 0; 1326 dm->data = mesh; 1327 1328 mesh->refct = 1; 1329 ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm), &mesh->coneSection);CHKERRQ(ierr); 1330 mesh->maxConeSize = 0; 1331 mesh->cones = NULL; 1332 mesh->coneOrientations = NULL; 1333 ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm), &mesh->supportSection);CHKERRQ(ierr); 1334 mesh->maxSupportSize = 0; 1335 mesh->supports = NULL; 1336 mesh->refinementUniform = PETSC_TRUE; 1337 mesh->refinementLimit = -1.0; 1338 1339 mesh->facesTmp = NULL; 1340 1341 mesh->tetgenOpts = NULL; 1342 mesh->triangleOpts = NULL; 1343 ierr = PetscPartitionerCreate(PetscObjectComm((PetscObject)dm), &mesh->partitioner);CHKERRQ(ierr); 1344 ierr = PetscPartitionerSetTypeFromOptions_Internal(mesh->partitioner);CHKERRQ(ierr); 1345 mesh->remeshBd = PETSC_FALSE; 1346 1347 mesh->subpointMap = NULL; 1348 1349 for (unit = 0; unit < NUM_PETSC_UNITS; ++unit) mesh->scale[unit] = 1.0; 1350 1351 mesh->regularRefinement = PETSC_FALSE; 1352 mesh->depthState = -1; 1353 mesh->globalVertexNumbers = NULL; 1354 mesh->globalCellNumbers = NULL; 1355 mesh->anchorSection = NULL; 1356 mesh->anchorIS = NULL; 1357 mesh->createanchors = NULL; 1358 mesh->computeanchormatrix = NULL; 1359 mesh->parentSection = NULL; 1360 mesh->parents = NULL; 1361 mesh->childIDs = NULL; 1362 mesh->childSection = NULL; 1363 mesh->children = NULL; 1364 mesh->referenceTree = NULL; 1365 mesh->getchildsymmetry = NULL; 1366 for (d = 0; d < 8; ++d) mesh->hybridPointMax[d] = PETSC_DETERMINE; 1367 mesh->vtkCellHeight = 0; 1368 mesh->useCone = PETSC_FALSE; 1369 mesh->useClosure = PETSC_TRUE; 1370 mesh->useAnchors = PETSC_FALSE; 1371 1372 mesh->maxProjectionHeight = 0; 1373 1374 mesh->printSetValues = PETSC_FALSE; 1375 mesh->printFEM = 0; 1376 mesh->printTol = 1.0e-10; 1377 1378 ierr = DMInitialize_Plex(dm);CHKERRQ(ierr); 1379 PetscFunctionReturn(0); 1380 } 1381 1382 /*@ 1383 DMPlexCreate - Creates a DMPlex object, which encapsulates an unstructured mesh, or CW complex, which can be expressed using a Hasse Diagram. 1384 1385 Collective on MPI_Comm 1386 1387 Input Parameter: 1388 . comm - The communicator for the DMPlex object 1389 1390 Output Parameter: 1391 . mesh - The DMPlex object 1392 1393 Level: beginner 1394 1395 .keywords: DMPlex, create 1396 @*/ 1397 PetscErrorCode DMPlexCreate(MPI_Comm comm, DM *mesh) 1398 { 1399 PetscErrorCode ierr; 1400 1401 PetscFunctionBegin; 1402 PetscValidPointer(mesh,2); 1403 ierr = DMCreate(comm, mesh);CHKERRQ(ierr); 1404 ierr = DMSetType(*mesh, DMPLEX);CHKERRQ(ierr); 1405 PetscFunctionReturn(0); 1406 } 1407 1408 /* 1409 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 1410 */ 1411 static PetscErrorCode DMPlexBuildFromCellList_Parallel_Private(DM dm, PetscInt numCells, PetscInt numVertices, PetscInt numCorners, const int cells[], PetscSF *sfVert) 1412 { 1413 PetscSF sfPoint; 1414 PetscLayout vLayout; 1415 PetscHashI vhash; 1416 PetscSFNode *remoteVerticesAdj, *vertexLocal, *vertexOwner, *remoteVertex; 1417 const PetscInt *vrange; 1418 PetscInt numVerticesAdj, off, *verticesAdj, numVerticesGhost = 0, *localVertex, *cone, c, p, v, g; 1419 PETSC_UNUSED PetscHashIIter ret, iter; 1420 PetscMPIInt rank, size; 1421 PetscErrorCode ierr; 1422 1423 PetscFunctionBegin; 1424 ierr = MPI_Comm_rank(PetscObjectComm((PetscObject) dm), &rank);CHKERRQ(ierr); 1425 ierr = MPI_Comm_size(PetscObjectComm((PetscObject) dm), &size);CHKERRQ(ierr); 1426 /* Partition vertices */ 1427 ierr = PetscLayoutCreate(PetscObjectComm((PetscObject) dm), &vLayout);CHKERRQ(ierr); 1428 ierr = PetscLayoutSetLocalSize(vLayout, numVertices);CHKERRQ(ierr); 1429 ierr = PetscLayoutSetBlockSize(vLayout, 1);CHKERRQ(ierr); 1430 ierr = PetscLayoutSetUp(vLayout);CHKERRQ(ierr); 1431 ierr = PetscLayoutGetRanges(vLayout, &vrange);CHKERRQ(ierr); 1432 /* Count vertices and map them to procs */ 1433 PetscHashICreate(vhash); 1434 for (c = 0; c < numCells; ++c) { 1435 for (p = 0; p < numCorners; ++p) { 1436 PetscHashIPut(vhash, cells[c*numCorners+p], ret, iter); 1437 } 1438 } 1439 PetscHashISize(vhash, numVerticesAdj); 1440 ierr = PetscMalloc1(numVerticesAdj, &verticesAdj);CHKERRQ(ierr); 1441 off = 0; ierr = PetscHashIGetKeys(vhash, &off, verticesAdj);CHKERRQ(ierr); 1442 if (off != numVerticesAdj) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid number of local vertices %D should be %D", off, numVerticesAdj); 1443 ierr = PetscSortInt(numVerticesAdj, verticesAdj);CHKERRQ(ierr); 1444 ierr = PetscMalloc1(numVerticesAdj, &remoteVerticesAdj);CHKERRQ(ierr); 1445 for (v = 0; v < numVerticesAdj; ++v) { 1446 const PetscInt gv = verticesAdj[v]; 1447 PetscInt vrank; 1448 1449 ierr = PetscFindInt(gv, size+1, vrange, &vrank);CHKERRQ(ierr); 1450 vrank = vrank < 0 ? -(vrank+2) : vrank; 1451 remoteVerticesAdj[v].index = gv - vrange[vrank]; 1452 remoteVerticesAdj[v].rank = vrank; 1453 } 1454 /* Create cones */ 1455 ierr = DMPlexSetChart(dm, 0, numCells+numVerticesAdj);CHKERRQ(ierr); 1456 for (c = 0; c < numCells; ++c) {ierr = DMPlexSetConeSize(dm, c, numCorners);CHKERRQ(ierr);} 1457 ierr = DMSetUp(dm);CHKERRQ(ierr); 1458 ierr = DMGetWorkArray(dm, numCorners, PETSC_INT, &cone);CHKERRQ(ierr); 1459 for (c = 0; c < numCells; ++c) { 1460 for (p = 0; p < numCorners; ++p) { 1461 const PetscInt gv = cells[c*numCorners+p]; 1462 PetscInt lv; 1463 1464 ierr = PetscFindInt(gv, numVerticesAdj, verticesAdj, &lv);CHKERRQ(ierr); 1465 if (lv < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Could not find global vertex %D in local connectivity", gv); 1466 cone[p] = lv+numCells; 1467 } 1468 ierr = DMPlexSetCone(dm, c, cone);CHKERRQ(ierr); 1469 } 1470 ierr = DMRestoreWorkArray(dm, numCorners, PETSC_INT, &cone);CHKERRQ(ierr); 1471 /* Create SF for vertices */ 1472 ierr = PetscSFCreate(PetscObjectComm((PetscObject)dm), sfVert);CHKERRQ(ierr); 1473 ierr = PetscObjectSetName((PetscObject) *sfVert, "Vertex Ownership SF");CHKERRQ(ierr); 1474 ierr = PetscSFSetFromOptions(*sfVert);CHKERRQ(ierr); 1475 ierr = PetscSFSetGraph(*sfVert, numVertices, numVerticesAdj, NULL, PETSC_OWN_POINTER, remoteVerticesAdj, PETSC_OWN_POINTER);CHKERRQ(ierr); 1476 ierr = PetscFree(verticesAdj);CHKERRQ(ierr); 1477 /* Build pointSF */ 1478 ierr = PetscMalloc2(numVerticesAdj, &vertexLocal, numVertices, &vertexOwner);CHKERRQ(ierr); 1479 for (v = 0; v < numVerticesAdj; ++v) {vertexLocal[v].index = v+numCells; vertexLocal[v].rank = rank;} 1480 for (v = 0; v < numVertices; ++v) {vertexOwner[v].index = -1; vertexOwner[v].rank = -1;} 1481 ierr = PetscSFReduceBegin(*sfVert, MPI_2INT, vertexLocal, vertexOwner, MPI_MAXLOC);CHKERRQ(ierr); 1482 ierr = PetscSFReduceEnd(*sfVert, MPI_2INT, vertexLocal, vertexOwner, MPI_MAXLOC);CHKERRQ(ierr); 1483 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); 1484 ierr = PetscSFBcastBegin(*sfVert, MPI_2INT, vertexOwner, vertexLocal);CHKERRQ(ierr); 1485 ierr = PetscSFBcastEnd(*sfVert, MPI_2INT, vertexOwner, vertexLocal);CHKERRQ(ierr); 1486 for (v = 0; v < numVerticesAdj; ++v) if (vertexLocal[v].rank != rank) ++numVerticesGhost; 1487 ierr = PetscMalloc1(numVerticesGhost, &localVertex);CHKERRQ(ierr); 1488 ierr = PetscMalloc1(numVerticesGhost, &remoteVertex);CHKERRQ(ierr); 1489 for (v = 0, g = 0; v < numVerticesAdj; ++v) { 1490 if (vertexLocal[v].rank != rank) { 1491 localVertex[g] = v+numCells; 1492 remoteVertex[g].index = vertexLocal[v].index; 1493 remoteVertex[g].rank = vertexLocal[v].rank; 1494 ++g; 1495 } 1496 } 1497 ierr = PetscFree2(vertexLocal, vertexOwner);CHKERRQ(ierr); 1498 if (g != numVerticesGhost) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid number of vertex ghosts %D should be %D", g, numVerticesGhost); 1499 ierr = DMGetPointSF(dm, &sfPoint);CHKERRQ(ierr); 1500 ierr = PetscObjectSetName((PetscObject) sfPoint, "point SF");CHKERRQ(ierr); 1501 ierr = PetscSFSetGraph(sfPoint, numCells+numVerticesAdj, numVerticesGhost, localVertex, PETSC_OWN_POINTER, remoteVertex, PETSC_OWN_POINTER);CHKERRQ(ierr); 1502 ierr = PetscLayoutDestroy(&vLayout);CHKERRQ(ierr); 1503 PetscHashIDestroy(vhash); 1504 /* Fill in the rest of the topology structure */ 1505 ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr); 1506 ierr = DMPlexStratify(dm);CHKERRQ(ierr); 1507 PetscFunctionReturn(0); 1508 } 1509 1510 /* 1511 This takes as input the coordinates for each owned vertex 1512 */ 1513 static PetscErrorCode DMPlexBuildCoordinates_Parallel_Private(DM dm, PetscInt spaceDim, PetscInt numCells, PetscInt numV, PetscSF sfVert, const PetscReal vertexCoords[]) 1514 { 1515 PetscSection coordSection; 1516 Vec coordinates; 1517 PetscScalar *coords; 1518 PetscInt numVertices, numVerticesAdj, coordSize, v; 1519 PetscErrorCode ierr; 1520 1521 PetscFunctionBegin; 1522 ierr = PetscSFGetGraph(sfVert, &numVertices, &numVerticesAdj, NULL, NULL);CHKERRQ(ierr); 1523 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 1524 ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr); 1525 ierr = PetscSectionSetFieldComponents(coordSection, 0, spaceDim);CHKERRQ(ierr); 1526 ierr = PetscSectionSetChart(coordSection, numCells, numCells + numVerticesAdj);CHKERRQ(ierr); 1527 for (v = numCells; v < numCells+numVerticesAdj; ++v) { 1528 ierr = PetscSectionSetDof(coordSection, v, spaceDim);CHKERRQ(ierr); 1529 ierr = PetscSectionSetFieldDof(coordSection, v, 0, spaceDim);CHKERRQ(ierr); 1530 } 1531 ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 1532 ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 1533 ierr = VecCreate(PetscObjectComm((PetscObject)dm), &coordinates);CHKERRQ(ierr); 1534 ierr = VecSetBlockSize(coordinates, spaceDim);CHKERRQ(ierr); 1535 ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr); 1536 ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 1537 ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr); 1538 ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 1539 { 1540 MPI_Datatype coordtype; 1541 1542 /* Need a temp buffer for coords if we have complex/single */ 1543 ierr = MPI_Type_contiguous(spaceDim, MPIU_SCALAR, &coordtype);CHKERRQ(ierr); 1544 ierr = MPI_Type_commit(&coordtype);CHKERRQ(ierr); 1545 #if defined(PETSC_USE_COMPLEX) 1546 { 1547 PetscScalar *svertexCoords; 1548 PetscInt i; 1549 ierr = PetscMalloc1(numV*spaceDim,&svertexCoords);CHKERRQ(ierr); 1550 for (i=0; i<numV*spaceDim; i++) svertexCoords[i] = vertexCoords[i]; 1551 ierr = PetscSFBcastBegin(sfVert, coordtype, svertexCoords, coords);CHKERRQ(ierr); 1552 ierr = PetscSFBcastEnd(sfVert, coordtype, svertexCoords, coords);CHKERRQ(ierr); 1553 ierr = PetscFree(svertexCoords);CHKERRQ(ierr); 1554 } 1555 #else 1556 ierr = PetscSFBcastBegin(sfVert, coordtype, vertexCoords, coords);CHKERRQ(ierr); 1557 ierr = PetscSFBcastEnd(sfVert, coordtype, vertexCoords, coords);CHKERRQ(ierr); 1558 #endif 1559 ierr = MPI_Type_free(&coordtype);CHKERRQ(ierr); 1560 } 1561 ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 1562 ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr); 1563 ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 1564 PetscFunctionReturn(0); 1565 } 1566 1567 /*@C 1568 DMPlexCreateFromCellListParallel - This takes as input common mesh generator output, a list of the vertices for each cell, and produces a DM 1569 1570 Input Parameters: 1571 + comm - The communicator 1572 . dim - The topological dimension of the mesh 1573 . numCells - The number of cells owned by this process 1574 . numVertices - The number of vertices owned by this process 1575 . numCorners - The number of vertices for each cell 1576 . interpolate - Flag indicating that intermediate mesh entities (faces, edges) should be created automatically 1577 . cells - An array of numCells*numCorners numbers, the global vertex numbers for each cell 1578 . spaceDim - The spatial dimension used for coordinates 1579 - vertexCoords - An array of numVertices*spaceDim numbers, the coordinates of each vertex 1580 1581 Output Parameter: 1582 + dm - The DM 1583 - vertexSF - Optional, SF describing complete vertex ownership 1584 1585 Note: Two triangles sharing a face 1586 $ 1587 $ 2 1588 $ / | \ 1589 $ / | \ 1590 $ / | \ 1591 $ 0 0 | 1 3 1592 $ \ | / 1593 $ \ | / 1594 $ \ | / 1595 $ 1 1596 would have input 1597 $ numCells = 2, numVertices = 4 1598 $ cells = [0 1 2 1 3 2] 1599 $ 1600 which would result in the DMPlex 1601 $ 1602 $ 4 1603 $ / | \ 1604 $ / | \ 1605 $ / | \ 1606 $ 2 0 | 1 5 1607 $ \ | / 1608 $ \ | / 1609 $ \ | / 1610 $ 3 1611 1612 Level: beginner 1613 1614 .seealso: DMPlexCreateFromCellList(), DMPlexCreateFromDAG(), DMPlexCreate() 1615 @*/ 1616 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) 1617 { 1618 PetscSF sfVert; 1619 PetscErrorCode ierr; 1620 1621 PetscFunctionBegin; 1622 ierr = DMCreate(comm, dm);CHKERRQ(ierr); 1623 ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 1624 PetscValidLogicalCollectiveInt(*dm, dim, 2); 1625 PetscValidLogicalCollectiveInt(*dm, spaceDim, 8); 1626 ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr); 1627 ierr = DMPlexBuildFromCellList_Parallel_Private(*dm, numCells, numVertices, numCorners, cells, &sfVert);CHKERRQ(ierr); 1628 if (interpolate) { 1629 DM idm = NULL; 1630 1631 ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr); 1632 ierr = DMDestroy(dm);CHKERRQ(ierr); 1633 *dm = idm; 1634 } 1635 ierr = DMPlexBuildCoordinates_Parallel_Private(*dm, spaceDim, numCells, numVertices,sfVert, vertexCoords);CHKERRQ(ierr); 1636 if (vertexSF) *vertexSF = sfVert; 1637 else {ierr = PetscSFDestroy(&sfVert);CHKERRQ(ierr);} 1638 PetscFunctionReturn(0); 1639 } 1640 1641 /* 1642 This takes as input the common mesh generator output, a list of the vertices for each cell 1643 */ 1644 static PetscErrorCode DMPlexBuildFromCellList_Private(DM dm, PetscInt numCells, PetscInt numVertices, PetscInt numCorners, const int cells[]) 1645 { 1646 PetscInt *cone, c, p; 1647 PetscErrorCode ierr; 1648 1649 PetscFunctionBegin; 1650 ierr = DMPlexSetChart(dm, 0, numCells+numVertices);CHKERRQ(ierr); 1651 for (c = 0; c < numCells; ++c) { 1652 ierr = DMPlexSetConeSize(dm, c, numCorners);CHKERRQ(ierr); 1653 } 1654 ierr = DMSetUp(dm);CHKERRQ(ierr); 1655 ierr = DMGetWorkArray(dm, numCorners, PETSC_INT, &cone);CHKERRQ(ierr); 1656 for (c = 0; c < numCells; ++c) { 1657 for (p = 0; p < numCorners; ++p) { 1658 cone[p] = cells[c*numCorners+p]+numCells; 1659 } 1660 ierr = DMPlexSetCone(dm, c, cone);CHKERRQ(ierr); 1661 } 1662 ierr = DMRestoreWorkArray(dm, numCorners, PETSC_INT, &cone);CHKERRQ(ierr); 1663 ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr); 1664 ierr = DMPlexStratify(dm);CHKERRQ(ierr); 1665 PetscFunctionReturn(0); 1666 } 1667 1668 /* 1669 This takes as input the coordinates for each vertex 1670 */ 1671 static PetscErrorCode DMPlexBuildCoordinates_Private(DM dm, PetscInt spaceDim, PetscInt numCells, PetscInt numVertices, const double vertexCoords[]) 1672 { 1673 PetscSection coordSection; 1674 Vec coordinates; 1675 DM cdm; 1676 PetscScalar *coords; 1677 PetscInt v, d; 1678 PetscErrorCode ierr; 1679 1680 PetscFunctionBegin; 1681 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 1682 ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr); 1683 ierr = PetscSectionSetFieldComponents(coordSection, 0, spaceDim);CHKERRQ(ierr); 1684 ierr = PetscSectionSetChart(coordSection, numCells, numCells + numVertices);CHKERRQ(ierr); 1685 for (v = numCells; v < numCells+numVertices; ++v) { 1686 ierr = PetscSectionSetDof(coordSection, v, spaceDim);CHKERRQ(ierr); 1687 ierr = PetscSectionSetFieldDof(coordSection, v, 0, spaceDim);CHKERRQ(ierr); 1688 } 1689 ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 1690 1691 ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr); 1692 ierr = DMCreateLocalVector(cdm, &coordinates);CHKERRQ(ierr); 1693 ierr = VecSetBlockSize(coordinates, spaceDim);CHKERRQ(ierr); 1694 ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr); 1695 ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 1696 for (v = 0; v < numVertices; ++v) { 1697 for (d = 0; d < spaceDim; ++d) { 1698 coords[v*spaceDim+d] = vertexCoords[v*spaceDim+d]; 1699 } 1700 } 1701 ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 1702 ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr); 1703 ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 1704 PetscFunctionReturn(0); 1705 } 1706 1707 /*@C 1708 DMPlexCreateFromCellList - This takes as input common mesh generator output, a list of the vertices for each cell, and produces a DM 1709 1710 Input Parameters: 1711 + comm - The communicator 1712 . dim - The topological dimension of the mesh 1713 . numCells - The number of cells 1714 . numVertices - The number of vertices 1715 . numCorners - The number of vertices for each cell 1716 . interpolate - Flag indicating that intermediate mesh entities (faces, edges) should be created automatically 1717 . cells - An array of numCells*numCorners numbers, the vertices for each cell 1718 . spaceDim - The spatial dimension used for coordinates 1719 - vertexCoords - An array of numVertices*spaceDim numbers, the coordinates of each vertex 1720 1721 Output Parameter: 1722 . dm - The DM 1723 1724 Note: Two triangles sharing a face 1725 $ 1726 $ 2 1727 $ / | \ 1728 $ / | \ 1729 $ / | \ 1730 $ 0 0 | 1 3 1731 $ \ | / 1732 $ \ | / 1733 $ \ | / 1734 $ 1 1735 would have input 1736 $ numCells = 2, numVertices = 4 1737 $ cells = [0 1 2 1 3 2] 1738 $ 1739 which would result in the DMPlex 1740 $ 1741 $ 4 1742 $ / | \ 1743 $ / | \ 1744 $ / | \ 1745 $ 2 0 | 1 5 1746 $ \ | / 1747 $ \ | / 1748 $ \ | / 1749 $ 3 1750 1751 Level: beginner 1752 1753 .seealso: DMPlexCreateFromDAG(), DMPlexCreate() 1754 @*/ 1755 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) 1756 { 1757 PetscErrorCode ierr; 1758 1759 PetscFunctionBegin; 1760 ierr = DMCreate(comm, dm);CHKERRQ(ierr); 1761 ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 1762 ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr); 1763 ierr = DMPlexBuildFromCellList_Private(*dm, numCells, numVertices, numCorners, cells);CHKERRQ(ierr); 1764 if (interpolate) { 1765 DM idm = NULL; 1766 1767 ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr); 1768 ierr = DMDestroy(dm);CHKERRQ(ierr); 1769 *dm = idm; 1770 } 1771 ierr = DMPlexBuildCoordinates_Private(*dm, spaceDim, numCells, numVertices, vertexCoords);CHKERRQ(ierr); 1772 PetscFunctionReturn(0); 1773 } 1774 1775 /*@ 1776 DMPlexCreateFromDAG - This takes as input the adjacency-list representation of the Directed Acyclic Graph (Hasse Diagram) encoding a mesh, and produces a DM 1777 1778 Input Parameters: 1779 + dm - The empty DM object, usually from DMCreate() and DMSetDimension() 1780 . depth - The depth of the DAG 1781 . numPoints - The number of points at each depth 1782 . coneSize - The cone size of each point 1783 . cones - The concatenation of the cone points for each point, the cone list must be oriented correctly for each point 1784 . coneOrientations - The orientation of each cone point 1785 - vertexCoords - An array of numVertices*dim numbers, the coordinates of each vertex 1786 1787 Output Parameter: 1788 . dm - The DM 1789 1790 Note: Two triangles sharing a face would have input 1791 $ depth = 1, numPoints = [4 2], coneSize = [3 3 0 0 0 0] 1792 $ cones = [2 3 4 3 5 4], coneOrientations = [0 0 0 0 0 0] 1793 $ vertexCoords = [-1.0 0.0 0.0 -1.0 0.0 1.0 1.0 0.0] 1794 $ 1795 which would result in the DMPlex 1796 $ 1797 $ 4 1798 $ / | \ 1799 $ / | \ 1800 $ / | \ 1801 $ 2 0 | 1 5 1802 $ \ | / 1803 $ \ | / 1804 $ \ | / 1805 $ 3 1806 $ 1807 $ Notice that all points are numbered consecutively, unlikely DMPlexCreateFromCellList() 1808 1809 Level: advanced 1810 1811 .seealso: DMPlexCreateFromCellList(), DMPlexCreate() 1812 @*/ 1813 PetscErrorCode DMPlexCreateFromDAG(DM dm, PetscInt depth, const PetscInt numPoints[], const PetscInt coneSize[], const PetscInt cones[], const PetscInt coneOrientations[], const PetscScalar vertexCoords[]) 1814 { 1815 Vec coordinates; 1816 PetscSection coordSection; 1817 PetscScalar *coords; 1818 PetscInt coordSize, firstVertex = -1, pStart = 0, pEnd = 0, p, v, dim, dimEmbed, d, off; 1819 PetscErrorCode ierr; 1820 1821 PetscFunctionBegin; 1822 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 1823 ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr); 1824 if (dimEmbed < dim) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Embedding dimension %d cannot be less than intrinsic dimension %d",dimEmbed,dim); 1825 for (d = 0; d <= depth; ++d) pEnd += numPoints[d]; 1826 ierr = DMPlexSetChart(dm, pStart, pEnd);CHKERRQ(ierr); 1827 for (p = pStart; p < pEnd; ++p) { 1828 ierr = DMPlexSetConeSize(dm, p, coneSize[p-pStart]);CHKERRQ(ierr); 1829 if (firstVertex < 0 && !coneSize[p - pStart]) { 1830 firstVertex = p - pStart; 1831 } 1832 } 1833 if (firstVertex < 0 && numPoints[0]) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Expected %d vertices but could not find any", numPoints[0]); 1834 ierr = DMSetUp(dm);CHKERRQ(ierr); /* Allocate space for cones */ 1835 for (p = pStart, off = 0; p < pEnd; off += coneSize[p-pStart], ++p) { 1836 ierr = DMPlexSetCone(dm, p, &cones[off]);CHKERRQ(ierr); 1837 ierr = DMPlexSetConeOrientation(dm, p, &coneOrientations[off]);CHKERRQ(ierr); 1838 } 1839 ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr); 1840 ierr = DMPlexStratify(dm);CHKERRQ(ierr); 1841 /* Build coordinates */ 1842 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 1843 ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr); 1844 ierr = PetscSectionSetFieldComponents(coordSection, 0, dimEmbed);CHKERRQ(ierr); 1845 ierr = PetscSectionSetChart(coordSection, firstVertex, firstVertex+numPoints[0]);CHKERRQ(ierr); 1846 for (v = firstVertex; v < firstVertex+numPoints[0]; ++v) { 1847 ierr = PetscSectionSetDof(coordSection, v, dimEmbed);CHKERRQ(ierr); 1848 ierr = PetscSectionSetFieldDof(coordSection, v, 0, dimEmbed);CHKERRQ(ierr); 1849 } 1850 ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 1851 ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 1852 ierr = VecCreate(PETSC_COMM_SELF, &coordinates);CHKERRQ(ierr); 1853 ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr); 1854 ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 1855 ierr = VecSetBlockSize(coordinates, dimEmbed);CHKERRQ(ierr); 1856 ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr); 1857 ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 1858 for (v = 0; v < numPoints[0]; ++v) { 1859 PetscInt off; 1860 1861 ierr = PetscSectionGetOffset(coordSection, v+firstVertex, &off);CHKERRQ(ierr); 1862 for (d = 0; d < dimEmbed; ++d) { 1863 coords[off+d] = vertexCoords[v*dimEmbed+d]; 1864 } 1865 } 1866 ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 1867 ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr); 1868 ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 1869 PetscFunctionReturn(0); 1870 } 1871 1872 /*@C 1873 DMPlexCreateFromFile - This takes a filename and produces a DM 1874 1875 Input Parameters: 1876 + comm - The communicator 1877 . filename - A file name 1878 - interpolate - Flag to create intermediate mesh pieces (edges, faces) 1879 1880 Output Parameter: 1881 . dm - The DM 1882 1883 Level: beginner 1884 1885 .seealso: DMPlexCreateFromDAG(), DMPlexCreateFromCellList(), DMPlexCreate() 1886 @*/ 1887 PetscErrorCode DMPlexCreateFromFile(MPI_Comm comm, const char filename[], PetscBool interpolate, DM *dm) 1888 { 1889 const char *extGmsh = ".msh"; 1890 const char *extCGNS = ".cgns"; 1891 const char *extExodus = ".exo"; 1892 const char *extFluent = ".cas"; 1893 const char *extHDF5 = ".h5"; 1894 const char *extMed = ".med"; 1895 size_t len; 1896 PetscBool isGmsh, isCGNS, isExodus, isFluent, isHDF5, isMed; 1897 PetscMPIInt rank; 1898 PetscErrorCode ierr; 1899 1900 PetscFunctionBegin; 1901 PetscValidPointer(filename, 2); 1902 PetscValidPointer(dm, 4); 1903 ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 1904 ierr = PetscStrlen(filename, &len);CHKERRQ(ierr); 1905 if (!len) SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Filename must be a valid path"); 1906 ierr = PetscStrncmp(&filename[PetscMax(0,len-4)], extGmsh, 4, &isGmsh);CHKERRQ(ierr); 1907 ierr = PetscStrncmp(&filename[PetscMax(0,len-5)], extCGNS, 5, &isCGNS);CHKERRQ(ierr); 1908 ierr = PetscStrncmp(&filename[PetscMax(0,len-4)], extExodus, 4, &isExodus);CHKERRQ(ierr); 1909 ierr = PetscStrncmp(&filename[PetscMax(0,len-4)], extFluent, 4, &isFluent);CHKERRQ(ierr); 1910 ierr = PetscStrncmp(&filename[PetscMax(0,len-3)], extHDF5, 3, &isHDF5);CHKERRQ(ierr); 1911 ierr = PetscStrncmp(&filename[PetscMax(0,len-4)], extMed, 4, &isMed);CHKERRQ(ierr); 1912 if (isGmsh) { 1913 ierr = DMPlexCreateGmshFromFile(comm, filename, interpolate, dm);CHKERRQ(ierr); 1914 } else if (isCGNS) { 1915 ierr = DMPlexCreateCGNSFromFile(comm, filename, interpolate, dm);CHKERRQ(ierr); 1916 } else if (isExodus) { 1917 ierr = DMPlexCreateExodusFromFile(comm, filename, interpolate, dm);CHKERRQ(ierr); 1918 } else if (isFluent) { 1919 ierr = DMPlexCreateFluentFromFile(comm, filename, interpolate, dm);CHKERRQ(ierr); 1920 } else if (isHDF5) { 1921 PetscViewer viewer; 1922 1923 ierr = PetscViewerCreate(comm, &viewer);CHKERRQ(ierr); 1924 ierr = PetscViewerSetType(viewer, PETSCVIEWERHDF5);CHKERRQ(ierr); 1925 ierr = PetscViewerFileSetMode(viewer, FILE_MODE_READ);CHKERRQ(ierr); 1926 ierr = PetscViewerFileSetName(viewer, filename);CHKERRQ(ierr); 1927 ierr = DMCreate(comm, dm);CHKERRQ(ierr); 1928 ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 1929 ierr = DMLoad(*dm, viewer);CHKERRQ(ierr); 1930 ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 1931 } else if (isMed) { 1932 ierr = DMPlexCreateMedFromFile(comm, filename, interpolate, dm);CHKERRQ(ierr); 1933 } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cannot load file %s: unrecognized extension", filename); 1934 PetscFunctionReturn(0); 1935 } 1936 1937 /*@ 1938 DMPlexCreateReferenceCell - Create a DMPLEX with the appropriate FEM reference cell 1939 1940 Collective on comm 1941 1942 Input Parameters: 1943 + comm - The communicator 1944 . dim - The spatial dimension 1945 - simplex - Flag for simplex, otherwise use a tensor-product cell 1946 1947 Output Parameter: 1948 . refdm - The reference cell 1949 1950 Level: intermediate 1951 1952 .keywords: reference cell 1953 .seealso: 1954 @*/ 1955 PetscErrorCode DMPlexCreateReferenceCell(MPI_Comm comm, PetscInt dim, PetscBool simplex, DM *refdm) 1956 { 1957 DM rdm; 1958 Vec coords; 1959 PetscErrorCode ierr; 1960 1961 PetscFunctionBeginUser; 1962 ierr = DMCreate(comm, &rdm);CHKERRQ(ierr); 1963 ierr = DMSetType(rdm, DMPLEX);CHKERRQ(ierr); 1964 ierr = DMSetDimension(rdm, dim);CHKERRQ(ierr); 1965 switch (dim) { 1966 case 0: 1967 { 1968 PetscInt numPoints[1] = {1}; 1969 PetscInt coneSize[1] = {0}; 1970 PetscInt cones[1] = {0}; 1971 PetscInt coneOrientations[1] = {0}; 1972 PetscScalar vertexCoords[1] = {0.0}; 1973 1974 ierr = DMPlexCreateFromDAG(rdm, 0, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr); 1975 } 1976 break; 1977 case 1: 1978 { 1979 PetscInt numPoints[2] = {2, 1}; 1980 PetscInt coneSize[3] = {2, 0, 0}; 1981 PetscInt cones[2] = {1, 2}; 1982 PetscInt coneOrientations[2] = {0, 0}; 1983 PetscScalar vertexCoords[2] = {-1.0, 1.0}; 1984 1985 ierr = DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr); 1986 } 1987 break; 1988 case 2: 1989 if (simplex) { 1990 PetscInt numPoints[2] = {3, 1}; 1991 PetscInt coneSize[4] = {3, 0, 0, 0}; 1992 PetscInt cones[3] = {1, 2, 3}; 1993 PetscInt coneOrientations[3] = {0, 0, 0}; 1994 PetscScalar vertexCoords[6] = {-1.0, -1.0, 1.0, -1.0, -1.0, 1.0}; 1995 1996 ierr = DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr); 1997 } else { 1998 PetscInt numPoints[2] = {4, 1}; 1999 PetscInt coneSize[5] = {4, 0, 0, 0, 0}; 2000 PetscInt cones[4] = {1, 2, 3, 4}; 2001 PetscInt coneOrientations[4] = {0, 0, 0, 0}; 2002 PetscScalar vertexCoords[8] = {-1.0, -1.0, 1.0, -1.0, 1.0, 1.0, -1.0, 1.0}; 2003 2004 ierr = DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr); 2005 } 2006 break; 2007 case 3: 2008 if (simplex) { 2009 PetscInt numPoints[2] = {4, 1}; 2010 PetscInt coneSize[5] = {4, 0, 0, 0, 0}; 2011 PetscInt cones[4] = {1, 3, 2, 4}; 2012 PetscInt coneOrientations[4] = {0, 0, 0, 0}; 2013 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}; 2014 2015 ierr = DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr); 2016 } else { 2017 PetscInt numPoints[2] = {8, 1}; 2018 PetscInt coneSize[9] = {8, 0, 0, 0, 0, 0, 0, 0, 0}; 2019 PetscInt cones[8] = {1, 4, 3, 2, 5, 6, 7, 8}; 2020 PetscInt coneOrientations[8] = {0, 0, 0, 0, 0, 0, 0, 0}; 2021 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, 2022 -1.0, -1.0, 1.0, 1.0, -1.0, 1.0, 1.0, 1.0, 1.0, -1.0, 1.0, 1.0}; 2023 2024 ierr = DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr); 2025 } 2026 break; 2027 default: 2028 SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Cannot create reference cell for dimension %d", dim); 2029 } 2030 *refdm = NULL; 2031 ierr = DMPlexInterpolate(rdm, refdm);CHKERRQ(ierr); 2032 if (rdm->coordinateDM) { 2033 DM ncdm; 2034 PetscSection cs; 2035 PetscInt pEnd = -1; 2036 2037 ierr = DMGetDefaultSection(rdm->coordinateDM, &cs);CHKERRQ(ierr); 2038 if (cs) {ierr = PetscSectionGetChart(cs, NULL, &pEnd);CHKERRQ(ierr);} 2039 if (pEnd >= 0) { 2040 ierr = DMClone(rdm->coordinateDM, &ncdm);CHKERRQ(ierr); 2041 ierr = DMSetDefaultSection(ncdm, cs);CHKERRQ(ierr); 2042 ierr = DMSetCoordinateDM(*refdm, ncdm);CHKERRQ(ierr); 2043 ierr = DMDestroy(&ncdm);CHKERRQ(ierr); 2044 } 2045 } 2046 ierr = DMGetCoordinatesLocal(rdm, &coords);CHKERRQ(ierr); 2047 if (coords) { 2048 ierr = DMSetCoordinatesLocal(*refdm, coords);CHKERRQ(ierr); 2049 } else { 2050 ierr = DMGetCoordinates(rdm, &coords);CHKERRQ(ierr); 2051 if (coords) {ierr = DMSetCoordinates(*refdm, coords);CHKERRQ(ierr);} 2052 } 2053 ierr = DMDestroy(&rdm);CHKERRQ(ierr); 2054 PetscFunctionReturn(0); 2055 } 2056