1 #define PETSCDM_DLL 2 #include <petsc-private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 3 #include <petscdmda.h> 4 #include <petscsf.h> 5 6 #undef __FUNCT__ 7 #define __FUNCT__ "DMPlexCreateDoublet" 8 /*@ 9 DMPlexCreateDoublet - Creates a mesh of two cells of the specified type, optionally with later refinement. 10 11 Collective on MPI_Comm 12 13 Input Parameters: 14 + comm - The communicator for the DM object 15 . dim - The spatial dimension 16 . simplex - Flag for simplicial cells, otherwise they are tensor product cells 17 . interpolate - Flag to create intermediate mesh pieces (edges, faces) 18 . refinementUniform - Flag for uniform parallel refinement 19 - refinementLimit - A nonzero number indicates the largest admissible volume for a refined cell 20 21 Output Parameter: 22 . dm - The DM object 23 24 Level: beginner 25 26 .keywords: DM, create 27 .seealso: DMSetType(), DMCreate() 28 @*/ 29 PetscErrorCode DMPlexCreateDoublet(MPI_Comm comm, PetscInt dim, PetscBool simplex, PetscBool interpolate, PetscBool refinementUniform, PetscReal refinementLimit, DM *newdm) 30 { 31 DM dm; 32 PetscInt p; 33 PetscMPIInt rank; 34 PetscErrorCode ierr; 35 36 PetscFunctionBegin; 37 ierr = DMCreate(comm, &dm);CHKERRQ(ierr); 38 ierr = DMSetType(dm, DMPLEX);CHKERRQ(ierr); 39 ierr = DMSetDimension(dm, dim);CHKERRQ(ierr); 40 ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 41 switch (dim) { 42 case 2: 43 if (simplex) {ierr = PetscObjectSetName((PetscObject) dm, "triangular");CHKERRQ(ierr);} 44 else {ierr = PetscObjectSetName((PetscObject) dm, "quadrilateral");CHKERRQ(ierr);} 45 break; 46 case 3: 47 if (simplex) {ierr = PetscObjectSetName((PetscObject) dm, "tetrahedral");CHKERRQ(ierr);} 48 else {ierr = PetscObjectSetName((PetscObject) dm, "hexahedral");CHKERRQ(ierr);} 49 break; 50 default: 51 SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Cannot make meshes for dimension %d", dim); 52 } 53 if (rank) { 54 PetscInt numPoints[2] = {0, 0}; 55 ierr = DMPlexCreateFromDAG(dm, 1, numPoints, NULL, NULL, NULL, NULL);CHKERRQ(ierr); 56 } else { 57 switch (dim) { 58 case 2: 59 if (simplex) { 60 PetscInt numPoints[2] = {4, 2}; 61 PetscInt coneSize[6] = {3, 3, 0, 0, 0, 0}; 62 PetscInt cones[6] = {2, 3, 4, 5, 4, 3}; 63 PetscInt coneOrientations[6] = {0, 0, 0, 0, 0, 0}; 64 PetscScalar vertexCoords[8] = {-0.5, 0.5, 0.0, 0.0, 0.0, 1.0, 0.5, 0.5}; 65 PetscInt markerPoints[8] = {2, 1, 3, 1, 4, 1, 5, 1}; 66 67 ierr = DMPlexCreateFromDAG(dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr); 68 for (p = 0; p < 4; ++p) {ierr = DMPlexSetLabelValue(dm, "marker", markerPoints[p*2], markerPoints[p*2+1]);CHKERRQ(ierr);} 69 } else { 70 PetscInt numPoints[2] = {6, 2}; 71 PetscInt coneSize[8] = {4, 4, 0, 0, 0, 0, 0, 0}; 72 PetscInt cones[8] = {2, 3, 4, 5, 3, 6, 7, 4}; 73 PetscInt coneOrientations[8] = {0, 0, 0, 0, 0, 0, 0, 0}; 74 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}; 75 76 ierr = DMPlexCreateFromDAG(dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr); 77 } 78 break; 79 case 3: 80 if (simplex) { 81 PetscInt numPoints[2] = {5, 2}; 82 PetscInt coneSize[7] = {4, 4, 0, 0, 0, 0, 0}; 83 PetscInt cones[8] = {4, 3, 5, 2, 5, 3, 4, 6}; 84 PetscInt coneOrientations[8] = {0, 0, 0, 0, 0, 0, 0, 0}; 85 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}; 86 PetscInt markerPoints[10] = {2, 1, 3, 1, 4, 1, 5, 1, 6, 1}; 87 88 ierr = DMPlexCreateFromDAG(dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr); 89 for (p = 0; p < 5; ++p) {ierr = DMPlexSetLabelValue(dm, "marker", markerPoints[p*2], markerPoints[p*2+1]);CHKERRQ(ierr);} 90 } else { 91 PetscInt numPoints[2] = {12, 2}; 92 PetscInt coneSize[14] = {8, 8, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0}; 93 PetscInt cones[16] = {2, 3, 4, 5, 6, 7, 8, 9, 5, 4, 10, 11, 7, 12, 13, 8}; 94 PetscInt coneOrientations[16] = {0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0}; 95 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, 96 -1.0, -0.5, 0.5, 0.0, -0.5, 0.5, 0.0, 0.5, 0.5, -1.0, 0.5, 0.5, 97 1.0, 0.5, -0.5, 1.0, -0.5, -0.5, 1.0, -0.5, 0.5, 1.0, 0.5, 0.5}; 98 99 ierr = DMPlexCreateFromDAG(dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr); 100 } 101 break; 102 default: 103 SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Cannot make meshes for dimension %d", dim); 104 } 105 } 106 *newdm = dm; 107 if (refinementLimit > 0.0) { 108 DM rdm; 109 const char *name; 110 111 ierr = DMPlexSetRefinementUniform(*newdm, PETSC_FALSE);CHKERRQ(ierr); 112 ierr = DMPlexSetRefinementLimit(*newdm, refinementLimit);CHKERRQ(ierr); 113 ierr = DMRefine(*newdm, comm, &rdm);CHKERRQ(ierr); 114 ierr = PetscObjectGetName((PetscObject) *newdm, &name);CHKERRQ(ierr); 115 ierr = PetscObjectSetName((PetscObject) rdm, name);CHKERRQ(ierr); 116 ierr = DMDestroy(newdm);CHKERRQ(ierr); 117 *newdm = rdm; 118 } 119 if (interpolate) { 120 DM idm = NULL; 121 const char *name; 122 123 ierr = DMPlexInterpolate(*newdm, &idm);CHKERRQ(ierr); 124 ierr = PetscObjectGetName((PetscObject) *newdm, &name);CHKERRQ(ierr); 125 ierr = PetscObjectSetName((PetscObject) idm, name);CHKERRQ(ierr); 126 ierr = DMPlexCopyCoordinates(*newdm, idm);CHKERRQ(ierr); 127 ierr = DMPlexCopyLabels(*newdm, idm);CHKERRQ(ierr); 128 ierr = DMDestroy(newdm);CHKERRQ(ierr); 129 *newdm = idm; 130 } 131 { 132 DM refinedMesh = NULL; 133 DM distributedMesh = NULL; 134 135 /* Distribute mesh over processes */ 136 ierr = DMPlexDistribute(*newdm, 0, NULL, &distributedMesh);CHKERRQ(ierr); 137 if (distributedMesh) { 138 ierr = DMDestroy(newdm);CHKERRQ(ierr); 139 *newdm = distributedMesh; 140 } 141 if (refinementUniform) { 142 ierr = DMPlexSetRefinementUniform(*newdm, refinementUniform);CHKERRQ(ierr); 143 ierr = DMRefine(*newdm, comm, &refinedMesh);CHKERRQ(ierr); 144 if (refinedMesh) { 145 ierr = DMDestroy(newdm);CHKERRQ(ierr); 146 *newdm = refinedMesh; 147 } 148 } 149 } 150 PetscFunctionReturn(0); 151 } 152 153 #undef __FUNCT__ 154 #define __FUNCT__ "DMPlexCreateSquareBoundary" 155 /*@ 156 DMPlexCreateSquareBoundary - Creates a 1D mesh the is the boundary of a square lattice. 157 158 Collective on MPI_Comm 159 160 Input Parameters: 161 + comm - The communicator for the DM object 162 . lower - The lower left corner coordinates 163 . upper - The upper right corner coordinates 164 - edges - The number of cells in each direction 165 166 Output Parameter: 167 . dm - The DM object 168 169 Note: Here is the numbering returned for 2 cells in each direction: 170 $ 18--5-17--4--16 171 $ | | | 172 $ 6 10 3 173 $ | | | 174 $ 19-11-20--9--15 175 $ | | | 176 $ 7 8 2 177 $ | | | 178 $ 12--0-13--1--14 179 180 Level: beginner 181 182 .keywords: DM, create 183 .seealso: DMPlexCreateBoxMesh(), DMPlexCreateCubeBoundary(), DMSetType(), DMCreate() 184 @*/ 185 PetscErrorCode DMPlexCreateSquareBoundary(DM dm, const PetscReal lower[], const PetscReal upper[], const PetscInt edges[]) 186 { 187 PetscInt numVertices = (edges[0]+1)*(edges[1]+1); 188 PetscInt numEdges = edges[0]*(edges[1]+1) + (edges[0]+1)*edges[1]; 189 PetscInt markerTop = 1; 190 PetscInt markerBottom = 1; 191 PetscInt markerRight = 1; 192 PetscInt markerLeft = 1; 193 PetscBool markerSeparate = PETSC_FALSE; 194 Vec coordinates; 195 PetscSection coordSection; 196 PetscScalar *coords; 197 PetscInt coordSize; 198 PetscMPIInt rank; 199 PetscInt v, vx, vy; 200 PetscErrorCode ierr; 201 202 PetscFunctionBegin; 203 ierr = PetscOptionsGetBool(((PetscObject) dm)->prefix, "-dm_plex_separate_marker", &markerSeparate, NULL);CHKERRQ(ierr); 204 if (markerSeparate) { 205 markerTop = 1; 206 markerBottom = 0; 207 markerRight = 0; 208 markerLeft = 0; 209 } 210 ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);CHKERRQ(ierr); 211 if (!rank) { 212 PetscInt e, ex, ey; 213 214 ierr = DMPlexSetChart(dm, 0, numEdges+numVertices);CHKERRQ(ierr); 215 for (e = 0; e < numEdges; ++e) { 216 ierr = DMPlexSetConeSize(dm, e, 2);CHKERRQ(ierr); 217 } 218 ierr = DMSetUp(dm);CHKERRQ(ierr); /* Allocate space for cones */ 219 for (vx = 0; vx <= edges[0]; vx++) { 220 for (ey = 0; ey < edges[1]; ey++) { 221 PetscInt edge = vx*edges[1] + ey + edges[0]*(edges[1]+1); 222 PetscInt vertex = ey*(edges[0]+1) + vx + numEdges; 223 PetscInt cone[2]; 224 225 cone[0] = vertex; cone[1] = vertex+edges[0]+1; 226 ierr = DMPlexSetCone(dm, edge, cone);CHKERRQ(ierr); 227 if (vx == edges[0]) { 228 ierr = DMPlexSetLabelValue(dm, "marker", edge, markerRight);CHKERRQ(ierr); 229 ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerRight);CHKERRQ(ierr); 230 if (ey == edges[1]-1) { 231 ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerRight);CHKERRQ(ierr); 232 } 233 } else if (vx == 0) { 234 ierr = DMPlexSetLabelValue(dm, "marker", edge, markerLeft);CHKERRQ(ierr); 235 ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerLeft);CHKERRQ(ierr); 236 if (ey == edges[1]-1) { 237 ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerLeft);CHKERRQ(ierr); 238 } 239 } 240 } 241 } 242 for (vy = 0; vy <= edges[1]; vy++) { 243 for (ex = 0; ex < edges[0]; ex++) { 244 PetscInt edge = vy*edges[0] + ex; 245 PetscInt vertex = vy*(edges[0]+1) + ex + numEdges; 246 PetscInt cone[2]; 247 248 cone[0] = vertex; cone[1] = vertex+1; 249 ierr = DMPlexSetCone(dm, edge, cone);CHKERRQ(ierr); 250 if (vy == edges[1]) { 251 ierr = DMPlexSetLabelValue(dm, "marker", edge, markerTop);CHKERRQ(ierr); 252 ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerTop);CHKERRQ(ierr); 253 if (ex == edges[0]-1) { 254 ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerTop);CHKERRQ(ierr); 255 } 256 } else if (vy == 0) { 257 ierr = DMPlexSetLabelValue(dm, "marker", edge, markerBottom);CHKERRQ(ierr); 258 ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerBottom);CHKERRQ(ierr); 259 if (ex == edges[0]-1) { 260 ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerBottom);CHKERRQ(ierr); 261 } 262 } 263 } 264 } 265 } 266 ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr); 267 ierr = DMPlexStratify(dm);CHKERRQ(ierr); 268 /* Build coordinates */ 269 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 270 ierr = PetscSectionSetChart(coordSection, numEdges, numEdges + numVertices);CHKERRQ(ierr); 271 for (v = numEdges; v < numEdges+numVertices; ++v) { 272 ierr = PetscSectionSetDof(coordSection, v, 2);CHKERRQ(ierr); 273 } 274 ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 275 ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 276 ierr = VecCreate(PetscObjectComm((PetscObject)dm), &coordinates);CHKERRQ(ierr); 277 ierr = VecSetBlockSize(coordinates, 2);CHKERRQ(ierr); 278 ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr); 279 ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 280 ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr); 281 ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 282 for (vy = 0; vy <= edges[1]; ++vy) { 283 for (vx = 0; vx <= edges[0]; ++vx) { 284 coords[(vy*(edges[0]+1)+vx)*2+0] = lower[0] + ((upper[0] - lower[0])/edges[0])*vx; 285 coords[(vy*(edges[0]+1)+vx)*2+1] = lower[1] + ((upper[1] - lower[1])/edges[1])*vy; 286 } 287 } 288 ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 289 ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr); 290 ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 291 PetscFunctionReturn(0); 292 } 293 294 #undef __FUNCT__ 295 #define __FUNCT__ "DMPlexCreateCubeBoundary" 296 /*@ 297 DMPlexCreateCubeBoundary - Creates a 2D mesh the is the boundary of a cubic lattice. 298 299 Collective on MPI_Comm 300 301 Input Parameters: 302 + comm - The communicator for the DM object 303 . lower - The lower left front corner coordinates 304 . upper - The upper right back corner coordinates 305 - edges - The number of cells in each direction 306 307 Output Parameter: 308 . dm - The DM object 309 310 Level: beginner 311 312 .keywords: DM, create 313 .seealso: DMPlexCreateBoxMesh(), DMPlexCreateSquareBoundary(), DMSetType(), DMCreate() 314 @*/ 315 PetscErrorCode DMPlexCreateCubeBoundary(DM dm, const PetscReal lower[], const PetscReal upper[], const PetscInt faces[]) 316 { 317 PetscInt vertices[3], numVertices; 318 PetscInt numFaces = 2*faces[0]*faces[1] + 2*faces[1]*faces[2] + 2*faces[0]*faces[2]; 319 Vec coordinates; 320 PetscSection coordSection; 321 PetscScalar *coords; 322 PetscInt coordSize; 323 PetscMPIInt rank; 324 PetscInt v, vx, vy, vz; 325 PetscInt voffset, iface=0, cone[4]; 326 PetscErrorCode ierr; 327 328 PetscFunctionBegin; 329 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"); 330 ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);CHKERRQ(ierr); 331 vertices[0] = faces[0]+1; vertices[1] = faces[1]+1; vertices[2] = faces[2]+1; 332 numVertices = vertices[0]*vertices[1]*vertices[2]; 333 if (!rank) { 334 PetscInt f; 335 336 ierr = DMPlexSetChart(dm, 0, numFaces+numVertices);CHKERRQ(ierr); 337 for (f = 0; f < numFaces; ++f) { 338 ierr = DMPlexSetConeSize(dm, f, 4);CHKERRQ(ierr); 339 } 340 ierr = DMSetUp(dm);CHKERRQ(ierr); /* Allocate space for cones */ 341 for (v = 0; v < numFaces+numVertices; ++v) { 342 ierr = DMPlexSetLabelValue(dm, "marker", v, 1);CHKERRQ(ierr); 343 } 344 345 /* Side 0 (Top) */ 346 for (vy = 0; vy < faces[1]; vy++) { 347 for (vx = 0; vx < faces[0]; vx++) { 348 voffset = numFaces + vertices[0]*vertices[1]*(vertices[2]-1) + vy*vertices[0] + vx; 349 cone[0] = voffset; cone[1] = voffset+1; cone[2] = voffset+vertices[0]+1; cone[3] = voffset+vertices[0]; 350 ierr = DMPlexSetCone(dm, iface, cone);CHKERRQ(ierr); 351 iface++; 352 } 353 } 354 355 /* Side 1 (Bottom) */ 356 for (vy = 0; vy < faces[1]; vy++) { 357 for (vx = 0; vx < faces[0]; vx++) { 358 voffset = numFaces + vy*(faces[0]+1) + vx; 359 cone[0] = voffset+1; cone[1] = voffset; cone[2] = voffset+vertices[0]; cone[3] = voffset+vertices[0]+1; 360 ierr = DMPlexSetCone(dm, iface, cone);CHKERRQ(ierr); 361 iface++; 362 } 363 } 364 365 /* Side 2 (Front) */ 366 for (vz = 0; vz < faces[2]; vz++) { 367 for (vx = 0; vx < faces[0]; vx++) { 368 voffset = numFaces + vz*vertices[0]*vertices[1] + vx; 369 cone[0] = voffset; cone[1] = voffset+1; cone[2] = voffset+vertices[0]*vertices[1]+1; cone[3] = voffset+vertices[0]*vertices[1]; 370 ierr = DMPlexSetCone(dm, iface, cone);CHKERRQ(ierr); 371 iface++; 372 } 373 } 374 375 /* Side 3 (Back) */ 376 for (vz = 0; vz < faces[2]; vz++) { 377 for (vx = 0; vx < faces[0]; vx++) { 378 voffset = numFaces + vz*vertices[0]*vertices[1] + vertices[0]*(vertices[1]-1) + vx; 379 cone[0] = voffset+vertices[0]*vertices[1]; cone[1] = voffset+vertices[0]*vertices[1]+1; 380 cone[2] = voffset+1; cone[3] = voffset; 381 ierr = DMPlexSetCone(dm, iface, cone);CHKERRQ(ierr); 382 iface++; 383 } 384 } 385 386 /* Side 4 (Left) */ 387 for (vz = 0; vz < faces[2]; vz++) { 388 for (vy = 0; vy < faces[1]; vy++) { 389 voffset = numFaces + vz*vertices[0]*vertices[1] + vy*vertices[0]; 390 cone[0] = voffset; cone[1] = voffset+vertices[0]*vertices[1]; 391 cone[2] = voffset+vertices[0]*vertices[1]+vertices[0]; cone[3] = voffset+vertices[0]; 392 ierr = DMPlexSetCone(dm, iface, cone);CHKERRQ(ierr); 393 iface++; 394 } 395 } 396 397 /* Side 5 (Right) */ 398 for (vz = 0; vz < faces[2]; vz++) { 399 for (vy = 0; vy < faces[1]; vy++) { 400 voffset = numFaces + vz*vertices[0]*vertices[1] + vy*vertices[0] + vx; 401 cone[0] = voffset+vertices[0]*vertices[1]; cone[1] = voffset; 402 cone[2] = voffset+vertices[0]; cone[3] = voffset+vertices[0]*vertices[1]+vertices[0]; 403 ierr = DMPlexSetCone(dm, iface, cone);CHKERRQ(ierr); 404 iface++; 405 } 406 } 407 } 408 ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr); 409 ierr = DMPlexStratify(dm);CHKERRQ(ierr); 410 /* Build coordinates */ 411 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 412 ierr = PetscSectionSetChart(coordSection, numFaces, numFaces + numVertices);CHKERRQ(ierr); 413 for (v = numFaces; v < numFaces+numVertices; ++v) { 414 ierr = PetscSectionSetDof(coordSection, v, 3);CHKERRQ(ierr); 415 } 416 ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 417 ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 418 ierr = VecCreate(PetscObjectComm((PetscObject)dm), &coordinates);CHKERRQ(ierr); 419 ierr = VecSetBlockSize(coordinates, 3);CHKERRQ(ierr); 420 ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr); 421 ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 422 ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr); 423 ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 424 for (vz = 0; vz <= faces[2]; ++vz) { 425 for (vy = 0; vy <= faces[1]; ++vy) { 426 for (vx = 0; vx <= faces[0]; ++vx) { 427 coords[((vz*(faces[1]+1)+vy)*(faces[0]+1)+vx)*3+0] = lower[0] + ((upper[0] - lower[0])/faces[0])*vx; 428 coords[((vz*(faces[1]+1)+vy)*(faces[0]+1)+vx)*3+1] = lower[1] + ((upper[1] - lower[1])/faces[1])*vy; 429 coords[((vz*(faces[1]+1)+vy)*(faces[0]+1)+vx)*3+2] = lower[2] + ((upper[2] - lower[2])/faces[2])*vz; 430 } 431 } 432 } 433 ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 434 ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr); 435 ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 436 PetscFunctionReturn(0); 437 } 438 439 #undef __FUNCT__ 440 #define __FUNCT__ "DMPlexCreateCubeMesh_Internal" 441 static PetscErrorCode DMPlexCreateCubeMesh_Internal(DM dm, const PetscReal lower[], const PetscReal upper[], const PetscInt edges[], DMBoundaryType bdX, DMBoundaryType bdY, DMBoundaryType bdZ) 442 { 443 PetscInt markerTop = 1; 444 PetscInt markerBottom = 1; 445 PetscInt markerFront = 1; 446 PetscInt markerBack = 1; 447 PetscInt markerRight = 1; 448 PetscInt markerLeft = 1; 449 PetscInt dim; 450 PetscBool markerSeparate = PETSC_FALSE; 451 PetscMPIInt rank; 452 PetscErrorCode ierr; 453 454 PetscFunctionBegin; 455 ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr); 456 ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);CHKERRQ(ierr); 457 ierr = PetscOptionsGetBool(((PetscObject) dm)->prefix, "-dm_plex_separate_marker", &markerSeparate, NULL);CHKERRQ(ierr); 458 if (markerSeparate) { 459 switch (dim) { 460 case 2: 461 markerTop = 3; 462 markerBottom = 1; 463 markerRight = 2; 464 markerLeft = 4; 465 break; 466 case 3: 467 markerBottom = 1; 468 markerTop = 2; 469 markerFront = 3; 470 markerBack = 4; 471 markerRight = 5; 472 markerLeft = 6; 473 break; 474 default: 475 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Dimension %d not supported",dim); 476 break; 477 } 478 } 479 { 480 const PetscInt numXEdges = !rank ? edges[0] : 0; 481 const PetscInt numYEdges = !rank ? edges[1] : 0; 482 const PetscInt numZEdges = !rank ? edges[2] : 0; 483 const PetscInt numXVertices = !rank ? (bdX == DM_BOUNDARY_PERIODIC || bdX == DM_BOUNDARY_TWIST ? edges[0] : edges[0]+1) : 0; 484 const PetscInt numYVertices = !rank ? (bdY == DM_BOUNDARY_PERIODIC || bdY == DM_BOUNDARY_TWIST ? edges[1] : edges[1]+1) : 0; 485 const PetscInt numZVertices = !rank ? (bdZ == DM_BOUNDARY_PERIODIC || bdY == DM_BOUNDARY_TWIST ? edges[2] : edges[2]+1) : 0; 486 const PetscInt numCells = numXEdges*numYEdges*numZEdges; 487 const PetscInt numXFaces = numYEdges*numZEdges; 488 const PetscInt numYFaces = numXEdges*numZEdges; 489 const PetscInt numZFaces = numXEdges*numYEdges; 490 const PetscInt numTotXFaces = numXVertices*numXFaces; 491 const PetscInt numTotYFaces = numYVertices*numYFaces; 492 const PetscInt numTotZFaces = numZVertices*numZFaces; 493 const PetscInt numFaces = numTotXFaces + numTotYFaces + numTotZFaces; 494 const PetscInt numTotXEdges = numXEdges*numYVertices*numZVertices; 495 const PetscInt numTotYEdges = numYEdges*numXVertices*numZVertices; 496 const PetscInt numTotZEdges = numZEdges*numXVertices*numYVertices; 497 const PetscInt numVertices = numXVertices*numYVertices*numZVertices; 498 const PetscInt numEdges = numTotXEdges + numTotYEdges + numTotZEdges; 499 const PetscInt firstVertex = (dim == 2) ? numFaces : numCells; 500 const PetscInt firstXFace = (dim == 2) ? 0 : numCells + numVertices; 501 const PetscInt firstYFace = firstXFace + numTotXFaces; 502 const PetscInt firstZFace = firstYFace + numTotYFaces; 503 const PetscInt firstXEdge = numCells + numFaces + numVertices; 504 const PetscInt firstYEdge = firstXEdge + numTotXEdges; 505 const PetscInt firstZEdge = firstYEdge + numTotYEdges; 506 Vec coordinates; 507 PetscSection coordSection; 508 PetscScalar *coords; 509 PetscInt coordSize; 510 PetscInt v, vx, vy, vz; 511 PetscInt c, f, fx, fy, fz, e, ex, ey, ez; 512 513 ierr = DMPlexSetChart(dm, 0, numCells+numFaces+numEdges+numVertices);CHKERRQ(ierr); 514 for (c = 0; c < numCells; c++) { 515 ierr = DMPlexSetConeSize(dm, c, 6);CHKERRQ(ierr); 516 } 517 for (f = firstXFace; f < firstXFace+numFaces; ++f) { 518 ierr = DMPlexSetConeSize(dm, f, 4);CHKERRQ(ierr); 519 } 520 for (e = firstXEdge; e < firstXEdge+numEdges; ++e) { 521 ierr = DMPlexSetConeSize(dm, e, 2);CHKERRQ(ierr); 522 } 523 ierr = DMSetUp(dm);CHKERRQ(ierr); /* Allocate space for cones */ 524 /* Build cells */ 525 for (fz = 0; fz < numZEdges; ++fz) { 526 for (fy = 0; fy < numYEdges; ++fy) { 527 for (fx = 0; fx < numXEdges; ++fx) { 528 PetscInt cell = (fz*numYEdges + fy)*numXEdges + fx; 529 PetscInt faceB = firstZFace + (fy*numXEdges+fx)*numZVertices + fz; 530 PetscInt faceT = firstZFace + (fy*numXEdges+fx)*numZVertices + ((fz+1)%numZVertices); 531 PetscInt faceF = firstYFace + (fz*numXEdges+fx)*numYVertices + fy; 532 PetscInt faceK = firstYFace + (fz*numXEdges+fx)*numYVertices + ((fy+1)%numYVertices); 533 PetscInt faceL = firstXFace + (fz*numYEdges+fy)*numXVertices + fx; 534 PetscInt faceR = firstXFace + (fz*numYEdges+fy)*numXVertices + ((fx+1)%numXVertices); 535 /* B, T, F, K, R, L */ 536 PetscInt ornt[8] = {-4, 0, 0, -1, 0, -4}; /* ??? */ 537 PetscInt cone[8]; 538 539 /* no boundary twisting in 3D */ 540 cone[0] = faceB; cone[1] = faceT; cone[2] = faceF; cone[3] = faceK; cone[4] = faceR; cone[5] = faceL; 541 ierr = DMPlexSetCone(dm, cell, cone);CHKERRQ(ierr); 542 ierr = DMPlexSetConeOrientation(dm, cell, ornt);CHKERRQ(ierr); 543 } 544 } 545 } 546 /* Build x faces */ 547 for (fz = 0; fz < numZEdges; ++fz) { 548 for (fy = 0; fy < numYEdges; ++fy) { 549 for (fx = 0; fx < numXVertices; ++fx) { 550 PetscInt face = firstXFace + (fz*numYEdges+fy)*numXVertices + fx; 551 PetscInt edgeL = firstZEdge + ( fy* numXVertices+fx)*numZEdges + fz; 552 PetscInt edgeR = firstZEdge + (((fy+1)%numYVertices)*numXVertices+fx)*numZEdges + fz; 553 PetscInt edgeB = firstYEdge + ( fz* numXVertices+fx)*numYEdges + fy; 554 PetscInt edgeT = firstYEdge + (((fz+1)%numZVertices)*numXVertices+fx)*numYEdges + fy; 555 PetscInt ornt[4] = {0, 0, -2, -2}; 556 PetscInt cone[4]; 557 558 if (dim == 3) { 559 /* markers */ 560 if (bdX != DM_BOUNDARY_PERIODIC) { 561 if (fx == numXVertices-1) { 562 ierr = DMPlexSetLabelValue(dm, "marker", face, markerRight);CHKERRQ(ierr); 563 } 564 else if (fx == 0) { 565 ierr = DMPlexSetLabelValue(dm, "marker", face, markerLeft);CHKERRQ(ierr); 566 } 567 } 568 } 569 cone[0] = edgeB; cone[1] = edgeR; cone[2] = edgeT; cone[3] = edgeL; 570 ierr = DMPlexSetCone(dm, face, cone);CHKERRQ(ierr); 571 ierr = DMPlexSetConeOrientation(dm, face, ornt);CHKERRQ(ierr); 572 } 573 } 574 } 575 /* Build y faces */ 576 for (fz = 0; fz < numZEdges; ++fz) { 577 for (fx = 0; fx < numYEdges; ++fx) { 578 for (fy = 0; fy < numYVertices; ++fy) { 579 PetscInt face = firstYFace + (fz*numXEdges+fx)*numYVertices + fy; 580 PetscInt edgeL = firstZEdge + (fy*numXVertices+ fx )*numZEdges + fz; 581 PetscInt edgeR = firstZEdge + (fy*numXVertices+((fx+1)%numXVertices))*numZEdges + fz; 582 PetscInt edgeB = firstXEdge + ( fz *numYVertices+fy)*numXEdges + fx; 583 PetscInt edgeT = firstXEdge + (((fz+1)%numZVertices)*numYVertices+fy)*numXEdges + fx; 584 PetscInt ornt[4] = {0, 0, -2, -2}; 585 PetscInt cone[4]; 586 587 if (dim == 3) { 588 /* markers */ 589 if (bdY != DM_BOUNDARY_PERIODIC) { 590 if (fy == numYVertices-1) { 591 ierr = DMPlexSetLabelValue(dm, "marker", face, markerBack);CHKERRQ(ierr); 592 } 593 else if (fy == 0) { 594 ierr = DMPlexSetLabelValue(dm, "marker", face, markerFront);CHKERRQ(ierr); 595 } 596 } 597 } 598 cone[0] = edgeB; cone[1] = edgeR; cone[2] = edgeT; cone[3] = edgeL; 599 ierr = DMPlexSetCone(dm, face, cone);CHKERRQ(ierr); 600 ierr = DMPlexSetConeOrientation(dm, face, ornt);CHKERRQ(ierr); 601 } 602 } 603 } 604 /* Build z faces */ 605 for (fy = 0; fy < numYEdges; ++fy) { 606 for (fx = 0; fx < numXEdges; ++fx) { 607 for (fz = 0; fz < numZVertices; fz++) { 608 PetscInt face = firstZFace + (fy*numXEdges+fx)*numZVertices + fz; 609 PetscInt edgeL = firstYEdge + (fz*numXVertices+ fx )*numYEdges + fy; 610 PetscInt edgeR = firstYEdge + (fz*numXVertices+((fx+1)%numXVertices))*numYEdges + fy; 611 PetscInt edgeB = firstXEdge + (fz*numYVertices+ fy )*numXEdges + fx; 612 PetscInt edgeT = firstXEdge + (fz*numYVertices+((fy+1)%numYVertices))*numXEdges + fx; 613 PetscInt ornt[4] = {0, 0, -2, -2}; 614 PetscInt cone[4]; 615 616 if (dim == 2) { 617 if (bdX == DM_BOUNDARY_TWIST && fx == numXEdges-1) {edgeR += numYEdges-1-2*fy; ornt[1] = -2;} 618 if (bdY == DM_BOUNDARY_TWIST && fy == numYEdges-1) {edgeT += numXEdges-1-2*fx; ornt[2] = 0;} 619 } 620 else { 621 /* markers */ 622 if (bdZ != DM_BOUNDARY_PERIODIC) { 623 if (fz == numZVertices-1) { 624 ierr = DMPlexSetLabelValue(dm, "marker", face, markerTop);CHKERRQ(ierr); 625 } 626 else if (fz == 0) { 627 ierr = DMPlexSetLabelValue(dm, "marker", face, markerBottom);CHKERRQ(ierr); 628 } 629 } 630 } 631 cone[0] = edgeB; cone[1] = edgeR; cone[2] = edgeT; cone[3] = edgeL; 632 ierr = DMPlexSetCone(dm, face, cone);CHKERRQ(ierr); 633 ierr = DMPlexSetConeOrientation(dm, face, ornt);CHKERRQ(ierr); 634 } 635 } 636 } 637 /* Build Z edges*/ 638 for (vy = 0; vy < numYVertices; vy++) { 639 for (vx = 0; vx < numXVertices; vx++) { 640 for (ez = 0; ez < numZEdges; ez++) { 641 const PetscInt edge = firstZEdge + (vy*numXVertices+vx)*numZEdges + ez; 642 const PetscInt vertexB = firstVertex + ( ez *numYVertices+vy)*numXVertices + vx; 643 const PetscInt vertexT = firstVertex + (((ez+1)%numZVertices)*numYVertices+vy)*numXVertices + vx; 644 PetscInt cone[2]; 645 646 if (dim == 3) { 647 if (bdX != DM_BOUNDARY_PERIODIC) { 648 if (vx == numXVertices-1) { 649 ierr = DMPlexSetLabelValue(dm, "marker", edge, markerRight);CHKERRQ(ierr); 650 } 651 else if (vx == 0) { 652 ierr = DMPlexSetLabelValue(dm, "marker", edge, markerLeft);CHKERRQ(ierr); 653 } 654 } 655 if (bdY != DM_BOUNDARY_PERIODIC) { 656 if (vy == numYVertices-1) { 657 ierr = DMPlexSetLabelValue(dm, "marker", edge, markerBack);CHKERRQ(ierr); 658 } 659 else if (vy == 0) { 660 ierr = DMPlexSetLabelValue(dm, "marker", edge, markerFront);CHKERRQ(ierr); 661 } 662 } 663 } 664 cone[0] = vertexB; cone[1] = vertexT; 665 ierr = DMPlexSetCone(dm, edge, cone);CHKERRQ(ierr); 666 } 667 } 668 } 669 /* Build Y edges*/ 670 for (vz = 0; vz < numZVertices; vz++) { 671 for (vx = 0; vx < numXVertices; vx++) { 672 for (ey = 0; ey < numYEdges; ey++) { 673 const PetscInt nextv = (dim == 2 && bdY == DM_BOUNDARY_TWIST && ey == numYEdges-1) ? (numXVertices-vx-1) : (vz*numYVertices+((ey+1)%numYVertices))*numXVertices + vx; 674 const PetscInt edge = firstYEdge + (vz*numXVertices+vx)*numYEdges + ey; 675 const PetscInt vertexF = firstVertex + (vz*numYVertices+ey)*numXVertices + vx; 676 const PetscInt vertexK = firstVertex + nextv; 677 PetscInt cone[2]; 678 679 cone[0] = vertexF; cone[1] = vertexK; 680 ierr = DMPlexSetCone(dm, edge, cone);CHKERRQ(ierr); 681 if (dim == 2) { 682 if ((bdX != DM_BOUNDARY_PERIODIC) && (bdX != DM_BOUNDARY_TWIST)) { 683 if (vx == numXVertices-1) { 684 ierr = DMPlexSetLabelValue(dm, "marker", edge, markerRight);CHKERRQ(ierr); 685 ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerRight);CHKERRQ(ierr); 686 if (ey == numYEdges-1) { 687 ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerRight);CHKERRQ(ierr); 688 } 689 } 690 else if (vx == 0) { 691 ierr = DMPlexSetLabelValue(dm, "marker", edge, markerLeft);CHKERRQ(ierr); 692 ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerLeft);CHKERRQ(ierr); 693 if (ey == numYEdges-1) { 694 ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerLeft);CHKERRQ(ierr); 695 } 696 } 697 } 698 } 699 else { 700 if (bdX != DM_BOUNDARY_PERIODIC) { 701 if (vx == numXVertices-1) { 702 ierr = DMPlexSetLabelValue(dm, "marker", edge, markerRight);CHKERRQ(ierr); 703 } 704 else if (vx == 0) { 705 ierr = DMPlexSetLabelValue(dm, "marker", edge, markerLeft);CHKERRQ(ierr); 706 } 707 } 708 if (bdZ != DM_BOUNDARY_PERIODIC) { 709 if (vz == numZVertices-1) { 710 ierr = DMPlexSetLabelValue(dm, "marker", edge, markerTop);CHKERRQ(ierr); 711 } 712 else if (vz == 0) { 713 ierr = DMPlexSetLabelValue(dm, "marker", edge, markerBottom);CHKERRQ(ierr); 714 } 715 } 716 } 717 } 718 } 719 } 720 /* Build X edges*/ 721 for (vz = 0; vz < numZVertices; vz++) { 722 for (vy = 0; vy < numYVertices; vy++) { 723 for (ex = 0; ex < numXEdges; ex++) { 724 const PetscInt nextv = (dim == 2 && bdX == DM_BOUNDARY_TWIST && ex == numXEdges-1) ? (numYVertices-vy-1)*numXVertices : (vz*numYVertices+vy)*numXVertices + (ex+1)%numXVertices; 725 const PetscInt edge = firstXEdge + (vz*numYVertices+vy)*numXEdges + ex; 726 const PetscInt vertexL = firstVertex + (vz*numYVertices+vy)*numXVertices + ex; 727 const PetscInt vertexR = firstVertex + nextv; 728 PetscInt cone[2]; 729 730 cone[0] = vertexL; cone[1] = vertexR; 731 ierr = DMPlexSetCone(dm, edge, cone);CHKERRQ(ierr); 732 if (dim == 2) { 733 if ((bdY != DM_BOUNDARY_PERIODIC) && (bdY != DM_BOUNDARY_TWIST)) { 734 if (vy == numYVertices-1) { 735 ierr = DMPlexSetLabelValue(dm, "marker", edge, markerTop);CHKERRQ(ierr); 736 ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerTop);CHKERRQ(ierr); 737 if (ex == numXEdges-1) { 738 ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerTop);CHKERRQ(ierr); 739 } 740 } 741 else if (vy == 0) { 742 ierr = DMPlexSetLabelValue(dm, "marker", edge, markerBottom);CHKERRQ(ierr); 743 ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerBottom);CHKERRQ(ierr); 744 if (ex == numXEdges-1) { 745 ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerBottom);CHKERRQ(ierr); 746 } 747 } 748 } 749 } 750 else { 751 if (bdY != DM_BOUNDARY_PERIODIC) { 752 if (vy == numYVertices-1) { 753 ierr = DMPlexSetLabelValue(dm, "marker", edge, markerBack);CHKERRQ(ierr); 754 } 755 else if (vy == 0) { 756 ierr = DMPlexSetLabelValue(dm, "marker", edge, markerFront);CHKERRQ(ierr); 757 } 758 } 759 if (bdZ != DM_BOUNDARY_PERIODIC) { 760 if (vz == numZVertices-1) { 761 ierr = DMPlexSetLabelValue(dm, "marker", edge, markerTop);CHKERRQ(ierr); 762 } 763 else if (vz == 0) { 764 ierr = DMPlexSetLabelValue(dm, "marker", edge, markerBottom);CHKERRQ(ierr); 765 } 766 } 767 } 768 } 769 } 770 } 771 ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr); 772 ierr = DMPlexStratify(dm);CHKERRQ(ierr); 773 /* Build coordinates */ 774 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 775 ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr); 776 ierr = PetscSectionSetFieldComponents(coordSection, 0, dim);CHKERRQ(ierr); 777 ierr = PetscSectionSetChart(coordSection, firstVertex, firstVertex+numVertices);CHKERRQ(ierr); 778 for (v = firstVertex; v < firstVertex+numVertices; ++v) { 779 ierr = PetscSectionSetDof(coordSection, v, dim);CHKERRQ(ierr); 780 ierr = PetscSectionSetFieldDof(coordSection, v, 0, dim);CHKERRQ(ierr); 781 } 782 ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 783 ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 784 ierr = VecCreate(PetscObjectComm((PetscObject)dm), &coordinates);CHKERRQ(ierr); 785 ierr = VecSetBlockSize(coordinates, dim);CHKERRQ(ierr); 786 ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr); 787 ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 788 ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr); 789 ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 790 for (vz = 0; vz < numZVertices; ++vz) { 791 for (vy = 0; vy < numYVertices; ++vy) { 792 for (vx = 0; vx < numXVertices; ++vx) { 793 coords[((vz*numYVertices+vy)*numXVertices+vx)*dim+0] = lower[0] + ((upper[0] - lower[0])/numXEdges)*vx; 794 coords[((vz*numYVertices+vy)*numXVertices+vx)*dim+1] = lower[1] + ((upper[1] - lower[1])/numYEdges)*vy; 795 if (dim == 3) { 796 coords[((vz*numYVertices+vy)*numXVertices+vx)*dim+2] = lower[2] + ((upper[2] - lower[2])/numZEdges)*vz; 797 } 798 } 799 } 800 } 801 ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 802 ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr); 803 ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 804 } 805 PetscFunctionReturn(0); 806 } 807 808 #undef __FUNCT__ 809 #define __FUNCT__ "DMPlexCreateSquareMesh" 810 /*@ 811 DMPlexCreateSquareMesh - Creates a 2D mesh for a square lattice. 812 813 Collective on MPI_Comm 814 815 Input Parameters: 816 + comm - The communicator for the DM object 817 . lower - The lower left corner coordinates 818 . upper - The upper right corner coordinates 819 . edges - The number of cells in each direction 820 . bdX - The boundary type for the X direction 821 - bdY - The boundary type for the Y direction 822 823 Output Parameter: 824 . dm - The DM object 825 826 Note: Here is the numbering returned for 2 cells in each direction: 827 $ 22--8-23--9--24 828 $ | | | 829 $ 13 2 14 3 15 830 $ | | | 831 $ 19--6-20--7--21 832 $ | | | 833 $ 10 0 11 1 12 834 $ | | | 835 $ 16--4-17--5--18 836 837 Level: beginner 838 839 .keywords: DM, create 840 .seealso: DMPlexCreateBoxMesh(), DMPlexCreateSquareBoundary(), DMPlexCreateCubeBoundary(), DMSetType(), DMCreate() 841 @*/ 842 PetscErrorCode DMPlexCreateSquareMesh(DM dm, const PetscReal lower[], const PetscReal upper[], const PetscInt edges[], DMBoundaryType bdX, DMBoundaryType bdY) 843 { 844 PetscReal lower3[3], upper3[3]; 845 PetscInt edges3[3]; 846 PetscErrorCode ierr; 847 848 PetscFunctionBegin; 849 lower3[0] = lower[0]; lower3[1] = lower[1]; lower3[2] = 0.; 850 upper3[0] = upper[0]; upper3[1] = upper[1]; upper3[2] = 0.; 851 edges3[0] = edges[0]; edges3[1] = edges[1]; edges3[2] = 0; 852 ierr = DMPlexCreateCubeMesh_Internal(dm, lower3, upper3, edges3, bdX, bdY, DM_BOUNDARY_NONE);CHKERRQ(ierr); 853 PetscFunctionReturn(0); 854 } 855 856 #undef __FUNCT__ 857 #define __FUNCT__ "DMPlexCreateBoxMesh" 858 /*@ 859 DMPlexCreateBoxMesh - Creates a mesh on the tensor product of unit intervals (box) using simplices. 860 861 Collective on MPI_Comm 862 863 Input Parameters: 864 + comm - The communicator for the DM object 865 . dim - The spatial dimension 866 - interpolate - Flag to create intermediate mesh pieces (edges, faces) 867 868 Output Parameter: 869 . dm - The DM object 870 871 Level: beginner 872 873 .keywords: DM, create 874 .seealso: DMPlexCreateHexBoxMesh(), DMSetType(), DMCreate() 875 @*/ 876 PetscErrorCode DMPlexCreateBoxMesh(MPI_Comm comm, PetscInt dim, PetscBool interpolate, DM *dm) 877 { 878 DM boundary; 879 PetscErrorCode ierr; 880 881 PetscFunctionBegin; 882 PetscValidPointer(dm, 4); 883 ierr = DMCreate(comm, &boundary);CHKERRQ(ierr); 884 PetscValidLogicalCollectiveInt(boundary,dim,2); 885 ierr = DMSetType(boundary, DMPLEX);CHKERRQ(ierr); 886 ierr = DMSetDimension(boundary, dim-1);CHKERRQ(ierr); 887 switch (dim) { 888 case 2: 889 { 890 PetscReal lower[2] = {0.0, 0.0}; 891 PetscReal upper[2] = {1.0, 1.0}; 892 PetscInt edges[2] = {2, 2}; 893 894 ierr = DMPlexCreateSquareBoundary(boundary, lower, upper, edges);CHKERRQ(ierr); 895 break; 896 } 897 case 3: 898 { 899 PetscReal lower[3] = {0.0, 0.0, 0.0}; 900 PetscReal upper[3] = {1.0, 1.0, 1.0}; 901 PetscInt faces[3] = {1, 1, 1}; 902 903 ierr = DMPlexCreateCubeBoundary(boundary, lower, upper, faces);CHKERRQ(ierr); 904 break; 905 } 906 default: 907 SETERRQ1(comm, PETSC_ERR_SUP, "Dimension not supported: %d", dim); 908 } 909 ierr = DMPlexGenerate(boundary, NULL, interpolate, dm);CHKERRQ(ierr); 910 ierr = DMDestroy(&boundary);CHKERRQ(ierr); 911 PetscFunctionReturn(0); 912 } 913 914 #undef __FUNCT__ 915 #define __FUNCT__ "DMPlexCreateHexBoxMesh" 916 /*@ 917 DMPlexCreateHexBoxMesh - Creates a mesh on the tensor product of unit intervals (box) using hexahedra. 918 919 Collective on MPI_Comm 920 921 Input Parameters: 922 + comm - The communicator for the DM object 923 . dim - The spatial dimension 924 . periodicX - The boundary type for the X direction 925 . periodicY - The boundary type for the Y direction 926 . periodicZ - The boundary type for the Z direction 927 - cells - The number of cells in each direction 928 929 Output Parameter: 930 . dm - The DM object 931 932 Level: beginner 933 934 .keywords: DM, create 935 .seealso: DMPlexCreateBoxMesh(), DMSetType(), DMCreate() 936 @*/ 937 PetscErrorCode DMPlexCreateHexBoxMesh(MPI_Comm comm, PetscInt dim, const PetscInt cells[], DMBoundaryType periodicX, DMBoundaryType periodicY, DMBoundaryType periodicZ, DM *dm) 938 { 939 PetscErrorCode ierr; 940 941 PetscFunctionBegin; 942 PetscValidPointer(dm, 4); 943 ierr = DMCreate(comm, dm);CHKERRQ(ierr); 944 PetscValidLogicalCollectiveInt(*dm,dim,2); 945 ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 946 ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr); 947 switch (dim) { 948 case 2: 949 { 950 PetscReal lower[2] = {0.0, 0.0}; 951 PetscReal upper[2] = {1.0, 1.0}; 952 953 ierr = DMPlexCreateSquareMesh(*dm, lower, upper, cells, periodicX, periodicY);CHKERRQ(ierr); 954 break; 955 } 956 case 3: 957 { 958 PetscReal lower[3] = {0.0, 0.0, 0.0}; 959 PetscReal upper[3] = {1.0, 1.0, 1.0}; 960 961 ierr = DMPlexCreateCubeMesh_Internal(*dm, lower, upper, cells, periodicX, periodicY, periodicZ);CHKERRQ(ierr); 962 break; 963 } 964 default: 965 SETERRQ1(comm, PETSC_ERR_SUP, "Dimension not supported: %d", dim); 966 } 967 PetscFunctionReturn(0); 968 } 969 970 /* External function declarations here */ 971 extern PetscErrorCode DMCreateInterpolation_Plex(DM dmCoarse, DM dmFine, Mat *interpolation, Vec *scaling); 972 extern PetscErrorCode DMCreateInjection_Plex(DM dmCoarse, DM dmFine, VecScatter *ctx); 973 extern PetscErrorCode DMCreateDefaultSection_Plex(DM dm); 974 extern PetscErrorCode DMCreateDefaultConstraints_Plex(DM dm); 975 extern PetscErrorCode DMCreateMatrix_Plex(DM dm, Mat *J); 976 extern PetscErrorCode DMCreateCoordinateDM_Plex(DM dm, DM *cdm); 977 extern PetscErrorCode DMRefine_Plex(DM dm, MPI_Comm comm, DM *dmRefined); 978 extern PetscErrorCode DMCoarsen_Plex(DM dm, MPI_Comm comm, DM *dmCoarsened); 979 extern PetscErrorCode DMRefineHierarchy_Plex(DM dm, PetscInt nlevels, DM dmRefined[]); 980 extern PetscErrorCode DMClone_Plex(DM dm, DM *newdm); 981 extern PetscErrorCode DMSetUp_Plex(DM dm); 982 extern PetscErrorCode DMDestroy_Plex(DM dm); 983 extern PetscErrorCode DMView_Plex(DM dm, PetscViewer viewer); 984 extern PetscErrorCode DMLoad_Plex(DM dm, PetscViewer viewer); 985 extern PetscErrorCode DMCreateSubDM_Plex(DM dm, PetscInt numFields, PetscInt fields[], IS *is, DM *subdm); 986 extern PetscErrorCode DMLocatePoints_Plex(DM dm, Vec v, IS *cellIS); 987 988 #undef __FUNCT__ 989 #define __FUNCT__ "DMPlexReplace_Static" 990 /* Replace dm with the contents of dmNew 991 - Share the DM_Plex structure 992 - Share the coordinates 993 - Share the SF 994 */ 995 static PetscErrorCode DMPlexReplace_Static(DM dm, DM dmNew) 996 { 997 PetscSF sf; 998 DM coordDM; 999 Vec coords; 1000 const PetscReal *maxCell, *L; 1001 PetscErrorCode ierr; 1002 1003 PetscFunctionBegin; 1004 ierr = DMGetPointSF(dmNew, &sf);CHKERRQ(ierr); 1005 ierr = DMSetPointSF(dm, sf);CHKERRQ(ierr); 1006 ierr = DMGetCoordinateDM(dmNew, &coordDM);CHKERRQ(ierr); 1007 ierr = DMGetCoordinatesLocal(dmNew, &coords);CHKERRQ(ierr); 1008 ierr = DMSetCoordinateDM(dm, coordDM);CHKERRQ(ierr); 1009 ierr = DMSetCoordinatesLocal(dm, coords);CHKERRQ(ierr); 1010 ierr = DMGetPeriodicity(dm, &maxCell, &L);CHKERRQ(ierr); 1011 if (L) {ierr = DMSetPeriodicity(dmNew, maxCell, L);CHKERRQ(ierr);} 1012 ierr = DMDestroy_Plex(dm);CHKERRQ(ierr); 1013 dm->data = dmNew->data; 1014 ((DM_Plex *) dmNew->data)->refct++; 1015 PetscFunctionReturn(0); 1016 } 1017 1018 #undef __FUNCT__ 1019 #define __FUNCT__ "DMPlexSwap_Static" 1020 /* Swap dm with the contents of dmNew 1021 - Swap the DM_Plex structure 1022 - Swap the coordinates 1023 */ 1024 static PetscErrorCode DMPlexSwap_Static(DM dmA, DM dmB) 1025 { 1026 DM coordDMA, coordDMB; 1027 Vec coordsA, coordsB; 1028 void *tmp; 1029 PetscErrorCode ierr; 1030 1031 PetscFunctionBegin; 1032 ierr = DMGetCoordinateDM(dmA, &coordDMA);CHKERRQ(ierr); 1033 ierr = DMGetCoordinateDM(dmB, &coordDMB);CHKERRQ(ierr); 1034 ierr = PetscObjectReference((PetscObject) coordDMA);CHKERRQ(ierr); 1035 ierr = DMSetCoordinateDM(dmA, coordDMB);CHKERRQ(ierr); 1036 ierr = DMSetCoordinateDM(dmB, coordDMA);CHKERRQ(ierr); 1037 ierr = PetscObjectDereference((PetscObject) coordDMA);CHKERRQ(ierr); 1038 1039 ierr = DMGetCoordinatesLocal(dmA, &coordsA);CHKERRQ(ierr); 1040 ierr = DMGetCoordinatesLocal(dmB, &coordsB);CHKERRQ(ierr); 1041 ierr = PetscObjectReference((PetscObject) coordsA);CHKERRQ(ierr); 1042 ierr = DMSetCoordinatesLocal(dmA, coordsB);CHKERRQ(ierr); 1043 ierr = DMSetCoordinatesLocal(dmB, coordsA);CHKERRQ(ierr); 1044 ierr = PetscObjectDereference((PetscObject) coordsA);CHKERRQ(ierr); 1045 tmp = dmA->data; 1046 dmA->data = dmB->data; 1047 dmB->data = tmp; 1048 PetscFunctionReturn(0); 1049 } 1050 1051 #undef __FUNCT__ 1052 #define __FUNCT__ "DMSetFromOptions_NonRefinement_Plex" 1053 PetscErrorCode DMSetFromOptions_NonRefinement_Plex(DM dm) 1054 { 1055 DM_Plex *mesh = (DM_Plex*) dm->data; 1056 DMBoundary b; 1057 PetscErrorCode ierr; 1058 1059 PetscFunctionBegin; 1060 /* Handle boundary conditions */ 1061 ierr = PetscOptionsBegin(PetscObjectComm((PetscObject) dm), NULL, "Boundary condition options", "");CHKERRQ(ierr); 1062 for (b = mesh->boundary; b; b = b->next) { 1063 char optname[1024]; 1064 PetscInt ids[1024], len = 1024, i; 1065 PetscBool flg; 1066 1067 ierr = PetscSNPrintf(optname, sizeof(optname), "-bc_%s", b->name);CHKERRQ(ierr); 1068 ierr = PetscMemzero(ids, sizeof(ids));CHKERRQ(ierr); 1069 ierr = PetscOptionsIntArray(optname, "List of boundary IDs", "", ids, &len, &flg);CHKERRQ(ierr); 1070 if (flg) { 1071 DMLabel label; 1072 1073 ierr = DMPlexGetLabel(dm, b->labelname, &label);CHKERRQ(ierr); 1074 for (i = 0; i < len; ++i) { 1075 PetscBool has; 1076 1077 ierr = DMLabelHasValue(label, ids[i], &has);CHKERRQ(ierr); 1078 if (!has) SETERRQ2(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Boundary id %D is not present in the label %s", ids[i], b->name); 1079 } 1080 b->numids = len; 1081 ierr = PetscFree(b->ids);CHKERRQ(ierr); 1082 ierr = PetscMalloc1(len, &b->ids);CHKERRQ(ierr); 1083 ierr = PetscMemcpy(b->ids, ids, len*sizeof(PetscInt));CHKERRQ(ierr); 1084 } 1085 } 1086 ierr = PetscOptionsEnd();CHKERRQ(ierr); 1087 /* Handle viewing */ 1088 ierr = PetscOptionsBool("-dm_plex_print_set_values", "Output all set values info", "DMView", PETSC_FALSE, &mesh->printSetValues, NULL);CHKERRQ(ierr); 1089 ierr = PetscOptionsInt("-dm_plex_print_fem", "Debug output level all fem computations", "DMView", 0, &mesh->printFEM, NULL);CHKERRQ(ierr); 1090 ierr = PetscOptionsReal("-dm_plex_print_tol", "Tolerance for FEM output", "DMView", mesh->printTol, &mesh->printTol, NULL);CHKERRQ(ierr); 1091 /* Projection behavior */ 1092 ierr = PetscOptionsInt("-dm_plex_max_projection_height", "Maxmimum mesh point height used to project locally", "DMPlexSetMaxProjectionHeight", 0, &mesh->maxProjectionHeight, NULL);CHKERRQ(ierr); 1093 PetscFunctionReturn(0); 1094 } 1095 1096 #undef __FUNCT__ 1097 #define __FUNCT__ "DMSetFromOptions_Plex" 1098 PetscErrorCode DMSetFromOptions_Plex(DM dm) 1099 { 1100 PetscInt refine = 0, r; 1101 PetscBool isHierarchy; 1102 PetscErrorCode ierr; 1103 1104 PetscFunctionBegin; 1105 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1106 ierr = PetscOptionsHead("DMPlex Options");CHKERRQ(ierr); 1107 /* Handle DMPlex refinement */ 1108 ierr = PetscOptionsInt("-dm_refine", "The number of uniform refinements", "DMCreate", refine, &refine, NULL);CHKERRQ(ierr); 1109 ierr = PetscOptionsInt("-dm_refine_hierarchy", "The number of uniform refinements", "DMCreate", refine, &refine, &isHierarchy);CHKERRQ(ierr); 1110 if (refine) {ierr = DMPlexSetRefinementUniform(dm, PETSC_TRUE);CHKERRQ(ierr);} 1111 if (refine && isHierarchy) { 1112 DM *dms; 1113 1114 ierr = PetscMalloc1(refine,&dms);CHKERRQ(ierr); 1115 ierr = DMRefineHierarchy(dm, refine, dms);CHKERRQ(ierr); 1116 /* Total hack since we do not pass in a pointer */ 1117 ierr = DMPlexSwap_Static(dm, dms[refine-1]);CHKERRQ(ierr); 1118 if (refine == 1) { 1119 ierr = DMPlexSetCoarseDM(dm, dms[0]);CHKERRQ(ierr); 1120 } else { 1121 ierr = DMPlexSetCoarseDM(dm, dms[refine-2]);CHKERRQ(ierr); 1122 ierr = DMPlexSetCoarseDM(dms[0], dms[refine-1]);CHKERRQ(ierr); 1123 } 1124 /* Free DMs */ 1125 for (r = 0; r < refine; ++r) { 1126 ierr = DMSetFromOptions_NonRefinement_Plex(dm);CHKERRQ(ierr); 1127 ierr = DMDestroy(&dms[r]);CHKERRQ(ierr); 1128 } 1129 ierr = PetscFree(dms);CHKERRQ(ierr); 1130 } else { 1131 for (r = 0; r < refine; ++r) { 1132 DM refinedMesh; 1133 1134 ierr = DMSetFromOptions_NonRefinement_Plex(dm);CHKERRQ(ierr); 1135 ierr = DMRefine(dm, PetscObjectComm((PetscObject) dm), &refinedMesh);CHKERRQ(ierr); 1136 /* Total hack since we do not pass in a pointer */ 1137 ierr = DMPlexReplace_Static(dm, refinedMesh);CHKERRQ(ierr); 1138 ierr = DMDestroy(&refinedMesh);CHKERRQ(ierr); 1139 } 1140 } 1141 ierr = DMSetFromOptions_NonRefinement_Plex(dm);CHKERRQ(ierr); 1142 ierr = PetscOptionsTail();CHKERRQ(ierr); 1143 PetscFunctionReturn(0); 1144 } 1145 1146 #undef __FUNCT__ 1147 #define __FUNCT__ "DMCreateGlobalVector_Plex" 1148 static PetscErrorCode DMCreateGlobalVector_Plex(DM dm,Vec *vec) 1149 { 1150 PetscErrorCode ierr; 1151 1152 PetscFunctionBegin; 1153 ierr = DMCreateGlobalVector_Section_Private(dm,vec);CHKERRQ(ierr); 1154 /* ierr = VecSetOperation(*vec, VECOP_DUPLICATE, (void(*)(void)) VecDuplicate_MPI_DM);CHKERRQ(ierr); */ 1155 ierr = VecSetOperation(*vec, VECOP_VIEW, (void (*)(void)) VecView_Plex);CHKERRQ(ierr); 1156 ierr = VecSetOperation(*vec, VECOP_LOAD, (void (*)(void)) VecLoad_Plex);CHKERRQ(ierr); 1157 PetscFunctionReturn(0); 1158 } 1159 1160 #undef __FUNCT__ 1161 #define __FUNCT__ "DMCreateLocalVector_Plex" 1162 static PetscErrorCode DMCreateLocalVector_Plex(DM dm,Vec *vec) 1163 { 1164 PetscErrorCode ierr; 1165 1166 PetscFunctionBegin; 1167 ierr = DMCreateLocalVector_Section_Private(dm,vec);CHKERRQ(ierr); 1168 ierr = VecSetOperation(*vec, VECOP_VIEW, (void (*)(void)) VecView_Plex_Local);CHKERRQ(ierr); 1169 ierr = VecSetOperation(*vec, VECOP_LOAD, (void (*)(void)) VecLoad_Plex_Local);CHKERRQ(ierr); 1170 PetscFunctionReturn(0); 1171 } 1172 1173 #undef __FUNCT__ 1174 #define __FUNCT__ "DMGetDimPoints_Plex" 1175 static PetscErrorCode DMGetDimPoints_Plex(DM dm, PetscInt dim, PetscInt *pStart, PetscInt *pEnd) 1176 { 1177 PetscInt depth, d; 1178 PetscErrorCode ierr; 1179 1180 PetscFunctionBegin; 1181 ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 1182 if (depth == 1) { 1183 ierr = DMGetDimension(dm, &d);CHKERRQ(ierr); 1184 if (dim == 0) {ierr = DMPlexGetDepthStratum(dm, dim, pStart, pEnd);CHKERRQ(ierr);} 1185 else if (dim == d) {ierr = DMPlexGetDepthStratum(dm, 1, pStart, pEnd);CHKERRQ(ierr);} 1186 else {*pStart = 0; *pEnd = 0;} 1187 } else { 1188 ierr = DMPlexGetDepthStratum(dm, dim, pStart, pEnd);CHKERRQ(ierr); 1189 } 1190 PetscFunctionReturn(0); 1191 } 1192 1193 #undef __FUNCT__ 1194 #define __FUNCT__ "DMInitialize_Plex" 1195 PetscErrorCode DMInitialize_Plex(DM dm) 1196 { 1197 PetscFunctionBegin; 1198 dm->ops->view = DMView_Plex; 1199 dm->ops->load = DMLoad_Plex; 1200 dm->ops->setfromoptions = DMSetFromOptions_Plex; 1201 dm->ops->clone = DMClone_Plex; 1202 dm->ops->setup = DMSetUp_Plex; 1203 dm->ops->createdefaultsection = DMCreateDefaultSection_Plex; 1204 dm->ops->createdefaultconstraints = DMCreateDefaultConstraints_Plex; 1205 dm->ops->createglobalvector = DMCreateGlobalVector_Plex; 1206 dm->ops->createlocalvector = DMCreateLocalVector_Plex; 1207 dm->ops->getlocaltoglobalmapping = NULL; 1208 dm->ops->createfieldis = NULL; 1209 dm->ops->createcoordinatedm = DMCreateCoordinateDM_Plex; 1210 dm->ops->getcoloring = NULL; 1211 dm->ops->creatematrix = DMCreateMatrix_Plex; 1212 dm->ops->createinterpolation = DMCreateInterpolation_Plex; 1213 dm->ops->getaggregates = NULL; 1214 dm->ops->getinjection = DMCreateInjection_Plex; 1215 dm->ops->refine = DMRefine_Plex; 1216 dm->ops->coarsen = DMCoarsen_Plex; 1217 dm->ops->refinehierarchy = DMRefineHierarchy_Plex; 1218 dm->ops->coarsenhierarchy = NULL; 1219 dm->ops->globaltolocalbegin = NULL; 1220 dm->ops->globaltolocalend = NULL; 1221 dm->ops->localtoglobalbegin = NULL; 1222 dm->ops->localtoglobalend = NULL; 1223 dm->ops->destroy = DMDestroy_Plex; 1224 dm->ops->createsubdm = DMCreateSubDM_Plex; 1225 dm->ops->getdimpoints = DMGetDimPoints_Plex; 1226 dm->ops->locatepoints = DMLocatePoints_Plex; 1227 PetscFunctionReturn(0); 1228 } 1229 1230 #undef __FUNCT__ 1231 #define __FUNCT__ "DMClone_Plex" 1232 PetscErrorCode DMClone_Plex(DM dm, DM *newdm) 1233 { 1234 DM_Plex *mesh = (DM_Plex *) dm->data; 1235 PetscErrorCode ierr; 1236 1237 PetscFunctionBegin; 1238 mesh->refct++; 1239 (*newdm)->data = mesh; 1240 ierr = PetscObjectChangeTypeName((PetscObject) *newdm, DMPLEX);CHKERRQ(ierr); 1241 ierr = DMInitialize_Plex(*newdm);CHKERRQ(ierr); 1242 PetscFunctionReturn(0); 1243 } 1244 1245 /*MC 1246 DMPLEX = "plex" - A DM object that encapsulates an unstructured mesh, or CW Complex, which can be expressed using a Hasse Diagram. 1247 In the local representation, Vecs contain all unknowns in the interior and shared boundary. This is 1248 specified by a PetscSection object. Ownership in the global representation is determined by 1249 ownership of the underlying DMPlex points. This is specified by another PetscSection object. 1250 1251 Level: intermediate 1252 1253 .seealso: DMType, DMPlexCreate(), DMCreate(), DMSetType() 1254 M*/ 1255 1256 #undef __FUNCT__ 1257 #define __FUNCT__ "DMCreate_Plex" 1258 PETSC_EXTERN PetscErrorCode DMCreate_Plex(DM dm) 1259 { 1260 DM_Plex *mesh; 1261 PetscInt unit, d; 1262 PetscErrorCode ierr; 1263 1264 PetscFunctionBegin; 1265 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1266 ierr = PetscNewLog(dm,&mesh);CHKERRQ(ierr); 1267 dm->dim = 0; 1268 dm->data = mesh; 1269 1270 mesh->refct = 1; 1271 ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm), &mesh->coneSection);CHKERRQ(ierr); 1272 mesh->maxConeSize = 0; 1273 mesh->cones = NULL; 1274 mesh->coneOrientations = NULL; 1275 ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm), &mesh->supportSection);CHKERRQ(ierr); 1276 mesh->maxSupportSize = 0; 1277 mesh->supports = NULL; 1278 mesh->refinementUniform = PETSC_TRUE; 1279 mesh->refinementLimit = -1.0; 1280 1281 mesh->facesTmp = NULL; 1282 1283 mesh->tetgenOpts = NULL; 1284 mesh->triangleOpts = NULL; 1285 ierr = PetscPartitionerCreate(PetscObjectComm((PetscObject)dm), &mesh->partitioner);CHKERRQ(ierr); 1286 ierr = PetscPartitionerSetTypeFromOptions_Internal(mesh->partitioner);CHKERRQ(ierr); 1287 1288 mesh->subpointMap = NULL; 1289 1290 for (unit = 0; unit < NUM_PETSC_UNITS; ++unit) mesh->scale[unit] = 1.0; 1291 1292 mesh->labels = NULL; 1293 mesh->depthLabel = NULL; 1294 mesh->globalVertexNumbers = NULL; 1295 mesh->globalCellNumbers = NULL; 1296 mesh->anchorSection = NULL; 1297 mesh->anchorIS = NULL; 1298 mesh->createanchors = NULL; 1299 mesh->computeanchormatrix = NULL; 1300 mesh->parentSection = NULL; 1301 mesh->parents = NULL; 1302 mesh->childIDs = NULL; 1303 mesh->childSection = NULL; 1304 mesh->children = NULL; 1305 mesh->referenceTree = NULL; 1306 mesh->getchildsymmetry = NULL; 1307 for (d = 0; d < 8; ++d) mesh->hybridPointMax[d] = PETSC_DETERMINE; 1308 mesh->vtkCellHeight = 0; 1309 mesh->useCone = PETSC_FALSE; 1310 mesh->useClosure = PETSC_TRUE; 1311 mesh->useAnchors = PETSC_FALSE; 1312 1313 mesh->maxProjectionHeight = 0; 1314 1315 mesh->printSetValues = PETSC_FALSE; 1316 mesh->printFEM = 0; 1317 mesh->printTol = 1.0e-10; 1318 1319 ierr = DMInitialize_Plex(dm);CHKERRQ(ierr); 1320 PetscFunctionReturn(0); 1321 } 1322 1323 #undef __FUNCT__ 1324 #define __FUNCT__ "DMPlexCreate" 1325 /*@ 1326 DMPlexCreate - Creates a DMPlex object, which encapsulates an unstructured mesh, or CW complex, which can be expressed using a Hasse Diagram. 1327 1328 Collective on MPI_Comm 1329 1330 Input Parameter: 1331 . comm - The communicator for the DMPlex object 1332 1333 Output Parameter: 1334 . mesh - The DMPlex object 1335 1336 Level: beginner 1337 1338 .keywords: DMPlex, create 1339 @*/ 1340 PetscErrorCode DMPlexCreate(MPI_Comm comm, DM *mesh) 1341 { 1342 PetscErrorCode ierr; 1343 1344 PetscFunctionBegin; 1345 PetscValidPointer(mesh,2); 1346 ierr = DMCreate(comm, mesh);CHKERRQ(ierr); 1347 ierr = DMSetType(*mesh, DMPLEX);CHKERRQ(ierr); 1348 PetscFunctionReturn(0); 1349 } 1350 1351 #undef __FUNCT__ 1352 #define __FUNCT__ "DMPlexBuildFromCellList_Private" 1353 /* 1354 This takes as input the common mesh generator output, a list of the vertices for each cell 1355 */ 1356 PetscErrorCode DMPlexBuildFromCellList_Private(DM dm, PetscInt numCells, PetscInt numVertices, PetscInt numCorners, const int cells[]) 1357 { 1358 PetscInt *cone, c, p; 1359 PetscErrorCode ierr; 1360 1361 PetscFunctionBegin; 1362 ierr = DMPlexSetChart(dm, 0, numCells+numVertices);CHKERRQ(ierr); 1363 for (c = 0; c < numCells; ++c) { 1364 ierr = DMPlexSetConeSize(dm, c, numCorners);CHKERRQ(ierr); 1365 } 1366 ierr = DMSetUp(dm);CHKERRQ(ierr); 1367 ierr = DMGetWorkArray(dm, numCorners, PETSC_INT, &cone);CHKERRQ(ierr); 1368 for (c = 0; c < numCells; ++c) { 1369 for (p = 0; p < numCorners; ++p) { 1370 cone[p] = cells[c*numCorners+p]+numCells; 1371 } 1372 ierr = DMPlexSetCone(dm, c, cone);CHKERRQ(ierr); 1373 } 1374 ierr = DMRestoreWorkArray(dm, numCorners, PETSC_INT, &cone);CHKERRQ(ierr); 1375 ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr); 1376 ierr = DMPlexStratify(dm);CHKERRQ(ierr); 1377 PetscFunctionReturn(0); 1378 } 1379 1380 #undef __FUNCT__ 1381 #define __FUNCT__ "DMPlexBuildCoordinates_Private" 1382 /* 1383 This takes as input the coordinates for each vertex 1384 */ 1385 PetscErrorCode DMPlexBuildCoordinates_Private(DM dm, PetscInt spaceDim, PetscInt numCells, PetscInt numVertices, const double vertexCoords[]) 1386 { 1387 PetscSection coordSection; 1388 Vec coordinates; 1389 PetscScalar *coords; 1390 PetscInt coordSize, v, d; 1391 PetscErrorCode ierr; 1392 1393 PetscFunctionBegin; 1394 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 1395 ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr); 1396 ierr = PetscSectionSetFieldComponents(coordSection, 0, spaceDim);CHKERRQ(ierr); 1397 ierr = PetscSectionSetChart(coordSection, numCells, numCells + numVertices);CHKERRQ(ierr); 1398 for (v = numCells; v < numCells+numVertices; ++v) { 1399 ierr = PetscSectionSetDof(coordSection, v, spaceDim);CHKERRQ(ierr); 1400 ierr = PetscSectionSetFieldDof(coordSection, v, 0, spaceDim);CHKERRQ(ierr); 1401 } 1402 ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 1403 ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 1404 ierr = VecCreate(PetscObjectComm((PetscObject)dm), &coordinates);CHKERRQ(ierr); 1405 ierr = VecSetBlockSize(coordinates, spaceDim);CHKERRQ(ierr); 1406 ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr); 1407 ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 1408 ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr); 1409 ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 1410 for (v = 0; v < numVertices; ++v) { 1411 for (d = 0; d < spaceDim; ++d) { 1412 coords[v*spaceDim+d] = vertexCoords[v*spaceDim+d]; 1413 } 1414 } 1415 ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 1416 ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr); 1417 ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 1418 PetscFunctionReturn(0); 1419 } 1420 1421 #undef __FUNCT__ 1422 #define __FUNCT__ "DMPlexCreateFromCellList" 1423 /*@C 1424 DMPlexCreateFromCellList - This takes as input common mesh generator output, a list of the vertices for each cell, and produces a DM 1425 1426 Input Parameters: 1427 + comm - The communicator 1428 . dim - The topological dimension of the mesh 1429 . numCells - The number of cells 1430 . numVertices - The number of vertices 1431 . numCorners - The number of vertices for each cell 1432 . interpolate - Flag indicating that intermediate mesh entities (faces, edges) should be created automatically 1433 . cells - An array of numCells*numCorners numbers, the vertices for each cell 1434 . spaceDim - The spatial dimension used for coordinates 1435 - vertexCoords - An array of numVertices*spaceDim numbers, the coordinates of each vertex 1436 1437 Output Parameter: 1438 . dm - The DM 1439 1440 Note: Two triangles sharing a face 1441 $ 1442 $ 2 1443 $ / | \ 1444 $ / | \ 1445 $ / | \ 1446 $ 0 0 | 1 3 1447 $ \ | / 1448 $ \ | / 1449 $ \ | / 1450 $ 1 1451 would have input 1452 $ numCells = 2, numVertices = 4 1453 $ cells = [0 1 2 1 3 2] 1454 $ 1455 which would result in the DMPlex 1456 $ 1457 $ 4 1458 $ / | \ 1459 $ / | \ 1460 $ / | \ 1461 $ 2 0 | 1 5 1462 $ \ | / 1463 $ \ | / 1464 $ \ | / 1465 $ 3 1466 1467 Level: beginner 1468 1469 .seealso: DMPlexCreateFromDAG(), DMPlexCreate() 1470 @*/ 1471 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) 1472 { 1473 PetscErrorCode ierr; 1474 1475 PetscFunctionBegin; 1476 ierr = DMCreate(comm, dm);CHKERRQ(ierr); 1477 ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 1478 ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr); 1479 ierr = DMPlexBuildFromCellList_Private(*dm, numCells, numVertices, numCorners, cells);CHKERRQ(ierr); 1480 if (interpolate) { 1481 DM idm = NULL; 1482 1483 ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr); 1484 ierr = DMDestroy(dm);CHKERRQ(ierr); 1485 *dm = idm; 1486 } 1487 ierr = DMPlexBuildCoordinates_Private(*dm, spaceDim, numCells, numVertices, vertexCoords);CHKERRQ(ierr); 1488 PetscFunctionReturn(0); 1489 } 1490 1491 #undef __FUNCT__ 1492 #define __FUNCT__ "DMPlexCreateFromDAG" 1493 /*@ 1494 DMPlexCreateFromDAG - This takes as input the adjacency-list representation of the Directed Acyclic Graph (Hasse Diagram) encoding a mesh, and produces a DM 1495 1496 Input Parameters: 1497 + dm - The empty DM object, usually from DMCreate() and DMSetDimension() 1498 . depth - The depth of the DAG 1499 . numPoints - The number of points at each depth 1500 . coneSize - The cone size of each point 1501 . cones - The concatenation of the cone points for each point, the cone list must be oriented correctly for each point 1502 . coneOrientations - The orientation of each cone point 1503 - vertexCoords - An array of numVertices*dim numbers, the coordinates of each vertex 1504 1505 Output Parameter: 1506 . dm - The DM 1507 1508 Note: Two triangles sharing a face would have input 1509 $ depth = 1, numPoints = [4 2], coneSize = [3 3 0 0 0 0] 1510 $ cones = [2 3 4 3 5 4], coneOrientations = [0 0 0 0 0 0] 1511 $ vertexCoords = [-1.0 0.0 0.0 -1.0 0.0 1.0 1.0 0.0] 1512 $ 1513 which would result in the DMPlex 1514 $ 1515 $ 4 1516 $ / | \ 1517 $ / | \ 1518 $ / | \ 1519 $ 2 0 | 1 5 1520 $ \ | / 1521 $ \ | / 1522 $ \ | / 1523 $ 3 1524 $ 1525 $ Notice that all points are numbered consecutively, unlikely DMPlexCreateFromCellList() 1526 1527 Level: advanced 1528 1529 .seealso: DMPlexCreateFromCellList(), DMPlexCreate() 1530 @*/ 1531 PetscErrorCode DMPlexCreateFromDAG(DM dm, PetscInt depth, const PetscInt numPoints[], const PetscInt coneSize[], const PetscInt cones[], const PetscInt coneOrientations[], const PetscScalar vertexCoords[]) 1532 { 1533 Vec coordinates; 1534 PetscSection coordSection; 1535 PetscScalar *coords; 1536 PetscInt coordSize, firstVertex = -1, pStart = 0, pEnd = 0, p, v, dim, dimEmbed, d, off; 1537 PetscErrorCode ierr; 1538 1539 PetscFunctionBegin; 1540 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 1541 ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr); 1542 if (dimEmbed < dim) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Embedding dimension %d cannot be less than intrinsic dimension %d",dimEmbed,dim); 1543 for (d = 0; d <= depth; ++d) pEnd += numPoints[d]; 1544 ierr = DMPlexSetChart(dm, pStart, pEnd);CHKERRQ(ierr); 1545 for (p = pStart; p < pEnd; ++p) { 1546 ierr = DMPlexSetConeSize(dm, p, coneSize[p-pStart]);CHKERRQ(ierr); 1547 if (firstVertex < 0 && !coneSize[p - pStart]) { 1548 firstVertex = p - pStart; 1549 } 1550 } 1551 if (firstVertex < 0 && numPoints[0]) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Expected %d vertices but could not find any", numPoints[0]); 1552 ierr = DMSetUp(dm);CHKERRQ(ierr); /* Allocate space for cones */ 1553 for (p = pStart, off = 0; p < pEnd; off += coneSize[p-pStart], ++p) { 1554 ierr = DMPlexSetCone(dm, p, &cones[off]);CHKERRQ(ierr); 1555 ierr = DMPlexSetConeOrientation(dm, p, &coneOrientations[off]);CHKERRQ(ierr); 1556 } 1557 ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr); 1558 ierr = DMPlexStratify(dm);CHKERRQ(ierr); 1559 /* Build coordinates */ 1560 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 1561 ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr); 1562 ierr = PetscSectionSetFieldComponents(coordSection, 0, dimEmbed);CHKERRQ(ierr); 1563 ierr = PetscSectionSetChart(coordSection, firstVertex, firstVertex+numPoints[0]);CHKERRQ(ierr); 1564 for (v = firstVertex; v < firstVertex+numPoints[0]; ++v) { 1565 ierr = PetscSectionSetDof(coordSection, v, dimEmbed);CHKERRQ(ierr); 1566 ierr = PetscSectionSetFieldDof(coordSection, v, 0, dimEmbed);CHKERRQ(ierr); 1567 } 1568 ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 1569 ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 1570 ierr = VecCreate(PetscObjectComm((PetscObject)dm), &coordinates);CHKERRQ(ierr); 1571 ierr = VecSetBlockSize(coordinates, dimEmbed);CHKERRQ(ierr); 1572 ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr); 1573 ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 1574 ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr); 1575 ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 1576 for (v = 0; v < numPoints[0]; ++v) { 1577 PetscInt off; 1578 1579 ierr = PetscSectionGetOffset(coordSection, v+firstVertex, &off);CHKERRQ(ierr); 1580 for (d = 0; d < dimEmbed; ++d) { 1581 coords[off+d] = vertexCoords[v*dimEmbed+d]; 1582 } 1583 } 1584 ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 1585 ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr); 1586 ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 1587 PetscFunctionReturn(0); 1588 } 1589 1590 #undef __FUNCT__ 1591 #define __FUNCT__ "DMPlexCreateFromFile" 1592 /*@C 1593 DMPlexCreateFromFile - This takes a filename and produces a DM 1594 1595 Input Parameters: 1596 + comm - The communicator 1597 . filename - A file name 1598 - interpolate - Flag to create intermediate mesh pieces (edges, faces) 1599 1600 Output Parameter: 1601 . dm - The DM 1602 1603 Level: beginner 1604 1605 .seealso: DMPlexCreateFromDAG(), DMPlexCreateFromCellList(), DMPlexCreate() 1606 @*/ 1607 PetscErrorCode DMPlexCreateFromFile(MPI_Comm comm, const char filename[], PetscBool interpolate, DM *dm) 1608 { 1609 const char *extGmsh = ".msh"; 1610 const char *extCGNS = ".cgns"; 1611 const char *extExodus = ".exo"; 1612 size_t len; 1613 PetscBool isGmsh, isCGNS, isExodus; 1614 PetscMPIInt rank; 1615 PetscErrorCode ierr; 1616 1617 PetscFunctionBegin; 1618 PetscValidPointer(filename, 2); 1619 PetscValidPointer(dm, 4); 1620 ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 1621 ierr = PetscStrlen(filename, &len);CHKERRQ(ierr); 1622 if (!len) SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Filename must be a valid path"); 1623 ierr = PetscStrncmp(&filename[PetscMax(0,len-4)], extGmsh, 4, &isGmsh);CHKERRQ(ierr); 1624 ierr = PetscStrncmp(&filename[PetscMax(0,len-5)], extCGNS, 5, &isCGNS);CHKERRQ(ierr); 1625 ierr = PetscStrncmp(&filename[PetscMax(0,len-4)], extExodus, 4, &isExodus);CHKERRQ(ierr); 1626 if (isGmsh) { 1627 PetscViewer viewer; 1628 1629 ierr = PetscViewerCreate(comm, &viewer);CHKERRQ(ierr); 1630 ierr = PetscViewerSetType(viewer, PETSCVIEWERASCII);CHKERRQ(ierr); 1631 ierr = PetscViewerFileSetMode(viewer, FILE_MODE_READ);CHKERRQ(ierr); 1632 ierr = PetscViewerFileSetName(viewer, filename);CHKERRQ(ierr); 1633 ierr = DMPlexCreateGmsh(comm, viewer, interpolate, dm);CHKERRQ(ierr); 1634 ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 1635 } else if (isCGNS) { 1636 ierr = DMPlexCreateCGNSFromFile(comm, filename, interpolate, dm);CHKERRQ(ierr); 1637 } else if (isExodus) { 1638 ierr = DMPlexCreateExodusFromFile(comm, filename, interpolate, dm);CHKERRQ(ierr); 1639 } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cannot load file %s: unrecognized extension", filename); 1640 PetscFunctionReturn(0); 1641 } 1642 1643 #undef __FUNCT__ 1644 #define __FUNCT__ "DMPlexCreateReferenceCell" 1645 /*@ 1646 DMPlexCreateReferenceCell - Create a DMPLEX with the appropriate FEM reference cell 1647 1648 Collective on comm 1649 1650 Input Parameters: 1651 + comm - The communicator 1652 . dim - The spatial dimension 1653 - simplex - Flag for simplex, otherwise use a tensor-product cell 1654 1655 Output Parameter: 1656 . refdm - The reference cell 1657 1658 Level: intermediate 1659 1660 .keywords: reference cell 1661 .seealso: 1662 @*/ 1663 PetscErrorCode DMPlexCreateReferenceCell(MPI_Comm comm, PetscInt dim, PetscBool simplex, DM *refdm) 1664 { 1665 DM rdm; 1666 PetscErrorCode ierr; 1667 1668 PetscFunctionBeginUser; 1669 ierr = DMCreate(comm, &rdm);CHKERRQ(ierr); 1670 ierr = DMSetType(rdm, DMPLEX);CHKERRQ(ierr); 1671 ierr = DMSetDimension(rdm, dim);CHKERRQ(ierr); 1672 switch (dim) { 1673 case 0: 1674 { 1675 PetscInt numPoints[1] = {1}; 1676 PetscInt coneSize[1] = {0}; 1677 PetscInt cones[1] = {0}; 1678 PetscInt coneOrientations[1] = {0}; 1679 PetscScalar vertexCoords[1] = {0.0}; 1680 1681 ierr = DMPlexCreateFromDAG(rdm, 0, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr); 1682 } 1683 break; 1684 case 1: 1685 { 1686 PetscInt numPoints[2] = {2, 1}; 1687 PetscInt coneSize[3] = {2, 0, 0}; 1688 PetscInt cones[2] = {1, 2}; 1689 PetscInt coneOrientations[2] = {0, 0}; 1690 PetscScalar vertexCoords[2] = {-1.0, 1.0}; 1691 1692 ierr = DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr); 1693 } 1694 break; 1695 case 2: 1696 if (simplex) { 1697 PetscInt numPoints[2] = {3, 1}; 1698 PetscInt coneSize[4] = {3, 0, 0, 0}; 1699 PetscInt cones[3] = {1, 2, 3}; 1700 PetscInt coneOrientations[3] = {0, 0, 0}; 1701 PetscScalar vertexCoords[6] = {-1.0, -1.0, 1.0, -1.0, -1.0, 1.0}; 1702 1703 ierr = DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr); 1704 } else { 1705 PetscInt numPoints[2] = {4, 1}; 1706 PetscInt coneSize[5] = {4, 0, 0, 0, 0}; 1707 PetscInt cones[4] = {1, 2, 3, 4}; 1708 PetscInt coneOrientations[4] = {0, 0, 0, 0}; 1709 PetscScalar vertexCoords[8] = {-1.0, -1.0, 1.0, -1.0, 1.0, 1.0, -1.0, 1.0}; 1710 1711 ierr = DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr); 1712 } 1713 break; 1714 case 3: 1715 if (simplex) { 1716 PetscInt numPoints[2] = {4, 1}; 1717 PetscInt coneSize[5] = {4, 0, 0, 0, 0}; 1718 PetscInt cones[4] = {1, 3, 2, 4}; 1719 PetscInt coneOrientations[4] = {0, 0, 0, 0}; 1720 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}; 1721 1722 ierr = DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr); 1723 } else { 1724 PetscInt numPoints[2] = {8, 1}; 1725 PetscInt coneSize[9] = {8, 0, 0, 0, 0, 0, 0, 0, 0}; 1726 PetscInt cones[8] = {1, 4, 3, 2, 5, 6, 7, 8}; 1727 PetscInt coneOrientations[8] = {0, 0, 0, 0, 0, 0, 0, 0}; 1728 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, 1729 -1.0, -1.0, 1.0, 1.0, -1.0, 1.0, 1.0, 1.0, 1.0, -1.0, 1.0, 1.0}; 1730 1731 ierr = DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr); 1732 } 1733 break; 1734 default: 1735 SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Cannot create reference cell for dimension %d", dim); 1736 } 1737 *refdm = NULL; 1738 ierr = DMPlexInterpolate(rdm, refdm);CHKERRQ(ierr); 1739 ierr = DMPlexCopyCoordinates(rdm, *refdm);CHKERRQ(ierr); 1740 ierr = DMDestroy(&rdm);CHKERRQ(ierr); 1741 PetscFunctionReturn(0); 1742 } 1743