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 434 Output Parameter: 435 . dm - The DM object 436 437 Note: Here is the numbering returned for 2 cells in each direction: 438 $ 22--8-23--9--24 439 $ | | | 440 $ 13 2 14 3 15 441 $ | | | 442 $ 19--6-20--7--21 443 $ | | | 444 $ 10 0 11 1 12 445 $ | | | 446 $ 16--4-17--5--18 447 448 Level: beginner 449 450 .keywords: DM, create 451 .seealso: DMPlexCreateBoxMesh(), DMPlexCreateSquareBoundary(), DMPlexCreateCubeBoundary(), DMSetType(), DMCreate() 452 @*/ 453 PetscErrorCode DMPlexCreateSquareMesh(DM dm, const PetscReal lower[], const PetscReal upper[], const PetscInt edges[]) 454 { 455 PetscInt markerTop = 1; 456 PetscInt markerBottom = 1; 457 PetscInt markerRight = 1; 458 PetscInt markerLeft = 1; 459 PetscBool markerSeparate = PETSC_FALSE; 460 PetscMPIInt rank; 461 PetscErrorCode ierr; 462 463 PetscFunctionBegin; 464 ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);CHKERRQ(ierr); 465 ierr = PetscOptionsGetBool(((PetscObject) dm)->prefix, "-dm_plex_separate_marker", &markerSeparate, NULL);CHKERRQ(ierr); 466 if (markerSeparate) { 467 markerTop = 3; 468 markerBottom = 1; 469 markerRight = 2; 470 markerLeft = 4; 471 } 472 { 473 const PetscInt numXEdges = !rank ? edges[0] : 0; 474 const PetscInt numYEdges = !rank ? edges[1] : 0; 475 const PetscInt numXVertices = !rank ? edges[0]+1 : 0; 476 const PetscInt numYVertices = !rank ? edges[1]+1 : 0; 477 const PetscInt numTotXEdges = numXEdges*numYVertices; 478 const PetscInt numTotYEdges = numYEdges*numXVertices; 479 const PetscInt numVertices = numXVertices*numYVertices; 480 const PetscInt numEdges = numTotXEdges + numTotYEdges; 481 const PetscInt numFaces = numXEdges*numYEdges; 482 const PetscInt firstVertex = numFaces; 483 const PetscInt firstXEdge = numFaces + numVertices; 484 const PetscInt firstYEdge = numFaces + numVertices + numTotXEdges; 485 Vec coordinates; 486 PetscSection coordSection; 487 PetscScalar *coords; 488 PetscInt coordSize; 489 PetscInt v, vx, vy; 490 PetscInt f, fx, fy, e, ex, ey; 491 492 ierr = DMPlexSetChart(dm, 0, numFaces+numEdges+numVertices);CHKERRQ(ierr); 493 for (f = 0; f < numFaces; ++f) { 494 ierr = DMPlexSetConeSize(dm, f, 4);CHKERRQ(ierr); 495 } 496 for (e = firstXEdge; e < firstXEdge+numEdges; ++e) { 497 ierr = DMPlexSetConeSize(dm, e, 2);CHKERRQ(ierr); 498 } 499 ierr = DMSetUp(dm);CHKERRQ(ierr); /* Allocate space for cones */ 500 /* Build faces */ 501 for (fy = 0; fy < numYEdges; fy++) { 502 for (fx = 0; fx < numXEdges; fx++) { 503 const PetscInt face = fy*numXEdges + fx; 504 const PetscInt edgeL = firstYEdge + fx*numYEdges + fy; 505 const PetscInt edgeB = firstXEdge + fy*numXEdges + fx; 506 const PetscInt ornt[4] = {0, 0, -2, -2}; 507 PetscInt cone[4]; 508 509 cone[0] = edgeB; cone[1] = edgeL+numYEdges; cone[2] = edgeB+numXEdges; cone[3] = edgeL; 510 ierr = DMPlexSetCone(dm, face, cone);CHKERRQ(ierr); 511 ierr = DMPlexSetConeOrientation(dm, face, ornt);CHKERRQ(ierr); 512 } 513 } 514 /* Build Y edges*/ 515 for (vx = 0; vx < numXVertices; vx++) { 516 for (ey = 0; ey < numYEdges; ey++) { 517 const PetscInt edge = firstYEdge + vx*numYEdges + ey; 518 const PetscInt vertex = firstVertex + ey*numXVertices + vx; 519 PetscInt cone[2]; 520 521 cone[0] = vertex; cone[1] = vertex+numXVertices; 522 ierr = DMPlexSetCone(dm, edge, cone);CHKERRQ(ierr); 523 if (vx == numXVertices-1) { 524 ierr = DMPlexSetLabelValue(dm, "marker", edge, markerRight);CHKERRQ(ierr); 525 ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerRight);CHKERRQ(ierr); 526 if (ey == numYEdges-1) { 527 ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerRight);CHKERRQ(ierr); 528 } 529 } else if (vx == 0) { 530 ierr = DMPlexSetLabelValue(dm, "marker", edge, markerLeft);CHKERRQ(ierr); 531 ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerLeft);CHKERRQ(ierr); 532 if (ey == numYEdges-1) { 533 ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerLeft);CHKERRQ(ierr); 534 } 535 } 536 } 537 } 538 /* Build X edges*/ 539 for (vy = 0; vy < numYVertices; vy++) { 540 for (ex = 0; ex < numXEdges; ex++) { 541 const PetscInt edge = firstXEdge + vy*numXEdges + ex; 542 const PetscInt vertex = firstVertex + vy*numXVertices + ex; 543 PetscInt cone[2]; 544 545 cone[0] = vertex; cone[1] = vertex+1; 546 ierr = DMPlexSetCone(dm, edge, cone);CHKERRQ(ierr); 547 if (vy == numYVertices-1) { 548 ierr = DMPlexSetLabelValue(dm, "marker", edge, markerTop);CHKERRQ(ierr); 549 ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerTop);CHKERRQ(ierr); 550 if (ex == numXEdges-1) { 551 ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerTop);CHKERRQ(ierr); 552 } 553 } else if (vy == 0) { 554 ierr = DMPlexSetLabelValue(dm, "marker", edge, markerBottom);CHKERRQ(ierr); 555 ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerBottom);CHKERRQ(ierr); 556 if (ex == numXEdges-1) { 557 ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerBottom);CHKERRQ(ierr); 558 } 559 } 560 } 561 } 562 ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr); 563 ierr = DMPlexStratify(dm);CHKERRQ(ierr); 564 /* Build coordinates */ 565 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 566 ierr = PetscSectionSetChart(coordSection, firstVertex, firstVertex+numVertices);CHKERRQ(ierr); 567 for (v = firstVertex; v < firstVertex+numVertices; ++v) { 568 ierr = PetscSectionSetDof(coordSection, v, 2);CHKERRQ(ierr); 569 } 570 ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 571 ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 572 ierr = VecCreate(PetscObjectComm((PetscObject)dm), &coordinates);CHKERRQ(ierr); 573 ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 574 ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr); 575 ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 576 for (vy = 0; vy < numYVertices; ++vy) { 577 for (vx = 0; vx < numXVertices; ++vx) { 578 coords[(vy*numXVertices+vx)*2+0] = lower[0] + ((upper[0] - lower[0])/numXEdges)*vx; 579 coords[(vy*numXVertices+vx)*2+1] = lower[1] + ((upper[1] - lower[1])/numYEdges)*vy; 580 } 581 } 582 ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 583 ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr); 584 ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 585 } 586 PetscFunctionReturn(0); 587 } 588 589 #undef __FUNCT__ 590 #define __FUNCT__ "DMPlexCreateBoxMesh" 591 /*@ 592 DMPlexCreateBoxMesh - Creates a mesh on the tensor product of unit intervals (box) using simplices. 593 594 Collective on MPI_Comm 595 596 Input Parameters: 597 + comm - The communicator for the DM object 598 . dim - The spatial dimension 599 - interpolate - Flag to create intermediate mesh pieces (edges, faces) 600 601 Output Parameter: 602 . dm - The DM object 603 604 Level: beginner 605 606 .keywords: DM, create 607 .seealso: DMPlexCreateHexBoxMesh(), DMSetType(), DMCreate() 608 @*/ 609 PetscErrorCode DMPlexCreateBoxMesh(MPI_Comm comm, PetscInt dim, PetscBool interpolate, DM *dm) 610 { 611 DM boundary; 612 PetscErrorCode ierr; 613 614 PetscFunctionBegin; 615 PetscValidPointer(dm, 4); 616 ierr = DMCreate(comm, &boundary);CHKERRQ(ierr); 617 PetscValidLogicalCollectiveInt(boundary,dim,2); 618 ierr = DMSetType(boundary, DMPLEX);CHKERRQ(ierr); 619 ierr = DMPlexSetDimension(boundary, dim-1);CHKERRQ(ierr); 620 switch (dim) { 621 case 2: 622 { 623 PetscReal lower[2] = {0.0, 0.0}; 624 PetscReal upper[2] = {1.0, 1.0}; 625 PetscInt edges[2] = {2, 2}; 626 627 ierr = DMPlexCreateSquareBoundary(boundary, lower, upper, edges);CHKERRQ(ierr); 628 break; 629 } 630 case 3: 631 { 632 PetscReal lower[3] = {0.0, 0.0, 0.0}; 633 PetscReal upper[3] = {1.0, 1.0, 1.0}; 634 PetscInt faces[3] = {1, 1, 1}; 635 636 ierr = DMPlexCreateCubeBoundary(boundary, lower, upper, faces);CHKERRQ(ierr); 637 break; 638 } 639 default: 640 SETERRQ1(comm, PETSC_ERR_SUP, "Dimension not supported: %d", dim); 641 } 642 ierr = DMPlexGenerate(boundary, NULL, interpolate, dm);CHKERRQ(ierr); 643 ierr = DMDestroy(&boundary);CHKERRQ(ierr); 644 PetscFunctionReturn(0); 645 } 646 647 #undef __FUNCT__ 648 #define __FUNCT__ "DMPlexCreateHexBoxMesh" 649 /*@ 650 DMPlexCreateHexBoxMesh - Creates a mesh on the tensor product of unit intervals (box) using hexahedra. 651 652 Collective on MPI_Comm 653 654 Input Parameters: 655 + comm - The communicator for the DM object 656 . dim - The spatial dimension 657 - cells - The number of cells in each direction 658 659 Output Parameter: 660 . dm - The DM object 661 662 Level: beginner 663 664 .keywords: DM, create 665 .seealso: DMPlexCreateBoxMesh(), DMSetType(), DMCreate() 666 @*/ 667 PetscErrorCode DMPlexCreateHexBoxMesh(MPI_Comm comm, PetscInt dim, const PetscInt cells[], DM *dm) 668 { 669 PetscErrorCode ierr; 670 671 PetscFunctionBegin; 672 PetscValidPointer(dm, 4); 673 ierr = DMCreate(comm, dm);CHKERRQ(ierr); 674 PetscValidLogicalCollectiveInt(*dm,dim,2); 675 ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 676 ierr = DMPlexSetDimension(*dm, dim);CHKERRQ(ierr); 677 switch (dim) { 678 case 2: 679 { 680 PetscReal lower[2] = {0.0, 0.0}; 681 PetscReal upper[2] = {1.0, 1.0}; 682 683 ierr = DMPlexCreateSquareMesh(*dm, lower, upper, cells);CHKERRQ(ierr); 684 break; 685 } 686 #if 0 687 case 3: 688 { 689 PetscReal lower[3] = {0.0, 0.0, 0.0}; 690 PetscReal upper[3] = {1.0, 1.0, 1.0}; 691 692 ierr = DMPlexCreateCubeMesh(boundary, lower, upper, cells);CHKERRQ(ierr); 693 break; 694 } 695 #endif 696 default: 697 SETERRQ1(comm, PETSC_ERR_SUP, "Dimension not supported: %d", dim); 698 } 699 PetscFunctionReturn(0); 700 } 701 702 /* External function declarations here */ 703 extern PetscErrorCode DMCreateInterpolation_Plex(DM dmCoarse, DM dmFine, Mat *interpolation, Vec *scaling); 704 extern PetscErrorCode DMCreateMatrix_Plex(DM dm, Mat *J); 705 extern PetscErrorCode DMCreateCoordinateDM_Plex(DM dm, DM *cdm); 706 extern PetscErrorCode DMRefine_Plex(DM dm, MPI_Comm comm, DM *dmRefined); 707 extern PetscErrorCode DMClone_Plex(DM dm, DM *newdm); 708 extern PetscErrorCode DMSetUp_Plex(DM dm); 709 extern PetscErrorCode DMDestroy_Plex(DM dm); 710 extern PetscErrorCode DMView_Plex(DM dm, PetscViewer viewer); 711 extern PetscErrorCode DMCreateSubDM_Plex(DM dm, PetscInt numFields, PetscInt fields[], IS *is, DM *subdm); 712 extern PetscErrorCode DMLocatePoints_Plex(DM dm, Vec v, IS *cellIS); 713 714 #undef __FUNCT__ 715 #define __FUNCT__ "DMCreateGlobalVector_Plex" 716 static PetscErrorCode DMCreateGlobalVector_Plex(DM dm,Vec *vec) 717 { 718 PetscErrorCode ierr; 719 720 PetscFunctionBegin; 721 ierr = DMCreateGlobalVector_Section_Private(dm,vec);CHKERRQ(ierr); 722 /* ierr = VecSetOperation(*vec, VECOP_DUPLICATE, (void(*)(void)) VecDuplicate_MPI_DM);CHKERRQ(ierr); */ 723 ierr = VecSetOperation(*vec, VECOP_VIEW, (void (*)(void))VecView_Plex);CHKERRQ(ierr); 724 PetscFunctionReturn(0); 725 } 726 727 #undef __FUNCT__ 728 #define __FUNCT__ "DMCreateLocalVector_Plex" 729 static PetscErrorCode DMCreateLocalVector_Plex(DM dm,Vec *vec) 730 { 731 PetscErrorCode ierr; 732 733 PetscFunctionBegin; 734 ierr = DMCreateLocalVector_Section_Private(dm,vec);CHKERRQ(ierr); 735 ierr = VecSetOperation(*vec, VECOP_VIEW, (void(*)(void)) VecView_Plex_Local);CHKERRQ(ierr); 736 PetscFunctionReturn(0); 737 } 738 739 #undef __FUNCT__ 740 #define __FUNCT__ "DMInitialize_Plex" 741 PetscErrorCode DMInitialize_Plex(DM dm) 742 { 743 PetscFunctionBegin; 744 dm->ops->view = DMView_Plex; 745 dm->ops->setfromoptions = DMSetFromOptions_Plex; 746 dm->ops->clone = DMClone_Plex; 747 dm->ops->setup = DMSetUp_Plex; 748 dm->ops->createglobalvector = DMCreateGlobalVector_Plex; 749 dm->ops->createlocalvector = DMCreateLocalVector_Plex; 750 dm->ops->getlocaltoglobalmapping = NULL; 751 dm->ops->getlocaltoglobalmappingblock = NULL; 752 dm->ops->createfieldis = NULL; 753 dm->ops->createcoordinatedm = DMCreateCoordinateDM_Plex; 754 dm->ops->getcoloring = 0; 755 dm->ops->creatematrix = DMCreateMatrix_Plex; 756 dm->ops->createinterpolation = 0; 757 dm->ops->getaggregates = 0; 758 dm->ops->getinjection = 0; 759 dm->ops->refine = DMRefine_Plex; 760 dm->ops->coarsen = 0; 761 dm->ops->refinehierarchy = 0; 762 dm->ops->coarsenhierarchy = 0; 763 dm->ops->globaltolocalbegin = NULL; 764 dm->ops->globaltolocalend = NULL; 765 dm->ops->localtoglobalbegin = NULL; 766 dm->ops->localtoglobalend = NULL; 767 dm->ops->destroy = DMDestroy_Plex; 768 dm->ops->createsubdm = DMCreateSubDM_Plex; 769 dm->ops->locatepoints = DMLocatePoints_Plex; 770 PetscFunctionReturn(0); 771 } 772 773 #undef __FUNCT__ 774 #define __FUNCT__ "DMClone_Plex" 775 PetscErrorCode DMClone_Plex(DM dm, DM *newdm) 776 { 777 DM_Plex *mesh = (DM_Plex *) dm->data; 778 PetscErrorCode ierr; 779 780 PetscFunctionBegin; 781 mesh->refct++; 782 (*newdm)->data = mesh; 783 ierr = PetscObjectChangeTypeName((PetscObject) *newdm, DMPLEX);CHKERRQ(ierr); 784 ierr = DMInitialize_Plex(*newdm);CHKERRQ(ierr); 785 PetscFunctionReturn(0); 786 } 787 788 /*MC 789 DMPLEX = "plex" - A DM object that encapsulates an unstructured mesh, or CW Complex, which can be expressed using a Hasse Diagram. 790 In the local representation, Vecs contain all unknowns in the interior and shared boundary. This is 791 specified by a PetscSection object. Ownership in the global representation is determined by 792 ownership of the underlying DMPlex points. This is specified by another PetscSection object. 793 794 Level: intermediate 795 796 .seealso: DMType, DMPlexCreate(), DMCreate(), DMSetType() 797 M*/ 798 799 #undef __FUNCT__ 800 #define __FUNCT__ "DMCreate_Plex" 801 PETSC_EXTERN PetscErrorCode DMCreate_Plex(DM dm) 802 { 803 DM_Plex *mesh; 804 PetscInt unit, d; 805 PetscErrorCode ierr; 806 807 PetscFunctionBegin; 808 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 809 ierr = PetscNewLog(dm,&mesh);CHKERRQ(ierr); 810 dm->data = mesh; 811 812 mesh->refct = 1; 813 mesh->dim = 0; 814 ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm), &mesh->coneSection);CHKERRQ(ierr); 815 mesh->maxConeSize = 0; 816 mesh->cones = NULL; 817 mesh->coneOrientations = NULL; 818 ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm), &mesh->supportSection);CHKERRQ(ierr); 819 mesh->maxSupportSize = 0; 820 mesh->supports = NULL; 821 mesh->refinementUniform = PETSC_TRUE; 822 mesh->refinementLimit = -1.0; 823 824 mesh->facesTmp = NULL; 825 826 mesh->subpointMap = NULL; 827 828 for (unit = 0; unit < NUM_PETSC_UNITS; ++unit) mesh->scale[unit] = 1.0; 829 830 mesh->labels = NULL; 831 mesh->depthLabel = NULL; 832 mesh->globalVertexNumbers = NULL; 833 mesh->globalCellNumbers = NULL; 834 for (d = 0; d < 8; ++d) mesh->hybridPointMax[d] = PETSC_DETERMINE; 835 mesh->vtkCellHeight = 0; 836 mesh->preallocCenterDim = -1; 837 838 mesh->printSetValues = PETSC_FALSE; 839 mesh->printFEM = 0; 840 mesh->printTol = 1.0e-10; 841 842 ierr = DMInitialize_Plex(dm);CHKERRQ(ierr); 843 PetscFunctionReturn(0); 844 } 845 846 #undef __FUNCT__ 847 #define __FUNCT__ "DMPlexCreate" 848 /*@ 849 DMPlexCreate - Creates a DMPlex object, which encapsulates an unstructured mesh, or CW complex, which can be expressed using a Hasse Diagram. 850 851 Collective on MPI_Comm 852 853 Input Parameter: 854 . comm - The communicator for the DMPlex object 855 856 Output Parameter: 857 . mesh - The DMPlex object 858 859 Level: beginner 860 861 .keywords: DMPlex, create 862 @*/ 863 PetscErrorCode DMPlexCreate(MPI_Comm comm, DM *mesh) 864 { 865 PetscErrorCode ierr; 866 867 PetscFunctionBegin; 868 PetscValidPointer(mesh,2); 869 ierr = DMCreate(comm, mesh);CHKERRQ(ierr); 870 ierr = DMSetType(*mesh, DMPLEX);CHKERRQ(ierr); 871 PetscFunctionReturn(0); 872 } 873 874 #undef __FUNCT__ 875 #define __FUNCT__ "DMPlexBuildFromCellList_Private" 876 /* 877 This takes as input the common mesh generator output, a list of the vertices for each cell 878 */ 879 PetscErrorCode DMPlexBuildFromCellList_Private(DM dm, PetscInt numCells, PetscInt numVertices, PetscInt numCorners, const int cells[]) 880 { 881 PetscInt *cone, c, p; 882 PetscErrorCode ierr; 883 884 PetscFunctionBegin; 885 ierr = DMPlexSetChart(dm, 0, numCells+numVertices);CHKERRQ(ierr); 886 for (c = 0; c < numCells; ++c) { 887 ierr = DMPlexSetConeSize(dm, c, numCorners);CHKERRQ(ierr); 888 } 889 ierr = DMSetUp(dm);CHKERRQ(ierr); 890 ierr = DMGetWorkArray(dm, numCorners, PETSC_INT, &cone);CHKERRQ(ierr); 891 for (c = 0; c < numCells; ++c) { 892 for (p = 0; p < numCorners; ++p) { 893 cone[p] = cells[c*numCorners+p]+numCells; 894 } 895 ierr = DMPlexSetCone(dm, c, cone);CHKERRQ(ierr); 896 } 897 ierr = DMRestoreWorkArray(dm, numCorners, PETSC_INT, &cone);CHKERRQ(ierr); 898 ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr); 899 ierr = DMPlexStratify(dm);CHKERRQ(ierr); 900 PetscFunctionReturn(0); 901 } 902 903 #undef __FUNCT__ 904 #define __FUNCT__ "DMPlexBuildCoordinates_Private" 905 /* 906 This takes as input the coordinates for each vertex 907 */ 908 PetscErrorCode DMPlexBuildCoordinates_Private(DM dm, PetscInt spaceDim, PetscInt numCells, PetscInt numVertices, const double vertexCoords[]) 909 { 910 PetscSection coordSection; 911 Vec coordinates; 912 PetscScalar *coords; 913 PetscInt coordSize, v, d; 914 PetscErrorCode ierr; 915 916 PetscFunctionBegin; 917 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 918 ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr); 919 ierr = PetscSectionSetFieldComponents(coordSection, 0, spaceDim);CHKERRQ(ierr); 920 ierr = PetscSectionSetChart(coordSection, numCells, numCells + numVertices);CHKERRQ(ierr); 921 for (v = numCells; v < numCells+numVertices; ++v) { 922 ierr = PetscSectionSetDof(coordSection, v, spaceDim);CHKERRQ(ierr); 923 ierr = PetscSectionSetFieldDof(coordSection, v, 0, spaceDim);CHKERRQ(ierr); 924 } 925 ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 926 ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 927 ierr = VecCreate(PetscObjectComm((PetscObject)dm), &coordinates);CHKERRQ(ierr); 928 ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr); 929 ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 930 ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr); 931 ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 932 for (v = 0; v < numVertices; ++v) { 933 for (d = 0; d < spaceDim; ++d) { 934 coords[v*spaceDim+d] = vertexCoords[v*spaceDim+d]; 935 } 936 } 937 ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 938 ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr); 939 ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 940 PetscFunctionReturn(0); 941 } 942 943 #undef __FUNCT__ 944 #define __FUNCT__ "DMPlexCreateFromCellList" 945 /*@C 946 DMPlexCreateFromCellList - This takes as input common mesh generator output, a list of the vertices for each cell, and produces a DM 947 948 Input Parameters: 949 + comm - The communicator 950 . dim - The topological dimension of the mesh 951 . numCells - The number of cells 952 . numVertices - The number of vertices 953 . numCorners - The number of vertices for each cell 954 . interpolate - Flag indicating that intermediate mesh entities (faces, edges) should be created automatically 955 . cells - An array of numCells*numCorners numbers, the vertices for each cell 956 . spaceDim - The spatial dimension used for coordinates 957 - vertexCoords - An array of numVertices*spaceDim numbers, the coordinates of each vertex 958 959 Output Parameter: 960 . dm - The DM 961 962 Note: Two triangles sharing a face 963 $ 964 $ 2 965 $ / | \ 966 $ / | \ 967 $ / | \ 968 $ 0 0 | 1 3 969 $ \ | / 970 $ \ | / 971 $ \ | / 972 $ 1 973 would have input 974 $ numCells = 2, numVertices = 4 975 $ cells = [0 1 2 1 3 2] 976 $ 977 which would result in the DMPlex 978 $ 979 $ 4 980 $ / | \ 981 $ / | \ 982 $ / | \ 983 $ 2 0 | 1 5 984 $ \ | / 985 $ \ | / 986 $ \ | / 987 $ 3 988 989 Level: beginner 990 991 .seealso: DMPlexCreateFromDAG(), DMPlexCreate() 992 @*/ 993 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) 994 { 995 PetscErrorCode ierr; 996 997 PetscFunctionBegin; 998 ierr = DMCreate(comm, dm);CHKERRQ(ierr); 999 ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 1000 ierr = DMPlexSetDimension(*dm, dim);CHKERRQ(ierr); 1001 ierr = DMPlexBuildFromCellList_Private(*dm, numCells, numVertices, numCorners, cells);CHKERRQ(ierr); 1002 if (interpolate) { 1003 DM idm; 1004 1005 ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr); 1006 ierr = DMDestroy(dm);CHKERRQ(ierr); 1007 *dm = idm; 1008 } 1009 ierr = DMPlexBuildCoordinates_Private(*dm, spaceDim, numCells, numVertices, vertexCoords);CHKERRQ(ierr); 1010 PetscFunctionReturn(0); 1011 } 1012 1013 #undef __FUNCT__ 1014 #define __FUNCT__ "DMPlexCreateFromDAG" 1015 /*@ 1016 DMPlexCreateFromDAG - This takes as input the adjacency-list representation of the Directed Acyclic Graph (Hasse Diagram) encoding a mesh, and produces a DM 1017 1018 Input Parameters: 1019 + dm - The empty DM object, usually from DMCreate() and DMPlexSetDimension() 1020 . depth - The depth of the DAG 1021 . numPoints - The number of points at each depth 1022 . coneSize - The cone size of each point 1023 . cones - The concatenation of the cone points for each point, the cone list must be oriented correctly for each point 1024 . coneOrientations - The orientation of each cone point 1025 - vertexCoords - An array of numVertices*dim numbers, the coordinates of each vertex 1026 1027 Output Parameter: 1028 . dm - The DM 1029 1030 Note: Two triangles sharing a face would have input 1031 $ depth = 1, numPoints = [4 2], coneSize = [3 3 0 0 0 0] 1032 $ cones = [2 3 4 3 5 4], coneOrientations = [0 0 0 0 0 0] 1033 $ vertexCoords = [-1.0 0.0 0.0 -1.0 0.0 1.0 1.0 0.0] 1034 $ 1035 which would result in the DMPlex 1036 $ 1037 $ 4 1038 $ / | \ 1039 $ / | \ 1040 $ / | \ 1041 $ 2 0 | 1 5 1042 $ \ | / 1043 $ \ | / 1044 $ \ | / 1045 $ 3 1046 $ 1047 $ Notice that all points are numbered consecutively, unlikely DMPlexCreateFromCellList() 1048 1049 Level: advanced 1050 1051 .seealso: DMPlexCreateFromCellList(), DMPlexCreate() 1052 @*/ 1053 PetscErrorCode DMPlexCreateFromDAG(DM dm, PetscInt depth, const PetscInt numPoints[], const PetscInt coneSize[], const PetscInt cones[], const PetscInt coneOrientations[], const PetscScalar vertexCoords[]) 1054 { 1055 Vec coordinates; 1056 PetscSection coordSection; 1057 PetscScalar *coords; 1058 PetscInt coordSize, firstVertex = numPoints[depth], pStart = 0, pEnd = 0, p, v, dim, d, off; 1059 PetscErrorCode ierr; 1060 1061 PetscFunctionBegin; 1062 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 1063 for (d = 0; d <= depth; ++d) pEnd += numPoints[d]; 1064 ierr = DMPlexSetChart(dm, pStart, pEnd);CHKERRQ(ierr); 1065 for (p = pStart; p < pEnd; ++p) { 1066 ierr = DMPlexSetConeSize(dm, p, coneSize[p-pStart]);CHKERRQ(ierr); 1067 } 1068 ierr = DMSetUp(dm);CHKERRQ(ierr); /* Allocate space for cones */ 1069 for (p = pStart, off = 0; p < pEnd; off += coneSize[p-pStart], ++p) { 1070 ierr = DMPlexSetCone(dm, p, &cones[off]);CHKERRQ(ierr); 1071 ierr = DMPlexSetConeOrientation(dm, p, &coneOrientations[off]);CHKERRQ(ierr); 1072 } 1073 ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr); 1074 ierr = DMPlexStratify(dm);CHKERRQ(ierr); 1075 /* Build coordinates */ 1076 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 1077 ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr); 1078 ierr = PetscSectionSetFieldComponents(coordSection, 0, dim);CHKERRQ(ierr); 1079 ierr = PetscSectionSetChart(coordSection, firstVertex, firstVertex+numPoints[0]);CHKERRQ(ierr); 1080 for (v = firstVertex; v < firstVertex+numPoints[0]; ++v) { 1081 ierr = PetscSectionSetDof(coordSection, v, dim);CHKERRQ(ierr); 1082 ierr = PetscSectionSetFieldDof(coordSection, v, 0, dim);CHKERRQ(ierr); 1083 } 1084 ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 1085 ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 1086 ierr = VecCreate(PetscObjectComm((PetscObject)dm), &coordinates);CHKERRQ(ierr); 1087 ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr); 1088 ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 1089 ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr); 1090 ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 1091 for (v = 0; v < numPoints[0]; ++v) { 1092 PetscInt off; 1093 1094 ierr = PetscSectionGetOffset(coordSection, v+firstVertex, &off);CHKERRQ(ierr); 1095 for (d = 0; d < dim; ++d) { 1096 coords[off+d] = vertexCoords[v*dim+d]; 1097 } 1098 } 1099 ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 1100 ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr); 1101 ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 1102 PetscFunctionReturn(0); 1103 } 1104