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 = DMSetLabelValue(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 = DMSetLabelValue(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 = DMCopyLabels(*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 const PetscInt numVertices = (edges[0]+1)*(edges[1]+1); 188 const 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)->options,((PetscObject) dm)->prefix, "-dm_plex_separate_marker", &markerSeparate, NULL);CHKERRQ(ierr); 204 if (markerSeparate) { 205 markerTop = 3; 206 markerBottom = 1; 207 markerRight = 2; 208 markerLeft = 4; 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 = DMSetLabelValue(dm, "marker", edge, markerRight);CHKERRQ(ierr); 229 ierr = DMSetLabelValue(dm, "marker", cone[0], markerRight);CHKERRQ(ierr); 230 if (ey == edges[1]-1) { 231 ierr = DMSetLabelValue(dm, "marker", cone[1], markerRight);CHKERRQ(ierr); 232 } 233 } else if (vx == 0) { 234 ierr = DMSetLabelValue(dm, "marker", edge, markerLeft);CHKERRQ(ierr); 235 ierr = DMSetLabelValue(dm, "marker", cone[0], markerLeft);CHKERRQ(ierr); 236 if (ey == edges[1]-1) { 237 ierr = DMSetLabelValue(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 = DMSetLabelValue(dm, "marker", edge, markerTop);CHKERRQ(ierr); 252 ierr = DMSetLabelValue(dm, "marker", cone[0], markerTop);CHKERRQ(ierr); 253 if (ex == edges[0]-1) { 254 ierr = DMSetLabelValue(dm, "marker", cone[1], markerTop);CHKERRQ(ierr); 255 } 256 } else if (vy == 0) { 257 ierr = DMSetLabelValue(dm, "marker", edge, markerBottom);CHKERRQ(ierr); 258 ierr = DMSetLabelValue(dm, "marker", cone[0], markerBottom);CHKERRQ(ierr); 259 if (ex == edges[0]-1) { 260 ierr = DMSetLabelValue(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 = DMSetLabelValue(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, faceMarkerTop = 1; 444 PetscInt markerBottom = 1, faceMarkerBottom = 1; 445 PetscInt markerFront = 1, faceMarkerFront = 1; 446 PetscInt markerBack = 1, faceMarkerBack = 1; 447 PetscInt markerRight = 1, faceMarkerRight = 1; 448 PetscInt markerLeft = 1, faceMarkerLeft = 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 switch (dim) { 458 case 2: 459 faceMarkerTop = 3; 460 faceMarkerBottom = 1; 461 faceMarkerRight = 2; 462 faceMarkerLeft = 4; 463 break; 464 case 3: 465 faceMarkerBottom = 1; 466 faceMarkerTop = 2; 467 faceMarkerFront = 3; 468 faceMarkerBack = 4; 469 faceMarkerRight = 5; 470 faceMarkerLeft = 6; 471 break; 472 default: 473 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Dimension %d not supported",dim); 474 break; 475 } 476 ierr = PetscOptionsGetBool(((PetscObject) dm)->options,((PetscObject) dm)->prefix, "-dm_plex_separate_marker", &markerSeparate, NULL);CHKERRQ(ierr); 477 if (markerSeparate) { 478 markerBottom = faceMarkerBottom; 479 markerTop = faceMarkerTop; 480 markerFront = faceMarkerFront; 481 markerBack = faceMarkerBack; 482 markerRight = faceMarkerRight; 483 markerLeft = faceMarkerLeft; 484 } 485 { 486 const PetscInt numXEdges = !rank ? edges[0] : 0; 487 const PetscInt numYEdges = !rank ? edges[1] : 0; 488 const PetscInt numZEdges = !rank ? edges[2] : 0; 489 const PetscInt numXVertices = !rank ? (bdX == DM_BOUNDARY_PERIODIC || bdX == DM_BOUNDARY_TWIST ? edges[0] : edges[0]+1) : 0; 490 const PetscInt numYVertices = !rank ? (bdY == DM_BOUNDARY_PERIODIC || bdY == DM_BOUNDARY_TWIST ? edges[1] : edges[1]+1) : 0; 491 const PetscInt numZVertices = !rank ? (bdZ == DM_BOUNDARY_PERIODIC || bdY == DM_BOUNDARY_TWIST ? edges[2] : edges[2]+1) : 0; 492 const PetscInt numCells = numXEdges*numYEdges*numZEdges; 493 const PetscInt numXFaces = numYEdges*numZEdges; 494 const PetscInt numYFaces = numXEdges*numZEdges; 495 const PetscInt numZFaces = numXEdges*numYEdges; 496 const PetscInt numTotXFaces = numXVertices*numXFaces; 497 const PetscInt numTotYFaces = numYVertices*numYFaces; 498 const PetscInt numTotZFaces = numZVertices*numZFaces; 499 const PetscInt numFaces = numTotXFaces + numTotYFaces + numTotZFaces; 500 const PetscInt numTotXEdges = numXEdges*numYVertices*numZVertices; 501 const PetscInt numTotYEdges = numYEdges*numXVertices*numZVertices; 502 const PetscInt numTotZEdges = numZEdges*numXVertices*numYVertices; 503 const PetscInt numVertices = numXVertices*numYVertices*numZVertices; 504 const PetscInt numEdges = numTotXEdges + numTotYEdges + numTotZEdges; 505 const PetscInt firstVertex = (dim == 2) ? numFaces : numCells; 506 const PetscInt firstXFace = (dim == 2) ? 0 : numCells + numVertices; 507 const PetscInt firstYFace = firstXFace + numTotXFaces; 508 const PetscInt firstZFace = firstYFace + numTotYFaces; 509 const PetscInt firstXEdge = numCells + numFaces + numVertices; 510 const PetscInt firstYEdge = firstXEdge + numTotXEdges; 511 const PetscInt firstZEdge = firstYEdge + numTotYEdges; 512 Vec coordinates; 513 PetscSection coordSection; 514 PetscScalar *coords; 515 PetscInt coordSize; 516 PetscInt v, vx, vy, vz; 517 PetscInt c, f, fx, fy, fz, e, ex, ey, ez; 518 519 ierr = DMPlexSetChart(dm, 0, numCells+numFaces+numEdges+numVertices);CHKERRQ(ierr); 520 for (c = 0; c < numCells; c++) { 521 ierr = DMPlexSetConeSize(dm, c, 6);CHKERRQ(ierr); 522 } 523 for (f = firstXFace; f < firstXFace+numFaces; ++f) { 524 ierr = DMPlexSetConeSize(dm, f, 4);CHKERRQ(ierr); 525 } 526 for (e = firstXEdge; e < firstXEdge+numEdges; ++e) { 527 ierr = DMPlexSetConeSize(dm, e, 2);CHKERRQ(ierr); 528 } 529 ierr = DMSetUp(dm);CHKERRQ(ierr); /* Allocate space for cones */ 530 /* Build cells */ 531 for (fz = 0; fz < numZEdges; ++fz) { 532 for (fy = 0; fy < numYEdges; ++fy) { 533 for (fx = 0; fx < numXEdges; ++fx) { 534 PetscInt cell = (fz*numYEdges + fy)*numXEdges + fx; 535 PetscInt faceB = firstZFace + (fy*numXEdges+fx)*numZVertices + fz; 536 PetscInt faceT = firstZFace + (fy*numXEdges+fx)*numZVertices + ((fz+1)%numZVertices); 537 PetscInt faceF = firstYFace + (fz*numXEdges+fx)*numYVertices + fy; 538 PetscInt faceK = firstYFace + (fz*numXEdges+fx)*numYVertices + ((fy+1)%numYVertices); 539 PetscInt faceL = firstXFace + (fz*numYEdges+fy)*numXVertices + fx; 540 PetscInt faceR = firstXFace + (fz*numYEdges+fy)*numXVertices + ((fx+1)%numXVertices); 541 /* B, T, F, K, R, L */ 542 PetscInt ornt[8] = {-4, 0, 0, -1, 0, -4}; /* ??? */ 543 PetscInt cone[8]; 544 545 /* no boundary twisting in 3D */ 546 cone[0] = faceB; cone[1] = faceT; cone[2] = faceF; cone[3] = faceK; cone[4] = faceR; cone[5] = faceL; 547 ierr = DMPlexSetCone(dm, cell, cone);CHKERRQ(ierr); 548 ierr = DMPlexSetConeOrientation(dm, cell, ornt);CHKERRQ(ierr); 549 } 550 } 551 } 552 /* Build x faces */ 553 for (fz = 0; fz < numZEdges; ++fz) { 554 for (fy = 0; fy < numYEdges; ++fy) { 555 for (fx = 0; fx < numXVertices; ++fx) { 556 PetscInt face = firstXFace + (fz*numYEdges+fy)*numXVertices + fx; 557 PetscInt edgeL = firstZEdge + ( fy* numXVertices+fx)*numZEdges + fz; 558 PetscInt edgeR = firstZEdge + (((fy+1)%numYVertices)*numXVertices+fx)*numZEdges + fz; 559 PetscInt edgeB = firstYEdge + ( fz* numXVertices+fx)*numYEdges + fy; 560 PetscInt edgeT = firstYEdge + (((fz+1)%numZVertices)*numXVertices+fx)*numYEdges + fy; 561 PetscInt ornt[4] = {0, 0, -2, -2}; 562 PetscInt cone[4]; 563 564 if (dim == 3) { 565 /* markers */ 566 if (bdX != DM_BOUNDARY_PERIODIC) { 567 if (fx == numXVertices-1) { 568 ierr = DMSetLabelValue(dm, "Face Sets", face, faceMarkerRight);CHKERRQ(ierr); 569 ierr = DMSetLabelValue(dm, "marker", face, markerRight);CHKERRQ(ierr); 570 } 571 else if (fx == 0) { 572 ierr = DMSetLabelValue(dm, "Face Sets", face, faceMarkerLeft);CHKERRQ(ierr); 573 ierr = DMSetLabelValue(dm, "marker", face, markerLeft);CHKERRQ(ierr); 574 } 575 } 576 } 577 cone[0] = edgeB; cone[1] = edgeR; cone[2] = edgeT; cone[3] = edgeL; 578 ierr = DMPlexSetCone(dm, face, cone);CHKERRQ(ierr); 579 ierr = DMPlexSetConeOrientation(dm, face, ornt);CHKERRQ(ierr); 580 } 581 } 582 } 583 /* Build y faces */ 584 for (fz = 0; fz < numZEdges; ++fz) { 585 for (fx = 0; fx < numYEdges; ++fx) { 586 for (fy = 0; fy < numYVertices; ++fy) { 587 PetscInt face = firstYFace + (fz*numXEdges+fx)*numYVertices + fy; 588 PetscInt edgeL = firstZEdge + (fy*numXVertices+ fx )*numZEdges + fz; 589 PetscInt edgeR = firstZEdge + (fy*numXVertices+((fx+1)%numXVertices))*numZEdges + fz; 590 PetscInt edgeB = firstXEdge + ( fz *numYVertices+fy)*numXEdges + fx; 591 PetscInt edgeT = firstXEdge + (((fz+1)%numZVertices)*numYVertices+fy)*numXEdges + fx; 592 PetscInt ornt[4] = {0, 0, -2, -2}; 593 PetscInt cone[4]; 594 595 if (dim == 3) { 596 /* markers */ 597 if (bdY != DM_BOUNDARY_PERIODIC) { 598 if (fy == numYVertices-1) { 599 ierr = DMSetLabelValue(dm, "Face Sets", face, faceMarkerBack);CHKERRQ(ierr); 600 ierr = DMSetLabelValue(dm, "marker", face, markerBack);CHKERRQ(ierr); 601 } 602 else if (fy == 0) { 603 ierr = DMSetLabelValue(dm, "Face Sets", face, faceMarkerFront);CHKERRQ(ierr); 604 ierr = DMSetLabelValue(dm, "marker", face, markerFront);CHKERRQ(ierr); 605 } 606 } 607 } 608 cone[0] = edgeB; cone[1] = edgeR; cone[2] = edgeT; cone[3] = edgeL; 609 ierr = DMPlexSetCone(dm, face, cone);CHKERRQ(ierr); 610 ierr = DMPlexSetConeOrientation(dm, face, ornt);CHKERRQ(ierr); 611 } 612 } 613 } 614 /* Build z faces */ 615 for (fy = 0; fy < numYEdges; ++fy) { 616 for (fx = 0; fx < numXEdges; ++fx) { 617 for (fz = 0; fz < numZVertices; fz++) { 618 PetscInt face = firstZFace + (fy*numXEdges+fx)*numZVertices + fz; 619 PetscInt edgeL = firstYEdge + (fz*numXVertices+ fx )*numYEdges + fy; 620 PetscInt edgeR = firstYEdge + (fz*numXVertices+((fx+1)%numXVertices))*numYEdges + fy; 621 PetscInt edgeB = firstXEdge + (fz*numYVertices+ fy )*numXEdges + fx; 622 PetscInt edgeT = firstXEdge + (fz*numYVertices+((fy+1)%numYVertices))*numXEdges + fx; 623 PetscInt ornt[4] = {0, 0, -2, -2}; 624 PetscInt cone[4]; 625 626 if (dim == 2) { 627 if (bdX == DM_BOUNDARY_TWIST && fx == numXEdges-1) {edgeR += numYEdges-1-2*fy; ornt[1] = -2;} 628 if (bdY == DM_BOUNDARY_TWIST && fy == numYEdges-1) {edgeT += numXEdges-1-2*fx; ornt[2] = 0;} 629 } 630 else { 631 /* markers */ 632 if (bdZ != DM_BOUNDARY_PERIODIC) { 633 if (fz == numZVertices-1) { 634 ierr = DMSetLabelValue(dm, "Face Sets", face, faceMarkerTop);CHKERRQ(ierr); 635 ierr = DMSetLabelValue(dm, "marker", face, markerTop);CHKERRQ(ierr); 636 } 637 else if (fz == 0) { 638 ierr = DMSetLabelValue(dm, "Face Sets", face, faceMarkerBottom);CHKERRQ(ierr); 639 ierr = DMSetLabelValue(dm, "marker", face, markerBottom);CHKERRQ(ierr); 640 } 641 } 642 } 643 cone[0] = edgeB; cone[1] = edgeR; cone[2] = edgeT; cone[3] = edgeL; 644 ierr = DMPlexSetCone(dm, face, cone);CHKERRQ(ierr); 645 ierr = DMPlexSetConeOrientation(dm, face, ornt);CHKERRQ(ierr); 646 } 647 } 648 } 649 /* Build Z edges*/ 650 for (vy = 0; vy < numYVertices; vy++) { 651 for (vx = 0; vx < numXVertices; vx++) { 652 for (ez = 0; ez < numZEdges; ez++) { 653 const PetscInt edge = firstZEdge + (vy*numXVertices+vx)*numZEdges + ez; 654 const PetscInt vertexB = firstVertex + ( ez *numYVertices+vy)*numXVertices + vx; 655 const PetscInt vertexT = firstVertex + (((ez+1)%numZVertices)*numYVertices+vy)*numXVertices + vx; 656 PetscInt cone[2]; 657 658 if (dim == 3) { 659 if (bdX != DM_BOUNDARY_PERIODIC) { 660 if (vx == numXVertices-1) { 661 ierr = DMSetLabelValue(dm, "marker", edge, markerRight);CHKERRQ(ierr); 662 } 663 else if (vx == 0) { 664 ierr = DMSetLabelValue(dm, "marker", edge, markerLeft);CHKERRQ(ierr); 665 } 666 } 667 if (bdY != DM_BOUNDARY_PERIODIC) { 668 if (vy == numYVertices-1) { 669 ierr = DMSetLabelValue(dm, "marker", edge, markerBack);CHKERRQ(ierr); 670 } 671 else if (vy == 0) { 672 ierr = DMSetLabelValue(dm, "marker", edge, markerFront);CHKERRQ(ierr); 673 } 674 } 675 } 676 cone[0] = vertexB; cone[1] = vertexT; 677 ierr = DMPlexSetCone(dm, edge, cone);CHKERRQ(ierr); 678 } 679 } 680 } 681 /* Build Y edges*/ 682 for (vz = 0; vz < numZVertices; vz++) { 683 for (vx = 0; vx < numXVertices; vx++) { 684 for (ey = 0; ey < numYEdges; ey++) { 685 const PetscInt nextv = (dim == 2 && bdY == DM_BOUNDARY_TWIST && ey == numYEdges-1) ? (numXVertices-vx-1) : (vz*numYVertices+((ey+1)%numYVertices))*numXVertices + vx; 686 const PetscInt edge = firstYEdge + (vz*numXVertices+vx)*numYEdges + ey; 687 const PetscInt vertexF = firstVertex + (vz*numYVertices+ey)*numXVertices + vx; 688 const PetscInt vertexK = firstVertex + nextv; 689 PetscInt cone[2]; 690 691 cone[0] = vertexF; cone[1] = vertexK; 692 ierr = DMPlexSetCone(dm, edge, cone);CHKERRQ(ierr); 693 if (dim == 2) { 694 if ((bdX != DM_BOUNDARY_PERIODIC) && (bdX != DM_BOUNDARY_TWIST)) { 695 if (vx == numXVertices-1) { 696 ierr = DMSetLabelValue(dm, "Face Sets", edge, faceMarkerRight);CHKERRQ(ierr); 697 ierr = DMSetLabelValue(dm, "marker", edge, markerRight);CHKERRQ(ierr); 698 ierr = DMSetLabelValue(dm, "marker", cone[0], markerRight);CHKERRQ(ierr); 699 if (ey == numYEdges-1) { 700 ierr = DMSetLabelValue(dm, "marker", cone[1], markerRight);CHKERRQ(ierr); 701 } 702 } 703 else if (vx == 0) { 704 ierr = DMSetLabelValue(dm, "Face Sets", edge, faceMarkerLeft);CHKERRQ(ierr); 705 ierr = DMSetLabelValue(dm, "marker", edge, markerLeft);CHKERRQ(ierr); 706 ierr = DMSetLabelValue(dm, "marker", cone[0], markerLeft);CHKERRQ(ierr); 707 if (ey == numYEdges-1) { 708 ierr = DMSetLabelValue(dm, "marker", cone[1], markerLeft);CHKERRQ(ierr); 709 } 710 } 711 } 712 } 713 else { 714 if (bdX != DM_BOUNDARY_PERIODIC) { 715 if (vx == numXVertices-1) { 716 ierr = DMSetLabelValue(dm, "marker", edge, markerRight);CHKERRQ(ierr); 717 } 718 else if (vx == 0) { 719 ierr = DMSetLabelValue(dm, "marker", edge, markerLeft);CHKERRQ(ierr); 720 } 721 } 722 if (bdZ != DM_BOUNDARY_PERIODIC) { 723 if (vz == numZVertices-1) { 724 ierr = DMSetLabelValue(dm, "marker", edge, markerTop);CHKERRQ(ierr); 725 } 726 else if (vz == 0) { 727 ierr = DMSetLabelValue(dm, "marker", edge, markerBottom);CHKERRQ(ierr); 728 } 729 } 730 } 731 } 732 } 733 } 734 /* Build X edges*/ 735 for (vz = 0; vz < numZVertices; vz++) { 736 for (vy = 0; vy < numYVertices; vy++) { 737 for (ex = 0; ex < numXEdges; ex++) { 738 const PetscInt nextv = (dim == 2 && bdX == DM_BOUNDARY_TWIST && ex == numXEdges-1) ? (numYVertices-vy-1)*numXVertices : (vz*numYVertices+vy)*numXVertices + (ex+1)%numXVertices; 739 const PetscInt edge = firstXEdge + (vz*numYVertices+vy)*numXEdges + ex; 740 const PetscInt vertexL = firstVertex + (vz*numYVertices+vy)*numXVertices + ex; 741 const PetscInt vertexR = firstVertex + nextv; 742 PetscInt cone[2]; 743 744 cone[0] = vertexL; cone[1] = vertexR; 745 ierr = DMPlexSetCone(dm, edge, cone);CHKERRQ(ierr); 746 if (dim == 2) { 747 if ((bdY != DM_BOUNDARY_PERIODIC) && (bdY != DM_BOUNDARY_TWIST)) { 748 if (vy == numYVertices-1) { 749 ierr = DMSetLabelValue(dm, "Face Sets", edge, faceMarkerTop);CHKERRQ(ierr); 750 ierr = DMSetLabelValue(dm, "marker", edge, markerTop);CHKERRQ(ierr); 751 ierr = DMSetLabelValue(dm, "marker", cone[0], markerTop);CHKERRQ(ierr); 752 if (ex == numXEdges-1) { 753 ierr = DMSetLabelValue(dm, "marker", cone[1], markerTop);CHKERRQ(ierr); 754 } 755 } 756 else if (vy == 0) { 757 ierr = DMSetLabelValue(dm, "Face Sets", edge, faceMarkerBottom);CHKERRQ(ierr); 758 ierr = DMSetLabelValue(dm, "marker", edge, markerBottom);CHKERRQ(ierr); 759 ierr = DMSetLabelValue(dm, "marker", cone[0], markerBottom);CHKERRQ(ierr); 760 if (ex == numXEdges-1) { 761 ierr = DMSetLabelValue(dm, "marker", cone[1], markerBottom);CHKERRQ(ierr); 762 } 763 } 764 } 765 } 766 else { 767 if (bdY != DM_BOUNDARY_PERIODIC) { 768 if (vy == numYVertices-1) { 769 ierr = DMSetLabelValue(dm, "marker", edge, markerBack);CHKERRQ(ierr); 770 } 771 else if (vy == 0) { 772 ierr = DMSetLabelValue(dm, "marker", edge, markerFront);CHKERRQ(ierr); 773 } 774 } 775 if (bdZ != DM_BOUNDARY_PERIODIC) { 776 if (vz == numZVertices-1) { 777 ierr = DMSetLabelValue(dm, "marker", edge, markerTop);CHKERRQ(ierr); 778 } 779 else if (vz == 0) { 780 ierr = DMSetLabelValue(dm, "marker", edge, markerBottom);CHKERRQ(ierr); 781 } 782 } 783 } 784 } 785 } 786 } 787 ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr); 788 ierr = DMPlexStratify(dm);CHKERRQ(ierr); 789 /* Build coordinates */ 790 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 791 ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr); 792 ierr = PetscSectionSetFieldComponents(coordSection, 0, dim);CHKERRQ(ierr); 793 ierr = PetscSectionSetChart(coordSection, firstVertex, firstVertex+numVertices);CHKERRQ(ierr); 794 for (v = firstVertex; v < firstVertex+numVertices; ++v) { 795 ierr = PetscSectionSetDof(coordSection, v, dim);CHKERRQ(ierr); 796 ierr = PetscSectionSetFieldDof(coordSection, v, 0, dim);CHKERRQ(ierr); 797 } 798 ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 799 ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 800 ierr = VecCreate(PetscObjectComm((PetscObject)dm), &coordinates);CHKERRQ(ierr); 801 ierr = VecSetBlockSize(coordinates, dim);CHKERRQ(ierr); 802 ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr); 803 ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 804 ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr); 805 ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 806 for (vz = 0; vz < numZVertices; ++vz) { 807 for (vy = 0; vy < numYVertices; ++vy) { 808 for (vx = 0; vx < numXVertices; ++vx) { 809 coords[((vz*numYVertices+vy)*numXVertices+vx)*dim+0] = lower[0] + ((upper[0] - lower[0])/numXEdges)*vx; 810 coords[((vz*numYVertices+vy)*numXVertices+vx)*dim+1] = lower[1] + ((upper[1] - lower[1])/numYEdges)*vy; 811 if (dim == 3) { 812 coords[((vz*numYVertices+vy)*numXVertices+vx)*dim+2] = lower[2] + ((upper[2] - lower[2])/numZEdges)*vz; 813 } 814 } 815 } 816 } 817 ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 818 ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr); 819 ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 820 } 821 PetscFunctionReturn(0); 822 } 823 824 #undef __FUNCT__ 825 #define __FUNCT__ "DMPlexCreateSquareMesh" 826 /*@ 827 DMPlexCreateSquareMesh - Creates a 2D mesh for a square lattice. 828 829 Collective on MPI_Comm 830 831 Input Parameters: 832 + comm - The communicator for the DM object 833 . lower - The lower left corner coordinates 834 . upper - The upper right corner coordinates 835 . edges - The number of cells in each direction 836 . bdX - The boundary type for the X direction 837 - bdY - The boundary type for the Y direction 838 839 Output Parameter: 840 . dm - The DM object 841 842 Note: Here is the numbering returned for 2 cells in each direction: 843 $ 22--8-23--9--24 844 $ | | | 845 $ 13 2 14 3 15 846 $ | | | 847 $ 19--6-20--7--21 848 $ | | | 849 $ 10 0 11 1 12 850 $ | | | 851 $ 16--4-17--5--18 852 853 Level: beginner 854 855 .keywords: DM, create 856 .seealso: DMPlexCreateBoxMesh(), DMPlexCreateSquareBoundary(), DMPlexCreateCubeBoundary(), DMSetType(), DMCreate() 857 @*/ 858 PetscErrorCode DMPlexCreateSquareMesh(DM dm, const PetscReal lower[], const PetscReal upper[], const PetscInt edges[], DMBoundaryType bdX, DMBoundaryType bdY) 859 { 860 PetscReal lower3[3], upper3[3]; 861 PetscInt edges3[3]; 862 PetscErrorCode ierr; 863 864 PetscFunctionBegin; 865 lower3[0] = lower[0]; lower3[1] = lower[1]; lower3[2] = 0.; 866 upper3[0] = upper[0]; upper3[1] = upper[1]; upper3[2] = 0.; 867 edges3[0] = edges[0]; edges3[1] = edges[1]; edges3[2] = 0; 868 ierr = DMPlexCreateCubeMesh_Internal(dm, lower3, upper3, edges3, bdX, bdY, DM_BOUNDARY_NONE);CHKERRQ(ierr); 869 PetscFunctionReturn(0); 870 } 871 872 #undef __FUNCT__ 873 #define __FUNCT__ "DMPlexCreateBoxMesh" 874 /*@ 875 DMPlexCreateBoxMesh - Creates a mesh on the tensor product of unit intervals (box) using simplices. 876 877 Collective on MPI_Comm 878 879 Input Parameters: 880 + comm - The communicator for the DM object 881 . dim - The spatial dimension 882 . numFaces - Number of faces per dimension 883 - interpolate - Flag to create intermediate mesh pieces (edges, faces) 884 885 Output Parameter: 886 . dm - The DM object 887 888 Level: beginner 889 890 .keywords: DM, create 891 .seealso: DMPlexCreateHexBoxMesh(), DMSetType(), DMCreate() 892 @*/ 893 PetscErrorCode DMPlexCreateBoxMesh(MPI_Comm comm, PetscInt dim, PetscInt numFaces, PetscBool interpolate, DM *dm) 894 { 895 DM boundary; 896 PetscErrorCode ierr; 897 898 PetscFunctionBegin; 899 PetscValidPointer(dm, 4); 900 ierr = DMCreate(comm, &boundary);CHKERRQ(ierr); 901 PetscValidLogicalCollectiveInt(boundary,dim,2); 902 ierr = DMSetType(boundary, DMPLEX);CHKERRQ(ierr); 903 ierr = DMSetDimension(boundary, dim-1);CHKERRQ(ierr); 904 switch (dim) { 905 case 2: 906 { 907 PetscReal lower[2] = {0.0, 0.0}; 908 PetscReal upper[2] = {1.0, 1.0}; 909 PetscInt edges[2]; 910 911 edges[0] = numFaces; edges[1] = numFaces; 912 ierr = DMPlexCreateSquareBoundary(boundary, lower, upper, edges);CHKERRQ(ierr); 913 break; 914 } 915 case 3: 916 { 917 PetscReal lower[3] = {0.0, 0.0, 0.0}; 918 PetscReal upper[3] = {1.0, 1.0, 1.0}; 919 PetscInt faces[3]; 920 921 faces[0] = numFaces; faces[1] = numFaces; faces[2] = numFaces; 922 ierr = DMPlexCreateCubeBoundary(boundary, lower, upper, faces);CHKERRQ(ierr); 923 break; 924 } 925 default: 926 SETERRQ1(comm, PETSC_ERR_SUP, "Dimension not supported: %d", dim); 927 } 928 ierr = DMPlexGenerate(boundary, NULL, interpolate, dm);CHKERRQ(ierr); 929 ierr = DMDestroy(&boundary);CHKERRQ(ierr); 930 PetscFunctionReturn(0); 931 } 932 933 #undef __FUNCT__ 934 #define __FUNCT__ "DMPlexCreateHexBoxMesh" 935 /*@ 936 DMPlexCreateHexBoxMesh - Creates a mesh on the tensor product of unit intervals (box) using hexahedra. 937 938 Collective on MPI_Comm 939 940 Input Parameters: 941 + comm - The communicator for the DM object 942 . dim - The spatial dimension 943 . periodicX - The boundary type for the X direction 944 . periodicY - The boundary type for the Y direction 945 . periodicZ - The boundary type for the Z direction 946 - cells - The number of cells in each direction 947 948 Output Parameter: 949 . dm - The DM object 950 951 Level: beginner 952 953 .keywords: DM, create 954 .seealso: DMPlexCreateBoxMesh(), DMSetType(), DMCreate() 955 @*/ 956 PetscErrorCode DMPlexCreateHexBoxMesh(MPI_Comm comm, PetscInt dim, const PetscInt cells[], DMBoundaryType periodicX, DMBoundaryType periodicY, DMBoundaryType periodicZ, DM *dm) 957 { 958 PetscErrorCode ierr; 959 960 PetscFunctionBegin; 961 PetscValidPointer(dm, 4); 962 ierr = DMCreate(comm, dm);CHKERRQ(ierr); 963 PetscValidLogicalCollectiveInt(*dm,dim,2); 964 ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 965 ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr); 966 switch (dim) { 967 case 2: 968 { 969 PetscReal lower[2] = {0.0, 0.0}; 970 PetscReal upper[2] = {1.0, 1.0}; 971 972 ierr = DMPlexCreateSquareMesh(*dm, lower, upper, cells, periodicX, periodicY);CHKERRQ(ierr); 973 break; 974 } 975 case 3: 976 { 977 PetscReal lower[3] = {0.0, 0.0, 0.0}; 978 PetscReal upper[3] = {1.0, 1.0, 1.0}; 979 980 ierr = DMPlexCreateCubeMesh_Internal(*dm, lower, upper, cells, periodicX, periodicY, periodicZ);CHKERRQ(ierr); 981 break; 982 } 983 default: 984 SETERRQ1(comm, PETSC_ERR_SUP, "Dimension not supported: %d", dim); 985 } 986 PetscFunctionReturn(0); 987 } 988 989 /* External function declarations here */ 990 extern PetscErrorCode DMCreateInterpolation_Plex(DM dmCoarse, DM dmFine, Mat *interpolation, Vec *scaling); 991 extern PetscErrorCode DMCreateInjection_Plex(DM dmCoarse, DM dmFine, Mat *mat); 992 extern PetscErrorCode DMCreateDefaultSection_Plex(DM dm); 993 extern PetscErrorCode DMCreateDefaultConstraints_Plex(DM dm); 994 extern PetscErrorCode DMCreateMatrix_Plex(DM dm, Mat *J); 995 extern PetscErrorCode DMCreateCoordinateDM_Plex(DM dm, DM *cdm); 996 extern PetscErrorCode DMRefine_Plex(DM dm, MPI_Comm comm, DM *dmRefined); 997 extern PetscErrorCode DMCoarsen_Plex(DM dm, MPI_Comm comm, DM *dmCoarsened); 998 extern PetscErrorCode DMRefineHierarchy_Plex(DM dm, PetscInt nlevels, DM dmRefined[]); 999 extern PetscErrorCode DMCoarsenHierarchy_Plex(DM dm, PetscInt nlevels, DM dmCoarsened[]); 1000 extern PetscErrorCode DMClone_Plex(DM dm, DM *newdm); 1001 extern PetscErrorCode DMSetUp_Plex(DM dm); 1002 extern PetscErrorCode DMDestroy_Plex(DM dm); 1003 extern PetscErrorCode DMView_Plex(DM dm, PetscViewer viewer); 1004 extern PetscErrorCode DMLoad_Plex(DM dm, PetscViewer viewer); 1005 extern PetscErrorCode DMCreateSubDM_Plex(DM dm, PetscInt numFields, PetscInt fields[], IS *is, DM *subdm); 1006 extern PetscErrorCode DMLocatePoints_Plex(DM dm, Vec v, PetscSF cellSF); 1007 extern PetscErrorCode DMProjectFunctionLocal_Plex(DM,PetscReal,PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal[],PetscInt,PetscScalar *,void *),void **,InsertMode,Vec); 1008 extern PetscErrorCode DMProjectFunctionLabelLocal_Plex(DM,PetscReal,DMLabel,PetscInt,const PetscInt[],PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal[],PetscInt,PetscScalar *,void *),void **,InsertMode,Vec); 1009 extern PetscErrorCode DMProjectFieldLocal_Plex(DM,Vec,void (**)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,const PetscReal[],PetscScalar[]),InsertMode,Vec); 1010 extern PetscErrorCode DMComputeL2Diff_Plex(DM,PetscReal,PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal[],PetscInt,PetscScalar *,void *),void **,Vec,PetscReal *); 1011 extern PetscErrorCode DMComputeL2GradientDiff_Plex(DM,PetscReal,PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal[], const PetscReal[],PetscInt,PetscScalar *,void *),void **,Vec,const PetscReal [],PetscReal *); 1012 extern PetscErrorCode DMComputeL2FieldDiff_Plex(DM,PetscReal,PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal[],PetscInt,PetscScalar *,void *),void **,Vec,PetscReal *); 1013 1014 #undef __FUNCT__ 1015 #define __FUNCT__ "DMPlexReplace_Static" 1016 /* Replace dm with the contents of dmNew 1017 - Share the DM_Plex structure 1018 - Share the coordinates 1019 - Share the SF 1020 */ 1021 static PetscErrorCode DMPlexReplace_Static(DM dm, DM dmNew) 1022 { 1023 PetscSF sf; 1024 DM coordDM, coarseDM; 1025 Vec coords; 1026 const PetscReal *maxCell, *L; 1027 const DMBoundaryType *bd; 1028 PetscErrorCode ierr; 1029 1030 PetscFunctionBegin; 1031 ierr = DMGetPointSF(dmNew, &sf);CHKERRQ(ierr); 1032 ierr = DMSetPointSF(dm, sf);CHKERRQ(ierr); 1033 ierr = DMGetCoordinateDM(dmNew, &coordDM);CHKERRQ(ierr); 1034 ierr = DMGetCoordinatesLocal(dmNew, &coords);CHKERRQ(ierr); 1035 ierr = DMSetCoordinateDM(dm, coordDM);CHKERRQ(ierr); 1036 ierr = DMSetCoordinatesLocal(dm, coords);CHKERRQ(ierr); 1037 ierr = DMGetPeriodicity(dm, &maxCell, &L, &bd);CHKERRQ(ierr); 1038 if (L) {ierr = DMSetPeriodicity(dmNew, maxCell, L, bd);CHKERRQ(ierr);} 1039 ierr = DMDestroy_Plex(dm);CHKERRQ(ierr); 1040 dm->data = dmNew->data; 1041 ((DM_Plex *) dmNew->data)->refct++; 1042 dmNew->labels->refct++; 1043 if (!--(dm->labels->refct)) { 1044 DMLabelLink next = dm->labels->next; 1045 1046 /* destroy the labels */ 1047 while (next) { 1048 DMLabelLink tmp = next->next; 1049 1050 ierr = DMLabelDestroy(&next->label);CHKERRQ(ierr); 1051 ierr = PetscFree(next);CHKERRQ(ierr); 1052 next = tmp; 1053 } 1054 ierr = PetscFree(dm->labels);CHKERRQ(ierr); 1055 } 1056 dm->labels = dmNew->labels; 1057 dm->depthLabel = dmNew->depthLabel; 1058 ierr = DMGetCoarseDM(dmNew,&coarseDM);CHKERRQ(ierr); 1059 ierr = DMSetCoarseDM(dm,coarseDM);CHKERRQ(ierr); 1060 PetscFunctionReturn(0); 1061 } 1062 1063 #undef __FUNCT__ 1064 #define __FUNCT__ "DMPlexSwap_Static" 1065 /* Swap dm with the contents of dmNew 1066 - Swap the DM_Plex structure 1067 - Swap the coordinates 1068 - Swap the point PetscSF 1069 */ 1070 static PetscErrorCode DMPlexSwap_Static(DM dmA, DM dmB) 1071 { 1072 DM coordDMA, coordDMB; 1073 Vec coordsA, coordsB; 1074 PetscSF sfA, sfB; 1075 void *tmp; 1076 DMLabelLinkList listTmp; 1077 DMLabel depthTmp; 1078 PetscErrorCode ierr; 1079 1080 PetscFunctionBegin; 1081 ierr = DMGetPointSF(dmA, &sfA);CHKERRQ(ierr); 1082 ierr = DMGetPointSF(dmB, &sfB);CHKERRQ(ierr); 1083 ierr = PetscObjectReference((PetscObject) sfA);CHKERRQ(ierr); 1084 ierr = DMSetPointSF(dmA, sfB);CHKERRQ(ierr); 1085 ierr = DMSetPointSF(dmB, sfA);CHKERRQ(ierr); 1086 ierr = PetscObjectDereference((PetscObject) sfA);CHKERRQ(ierr); 1087 1088 ierr = DMGetCoordinateDM(dmA, &coordDMA);CHKERRQ(ierr); 1089 ierr = DMGetCoordinateDM(dmB, &coordDMB);CHKERRQ(ierr); 1090 ierr = PetscObjectReference((PetscObject) coordDMA);CHKERRQ(ierr); 1091 ierr = DMSetCoordinateDM(dmA, coordDMB);CHKERRQ(ierr); 1092 ierr = DMSetCoordinateDM(dmB, coordDMA);CHKERRQ(ierr); 1093 ierr = PetscObjectDereference((PetscObject) coordDMA);CHKERRQ(ierr); 1094 1095 ierr = DMGetCoordinatesLocal(dmA, &coordsA);CHKERRQ(ierr); 1096 ierr = DMGetCoordinatesLocal(dmB, &coordsB);CHKERRQ(ierr); 1097 ierr = PetscObjectReference((PetscObject) coordsA);CHKERRQ(ierr); 1098 ierr = DMSetCoordinatesLocal(dmA, coordsB);CHKERRQ(ierr); 1099 ierr = DMSetCoordinatesLocal(dmB, coordsA);CHKERRQ(ierr); 1100 ierr = PetscObjectDereference((PetscObject) coordsA);CHKERRQ(ierr); 1101 1102 tmp = dmA->data; 1103 dmA->data = dmB->data; 1104 dmB->data = tmp; 1105 listTmp = dmA->labels; 1106 dmA->labels = dmB->labels; 1107 dmB->labels = listTmp; 1108 depthTmp = dmA->depthLabel; 1109 dmA->depthLabel = dmB->depthLabel; 1110 dmB->depthLabel = depthTmp; 1111 PetscFunctionReturn(0); 1112 } 1113 1114 #undef __FUNCT__ 1115 #define __FUNCT__ "DMSetFromOptions_NonRefinement_Plex" 1116 PetscErrorCode DMSetFromOptions_NonRefinement_Plex(PetscOptionItems *PetscOptionsObject,DM dm) 1117 { 1118 DM_Plex *mesh = (DM_Plex*) dm->data; 1119 DMBoundary b; 1120 PetscErrorCode ierr; 1121 1122 PetscFunctionBegin; 1123 /* Handle boundary conditions */ 1124 ierr = PetscOptionsBegin(PetscObjectComm((PetscObject) dm), NULL, "Boundary condition options", "");CHKERRQ(ierr); 1125 for (b = dm->boundary->next; b; b = b->next) { 1126 char optname[1024]; 1127 PetscInt ids[1024], len = 1024, i; 1128 PetscBool flg; 1129 1130 ierr = PetscSNPrintf(optname, sizeof(optname), "-bc_%s", b->name);CHKERRQ(ierr); 1131 ierr = PetscMemzero(ids, sizeof(ids));CHKERRQ(ierr); 1132 ierr = PetscOptionsIntArray(optname, "List of boundary IDs", "", ids, &len, &flg);CHKERRQ(ierr); 1133 if (flg) { 1134 DMLabel label; 1135 1136 ierr = DMGetLabel(dm, b->labelname, &label);CHKERRQ(ierr); 1137 for (i = 0; i < len; ++i) { 1138 PetscBool has; 1139 1140 ierr = DMLabelHasValue(label, ids[i], &has);CHKERRQ(ierr); 1141 if (!has) SETERRQ2(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Boundary id %D is not present in the label %s", ids[i], b->name); 1142 } 1143 b->numids = len; 1144 ierr = PetscFree(b->ids);CHKERRQ(ierr); 1145 ierr = PetscMalloc1(len, &b->ids);CHKERRQ(ierr); 1146 ierr = PetscMemcpy(b->ids, ids, len*sizeof(PetscInt));CHKERRQ(ierr); 1147 } 1148 ierr = PetscSNPrintf(optname, sizeof(optname), "-bc_%s_comp", b->name);CHKERRQ(ierr); 1149 ierr = PetscMemzero(ids, sizeof(ids));CHKERRQ(ierr); 1150 ierr = PetscOptionsIntArray(optname, "List of boundary field components", "", ids, &len, &flg);CHKERRQ(ierr); 1151 if (flg) { 1152 b->numcomps = len; 1153 ierr = PetscFree(b->comps);CHKERRQ(ierr); 1154 ierr = PetscMalloc1(len, &b->comps);CHKERRQ(ierr); 1155 ierr = PetscMemcpy(b->comps, ids, len*sizeof(PetscInt));CHKERRQ(ierr); 1156 } 1157 } 1158 ierr = PetscOptionsEnd();CHKERRQ(ierr); 1159 /* Handle viewing */ 1160 ierr = PetscOptionsBool("-dm_plex_print_set_values", "Output all set values info", "DMView", PETSC_FALSE, &mesh->printSetValues, NULL);CHKERRQ(ierr); 1161 ierr = PetscOptionsInt("-dm_plex_print_fem", "Debug output level all fem computations", "DMView", 0, &mesh->printFEM, NULL);CHKERRQ(ierr); 1162 ierr = PetscOptionsReal("-dm_plex_print_tol", "Tolerance for FEM output", "DMView", mesh->printTol, &mesh->printTol, NULL);CHKERRQ(ierr); 1163 /* Point Location */ 1164 ierr = PetscOptionsBool("-dm_plex_hash_location", "Use grid hashing for point location", "DMView", PETSC_FALSE, &mesh->useHashLocation, NULL);CHKERRQ(ierr); 1165 /* Generation and remeshing */ 1166 ierr = PetscOptionsBool("-dm_plex_remesh_bd", "Allow changes to the boundary on remeshing", "DMView", PETSC_FALSE, &mesh->remeshBd, NULL);CHKERRQ(ierr); 1167 /* Projection behavior */ 1168 ierr = PetscOptionsInt("-dm_plex_max_projection_height", "Maxmimum mesh point height used to project locally", "DMPlexSetMaxProjectionHeight", 0, &mesh->maxProjectionHeight, NULL);CHKERRQ(ierr); 1169 ierr = PetscOptionsBool("-dm_plex_regular_refinement", "Use special nested projection algorithm for regular refinement", "DMPlexSetRegularRefinement", mesh->regularRefinement, &mesh->regularRefinement, NULL);CHKERRQ(ierr); 1170 PetscFunctionReturn(0); 1171 } 1172 1173 #undef __FUNCT__ 1174 #define __FUNCT__ "DMSetFromOptions_Plex" 1175 PetscErrorCode DMSetFromOptions_Plex(PetscOptionItems *PetscOptionsObject,DM dm) 1176 { 1177 PetscInt refine = 0, coarsen = 0, r; 1178 PetscBool isHierarchy; 1179 PetscErrorCode ierr; 1180 1181 PetscFunctionBegin; 1182 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1183 ierr = PetscOptionsHead(PetscOptionsObject,"DMPlex Options");CHKERRQ(ierr); 1184 /* Handle DMPlex refinement */ 1185 ierr = PetscOptionsInt("-dm_refine", "The number of uniform refinements", "DMCreate", refine, &refine, NULL);CHKERRQ(ierr); 1186 ierr = PetscOptionsInt("-dm_refine_hierarchy", "The number of uniform refinements", "DMCreate", refine, &refine, &isHierarchy);CHKERRQ(ierr); 1187 if (refine) {ierr = DMPlexSetRefinementUniform(dm, PETSC_TRUE);CHKERRQ(ierr);} 1188 if (refine && isHierarchy) { 1189 DM *dms, coarseDM; 1190 1191 ierr = DMGetCoarseDM(dm, &coarseDM);CHKERRQ(ierr); 1192 ierr = PetscObjectReference((PetscObject)coarseDM);CHKERRQ(ierr); 1193 ierr = PetscMalloc1(refine,&dms);CHKERRQ(ierr); 1194 ierr = DMRefineHierarchy(dm, refine, dms);CHKERRQ(ierr); 1195 /* Total hack since we do not pass in a pointer */ 1196 ierr = DMPlexSwap_Static(dm, dms[refine-1]);CHKERRQ(ierr); 1197 if (refine == 1) { 1198 ierr = DMSetCoarseDM(dm, dms[0]);CHKERRQ(ierr); 1199 ierr = DMPlexSetRegularRefinement(dm, PETSC_TRUE);CHKERRQ(ierr); 1200 } else { 1201 ierr = DMSetCoarseDM(dm, dms[refine-2]);CHKERRQ(ierr); 1202 ierr = DMPlexSetRegularRefinement(dm, PETSC_TRUE);CHKERRQ(ierr); 1203 ierr = DMSetCoarseDM(dms[0], dms[refine-1]);CHKERRQ(ierr); 1204 ierr = DMPlexSetRegularRefinement(dms[0], PETSC_TRUE);CHKERRQ(ierr); 1205 } 1206 ierr = DMSetCoarseDM(dms[refine-1], coarseDM);CHKERRQ(ierr); 1207 ierr = PetscObjectDereference((PetscObject)coarseDM);CHKERRQ(ierr); 1208 /* Free DMs */ 1209 for (r = 0; r < refine; ++r) { 1210 ierr = DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dms[r]);CHKERRQ(ierr); 1211 ierr = DMDestroy(&dms[r]);CHKERRQ(ierr); 1212 } 1213 ierr = PetscFree(dms);CHKERRQ(ierr); 1214 } else { 1215 for (r = 0; r < refine; ++r) { 1216 DM refinedMesh; 1217 1218 ierr = DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dm);CHKERRQ(ierr); 1219 ierr = DMRefine(dm, PetscObjectComm((PetscObject) dm), &refinedMesh);CHKERRQ(ierr); 1220 /* Total hack since we do not pass in a pointer */ 1221 ierr = DMPlexReplace_Static(dm, refinedMesh);CHKERRQ(ierr); 1222 ierr = DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dm);CHKERRQ(ierr); 1223 ierr = DMDestroy(&refinedMesh);CHKERRQ(ierr); 1224 } 1225 } 1226 /* Handle DMPlex coarsening */ 1227 ierr = PetscOptionsInt("-dm_coarsen", "Coarsen the mesh", "DMCreate", coarsen, &coarsen, NULL);CHKERRQ(ierr); 1228 ierr = PetscOptionsInt("-dm_coarsen_hierarchy", "The number of coarsenings", "DMCreate", coarsen, &coarsen, &isHierarchy);CHKERRQ(ierr); 1229 if (coarsen && isHierarchy) { 1230 DM *dms; 1231 1232 ierr = PetscMalloc1(coarsen, &dms);CHKERRQ(ierr); 1233 ierr = DMCoarsenHierarchy(dm, coarsen, dms);CHKERRQ(ierr); 1234 /* Free DMs */ 1235 for (r = 0; r < coarsen; ++r) { 1236 ierr = DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dms[r]);CHKERRQ(ierr); 1237 ierr = DMDestroy(&dms[r]);CHKERRQ(ierr); 1238 } 1239 ierr = PetscFree(dms);CHKERRQ(ierr); 1240 } else { 1241 for (r = 0; r < coarsen; ++r) { 1242 DM coarseMesh; 1243 1244 ierr = DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dm);CHKERRQ(ierr); 1245 ierr = DMCoarsen(dm, PetscObjectComm((PetscObject) dm), &coarseMesh);CHKERRQ(ierr); 1246 /* Total hack since we do not pass in a pointer */ 1247 ierr = DMPlexReplace_Static(dm, coarseMesh);CHKERRQ(ierr); 1248 ierr = DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dm);CHKERRQ(ierr); 1249 ierr = DMDestroy(&coarseMesh);CHKERRQ(ierr); 1250 } 1251 } 1252 /* Handle */ 1253 ierr = DMSetFromOptions_NonRefinement_Plex(PetscOptionsObject, dm);CHKERRQ(ierr); 1254 ierr = PetscOptionsTail();CHKERRQ(ierr); 1255 PetscFunctionReturn(0); 1256 } 1257 1258 #undef __FUNCT__ 1259 #define __FUNCT__ "DMCreateGlobalVector_Plex" 1260 static PetscErrorCode DMCreateGlobalVector_Plex(DM dm,Vec *vec) 1261 { 1262 PetscErrorCode ierr; 1263 1264 PetscFunctionBegin; 1265 ierr = DMCreateGlobalVector_Section_Private(dm,vec);CHKERRQ(ierr); 1266 /* ierr = VecSetOperation(*vec, VECOP_DUPLICATE, (void(*)(void)) VecDuplicate_MPI_DM);CHKERRQ(ierr); */ 1267 ierr = VecSetOperation(*vec, VECOP_VIEW, (void (*)(void)) VecView_Plex);CHKERRQ(ierr); 1268 ierr = VecSetOperation(*vec, VECOP_VIEWNATIVE, (void (*)(void)) VecView_Plex_Native);CHKERRQ(ierr); 1269 ierr = VecSetOperation(*vec, VECOP_LOAD, (void (*)(void)) VecLoad_Plex);CHKERRQ(ierr); 1270 ierr = VecSetOperation(*vec, VECOP_LOADNATIVE, (void (*)(void)) VecLoad_Plex_Native);CHKERRQ(ierr); 1271 PetscFunctionReturn(0); 1272 } 1273 1274 #undef __FUNCT__ 1275 #define __FUNCT__ "DMCreateLocalVector_Plex" 1276 static PetscErrorCode DMCreateLocalVector_Plex(DM dm,Vec *vec) 1277 { 1278 PetscErrorCode ierr; 1279 1280 PetscFunctionBegin; 1281 ierr = DMCreateLocalVector_Section_Private(dm,vec);CHKERRQ(ierr); 1282 ierr = VecSetOperation(*vec, VECOP_VIEW, (void (*)(void)) VecView_Plex_Local);CHKERRQ(ierr); 1283 ierr = VecSetOperation(*vec, VECOP_LOAD, (void (*)(void)) VecLoad_Plex_Local);CHKERRQ(ierr); 1284 PetscFunctionReturn(0); 1285 } 1286 1287 #undef __FUNCT__ 1288 #define __FUNCT__ "DMGetDimPoints_Plex" 1289 static PetscErrorCode DMGetDimPoints_Plex(DM dm, PetscInt dim, PetscInt *pStart, PetscInt *pEnd) 1290 { 1291 PetscInt depth, d; 1292 PetscErrorCode ierr; 1293 1294 PetscFunctionBegin; 1295 ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 1296 if (depth == 1) { 1297 ierr = DMGetDimension(dm, &d);CHKERRQ(ierr); 1298 if (dim == 0) {ierr = DMPlexGetDepthStratum(dm, dim, pStart, pEnd);CHKERRQ(ierr);} 1299 else if (dim == d) {ierr = DMPlexGetDepthStratum(dm, 1, pStart, pEnd);CHKERRQ(ierr);} 1300 else {*pStart = 0; *pEnd = 0;} 1301 } else { 1302 ierr = DMPlexGetDepthStratum(dm, dim, pStart, pEnd);CHKERRQ(ierr); 1303 } 1304 PetscFunctionReturn(0); 1305 } 1306 1307 #undef __FUNCT__ 1308 #define __FUNCT__ "DMInitialize_Plex" 1309 PetscErrorCode DMInitialize_Plex(DM dm) 1310 { 1311 PetscFunctionBegin; 1312 dm->ops->view = DMView_Plex; 1313 dm->ops->load = DMLoad_Plex; 1314 dm->ops->setfromoptions = DMSetFromOptions_Plex; 1315 dm->ops->clone = DMClone_Plex; 1316 dm->ops->setup = DMSetUp_Plex; 1317 dm->ops->createdefaultsection = DMCreateDefaultSection_Plex; 1318 dm->ops->createdefaultconstraints = DMCreateDefaultConstraints_Plex; 1319 dm->ops->createglobalvector = DMCreateGlobalVector_Plex; 1320 dm->ops->createlocalvector = DMCreateLocalVector_Plex; 1321 dm->ops->getlocaltoglobalmapping = NULL; 1322 dm->ops->createfieldis = NULL; 1323 dm->ops->createcoordinatedm = DMCreateCoordinateDM_Plex; 1324 dm->ops->getcoloring = NULL; 1325 dm->ops->creatematrix = DMCreateMatrix_Plex; 1326 dm->ops->createinterpolation = DMCreateInterpolation_Plex; 1327 dm->ops->getaggregates = NULL; 1328 dm->ops->getinjection = DMCreateInjection_Plex; 1329 dm->ops->refine = DMRefine_Plex; 1330 dm->ops->coarsen = DMCoarsen_Plex; 1331 dm->ops->refinehierarchy = DMRefineHierarchy_Plex; 1332 dm->ops->coarsenhierarchy = DMCoarsenHierarchy_Plex; 1333 dm->ops->globaltolocalbegin = NULL; 1334 dm->ops->globaltolocalend = NULL; 1335 dm->ops->localtoglobalbegin = NULL; 1336 dm->ops->localtoglobalend = NULL; 1337 dm->ops->destroy = DMDestroy_Plex; 1338 dm->ops->createsubdm = DMCreateSubDM_Plex; 1339 dm->ops->getdimpoints = DMGetDimPoints_Plex; 1340 dm->ops->locatepoints = DMLocatePoints_Plex; 1341 dm->ops->projectfunctionlocal = DMProjectFunctionLocal_Plex; 1342 dm->ops->projectfunctionlabellocal = DMProjectFunctionLabelLocal_Plex; 1343 dm->ops->projectfieldlocal = DMProjectFieldLocal_Plex; 1344 dm->ops->computel2diff = DMComputeL2Diff_Plex; 1345 dm->ops->computel2gradientdiff = DMComputeL2GradientDiff_Plex; 1346 dm->ops->computel2fielddiff = DMComputeL2FieldDiff_Plex; 1347 PetscFunctionReturn(0); 1348 } 1349 1350 #undef __FUNCT__ 1351 #define __FUNCT__ "DMClone_Plex" 1352 PetscErrorCode DMClone_Plex(DM dm, DM *newdm) 1353 { 1354 DM_Plex *mesh = (DM_Plex *) dm->data; 1355 PetscErrorCode ierr; 1356 1357 PetscFunctionBegin; 1358 mesh->refct++; 1359 (*newdm)->data = mesh; 1360 ierr = PetscObjectChangeTypeName((PetscObject) *newdm, DMPLEX);CHKERRQ(ierr); 1361 ierr = DMInitialize_Plex(*newdm);CHKERRQ(ierr); 1362 PetscFunctionReturn(0); 1363 } 1364 1365 /*MC 1366 DMPLEX = "plex" - A DM object that encapsulates an unstructured mesh, or CW Complex, which can be expressed using a Hasse Diagram. 1367 In the local representation, Vecs contain all unknowns in the interior and shared boundary. This is 1368 specified by a PetscSection object. Ownership in the global representation is determined by 1369 ownership of the underlying DMPlex points. This is specified by another PetscSection object. 1370 1371 Level: intermediate 1372 1373 .seealso: DMType, DMPlexCreate(), DMCreate(), DMSetType() 1374 M*/ 1375 1376 #undef __FUNCT__ 1377 #define __FUNCT__ "DMCreate_Plex" 1378 PETSC_EXTERN PetscErrorCode DMCreate_Plex(DM dm) 1379 { 1380 DM_Plex *mesh; 1381 PetscInt unit, d; 1382 PetscErrorCode ierr; 1383 1384 PetscFunctionBegin; 1385 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1386 ierr = PetscNewLog(dm,&mesh);CHKERRQ(ierr); 1387 dm->dim = 0; 1388 dm->data = mesh; 1389 1390 mesh->refct = 1; 1391 ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm), &mesh->coneSection);CHKERRQ(ierr); 1392 mesh->maxConeSize = 0; 1393 mesh->cones = NULL; 1394 mesh->coneOrientations = NULL; 1395 ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm), &mesh->supportSection);CHKERRQ(ierr); 1396 mesh->maxSupportSize = 0; 1397 mesh->supports = NULL; 1398 mesh->refinementUniform = PETSC_TRUE; 1399 mesh->refinementLimit = -1.0; 1400 1401 mesh->facesTmp = NULL; 1402 1403 mesh->tetgenOpts = NULL; 1404 mesh->triangleOpts = NULL; 1405 ierr = PetscPartitionerCreate(PetscObjectComm((PetscObject)dm), &mesh->partitioner);CHKERRQ(ierr); 1406 ierr = PetscPartitionerSetTypeFromOptions_Internal(mesh->partitioner);CHKERRQ(ierr); 1407 mesh->remeshBd = PETSC_FALSE; 1408 1409 mesh->subpointMap = NULL; 1410 1411 for (unit = 0; unit < NUM_PETSC_UNITS; ++unit) mesh->scale[unit] = 1.0; 1412 1413 mesh->regularRefinement = PETSC_FALSE; 1414 mesh->depthState = -1; 1415 mesh->globalVertexNumbers = NULL; 1416 mesh->globalCellNumbers = NULL; 1417 mesh->anchorSection = NULL; 1418 mesh->anchorIS = NULL; 1419 mesh->createanchors = NULL; 1420 mesh->computeanchormatrix = NULL; 1421 mesh->parentSection = NULL; 1422 mesh->parents = NULL; 1423 mesh->childIDs = NULL; 1424 mesh->childSection = NULL; 1425 mesh->children = NULL; 1426 mesh->referenceTree = NULL; 1427 mesh->getchildsymmetry = NULL; 1428 for (d = 0; d < 8; ++d) mesh->hybridPointMax[d] = PETSC_DETERMINE; 1429 mesh->vtkCellHeight = 0; 1430 mesh->useCone = PETSC_FALSE; 1431 mesh->useClosure = PETSC_TRUE; 1432 mesh->useAnchors = PETSC_FALSE; 1433 1434 mesh->maxProjectionHeight = 0; 1435 1436 mesh->printSetValues = PETSC_FALSE; 1437 mesh->printFEM = 0; 1438 mesh->printTol = 1.0e-10; 1439 1440 ierr = DMInitialize_Plex(dm);CHKERRQ(ierr); 1441 PetscFunctionReturn(0); 1442 } 1443 1444 #undef __FUNCT__ 1445 #define __FUNCT__ "DMPlexCreate" 1446 /*@ 1447 DMPlexCreate - Creates a DMPlex object, which encapsulates an unstructured mesh, or CW complex, which can be expressed using a Hasse Diagram. 1448 1449 Collective on MPI_Comm 1450 1451 Input Parameter: 1452 . comm - The communicator for the DMPlex object 1453 1454 Output Parameter: 1455 . mesh - The DMPlex object 1456 1457 Level: beginner 1458 1459 .keywords: DMPlex, create 1460 @*/ 1461 PetscErrorCode DMPlexCreate(MPI_Comm comm, DM *mesh) 1462 { 1463 PetscErrorCode ierr; 1464 1465 PetscFunctionBegin; 1466 PetscValidPointer(mesh,2); 1467 ierr = DMCreate(comm, mesh);CHKERRQ(ierr); 1468 ierr = DMSetType(*mesh, DMPLEX);CHKERRQ(ierr); 1469 PetscFunctionReturn(0); 1470 } 1471 1472 #undef __FUNCT__ 1473 #define __FUNCT__ "DMPlexBuildFromCellList_Parallel_Private" 1474 /* 1475 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 1476 */ 1477 static PetscErrorCode DMPlexBuildFromCellList_Parallel_Private(DM dm, PetscInt numCells, PetscInt numVertices, PetscInt numCorners, const int cells[], PetscSF *sfVert) 1478 { 1479 PetscSF sfPoint; 1480 PetscLayout vLayout; 1481 PetscHashI vhash; 1482 PetscSFNode *remoteVerticesAdj, *vertexLocal, *vertexOwner, *remoteVertex; 1483 const PetscInt *vrange; 1484 PetscInt numVerticesAdj, off, *verticesAdj, numVerticesGhost = 0, *localVertex, *cone, c, p, v, g; 1485 PETSC_UNUSED PetscHashIIter ret, iter; 1486 PetscMPIInt rank, numProcs; 1487 PetscErrorCode ierr; 1488 1489 PetscFunctionBegin; 1490 ierr = MPI_Comm_rank(PetscObjectComm((PetscObject) dm), &rank);CHKERRQ(ierr); 1491 ierr = MPI_Comm_size(PetscObjectComm((PetscObject) dm), &numProcs);CHKERRQ(ierr); 1492 /* Partition vertices */ 1493 ierr = PetscLayoutCreate(PetscObjectComm((PetscObject) dm), &vLayout);CHKERRQ(ierr); 1494 ierr = PetscLayoutSetLocalSize(vLayout, numVertices);CHKERRQ(ierr); 1495 ierr = PetscLayoutSetBlockSize(vLayout, 1);CHKERRQ(ierr); 1496 ierr = PetscLayoutSetUp(vLayout);CHKERRQ(ierr); 1497 ierr = PetscLayoutGetRanges(vLayout, &vrange);CHKERRQ(ierr); 1498 /* Count vertices and map them to procs */ 1499 PetscHashICreate(vhash); 1500 for (c = 0; c < numCells; ++c) { 1501 for (p = 0; p < numCorners; ++p) { 1502 PetscHashIPut(vhash, cells[c*numCorners+p], ret, iter); 1503 } 1504 } 1505 PetscHashISize(vhash, numVerticesAdj); 1506 ierr = PetscMalloc1(numVerticesAdj, &verticesAdj);CHKERRQ(ierr); 1507 off = 0; ierr = PetscHashIGetKeys(vhash, &off, verticesAdj);CHKERRQ(ierr); 1508 if (off != numVerticesAdj) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid number of local vertices %D should be %D", off, numVerticesAdj); 1509 ierr = PetscMalloc1(numVerticesAdj, &remoteVerticesAdj);CHKERRQ(ierr); 1510 for (v = 0; v < numVerticesAdj; ++v) { 1511 const PetscInt gv = verticesAdj[v]; 1512 PetscInt vrank; 1513 1514 ierr = PetscFindInt(gv, numProcs+1, vrange, &vrank);CHKERRQ(ierr); 1515 vrank = vrank < 0 ? -(vrank+2) : vrank; 1516 remoteVerticesAdj[v].index = gv - vrange[vrank]; 1517 remoteVerticesAdj[v].rank = vrank; 1518 } 1519 /* Create cones */ 1520 ierr = DMPlexSetChart(dm, 0, numCells+numVerticesAdj);CHKERRQ(ierr); 1521 for (c = 0; c < numCells; ++c) {ierr = DMPlexSetConeSize(dm, c, numCorners);CHKERRQ(ierr);} 1522 ierr = DMSetUp(dm);CHKERRQ(ierr); 1523 ierr = DMGetWorkArray(dm, numCorners, PETSC_INT, &cone);CHKERRQ(ierr); 1524 for (c = 0; c < numCells; ++c) { 1525 for (p = 0; p < numCorners; ++p) { 1526 const PetscInt gv = cells[c*numCorners+p]; 1527 PetscInt lv; 1528 1529 ierr = PetscFindInt(gv, numVerticesAdj, verticesAdj, &lv);CHKERRQ(ierr); 1530 if (lv < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Could not find global vertex %D in local connectivity", gv); 1531 cone[p] = lv+numCells; 1532 } 1533 ierr = DMPlexSetCone(dm, c, cone);CHKERRQ(ierr); 1534 } 1535 ierr = DMRestoreWorkArray(dm, numCorners, PETSC_INT, &cone);CHKERRQ(ierr); 1536 /* Create SF for vertices */ 1537 ierr = PetscSFCreate(PetscObjectComm((PetscObject)dm), sfVert);CHKERRQ(ierr); 1538 ierr = PetscObjectSetName((PetscObject) *sfVert, "Vertex Ownership SF");CHKERRQ(ierr); 1539 ierr = PetscSFSetFromOptions(*sfVert);CHKERRQ(ierr); 1540 ierr = PetscSFSetGraph(*sfVert, numVertices, numVerticesAdj, NULL, PETSC_OWN_POINTER, remoteVerticesAdj, PETSC_OWN_POINTER);CHKERRQ(ierr); 1541 ierr = PetscFree(verticesAdj);CHKERRQ(ierr); 1542 /* Build pointSF */ 1543 ierr = PetscMalloc2(numVerticesAdj, &vertexLocal, numVertices, &vertexOwner);CHKERRQ(ierr); 1544 for (v = 0; v < numVerticesAdj; ++v) {vertexLocal[v].index = v+numCells; vertexLocal[v].rank = rank;} 1545 for (v = 0; v < numVertices; ++v) {vertexOwner[v].index = -1; vertexOwner[v].rank = -1;} 1546 ierr = PetscSFReduceBegin(*sfVert, MPI_2INT, vertexLocal, vertexOwner, MPI_MAXLOC);CHKERRQ(ierr); 1547 ierr = PetscSFReduceEnd(*sfVert, MPI_2INT, vertexLocal, vertexOwner, MPI_MAXLOC);CHKERRQ(ierr); 1548 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); 1549 ierr = PetscSFBcastBegin(*sfVert, MPI_2INT, vertexOwner, vertexLocal);CHKERRQ(ierr); 1550 ierr = PetscSFBcastEnd(*sfVert, MPI_2INT, vertexOwner, vertexLocal);CHKERRQ(ierr); 1551 for (v = 0; v < numVerticesAdj; ++v) if (vertexLocal[v].rank != rank) ++numVerticesGhost; 1552 ierr = PetscMalloc1(numVerticesGhost, &localVertex);CHKERRQ(ierr); 1553 ierr = PetscMalloc1(numVerticesGhost, &remoteVertex);CHKERRQ(ierr); 1554 for (v = 0, g = 0; v < numVerticesAdj; ++v) { 1555 if (vertexLocal[v].rank != rank) { 1556 localVertex[g] = v+numCells; 1557 remoteVertex[g].index = vertexLocal[v].index; 1558 remoteVertex[g].rank = vertexLocal[v].rank; 1559 ++g; 1560 } 1561 } 1562 ierr = PetscFree2(vertexLocal, vertexOwner);CHKERRQ(ierr); 1563 if (g != numVerticesGhost) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid number of vertex ghosts %D should be %D", g, numVerticesGhost); 1564 ierr = DMGetPointSF(dm, &sfPoint);CHKERRQ(ierr); 1565 ierr = PetscObjectSetName((PetscObject) sfPoint, "point SF");CHKERRQ(ierr); 1566 ierr = PetscSFSetGraph(sfPoint, numCells+numVerticesAdj, numVerticesGhost, localVertex, PETSC_OWN_POINTER, remoteVertex, PETSC_OWN_POINTER);CHKERRQ(ierr); 1567 ierr = PetscLayoutDestroy(&vLayout);CHKERRQ(ierr); 1568 PetscHashIDestroy(vhash); 1569 /* Fill in the rest of the topology structure */ 1570 ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr); 1571 ierr = DMPlexStratify(dm);CHKERRQ(ierr); 1572 PetscFunctionReturn(0); 1573 } 1574 1575 #undef __FUNCT__ 1576 #define __FUNCT__ "DMPlexBuildCoordinates_Parallel_Private" 1577 /* 1578 This takes as input the coordinates for each owned vertex 1579 */ 1580 static PetscErrorCode DMPlexBuildCoordinates_Parallel_Private(DM dm, PetscInt spaceDim, PetscInt numCells, PetscSF sfVert, const double vertexCoords[]) 1581 { 1582 PetscSection coordSection; 1583 Vec coordinates; 1584 PetscScalar *coords; 1585 PetscInt numVertices, numVerticesAdj, coordSize, v; 1586 PetscErrorCode ierr; 1587 1588 PetscFunctionBegin; 1589 ierr = PetscSFGetGraph(sfVert, &numVertices, &numVerticesAdj, NULL, NULL);CHKERRQ(ierr); 1590 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 1591 ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr); 1592 ierr = PetscSectionSetFieldComponents(coordSection, 0, spaceDim);CHKERRQ(ierr); 1593 ierr = PetscSectionSetChart(coordSection, numCells, numCells + numVerticesAdj);CHKERRQ(ierr); 1594 for (v = numCells; v < numCells+numVerticesAdj; ++v) { 1595 ierr = PetscSectionSetDof(coordSection, v, spaceDim);CHKERRQ(ierr); 1596 ierr = PetscSectionSetFieldDof(coordSection, v, 0, spaceDim);CHKERRQ(ierr); 1597 } 1598 ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 1599 ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 1600 ierr = VecCreate(PetscObjectComm((PetscObject)dm), &coordinates);CHKERRQ(ierr); 1601 ierr = VecSetBlockSize(coordinates, spaceDim);CHKERRQ(ierr); 1602 ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr); 1603 ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 1604 ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr); 1605 ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 1606 { 1607 MPI_Datatype coordtype; 1608 1609 /* Need a temp buffer for coords if we have complex/single */ 1610 ierr = MPI_Type_contiguous(spaceDim, MPI_DOUBLE, &coordtype);CHKERRQ(ierr); 1611 ierr = MPI_Type_commit(&coordtype);CHKERRQ(ierr); 1612 ierr = PetscSFBcastBegin(sfVert, coordtype, vertexCoords, coords);CHKERRQ(ierr); 1613 ierr = PetscSFBcastEnd(sfVert, coordtype, vertexCoords, coords);CHKERRQ(ierr); 1614 ierr = MPI_Type_free(&coordtype);CHKERRQ(ierr); 1615 } 1616 ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 1617 ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr); 1618 ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 1619 PetscFunctionReturn(0); 1620 } 1621 1622 #undef __FUNCT__ 1623 #define __FUNCT__ "DMPlexCreateFromCellListParallel" 1624 /*@C 1625 DMPlexCreateFromCellListParallel - This takes as input common mesh generator output, a list of the vertices for each cell, and produces a DM 1626 1627 Input Parameters: 1628 + comm - The communicator 1629 . dim - The topological dimension of the mesh 1630 . numCells - The number of cells owned by this process 1631 . numVertices - The number of vertices owned by this process 1632 . numCorners - The number of vertices for each cell 1633 . interpolate - Flag indicating that intermediate mesh entities (faces, edges) should be created automatically 1634 . cells - An array of numCells*numCorners numbers, the global vertex numbers for each cell 1635 . spaceDim - The spatial dimension used for coordinates 1636 - vertexCoords - An array of numVertices*spaceDim numbers, the coordinates of each vertex 1637 1638 Output Parameter: 1639 . dm - The DM 1640 1641 Note: Two triangles sharing a face 1642 $ 1643 $ 2 1644 $ / | \ 1645 $ / | \ 1646 $ / | \ 1647 $ 0 0 | 1 3 1648 $ \ | / 1649 $ \ | / 1650 $ \ | / 1651 $ 1 1652 would have input 1653 $ numCells = 2, numVertices = 4 1654 $ cells = [0 1 2 1 3 2] 1655 $ 1656 which would result in the DMPlex 1657 $ 1658 $ 4 1659 $ / | \ 1660 $ / | \ 1661 $ / | \ 1662 $ 2 0 | 1 5 1663 $ \ | / 1664 $ \ | / 1665 $ \ | / 1666 $ 3 1667 1668 Level: beginner 1669 1670 .seealso: DMPlexCreateFromCellList(), DMPlexCreateFromDAG(), DMPlexCreate() 1671 @*/ 1672 PetscErrorCode DMPlexCreateFromCellListParallel(MPI_Comm comm, PetscInt dim, PetscInt numCells, PetscInt numVertices, PetscInt numCorners, PetscBool interpolate, const int cells[], PetscInt spaceDim, const double vertexCoords[], DM *dm) 1673 { 1674 PetscSF sfVert; 1675 PetscErrorCode ierr; 1676 1677 PetscFunctionBegin; 1678 ierr = DMCreate(comm, dm);CHKERRQ(ierr); 1679 ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 1680 PetscValidLogicalCollectiveInt(*dm, dim, 2); 1681 PetscValidLogicalCollectiveInt(*dm, spaceDim, 8); 1682 ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr); 1683 ierr = DMPlexBuildFromCellList_Parallel_Private(*dm, numCells, numVertices, numCorners, cells, &sfVert);CHKERRQ(ierr); 1684 if (interpolate) { 1685 DM idm = NULL; 1686 1687 ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr); 1688 ierr = DMDestroy(dm);CHKERRQ(ierr); 1689 *dm = idm; 1690 } 1691 ierr = DMPlexBuildCoordinates_Parallel_Private(*dm, spaceDim, numCells, sfVert, vertexCoords);CHKERRQ(ierr); 1692 ierr = PetscSFDestroy(&sfVert);CHKERRQ(ierr); 1693 PetscFunctionReturn(0); 1694 } 1695 1696 #undef __FUNCT__ 1697 #define __FUNCT__ "DMPlexBuildFromCellList_Private" 1698 /* 1699 This takes as input the common mesh generator output, a list of the vertices for each cell 1700 */ 1701 static PetscErrorCode DMPlexBuildFromCellList_Private(DM dm, PetscInt numCells, PetscInt numVertices, PetscInt numCorners, const int cells[]) 1702 { 1703 PetscInt *cone, c, p; 1704 PetscErrorCode ierr; 1705 1706 PetscFunctionBegin; 1707 ierr = DMPlexSetChart(dm, 0, numCells+numVertices);CHKERRQ(ierr); 1708 for (c = 0; c < numCells; ++c) { 1709 ierr = DMPlexSetConeSize(dm, c, numCorners);CHKERRQ(ierr); 1710 } 1711 ierr = DMSetUp(dm);CHKERRQ(ierr); 1712 ierr = DMGetWorkArray(dm, numCorners, PETSC_INT, &cone);CHKERRQ(ierr); 1713 for (c = 0; c < numCells; ++c) { 1714 for (p = 0; p < numCorners; ++p) { 1715 cone[p] = cells[c*numCorners+p]+numCells; 1716 } 1717 ierr = DMPlexSetCone(dm, c, cone);CHKERRQ(ierr); 1718 } 1719 ierr = DMRestoreWorkArray(dm, numCorners, PETSC_INT, &cone);CHKERRQ(ierr); 1720 ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr); 1721 ierr = DMPlexStratify(dm);CHKERRQ(ierr); 1722 PetscFunctionReturn(0); 1723 } 1724 1725 #undef __FUNCT__ 1726 #define __FUNCT__ "DMPlexBuildCoordinates_Private" 1727 /* 1728 This takes as input the coordinates for each vertex 1729 */ 1730 static PetscErrorCode DMPlexBuildCoordinates_Private(DM dm, PetscInt spaceDim, PetscInt numCells, PetscInt numVertices, const double vertexCoords[]) 1731 { 1732 PetscSection coordSection; 1733 Vec coordinates; 1734 PetscScalar *coords; 1735 PetscInt coordSize, v, d; 1736 PetscErrorCode ierr; 1737 1738 PetscFunctionBegin; 1739 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 1740 ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr); 1741 ierr = PetscSectionSetFieldComponents(coordSection, 0, spaceDim);CHKERRQ(ierr); 1742 ierr = PetscSectionSetChart(coordSection, numCells, numCells + numVertices);CHKERRQ(ierr); 1743 for (v = numCells; v < numCells+numVertices; ++v) { 1744 ierr = PetscSectionSetDof(coordSection, v, spaceDim);CHKERRQ(ierr); 1745 ierr = PetscSectionSetFieldDof(coordSection, v, 0, spaceDim);CHKERRQ(ierr); 1746 } 1747 ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 1748 ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 1749 ierr = VecCreate(PetscObjectComm((PetscObject)dm), &coordinates);CHKERRQ(ierr); 1750 ierr = VecSetBlockSize(coordinates, spaceDim);CHKERRQ(ierr); 1751 ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr); 1752 ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 1753 ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr); 1754 ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 1755 for (v = 0; v < numVertices; ++v) { 1756 for (d = 0; d < spaceDim; ++d) { 1757 coords[v*spaceDim+d] = vertexCoords[v*spaceDim+d]; 1758 } 1759 } 1760 ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 1761 ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr); 1762 ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 1763 PetscFunctionReturn(0); 1764 } 1765 1766 #undef __FUNCT__ 1767 #define __FUNCT__ "DMPlexCreateFromCellList" 1768 /*@C 1769 DMPlexCreateFromCellList - This takes as input common mesh generator output, a list of the vertices for each cell, and produces a DM 1770 1771 Input Parameters: 1772 + comm - The communicator 1773 . dim - The topological dimension of the mesh 1774 . numCells - The number of cells 1775 . numVertices - The number of vertices 1776 . numCorners - The number of vertices for each cell 1777 . interpolate - Flag indicating that intermediate mesh entities (faces, edges) should be created automatically 1778 . cells - An array of numCells*numCorners numbers, the vertices for each cell 1779 . spaceDim - The spatial dimension used for coordinates 1780 - vertexCoords - An array of numVertices*spaceDim numbers, the coordinates of each vertex 1781 1782 Output Parameter: 1783 . dm - The DM 1784 1785 Note: Two triangles sharing a face 1786 $ 1787 $ 2 1788 $ / | \ 1789 $ / | \ 1790 $ / | \ 1791 $ 0 0 | 1 3 1792 $ \ | / 1793 $ \ | / 1794 $ \ | / 1795 $ 1 1796 would have input 1797 $ numCells = 2, numVertices = 4 1798 $ cells = [0 1 2 1 3 2] 1799 $ 1800 which would result in the DMPlex 1801 $ 1802 $ 4 1803 $ / | \ 1804 $ / | \ 1805 $ / | \ 1806 $ 2 0 | 1 5 1807 $ \ | / 1808 $ \ | / 1809 $ \ | / 1810 $ 3 1811 1812 Level: beginner 1813 1814 .seealso: DMPlexCreateFromDAG(), DMPlexCreate() 1815 @*/ 1816 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) 1817 { 1818 PetscErrorCode ierr; 1819 1820 PetscFunctionBegin; 1821 ierr = DMCreate(comm, dm);CHKERRQ(ierr); 1822 ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 1823 ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr); 1824 ierr = DMPlexBuildFromCellList_Private(*dm, numCells, numVertices, numCorners, cells);CHKERRQ(ierr); 1825 if (interpolate) { 1826 DM idm = NULL; 1827 1828 ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr); 1829 ierr = DMDestroy(dm);CHKERRQ(ierr); 1830 *dm = idm; 1831 } 1832 ierr = DMPlexBuildCoordinates_Private(*dm, spaceDim, numCells, numVertices, vertexCoords);CHKERRQ(ierr); 1833 PetscFunctionReturn(0); 1834 } 1835 1836 #undef __FUNCT__ 1837 #define __FUNCT__ "DMPlexCreateFromDAG" 1838 /*@ 1839 DMPlexCreateFromDAG - This takes as input the adjacency-list representation of the Directed Acyclic Graph (Hasse Diagram) encoding a mesh, and produces a DM 1840 1841 Input Parameters: 1842 + dm - The empty DM object, usually from DMCreate() and DMSetDimension() 1843 . depth - The depth of the DAG 1844 . numPoints - The number of points at each depth 1845 . coneSize - The cone size of each point 1846 . cones - The concatenation of the cone points for each point, the cone list must be oriented correctly for each point 1847 . coneOrientations - The orientation of each cone point 1848 - vertexCoords - An array of numVertices*dim numbers, the coordinates of each vertex 1849 1850 Output Parameter: 1851 . dm - The DM 1852 1853 Note: Two triangles sharing a face would have input 1854 $ depth = 1, numPoints = [4 2], coneSize = [3 3 0 0 0 0] 1855 $ cones = [2 3 4 3 5 4], coneOrientations = [0 0 0 0 0 0] 1856 $ vertexCoords = [-1.0 0.0 0.0 -1.0 0.0 1.0 1.0 0.0] 1857 $ 1858 which would result in the DMPlex 1859 $ 1860 $ 4 1861 $ / | \ 1862 $ / | \ 1863 $ / | \ 1864 $ 2 0 | 1 5 1865 $ \ | / 1866 $ \ | / 1867 $ \ | / 1868 $ 3 1869 $ 1870 $ Notice that all points are numbered consecutively, unlikely DMPlexCreateFromCellList() 1871 1872 Level: advanced 1873 1874 .seealso: DMPlexCreateFromCellList(), DMPlexCreate() 1875 @*/ 1876 PetscErrorCode DMPlexCreateFromDAG(DM dm, PetscInt depth, const PetscInt numPoints[], const PetscInt coneSize[], const PetscInt cones[], const PetscInt coneOrientations[], const PetscScalar vertexCoords[]) 1877 { 1878 Vec coordinates; 1879 PetscSection coordSection; 1880 PetscScalar *coords; 1881 PetscInt coordSize, firstVertex = -1, pStart = 0, pEnd = 0, p, v, dim, dimEmbed, d, off; 1882 PetscErrorCode ierr; 1883 1884 PetscFunctionBegin; 1885 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 1886 ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr); 1887 if (dimEmbed < dim) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Embedding dimension %d cannot be less than intrinsic dimension %d",dimEmbed,dim); 1888 for (d = 0; d <= depth; ++d) pEnd += numPoints[d]; 1889 ierr = DMPlexSetChart(dm, pStart, pEnd);CHKERRQ(ierr); 1890 for (p = pStart; p < pEnd; ++p) { 1891 ierr = DMPlexSetConeSize(dm, p, coneSize[p-pStart]);CHKERRQ(ierr); 1892 if (firstVertex < 0 && !coneSize[p - pStart]) { 1893 firstVertex = p - pStart; 1894 } 1895 } 1896 if (firstVertex < 0 && numPoints[0]) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Expected %d vertices but could not find any", numPoints[0]); 1897 ierr = DMSetUp(dm);CHKERRQ(ierr); /* Allocate space for cones */ 1898 for (p = pStart, off = 0; p < pEnd; off += coneSize[p-pStart], ++p) { 1899 ierr = DMPlexSetCone(dm, p, &cones[off]);CHKERRQ(ierr); 1900 ierr = DMPlexSetConeOrientation(dm, p, &coneOrientations[off]);CHKERRQ(ierr); 1901 } 1902 ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr); 1903 ierr = DMPlexStratify(dm);CHKERRQ(ierr); 1904 /* Build coordinates */ 1905 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 1906 ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr); 1907 ierr = PetscSectionSetFieldComponents(coordSection, 0, dimEmbed);CHKERRQ(ierr); 1908 ierr = PetscSectionSetChart(coordSection, firstVertex, firstVertex+numPoints[0]);CHKERRQ(ierr); 1909 for (v = firstVertex; v < firstVertex+numPoints[0]; ++v) { 1910 ierr = PetscSectionSetDof(coordSection, v, dimEmbed);CHKERRQ(ierr); 1911 ierr = PetscSectionSetFieldDof(coordSection, v, 0, dimEmbed);CHKERRQ(ierr); 1912 } 1913 ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 1914 ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 1915 ierr = VecCreate(PetscObjectComm((PetscObject)dm), &coordinates);CHKERRQ(ierr); 1916 ierr = VecSetBlockSize(coordinates, dimEmbed);CHKERRQ(ierr); 1917 ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr); 1918 ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 1919 ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr); 1920 ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 1921 for (v = 0; v < numPoints[0]; ++v) { 1922 PetscInt off; 1923 1924 ierr = PetscSectionGetOffset(coordSection, v+firstVertex, &off);CHKERRQ(ierr); 1925 for (d = 0; d < dimEmbed; ++d) { 1926 coords[off+d] = vertexCoords[v*dimEmbed+d]; 1927 } 1928 } 1929 ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 1930 ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr); 1931 ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 1932 PetscFunctionReturn(0); 1933 } 1934 1935 #undef __FUNCT__ 1936 #define __FUNCT__ "DMPlexCreateFromFile" 1937 /*@C 1938 DMPlexCreateFromFile - This takes a filename and produces a DM 1939 1940 Input Parameters: 1941 + comm - The communicator 1942 . filename - A file name 1943 - interpolate - Flag to create intermediate mesh pieces (edges, faces) 1944 1945 Output Parameter: 1946 . dm - The DM 1947 1948 Level: beginner 1949 1950 .seealso: DMPlexCreateFromDAG(), DMPlexCreateFromCellList(), DMPlexCreate() 1951 @*/ 1952 PetscErrorCode DMPlexCreateFromFile(MPI_Comm comm, const char filename[], PetscBool interpolate, DM *dm) 1953 { 1954 const char *extGmsh = ".msh"; 1955 const char *extCGNS = ".cgns"; 1956 const char *extExodus = ".exo"; 1957 const char *extFluent = ".cas"; 1958 const char *extHDF5 = ".h5"; 1959 size_t len; 1960 PetscBool isGmsh, isCGNS, isExodus, isFluent, isHDF5; 1961 PetscMPIInt rank; 1962 PetscErrorCode ierr; 1963 1964 PetscFunctionBegin; 1965 PetscValidPointer(filename, 2); 1966 PetscValidPointer(dm, 4); 1967 ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 1968 ierr = PetscStrlen(filename, &len);CHKERRQ(ierr); 1969 if (!len) SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Filename must be a valid path"); 1970 ierr = PetscStrncmp(&filename[PetscMax(0,len-4)], extGmsh, 4, &isGmsh);CHKERRQ(ierr); 1971 ierr = PetscStrncmp(&filename[PetscMax(0,len-5)], extCGNS, 5, &isCGNS);CHKERRQ(ierr); 1972 ierr = PetscStrncmp(&filename[PetscMax(0,len-4)], extExodus, 4, &isExodus);CHKERRQ(ierr); 1973 ierr = PetscStrncmp(&filename[PetscMax(0,len-4)], extFluent, 4, &isFluent);CHKERRQ(ierr); 1974 ierr = PetscStrncmp(&filename[PetscMax(0,len-3)], extHDF5, 3, &isHDF5);CHKERRQ(ierr); 1975 if (isGmsh) { 1976 ierr = DMPlexCreateGmshFromFile(comm, filename, interpolate, dm);CHKERRQ(ierr); 1977 } else if (isCGNS) { 1978 ierr = DMPlexCreateCGNSFromFile(comm, filename, interpolate, dm);CHKERRQ(ierr); 1979 } else if (isExodus) { 1980 ierr = DMPlexCreateExodusFromFile(comm, filename, interpolate, dm);CHKERRQ(ierr); 1981 } else if (isFluent) { 1982 ierr = DMPlexCreateFluentFromFile(comm, filename, interpolate, dm);CHKERRQ(ierr); 1983 } else if (isHDF5) { 1984 PetscViewer viewer; 1985 1986 ierr = PetscViewerCreate(comm, &viewer);CHKERRQ(ierr); 1987 ierr = PetscViewerSetType(viewer, PETSCVIEWERHDF5);CHKERRQ(ierr); 1988 ierr = PetscViewerFileSetMode(viewer, FILE_MODE_READ);CHKERRQ(ierr); 1989 ierr = PetscViewerFileSetName(viewer, filename);CHKERRQ(ierr); 1990 ierr = DMCreate(comm, dm);CHKERRQ(ierr); 1991 ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 1992 ierr = DMLoad(*dm, viewer);CHKERRQ(ierr); 1993 ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 1994 } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cannot load file %s: unrecognized extension", filename); 1995 PetscFunctionReturn(0); 1996 } 1997 1998 #undef __FUNCT__ 1999 #define __FUNCT__ "DMPlexCreateReferenceCell" 2000 /*@ 2001 DMPlexCreateReferenceCell - Create a DMPLEX with the appropriate FEM reference cell 2002 2003 Collective on comm 2004 2005 Input Parameters: 2006 + comm - The communicator 2007 . dim - The spatial dimension 2008 - simplex - Flag for simplex, otherwise use a tensor-product cell 2009 2010 Output Parameter: 2011 . refdm - The reference cell 2012 2013 Level: intermediate 2014 2015 .keywords: reference cell 2016 .seealso: 2017 @*/ 2018 PetscErrorCode DMPlexCreateReferenceCell(MPI_Comm comm, PetscInt dim, PetscBool simplex, DM *refdm) 2019 { 2020 DM rdm; 2021 Vec coords; 2022 PetscErrorCode ierr; 2023 2024 PetscFunctionBeginUser; 2025 ierr = DMCreate(comm, &rdm);CHKERRQ(ierr); 2026 ierr = DMSetType(rdm, DMPLEX);CHKERRQ(ierr); 2027 ierr = DMSetDimension(rdm, dim);CHKERRQ(ierr); 2028 switch (dim) { 2029 case 0: 2030 { 2031 PetscInt numPoints[1] = {1}; 2032 PetscInt coneSize[1] = {0}; 2033 PetscInt cones[1] = {0}; 2034 PetscInt coneOrientations[1] = {0}; 2035 PetscScalar vertexCoords[1] = {0.0}; 2036 2037 ierr = DMPlexCreateFromDAG(rdm, 0, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr); 2038 } 2039 break; 2040 case 1: 2041 { 2042 PetscInt numPoints[2] = {2, 1}; 2043 PetscInt coneSize[3] = {2, 0, 0}; 2044 PetscInt cones[2] = {1, 2}; 2045 PetscInt coneOrientations[2] = {0, 0}; 2046 PetscScalar vertexCoords[2] = {-1.0, 1.0}; 2047 2048 ierr = DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr); 2049 } 2050 break; 2051 case 2: 2052 if (simplex) { 2053 PetscInt numPoints[2] = {3, 1}; 2054 PetscInt coneSize[4] = {3, 0, 0, 0}; 2055 PetscInt cones[3] = {1, 2, 3}; 2056 PetscInt coneOrientations[3] = {0, 0, 0}; 2057 PetscScalar vertexCoords[6] = {-1.0, -1.0, 1.0, -1.0, -1.0, 1.0}; 2058 2059 ierr = DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr); 2060 } else { 2061 PetscInt numPoints[2] = {4, 1}; 2062 PetscInt coneSize[5] = {4, 0, 0, 0, 0}; 2063 PetscInt cones[4] = {1, 2, 3, 4}; 2064 PetscInt coneOrientations[4] = {0, 0, 0, 0}; 2065 PetscScalar vertexCoords[8] = {-1.0, -1.0, 1.0, -1.0, 1.0, 1.0, -1.0, 1.0}; 2066 2067 ierr = DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr); 2068 } 2069 break; 2070 case 3: 2071 if (simplex) { 2072 PetscInt numPoints[2] = {4, 1}; 2073 PetscInt coneSize[5] = {4, 0, 0, 0, 0}; 2074 PetscInt cones[4] = {1, 3, 2, 4}; 2075 PetscInt coneOrientations[4] = {0, 0, 0, 0}; 2076 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}; 2077 2078 ierr = DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr); 2079 } else { 2080 PetscInt numPoints[2] = {8, 1}; 2081 PetscInt coneSize[9] = {8, 0, 0, 0, 0, 0, 0, 0, 0}; 2082 PetscInt cones[8] = {1, 4, 3, 2, 5, 6, 7, 8}; 2083 PetscInt coneOrientations[8] = {0, 0, 0, 0, 0, 0, 0, 0}; 2084 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, 2085 -1.0, -1.0, 1.0, 1.0, -1.0, 1.0, 1.0, 1.0, 1.0, -1.0, 1.0, 1.0}; 2086 2087 ierr = DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr); 2088 } 2089 break; 2090 default: 2091 SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Cannot create reference cell for dimension %d", dim); 2092 } 2093 *refdm = NULL; 2094 ierr = DMPlexInterpolate(rdm, refdm);CHKERRQ(ierr); 2095 if (rdm->coordinateDM) { 2096 DM ncdm; 2097 PetscSection cs; 2098 PetscInt pEnd = -1; 2099 2100 ierr = DMGetDefaultSection(rdm->coordinateDM, &cs);CHKERRQ(ierr); 2101 if (cs) {ierr = PetscSectionGetChart(cs, NULL, &pEnd);CHKERRQ(ierr);} 2102 if (pEnd >= 0) { 2103 ierr = DMClone(rdm->coordinateDM, &ncdm);CHKERRQ(ierr); 2104 ierr = DMSetDefaultSection(ncdm, cs);CHKERRQ(ierr); 2105 ierr = DMSetCoordinateDM(*refdm, ncdm);CHKERRQ(ierr); 2106 ierr = DMDestroy(&ncdm);CHKERRQ(ierr); 2107 } 2108 } 2109 ierr = DMGetCoordinatesLocal(rdm, &coords);CHKERRQ(ierr); 2110 if (coords) { 2111 ierr = DMSetCoordinatesLocal(*refdm, coords);CHKERRQ(ierr); 2112 } else { 2113 ierr = DMGetCoordinates(rdm, &coords);CHKERRQ(ierr); 2114 if (coords) {ierr = DMSetCoordinates(*refdm, coords);CHKERRQ(ierr);} 2115 } 2116 ierr = DMDestroy(&rdm);CHKERRQ(ierr); 2117 PetscFunctionReturn(0); 2118 } 2119