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__ "DMSetFromOptions_Plex" 8 PetscErrorCode DMSetFromOptions_Plex(DM dm) 9 { 10 DM_Plex *mesh = (DM_Plex*) dm->data; 11 PetscErrorCode ierr; 12 13 PetscFunctionBegin; 14 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 15 ierr = PetscOptionsHead("DMPlex Options");CHKERRQ(ierr); 16 /* Handle DMPlex refinement */ 17 /* Handle associated vectors */ 18 /* Handle viewing */ 19 ierr = PetscOptionsBool("-dm_plex_print_set_values", "Output all set values info", "DMView", PETSC_FALSE, &mesh->printSetValues, NULL);CHKERRQ(ierr); 20 ierr = PetscOptionsInt("-dm_plex_print_fem", "Debug output level all fem computations", "DMView", 0, &mesh->printFEM, NULL);CHKERRQ(ierr); 21 ierr = PetscOptionsReal("-dm_plex_print_tol", "Tolerance for FEM output", "DMView", PETSC_FALSE, &mesh->printTol, NULL);CHKERRQ(ierr); 22 ierr = PetscOptionsTail();CHKERRQ(ierr); 23 PetscFunctionReturn(0); 24 } 25 26 #undef __FUNCT__ 27 #define __FUNCT__ "DMPlexCreateDoublet" 28 /*@ 29 DMPlexCreateDoublet - Creates a mesh of two cells of the specified type, optionally with later refinement. 30 31 Collective on MPI_Comm 32 33 Input Parameters: 34 + comm - The communicator for the DM object 35 . dim - The spatial dimension 36 . simplex - Flag for simplicial cells, otherwise they are tensor product cells 37 . interpolate - Flag to create intermediate mesh pieces (edges, faces) 38 . refinementUniform - Flag for uniform parallel refinement 39 - refinementLimit - A nonzero number indicates the largest admissible volume for a refined cell 40 41 Output Parameter: 42 . dm - The DM object 43 44 Level: beginner 45 46 .keywords: DM, create 47 .seealso: DMSetType(), DMCreate() 48 @*/ 49 PetscErrorCode DMPlexCreateDoublet(MPI_Comm comm, PetscInt dim, PetscBool simplex, PetscBool interpolate, PetscBool refinementUniform, PetscReal refinementLimit, DM *newdm) 50 { 51 DM dm; 52 PetscInt p; 53 PetscMPIInt rank; 54 PetscErrorCode ierr; 55 56 PetscFunctionBegin; 57 ierr = DMCreate(comm, &dm);CHKERRQ(ierr); 58 ierr = DMSetType(dm, DMPLEX);CHKERRQ(ierr); 59 ierr = DMPlexSetDimension(dm, dim);CHKERRQ(ierr); 60 ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 61 switch (dim) { 62 case 2: 63 if (simplex) {ierr = PetscObjectSetName((PetscObject) dm, "triangular");CHKERRQ(ierr);} 64 else {ierr = PetscObjectSetName((PetscObject) dm, "quadrilateral");CHKERRQ(ierr);} 65 break; 66 case 3: 67 if (simplex) {ierr = PetscObjectSetName((PetscObject) dm, "tetrahedral");CHKERRQ(ierr);} 68 else {ierr = PetscObjectSetName((PetscObject) dm, "hexahedral");CHKERRQ(ierr);} 69 break; 70 default: 71 SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Cannot make meshes for dimension %d", dim); 72 } 73 if (rank) { 74 PetscInt numPoints[2] = {0, 0}; 75 ierr = DMPlexCreateFromDAG(dm, 1, numPoints, NULL, NULL, NULL, NULL);CHKERRQ(ierr); 76 } else { 77 switch (dim) { 78 case 2: 79 if (simplex) { 80 PetscInt numPoints[2] = {4, 2}; 81 PetscInt coneSize[6] = {3, 3, 0, 0, 0, 0}; 82 PetscInt cones[6] = {2, 3, 4, 5, 4, 3}; 83 PetscInt coneOrientations[6] = {0, 0, 0, 0, 0, 0}; 84 PetscScalar vertexCoords[8] = {-0.5, 0.5, 0.0, 0.0, 0.0, 1.0, 0.5, 0.5}; 85 PetscInt markerPoints[8] = {2, 1, 3, 1, 4, 1, 5, 1}; 86 87 ierr = DMPlexCreateFromDAG(dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr); 88 for (p = 0; p < 4; ++p) {ierr = DMPlexSetLabelValue(dm, "marker", markerPoints[p*2], markerPoints[p*2+1]);CHKERRQ(ierr);} 89 } else { 90 PetscInt numPoints[2] = {6, 2}; 91 PetscInt coneSize[8] = {4, 4, 0, 0, 0, 0, 0, 0}; 92 PetscInt cones[8] = {2, 3, 4, 5, 3, 6, 7, 4}; 93 PetscInt coneOrientations[8] = {0, 0, 0, 0, 0, 0, 0, 0}; 94 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}; 95 96 ierr = DMPlexCreateFromDAG(dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr); 97 } 98 break; 99 case 3: 100 if (simplex) { 101 PetscInt numPoints[2] = {5, 2}; 102 PetscInt coneSize[7] = {4, 4, 0, 0, 0, 0, 0}; 103 PetscInt cones[8] = {4, 3, 5, 2, 5, 3, 4, 6}; 104 PetscInt coneOrientations[8] = {0, 0, 0, 0, 0, 0, 0, 0}; 105 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}; 106 PetscInt markerPoints[10] = {2, 1, 3, 1, 4, 1, 5, 1, 6, 1}; 107 108 ierr = DMPlexCreateFromDAG(dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr); 109 for (p = 0; p < 5; ++p) {ierr = DMPlexSetLabelValue(dm, "marker", markerPoints[p*2], markerPoints[p*2+1]);CHKERRQ(ierr);} 110 } else { 111 PetscInt numPoints[2] = {12, 2}; 112 PetscInt coneSize[14] = {8, 8, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0}; 113 PetscInt cones[16] = {2, 3, 4, 5, 6, 7, 8, 9, 5, 4, 10, 11, 7, 12, 13, 8}; 114 PetscInt coneOrientations[16] = {0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0}; 115 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, 116 -1.0, -0.5, 0.5, 0.0, -0.5, 0.5, 0.0, 0.5, 0.5, -1.0, 0.5, 0.5, 117 1.0, 0.5, -0.5, 1.0, -0.5, -0.5, 1.0, -0.5, 0.5, 1.0, 0.5, 0.5}; 118 119 ierr = DMPlexCreateFromDAG(dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr); 120 } 121 break; 122 default: 123 SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Cannot make meshes for dimension %d", dim); 124 } 125 } 126 *newdm = dm; 127 if (refinementLimit > 0.0) { 128 DM rdm; 129 const char *name; 130 131 ierr = DMPlexSetRefinementUniform(*newdm, PETSC_FALSE);CHKERRQ(ierr); 132 ierr = DMPlexSetRefinementLimit(*newdm, refinementLimit);CHKERRQ(ierr); 133 ierr = DMRefine(*newdm, comm, &rdm);CHKERRQ(ierr); 134 ierr = PetscObjectGetName((PetscObject) *newdm, &name);CHKERRQ(ierr); 135 ierr = PetscObjectSetName((PetscObject) rdm, name);CHKERRQ(ierr); 136 ierr = DMDestroy(newdm);CHKERRQ(ierr); 137 *newdm = rdm; 138 } 139 if (interpolate) { 140 DM idm; 141 const char *name; 142 143 ierr = DMPlexInterpolate(*newdm, &idm);CHKERRQ(ierr); 144 ierr = PetscObjectGetName((PetscObject) *newdm, &name);CHKERRQ(ierr); 145 ierr = PetscObjectSetName((PetscObject) idm, name);CHKERRQ(ierr); 146 ierr = DMPlexCopyCoordinates(*newdm, idm);CHKERRQ(ierr); 147 ierr = DMPlexCopyLabels(*newdm, idm);CHKERRQ(ierr); 148 ierr = DMDestroy(newdm);CHKERRQ(ierr); 149 *newdm = idm; 150 } 151 { 152 DM refinedMesh = NULL; 153 DM distributedMesh = NULL; 154 155 /* Distribute mesh over processes */ 156 ierr = DMPlexDistribute(*newdm, NULL, 0, NULL, &distributedMesh);CHKERRQ(ierr); 157 if (distributedMesh) { 158 ierr = DMDestroy(newdm);CHKERRQ(ierr); 159 *newdm = distributedMesh; 160 } 161 if (refinementUniform) { 162 ierr = DMPlexSetRefinementUniform(*newdm, refinementUniform);CHKERRQ(ierr); 163 ierr = DMRefine(*newdm, comm, &refinedMesh);CHKERRQ(ierr); 164 if (refinedMesh) { 165 ierr = DMDestroy(newdm);CHKERRQ(ierr); 166 *newdm = refinedMesh; 167 } 168 } 169 } 170 PetscFunctionReturn(0); 171 } 172 173 #undef __FUNCT__ 174 #define __FUNCT__ "DMPlexCreateSquareBoundary" 175 /*@ 176 DMPlexCreateSquareBoundary - Creates a 1D mesh the is the boundary of a square lattice. 177 178 Collective on MPI_Comm 179 180 Input Parameters: 181 + comm - The communicator for the DM object 182 . lower - The lower left corner coordinates 183 . upper - The upper right corner coordinates 184 - edges - The number of cells in each direction 185 186 Output Parameter: 187 . dm - The DM object 188 189 Note: Here is the numbering returned for 2 cells in each direction: 190 $ 18--5-17--4--16 191 $ | | | 192 $ 6 10 3 193 $ | | | 194 $ 19-11-20--9--15 195 $ | | | 196 $ 7 8 2 197 $ | | | 198 $ 12--0-13--1--14 199 200 Level: beginner 201 202 .keywords: DM, create 203 .seealso: DMPlexCreateBoxMesh(), DMPlexCreateCubeBoundary(), DMSetType(), DMCreate() 204 @*/ 205 PetscErrorCode DMPlexCreateSquareBoundary(DM dm, const PetscReal lower[], const PetscReal upper[], const PetscInt edges[]) 206 { 207 PetscInt numVertices = (edges[0]+1)*(edges[1]+1); 208 PetscInt numEdges = edges[0]*(edges[1]+1) + (edges[0]+1)*edges[1]; 209 PetscInt markerTop = 1; 210 PetscInt markerBottom = 1; 211 PetscInt markerRight = 1; 212 PetscInt markerLeft = 1; 213 PetscBool markerSeparate = PETSC_FALSE; 214 Vec coordinates; 215 PetscSection coordSection; 216 PetscScalar *coords; 217 PetscInt coordSize; 218 PetscMPIInt rank; 219 PetscInt v, vx, vy; 220 PetscErrorCode ierr; 221 222 PetscFunctionBegin; 223 ierr = PetscOptionsGetBool(((PetscObject) dm)->prefix, "-dm_plex_separate_marker", &markerSeparate, NULL);CHKERRQ(ierr); 224 if (markerSeparate) { 225 markerTop = 1; 226 markerBottom = 0; 227 markerRight = 0; 228 markerLeft = 0; 229 } 230 ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);CHKERRQ(ierr); 231 if (!rank) { 232 PetscInt e, ex, ey; 233 234 ierr = DMPlexSetChart(dm, 0, numEdges+numVertices);CHKERRQ(ierr); 235 for (e = 0; e < numEdges; ++e) { 236 ierr = DMPlexSetConeSize(dm, e, 2);CHKERRQ(ierr); 237 } 238 ierr = DMSetUp(dm);CHKERRQ(ierr); /* Allocate space for cones */ 239 for (vx = 0; vx <= edges[0]; vx++) { 240 for (ey = 0; ey < edges[1]; ey++) { 241 PetscInt edge = vx*edges[1] + ey + edges[0]*(edges[1]+1); 242 PetscInt vertex = ey*(edges[0]+1) + vx + numEdges; 243 PetscInt cone[2]; 244 245 cone[0] = vertex; cone[1] = vertex+edges[0]+1; 246 ierr = DMPlexSetCone(dm, edge, cone);CHKERRQ(ierr); 247 if (vx == edges[0]) { 248 ierr = DMPlexSetLabelValue(dm, "marker", edge, markerRight);CHKERRQ(ierr); 249 ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerRight);CHKERRQ(ierr); 250 if (ey == edges[1]-1) { 251 ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerRight);CHKERRQ(ierr); 252 } 253 } else if (vx == 0) { 254 ierr = DMPlexSetLabelValue(dm, "marker", edge, markerLeft);CHKERRQ(ierr); 255 ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerLeft);CHKERRQ(ierr); 256 if (ey == edges[1]-1) { 257 ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerLeft);CHKERRQ(ierr); 258 } 259 } 260 } 261 } 262 for (vy = 0; vy <= edges[1]; vy++) { 263 for (ex = 0; ex < edges[0]; ex++) { 264 PetscInt edge = vy*edges[0] + ex; 265 PetscInt vertex = vy*(edges[0]+1) + ex + numEdges; 266 PetscInt cone[2]; 267 268 cone[0] = vertex; cone[1] = vertex+1; 269 ierr = DMPlexSetCone(dm, edge, cone);CHKERRQ(ierr); 270 if (vy == edges[1]) { 271 ierr = DMPlexSetLabelValue(dm, "marker", edge, markerTop);CHKERRQ(ierr); 272 ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerTop);CHKERRQ(ierr); 273 if (ex == edges[0]-1) { 274 ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerTop);CHKERRQ(ierr); 275 } 276 } else if (vy == 0) { 277 ierr = DMPlexSetLabelValue(dm, "marker", edge, markerBottom);CHKERRQ(ierr); 278 ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerBottom);CHKERRQ(ierr); 279 if (ex == edges[0]-1) { 280 ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerBottom);CHKERRQ(ierr); 281 } 282 } 283 } 284 } 285 } 286 ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr); 287 ierr = DMPlexStratify(dm);CHKERRQ(ierr); 288 /* Build coordinates */ 289 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 290 ierr = PetscSectionSetChart(coordSection, numEdges, numEdges + numVertices);CHKERRQ(ierr); 291 for (v = numEdges; v < numEdges+numVertices; ++v) { 292 ierr = PetscSectionSetDof(coordSection, v, 2);CHKERRQ(ierr); 293 } 294 ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 295 ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 296 ierr = VecCreate(PetscObjectComm((PetscObject)dm), &coordinates);CHKERRQ(ierr); 297 ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 298 ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr); 299 ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 300 for (vy = 0; vy <= edges[1]; ++vy) { 301 for (vx = 0; vx <= edges[0]; ++vx) { 302 coords[(vy*(edges[0]+1)+vx)*2+0] = lower[0] + ((upper[0] - lower[0])/edges[0])*vx; 303 coords[(vy*(edges[0]+1)+vx)*2+1] = lower[1] + ((upper[1] - lower[1])/edges[1])*vy; 304 } 305 } 306 ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 307 ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr); 308 ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 309 PetscFunctionReturn(0); 310 } 311 312 #undef __FUNCT__ 313 #define __FUNCT__ "DMPlexCreateCubeBoundary" 314 /*@ 315 DMPlexCreateCubeBoundary - Creates a 2D mesh the is the boundary of a cubic lattice. 316 317 Collective on MPI_Comm 318 319 Input Parameters: 320 + comm - The communicator for the DM object 321 . lower - The lower left front corner coordinates 322 . upper - The upper right back corner coordinates 323 - edges - The number of cells in each direction 324 325 Output Parameter: 326 . dm - The DM object 327 328 Level: beginner 329 330 .keywords: DM, create 331 .seealso: DMPlexCreateBoxMesh(), DMPlexCreateSquareBoundary(), DMSetType(), DMCreate() 332 @*/ 333 PetscErrorCode DMPlexCreateCubeBoundary(DM dm, const PetscReal lower[], const PetscReal upper[], const PetscInt faces[]) 334 { 335 PetscInt numVertices = (faces[0]+1)*(faces[1]+1)*(faces[2]+1); 336 PetscInt numFaces = 6; 337 Vec coordinates; 338 PetscSection coordSection; 339 PetscScalar *coords; 340 PetscInt coordSize; 341 PetscMPIInt rank; 342 PetscInt v, vx, vy, vz; 343 PetscErrorCode ierr; 344 345 PetscFunctionBegin; 346 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"); 347 if ((faces[0] > 1) || (faces[1] > 1) || (faces[2] > 1)) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Currently can't handle more than 1 face per side"); 348 ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);CHKERRQ(ierr); 349 if (!rank) { 350 PetscInt f; 351 352 ierr = DMPlexSetChart(dm, 0, numFaces+numVertices);CHKERRQ(ierr); 353 for (f = 0; f < numFaces; ++f) { 354 ierr = DMPlexSetConeSize(dm, f, 4);CHKERRQ(ierr); 355 } 356 ierr = DMSetUp(dm);CHKERRQ(ierr); /* Allocate space for cones */ 357 for (v = 0; v < numFaces+numVertices; ++v) { 358 ierr = DMPlexSetLabelValue(dm, "marker", v, 1);CHKERRQ(ierr); 359 } 360 { /* Side 0 (Front) */ 361 PetscInt cone[4]; 362 cone[0] = numFaces+4; cone[1] = numFaces+5; cone[2] = numFaces+7; cone[3] = numFaces+6; 363 ierr = DMPlexSetCone(dm, 0, cone);CHKERRQ(ierr); 364 } 365 { /* Side 1 (Back) */ 366 PetscInt cone[4]; 367 cone[0] = numFaces+1; cone[1] = numFaces+0; cone[2] = numFaces+2; cone[3] = numFaces+3; 368 ierr = DMPlexSetCone(dm, 1, cone);CHKERRQ(ierr); 369 } 370 { /* Side 2 (Bottom) */ 371 PetscInt cone[4]; 372 cone[0] = numFaces+0; cone[1] = numFaces+1; cone[2] = numFaces+5; cone[3] = numFaces+4; 373 ierr = DMPlexSetCone(dm, 2, cone);CHKERRQ(ierr); 374 } 375 { /* Side 3 (Top) */ 376 PetscInt cone[4]; 377 cone[0] = numFaces+6; cone[1] = numFaces+7; cone[2] = numFaces+3; cone[3] = numFaces+2; 378 ierr = DMPlexSetCone(dm, 3, cone);CHKERRQ(ierr); 379 } 380 { /* Side 4 (Left) */ 381 PetscInt cone[4]; 382 cone[0] = numFaces+0; cone[1] = numFaces+4; cone[2] = numFaces+6; cone[3] = numFaces+2; 383 ierr = DMPlexSetCone(dm, 4, cone);CHKERRQ(ierr); 384 } 385 { /* Side 5 (Right) */ 386 PetscInt cone[4]; 387 cone[0] = numFaces+5; cone[1] = numFaces+1; cone[2] = numFaces+3; cone[3] = numFaces+7; 388 ierr = DMPlexSetCone(dm, 5, cone);CHKERRQ(ierr); 389 } 390 } 391 ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr); 392 ierr = DMPlexStratify(dm);CHKERRQ(ierr); 393 /* Build coordinates */ 394 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 395 ierr = PetscSectionSetChart(coordSection, numFaces, numFaces + numVertices);CHKERRQ(ierr); 396 for (v = numFaces; v < numFaces+numVertices; ++v) { 397 ierr = PetscSectionSetDof(coordSection, v, 3);CHKERRQ(ierr); 398 } 399 ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 400 ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 401 ierr = VecCreate(PetscObjectComm((PetscObject)dm), &coordinates);CHKERRQ(ierr); 402 ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr); 403 ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 404 ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr); 405 ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 406 for (vz = 0; vz <= faces[2]; ++vz) { 407 for (vy = 0; vy <= faces[1]; ++vy) { 408 for (vx = 0; vx <= faces[0]; ++vx) { 409 coords[((vz*(faces[1]+1)+vy)*(faces[0]+1)+vx)*3+0] = lower[0] + ((upper[0] - lower[0])/faces[0])*vx; 410 coords[((vz*(faces[1]+1)+vy)*(faces[0]+1)+vx)*3+1] = lower[1] + ((upper[1] - lower[1])/faces[1])*vy; 411 coords[((vz*(faces[1]+1)+vy)*(faces[0]+1)+vx)*3+2] = lower[2] + ((upper[2] - lower[2])/faces[2])*vz; 412 } 413 } 414 } 415 ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 416 ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr); 417 ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 418 PetscFunctionReturn(0); 419 } 420 421 #undef __FUNCT__ 422 #define __FUNCT__ "DMPlexCreateSquareMesh" 423 /*@ 424 DMPlexCreateSquareMesh - Creates a 2D mesh for a square lattice. 425 426 Collective on MPI_Comm 427 428 Input Parameters: 429 + comm - The communicator for the DM object 430 . lower - The lower left corner coordinates 431 . upper - The upper right corner coordinates 432 . edges - The number of cells in each direction 433 . bdX - The boundary type for the X direction 434 - bdY - The boundary type for the Y direction 435 436 Output Parameter: 437 . dm - The DM object 438 439 Note: Here is the numbering returned for 2 cells in each direction: 440 $ 22--8-23--9--24 441 $ | | | 442 $ 13 2 14 3 15 443 $ | | | 444 $ 19--6-20--7--21 445 $ | | | 446 $ 10 0 11 1 12 447 $ | | | 448 $ 16--4-17--5--18 449 450 Level: beginner 451 452 .keywords: DM, create 453 .seealso: DMPlexCreateBoxMesh(), DMPlexCreateSquareBoundary(), DMPlexCreateCubeBoundary(), DMSetType(), DMCreate() 454 @*/ 455 PetscErrorCode DMPlexCreateSquareMesh(DM dm, const PetscReal lower[], const PetscReal upper[], const PetscInt edges[], DMBoundaryType bdX, DMBoundaryType bdY) 456 { 457 PetscInt markerTop = 1; 458 PetscInt markerBottom = 1; 459 PetscInt markerRight = 1; 460 PetscInt markerLeft = 1; 461 PetscBool markerSeparate = PETSC_FALSE; 462 PetscMPIInt rank; 463 PetscErrorCode ierr; 464 465 PetscFunctionBegin; 466 ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);CHKERRQ(ierr); 467 ierr = PetscOptionsGetBool(((PetscObject) dm)->prefix, "-dm_plex_separate_marker", &markerSeparate, NULL);CHKERRQ(ierr); 468 if (markerSeparate) { 469 markerTop = 3; 470 markerBottom = 1; 471 markerRight = 2; 472 markerLeft = 4; 473 } 474 { 475 const PetscInt numXEdges = !rank ? edges[0] : 0; 476 const PetscInt numYEdges = !rank ? edges[1] : 0; 477 const PetscInt numXVertices = !rank ? (bdX == DM_BOUNDARY_PERIODIC || bdX == DM_BOUNDARY_TWIST ? edges[0] : edges[0]+1) : 0; 478 const PetscInt numYVertices = !rank ? (bdY == DM_BOUNDARY_PERIODIC || bdY == DM_BOUNDARY_TWIST ? edges[1] : edges[1]+1) : 0; 479 const PetscInt numTotXEdges = numXEdges*numYVertices; 480 const PetscInt numTotYEdges = numYEdges*numXVertices; 481 const PetscInt numVertices = numXVertices*numYVertices; 482 const PetscInt numEdges = numTotXEdges + numTotYEdges; 483 const PetscInt numFaces = numXEdges*numYEdges; 484 const PetscInt firstVertex = numFaces; 485 const PetscInt firstXEdge = numFaces + numVertices; 486 const PetscInt firstYEdge = numFaces + numVertices + numTotXEdges; 487 Vec coordinates; 488 PetscSection coordSection; 489 PetscScalar *coords; 490 PetscInt coordSize; 491 PetscInt v, vx, vy; 492 PetscInt f, fx, fy, e, ex, ey; 493 494 ierr = DMPlexSetChart(dm, 0, numFaces+numEdges+numVertices);CHKERRQ(ierr); 495 for (f = 0; f < numFaces; ++f) { 496 ierr = DMPlexSetConeSize(dm, f, 4);CHKERRQ(ierr); 497 } 498 for (e = firstXEdge; e < firstXEdge+numEdges; ++e) { 499 ierr = DMPlexSetConeSize(dm, e, 2);CHKERRQ(ierr); 500 } 501 ierr = DMSetUp(dm);CHKERRQ(ierr); /* Allocate space for cones */ 502 /* Build faces */ 503 for (fy = 0; fy < numYEdges; fy++) { 504 for (fx = 0; fx < numXEdges; fx++) { 505 const PetscInt face = fy*numXEdges + fx; 506 const PetscInt edgeL = firstYEdge + fx *numYEdges + fy; 507 const PetscInt edgeR = firstYEdge + ((fx+1)%numXVertices)*numYEdges + fy; 508 const PetscInt edgeB = firstXEdge + fy *numXEdges + fx; 509 const PetscInt edgeT = firstXEdge + ((fy+1)%numYVertices)*numXEdges + fx; 510 const PetscInt ornt[4] = {0, 0, -2, -2}; 511 PetscInt cone[4]; 512 513 cone[0] = edgeB; cone[1] = edgeR; cone[2] = edgeT; cone[3] = edgeL; 514 ierr = DMPlexSetCone(dm, face, cone);CHKERRQ(ierr); 515 ierr = DMPlexSetConeOrientation(dm, face, ornt);CHKERRQ(ierr); 516 } 517 } 518 /* Build Y edges*/ 519 for (vx = 0; vx < numXVertices; vx++) { 520 for (ey = 0; ey < numYEdges; ey++) { 521 const PetscInt nextv = (bdY == DM_BOUNDARY_TWIST && ey == numYEdges-1) ? (numXVertices-vx-1) : ((ey+1)%numYVertices)*numXVertices + vx; 522 const PetscInt edge = firstYEdge + vx*numYEdges + ey; 523 const PetscInt vertexB = firstVertex + ey*numXVertices + vx; 524 const PetscInt vertexT = firstVertex + nextv; 525 PetscInt cone[2]; 526 527 cone[0] = vertexB; cone[1] = vertexT; 528 ierr = DMPlexSetCone(dm, edge, cone);CHKERRQ(ierr); 529 if ((bdX != DM_BOUNDARY_PERIODIC) && (bdX != DM_BOUNDARY_TWIST)) { 530 if (vx == numXVertices-1) { 531 ierr = DMPlexSetLabelValue(dm, "marker", edge, markerRight);CHKERRQ(ierr); 532 ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerRight);CHKERRQ(ierr); 533 if (ey == numYEdges-1) { 534 ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerRight);CHKERRQ(ierr); 535 } 536 } else if (vx == 0) { 537 ierr = DMPlexSetLabelValue(dm, "marker", edge, markerLeft);CHKERRQ(ierr); 538 ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerLeft);CHKERRQ(ierr); 539 if (ey == numYEdges-1) { 540 ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerLeft);CHKERRQ(ierr); 541 } 542 } 543 } 544 } 545 } 546 /* Build X edges*/ 547 for (vy = 0; vy < numYVertices; vy++) { 548 for (ex = 0; ex < numXEdges; ex++) { 549 const PetscInt nextv = (bdX == DM_BOUNDARY_TWIST && ex == numXEdges-1) ? (numYVertices-vy-1)*numXVertices : vy*numXVertices + (ex+1)%numXVertices; 550 const PetscInt edge = firstXEdge + vy*numXEdges + ex; 551 const PetscInt vertexL = firstVertex + vy*numXVertices + ex; 552 const PetscInt vertexR = firstVertex + nextv; 553 PetscInt cone[2]; 554 555 cone[0] = vertexL; cone[1] = vertexR; 556 ierr = DMPlexSetCone(dm, edge, cone);CHKERRQ(ierr); 557 if ((bdY != DM_BOUNDARY_PERIODIC) && (bdY != DM_BOUNDARY_TWIST)) { 558 if (vy == numYVertices-1) { 559 ierr = DMPlexSetLabelValue(dm, "marker", edge, markerTop);CHKERRQ(ierr); 560 ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerTop);CHKERRQ(ierr); 561 if (ex == numXEdges-1) { 562 ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerTop);CHKERRQ(ierr); 563 } 564 } else if (vy == 0) { 565 ierr = DMPlexSetLabelValue(dm, "marker", edge, markerBottom);CHKERRQ(ierr); 566 ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerBottom);CHKERRQ(ierr); 567 if (ex == numXEdges-1) { 568 ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerBottom);CHKERRQ(ierr); 569 } 570 } 571 } 572 } 573 } 574 ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr); 575 ierr = DMPlexStratify(dm);CHKERRQ(ierr); 576 /* Build coordinates */ 577 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 578 ierr = PetscSectionSetChart(coordSection, firstVertex, firstVertex+numVertices);CHKERRQ(ierr); 579 for (v = firstVertex; v < firstVertex+numVertices; ++v) { 580 ierr = PetscSectionSetDof(coordSection, v, 2);CHKERRQ(ierr); 581 } 582 ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 583 ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 584 ierr = VecCreate(PetscObjectComm((PetscObject)dm), &coordinates);CHKERRQ(ierr); 585 ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 586 ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr); 587 ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 588 for (vy = 0; vy < numYVertices; ++vy) { 589 for (vx = 0; vx < numXVertices; ++vx) { 590 coords[(vy*numXVertices+vx)*2+0] = lower[0] + ((upper[0] - lower[0])/numXEdges)*vx; 591 coords[(vy*numXVertices+vx)*2+1] = lower[1] + ((upper[1] - lower[1])/numYEdges)*vy; 592 } 593 } 594 ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 595 ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr); 596 ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 597 } 598 PetscFunctionReturn(0); 599 } 600 601 #undef __FUNCT__ 602 #define __FUNCT__ "DMPlexCreateBoxMesh" 603 /*@ 604 DMPlexCreateBoxMesh - Creates a mesh on the tensor product of unit intervals (box) using simplices. 605 606 Collective on MPI_Comm 607 608 Input Parameters: 609 + comm - The communicator for the DM object 610 . dim - The spatial dimension 611 - interpolate - Flag to create intermediate mesh pieces (edges, faces) 612 613 Output Parameter: 614 . dm - The DM object 615 616 Level: beginner 617 618 .keywords: DM, create 619 .seealso: DMPlexCreateHexBoxMesh(), DMSetType(), DMCreate() 620 @*/ 621 PetscErrorCode DMPlexCreateBoxMesh(MPI_Comm comm, PetscInt dim, PetscBool interpolate, DM *dm) 622 { 623 DM boundary; 624 PetscErrorCode ierr; 625 626 PetscFunctionBegin; 627 PetscValidPointer(dm, 4); 628 ierr = DMCreate(comm, &boundary);CHKERRQ(ierr); 629 PetscValidLogicalCollectiveInt(boundary,dim,2); 630 ierr = DMSetType(boundary, DMPLEX);CHKERRQ(ierr); 631 ierr = DMPlexSetDimension(boundary, dim-1);CHKERRQ(ierr); 632 switch (dim) { 633 case 2: 634 { 635 PetscReal lower[2] = {0.0, 0.0}; 636 PetscReal upper[2] = {1.0, 1.0}; 637 PetscInt edges[2] = {2, 2}; 638 639 ierr = DMPlexCreateSquareBoundary(boundary, lower, upper, edges);CHKERRQ(ierr); 640 break; 641 } 642 case 3: 643 { 644 PetscReal lower[3] = {0.0, 0.0, 0.0}; 645 PetscReal upper[3] = {1.0, 1.0, 1.0}; 646 PetscInt faces[3] = {1, 1, 1}; 647 648 ierr = DMPlexCreateCubeBoundary(boundary, lower, upper, faces);CHKERRQ(ierr); 649 break; 650 } 651 default: 652 SETERRQ1(comm, PETSC_ERR_SUP, "Dimension not supported: %d", dim); 653 } 654 ierr = DMPlexGenerate(boundary, NULL, interpolate, dm);CHKERRQ(ierr); 655 ierr = DMDestroy(&boundary);CHKERRQ(ierr); 656 PetscFunctionReturn(0); 657 } 658 659 #undef __FUNCT__ 660 #define __FUNCT__ "DMPlexCreateHexBoxMesh" 661 /*@ 662 DMPlexCreateHexBoxMesh - Creates a mesh on the tensor product of unit intervals (box) using hexahedra. 663 664 Collective on MPI_Comm 665 666 Input Parameters: 667 + comm - The communicator for the DM object 668 . dim - The spatial dimension 669 . periodicX - The boundary type for the X direction 670 . periodicY - The boundary type for the Y direction 671 . periodicZ - The boundary type for the Z direction 672 - cells - The number of cells in each direction 673 674 Output Parameter: 675 . dm - The DM object 676 677 Level: beginner 678 679 .keywords: DM, create 680 .seealso: DMPlexCreateBoxMesh(), DMSetType(), DMCreate() 681 @*/ 682 PetscErrorCode DMPlexCreateHexBoxMesh(MPI_Comm comm, PetscInt dim, const PetscInt cells[], DMBoundaryType periodicX, DMBoundaryType periodicY, DMBoundaryType periodicZ, DM *dm) 683 { 684 PetscErrorCode ierr; 685 686 PetscFunctionBegin; 687 PetscValidPointer(dm, 4); 688 ierr = DMCreate(comm, dm);CHKERRQ(ierr); 689 PetscValidLogicalCollectiveInt(*dm,dim,2); 690 ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 691 ierr = DMPlexSetDimension(*dm, dim);CHKERRQ(ierr); 692 switch (dim) { 693 case 2: 694 { 695 PetscReal lower[2] = {0.0, 0.0}; 696 PetscReal upper[2] = {1.0, 1.0}; 697 698 ierr = DMPlexCreateSquareMesh(*dm, lower, upper, cells, periodicX, periodicY);CHKERRQ(ierr); 699 break; 700 } 701 #if 0 702 case 3: 703 { 704 PetscReal lower[3] = {0.0, 0.0, 0.0}; 705 PetscReal upper[3] = {1.0, 1.0, 1.0}; 706 707 ierr = DMPlexCreateCubeMesh(boundary, lower, upper, cells, periodicX, periodicY, periodicZ);CHKERRQ(ierr); 708 break; 709 } 710 #endif 711 default: 712 SETERRQ1(comm, PETSC_ERR_SUP, "Dimension not supported: %d", dim); 713 } 714 PetscFunctionReturn(0); 715 } 716 717 /* External function declarations here */ 718 extern PetscErrorCode DMCreateInterpolation_Plex(DM dmCoarse, DM dmFine, Mat *interpolation, Vec *scaling); 719 extern PetscErrorCode DMCreateMatrix_Plex(DM dm, Mat *J); 720 extern PetscErrorCode DMCreateCoordinateDM_Plex(DM dm, DM *cdm); 721 extern PetscErrorCode DMRefine_Plex(DM dm, MPI_Comm comm, DM *dmRefined); 722 extern PetscErrorCode DMClone_Plex(DM dm, DM *newdm); 723 extern PetscErrorCode DMSetUp_Plex(DM dm); 724 extern PetscErrorCode DMDestroy_Plex(DM dm); 725 extern PetscErrorCode DMView_Plex(DM dm, PetscViewer viewer); 726 extern PetscErrorCode DMCreateSubDM_Plex(DM dm, PetscInt numFields, PetscInt fields[], IS *is, DM *subdm); 727 extern PetscErrorCode DMLocatePoints_Plex(DM dm, Vec v, IS *cellIS); 728 729 #undef __FUNCT__ 730 #define __FUNCT__ "DMCreateGlobalVector_Plex" 731 static PetscErrorCode DMCreateGlobalVector_Plex(DM dm,Vec *vec) 732 { 733 PetscErrorCode ierr; 734 735 PetscFunctionBegin; 736 ierr = DMCreateGlobalVector_Section_Private(dm,vec);CHKERRQ(ierr); 737 /* ierr = VecSetOperation(*vec, VECOP_DUPLICATE, (void(*)(void)) VecDuplicate_MPI_DM);CHKERRQ(ierr); */ 738 ierr = VecSetOperation(*vec, VECOP_VIEW, (void (*)(void))VecView_Plex);CHKERRQ(ierr); 739 PetscFunctionReturn(0); 740 } 741 742 #undef __FUNCT__ 743 #define __FUNCT__ "DMCreateLocalVector_Plex" 744 static PetscErrorCode DMCreateLocalVector_Plex(DM dm,Vec *vec) 745 { 746 PetscErrorCode ierr; 747 748 PetscFunctionBegin; 749 ierr = DMCreateLocalVector_Section_Private(dm,vec);CHKERRQ(ierr); 750 ierr = VecSetOperation(*vec, VECOP_VIEW, (void(*)(void)) VecView_Plex_Local);CHKERRQ(ierr); 751 PetscFunctionReturn(0); 752 } 753 754 #undef __FUNCT__ 755 #define __FUNCT__ "DMInitialize_Plex" 756 PetscErrorCode DMInitialize_Plex(DM dm) 757 { 758 PetscFunctionBegin; 759 dm->ops->view = DMView_Plex; 760 dm->ops->setfromoptions = DMSetFromOptions_Plex; 761 dm->ops->clone = DMClone_Plex; 762 dm->ops->setup = DMSetUp_Plex; 763 dm->ops->createglobalvector = DMCreateGlobalVector_Plex; 764 dm->ops->createlocalvector = DMCreateLocalVector_Plex; 765 dm->ops->getlocaltoglobalmapping = NULL; 766 dm->ops->getlocaltoglobalmappingblock = NULL; 767 dm->ops->createfieldis = NULL; 768 dm->ops->createcoordinatedm = DMCreateCoordinateDM_Plex; 769 dm->ops->getcoloring = 0; 770 dm->ops->creatematrix = DMCreateMatrix_Plex; 771 dm->ops->createinterpolation = 0; 772 dm->ops->getaggregates = 0; 773 dm->ops->getinjection = 0; 774 dm->ops->refine = DMRefine_Plex; 775 dm->ops->coarsen = 0; 776 dm->ops->refinehierarchy = 0; 777 dm->ops->coarsenhierarchy = 0; 778 dm->ops->globaltolocalbegin = NULL; 779 dm->ops->globaltolocalend = NULL; 780 dm->ops->localtoglobalbegin = NULL; 781 dm->ops->localtoglobalend = NULL; 782 dm->ops->destroy = DMDestroy_Plex; 783 dm->ops->createsubdm = DMCreateSubDM_Plex; 784 dm->ops->locatepoints = DMLocatePoints_Plex; 785 PetscFunctionReturn(0); 786 } 787 788 #undef __FUNCT__ 789 #define __FUNCT__ "DMClone_Plex" 790 PetscErrorCode DMClone_Plex(DM dm, DM *newdm) 791 { 792 DM_Plex *mesh = (DM_Plex *) dm->data; 793 PetscErrorCode ierr; 794 795 PetscFunctionBegin; 796 mesh->refct++; 797 (*newdm)->data = mesh; 798 ierr = PetscObjectChangeTypeName((PetscObject) *newdm, DMPLEX);CHKERRQ(ierr); 799 ierr = DMInitialize_Plex(*newdm);CHKERRQ(ierr); 800 PetscFunctionReturn(0); 801 } 802 803 /*MC 804 DMPLEX = "plex" - A DM object that encapsulates an unstructured mesh, or CW Complex, which can be expressed using a Hasse Diagram. 805 In the local representation, Vecs contain all unknowns in the interior and shared boundary. This is 806 specified by a PetscSection object. Ownership in the global representation is determined by 807 ownership of the underlying DMPlex points. This is specified by another PetscSection object. 808 809 Level: intermediate 810 811 .seealso: DMType, DMPlexCreate(), DMCreate(), DMSetType() 812 M*/ 813 814 #undef __FUNCT__ 815 #define __FUNCT__ "DMCreate_Plex" 816 PETSC_EXTERN PetscErrorCode DMCreate_Plex(DM dm) 817 { 818 DM_Plex *mesh; 819 PetscInt unit, d; 820 PetscErrorCode ierr; 821 822 PetscFunctionBegin; 823 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 824 ierr = PetscNewLog(dm,&mesh);CHKERRQ(ierr); 825 dm->data = mesh; 826 827 mesh->refct = 1; 828 mesh->dim = 0; 829 ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm), &mesh->coneSection);CHKERRQ(ierr); 830 mesh->maxConeSize = 0; 831 mesh->cones = NULL; 832 mesh->coneOrientations = NULL; 833 ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm), &mesh->supportSection);CHKERRQ(ierr); 834 mesh->maxSupportSize = 0; 835 mesh->supports = NULL; 836 mesh->refinementUniform = PETSC_TRUE; 837 mesh->refinementLimit = -1.0; 838 839 mesh->facesTmp = NULL; 840 841 mesh->subpointMap = NULL; 842 843 for (unit = 0; unit < NUM_PETSC_UNITS; ++unit) mesh->scale[unit] = 1.0; 844 845 mesh->labels = NULL; 846 mesh->depthLabel = NULL; 847 mesh->globalVertexNumbers = NULL; 848 mesh->globalCellNumbers = NULL; 849 for (d = 0; d < 8; ++d) mesh->hybridPointMax[d] = PETSC_DETERMINE; 850 mesh->vtkCellHeight = 0; 851 mesh->preallocCenterDim = -1; 852 853 mesh->printSetValues = PETSC_FALSE; 854 mesh->printFEM = 0; 855 mesh->printTol = 1.0e-10; 856 857 ierr = DMInitialize_Plex(dm);CHKERRQ(ierr); 858 PetscFunctionReturn(0); 859 } 860 861 #undef __FUNCT__ 862 #define __FUNCT__ "DMPlexCreate" 863 /*@ 864 DMPlexCreate - Creates a DMPlex object, which encapsulates an unstructured mesh, or CW complex, which can be expressed using a Hasse Diagram. 865 866 Collective on MPI_Comm 867 868 Input Parameter: 869 . comm - The communicator for the DMPlex object 870 871 Output Parameter: 872 . mesh - The DMPlex object 873 874 Level: beginner 875 876 .keywords: DMPlex, create 877 @*/ 878 PetscErrorCode DMPlexCreate(MPI_Comm comm, DM *mesh) 879 { 880 PetscErrorCode ierr; 881 882 PetscFunctionBegin; 883 PetscValidPointer(mesh,2); 884 ierr = DMCreate(comm, mesh);CHKERRQ(ierr); 885 ierr = DMSetType(*mesh, DMPLEX);CHKERRQ(ierr); 886 PetscFunctionReturn(0); 887 } 888 889 #undef __FUNCT__ 890 #define __FUNCT__ "DMPlexBuildFromCellList_Private" 891 /* 892 This takes as input the common mesh generator output, a list of the vertices for each cell 893 */ 894 PetscErrorCode DMPlexBuildFromCellList_Private(DM dm, PetscInt numCells, PetscInt numVertices, PetscInt numCorners, const int cells[]) 895 { 896 PetscInt *cone, c, p; 897 PetscErrorCode ierr; 898 899 PetscFunctionBegin; 900 ierr = DMPlexSetChart(dm, 0, numCells+numVertices);CHKERRQ(ierr); 901 for (c = 0; c < numCells; ++c) { 902 ierr = DMPlexSetConeSize(dm, c, numCorners);CHKERRQ(ierr); 903 } 904 ierr = DMSetUp(dm);CHKERRQ(ierr); 905 ierr = DMGetWorkArray(dm, numCorners, PETSC_INT, &cone);CHKERRQ(ierr); 906 for (c = 0; c < numCells; ++c) { 907 for (p = 0; p < numCorners; ++p) { 908 cone[p] = cells[c*numCorners+p]+numCells; 909 } 910 ierr = DMPlexSetCone(dm, c, cone);CHKERRQ(ierr); 911 } 912 ierr = DMRestoreWorkArray(dm, numCorners, PETSC_INT, &cone);CHKERRQ(ierr); 913 ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr); 914 ierr = DMPlexStratify(dm);CHKERRQ(ierr); 915 PetscFunctionReturn(0); 916 } 917 918 #undef __FUNCT__ 919 #define __FUNCT__ "DMPlexBuildCoordinates_Private" 920 /* 921 This takes as input the coordinates for each vertex 922 */ 923 PetscErrorCode DMPlexBuildCoordinates_Private(DM dm, PetscInt spaceDim, PetscInt numCells, PetscInt numVertices, const double vertexCoords[]) 924 { 925 PetscSection coordSection; 926 Vec coordinates; 927 PetscScalar *coords; 928 PetscInt coordSize, v, d; 929 PetscErrorCode ierr; 930 931 PetscFunctionBegin; 932 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 933 ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr); 934 ierr = PetscSectionSetFieldComponents(coordSection, 0, spaceDim);CHKERRQ(ierr); 935 ierr = PetscSectionSetChart(coordSection, numCells, numCells + numVertices);CHKERRQ(ierr); 936 for (v = numCells; v < numCells+numVertices; ++v) { 937 ierr = PetscSectionSetDof(coordSection, v, spaceDim);CHKERRQ(ierr); 938 ierr = PetscSectionSetFieldDof(coordSection, v, 0, spaceDim);CHKERRQ(ierr); 939 } 940 ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 941 ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 942 ierr = VecCreate(PetscObjectComm((PetscObject)dm), &coordinates);CHKERRQ(ierr); 943 ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr); 944 ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 945 ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr); 946 ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 947 for (v = 0; v < numVertices; ++v) { 948 for (d = 0; d < spaceDim; ++d) { 949 coords[v*spaceDim+d] = vertexCoords[v*spaceDim+d]; 950 } 951 } 952 ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 953 ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr); 954 ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 955 PetscFunctionReturn(0); 956 } 957 958 #undef __FUNCT__ 959 #define __FUNCT__ "DMPlexCreateFromCellList" 960 /*@C 961 DMPlexCreateFromCellList - This takes as input common mesh generator output, a list of the vertices for each cell, and produces a DM 962 963 Input Parameters: 964 + comm - The communicator 965 . dim - The topological dimension of the mesh 966 . numCells - The number of cells 967 . numVertices - The number of vertices 968 . numCorners - The number of vertices for each cell 969 . interpolate - Flag indicating that intermediate mesh entities (faces, edges) should be created automatically 970 . cells - An array of numCells*numCorners numbers, the vertices for each cell 971 . spaceDim - The spatial dimension used for coordinates 972 - vertexCoords - An array of numVertices*spaceDim numbers, the coordinates of each vertex 973 974 Output Parameter: 975 . dm - The DM 976 977 Note: Two triangles sharing a face 978 $ 979 $ 2 980 $ / | \ 981 $ / | \ 982 $ / | \ 983 $ 0 0 | 1 3 984 $ \ | / 985 $ \ | / 986 $ \ | / 987 $ 1 988 would have input 989 $ numCells = 2, numVertices = 4 990 $ cells = [0 1 2 1 3 2] 991 $ 992 which would result in the DMPlex 993 $ 994 $ 4 995 $ / | \ 996 $ / | \ 997 $ / | \ 998 $ 2 0 | 1 5 999 $ \ | / 1000 $ \ | / 1001 $ \ | / 1002 $ 3 1003 1004 Level: beginner 1005 1006 .seealso: DMPlexCreateFromDAG(), DMPlexCreate() 1007 @*/ 1008 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) 1009 { 1010 PetscErrorCode ierr; 1011 1012 PetscFunctionBegin; 1013 ierr = DMCreate(comm, dm);CHKERRQ(ierr); 1014 ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 1015 ierr = DMPlexSetDimension(*dm, dim);CHKERRQ(ierr); 1016 ierr = DMPlexBuildFromCellList_Private(*dm, numCells, numVertices, numCorners, cells);CHKERRQ(ierr); 1017 if (interpolate) { 1018 DM idm; 1019 1020 ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr); 1021 ierr = DMDestroy(dm);CHKERRQ(ierr); 1022 *dm = idm; 1023 } 1024 ierr = DMPlexBuildCoordinates_Private(*dm, spaceDim, numCells, numVertices, vertexCoords);CHKERRQ(ierr); 1025 PetscFunctionReturn(0); 1026 } 1027 1028 #undef __FUNCT__ 1029 #define __FUNCT__ "DMPlexCreateFromDAG" 1030 /*@ 1031 DMPlexCreateFromDAG - This takes as input the adjacency-list representation of the Directed Acyclic Graph (Hasse Diagram) encoding a mesh, and produces a DM 1032 1033 Input Parameters: 1034 + dm - The empty DM object, usually from DMCreate() and DMPlexSetDimension() 1035 . depth - The depth of the DAG 1036 . numPoints - The number of points at each depth 1037 . coneSize - The cone size of each point 1038 . cones - The concatenation of the cone points for each point, the cone list must be oriented correctly for each point 1039 . coneOrientations - The orientation of each cone point 1040 - vertexCoords - An array of numVertices*dim numbers, the coordinates of each vertex 1041 1042 Output Parameter: 1043 . dm - The DM 1044 1045 Note: Two triangles sharing a face would have input 1046 $ depth = 1, numPoints = [4 2], coneSize = [3 3 0 0 0 0] 1047 $ cones = [2 3 4 3 5 4], coneOrientations = [0 0 0 0 0 0] 1048 $ vertexCoords = [-1.0 0.0 0.0 -1.0 0.0 1.0 1.0 0.0] 1049 $ 1050 which would result in the DMPlex 1051 $ 1052 $ 4 1053 $ / | \ 1054 $ / | \ 1055 $ / | \ 1056 $ 2 0 | 1 5 1057 $ \ | / 1058 $ \ | / 1059 $ \ | / 1060 $ 3 1061 $ 1062 $ Notice that all points are numbered consecutively, unlikely DMPlexCreateFromCellList() 1063 1064 Level: advanced 1065 1066 .seealso: DMPlexCreateFromCellList(), DMPlexCreate() 1067 @*/ 1068 PetscErrorCode DMPlexCreateFromDAG(DM dm, PetscInt depth, const PetscInt numPoints[], const PetscInt coneSize[], const PetscInt cones[], const PetscInt coneOrientations[], const PetscScalar vertexCoords[]) 1069 { 1070 Vec coordinates; 1071 PetscSection coordSection; 1072 PetscScalar *coords; 1073 PetscInt coordSize, firstVertex = numPoints[depth], pStart = 0, pEnd = 0, p, v, dim, d, off; 1074 PetscErrorCode ierr; 1075 1076 PetscFunctionBegin; 1077 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 1078 for (d = 0; d <= depth; ++d) pEnd += numPoints[d]; 1079 ierr = DMPlexSetChart(dm, pStart, pEnd);CHKERRQ(ierr); 1080 for (p = pStart; p < pEnd; ++p) { 1081 ierr = DMPlexSetConeSize(dm, p, coneSize[p-pStart]);CHKERRQ(ierr); 1082 } 1083 ierr = DMSetUp(dm);CHKERRQ(ierr); /* Allocate space for cones */ 1084 for (p = pStart, off = 0; p < pEnd; off += coneSize[p-pStart], ++p) { 1085 ierr = DMPlexSetCone(dm, p, &cones[off]);CHKERRQ(ierr); 1086 ierr = DMPlexSetConeOrientation(dm, p, &coneOrientations[off]);CHKERRQ(ierr); 1087 } 1088 ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr); 1089 ierr = DMPlexStratify(dm);CHKERRQ(ierr); 1090 /* Build coordinates */ 1091 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 1092 ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr); 1093 ierr = PetscSectionSetFieldComponents(coordSection, 0, dim);CHKERRQ(ierr); 1094 ierr = PetscSectionSetChart(coordSection, firstVertex, firstVertex+numPoints[0]);CHKERRQ(ierr); 1095 for (v = firstVertex; v < firstVertex+numPoints[0]; ++v) { 1096 ierr = PetscSectionSetDof(coordSection, v, dim);CHKERRQ(ierr); 1097 ierr = PetscSectionSetFieldDof(coordSection, v, 0, dim);CHKERRQ(ierr); 1098 } 1099 ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 1100 ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 1101 ierr = VecCreate(PetscObjectComm((PetscObject)dm), &coordinates);CHKERRQ(ierr); 1102 ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr); 1103 ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 1104 ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr); 1105 ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 1106 for (v = 0; v < numPoints[0]; ++v) { 1107 PetscInt off; 1108 1109 ierr = PetscSectionGetOffset(coordSection, v+firstVertex, &off);CHKERRQ(ierr); 1110 for (d = 0; d < dim; ++d) { 1111 coords[off+d] = vertexCoords[v*dim+d]; 1112 } 1113 } 1114 ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 1115 ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr); 1116 ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 1117 PetscFunctionReturn(0); 1118 } 1119