1 #define PETSCDM_DLL 2 #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 3 #include <petsc/private/hashseti.h> /*I "petscdmplex.h" I*/ 4 #include <petscsf.h> 5 #include <petscdmplextransform.h> 6 #include <petscdmlabelephemeral.h> 7 #include <petsc/private/kernels/blockmatmult.h> 8 #include <petsc/private/kernels/blockinvert.h> 9 10 PetscLogEvent DMPLEX_CreateFromFile, DMPLEX_CreateFromOptions, DMPLEX_BuildFromCellList, DMPLEX_BuildCoordinatesFromCellList; 11 12 /* External function declarations here */ 13 static PetscErrorCode DMInitialize_Plex(DM dm); 14 15 /* This copies internal things in the Plex structure that we generally want when making a new, related Plex */ 16 PetscErrorCode DMPlexCopy_Internal(DM dmin, PetscBool copyPeriodicity, PetscBool copyOverlap, DM dmout) 17 { 18 const PetscReal *maxCell, *Lstart, *L; 19 PetscBool dist, useCeed; 20 DMPlexReorderDefaultFlag reorder; 21 22 PetscFunctionBegin; 23 if (copyPeriodicity) { 24 PetscCall(DMGetPeriodicity(dmin, &maxCell, &Lstart, &L)); 25 PetscCall(DMSetPeriodicity(dmout, maxCell, Lstart, L)); 26 } 27 PetscCall(DMPlexDistributeGetDefault(dmin, &dist)); 28 PetscCall(DMPlexDistributeSetDefault(dmout, dist)); 29 PetscCall(DMPlexReorderGetDefault(dmin, &reorder)); 30 PetscCall(DMPlexReorderSetDefault(dmout, reorder)); 31 PetscCall(DMPlexGetUseCeed(dmin, &useCeed)); 32 PetscCall(DMPlexSetUseCeed(dmout, useCeed)); 33 ((DM_Plex *)dmout->data)->useHashLocation = ((DM_Plex *)dmin->data)->useHashLocation; 34 ((DM_Plex *)dmout->data)->printSetValues = ((DM_Plex *)dmin->data)->printSetValues; 35 ((DM_Plex *)dmout->data)->printFEM = ((DM_Plex *)dmin->data)->printFEM; 36 ((DM_Plex *)dmout->data)->printFVM = ((DM_Plex *)dmin->data)->printFVM; 37 ((DM_Plex *)dmout->data)->printL2 = ((DM_Plex *)dmin->data)->printL2; 38 ((DM_Plex *)dmout->data)->printLocate = ((DM_Plex *)dmin->data)->printLocate; 39 ((DM_Plex *)dmout->data)->printTol = ((DM_Plex *)dmin->data)->printTol; 40 if (copyOverlap) PetscCall(DMPlexSetOverlap_Plex(dmout, dmin, 0)); 41 PetscFunctionReturn(PETSC_SUCCESS); 42 } 43 44 /* Replace dm with the contents of ndm, and then destroy ndm 45 - Share the DM_Plex structure 46 - Share the coordinates 47 - Share the SF 48 */ 49 PetscErrorCode DMPlexReplace_Internal(DM dm, DM *ndm) 50 { 51 PetscSF sf; 52 DM dmNew = *ndm, coordDM, coarseDM; 53 Vec coords; 54 const PetscReal *maxCell, *Lstart, *L; 55 PetscInt dim, cdim; 56 57 PetscFunctionBegin; 58 if (dm == dmNew) { 59 PetscCall(DMDestroy(ndm)); 60 PetscFunctionReturn(PETSC_SUCCESS); 61 } 62 dm->setupcalled = dmNew->setupcalled; 63 PetscCall(DMGetDimension(dmNew, &dim)); 64 PetscCall(DMSetDimension(dm, dim)); 65 PetscCall(DMGetCoordinateDim(dmNew, &cdim)); 66 PetscCall(DMSetCoordinateDim(dm, cdim)); 67 PetscCall(DMGetPointSF(dmNew, &sf)); 68 PetscCall(DMSetPointSF(dm, sf)); 69 PetscCall(DMGetCoordinateDM(dmNew, &coordDM)); 70 PetscCall(DMGetCoordinatesLocal(dmNew, &coords)); 71 PetscCall(DMSetCoordinateDM(dm, coordDM)); 72 PetscCall(DMSetCoordinatesLocal(dm, coords)); 73 PetscCall(DMGetCellCoordinateDM(dmNew, &coordDM)); 74 PetscCall(DMGetCellCoordinatesLocal(dmNew, &coords)); 75 PetscCall(DMSetCellCoordinateDM(dm, coordDM)); 76 PetscCall(DMSetCellCoordinatesLocal(dm, coords)); 77 /* Do not want to create the coordinate field if it does not already exist, so do not call DMGetCoordinateField() */ 78 PetscCall(DMFieldDestroy(&dm->coordinates[0].field)); 79 dm->coordinates[0].field = dmNew->coordinates[0].field; 80 ((DM_Plex *)dmNew->data)->coordFunc = ((DM_Plex *)dm->data)->coordFunc; 81 PetscCall(DMGetPeriodicity(dmNew, &maxCell, &Lstart, &L)); 82 PetscCall(DMSetPeriodicity(dm, maxCell, Lstart, L)); 83 PetscCall(DMPlexGetGlobalToNaturalSF(dmNew, &sf)); 84 PetscCall(DMPlexSetGlobalToNaturalSF(dm, sf)); 85 PetscCall(DMDestroy_Plex(dm)); 86 PetscCall(DMInitialize_Plex(dm)); 87 dm->data = dmNew->data; 88 ((DM_Plex *)dmNew->data)->refct++; 89 PetscCall(DMPlexGetIsoperiodicFaceSF(dm, &sf)); 90 PetscCall(DMPlexSetIsoperiodicFaceSF(dm, sf)); // for the compose function effect on dm 91 PetscCall(DMDestroyLabelLinkList_Internal(dm)); 92 PetscCall(DMCopyLabels(dmNew, dm, PETSC_OWN_POINTER, PETSC_TRUE, DM_COPY_LABELS_FAIL)); 93 PetscCall(DMGetCoarseDM(dmNew, &coarseDM)); 94 PetscCall(DMSetCoarseDM(dm, coarseDM)); 95 PetscCall(DMDestroy(ndm)); 96 PetscFunctionReturn(PETSC_SUCCESS); 97 } 98 99 /* Swap dm with the contents of dmNew 100 - Swap the DM_Plex structure 101 - Swap the coordinates 102 - Swap the point PetscSF 103 */ 104 static PetscErrorCode DMPlexSwap_Static(DM dmA, DM dmB) 105 { 106 DM coordDMA, coordDMB; 107 Vec coordsA, coordsB; 108 PetscSF sfA, sfB; 109 DMField fieldTmp; 110 void *tmp; 111 DMLabelLink listTmp; 112 DMLabel depthTmp; 113 PetscInt tmpI; 114 115 PetscFunctionBegin; 116 if (dmA == dmB) PetscFunctionReturn(PETSC_SUCCESS); 117 PetscCall(DMGetPointSF(dmA, &sfA)); 118 PetscCall(DMGetPointSF(dmB, &sfB)); 119 PetscCall(PetscObjectReference((PetscObject)sfA)); 120 PetscCall(DMSetPointSF(dmA, sfB)); 121 PetscCall(DMSetPointSF(dmB, sfA)); 122 PetscCall(PetscObjectDereference((PetscObject)sfA)); 123 124 PetscCall(DMGetCoordinateDM(dmA, &coordDMA)); 125 PetscCall(DMGetCoordinateDM(dmB, &coordDMB)); 126 PetscCall(PetscObjectReference((PetscObject)coordDMA)); 127 PetscCall(DMSetCoordinateDM(dmA, coordDMB)); 128 PetscCall(DMSetCoordinateDM(dmB, coordDMA)); 129 PetscCall(PetscObjectDereference((PetscObject)coordDMA)); 130 131 PetscCall(DMGetCoordinatesLocal(dmA, &coordsA)); 132 PetscCall(DMGetCoordinatesLocal(dmB, &coordsB)); 133 PetscCall(PetscObjectReference((PetscObject)coordsA)); 134 PetscCall(DMSetCoordinatesLocal(dmA, coordsB)); 135 PetscCall(DMSetCoordinatesLocal(dmB, coordsA)); 136 PetscCall(PetscObjectDereference((PetscObject)coordsA)); 137 138 PetscCall(DMGetCellCoordinateDM(dmA, &coordDMA)); 139 PetscCall(DMGetCellCoordinateDM(dmB, &coordDMB)); 140 PetscCall(PetscObjectReference((PetscObject)coordDMA)); 141 PetscCall(DMSetCellCoordinateDM(dmA, coordDMB)); 142 PetscCall(DMSetCellCoordinateDM(dmB, coordDMA)); 143 PetscCall(PetscObjectDereference((PetscObject)coordDMA)); 144 145 PetscCall(DMGetCellCoordinatesLocal(dmA, &coordsA)); 146 PetscCall(DMGetCellCoordinatesLocal(dmB, &coordsB)); 147 PetscCall(PetscObjectReference((PetscObject)coordsA)); 148 PetscCall(DMSetCellCoordinatesLocal(dmA, coordsB)); 149 PetscCall(DMSetCellCoordinatesLocal(dmB, coordsA)); 150 PetscCall(PetscObjectDereference((PetscObject)coordsA)); 151 152 fieldTmp = dmA->coordinates[0].field; 153 dmA->coordinates[0].field = dmB->coordinates[0].field; 154 dmB->coordinates[0].field = fieldTmp; 155 fieldTmp = dmA->coordinates[1].field; 156 dmA->coordinates[1].field = dmB->coordinates[1].field; 157 dmB->coordinates[1].field = fieldTmp; 158 tmp = dmA->data; 159 dmA->data = dmB->data; 160 dmB->data = tmp; 161 listTmp = dmA->labels; 162 dmA->labels = dmB->labels; 163 dmB->labels = listTmp; 164 depthTmp = dmA->depthLabel; 165 dmA->depthLabel = dmB->depthLabel; 166 dmB->depthLabel = depthTmp; 167 depthTmp = dmA->celltypeLabel; 168 dmA->celltypeLabel = dmB->celltypeLabel; 169 dmB->celltypeLabel = depthTmp; 170 tmpI = dmA->levelup; 171 dmA->levelup = dmB->levelup; 172 dmB->levelup = tmpI; 173 PetscFunctionReturn(PETSC_SUCCESS); 174 } 175 176 PetscErrorCode DMPlexInterpolateInPlace_Internal(DM dm) 177 { 178 DM idm; 179 180 PetscFunctionBegin; 181 PetscCall(DMPlexInterpolate(dm, &idm)); 182 PetscCall(DMPlexCopyCoordinates(dm, idm)); 183 PetscCall(DMPlexReplace_Internal(dm, &idm)); 184 PetscFunctionReturn(PETSC_SUCCESS); 185 } 186 187 /*@C 188 DMPlexCreateCoordinateSpace - Creates a finite element space for the coordinates 189 190 Collective 191 192 Input Parameters: 193 + dm - The `DMPLEX` 194 . degree - The degree of the finite element or `PETSC_DECIDE` 195 . project - Flag to project current coordinates into the space 196 - coordFunc - An optional function to map new points from refinement to the surface 197 198 Level: advanced 199 200 .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `PetscPointFunc`, `PetscFECreateLagrange()`, `DMGetCoordinateDM()` 201 @*/ 202 PetscErrorCode DMPlexCreateCoordinateSpace(DM dm, PetscInt degree, PetscBool project, PetscPointFunc coordFunc) 203 { 204 DM_Plex *mesh = (DM_Plex *)dm->data; 205 PetscFE fe = NULL; 206 DM cdm; 207 PetscInt dim, dE, qorder; 208 209 PetscFunctionBegin; 210 PetscCall(DMGetDimension(dm, &dim)); 211 PetscCall(DMGetCoordinateDim(dm, &dE)); 212 qorder = degree; 213 PetscCall(DMGetCoordinateDM(dm, &cdm)); 214 PetscObjectOptionsBegin((PetscObject)cdm); 215 PetscCall(PetscOptionsBoundedInt("-default_quadrature_order", "Quadrature order is one less than quadrature points per edge", "DMPlexCreateCoordinateSpace", qorder, &qorder, NULL, 0)); 216 PetscOptionsEnd(); 217 if (degree >= 0) { 218 DMPolytopeType ct = DM_POLYTOPE_UNKNOWN; 219 PetscInt cStart, cEnd, gct; 220 221 PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 222 if (cEnd > cStart) PetscCall(DMPlexGetCellType(dm, cStart, &ct)); 223 gct = (PetscInt)ct; 224 PetscCall(MPIU_Allreduce(MPI_IN_PLACE, &gct, 1, MPIU_INT, MPI_MIN, PetscObjectComm((PetscObject)dm))); 225 ct = (DMPolytopeType)gct; 226 // Work around current bug in PetscDualSpaceSetUp_Lagrange() 227 // Can be seen in plex_tutorials-ex10_1 228 if (ct != DM_POLYTOPE_SEG_PRISM_TENSOR && ct != DM_POLYTOPE_TRI_PRISM_TENSOR && ct != DM_POLYTOPE_QUAD_PRISM_TENSOR) PetscCall(PetscFECreateLagrangeByCell(PETSC_COMM_SELF, dim, dE, ct, degree, qorder, &fe)); 229 } 230 PetscCall(DMSetCoordinateDisc(dm, fe, project)); 231 PetscCall(PetscFEDestroy(&fe)); 232 mesh->coordFunc = coordFunc; 233 PetscFunctionReturn(PETSC_SUCCESS); 234 } 235 236 /*@ 237 DMPlexCreateDoublet - Creates a mesh of two cells of the specified type, optionally with later refinement. 238 239 Collective 240 241 Input Parameters: 242 + comm - The communicator for the `DM` object 243 . dim - The spatial dimension 244 . simplex - Flag for simplicial cells, otherwise they are tensor product cells 245 . interpolate - Flag to create intermediate mesh pieces (edges, faces) 246 - refinementLimit - A nonzero number indicates the largest admissible volume for a refined cell 247 248 Output Parameter: 249 . newdm - The `DM` object 250 251 Level: beginner 252 253 .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMSetType()`, `DMCreate()` 254 @*/ 255 PetscErrorCode DMPlexCreateDoublet(MPI_Comm comm, PetscInt dim, PetscBool simplex, PetscBool interpolate, PetscReal refinementLimit, DM *newdm) 256 { 257 DM dm; 258 PetscMPIInt rank; 259 260 PetscFunctionBegin; 261 PetscCall(DMCreate(comm, &dm)); 262 PetscCall(DMSetType(dm, DMPLEX)); 263 PetscCall(DMSetDimension(dm, dim)); 264 PetscCall(PetscLogEventBegin(DMPLEX_Generate, dm, 0, 0, 0)); 265 PetscCallMPI(MPI_Comm_rank(comm, &rank)); 266 switch (dim) { 267 case 2: 268 if (simplex) PetscCall(PetscObjectSetName((PetscObject)dm, "triangular")); 269 else PetscCall(PetscObjectSetName((PetscObject)dm, "quadrilateral")); 270 break; 271 case 3: 272 if (simplex) PetscCall(PetscObjectSetName((PetscObject)dm, "tetrahedral")); 273 else PetscCall(PetscObjectSetName((PetscObject)dm, "hexahedral")); 274 break; 275 default: 276 SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "Cannot make meshes for dimension %" PetscInt_FMT, dim); 277 } 278 if (rank) { 279 PetscInt numPoints[2] = {0, 0}; 280 PetscCall(DMPlexCreateFromDAG(dm, 1, numPoints, NULL, NULL, NULL, NULL)); 281 } else { 282 switch (dim) { 283 case 2: 284 if (simplex) { 285 PetscInt numPoints[2] = {4, 2}; 286 PetscInt coneSize[6] = {3, 3, 0, 0, 0, 0}; 287 PetscInt cones[6] = {2, 3, 4, 5, 4, 3}; 288 PetscInt coneOrientations[6] = {0, 0, 0, 0, 0, 0}; 289 PetscScalar vertexCoords[8] = {-0.5, 0.5, 0.0, 0.0, 0.0, 1.0, 0.5, 0.5}; 290 291 PetscCall(DMPlexCreateFromDAG(dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords)); 292 } else { 293 PetscInt numPoints[2] = {6, 2}; 294 PetscInt coneSize[8] = {4, 4, 0, 0, 0, 0, 0, 0}; 295 PetscInt cones[8] = {2, 3, 4, 5, 3, 6, 7, 4}; 296 PetscInt coneOrientations[8] = {0, 0, 0, 0, 0, 0, 0, 0}; 297 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}; 298 299 PetscCall(DMPlexCreateFromDAG(dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords)); 300 } 301 break; 302 case 3: 303 if (simplex) { 304 PetscInt numPoints[2] = {5, 2}; 305 PetscInt coneSize[7] = {4, 4, 0, 0, 0, 0, 0}; 306 PetscInt cones[8] = {4, 3, 5, 2, 5, 3, 4, 6}; 307 PetscInt coneOrientations[8] = {0, 0, 0, 0, 0, 0, 0, 0}; 308 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}; 309 310 PetscCall(DMPlexCreateFromDAG(dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords)); 311 } else { 312 PetscInt numPoints[2] = {12, 2}; 313 PetscInt coneSize[14] = {8, 8, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0}; 314 PetscInt cones[16] = {2, 3, 4, 5, 6, 7, 8, 9, 5, 4, 10, 11, 7, 12, 13, 8}; 315 PetscInt coneOrientations[16] = {0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0}; 316 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, -1.0, -0.5, 0.5, 0.0, -0.5, 0.5, 0.0, 0.5, 0.5, -1.0, 0.5, 0.5, 1.0, 0.5, -0.5, 1.0, -0.5, -0.5, 1.0, -0.5, 0.5, 1.0, 0.5, 0.5}; 317 318 PetscCall(DMPlexCreateFromDAG(dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords)); 319 } 320 break; 321 default: 322 SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "Cannot make meshes for dimension %" PetscInt_FMT, dim); 323 } 324 } 325 PetscCall(PetscLogEventEnd(DMPLEX_Generate, dm, 0, 0, 0)); 326 *newdm = dm; 327 if (refinementLimit > 0.0) { 328 DM rdm; 329 const char *name; 330 331 PetscCall(DMPlexSetRefinementUniform(*newdm, PETSC_FALSE)); 332 PetscCall(DMPlexSetRefinementLimit(*newdm, refinementLimit)); 333 PetscCall(DMRefine(*newdm, comm, &rdm)); 334 PetscCall(PetscObjectGetName((PetscObject)*newdm, &name)); 335 PetscCall(PetscObjectSetName((PetscObject)rdm, name)); 336 PetscCall(DMDestroy(newdm)); 337 *newdm = rdm; 338 } 339 if (interpolate) { 340 DM idm; 341 342 PetscCall(DMPlexInterpolate(*newdm, &idm)); 343 PetscCall(DMDestroy(newdm)); 344 *newdm = idm; 345 } 346 PetscFunctionReturn(PETSC_SUCCESS); 347 } 348 349 static PetscErrorCode DMPlexCreateBoxSurfaceMesh_Tensor_1D_Internal(DM dm, const PetscReal lower[], const PetscReal upper[], const PetscInt edges[]) 350 { 351 const PetscInt numVertices = 2; 352 PetscInt markerRight = 1; 353 PetscInt markerLeft = 1; 354 PetscBool markerSeparate = PETSC_FALSE; 355 Vec coordinates; 356 PetscSection coordSection; 357 PetscScalar *coords; 358 PetscInt coordSize; 359 PetscMPIInt rank; 360 PetscInt cdim = 1, v; 361 362 PetscFunctionBegin; 363 PetscCall(PetscOptionsGetBool(((PetscObject)dm)->options, ((PetscObject)dm)->prefix, "-dm_plex_separate_marker", &markerSeparate, NULL)); 364 if (markerSeparate) { 365 markerRight = 2; 366 markerLeft = 1; 367 } 368 PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank)); 369 if (rank == 0) { 370 PetscCall(DMPlexSetChart(dm, 0, numVertices)); 371 PetscCall(DMSetUp(dm)); /* Allocate space for cones */ 372 PetscCall(DMSetLabelValue(dm, "marker", 0, markerLeft)); 373 PetscCall(DMSetLabelValue(dm, "marker", 1, markerRight)); 374 } 375 PetscCall(DMPlexSymmetrize(dm)); 376 PetscCall(DMPlexStratify(dm)); 377 /* Build coordinates */ 378 PetscCall(DMSetCoordinateDim(dm, cdim)); 379 PetscCall(DMGetCoordinateSection(dm, &coordSection)); 380 PetscCall(PetscSectionSetNumFields(coordSection, 1)); 381 PetscCall(PetscSectionSetChart(coordSection, 0, numVertices)); 382 PetscCall(PetscSectionSetFieldComponents(coordSection, 0, cdim)); 383 for (v = 0; v < numVertices; ++v) { 384 PetscCall(PetscSectionSetDof(coordSection, v, cdim)); 385 PetscCall(PetscSectionSetFieldDof(coordSection, v, 0, cdim)); 386 } 387 PetscCall(PetscSectionSetUp(coordSection)); 388 PetscCall(PetscSectionGetStorageSize(coordSection, &coordSize)); 389 PetscCall(VecCreate(PETSC_COMM_SELF, &coordinates)); 390 PetscCall(PetscObjectSetName((PetscObject)coordinates, "coordinates")); 391 PetscCall(VecSetSizes(coordinates, coordSize, PETSC_DETERMINE)); 392 PetscCall(VecSetBlockSize(coordinates, cdim)); 393 PetscCall(VecSetType(coordinates, VECSTANDARD)); 394 PetscCall(VecGetArray(coordinates, &coords)); 395 coords[0] = lower[0]; 396 coords[1] = upper[0]; 397 PetscCall(VecRestoreArray(coordinates, &coords)); 398 PetscCall(DMSetCoordinatesLocal(dm, coordinates)); 399 PetscCall(VecDestroy(&coordinates)); 400 PetscFunctionReturn(PETSC_SUCCESS); 401 } 402 403 static PetscErrorCode DMPlexCreateBoxSurfaceMesh_Tensor_2D_Internal(DM dm, const PetscReal lower[], const PetscReal upper[], const PetscInt edges[]) 404 { 405 const PetscInt numVertices = (edges[0] + 1) * (edges[1] + 1); 406 const PetscInt numEdges = edges[0] * (edges[1] + 1) + (edges[0] + 1) * edges[1]; 407 PetscInt markerTop = 1; 408 PetscInt markerBottom = 1; 409 PetscInt markerRight = 1; 410 PetscInt markerLeft = 1; 411 PetscBool markerSeparate = PETSC_FALSE; 412 Vec coordinates; 413 PetscSection coordSection; 414 PetscScalar *coords; 415 PetscInt coordSize; 416 PetscMPIInt rank; 417 PetscInt v, vx, vy; 418 419 PetscFunctionBegin; 420 PetscCall(PetscOptionsGetBool(((PetscObject)dm)->options, ((PetscObject)dm)->prefix, "-dm_plex_separate_marker", &markerSeparate, NULL)); 421 if (markerSeparate) { 422 markerTop = 3; 423 markerBottom = 1; 424 markerRight = 2; 425 markerLeft = 4; 426 } 427 PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank)); 428 if (rank == 0) { 429 PetscInt e, ex, ey; 430 431 PetscCall(DMPlexSetChart(dm, 0, numEdges + numVertices)); 432 for (e = 0; e < numEdges; ++e) PetscCall(DMPlexSetConeSize(dm, e, 2)); 433 PetscCall(DMSetUp(dm)); /* Allocate space for cones */ 434 for (vx = 0; vx <= edges[0]; vx++) { 435 for (ey = 0; ey < edges[1]; ey++) { 436 PetscInt edge = vx * edges[1] + ey + edges[0] * (edges[1] + 1); 437 PetscInt vertex = ey * (edges[0] + 1) + vx + numEdges; 438 PetscInt cone[2]; 439 440 cone[0] = vertex; 441 cone[1] = vertex + edges[0] + 1; 442 PetscCall(DMPlexSetCone(dm, edge, cone)); 443 if (vx == edges[0]) { 444 PetscCall(DMSetLabelValue(dm, "marker", edge, markerRight)); 445 PetscCall(DMSetLabelValue(dm, "marker", cone[0], markerRight)); 446 if (ey == edges[1] - 1) { 447 PetscCall(DMSetLabelValue(dm, "marker", cone[1], markerRight)); 448 PetscCall(DMSetLabelValue(dm, "Face Sets", cone[1], markerRight)); 449 } 450 } else if (vx == 0) { 451 PetscCall(DMSetLabelValue(dm, "marker", edge, markerLeft)); 452 PetscCall(DMSetLabelValue(dm, "marker", cone[0], markerLeft)); 453 if (ey == edges[1] - 1) { 454 PetscCall(DMSetLabelValue(dm, "marker", cone[1], markerLeft)); 455 PetscCall(DMSetLabelValue(dm, "Face Sets", cone[1], markerLeft)); 456 } 457 } 458 } 459 } 460 for (vy = 0; vy <= edges[1]; vy++) { 461 for (ex = 0; ex < edges[0]; ex++) { 462 PetscInt edge = vy * edges[0] + ex; 463 PetscInt vertex = vy * (edges[0] + 1) + ex + numEdges; 464 PetscInt cone[2]; 465 466 cone[0] = vertex; 467 cone[1] = vertex + 1; 468 PetscCall(DMPlexSetCone(dm, edge, cone)); 469 if (vy == edges[1]) { 470 PetscCall(DMSetLabelValue(dm, "marker", edge, markerTop)); 471 PetscCall(DMSetLabelValue(dm, "marker", cone[0], markerTop)); 472 if (ex == edges[0] - 1) { 473 PetscCall(DMSetLabelValue(dm, "marker", cone[1], markerTop)); 474 PetscCall(DMSetLabelValue(dm, "Face Sets", cone[1], markerTop)); 475 } 476 } else if (vy == 0) { 477 PetscCall(DMSetLabelValue(dm, "marker", edge, markerBottom)); 478 PetscCall(DMSetLabelValue(dm, "marker", cone[0], markerBottom)); 479 if (ex == edges[0] - 1) { 480 PetscCall(DMSetLabelValue(dm, "marker", cone[1], markerBottom)); 481 PetscCall(DMSetLabelValue(dm, "Face Sets", cone[1], markerBottom)); 482 } 483 } 484 } 485 } 486 } 487 PetscCall(DMPlexSymmetrize(dm)); 488 PetscCall(DMPlexStratify(dm)); 489 /* Build coordinates */ 490 PetscCall(DMSetCoordinateDim(dm, 2)); 491 PetscCall(DMGetCoordinateSection(dm, &coordSection)); 492 PetscCall(PetscSectionSetNumFields(coordSection, 1)); 493 PetscCall(PetscSectionSetChart(coordSection, numEdges, numEdges + numVertices)); 494 PetscCall(PetscSectionSetFieldComponents(coordSection, 0, 2)); 495 for (v = numEdges; v < numEdges + numVertices; ++v) { 496 PetscCall(PetscSectionSetDof(coordSection, v, 2)); 497 PetscCall(PetscSectionSetFieldDof(coordSection, v, 0, 2)); 498 } 499 PetscCall(PetscSectionSetUp(coordSection)); 500 PetscCall(PetscSectionGetStorageSize(coordSection, &coordSize)); 501 PetscCall(VecCreate(PETSC_COMM_SELF, &coordinates)); 502 PetscCall(PetscObjectSetName((PetscObject)coordinates, "coordinates")); 503 PetscCall(VecSetSizes(coordinates, coordSize, PETSC_DETERMINE)); 504 PetscCall(VecSetBlockSize(coordinates, 2)); 505 PetscCall(VecSetType(coordinates, VECSTANDARD)); 506 PetscCall(VecGetArray(coordinates, &coords)); 507 for (vy = 0; vy <= edges[1]; ++vy) { 508 for (vx = 0; vx <= edges[0]; ++vx) { 509 coords[(vy * (edges[0] + 1) + vx) * 2 + 0] = lower[0] + ((upper[0] - lower[0]) / edges[0]) * vx; 510 coords[(vy * (edges[0] + 1) + vx) * 2 + 1] = lower[1] + ((upper[1] - lower[1]) / edges[1]) * vy; 511 } 512 } 513 PetscCall(VecRestoreArray(coordinates, &coords)); 514 PetscCall(DMSetCoordinatesLocal(dm, coordinates)); 515 PetscCall(VecDestroy(&coordinates)); 516 PetscFunctionReturn(PETSC_SUCCESS); 517 } 518 519 static PetscErrorCode DMPlexCreateBoxSurfaceMesh_Tensor_3D_Internal(DM dm, const PetscReal lower[], const PetscReal upper[], const PetscInt faces[]) 520 { 521 PetscInt vertices[3], numVertices; 522 PetscInt numFaces = 2 * faces[0] * faces[1] + 2 * faces[1] * faces[2] + 2 * faces[0] * faces[2]; 523 PetscInt markerTop = 1; 524 PetscInt markerBottom = 1; 525 PetscInt markerFront = 1; 526 PetscInt markerBack = 1; 527 PetscInt markerRight = 1; 528 PetscInt markerLeft = 1; 529 PetscBool markerSeparate = PETSC_FALSE; 530 Vec coordinates; 531 PetscSection coordSection; 532 PetscScalar *coords; 533 PetscInt coordSize; 534 PetscMPIInt rank; 535 PetscInt v, vx, vy, vz; 536 PetscInt voffset, iface = 0, cone[4]; 537 538 PetscFunctionBegin; 539 PetscCheck(faces[0] >= 1 && faces[1] >= 1 && faces[2] >= 1, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Must have at least 1 face per side"); 540 PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank)); 541 PetscCall(PetscOptionsGetBool(((PetscObject)dm)->options, ((PetscObject)dm)->prefix, "-dm_plex_separate_marker", &markerSeparate, NULL)); 542 if (markerSeparate) { 543 markerBottom = 1; 544 markerTop = 2; 545 markerFront = 3; 546 markerBack = 4; 547 markerRight = 5; 548 markerLeft = 6; 549 } 550 vertices[0] = faces[0] + 1; 551 vertices[1] = faces[1] + 1; 552 vertices[2] = faces[2] + 1; 553 numVertices = vertices[0] * vertices[1] * vertices[2]; 554 if (rank == 0) { 555 PetscInt f; 556 557 PetscCall(DMPlexSetChart(dm, 0, numFaces + numVertices)); 558 for (f = 0; f < numFaces; ++f) PetscCall(DMPlexSetConeSize(dm, f, 4)); 559 PetscCall(DMSetUp(dm)); /* Allocate space for cones */ 560 561 /* Side 0 (Top) */ 562 for (vy = 0; vy < faces[1]; vy++) { 563 for (vx = 0; vx < faces[0]; vx++) { 564 voffset = numFaces + vertices[0] * vertices[1] * (vertices[2] - 1) + vy * vertices[0] + vx; 565 cone[0] = voffset; 566 cone[1] = voffset + 1; 567 cone[2] = voffset + vertices[0] + 1; 568 cone[3] = voffset + vertices[0]; 569 PetscCall(DMPlexSetCone(dm, iface, cone)); 570 PetscCall(DMSetLabelValue(dm, "marker", iface, markerTop)); 571 PetscCall(DMSetLabelValue(dm, "marker", voffset + 0, markerTop)); 572 PetscCall(DMSetLabelValue(dm, "marker", voffset + 1, markerTop)); 573 PetscCall(DMSetLabelValue(dm, "marker", voffset + vertices[0] + 0, markerTop)); 574 PetscCall(DMSetLabelValue(dm, "marker", voffset + vertices[0] + 1, markerTop)); 575 iface++; 576 } 577 } 578 579 /* Side 1 (Bottom) */ 580 for (vy = 0; vy < faces[1]; vy++) { 581 for (vx = 0; vx < faces[0]; vx++) { 582 voffset = numFaces + vy * (faces[0] + 1) + vx; 583 cone[0] = voffset + 1; 584 cone[1] = voffset; 585 cone[2] = voffset + vertices[0]; 586 cone[3] = voffset + vertices[0] + 1; 587 PetscCall(DMPlexSetCone(dm, iface, cone)); 588 PetscCall(DMSetLabelValue(dm, "marker", iface, markerBottom)); 589 PetscCall(DMSetLabelValue(dm, "marker", voffset + 0, markerBottom)); 590 PetscCall(DMSetLabelValue(dm, "marker", voffset + 1, markerBottom)); 591 PetscCall(DMSetLabelValue(dm, "marker", voffset + vertices[0] + 0, markerBottom)); 592 PetscCall(DMSetLabelValue(dm, "marker", voffset + vertices[0] + 1, markerBottom)); 593 iface++; 594 } 595 } 596 597 /* Side 2 (Front) */ 598 for (vz = 0; vz < faces[2]; vz++) { 599 for (vx = 0; vx < faces[0]; vx++) { 600 voffset = numFaces + vz * vertices[0] * vertices[1] + vx; 601 cone[0] = voffset; 602 cone[1] = voffset + 1; 603 cone[2] = voffset + vertices[0] * vertices[1] + 1; 604 cone[3] = voffset + vertices[0] * vertices[1]; 605 PetscCall(DMPlexSetCone(dm, iface, cone)); 606 PetscCall(DMSetLabelValue(dm, "marker", iface, markerFront)); 607 PetscCall(DMSetLabelValue(dm, "marker", voffset + 0, markerFront)); 608 PetscCall(DMSetLabelValue(dm, "marker", voffset + 1, markerFront)); 609 PetscCall(DMSetLabelValue(dm, "marker", voffset + vertices[0] * vertices[1] + 0, markerFront)); 610 PetscCall(DMSetLabelValue(dm, "marker", voffset + vertices[0] * vertices[1] + 1, markerFront)); 611 iface++; 612 } 613 } 614 615 /* Side 3 (Back) */ 616 for (vz = 0; vz < faces[2]; vz++) { 617 for (vx = 0; vx < faces[0]; vx++) { 618 voffset = numFaces + vz * vertices[0] * vertices[1] + vertices[0] * (vertices[1] - 1) + vx; 619 cone[0] = voffset + vertices[0] * vertices[1]; 620 cone[1] = voffset + vertices[0] * vertices[1] + 1; 621 cone[2] = voffset + 1; 622 cone[3] = voffset; 623 PetscCall(DMPlexSetCone(dm, iface, cone)); 624 PetscCall(DMSetLabelValue(dm, "marker", iface, markerBack)); 625 PetscCall(DMSetLabelValue(dm, "marker", voffset + 0, markerBack)); 626 PetscCall(DMSetLabelValue(dm, "marker", voffset + 1, markerBack)); 627 PetscCall(DMSetLabelValue(dm, "marker", voffset + vertices[0] * vertices[1] + 0, markerBack)); 628 PetscCall(DMSetLabelValue(dm, "marker", voffset + vertices[0] * vertices[1] + 1, markerBack)); 629 iface++; 630 } 631 } 632 633 /* Side 4 (Left) */ 634 for (vz = 0; vz < faces[2]; vz++) { 635 for (vy = 0; vy < faces[1]; vy++) { 636 voffset = numFaces + vz * vertices[0] * vertices[1] + vy * vertices[0]; 637 cone[0] = voffset; 638 cone[1] = voffset + vertices[0] * vertices[1]; 639 cone[2] = voffset + vertices[0] * vertices[1] + vertices[0]; 640 cone[3] = voffset + vertices[0]; 641 PetscCall(DMPlexSetCone(dm, iface, cone)); 642 PetscCall(DMSetLabelValue(dm, "marker", iface, markerLeft)); 643 PetscCall(DMSetLabelValue(dm, "marker", voffset + 0, markerLeft)); 644 PetscCall(DMSetLabelValue(dm, "marker", voffset + vertices[0] + 0, markerLeft)); 645 PetscCall(DMSetLabelValue(dm, "marker", voffset + vertices[1] + 0, markerLeft)); 646 PetscCall(DMSetLabelValue(dm, "marker", voffset + vertices[0] * vertices[1] + vertices[0], markerLeft)); 647 iface++; 648 } 649 } 650 651 /* Side 5 (Right) */ 652 for (vz = 0; vz < faces[2]; vz++) { 653 for (vy = 0; vy < faces[1]; vy++) { 654 voffset = numFaces + vz * vertices[0] * vertices[1] + vy * vertices[0] + faces[0]; 655 cone[0] = voffset + vertices[0] * vertices[1]; 656 cone[1] = voffset; 657 cone[2] = voffset + vertices[0]; 658 cone[3] = voffset + vertices[0] * vertices[1] + vertices[0]; 659 PetscCall(DMPlexSetCone(dm, iface, cone)); 660 PetscCall(DMSetLabelValue(dm, "marker", iface, markerRight)); 661 PetscCall(DMSetLabelValue(dm, "marker", voffset + 0, markerRight)); 662 PetscCall(DMSetLabelValue(dm, "marker", voffset + vertices[0] + 0, markerRight)); 663 PetscCall(DMSetLabelValue(dm, "marker", voffset + vertices[0] * vertices[1] + 0, markerRight)); 664 PetscCall(DMSetLabelValue(dm, "marker", voffset + vertices[0] * vertices[1] + vertices[0], markerRight)); 665 iface++; 666 } 667 } 668 } 669 PetscCall(DMPlexSymmetrize(dm)); 670 PetscCall(DMPlexStratify(dm)); 671 /* Build coordinates */ 672 PetscCall(DMSetCoordinateDim(dm, 3)); 673 PetscCall(DMGetCoordinateSection(dm, &coordSection)); 674 PetscCall(PetscSectionSetNumFields(coordSection, 1)); 675 PetscCall(PetscSectionSetChart(coordSection, numFaces, numFaces + numVertices)); 676 PetscCall(PetscSectionSetFieldComponents(coordSection, 0, 3)); 677 for (v = numFaces; v < numFaces + numVertices; ++v) { 678 PetscCall(PetscSectionSetDof(coordSection, v, 3)); 679 PetscCall(PetscSectionSetFieldDof(coordSection, v, 0, 3)); 680 } 681 PetscCall(PetscSectionSetUp(coordSection)); 682 PetscCall(PetscSectionGetStorageSize(coordSection, &coordSize)); 683 PetscCall(VecCreate(PETSC_COMM_SELF, &coordinates)); 684 PetscCall(PetscObjectSetName((PetscObject)coordinates, "coordinates")); 685 PetscCall(VecSetSizes(coordinates, coordSize, PETSC_DETERMINE)); 686 PetscCall(VecSetBlockSize(coordinates, 3)); 687 PetscCall(VecSetType(coordinates, VECSTANDARD)); 688 PetscCall(VecGetArray(coordinates, &coords)); 689 for (vz = 0; vz <= faces[2]; ++vz) { 690 for (vy = 0; vy <= faces[1]; ++vy) { 691 for (vx = 0; vx <= faces[0]; ++vx) { 692 coords[((vz * (faces[1] + 1) + vy) * (faces[0] + 1) + vx) * 3 + 0] = lower[0] + ((upper[0] - lower[0]) / faces[0]) * vx; 693 coords[((vz * (faces[1] + 1) + vy) * (faces[0] + 1) + vx) * 3 + 1] = lower[1] + ((upper[1] - lower[1]) / faces[1]) * vy; 694 coords[((vz * (faces[1] + 1) + vy) * (faces[0] + 1) + vx) * 3 + 2] = lower[2] + ((upper[2] - lower[2]) / faces[2]) * vz; 695 } 696 } 697 } 698 PetscCall(VecRestoreArray(coordinates, &coords)); 699 PetscCall(DMSetCoordinatesLocal(dm, coordinates)); 700 PetscCall(VecDestroy(&coordinates)); 701 PetscFunctionReturn(PETSC_SUCCESS); 702 } 703 704 static PetscErrorCode DMPlexCreateBoxSurfaceMesh_Internal(DM dm, PetscInt dim, const PetscInt faces[], const PetscReal lower[], const PetscReal upper[], PetscBool interpolate) 705 { 706 PetscFunctionBegin; 707 PetscValidLogicalCollectiveInt(dm, dim, 2); 708 PetscCall(PetscLogEventBegin(DMPLEX_Generate, dm, 0, 0, 0)); 709 PetscCall(DMSetDimension(dm, dim - 1)); 710 PetscCall(DMSetCoordinateDim(dm, dim)); 711 switch (dim) { 712 case 1: 713 PetscCall(DMPlexCreateBoxSurfaceMesh_Tensor_1D_Internal(dm, lower, upper, faces)); 714 break; 715 case 2: 716 PetscCall(DMPlexCreateBoxSurfaceMesh_Tensor_2D_Internal(dm, lower, upper, faces)); 717 break; 718 case 3: 719 PetscCall(DMPlexCreateBoxSurfaceMesh_Tensor_3D_Internal(dm, lower, upper, faces)); 720 break; 721 default: 722 SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Dimension not supported: %" PetscInt_FMT, dim); 723 } 724 PetscCall(PetscLogEventEnd(DMPLEX_Generate, dm, 0, 0, 0)); 725 if (interpolate) PetscCall(DMPlexInterpolateInPlace_Internal(dm)); 726 PetscFunctionReturn(PETSC_SUCCESS); 727 } 728 729 /*@C 730 DMPlexCreateBoxSurfaceMesh - Creates a mesh on the surface of the tensor product of unit intervals (box) using tensor cells (hexahedra). 731 732 Collective 733 734 Input Parameters: 735 + comm - The communicator for the `DM` object 736 . dim - The spatial dimension of the box, so the resulting mesh is has dimension `dim`-1 737 . faces - Number of faces per dimension, or `NULL` for (1,) in 1D and (2, 2) in 2D and (1, 1, 1) in 3D 738 . lower - The lower left corner, or `NULL` for (0, 0, 0) 739 . upper - The upper right corner, or `NULL` for (1, 1, 1) 740 - interpolate - Flag to create intermediate mesh pieces (edges, faces) 741 742 Output Parameter: 743 . dm - The `DM` object 744 745 Level: beginner 746 747 .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMSetFromOptions()`, `DMPlexCreateBoxMesh()`, `DMPlexCreateFromFile()`, `DMSetType()`, `DMCreate()` 748 @*/ 749 PetscErrorCode DMPlexCreateBoxSurfaceMesh(MPI_Comm comm, PetscInt dim, const PetscInt faces[], const PetscReal lower[], const PetscReal upper[], PetscBool interpolate, DM *dm) 750 { 751 PetscInt fac[3] = {1, 1, 1}; 752 PetscReal low[3] = {0, 0, 0}; 753 PetscReal upp[3] = {1, 1, 1}; 754 755 PetscFunctionBegin; 756 PetscCall(DMCreate(comm, dm)); 757 PetscCall(DMSetType(*dm, DMPLEX)); 758 PetscCall(DMPlexCreateBoxSurfaceMesh_Internal(*dm, dim, faces ? faces : fac, lower ? lower : low, upper ? upper : upp, interpolate)); 759 PetscFunctionReturn(PETSC_SUCCESS); 760 } 761 762 static PetscErrorCode DMPlexCreateLineMesh_Internal(DM dm, PetscInt segments, PetscReal lower, PetscReal upper, DMBoundaryType bd) 763 { 764 PetscInt i, fStart, fEnd, numCells = 0, numVerts = 0; 765 PetscInt numPoints[2], *coneSize, *cones, *coneOrientations; 766 PetscScalar *vertexCoords; 767 PetscReal L, maxCell; 768 PetscBool markerSeparate = PETSC_FALSE; 769 PetscInt markerLeft = 1, faceMarkerLeft = 1; 770 PetscInt markerRight = 1, faceMarkerRight = 2; 771 PetscBool wrap = (bd == DM_BOUNDARY_PERIODIC || bd == DM_BOUNDARY_TWIST) ? PETSC_TRUE : PETSC_FALSE; 772 PetscMPIInt rank; 773 774 PetscFunctionBegin; 775 PetscAssertPointer(dm, 1); 776 777 PetscCall(DMSetDimension(dm, 1)); 778 PetscCall(DMCreateLabel(dm, "marker")); 779 PetscCall(DMCreateLabel(dm, "Face Sets")); 780 781 PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank)); 782 if (rank == 0) numCells = segments; 783 if (rank == 0) numVerts = segments + (wrap ? 0 : 1); 784 785 numPoints[0] = numVerts; 786 numPoints[1] = numCells; 787 PetscCall(PetscMalloc4(numCells + numVerts, &coneSize, numCells * 2, &cones, numCells + numVerts, &coneOrientations, numVerts, &vertexCoords)); 788 PetscCall(PetscArrayzero(coneOrientations, numCells + numVerts)); 789 for (i = 0; i < numCells; ++i) coneSize[i] = 2; 790 for (i = 0; i < numVerts; ++i) coneSize[numCells + i] = 0; 791 for (i = 0; i < numCells; ++i) { 792 cones[2 * i] = numCells + i % numVerts; 793 cones[2 * i + 1] = numCells + (i + 1) % numVerts; 794 } 795 for (i = 0; i < numVerts; ++i) vertexCoords[i] = lower + (upper - lower) * ((PetscReal)i / (PetscReal)numCells); 796 PetscCall(DMPlexCreateFromDAG(dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords)); 797 PetscCall(PetscFree4(coneSize, cones, coneOrientations, vertexCoords)); 798 799 PetscCall(PetscOptionsGetBool(((PetscObject)dm)->options, ((PetscObject)dm)->prefix, "-dm_plex_separate_marker", &markerSeparate, NULL)); 800 if (markerSeparate) { 801 markerLeft = faceMarkerLeft; 802 markerRight = faceMarkerRight; 803 } 804 if (!wrap && rank == 0) { 805 PetscCall(DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd)); 806 PetscCall(DMSetLabelValue(dm, "marker", fStart, markerLeft)); 807 PetscCall(DMSetLabelValue(dm, "marker", fEnd - 1, markerRight)); 808 PetscCall(DMSetLabelValue(dm, "Face Sets", fStart, faceMarkerLeft)); 809 PetscCall(DMSetLabelValue(dm, "Face Sets", fEnd - 1, faceMarkerRight)); 810 } 811 if (wrap) { 812 L = upper - lower; 813 maxCell = (PetscReal)1.1 * (L / (PetscReal)PetscMax(1, segments)); 814 PetscCall(DMSetPeriodicity(dm, &maxCell, &lower, &L)); 815 } 816 PetscCall(DMPlexSetRefinementUniform(dm, PETSC_TRUE)); 817 PetscFunctionReturn(PETSC_SUCCESS); 818 } 819 820 static PetscErrorCode DMPlexCreateBoxMesh_Simplex_Internal(DM dm, PetscInt dim, const PetscInt faces[], const PetscReal lower[], const PetscReal upper[], const DMBoundaryType periodicity[], PetscBool interpolate) 821 { 822 DM boundary, vol; 823 DMLabel bdlabel; 824 825 PetscFunctionBegin; 826 PetscAssertPointer(dm, 1); 827 for (PetscInt i = 0; i < dim; ++i) PetscCheck(periodicity[i] == DM_BOUNDARY_NONE, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Periodicity is not supported for simplex meshes"); 828 PetscCall(DMCreate(PetscObjectComm((PetscObject)dm), &boundary)); 829 PetscCall(DMSetType(boundary, DMPLEX)); 830 PetscCall(DMPlexCreateBoxSurfaceMesh_Internal(boundary, dim, faces, lower, upper, PETSC_FALSE)); 831 PetscCall(DMPlexGenerate(boundary, NULL, interpolate, &vol)); 832 PetscCall(DMGetLabel(vol, "marker", &bdlabel)); 833 if (bdlabel) PetscCall(DMPlexLabelComplete(vol, bdlabel)); 834 PetscCall(DMPlexCopy_Internal(dm, PETSC_TRUE, PETSC_FALSE, vol)); 835 PetscCall(DMPlexReplace_Internal(dm, &vol)); 836 PetscCall(DMDestroy(&boundary)); 837 PetscFunctionReturn(PETSC_SUCCESS); 838 } 839 840 static PetscErrorCode DMPlexCreateCubeMesh_Internal(DM dm, const PetscReal lower[], const PetscReal upper[], const PetscInt edges[], DMBoundaryType bdX, DMBoundaryType bdY, DMBoundaryType bdZ) 841 { 842 DMLabel cutLabel = NULL; 843 PetscInt markerTop = 1, faceMarkerTop = 1; 844 PetscInt markerBottom = 1, faceMarkerBottom = 1; 845 PetscInt markerFront = 1, faceMarkerFront = 1; 846 PetscInt markerBack = 1, faceMarkerBack = 1; 847 PetscInt markerRight = 1, faceMarkerRight = 1; 848 PetscInt markerLeft = 1, faceMarkerLeft = 1; 849 PetscInt dim; 850 PetscBool markerSeparate = PETSC_FALSE, cutMarker = PETSC_FALSE; 851 PetscMPIInt rank; 852 853 PetscFunctionBegin; 854 PetscCall(DMGetDimension(dm, &dim)); 855 PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank)); 856 PetscCall(DMCreateLabel(dm, "marker")); 857 PetscCall(DMCreateLabel(dm, "Face Sets")); 858 PetscCall(PetscOptionsGetBool(((PetscObject)dm)->options, ((PetscObject)dm)->prefix, "-dm_plex_periodic_cut", &cutMarker, NULL)); 859 if (bdX == DM_BOUNDARY_PERIODIC || bdX == DM_BOUNDARY_TWIST || bdY == DM_BOUNDARY_PERIODIC || bdY == DM_BOUNDARY_TWIST || bdZ == DM_BOUNDARY_PERIODIC || bdZ == DM_BOUNDARY_TWIST) { 860 if (cutMarker) { 861 PetscCall(DMCreateLabel(dm, "periodic_cut")); 862 PetscCall(DMGetLabel(dm, "periodic_cut", &cutLabel)); 863 } 864 } 865 switch (dim) { 866 case 2: 867 faceMarkerTop = 3; 868 faceMarkerBottom = 1; 869 faceMarkerRight = 2; 870 faceMarkerLeft = 4; 871 break; 872 case 3: 873 faceMarkerBottom = 1; 874 faceMarkerTop = 2; 875 faceMarkerFront = 3; 876 faceMarkerBack = 4; 877 faceMarkerRight = 5; 878 faceMarkerLeft = 6; 879 break; 880 default: 881 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Dimension %" PetscInt_FMT " not supported", dim); 882 } 883 PetscCall(PetscOptionsGetBool(((PetscObject)dm)->options, ((PetscObject)dm)->prefix, "-dm_plex_separate_marker", &markerSeparate, NULL)); 884 if (markerSeparate) { 885 markerBottom = faceMarkerBottom; 886 markerTop = faceMarkerTop; 887 markerFront = faceMarkerFront; 888 markerBack = faceMarkerBack; 889 markerRight = faceMarkerRight; 890 markerLeft = faceMarkerLeft; 891 } 892 { 893 const PetscInt numXEdges = rank == 0 ? edges[0] : 0; 894 const PetscInt numYEdges = rank == 0 ? edges[1] : 0; 895 const PetscInt numZEdges = rank == 0 ? edges[2] : 0; 896 const PetscInt numXVertices = rank == 0 ? (bdX == DM_BOUNDARY_PERIODIC || bdX == DM_BOUNDARY_TWIST ? edges[0] : edges[0] + 1) : 0; 897 const PetscInt numYVertices = rank == 0 ? (bdY == DM_BOUNDARY_PERIODIC || bdY == DM_BOUNDARY_TWIST ? edges[1] : edges[1] + 1) : 0; 898 const PetscInt numZVertices = rank == 0 ? (bdZ == DM_BOUNDARY_PERIODIC || bdZ == DM_BOUNDARY_TWIST ? edges[2] : edges[2] + 1) : 0; 899 const PetscInt numCells = numXEdges * numYEdges * numZEdges; 900 const PetscInt numXFaces = numYEdges * numZEdges; 901 const PetscInt numYFaces = numXEdges * numZEdges; 902 const PetscInt numZFaces = numXEdges * numYEdges; 903 const PetscInt numTotXFaces = numXVertices * numXFaces; 904 const PetscInt numTotYFaces = numYVertices * numYFaces; 905 const PetscInt numTotZFaces = numZVertices * numZFaces; 906 const PetscInt numFaces = numTotXFaces + numTotYFaces + numTotZFaces; 907 const PetscInt numTotXEdges = numXEdges * numYVertices * numZVertices; 908 const PetscInt numTotYEdges = numYEdges * numXVertices * numZVertices; 909 const PetscInt numTotZEdges = numZEdges * numXVertices * numYVertices; 910 const PetscInt numVertices = numXVertices * numYVertices * numZVertices; 911 const PetscInt numEdges = numTotXEdges + numTotYEdges + numTotZEdges; 912 const PetscInt firstVertex = (dim == 2) ? numFaces : numCells; 913 const PetscInt firstXFace = (dim == 2) ? 0 : numCells + numVertices; 914 const PetscInt firstYFace = firstXFace + numTotXFaces; 915 const PetscInt firstZFace = firstYFace + numTotYFaces; 916 const PetscInt firstXEdge = numCells + numFaces + numVertices; 917 const PetscInt firstYEdge = firstXEdge + numTotXEdges; 918 const PetscInt firstZEdge = firstYEdge + numTotYEdges; 919 Vec coordinates; 920 PetscSection coordSection; 921 PetscScalar *coords; 922 PetscInt coordSize; 923 PetscInt v, vx, vy, vz; 924 PetscInt c, f, fx, fy, fz, e, ex, ey, ez; 925 926 PetscCall(DMPlexSetChart(dm, 0, numCells + numFaces + numEdges + numVertices)); 927 for (c = 0; c < numCells; c++) PetscCall(DMPlexSetConeSize(dm, c, 6)); 928 for (f = firstXFace; f < firstXFace + numFaces; ++f) PetscCall(DMPlexSetConeSize(dm, f, 4)); 929 for (e = firstXEdge; e < firstXEdge + numEdges; ++e) PetscCall(DMPlexSetConeSize(dm, e, 2)); 930 PetscCall(DMSetUp(dm)); /* Allocate space for cones */ 931 /* Build cells */ 932 for (fz = 0; fz < numZEdges; ++fz) { 933 for (fy = 0; fy < numYEdges; ++fy) { 934 for (fx = 0; fx < numXEdges; ++fx) { 935 PetscInt cell = (fz * numYEdges + fy) * numXEdges + fx; 936 PetscInt faceB = firstZFace + (fy * numXEdges + fx) * numZVertices + fz; 937 PetscInt faceT = firstZFace + (fy * numXEdges + fx) * numZVertices + ((fz + 1) % numZVertices); 938 PetscInt faceF = firstYFace + (fz * numXEdges + fx) * numYVertices + fy; 939 PetscInt faceK = firstYFace + (fz * numXEdges + fx) * numYVertices + ((fy + 1) % numYVertices); 940 PetscInt faceL = firstXFace + (fz * numYEdges + fy) * numXVertices + fx; 941 PetscInt faceR = firstXFace + (fz * numYEdges + fy) * numXVertices + ((fx + 1) % numXVertices); 942 /* B, T, F, K, R, L */ 943 PetscInt ornt[6] = {-2, 0, 0, -3, 0, -2}; /* ??? */ 944 PetscInt cone[6]; 945 946 /* no boundary twisting in 3D */ 947 cone[0] = faceB; 948 cone[1] = faceT; 949 cone[2] = faceF; 950 cone[3] = faceK; 951 cone[4] = faceR; 952 cone[5] = faceL; 953 PetscCall(DMPlexSetCone(dm, cell, cone)); 954 PetscCall(DMPlexSetConeOrientation(dm, cell, ornt)); 955 if (bdX != DM_BOUNDARY_NONE && fx == numXEdges - 1 && cutLabel) PetscCall(DMLabelSetValue(cutLabel, cell, 2)); 956 if (bdY != DM_BOUNDARY_NONE && fy == numYEdges - 1 && cutLabel) PetscCall(DMLabelSetValue(cutLabel, cell, 2)); 957 if (bdZ != DM_BOUNDARY_NONE && fz == numZEdges - 1 && cutLabel) PetscCall(DMLabelSetValue(cutLabel, cell, 2)); 958 } 959 } 960 } 961 /* Build x faces */ 962 for (fz = 0; fz < numZEdges; ++fz) { 963 for (fy = 0; fy < numYEdges; ++fy) { 964 for (fx = 0; fx < numXVertices; ++fx) { 965 PetscInt face = firstXFace + (fz * numYEdges + fy) * numXVertices + fx; 966 PetscInt edgeL = firstZEdge + (fy * numXVertices + fx) * numZEdges + fz; 967 PetscInt edgeR = firstZEdge + (((fy + 1) % numYVertices) * numXVertices + fx) * numZEdges + fz; 968 PetscInt edgeB = firstYEdge + (fz * numXVertices + fx) * numYEdges + fy; 969 PetscInt edgeT = firstYEdge + (((fz + 1) % numZVertices) * numXVertices + fx) * numYEdges + fy; 970 PetscInt ornt[4] = {0, 0, -1, -1}; 971 PetscInt cone[4]; 972 973 if (dim == 3) { 974 /* markers */ 975 if (bdX != DM_BOUNDARY_PERIODIC) { 976 if (fx == numXVertices - 1) { 977 PetscCall(DMSetLabelValue(dm, "Face Sets", face, faceMarkerRight)); 978 PetscCall(DMSetLabelValue(dm, "marker", face, markerRight)); 979 } else if (fx == 0) { 980 PetscCall(DMSetLabelValue(dm, "Face Sets", face, faceMarkerLeft)); 981 PetscCall(DMSetLabelValue(dm, "marker", face, markerLeft)); 982 } 983 } 984 } 985 cone[0] = edgeB; 986 cone[1] = edgeR; 987 cone[2] = edgeT; 988 cone[3] = edgeL; 989 PetscCall(DMPlexSetCone(dm, face, cone)); 990 PetscCall(DMPlexSetConeOrientation(dm, face, ornt)); 991 } 992 } 993 } 994 /* Build y faces */ 995 for (fz = 0; fz < numZEdges; ++fz) { 996 for (fx = 0; fx < numXEdges; ++fx) { 997 for (fy = 0; fy < numYVertices; ++fy) { 998 PetscInt face = firstYFace + (fz * numXEdges + fx) * numYVertices + fy; 999 PetscInt edgeL = firstZEdge + (fy * numXVertices + fx) * numZEdges + fz; 1000 PetscInt edgeR = firstZEdge + (fy * numXVertices + ((fx + 1) % numXVertices)) * numZEdges + fz; 1001 PetscInt edgeB = firstXEdge + (fz * numYVertices + fy) * numXEdges + fx; 1002 PetscInt edgeT = firstXEdge + (((fz + 1) % numZVertices) * numYVertices + fy) * numXEdges + fx; 1003 PetscInt ornt[4] = {0, 0, -1, -1}; 1004 PetscInt cone[4]; 1005 1006 if (dim == 3) { 1007 /* markers */ 1008 if (bdY != DM_BOUNDARY_PERIODIC) { 1009 if (fy == numYVertices - 1) { 1010 PetscCall(DMSetLabelValue(dm, "Face Sets", face, faceMarkerBack)); 1011 PetscCall(DMSetLabelValue(dm, "marker", face, markerBack)); 1012 } else if (fy == 0) { 1013 PetscCall(DMSetLabelValue(dm, "Face Sets", face, faceMarkerFront)); 1014 PetscCall(DMSetLabelValue(dm, "marker", face, markerFront)); 1015 } 1016 } 1017 } 1018 cone[0] = edgeB; 1019 cone[1] = edgeR; 1020 cone[2] = edgeT; 1021 cone[3] = edgeL; 1022 PetscCall(DMPlexSetCone(dm, face, cone)); 1023 PetscCall(DMPlexSetConeOrientation(dm, face, ornt)); 1024 } 1025 } 1026 } 1027 /* Build z faces */ 1028 for (fy = 0; fy < numYEdges; ++fy) { 1029 for (fx = 0; fx < numXEdges; ++fx) { 1030 for (fz = 0; fz < numZVertices; fz++) { 1031 PetscInt face = firstZFace + (fy * numXEdges + fx) * numZVertices + fz; 1032 PetscInt edgeL = firstYEdge + (fz * numXVertices + fx) * numYEdges + fy; 1033 PetscInt edgeR = firstYEdge + (fz * numXVertices + ((fx + 1) % numXVertices)) * numYEdges + fy; 1034 PetscInt edgeB = firstXEdge + (fz * numYVertices + fy) * numXEdges + fx; 1035 PetscInt edgeT = firstXEdge + (fz * numYVertices + ((fy + 1) % numYVertices)) * numXEdges + fx; 1036 PetscInt ornt[4] = {0, 0, -1, -1}; 1037 PetscInt cone[4]; 1038 1039 if (dim == 2) { 1040 if (bdX == DM_BOUNDARY_TWIST && fx == numXEdges - 1) { 1041 edgeR += numYEdges - 1 - 2 * fy; 1042 ornt[1] = -1; 1043 } 1044 if (bdY == DM_BOUNDARY_TWIST && fy == numYEdges - 1) { 1045 edgeT += numXEdges - 1 - 2 * fx; 1046 ornt[2] = 0; 1047 } 1048 if (bdX != DM_BOUNDARY_NONE && fx == numXEdges - 1 && cutLabel) PetscCall(DMLabelSetValue(cutLabel, face, 2)); 1049 if (bdY != DM_BOUNDARY_NONE && fy == numYEdges - 1 && cutLabel) PetscCall(DMLabelSetValue(cutLabel, face, 2)); 1050 } else { 1051 /* markers */ 1052 if (bdZ != DM_BOUNDARY_PERIODIC) { 1053 if (fz == numZVertices - 1) { 1054 PetscCall(DMSetLabelValue(dm, "Face Sets", face, faceMarkerTop)); 1055 PetscCall(DMSetLabelValue(dm, "marker", face, markerTop)); 1056 } else if (fz == 0) { 1057 PetscCall(DMSetLabelValue(dm, "Face Sets", face, faceMarkerBottom)); 1058 PetscCall(DMSetLabelValue(dm, "marker", face, markerBottom)); 1059 } 1060 } 1061 } 1062 cone[0] = edgeB; 1063 cone[1] = edgeR; 1064 cone[2] = edgeT; 1065 cone[3] = edgeL; 1066 PetscCall(DMPlexSetCone(dm, face, cone)); 1067 PetscCall(DMPlexSetConeOrientation(dm, face, ornt)); 1068 } 1069 } 1070 } 1071 /* Build Z edges*/ 1072 for (vy = 0; vy < numYVertices; vy++) { 1073 for (vx = 0; vx < numXVertices; vx++) { 1074 for (ez = 0; ez < numZEdges; ez++) { 1075 const PetscInt edge = firstZEdge + (vy * numXVertices + vx) * numZEdges + ez; 1076 const PetscInt vertexB = firstVertex + (ez * numYVertices + vy) * numXVertices + vx; 1077 const PetscInt vertexT = firstVertex + (((ez + 1) % numZVertices) * numYVertices + vy) * numXVertices + vx; 1078 PetscInt cone[2]; 1079 1080 cone[0] = vertexB; 1081 cone[1] = vertexT; 1082 PetscCall(DMPlexSetCone(dm, edge, cone)); 1083 if (dim == 3) { 1084 if (bdX != DM_BOUNDARY_PERIODIC) { 1085 if (vx == numXVertices - 1) { 1086 PetscCall(DMSetLabelValue(dm, "marker", edge, markerRight)); 1087 PetscCall(DMSetLabelValue(dm, "marker", cone[0], markerRight)); 1088 if (ez == numZEdges - 1) PetscCall(DMSetLabelValue(dm, "marker", cone[1], markerRight)); 1089 } else if (vx == 0) { 1090 PetscCall(DMSetLabelValue(dm, "marker", edge, markerLeft)); 1091 PetscCall(DMSetLabelValue(dm, "marker", cone[0], markerLeft)); 1092 if (ez == numZEdges - 1) PetscCall(DMSetLabelValue(dm, "marker", cone[1], markerLeft)); 1093 } 1094 } 1095 if (bdY != DM_BOUNDARY_PERIODIC) { 1096 if (vy == numYVertices - 1) { 1097 PetscCall(DMSetLabelValue(dm, "marker", edge, markerBack)); 1098 PetscCall(DMSetLabelValue(dm, "marker", cone[0], markerBack)); 1099 if (ez == numZEdges - 1) PetscCall(DMSetLabelValue(dm, "marker", cone[1], markerBack)); 1100 } else if (vy == 0) { 1101 PetscCall(DMSetLabelValue(dm, "marker", edge, markerFront)); 1102 PetscCall(DMSetLabelValue(dm, "marker", cone[0], markerFront)); 1103 if (ez == numZEdges - 1) PetscCall(DMSetLabelValue(dm, "marker", cone[1], markerFront)); 1104 } 1105 } 1106 } 1107 } 1108 } 1109 } 1110 /* Build Y edges*/ 1111 for (vz = 0; vz < numZVertices; vz++) { 1112 for (vx = 0; vx < numXVertices; vx++) { 1113 for (ey = 0; ey < numYEdges; ey++) { 1114 const PetscInt nextv = (dim == 2 && bdY == DM_BOUNDARY_TWIST && ey == numYEdges - 1) ? (numXVertices - vx - 1) : (vz * numYVertices + ((ey + 1) % numYVertices)) * numXVertices + vx; 1115 const PetscInt edge = firstYEdge + (vz * numXVertices + vx) * numYEdges + ey; 1116 const PetscInt vertexF = firstVertex + (vz * numYVertices + ey) * numXVertices + vx; 1117 const PetscInt vertexK = firstVertex + nextv; 1118 PetscInt cone[2]; 1119 1120 cone[0] = vertexF; 1121 cone[1] = vertexK; 1122 PetscCall(DMPlexSetCone(dm, edge, cone)); 1123 if (dim == 2) { 1124 if ((bdX != DM_BOUNDARY_PERIODIC) && (bdX != DM_BOUNDARY_TWIST)) { 1125 if (vx == numXVertices - 1) { 1126 PetscCall(DMSetLabelValue(dm, "Face Sets", edge, faceMarkerRight)); 1127 PetscCall(DMSetLabelValue(dm, "marker", edge, markerRight)); 1128 PetscCall(DMSetLabelValue(dm, "marker", cone[0], markerRight)); 1129 if (ey == numYEdges - 1) PetscCall(DMSetLabelValue(dm, "marker", cone[1], markerRight)); 1130 } else if (vx == 0) { 1131 PetscCall(DMSetLabelValue(dm, "Face Sets", edge, faceMarkerLeft)); 1132 PetscCall(DMSetLabelValue(dm, "marker", edge, markerLeft)); 1133 PetscCall(DMSetLabelValue(dm, "marker", cone[0], markerLeft)); 1134 if (ey == numYEdges - 1) PetscCall(DMSetLabelValue(dm, "marker", cone[1], markerLeft)); 1135 } 1136 } else { 1137 if (vx == 0 && cutLabel) { 1138 PetscCall(DMLabelSetValue(cutLabel, edge, 1)); 1139 PetscCall(DMLabelSetValue(cutLabel, cone[0], 1)); 1140 if (ey == numYEdges - 1) PetscCall(DMLabelSetValue(cutLabel, cone[1], 1)); 1141 } 1142 } 1143 } else { 1144 if (bdX != DM_BOUNDARY_PERIODIC) { 1145 if (vx == numXVertices - 1) { 1146 PetscCall(DMSetLabelValue(dm, "marker", edge, markerRight)); 1147 PetscCall(DMSetLabelValue(dm, "marker", cone[0], markerRight)); 1148 if (ey == numYEdges - 1) PetscCall(DMSetLabelValue(dm, "marker", cone[1], markerRight)); 1149 } else if (vx == 0) { 1150 PetscCall(DMSetLabelValue(dm, "marker", edge, markerLeft)); 1151 PetscCall(DMSetLabelValue(dm, "marker", cone[0], markerLeft)); 1152 if (ey == numYEdges - 1) PetscCall(DMSetLabelValue(dm, "marker", cone[1], markerLeft)); 1153 } 1154 } 1155 if (bdZ != DM_BOUNDARY_PERIODIC) { 1156 if (vz == numZVertices - 1) { 1157 PetscCall(DMSetLabelValue(dm, "marker", edge, markerTop)); 1158 PetscCall(DMSetLabelValue(dm, "marker", cone[0], markerTop)); 1159 if (ey == numYEdges - 1) PetscCall(DMSetLabelValue(dm, "marker", cone[1], markerTop)); 1160 } else if (vz == 0) { 1161 PetscCall(DMSetLabelValue(dm, "marker", edge, markerBottom)); 1162 PetscCall(DMSetLabelValue(dm, "marker", cone[0], markerBottom)); 1163 if (ey == numYEdges - 1) PetscCall(DMSetLabelValue(dm, "marker", cone[1], markerBottom)); 1164 } 1165 } 1166 } 1167 } 1168 } 1169 } 1170 /* Build X edges*/ 1171 for (vz = 0; vz < numZVertices; vz++) { 1172 for (vy = 0; vy < numYVertices; vy++) { 1173 for (ex = 0; ex < numXEdges; ex++) { 1174 const PetscInt nextv = (dim == 2 && bdX == DM_BOUNDARY_TWIST && ex == numXEdges - 1) ? (numYVertices - vy - 1) * numXVertices : (vz * numYVertices + vy) * numXVertices + (ex + 1) % numXVertices; 1175 const PetscInt edge = firstXEdge + (vz * numYVertices + vy) * numXEdges + ex; 1176 const PetscInt vertexL = firstVertex + (vz * numYVertices + vy) * numXVertices + ex; 1177 const PetscInt vertexR = firstVertex + nextv; 1178 PetscInt cone[2]; 1179 1180 cone[0] = vertexL; 1181 cone[1] = vertexR; 1182 PetscCall(DMPlexSetCone(dm, edge, cone)); 1183 if (dim == 2) { 1184 if ((bdY != DM_BOUNDARY_PERIODIC) && (bdY != DM_BOUNDARY_TWIST)) { 1185 if (vy == numYVertices - 1) { 1186 PetscCall(DMSetLabelValue(dm, "Face Sets", edge, faceMarkerTop)); 1187 PetscCall(DMSetLabelValue(dm, "marker", edge, markerTop)); 1188 PetscCall(DMSetLabelValue(dm, "marker", cone[0], markerTop)); 1189 if (ex == numXEdges - 1) PetscCall(DMSetLabelValue(dm, "marker", cone[1], markerTop)); 1190 } else if (vy == 0) { 1191 PetscCall(DMSetLabelValue(dm, "Face Sets", edge, faceMarkerBottom)); 1192 PetscCall(DMSetLabelValue(dm, "marker", edge, markerBottom)); 1193 PetscCall(DMSetLabelValue(dm, "marker", cone[0], markerBottom)); 1194 if (ex == numXEdges - 1) PetscCall(DMSetLabelValue(dm, "marker", cone[1], markerBottom)); 1195 } 1196 } else { 1197 if (vy == 0 && cutLabel) { 1198 PetscCall(DMLabelSetValue(cutLabel, edge, 1)); 1199 PetscCall(DMLabelSetValue(cutLabel, cone[0], 1)); 1200 if (ex == numXEdges - 1) PetscCall(DMLabelSetValue(cutLabel, cone[1], 1)); 1201 } 1202 } 1203 } else { 1204 if (bdY != DM_BOUNDARY_PERIODIC) { 1205 if (vy == numYVertices - 1) { 1206 PetscCall(DMSetLabelValue(dm, "marker", edge, markerBack)); 1207 PetscCall(DMSetLabelValue(dm, "marker", cone[0], markerBack)); 1208 if (ex == numXEdges - 1) PetscCall(DMSetLabelValue(dm, "marker", cone[1], markerBack)); 1209 } else if (vy == 0) { 1210 PetscCall(DMSetLabelValue(dm, "marker", edge, markerFront)); 1211 PetscCall(DMSetLabelValue(dm, "marker", cone[0], markerFront)); 1212 if (ex == numXEdges - 1) PetscCall(DMSetLabelValue(dm, "marker", cone[1], markerFront)); 1213 } 1214 } 1215 if (bdZ != DM_BOUNDARY_PERIODIC) { 1216 if (vz == numZVertices - 1) { 1217 PetscCall(DMSetLabelValue(dm, "marker", edge, markerTop)); 1218 PetscCall(DMSetLabelValue(dm, "marker", cone[0], markerTop)); 1219 if (ex == numXEdges - 1) PetscCall(DMSetLabelValue(dm, "marker", cone[1], markerTop)); 1220 } else if (vz == 0) { 1221 PetscCall(DMSetLabelValue(dm, "marker", edge, markerBottom)); 1222 PetscCall(DMSetLabelValue(dm, "marker", cone[0], markerBottom)); 1223 if (ex == numXEdges - 1) PetscCall(DMSetLabelValue(dm, "marker", cone[1], markerBottom)); 1224 } 1225 } 1226 } 1227 } 1228 } 1229 } 1230 PetscCall(DMPlexSymmetrize(dm)); 1231 PetscCall(DMPlexStratify(dm)); 1232 /* Build coordinates */ 1233 PetscCall(DMGetCoordinateSection(dm, &coordSection)); 1234 PetscCall(PetscSectionSetNumFields(coordSection, 1)); 1235 PetscCall(PetscSectionSetFieldComponents(coordSection, 0, dim)); 1236 PetscCall(PetscSectionSetChart(coordSection, firstVertex, firstVertex + numVertices)); 1237 for (v = firstVertex; v < firstVertex + numVertices; ++v) { 1238 PetscCall(PetscSectionSetDof(coordSection, v, dim)); 1239 PetscCall(PetscSectionSetFieldDof(coordSection, v, 0, dim)); 1240 } 1241 PetscCall(PetscSectionSetUp(coordSection)); 1242 PetscCall(PetscSectionGetStorageSize(coordSection, &coordSize)); 1243 PetscCall(VecCreate(PETSC_COMM_SELF, &coordinates)); 1244 PetscCall(PetscObjectSetName((PetscObject)coordinates, "coordinates")); 1245 PetscCall(VecSetSizes(coordinates, coordSize, PETSC_DETERMINE)); 1246 PetscCall(VecSetBlockSize(coordinates, dim)); 1247 PetscCall(VecSetType(coordinates, VECSTANDARD)); 1248 PetscCall(VecGetArray(coordinates, &coords)); 1249 for (vz = 0; vz < numZVertices; ++vz) { 1250 for (vy = 0; vy < numYVertices; ++vy) { 1251 for (vx = 0; vx < numXVertices; ++vx) { 1252 coords[((vz * numYVertices + vy) * numXVertices + vx) * dim + 0] = lower[0] + ((upper[0] - lower[0]) / numXEdges) * vx; 1253 coords[((vz * numYVertices + vy) * numXVertices + vx) * dim + 1] = lower[1] + ((upper[1] - lower[1]) / numYEdges) * vy; 1254 if (dim == 3) coords[((vz * numYVertices + vy) * numXVertices + vx) * dim + 2] = lower[2] + ((upper[2] - lower[2]) / numZEdges) * vz; 1255 } 1256 } 1257 } 1258 PetscCall(VecRestoreArray(coordinates, &coords)); 1259 PetscCall(DMSetCoordinatesLocal(dm, coordinates)); 1260 PetscCall(VecDestroy(&coordinates)); 1261 } 1262 PetscFunctionReturn(PETSC_SUCCESS); 1263 } 1264 1265 static PetscErrorCode DMPlexCreateBoxMesh_Tensor_Internal(DM dm, PetscInt dim, const PetscInt faces[], const PetscReal lower[], const PetscReal upper[], const DMBoundaryType periodicity[]) 1266 { 1267 DMBoundaryType bdt[3] = {DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE}; 1268 PetscInt fac[3] = {0, 0, 0}, d; 1269 1270 PetscFunctionBegin; 1271 PetscAssertPointer(dm, 1); 1272 PetscValidLogicalCollectiveInt(dm, dim, 2); 1273 PetscCall(DMSetDimension(dm, dim)); 1274 for (d = 0; d < dim; ++d) { 1275 fac[d] = faces[d]; 1276 bdt[d] = periodicity[d]; 1277 } 1278 PetscCall(DMPlexCreateCubeMesh_Internal(dm, lower, upper, fac, bdt[0], bdt[1], bdt[2])); 1279 if (periodicity[0] == DM_BOUNDARY_PERIODIC || periodicity[0] == DM_BOUNDARY_TWIST || periodicity[1] == DM_BOUNDARY_PERIODIC || periodicity[1] == DM_BOUNDARY_TWIST || (dim > 2 && (periodicity[2] == DM_BOUNDARY_PERIODIC || periodicity[2] == DM_BOUNDARY_TWIST))) { 1280 PetscReal L[3] = {-1., -1., 0.}; 1281 PetscReal maxCell[3] = {-1., -1., 0.}; 1282 1283 for (d = 0; d < dim; ++d) { 1284 if (periodicity[d] != DM_BOUNDARY_NONE) { 1285 L[d] = upper[d] - lower[d]; 1286 maxCell[d] = 1.1 * (L[d] / PetscMax(1, faces[d])); 1287 } 1288 } 1289 PetscCall(DMSetPeriodicity(dm, maxCell, lower, L)); 1290 } 1291 PetscCall(DMPlexSetRefinementUniform(dm, PETSC_TRUE)); 1292 PetscFunctionReturn(PETSC_SUCCESS); 1293 } 1294 1295 static PetscErrorCode DMPlexCreateBoxMesh_Internal(DM dm, DMPlexShape shape, PetscInt dim, PetscBool simplex, const PetscInt faces[], const PetscReal lower[], const PetscReal upper[], const DMBoundaryType periodicity[], PetscBool interpolate) 1296 { 1297 PetscFunctionBegin; 1298 PetscCall(PetscLogEventBegin(DMPLEX_Generate, dm, 0, 0, 0)); 1299 if (shape == DM_SHAPE_ZBOX) PetscCall(DMPlexCreateBoxMesh_Tensor_SFC_Internal(dm, dim, faces, lower, upper, periodicity, interpolate)); 1300 else if (dim == 1) PetscCall(DMPlexCreateLineMesh_Internal(dm, faces[0], lower[0], upper[0], periodicity[0])); 1301 else if (simplex) PetscCall(DMPlexCreateBoxMesh_Simplex_Internal(dm, dim, faces, lower, upper, periodicity, interpolate)); 1302 else PetscCall(DMPlexCreateBoxMesh_Tensor_Internal(dm, dim, faces, lower, upper, periodicity)); 1303 if (!interpolate && dim > 1 && !simplex) { 1304 DM udm; 1305 1306 PetscCall(DMPlexUninterpolate(dm, &udm)); 1307 PetscCall(DMPlexCopyCoordinates(dm, udm)); 1308 PetscCall(DMPlexReplace_Internal(dm, &udm)); 1309 } 1310 PetscCall(PetscLogEventEnd(DMPLEX_Generate, dm, 0, 0, 0)); 1311 PetscFunctionReturn(PETSC_SUCCESS); 1312 } 1313 1314 /*@C 1315 DMPlexCreateBoxMesh - Creates a mesh on the tensor product of unit intervals (box) using simplices or tensor cells (hexahedra). 1316 1317 Collective 1318 1319 Input Parameters: 1320 + comm - The communicator for the `DM` object 1321 . dim - The spatial dimension 1322 . simplex - `PETSC_TRUE` for simplices, `PETSC_FALSE` for tensor cells 1323 . faces - Number of faces per dimension, or `NULL` for (1,) in 1D and (2, 2) in 2D and (1, 1, 1) in 3D 1324 . lower - The lower left corner, or `NULL` for (0, 0, 0) 1325 . upper - The upper right corner, or `NULL` for (1, 1, 1) 1326 . periodicity - The boundary type for the X,Y,Z direction, or `NULL` for `DM_BOUNDARY_NONE` 1327 - interpolate - Flag to create intermediate mesh pieces (edges, faces) 1328 1329 Output Parameter: 1330 . dm - The `DM` object 1331 1332 Level: beginner 1333 1334 Note: 1335 To customize this mesh using options, use 1336 .vb 1337 DMCreate(comm, &dm); 1338 DMSetType(dm, DMPLEX); 1339 DMSetFromOptions(dm); 1340 .ve 1341 and use the options in `DMSetFromOptions()`. 1342 1343 Here is the numbering returned for 2 faces in each direction for tensor cells\: 1344 .vb 1345 10---17---11---18----12 1346 | | | 1347 | | | 1348 20 2 22 3 24 1349 | | | 1350 | | | 1351 7---15----8---16----9 1352 | | | 1353 | | | 1354 19 0 21 1 23 1355 | | | 1356 | | | 1357 4---13----5---14----6 1358 .ve 1359 and for simplicial cells 1360 .vb 1361 14----8---15----9----16 1362 |\ 5 |\ 7 | 1363 | \ | \ | 1364 13 2 14 3 15 1365 | 4 \ | 6 \ | 1366 | \ | \ | 1367 11----6---12----7----13 1368 |\ |\ | 1369 | \ 1 | \ 3 | 1370 10 0 11 1 12 1371 | 0 \ | 2 \ | 1372 | \ | \ | 1373 8----4----9----5----10 1374 .ve 1375 1376 .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMSetFromOptions()`, `DMPlexCreateFromFile()`, `DMPlexCreateHexCylinderMesh()`, `DMSetType()`, `DMCreate()` 1377 @*/ 1378 PetscErrorCode DMPlexCreateBoxMesh(MPI_Comm comm, PetscInt dim, PetscBool simplex, const PetscInt faces[], const PetscReal lower[], const PetscReal upper[], const DMBoundaryType periodicity[], PetscBool interpolate, DM *dm) 1379 { 1380 PetscInt fac[3] = {1, 1, 1}; 1381 PetscReal low[3] = {0, 0, 0}; 1382 PetscReal upp[3] = {1, 1, 1}; 1383 DMBoundaryType bdt[3] = {DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE}; 1384 1385 PetscFunctionBegin; 1386 PetscCall(DMCreate(comm, dm)); 1387 PetscCall(DMSetType(*dm, DMPLEX)); 1388 PetscCall(DMPlexCreateBoxMesh_Internal(*dm, DM_SHAPE_BOX, dim, simplex, faces ? faces : fac, lower ? lower : low, upper ? upper : upp, periodicity ? periodicity : bdt, interpolate)); 1389 if (periodicity) PetscCall(DMLocalizeCoordinates(*dm)); 1390 PetscFunctionReturn(PETSC_SUCCESS); 1391 } 1392 1393 static PetscErrorCode DMPlexCreateWedgeBoxMesh_Internal(DM dm, const PetscInt faces[], const PetscReal lower[], const PetscReal upper[], const DMBoundaryType periodicity[]) 1394 { 1395 DM bdm, vol; 1396 PetscInt i; 1397 1398 PetscFunctionBegin; 1399 // TODO Now we can support periodicity 1400 for (i = 0; i < 3; ++i) PetscCheck(periodicity[i] == DM_BOUNDARY_NONE, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Periodicity not yet supported"); 1401 PetscCall(DMCreate(PetscObjectComm((PetscObject)dm), &bdm)); 1402 PetscCall(DMSetType(bdm, DMPLEX)); 1403 PetscCall(DMSetDimension(bdm, 2)); 1404 PetscCall(PetscLogEventBegin(DMPLEX_Generate, bdm, 0, 0, 0)); 1405 PetscCall(DMPlexCreateBoxMesh_Simplex_Internal(bdm, 2, faces, lower, upper, periodicity, PETSC_TRUE)); 1406 PetscCall(DMPlexExtrude(bdm, faces[2], upper[2] - lower[2], PETSC_TRUE, PETSC_FALSE, PETSC_FALSE, NULL, NULL, &vol)); 1407 PetscCall(PetscLogEventEnd(DMPLEX_Generate, bdm, 0, 0, 0)); 1408 PetscCall(DMDestroy(&bdm)); 1409 PetscCall(DMPlexReplace_Internal(dm, &vol)); 1410 if (lower[2] != 0.0) { 1411 Vec v; 1412 PetscScalar *x; 1413 PetscInt cDim, n; 1414 1415 PetscCall(DMGetCoordinatesLocal(dm, &v)); 1416 PetscCall(VecGetBlockSize(v, &cDim)); 1417 PetscCall(VecGetLocalSize(v, &n)); 1418 PetscCall(VecGetArray(v, &x)); 1419 x += cDim; 1420 for (i = 0; i < n; i += cDim) x[i] += lower[2]; 1421 PetscCall(VecRestoreArray(v, &x)); 1422 PetscCall(DMSetCoordinatesLocal(dm, v)); 1423 } 1424 PetscFunctionReturn(PETSC_SUCCESS); 1425 } 1426 1427 /*@ 1428 DMPlexCreateWedgeBoxMesh - Creates a 3-D mesh tesselating the (x,y) plane and extruding in the third direction using wedge cells. 1429 1430 Collective 1431 1432 Input Parameters: 1433 + comm - The communicator for the `DM` object 1434 . faces - Number of faces per dimension, or `NULL` for (1, 1, 1) 1435 . lower - The lower left corner, or `NULL` for (0, 0, 0) 1436 . upper - The upper right corner, or `NULL` for (1, 1, 1) 1437 . periodicity - The boundary type for the X,Y,Z direction, or `NULL` for `DM_BOUNDARY_NONE` 1438 . orderHeight - If `PETSC_TRUE`, orders the extruded cells in the height first. Otherwise, orders the cell on the layers first 1439 - interpolate - Flag to create intermediate mesh pieces (edges, faces) 1440 1441 Output Parameter: 1442 . dm - The `DM` object 1443 1444 Level: beginner 1445 1446 .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexCreateHexCylinderMesh()`, `DMPlexCreateWedgeCylinderMesh()`, `DMExtrude()`, `DMPlexCreateBoxMesh()`, `DMSetType()`, `DMCreate()` 1447 @*/ 1448 PetscErrorCode DMPlexCreateWedgeBoxMesh(MPI_Comm comm, const PetscInt faces[], const PetscReal lower[], const PetscReal upper[], const DMBoundaryType periodicity[], PetscBool orderHeight, PetscBool interpolate, DM *dm) 1449 { 1450 PetscInt fac[3] = {1, 1, 1}; 1451 PetscReal low[3] = {0, 0, 0}; 1452 PetscReal upp[3] = {1, 1, 1}; 1453 DMBoundaryType bdt[3] = {DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE}; 1454 1455 PetscFunctionBegin; 1456 PetscCall(DMCreate(comm, dm)); 1457 PetscCall(DMSetType(*dm, DMPLEX)); 1458 PetscCall(DMPlexCreateWedgeBoxMesh_Internal(*dm, faces ? faces : fac, lower ? lower : low, upper ? upper : upp, periodicity ? periodicity : bdt)); 1459 if (!interpolate) { 1460 DM udm; 1461 1462 PetscCall(DMPlexUninterpolate(*dm, &udm)); 1463 PetscCall(DMPlexReplace_Internal(*dm, &udm)); 1464 } 1465 if (periodicity) PetscCall(DMLocalizeCoordinates(*dm)); 1466 PetscFunctionReturn(PETSC_SUCCESS); 1467 } 1468 1469 /* 1470 DMPlexTensorPointLexicographic_Private - Returns all tuples of size 'len' with nonnegative integers that are all less than or equal to 'max' for that dimension. 1471 1472 Input Parameters: 1473 + len - The length of the tuple 1474 . max - The maximum for each dimension, so values are in [0, max) 1475 - tup - A tuple of length len+1: tup[len] > 0 indicates a stopping condition 1476 1477 Output Parameter: 1478 . tup - A tuple of `len` integers whose entries are at most `max` 1479 1480 Level: developer 1481 1482 Note: 1483 Ordering is lexicographic with lowest index as least significant in ordering. 1484 e.g. for len == 2 and max == 2, this will return, in order, {0,0}, {1,0}, {2,0}, {0,1}, {1,1}, {2,1}, {0,2}, {1,2}, {2,2}. 1485 1486 .seealso: PetscDualSpaceTensorPointLexicographic_Internal(), PetscDualSpaceLatticePointLexicographic_Internal() 1487 */ 1488 static PetscErrorCode DMPlexTensorPointLexicographic_Private(PetscInt len, const PetscInt max[], PetscInt tup[]) 1489 { 1490 PetscInt i; 1491 1492 PetscFunctionBegin; 1493 for (i = 0; i < len; ++i) { 1494 if (tup[i] < max[i] - 1) { 1495 break; 1496 } else { 1497 tup[i] = 0; 1498 } 1499 } 1500 if (i == len) tup[i - 1] = max[i - 1]; 1501 else ++tup[i]; 1502 PetscFunctionReturn(PETSC_SUCCESS); 1503 } 1504 1505 static PetscInt TupleToIndex_Private(PetscInt len, const PetscInt max[], const PetscInt tup[]) 1506 { 1507 PetscInt i, idx = tup[len - 1]; 1508 1509 for (i = len - 2; i >= 0; --i) { 1510 idx *= max[i]; 1511 idx += tup[i]; 1512 } 1513 return idx; 1514 } 1515 1516 static PetscErrorCode DestroyExtent_Private(void *extent) 1517 { 1518 return PetscFree(extent); 1519 } 1520 1521 static PetscErrorCode DMPlexCreateHypercubicMesh_Internal(DM dm, PetscInt dim, const PetscReal lower[], const PetscReal upper[], const PetscInt edges[], const DMBoundaryType bd[]) 1522 { 1523 Vec coordinates; 1524 PetscSection coordSection; 1525 DMLabel cutLabel = NULL; 1526 PetscBool cutMarker = PETSC_FALSE; 1527 PetscBool periodic = PETSC_FALSE; 1528 PetscInt numCells = 1, c; 1529 PetscInt numVertices = 1, v; 1530 PetscScalar *coords; 1531 PetscInt *vertices, *vert, *vtmp, *supp, cone[2]; 1532 PetscInt d, e, cell = 0, coordSize; 1533 PetscMPIInt rank; 1534 1535 PetscFunctionBegin; 1536 PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank)); 1537 PetscCall(DMSetDimension(dm, dim)); 1538 PetscCall(PetscCalloc4(dim, &vertices, dim, &vert, dim, &vtmp, 2 * dim, &supp)); 1539 PetscCall(DMCreateLabel(dm, "marker")); 1540 PetscCall(PetscOptionsGetBool(((PetscObject)dm)->options, ((PetscObject)dm)->prefix, "-dm_plex_periodic_cut", &cutMarker, NULL)); 1541 for (d = 0; d < dim; ++d) periodic = (periodic || bd[d] == DM_BOUNDARY_PERIODIC) ? PETSC_TRUE : PETSC_FALSE; 1542 if (periodic && cutMarker) { 1543 PetscCall(DMCreateLabel(dm, "periodic_cut")); 1544 PetscCall(DMGetLabel(dm, "periodic_cut", &cutLabel)); 1545 } 1546 for (d = 0; d < dim; ++d) PetscCheck(bd[d] == DM_BOUNDARY_PERIODIC, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Hypercubic mesh must be periodic now"); 1547 for (d = 0; d < dim; ++d) { 1548 vertices[d] = edges[d]; 1549 numVertices *= vertices[d]; 1550 } 1551 numCells = numVertices * dim; 1552 PetscCall(DMPlexSetChart(dm, 0, numCells + numVertices)); 1553 for (c = 0; c < numCells; ++c) PetscCall(DMPlexSetConeSize(dm, c, 2)); 1554 for (v = numCells; v < numCells + numVertices; ++v) PetscCall(DMPlexSetSupportSize(dm, v, 2 * dim)); 1555 /* TODO Loop over boundary and reset support sizes */ 1556 PetscCall(DMSetUp(dm)); /* Allocate space for cones and supports */ 1557 /* Build cell cones and vertex supports */ 1558 PetscCall(DMCreateLabel(dm, "celltype")); 1559 while (vert[dim - 1] < vertices[dim - 1]) { 1560 const PetscInt vertex = TupleToIndex_Private(dim, vertices, vert) + numCells; 1561 PetscInt s = 0; 1562 1563 PetscCall(PetscPrintf(PETSC_COMM_SELF, "Vertex %" PetscInt_FMT ":", vertex)); 1564 for (d = 0; d < dim; ++d) PetscCall(PetscPrintf(PETSC_COMM_SELF, " %" PetscInt_FMT, vert[d])); 1565 PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n")); 1566 PetscCall(DMPlexSetCellType(dm, vertex, DM_POLYTOPE_POINT)); 1567 for (d = 0; d < dim; ++d) { 1568 for (e = 0; e < dim; ++e) vtmp[e] = vert[e]; 1569 vtmp[d] = (vert[d] + 1) % vertices[d]; 1570 cone[0] = vertex; 1571 cone[1] = TupleToIndex_Private(dim, vertices, vtmp) + numCells; 1572 PetscCall(PetscPrintf(PETSC_COMM_SELF, " Vertex %" PetscInt_FMT ":", cone[1])); 1573 for (e = 0; e < dim; ++e) PetscCall(PetscPrintf(PETSC_COMM_SELF, " %" PetscInt_FMT, vtmp[e])); 1574 PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n")); 1575 PetscCall(DMPlexSetCone(dm, cell, cone)); 1576 PetscCall(DMPlexSetCellType(dm, cell, DM_POLYTOPE_SEGMENT)); 1577 PetscCall(PetscPrintf(PETSC_COMM_SELF, " Edge %" PetscInt_FMT " (%" PetscInt_FMT " %" PetscInt_FMT ")\n", cell, cone[0], cone[1])); 1578 ++cell; 1579 } 1580 for (d = 0; d < dim; ++d) { 1581 for (e = 0; e < dim; ++e) vtmp[e] = vert[e]; 1582 vtmp[d] = (vert[d] + vertices[d] - 1) % vertices[d]; 1583 supp[s++] = TupleToIndex_Private(dim, vertices, vtmp) * dim + d; 1584 supp[s++] = (vertex - numCells) * dim + d; 1585 PetscCall(DMPlexSetSupport(dm, vertex, supp)); 1586 } 1587 PetscCall(DMPlexTensorPointLexicographic_Private(dim, vertices, vert)); 1588 } 1589 PetscCall(DMPlexStratify(dm)); 1590 /* Build coordinates */ 1591 PetscCall(DMGetCoordinateSection(dm, &coordSection)); 1592 PetscCall(PetscSectionSetNumFields(coordSection, 1)); 1593 PetscCall(PetscSectionSetFieldComponents(coordSection, 0, dim)); 1594 PetscCall(PetscSectionSetChart(coordSection, numCells, numCells + numVertices)); 1595 for (v = numCells; v < numCells + numVertices; ++v) { 1596 PetscCall(PetscSectionSetDof(coordSection, v, dim)); 1597 PetscCall(PetscSectionSetFieldDof(coordSection, v, 0, dim)); 1598 } 1599 PetscCall(PetscSectionSetUp(coordSection)); 1600 PetscCall(PetscSectionGetStorageSize(coordSection, &coordSize)); 1601 PetscCall(VecCreate(PETSC_COMM_SELF, &coordinates)); 1602 PetscCall(PetscObjectSetName((PetscObject)coordinates, "coordinates")); 1603 PetscCall(VecSetSizes(coordinates, coordSize, PETSC_DETERMINE)); 1604 PetscCall(VecSetBlockSize(coordinates, dim)); 1605 PetscCall(VecSetType(coordinates, VECSTANDARD)); 1606 PetscCall(VecGetArray(coordinates, &coords)); 1607 for (d = 0; d < dim; ++d) vert[d] = 0; 1608 while (vert[dim - 1] < vertices[dim - 1]) { 1609 const PetscInt vertex = TupleToIndex_Private(dim, vertices, vert); 1610 1611 for (d = 0; d < dim; ++d) coords[vertex * dim + d] = lower[d] + ((upper[d] - lower[d]) / vertices[d]) * vert[d]; 1612 PetscCall(PetscPrintf(PETSC_COMM_SELF, "Vertex %" PetscInt_FMT ":", vertex)); 1613 for (d = 0; d < dim; ++d) PetscCall(PetscPrintf(PETSC_COMM_SELF, " %" PetscInt_FMT, vert[d])); 1614 for (d = 0; d < dim; ++d) PetscCall(PetscPrintf(PETSC_COMM_SELF, " %g", (double)PetscRealPart(coords[vertex * dim + d]))); 1615 PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n")); 1616 PetscCall(DMPlexTensorPointLexicographic_Private(dim, vertices, vert)); 1617 } 1618 PetscCall(VecRestoreArray(coordinates, &coords)); 1619 PetscCall(DMSetCoordinatesLocal(dm, coordinates)); 1620 PetscCall(VecDestroy(&coordinates)); 1621 PetscCall(PetscFree4(vertices, vert, vtmp, supp)); 1622 //PetscCall(DMSetPeriodicity(dm, NULL, lower, upper)); 1623 // Attach the extent 1624 { 1625 PetscContainer c; 1626 PetscInt *extent; 1627 1628 PetscCall(PetscMalloc1(dim, &extent)); 1629 for (PetscInt d = 0; d < dim; ++d) extent[d] = edges[d]; 1630 PetscCall(PetscContainerCreate(PETSC_COMM_SELF, &c)); 1631 PetscCall(PetscContainerSetUserDestroy(c, DestroyExtent_Private)); 1632 PetscCall(PetscContainerSetPointer(c, extent)); 1633 PetscCall(PetscObjectCompose((PetscObject)dm, "_extent", (PetscObject)c)); 1634 PetscCall(PetscContainerDestroy(&c)); 1635 } 1636 PetscFunctionReturn(PETSC_SUCCESS); 1637 } 1638 1639 /*@C 1640 DMPlexCreateHypercubicMesh - Creates a periodic mesh on the tensor product of unit intervals using only vertices and edges. 1641 1642 Collective 1643 1644 Input Parameters: 1645 + comm - The communicator for the DM object 1646 . dim - The spatial dimension 1647 . edges - Number of edges per dimension, or `NULL` for (1,) in 1D and (2, 2) in 2D and (1, 1, 1) in 3D 1648 . lower - The lower left corner, or `NULL` for (0, 0, 0) 1649 - upper - The upper right corner, or `NULL` for (1, 1, 1) 1650 1651 Output Parameter: 1652 . dm - The DM object 1653 1654 Level: beginner 1655 1656 Note: 1657 If you want to customize this mesh using options, you just need to 1658 .vb 1659 DMCreate(comm, &dm); 1660 DMSetType(dm, DMPLEX); 1661 DMSetFromOptions(dm); 1662 .ve 1663 and use the options on the `DMSetFromOptions()` page. 1664 1665 The vertices are numbered is lexicographic order, and the dim edges exiting a vertex in the positive orthant are number consecutively, 1666 .vb 1667 18--0-19--2-20--4-18 1668 | | | | 1669 13 15 17 13 1670 | | | | 1671 24-12-25-14-26-16-24 1672 | | | | 1673 7 9 11 7 1674 | | | | 1675 21--6-22--8-23-10-21 1676 | | | | 1677 1 3 5 1 1678 | | | | 1679 18--0-19--2-20--4-18 1680 .ve 1681 1682 .seealso: `DMSetFromOptions()`, `DMPlexCreateFromFile()`, `DMPlexCreateHexCylinderMesh()`, `DMSetType()`, `DMCreate()` 1683 @*/ 1684 PetscErrorCode DMPlexCreateHypercubicMesh(MPI_Comm comm, PetscInt dim, const PetscInt edges[], const PetscReal lower[], const PetscReal upper[], DM *dm) 1685 { 1686 PetscInt *edg; 1687 PetscReal *low, *upp; 1688 DMBoundaryType *bdt; 1689 PetscInt d; 1690 1691 PetscFunctionBegin; 1692 PetscCall(DMCreate(comm, dm)); 1693 PetscCall(DMSetType(*dm, DMPLEX)); 1694 PetscCall(PetscMalloc4(dim, &edg, dim, &low, dim, &upp, dim, &bdt)); 1695 for (d = 0; d < dim; ++d) { 1696 edg[d] = edges ? edges[d] : 1; 1697 low[d] = lower ? lower[d] : 0.; 1698 upp[d] = upper ? upper[d] : 1.; 1699 bdt[d] = DM_BOUNDARY_PERIODIC; 1700 } 1701 PetscCall(DMPlexCreateHypercubicMesh_Internal(*dm, dim, low, upp, edg, bdt)); 1702 PetscCall(PetscFree4(edg, low, upp, bdt)); 1703 PetscFunctionReturn(PETSC_SUCCESS); 1704 } 1705 1706 /*@C 1707 DMPlexSetOptionsPrefix - Sets the prefix used for searching for all `DM` options in the database. 1708 1709 Logically Collective 1710 1711 Input Parameters: 1712 + dm - the `DM` context 1713 - prefix - the prefix to prepend to all option names 1714 1715 Level: advanced 1716 1717 Note: 1718 A hyphen (-) must NOT be given at the beginning of the prefix name. 1719 The first character of all runtime options is AUTOMATICALLY the hyphen. 1720 1721 .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `SNESSetFromOptions()` 1722 @*/ 1723 PetscErrorCode DMPlexSetOptionsPrefix(DM dm, const char prefix[]) 1724 { 1725 DM_Plex *mesh = (DM_Plex *)dm->data; 1726 1727 PetscFunctionBegin; 1728 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1729 PetscCall(PetscObjectSetOptionsPrefix((PetscObject)dm, prefix)); 1730 PetscCall(PetscObjectSetOptionsPrefix((PetscObject)mesh->partitioner, prefix)); 1731 PetscFunctionReturn(PETSC_SUCCESS); 1732 } 1733 1734 /* Remap geometry to cylinder 1735 TODO: This only works for a single refinement, then it is broken 1736 1737 Interior square: Linear interpolation is correct 1738 The other cells all have vertices on rays from the origin. We want to uniformly expand the spacing 1739 such that the last vertex is on the unit circle. So the closest and farthest vertices are at distance 1740 1741 phi = arctan(y/x) 1742 d_close = sqrt(1/8 + 1/4 sin^2(phi)) 1743 d_far = sqrt(1/2 + sin^2(phi)) 1744 1745 so we remap them using 1746 1747 x_new = x_close + (x - x_close) (1 - d_close) / (d_far - d_close) 1748 y_new = y_close + (y - y_close) (1 - d_close) / (d_far - d_close) 1749 1750 If pi/4 < phi < 3pi/4 or -3pi/4 < phi < -pi/4, then we switch x and y. 1751 */ 1752 static void snapToCylinder(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 1753 { 1754 const PetscReal dis = 1.0 / PetscSqrtReal(2.0); 1755 const PetscReal ds2 = 0.5 * dis; 1756 1757 if ((PetscAbsScalar(u[0]) <= ds2) && (PetscAbsScalar(u[1]) <= ds2)) { 1758 f0[0] = u[0]; 1759 f0[1] = u[1]; 1760 } else { 1761 PetscReal phi, sinp, cosp, dc, df, x, y, xc, yc; 1762 1763 x = PetscRealPart(u[0]); 1764 y = PetscRealPart(u[1]); 1765 phi = PetscAtan2Real(y, x); 1766 sinp = PetscSinReal(phi); 1767 cosp = PetscCosReal(phi); 1768 if ((PetscAbsReal(phi) > PETSC_PI / 4.0) && (PetscAbsReal(phi) < 3.0 * PETSC_PI / 4.0)) { 1769 dc = PetscAbsReal(ds2 / sinp); 1770 df = PetscAbsReal(dis / sinp); 1771 xc = ds2 * x / PetscAbsReal(y); 1772 yc = ds2 * PetscSignReal(y); 1773 } else { 1774 dc = PetscAbsReal(ds2 / cosp); 1775 df = PetscAbsReal(dis / cosp); 1776 xc = ds2 * PetscSignReal(x); 1777 yc = ds2 * y / PetscAbsReal(x); 1778 } 1779 f0[0] = xc + (u[0] - xc) * (1.0 - dc) / (df - dc); 1780 f0[1] = yc + (u[1] - yc) * (1.0 - dc) / (df - dc); 1781 } 1782 f0[2] = u[2]; 1783 } 1784 1785 static PetscErrorCode DMPlexCreateHexCylinderMesh_Internal(DM dm, DMBoundaryType periodicZ) 1786 { 1787 const PetscInt dim = 3; 1788 PetscInt numCells, numVertices; 1789 PetscMPIInt rank; 1790 1791 PetscFunctionBegin; 1792 PetscCall(PetscLogEventBegin(DMPLEX_Generate, dm, 0, 0, 0)); 1793 PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank)); 1794 PetscCall(DMSetDimension(dm, dim)); 1795 /* Create topology */ 1796 { 1797 PetscInt cone[8], c; 1798 1799 numCells = rank == 0 ? 5 : 0; 1800 numVertices = rank == 0 ? 16 : 0; 1801 if (periodicZ == DM_BOUNDARY_PERIODIC) { 1802 numCells *= 3; 1803 numVertices = rank == 0 ? 24 : 0; 1804 } 1805 PetscCall(DMPlexSetChart(dm, 0, numCells + numVertices)); 1806 for (c = 0; c < numCells; c++) PetscCall(DMPlexSetConeSize(dm, c, 8)); 1807 PetscCall(DMSetUp(dm)); 1808 if (rank == 0) { 1809 if (periodicZ == DM_BOUNDARY_PERIODIC) { 1810 cone[0] = 15; 1811 cone[1] = 18; 1812 cone[2] = 17; 1813 cone[3] = 16; 1814 cone[4] = 31; 1815 cone[5] = 32; 1816 cone[6] = 33; 1817 cone[7] = 34; 1818 PetscCall(DMPlexSetCone(dm, 0, cone)); 1819 cone[0] = 16; 1820 cone[1] = 17; 1821 cone[2] = 24; 1822 cone[3] = 23; 1823 cone[4] = 32; 1824 cone[5] = 36; 1825 cone[6] = 37; 1826 cone[7] = 33; /* 22 25 26 21 */ 1827 PetscCall(DMPlexSetCone(dm, 1, cone)); 1828 cone[0] = 18; 1829 cone[1] = 27; 1830 cone[2] = 24; 1831 cone[3] = 17; 1832 cone[4] = 34; 1833 cone[5] = 33; 1834 cone[6] = 37; 1835 cone[7] = 38; 1836 PetscCall(DMPlexSetCone(dm, 2, cone)); 1837 cone[0] = 29; 1838 cone[1] = 27; 1839 cone[2] = 18; 1840 cone[3] = 15; 1841 cone[4] = 35; 1842 cone[5] = 31; 1843 cone[6] = 34; 1844 cone[7] = 38; 1845 PetscCall(DMPlexSetCone(dm, 3, cone)); 1846 cone[0] = 29; 1847 cone[1] = 15; 1848 cone[2] = 16; 1849 cone[3] = 23; 1850 cone[4] = 35; 1851 cone[5] = 36; 1852 cone[6] = 32; 1853 cone[7] = 31; 1854 PetscCall(DMPlexSetCone(dm, 4, cone)); 1855 1856 cone[0] = 31; 1857 cone[1] = 34; 1858 cone[2] = 33; 1859 cone[3] = 32; 1860 cone[4] = 19; 1861 cone[5] = 22; 1862 cone[6] = 21; 1863 cone[7] = 20; 1864 PetscCall(DMPlexSetCone(dm, 5, cone)); 1865 cone[0] = 32; 1866 cone[1] = 33; 1867 cone[2] = 37; 1868 cone[3] = 36; 1869 cone[4] = 22; 1870 cone[5] = 25; 1871 cone[6] = 26; 1872 cone[7] = 21; 1873 PetscCall(DMPlexSetCone(dm, 6, cone)); 1874 cone[0] = 34; 1875 cone[1] = 38; 1876 cone[2] = 37; 1877 cone[3] = 33; 1878 cone[4] = 20; 1879 cone[5] = 21; 1880 cone[6] = 26; 1881 cone[7] = 28; 1882 PetscCall(DMPlexSetCone(dm, 7, cone)); 1883 cone[0] = 35; 1884 cone[1] = 38; 1885 cone[2] = 34; 1886 cone[3] = 31; 1887 cone[4] = 30; 1888 cone[5] = 19; 1889 cone[6] = 20; 1890 cone[7] = 28; 1891 PetscCall(DMPlexSetCone(dm, 8, cone)); 1892 cone[0] = 35; 1893 cone[1] = 31; 1894 cone[2] = 32; 1895 cone[3] = 36; 1896 cone[4] = 30; 1897 cone[5] = 25; 1898 cone[6] = 22; 1899 cone[7] = 19; 1900 PetscCall(DMPlexSetCone(dm, 9, cone)); 1901 1902 cone[0] = 19; 1903 cone[1] = 20; 1904 cone[2] = 21; 1905 cone[3] = 22; 1906 cone[4] = 15; 1907 cone[5] = 16; 1908 cone[6] = 17; 1909 cone[7] = 18; 1910 PetscCall(DMPlexSetCone(dm, 10, cone)); 1911 cone[0] = 22; 1912 cone[1] = 21; 1913 cone[2] = 26; 1914 cone[3] = 25; 1915 cone[4] = 16; 1916 cone[5] = 23; 1917 cone[6] = 24; 1918 cone[7] = 17; 1919 PetscCall(DMPlexSetCone(dm, 11, cone)); 1920 cone[0] = 20; 1921 cone[1] = 28; 1922 cone[2] = 26; 1923 cone[3] = 21; 1924 cone[4] = 18; 1925 cone[5] = 17; 1926 cone[6] = 24; 1927 cone[7] = 27; 1928 PetscCall(DMPlexSetCone(dm, 12, cone)); 1929 cone[0] = 30; 1930 cone[1] = 28; 1931 cone[2] = 20; 1932 cone[3] = 19; 1933 cone[4] = 29; 1934 cone[5] = 15; 1935 cone[6] = 18; 1936 cone[7] = 27; 1937 PetscCall(DMPlexSetCone(dm, 13, cone)); 1938 cone[0] = 30; 1939 cone[1] = 19; 1940 cone[2] = 22; 1941 cone[3] = 25; 1942 cone[4] = 29; 1943 cone[5] = 23; 1944 cone[6] = 16; 1945 cone[7] = 15; 1946 PetscCall(DMPlexSetCone(dm, 14, cone)); 1947 } else { 1948 cone[0] = 5; 1949 cone[1] = 8; 1950 cone[2] = 7; 1951 cone[3] = 6; 1952 cone[4] = 9; 1953 cone[5] = 12; 1954 cone[6] = 11; 1955 cone[7] = 10; 1956 PetscCall(DMPlexSetCone(dm, 0, cone)); 1957 cone[0] = 6; 1958 cone[1] = 7; 1959 cone[2] = 14; 1960 cone[3] = 13; 1961 cone[4] = 12; 1962 cone[5] = 15; 1963 cone[6] = 16; 1964 cone[7] = 11; 1965 PetscCall(DMPlexSetCone(dm, 1, cone)); 1966 cone[0] = 8; 1967 cone[1] = 17; 1968 cone[2] = 14; 1969 cone[3] = 7; 1970 cone[4] = 10; 1971 cone[5] = 11; 1972 cone[6] = 16; 1973 cone[7] = 18; 1974 PetscCall(DMPlexSetCone(dm, 2, cone)); 1975 cone[0] = 19; 1976 cone[1] = 17; 1977 cone[2] = 8; 1978 cone[3] = 5; 1979 cone[4] = 20; 1980 cone[5] = 9; 1981 cone[6] = 10; 1982 cone[7] = 18; 1983 PetscCall(DMPlexSetCone(dm, 3, cone)); 1984 cone[0] = 19; 1985 cone[1] = 5; 1986 cone[2] = 6; 1987 cone[3] = 13; 1988 cone[4] = 20; 1989 cone[5] = 15; 1990 cone[6] = 12; 1991 cone[7] = 9; 1992 PetscCall(DMPlexSetCone(dm, 4, cone)); 1993 } 1994 } 1995 PetscCall(DMPlexSymmetrize(dm)); 1996 PetscCall(DMPlexStratify(dm)); 1997 } 1998 /* Create cube geometry */ 1999 { 2000 Vec coordinates; 2001 PetscSection coordSection; 2002 PetscScalar *coords; 2003 PetscInt coordSize, v; 2004 const PetscReal dis = 1.0 / PetscSqrtReal(2.0); 2005 const PetscReal ds2 = dis / 2.0; 2006 2007 /* Build coordinates */ 2008 PetscCall(DMGetCoordinateSection(dm, &coordSection)); 2009 PetscCall(PetscSectionSetNumFields(coordSection, 1)); 2010 PetscCall(PetscSectionSetFieldComponents(coordSection, 0, dim)); 2011 PetscCall(PetscSectionSetChart(coordSection, numCells, numCells + numVertices)); 2012 for (v = numCells; v < numCells + numVertices; ++v) { 2013 PetscCall(PetscSectionSetDof(coordSection, v, dim)); 2014 PetscCall(PetscSectionSetFieldDof(coordSection, v, 0, dim)); 2015 } 2016 PetscCall(PetscSectionSetUp(coordSection)); 2017 PetscCall(PetscSectionGetStorageSize(coordSection, &coordSize)); 2018 PetscCall(VecCreate(PETSC_COMM_SELF, &coordinates)); 2019 PetscCall(PetscObjectSetName((PetscObject)coordinates, "coordinates")); 2020 PetscCall(VecSetSizes(coordinates, coordSize, PETSC_DETERMINE)); 2021 PetscCall(VecSetBlockSize(coordinates, dim)); 2022 PetscCall(VecSetType(coordinates, VECSTANDARD)); 2023 PetscCall(VecGetArray(coordinates, &coords)); 2024 if (rank == 0) { 2025 coords[0 * dim + 0] = -ds2; 2026 coords[0 * dim + 1] = -ds2; 2027 coords[0 * dim + 2] = 0.0; 2028 coords[1 * dim + 0] = ds2; 2029 coords[1 * dim + 1] = -ds2; 2030 coords[1 * dim + 2] = 0.0; 2031 coords[2 * dim + 0] = ds2; 2032 coords[2 * dim + 1] = ds2; 2033 coords[2 * dim + 2] = 0.0; 2034 coords[3 * dim + 0] = -ds2; 2035 coords[3 * dim + 1] = ds2; 2036 coords[3 * dim + 2] = 0.0; 2037 coords[4 * dim + 0] = -ds2; 2038 coords[4 * dim + 1] = -ds2; 2039 coords[4 * dim + 2] = 1.0; 2040 coords[5 * dim + 0] = -ds2; 2041 coords[5 * dim + 1] = ds2; 2042 coords[5 * dim + 2] = 1.0; 2043 coords[6 * dim + 0] = ds2; 2044 coords[6 * dim + 1] = ds2; 2045 coords[6 * dim + 2] = 1.0; 2046 coords[7 * dim + 0] = ds2; 2047 coords[7 * dim + 1] = -ds2; 2048 coords[7 * dim + 2] = 1.0; 2049 coords[8 * dim + 0] = dis; 2050 coords[8 * dim + 1] = -dis; 2051 coords[8 * dim + 2] = 0.0; 2052 coords[9 * dim + 0] = dis; 2053 coords[9 * dim + 1] = dis; 2054 coords[9 * dim + 2] = 0.0; 2055 coords[10 * dim + 0] = dis; 2056 coords[10 * dim + 1] = -dis; 2057 coords[10 * dim + 2] = 1.0; 2058 coords[11 * dim + 0] = dis; 2059 coords[11 * dim + 1] = dis; 2060 coords[11 * dim + 2] = 1.0; 2061 coords[12 * dim + 0] = -dis; 2062 coords[12 * dim + 1] = dis; 2063 coords[12 * dim + 2] = 0.0; 2064 coords[13 * dim + 0] = -dis; 2065 coords[13 * dim + 1] = dis; 2066 coords[13 * dim + 2] = 1.0; 2067 coords[14 * dim + 0] = -dis; 2068 coords[14 * dim + 1] = -dis; 2069 coords[14 * dim + 2] = 0.0; 2070 coords[15 * dim + 0] = -dis; 2071 coords[15 * dim + 1] = -dis; 2072 coords[15 * dim + 2] = 1.0; 2073 if (periodicZ == DM_BOUNDARY_PERIODIC) { 2074 /* 15 31 19 */ coords[16 * dim + 0] = -ds2; 2075 coords[16 * dim + 1] = -ds2; 2076 coords[16 * dim + 2] = 0.5; 2077 /* 16 32 22 */ coords[17 * dim + 0] = ds2; 2078 coords[17 * dim + 1] = -ds2; 2079 coords[17 * dim + 2] = 0.5; 2080 /* 17 33 21 */ coords[18 * dim + 0] = ds2; 2081 coords[18 * dim + 1] = ds2; 2082 coords[18 * dim + 2] = 0.5; 2083 /* 18 34 20 */ coords[19 * dim + 0] = -ds2; 2084 coords[19 * dim + 1] = ds2; 2085 coords[19 * dim + 2] = 0.5; 2086 /* 29 35 30 */ coords[20 * dim + 0] = -dis; 2087 coords[20 * dim + 1] = -dis; 2088 coords[20 * dim + 2] = 0.5; 2089 /* 23 36 25 */ coords[21 * dim + 0] = dis; 2090 coords[21 * dim + 1] = -dis; 2091 coords[21 * dim + 2] = 0.5; 2092 /* 24 37 26 */ coords[22 * dim + 0] = dis; 2093 coords[22 * dim + 1] = dis; 2094 coords[22 * dim + 2] = 0.5; 2095 /* 27 38 28 */ coords[23 * dim + 0] = -dis; 2096 coords[23 * dim + 1] = dis; 2097 coords[23 * dim + 2] = 0.5; 2098 } 2099 } 2100 PetscCall(VecRestoreArray(coordinates, &coords)); 2101 PetscCall(DMSetCoordinatesLocal(dm, coordinates)); 2102 PetscCall(VecDestroy(&coordinates)); 2103 } 2104 /* Create periodicity */ 2105 if (periodicZ == DM_BOUNDARY_PERIODIC || periodicZ == DM_BOUNDARY_TWIST) { 2106 PetscReal L[3] = {-1., -1., 0.}; 2107 PetscReal maxCell[3] = {-1., -1., 0.}; 2108 PetscReal lower[3] = {0.0, 0.0, 0.0}; 2109 PetscReal upper[3] = {1.0, 1.0, 1.5}; 2110 PetscInt numZCells = 3; 2111 2112 L[2] = upper[2] - lower[2]; 2113 maxCell[2] = 1.1 * (L[2] / numZCells); 2114 PetscCall(DMSetPeriodicity(dm, maxCell, lower, L)); 2115 } 2116 { 2117 DM cdm; 2118 PetscDS cds; 2119 PetscScalar c[2] = {1.0, 1.0}; 2120 2121 PetscCall(DMPlexCreateCoordinateSpace(dm, 1, PETSC_TRUE, snapToCylinder)); 2122 PetscCall(DMGetCoordinateDM(dm, &cdm)); 2123 PetscCall(DMGetDS(cdm, &cds)); 2124 PetscCall(PetscDSSetConstants(cds, 2, c)); 2125 } 2126 PetscCall(PetscLogEventEnd(DMPLEX_Generate, dm, 0, 0, 0)); 2127 2128 /* Wait for coordinate creation before doing in-place modification */ 2129 PetscCall(DMPlexInterpolateInPlace_Internal(dm)); 2130 PetscFunctionReturn(PETSC_SUCCESS); 2131 } 2132 2133 /*@ 2134 DMPlexCreateHexCylinderMesh - Creates a mesh on the tensor product of the unit interval with the circle (cylinder) using hexahedra. 2135 2136 Collective 2137 2138 Input Parameters: 2139 + comm - The communicator for the `DM` object 2140 - periodicZ - The boundary type for the Z direction 2141 2142 Output Parameter: 2143 . dm - The `DM` object 2144 2145 Level: beginner 2146 2147 Note: 2148 Here is the output numbering looking from the bottom of the cylinder\: 2149 .vb 2150 17-----14 2151 | | 2152 | 2 | 2153 | | 2154 17-----8-----7-----14 2155 | | | | 2156 | 3 | 0 | 1 | 2157 | | | | 2158 19-----5-----6-----13 2159 | | 2160 | 4 | 2161 | | 2162 19-----13 2163 2164 and up through the top 2165 2166 18-----16 2167 | | 2168 | 2 | 2169 | | 2170 18----10----11-----16 2171 | | | | 2172 | 3 | 0 | 1 | 2173 | | | | 2174 20-----9----12-----15 2175 | | 2176 | 4 | 2177 | | 2178 20-----15 2179 .ve 2180 2181 .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexCreateBoxMesh()`, `DMSetType()`, `DMCreate()` 2182 @*/ 2183 PetscErrorCode DMPlexCreateHexCylinderMesh(MPI_Comm comm, DMBoundaryType periodicZ, DM *dm) 2184 { 2185 PetscFunctionBegin; 2186 PetscAssertPointer(dm, 3); 2187 PetscCall(DMCreate(comm, dm)); 2188 PetscCall(DMSetType(*dm, DMPLEX)); 2189 PetscCall(DMPlexCreateHexCylinderMesh_Internal(*dm, periodicZ)); 2190 PetscFunctionReturn(PETSC_SUCCESS); 2191 } 2192 2193 static PetscErrorCode DMPlexCreateWedgeCylinderMesh_Internal(DM dm, PetscInt n, PetscBool interpolate) 2194 { 2195 const PetscInt dim = 3; 2196 PetscInt numCells, numVertices, v; 2197 PetscMPIInt rank; 2198 2199 PetscFunctionBegin; 2200 PetscCheck(n >= 0, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Number of wedges %" PetscInt_FMT " cannot be negative", n); 2201 PetscCall(PetscLogEventBegin(DMPLEX_Generate, dm, 0, 0, 0)); 2202 PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank)); 2203 PetscCall(DMSetDimension(dm, dim)); 2204 /* Must create the celltype label here so that we do not automatically try to compute the types */ 2205 PetscCall(DMCreateLabel(dm, "celltype")); 2206 /* Create topology */ 2207 { 2208 PetscInt cone[6], c; 2209 2210 numCells = rank == 0 ? n : 0; 2211 numVertices = rank == 0 ? 2 * (n + 1) : 0; 2212 PetscCall(DMPlexSetChart(dm, 0, numCells + numVertices)); 2213 for (c = 0; c < numCells; c++) PetscCall(DMPlexSetConeSize(dm, c, 6)); 2214 PetscCall(DMSetUp(dm)); 2215 for (c = 0; c < numCells; c++) { 2216 cone[0] = c + n * 1; 2217 cone[1] = (c + 1) % n + n * 1; 2218 cone[2] = 0 + 3 * n; 2219 cone[3] = c + n * 2; 2220 cone[4] = (c + 1) % n + n * 2; 2221 cone[5] = 1 + 3 * n; 2222 PetscCall(DMPlexSetCone(dm, c, cone)); 2223 PetscCall(DMPlexSetCellType(dm, c, DM_POLYTOPE_TRI_PRISM_TENSOR)); 2224 } 2225 PetscCall(DMPlexSymmetrize(dm)); 2226 PetscCall(DMPlexStratify(dm)); 2227 } 2228 for (v = numCells; v < numCells + numVertices; ++v) PetscCall(DMPlexSetCellType(dm, v, DM_POLYTOPE_POINT)); 2229 /* Create cylinder geometry */ 2230 { 2231 Vec coordinates; 2232 PetscSection coordSection; 2233 PetscScalar *coords; 2234 PetscInt coordSize, c; 2235 2236 /* Build coordinates */ 2237 PetscCall(DMGetCoordinateSection(dm, &coordSection)); 2238 PetscCall(PetscSectionSetNumFields(coordSection, 1)); 2239 PetscCall(PetscSectionSetFieldComponents(coordSection, 0, dim)); 2240 PetscCall(PetscSectionSetChart(coordSection, numCells, numCells + numVertices)); 2241 for (v = numCells; v < numCells + numVertices; ++v) { 2242 PetscCall(PetscSectionSetDof(coordSection, v, dim)); 2243 PetscCall(PetscSectionSetFieldDof(coordSection, v, 0, dim)); 2244 } 2245 PetscCall(PetscSectionSetUp(coordSection)); 2246 PetscCall(PetscSectionGetStorageSize(coordSection, &coordSize)); 2247 PetscCall(VecCreate(PETSC_COMM_SELF, &coordinates)); 2248 PetscCall(PetscObjectSetName((PetscObject)coordinates, "coordinates")); 2249 PetscCall(VecSetSizes(coordinates, coordSize, PETSC_DETERMINE)); 2250 PetscCall(VecSetBlockSize(coordinates, dim)); 2251 PetscCall(VecSetType(coordinates, VECSTANDARD)); 2252 PetscCall(VecGetArray(coordinates, &coords)); 2253 for (c = 0; c < numCells; c++) { 2254 coords[(c + 0 * n) * dim + 0] = PetscCosReal(2.0 * c * PETSC_PI / n); 2255 coords[(c + 0 * n) * dim + 1] = PetscSinReal(2.0 * c * PETSC_PI / n); 2256 coords[(c + 0 * n) * dim + 2] = 1.0; 2257 coords[(c + 1 * n) * dim + 0] = PetscCosReal(2.0 * c * PETSC_PI / n); 2258 coords[(c + 1 * n) * dim + 1] = PetscSinReal(2.0 * c * PETSC_PI / n); 2259 coords[(c + 1 * n) * dim + 2] = 0.0; 2260 } 2261 if (rank == 0) { 2262 coords[(2 * n + 0) * dim + 0] = 0.0; 2263 coords[(2 * n + 0) * dim + 1] = 0.0; 2264 coords[(2 * n + 0) * dim + 2] = 1.0; 2265 coords[(2 * n + 1) * dim + 0] = 0.0; 2266 coords[(2 * n + 1) * dim + 1] = 0.0; 2267 coords[(2 * n + 1) * dim + 2] = 0.0; 2268 } 2269 PetscCall(VecRestoreArray(coordinates, &coords)); 2270 PetscCall(DMSetCoordinatesLocal(dm, coordinates)); 2271 PetscCall(VecDestroy(&coordinates)); 2272 } 2273 PetscCall(PetscLogEventEnd(DMPLEX_Generate, dm, 0, 0, 0)); 2274 /* Interpolate */ 2275 if (interpolate) PetscCall(DMPlexInterpolateInPlace_Internal(dm)); 2276 PetscFunctionReturn(PETSC_SUCCESS); 2277 } 2278 2279 /*@ 2280 DMPlexCreateWedgeCylinderMesh - Creates a mesh on the tensor product of the unit interval with the circle (cylinder) using wedges. 2281 2282 Collective 2283 2284 Input Parameters: 2285 + comm - The communicator for the `DM` object 2286 . n - The number of wedges around the origin 2287 - interpolate - Create edges and faces 2288 2289 Output Parameter: 2290 . dm - The `DM` object 2291 2292 Level: beginner 2293 2294 .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexCreateHexCylinderMesh()`, `DMPlexCreateBoxMesh()`, `DMSetType()`, `DMCreate()` 2295 @*/ 2296 PetscErrorCode DMPlexCreateWedgeCylinderMesh(MPI_Comm comm, PetscInt n, PetscBool interpolate, DM *dm) 2297 { 2298 PetscFunctionBegin; 2299 PetscAssertPointer(dm, 4); 2300 PetscCall(DMCreate(comm, dm)); 2301 PetscCall(DMSetType(*dm, DMPLEX)); 2302 PetscCall(DMPlexCreateWedgeCylinderMesh_Internal(*dm, n, interpolate)); 2303 PetscFunctionReturn(PETSC_SUCCESS); 2304 } 2305 2306 static inline PetscReal DiffNormReal(PetscInt dim, const PetscReal x[], const PetscReal y[]) 2307 { 2308 PetscReal prod = 0.0; 2309 PetscInt i; 2310 for (i = 0; i < dim; ++i) prod += PetscSqr(x[i] - y[i]); 2311 return PetscSqrtReal(prod); 2312 } 2313 static inline PetscReal DotReal(PetscInt dim, const PetscReal x[], const PetscReal y[]) 2314 { 2315 PetscReal prod = 0.0; 2316 PetscInt i; 2317 for (i = 0; i < dim; ++i) prod += x[i] * y[i]; 2318 return prod; 2319 } 2320 2321 /* The first constant is the sphere radius */ 2322 static void snapToSphere(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 2323 { 2324 PetscReal r = PetscRealPart(constants[0]); 2325 PetscReal norm2 = 0.0, fac; 2326 PetscInt n = uOff[1] - uOff[0], d; 2327 2328 for (d = 0; d < n; ++d) norm2 += PetscSqr(PetscRealPart(u[d])); 2329 fac = r / PetscSqrtReal(norm2); 2330 for (d = 0; d < n; ++d) f0[d] = u[d] * fac; 2331 } 2332 2333 static PetscErrorCode DMPlexCreateSphereMesh_Internal(DM dm, PetscInt dim, PetscBool simplex, PetscReal R) 2334 { 2335 const PetscInt embedDim = dim + 1; 2336 PetscSection coordSection; 2337 Vec coordinates; 2338 PetscScalar *coords; 2339 PetscReal *coordsIn; 2340 PetscInt numCells, numEdges, numVerts = 0, firstVertex = 0, v, firstEdge, coordSize, d, e; 2341 PetscMPIInt rank; 2342 2343 PetscFunctionBegin; 2344 PetscValidLogicalCollectiveBool(dm, simplex, 3); 2345 PetscCall(PetscLogEventBegin(DMPLEX_Generate, dm, 0, 0, 0)); 2346 PetscCall(DMSetDimension(dm, dim)); 2347 PetscCall(DMSetCoordinateDim(dm, dim + 1)); 2348 PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank)); 2349 switch (dim) { 2350 case 1: 2351 numCells = 16; 2352 numVerts = numCells; 2353 2354 // Build Topology 2355 PetscCall(DMPlexSetChart(dm, 0, numCells + numVerts)); 2356 for (PetscInt c = 0; c < numCells; c++) PetscCall(DMPlexSetConeSize(dm, c, embedDim)); 2357 PetscCall(DMSetUp(dm)); 2358 for (PetscInt c = 0; c < numCells; ++c) { 2359 PetscInt cone[2]; 2360 2361 cone[0] = c + numCells; 2362 cone[1] = (c + 1) % numVerts + numCells; 2363 PetscCall(DMPlexSetCone(dm, c, cone)); 2364 } 2365 PetscCall(DMPlexSymmetrize(dm)); 2366 PetscCall(DMPlexStratify(dm)); 2367 PetscCall(PetscMalloc1(numVerts * embedDim, &coordsIn)); 2368 for (PetscInt v = 0; v < numVerts; ++v) { 2369 const PetscReal rad = 2. * PETSC_PI * v / numVerts; 2370 2371 coordsIn[v * embedDim + 0] = PetscCosReal(rad); 2372 coordsIn[v * embedDim + 1] = PetscSinReal(rad); 2373 } 2374 break; 2375 case 2: 2376 if (simplex) { 2377 const PetscReal radius = PetscSqrtReal(1 + PETSC_PHI * PETSC_PHI) / (1.0 + PETSC_PHI); 2378 const PetscReal edgeLen = 2.0 / (1.0 + PETSC_PHI) * (R / radius); 2379 const PetscInt degree = 5; 2380 PetscReal vertex[3] = {0.0, 1.0 / (1.0 + PETSC_PHI), PETSC_PHI / (1.0 + PETSC_PHI)}; 2381 PetscInt s[3] = {1, 1, 1}; 2382 PetscInt cone[3]; 2383 PetscInt *graph; 2384 2385 vertex[0] *= R / radius; 2386 vertex[1] *= R / radius; 2387 vertex[2] *= R / radius; 2388 numCells = rank == 0 ? 20 : 0; 2389 numVerts = rank == 0 ? 12 : 0; 2390 firstVertex = numCells; 2391 /* Use icosahedron, which for a R-sphere has coordinates which are all cyclic permutations of 2392 2393 (0, \pm 1/\phi+1, \pm \phi/\phi+1) 2394 2395 where \phi^2 - \phi - 1 = 0, meaning \phi is the golden ratio \frac{1 + \sqrt{5}}{2}. The edge 2396 length is then given by 2/(1+\phi) = 2 * 0.38197 = 0.76393. 2397 */ 2398 /* Construct vertices */ 2399 PetscCall(PetscCalloc1(numVerts * embedDim, &coordsIn)); 2400 if (rank == 0) { 2401 for (PetscInt p = 0, i = 0; p < embedDim; ++p) { 2402 for (s[1] = -1; s[1] < 2; s[1] += 2) { 2403 for (s[2] = -1; s[2] < 2; s[2] += 2) { 2404 for (d = 0; d < embedDim; ++d) coordsIn[i * embedDim + d] = s[(d + p) % embedDim] * vertex[(d + p) % embedDim]; 2405 ++i; 2406 } 2407 } 2408 } 2409 } 2410 /* Construct graph */ 2411 PetscCall(PetscCalloc1(numVerts * numVerts, &graph)); 2412 for (PetscInt i = 0; i < numVerts; ++i) { 2413 PetscInt k = 0; 2414 for (PetscInt j = 0; j < numVerts; ++j) { 2415 if (PetscAbsReal(DiffNormReal(embedDim, &coordsIn[i * embedDim], &coordsIn[j * embedDim]) - edgeLen) < PETSC_SMALL) { 2416 graph[i * numVerts + j] = 1; 2417 ++k; 2418 } 2419 } 2420 PetscCheck(k == degree, PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Invalid icosahedron, vertex %" PetscInt_FMT " degree %" PetscInt_FMT " != %" PetscInt_FMT, i, k, degree); 2421 } 2422 /* Build Topology */ 2423 PetscCall(DMPlexSetChart(dm, 0, numCells + numVerts)); 2424 for (PetscInt c = 0; c < numCells; c++) PetscCall(DMPlexSetConeSize(dm, c, embedDim)); 2425 PetscCall(DMSetUp(dm)); /* Allocate space for cones */ 2426 /* Cells */ 2427 for (PetscInt i = 0, c = 0; i < numVerts; ++i) { 2428 for (PetscInt j = 0; j < i; ++j) { 2429 for (PetscInt k = 0; k < j; ++k) { 2430 if (graph[i * numVerts + j] && graph[j * numVerts + k] && graph[k * numVerts + i]) { 2431 cone[0] = firstVertex + i; 2432 cone[1] = firstVertex + j; 2433 cone[2] = firstVertex + k; 2434 /* Check orientation */ 2435 { 2436 const PetscInt epsilon[3][3][3] = { 2437 {{0, 0, 0}, {0, 0, 1}, {0, -1, 0}}, 2438 {{0, 0, -1}, {0, 0, 0}, {1, 0, 0} }, 2439 {{0, 1, 0}, {-1, 0, 0}, {0, 0, 0} } 2440 }; 2441 PetscReal normal[3]; 2442 PetscInt e, f; 2443 2444 for (d = 0; d < embedDim; ++d) { 2445 normal[d] = 0.0; 2446 for (e = 0; e < embedDim; ++e) { 2447 for (f = 0; f < embedDim; ++f) normal[d] += epsilon[d][e][f] * (coordsIn[j * embedDim + e] - coordsIn[i * embedDim + e]) * (coordsIn[k * embedDim + f] - coordsIn[i * embedDim + f]); 2448 } 2449 } 2450 if (DotReal(embedDim, normal, &coordsIn[i * embedDim]) < 0) { 2451 PetscInt tmp = cone[1]; 2452 cone[1] = cone[2]; 2453 cone[2] = tmp; 2454 } 2455 } 2456 PetscCall(DMPlexSetCone(dm, c++, cone)); 2457 } 2458 } 2459 } 2460 } 2461 PetscCall(DMPlexSymmetrize(dm)); 2462 PetscCall(DMPlexStratify(dm)); 2463 PetscCall(PetscFree(graph)); 2464 } else { 2465 /* 2466 12-21--13 2467 | | 2468 25 4 24 2469 | | 2470 12-25--9-16--8-24--13 2471 | | | | 2472 23 5 17 0 15 3 22 2473 | | | | 2474 10-20--6-14--7-19--11 2475 | | 2476 20 1 19 2477 | | 2478 10-18--11 2479 | | 2480 23 2 22 2481 | | 2482 12-21--13 2483 */ 2484 PetscInt cone[4], ornt[4]; 2485 2486 numCells = rank == 0 ? 6 : 0; 2487 numEdges = rank == 0 ? 12 : 0; 2488 numVerts = rank == 0 ? 8 : 0; 2489 firstVertex = numCells; 2490 firstEdge = numCells + numVerts; 2491 /* Build Topology */ 2492 PetscCall(DMPlexSetChart(dm, 0, numCells + numEdges + numVerts)); 2493 for (PetscInt c = 0; c < numCells; c++) PetscCall(DMPlexSetConeSize(dm, c, 4)); 2494 for (e = firstEdge; e < firstEdge + numEdges; ++e) PetscCall(DMPlexSetConeSize(dm, e, 2)); 2495 PetscCall(DMSetUp(dm)); /* Allocate space for cones */ 2496 if (rank == 0) { 2497 /* Cell 0 */ 2498 cone[0] = 14; 2499 cone[1] = 15; 2500 cone[2] = 16; 2501 cone[3] = 17; 2502 PetscCall(DMPlexSetCone(dm, 0, cone)); 2503 ornt[0] = 0; 2504 ornt[1] = 0; 2505 ornt[2] = 0; 2506 ornt[3] = 0; 2507 PetscCall(DMPlexSetConeOrientation(dm, 0, ornt)); 2508 /* Cell 1 */ 2509 cone[0] = 18; 2510 cone[1] = 19; 2511 cone[2] = 14; 2512 cone[3] = 20; 2513 PetscCall(DMPlexSetCone(dm, 1, cone)); 2514 ornt[0] = 0; 2515 ornt[1] = 0; 2516 ornt[2] = -1; 2517 ornt[3] = 0; 2518 PetscCall(DMPlexSetConeOrientation(dm, 1, ornt)); 2519 /* Cell 2 */ 2520 cone[0] = 21; 2521 cone[1] = 22; 2522 cone[2] = 18; 2523 cone[3] = 23; 2524 PetscCall(DMPlexSetCone(dm, 2, cone)); 2525 ornt[0] = 0; 2526 ornt[1] = 0; 2527 ornt[2] = -1; 2528 ornt[3] = 0; 2529 PetscCall(DMPlexSetConeOrientation(dm, 2, ornt)); 2530 /* Cell 3 */ 2531 cone[0] = 19; 2532 cone[1] = 22; 2533 cone[2] = 24; 2534 cone[3] = 15; 2535 PetscCall(DMPlexSetCone(dm, 3, cone)); 2536 ornt[0] = -1; 2537 ornt[1] = -1; 2538 ornt[2] = 0; 2539 ornt[3] = -1; 2540 PetscCall(DMPlexSetConeOrientation(dm, 3, ornt)); 2541 /* Cell 4 */ 2542 cone[0] = 16; 2543 cone[1] = 24; 2544 cone[2] = 21; 2545 cone[3] = 25; 2546 PetscCall(DMPlexSetCone(dm, 4, cone)); 2547 ornt[0] = -1; 2548 ornt[1] = -1; 2549 ornt[2] = -1; 2550 ornt[3] = 0; 2551 PetscCall(DMPlexSetConeOrientation(dm, 4, ornt)); 2552 /* Cell 5 */ 2553 cone[0] = 20; 2554 cone[1] = 17; 2555 cone[2] = 25; 2556 cone[3] = 23; 2557 PetscCall(DMPlexSetCone(dm, 5, cone)); 2558 ornt[0] = -1; 2559 ornt[1] = -1; 2560 ornt[2] = -1; 2561 ornt[3] = -1; 2562 PetscCall(DMPlexSetConeOrientation(dm, 5, ornt)); 2563 /* Edges */ 2564 cone[0] = 6; 2565 cone[1] = 7; 2566 PetscCall(DMPlexSetCone(dm, 14, cone)); 2567 cone[0] = 7; 2568 cone[1] = 8; 2569 PetscCall(DMPlexSetCone(dm, 15, cone)); 2570 cone[0] = 8; 2571 cone[1] = 9; 2572 PetscCall(DMPlexSetCone(dm, 16, cone)); 2573 cone[0] = 9; 2574 cone[1] = 6; 2575 PetscCall(DMPlexSetCone(dm, 17, cone)); 2576 cone[0] = 10; 2577 cone[1] = 11; 2578 PetscCall(DMPlexSetCone(dm, 18, cone)); 2579 cone[0] = 11; 2580 cone[1] = 7; 2581 PetscCall(DMPlexSetCone(dm, 19, cone)); 2582 cone[0] = 6; 2583 cone[1] = 10; 2584 PetscCall(DMPlexSetCone(dm, 20, cone)); 2585 cone[0] = 12; 2586 cone[1] = 13; 2587 PetscCall(DMPlexSetCone(dm, 21, cone)); 2588 cone[0] = 13; 2589 cone[1] = 11; 2590 PetscCall(DMPlexSetCone(dm, 22, cone)); 2591 cone[0] = 10; 2592 cone[1] = 12; 2593 PetscCall(DMPlexSetCone(dm, 23, cone)); 2594 cone[0] = 13; 2595 cone[1] = 8; 2596 PetscCall(DMPlexSetCone(dm, 24, cone)); 2597 cone[0] = 12; 2598 cone[1] = 9; 2599 PetscCall(DMPlexSetCone(dm, 25, cone)); 2600 } 2601 PetscCall(DMPlexSymmetrize(dm)); 2602 PetscCall(DMPlexStratify(dm)); 2603 /* Build coordinates */ 2604 PetscCall(PetscCalloc1(numVerts * embedDim, &coordsIn)); 2605 if (rank == 0) { 2606 coordsIn[0 * embedDim + 0] = -R; 2607 coordsIn[0 * embedDim + 1] = R; 2608 coordsIn[0 * embedDim + 2] = -R; 2609 coordsIn[1 * embedDim + 0] = R; 2610 coordsIn[1 * embedDim + 1] = R; 2611 coordsIn[1 * embedDim + 2] = -R; 2612 coordsIn[2 * embedDim + 0] = R; 2613 coordsIn[2 * embedDim + 1] = -R; 2614 coordsIn[2 * embedDim + 2] = -R; 2615 coordsIn[3 * embedDim + 0] = -R; 2616 coordsIn[3 * embedDim + 1] = -R; 2617 coordsIn[3 * embedDim + 2] = -R; 2618 coordsIn[4 * embedDim + 0] = -R; 2619 coordsIn[4 * embedDim + 1] = R; 2620 coordsIn[4 * embedDim + 2] = R; 2621 coordsIn[5 * embedDim + 0] = R; 2622 coordsIn[5 * embedDim + 1] = R; 2623 coordsIn[5 * embedDim + 2] = R; 2624 coordsIn[6 * embedDim + 0] = -R; 2625 coordsIn[6 * embedDim + 1] = -R; 2626 coordsIn[6 * embedDim + 2] = R; 2627 coordsIn[7 * embedDim + 0] = R; 2628 coordsIn[7 * embedDim + 1] = -R; 2629 coordsIn[7 * embedDim + 2] = R; 2630 } 2631 } 2632 break; 2633 case 3: 2634 if (simplex) { 2635 const PetscReal edgeLen = 1.0 / PETSC_PHI; 2636 PetscReal vertexA[4] = {0.5, 0.5, 0.5, 0.5}; 2637 PetscReal vertexB[4] = {1.0, 0.0, 0.0, 0.0}; 2638 PetscReal vertexC[4] = {0.5, 0.5 * PETSC_PHI, 0.5 / PETSC_PHI, 0.0}; 2639 const PetscInt degree = 12; 2640 PetscInt s[4] = {1, 1, 1}; 2641 PetscInt evenPerm[12][4] = { 2642 {0, 1, 2, 3}, 2643 {0, 2, 3, 1}, 2644 {0, 3, 1, 2}, 2645 {1, 0, 3, 2}, 2646 {1, 2, 0, 3}, 2647 {1, 3, 2, 0}, 2648 {2, 0, 1, 3}, 2649 {2, 1, 3, 0}, 2650 {2, 3, 0, 1}, 2651 {3, 0, 2, 1}, 2652 {3, 1, 0, 2}, 2653 {3, 2, 1, 0} 2654 }; 2655 PetscInt cone[4]; 2656 PetscInt *graph, p, i, j, k, l; 2657 2658 vertexA[0] *= R; 2659 vertexA[1] *= R; 2660 vertexA[2] *= R; 2661 vertexA[3] *= R; 2662 vertexB[0] *= R; 2663 vertexB[1] *= R; 2664 vertexB[2] *= R; 2665 vertexB[3] *= R; 2666 vertexC[0] *= R; 2667 vertexC[1] *= R; 2668 vertexC[2] *= R; 2669 vertexC[3] *= R; 2670 numCells = rank == 0 ? 600 : 0; 2671 numVerts = rank == 0 ? 120 : 0; 2672 firstVertex = numCells; 2673 /* Use the 600-cell, which for a unit sphere has coordinates which are 2674 2675 1/2 (\pm 1, \pm 1, \pm 1, \pm 1) 16 2676 (\pm 1, 0, 0, 0) all cyclic permutations 8 2677 1/2 (\pm 1, \pm phi, \pm 1/phi, 0) all even permutations 96 2678 2679 where \phi^2 - \phi - 1 = 0, meaning \phi is the golden ratio \frac{1 + \sqrt{5}}{2}. The edge 2680 length is then given by 1/\phi = 0.61803. 2681 2682 http://buzzard.pugetsound.edu/sage-practice/ch03s03.html 2683 http://mathworld.wolfram.com/600-Cell.html 2684 */ 2685 /* Construct vertices */ 2686 PetscCall(PetscCalloc1(numVerts * embedDim, &coordsIn)); 2687 i = 0; 2688 if (rank == 0) { 2689 for (s[0] = -1; s[0] < 2; s[0] += 2) { 2690 for (s[1] = -1; s[1] < 2; s[1] += 2) { 2691 for (s[2] = -1; s[2] < 2; s[2] += 2) { 2692 for (s[3] = -1; s[3] < 2; s[3] += 2) { 2693 for (d = 0; d < embedDim; ++d) coordsIn[i * embedDim + d] = s[d] * vertexA[d]; 2694 ++i; 2695 } 2696 } 2697 } 2698 } 2699 for (p = 0; p < embedDim; ++p) { 2700 s[1] = s[2] = s[3] = 1; 2701 for (s[0] = -1; s[0] < 2; s[0] += 2) { 2702 for (d = 0; d < embedDim; ++d) coordsIn[i * embedDim + d] = s[(d + p) % embedDim] * vertexB[(d + p) % embedDim]; 2703 ++i; 2704 } 2705 } 2706 for (p = 0; p < 12; ++p) { 2707 s[3] = 1; 2708 for (s[0] = -1; s[0] < 2; s[0] += 2) { 2709 for (s[1] = -1; s[1] < 2; s[1] += 2) { 2710 for (s[2] = -1; s[2] < 2; s[2] += 2) { 2711 for (d = 0; d < embedDim; ++d) coordsIn[i * embedDim + d] = s[evenPerm[p][d]] * vertexC[evenPerm[p][d]]; 2712 ++i; 2713 } 2714 } 2715 } 2716 } 2717 } 2718 PetscCheck(i == numVerts, PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Invalid 600-cell, vertices %" PetscInt_FMT " != %" PetscInt_FMT, i, numVerts); 2719 /* Construct graph */ 2720 PetscCall(PetscCalloc1(numVerts * numVerts, &graph)); 2721 for (i = 0; i < numVerts; ++i) { 2722 for (j = 0, k = 0; j < numVerts; ++j) { 2723 if (PetscAbsReal(DiffNormReal(embedDim, &coordsIn[i * embedDim], &coordsIn[j * embedDim]) - edgeLen) < PETSC_SMALL) { 2724 graph[i * numVerts + j] = 1; 2725 ++k; 2726 } 2727 } 2728 PetscCheck(k == degree, PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Invalid 600-cell, vertex %" PetscInt_FMT " degree %" PetscInt_FMT " != %" PetscInt_FMT, i, k, degree); 2729 } 2730 /* Build Topology */ 2731 PetscCall(DMPlexSetChart(dm, 0, numCells + numVerts)); 2732 for (PetscInt c = 0; c < numCells; c++) PetscCall(DMPlexSetConeSize(dm, c, embedDim)); 2733 PetscCall(DMSetUp(dm)); /* Allocate space for cones */ 2734 /* Cells */ 2735 if (rank == 0) { 2736 for (PetscInt i = 0, c = 0; i < numVerts; ++i) { 2737 for (j = 0; j < i; ++j) { 2738 for (k = 0; k < j; ++k) { 2739 for (l = 0; l < k; ++l) { 2740 if (graph[i * numVerts + j] && graph[j * numVerts + k] && graph[k * numVerts + i] && graph[l * numVerts + i] && graph[l * numVerts + j] && graph[l * numVerts + k]) { 2741 cone[0] = firstVertex + i; 2742 cone[1] = firstVertex + j; 2743 cone[2] = firstVertex + k; 2744 cone[3] = firstVertex + l; 2745 /* Check orientation: https://ef.gy/linear-algebra:normal-vectors-in-higher-dimensional-spaces */ 2746 { 2747 const PetscInt epsilon[4][4][4][4] = { 2748 {{{0, 0, 0, 0}, {0, 0, 0, 0}, {0, 0, 0, 0}, {0, 0, 0, 0}}, {{0, 0, 0, 0}, {0, 0, 0, 0}, {0, 0, 0, 1}, {0, 0, -1, 0}}, {{0, 0, 0, 0}, {0, 0, 0, -1}, {0, 0, 0, 0}, {0, 1, 0, 0}}, {{0, 0, 0, 0}, {0, 0, 1, 0}, {0, -1, 0, 0}, {0, 0, 0, 0}}}, 2749 2750 {{{0, 0, 0, 0}, {0, 0, 0, 0}, {0, 0, 0, -1}, {0, 0, 1, 0}}, {{0, 0, 0, 0}, {0, 0, 0, 0}, {0, 0, 0, 0}, {0, 0, 0, 0}}, {{0, 0, 0, 1}, {0, 0, 0, 0}, {0, 0, 0, 0}, {-1, 0, 0, 0}}, {{0, 0, -1, 0}, {0, 0, 0, 0}, {1, 0, 0, 0}, {0, 0, 0, 0}}}, 2751 2752 {{{0, 0, 0, 0}, {0, 0, 0, 1}, {0, 0, 0, 0}, {0, -1, 0, 0}}, {{0, 0, 0, -1}, {0, 0, 0, 0}, {0, 0, 0, 0}, {1, 0, 0, 0}}, {{0, 0, 0, 0}, {0, 0, 0, 0}, {0, 0, 0, 0}, {0, 0, 0, 0}}, {{0, 1, 0, 0}, {-1, 0, 0, 0}, {0, 0, 0, 0}, {0, 0, 0, 0}}}, 2753 2754 {{{0, 0, 0, 0}, {0, 0, -1, 0}, {0, 1, 0, 0}, {0, 0, 0, 0}}, {{0, 0, 1, 0}, {0, 0, 0, 0}, {-1, 0, 0, 0}, {0, 0, 0, 0}}, {{0, -1, 0, 0}, {1, 0, 0, 0}, {0, 0, 0, 0}, {0, 0, 0, 0}}, {{0, 0, 0, 0}, {0, 0, 0, 0}, {0, 0, 0, 0}, {0, 0, 0, 0}} } 2755 }; 2756 PetscReal normal[4]; 2757 PetscInt e, f, g; 2758 2759 for (d = 0; d < embedDim; ++d) { 2760 normal[d] = 0.0; 2761 for (e = 0; e < embedDim; ++e) { 2762 for (f = 0; f < embedDim; ++f) { 2763 for (g = 0; g < embedDim; ++g) { 2764 normal[d] += epsilon[d][e][f][g] * (coordsIn[j * embedDim + e] - coordsIn[i * embedDim + e]) * (coordsIn[k * embedDim + f] - coordsIn[i * embedDim + f]) * (coordsIn[l * embedDim + f] - coordsIn[i * embedDim + f]); 2765 } 2766 } 2767 } 2768 } 2769 if (DotReal(embedDim, normal, &coordsIn[i * embedDim]) < 0) { 2770 PetscInt tmp = cone[1]; 2771 cone[1] = cone[2]; 2772 cone[2] = tmp; 2773 } 2774 } 2775 PetscCall(DMPlexSetCone(dm, c++, cone)); 2776 } 2777 } 2778 } 2779 } 2780 } 2781 } 2782 PetscCall(DMPlexSymmetrize(dm)); 2783 PetscCall(DMPlexStratify(dm)); 2784 PetscCall(PetscFree(graph)); 2785 } 2786 break; 2787 default: 2788 SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported dimension for sphere: %" PetscInt_FMT, dim); 2789 } 2790 /* Create coordinates */ 2791 PetscCall(DMGetCoordinateSection(dm, &coordSection)); 2792 PetscCall(PetscSectionSetNumFields(coordSection, 1)); 2793 PetscCall(PetscSectionSetFieldComponents(coordSection, 0, embedDim)); 2794 PetscCall(PetscSectionSetChart(coordSection, firstVertex, firstVertex + numVerts)); 2795 for (v = firstVertex; v < firstVertex + numVerts; ++v) { 2796 PetscCall(PetscSectionSetDof(coordSection, v, embedDim)); 2797 PetscCall(PetscSectionSetFieldDof(coordSection, v, 0, embedDim)); 2798 } 2799 PetscCall(PetscSectionSetUp(coordSection)); 2800 PetscCall(PetscSectionGetStorageSize(coordSection, &coordSize)); 2801 PetscCall(VecCreate(PETSC_COMM_SELF, &coordinates)); 2802 PetscCall(VecSetBlockSize(coordinates, embedDim)); 2803 PetscCall(PetscObjectSetName((PetscObject)coordinates, "coordinates")); 2804 PetscCall(VecSetSizes(coordinates, coordSize, PETSC_DETERMINE)); 2805 PetscCall(VecSetType(coordinates, VECSTANDARD)); 2806 PetscCall(VecGetArray(coordinates, &coords)); 2807 for (v = 0; v < numVerts; ++v) 2808 for (d = 0; d < embedDim; ++d) coords[v * embedDim + d] = coordsIn[v * embedDim + d]; 2809 PetscCall(VecRestoreArray(coordinates, &coords)); 2810 PetscCall(DMSetCoordinatesLocal(dm, coordinates)); 2811 PetscCall(VecDestroy(&coordinates)); 2812 PetscCall(PetscFree(coordsIn)); 2813 { 2814 DM cdm; 2815 PetscDS cds; 2816 PetscScalar c = R; 2817 2818 PetscCall(DMPlexCreateCoordinateSpace(dm, 1, PETSC_TRUE, snapToSphere)); 2819 PetscCall(DMGetCoordinateDM(dm, &cdm)); 2820 PetscCall(DMGetDS(cdm, &cds)); 2821 PetscCall(PetscDSSetConstants(cds, 1, &c)); 2822 } 2823 PetscCall(PetscLogEventEnd(DMPLEX_Generate, dm, 0, 0, 0)); 2824 /* Wait for coordinate creation before doing in-place modification */ 2825 if (simplex) PetscCall(DMPlexInterpolateInPlace_Internal(dm)); 2826 PetscFunctionReturn(PETSC_SUCCESS); 2827 } 2828 2829 typedef void (*TPSEvaluateFunc)(const PetscReal[], PetscReal *, PetscReal[], PetscReal (*)[3]); 2830 2831 /* 2832 The Schwarz P implicit surface is 2833 2834 f(x) = cos(x0) + cos(x1) + cos(x2) = 0 2835 */ 2836 static void TPSEvaluate_SchwarzP(const PetscReal y[3], PetscReal *f, PetscReal grad[], PetscReal (*hess)[3]) 2837 { 2838 PetscReal c[3] = {PetscCosReal(y[0] * PETSC_PI), PetscCosReal(y[1] * PETSC_PI), PetscCosReal(y[2] * PETSC_PI)}; 2839 PetscReal g[3] = {-PetscSinReal(y[0] * PETSC_PI), -PetscSinReal(y[1] * PETSC_PI), -PetscSinReal(y[2] * PETSC_PI)}; 2840 f[0] = c[0] + c[1] + c[2]; 2841 for (PetscInt i = 0; i < 3; i++) { 2842 grad[i] = PETSC_PI * g[i]; 2843 for (PetscInt j = 0; j < 3; j++) hess[i][j] = (i == j) ? -PetscSqr(PETSC_PI) * c[i] : 0.; 2844 } 2845 } 2846 2847 // u[] is a tentative normal on input. Replace with the implicit function gradient in the same direction 2848 static PetscErrorCode TPSExtrudeNormalFunc_SchwarzP(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt r, PetscScalar u[], void *ctx) 2849 { 2850 for (PetscInt i = 0; i < 3; i++) u[i] = -PETSC_PI * PetscSinReal(x[i] * PETSC_PI); 2851 return PETSC_SUCCESS; 2852 } 2853 2854 /* 2855 The Gyroid implicit surface is 2856 2857 f(x,y,z) = sin(pi * x) * cos (pi * (y + 1/2)) + sin(pi * (y + 1/2)) * cos(pi * (z + 1/4)) + sin(pi * (z + 1/4)) * cos(pi * x) 2858 2859 */ 2860 static void TPSEvaluate_Gyroid(const PetscReal y[3], PetscReal *f, PetscReal grad[], PetscReal (*hess)[3]) 2861 { 2862 PetscReal s[3] = {PetscSinReal(PETSC_PI * y[0]), PetscSinReal(PETSC_PI * (y[1] + .5)), PetscSinReal(PETSC_PI * (y[2] + .25))}; 2863 PetscReal c[3] = {PetscCosReal(PETSC_PI * y[0]), PetscCosReal(PETSC_PI * (y[1] + .5)), PetscCosReal(PETSC_PI * (y[2] + .25))}; 2864 f[0] = s[0] * c[1] + s[1] * c[2] + s[2] * c[0]; 2865 grad[0] = PETSC_PI * (c[0] * c[1] - s[2] * s[0]); 2866 grad[1] = PETSC_PI * (c[1] * c[2] - s[0] * s[1]); 2867 grad[2] = PETSC_PI * (c[2] * c[0] - s[1] * s[2]); 2868 hess[0][0] = -PetscSqr(PETSC_PI) * (s[0] * c[1] + s[2] * c[0]); 2869 hess[0][1] = -PetscSqr(PETSC_PI) * (c[0] * s[1]); 2870 hess[0][2] = -PetscSqr(PETSC_PI) * (c[2] * s[0]); 2871 hess[1][0] = -PetscSqr(PETSC_PI) * (s[1] * c[2] + s[0] * c[1]); 2872 hess[1][1] = -PetscSqr(PETSC_PI) * (c[1] * s[2]); 2873 hess[2][2] = -PetscSqr(PETSC_PI) * (c[0] * s[1]); 2874 hess[2][0] = -PetscSqr(PETSC_PI) * (s[2] * c[0] + s[1] * c[2]); 2875 hess[2][1] = -PetscSqr(PETSC_PI) * (c[2] * s[0]); 2876 hess[2][2] = -PetscSqr(PETSC_PI) * (c[1] * s[2]); 2877 } 2878 2879 // u[] is a tentative normal on input. Replace with the implicit function gradient in the same direction 2880 static PetscErrorCode TPSExtrudeNormalFunc_Gyroid(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt r, PetscScalar u[], void *ctx) 2881 { 2882 PetscReal s[3] = {PetscSinReal(PETSC_PI * x[0]), PetscSinReal(PETSC_PI * (x[1] + .5)), PetscSinReal(PETSC_PI * (x[2] + .25))}; 2883 PetscReal c[3] = {PetscCosReal(PETSC_PI * x[0]), PetscCosReal(PETSC_PI * (x[1] + .5)), PetscCosReal(PETSC_PI * (x[2] + .25))}; 2884 u[0] = PETSC_PI * (c[0] * c[1] - s[2] * s[0]); 2885 u[1] = PETSC_PI * (c[1] * c[2] - s[0] * s[1]); 2886 u[2] = PETSC_PI * (c[2] * c[0] - s[1] * s[2]); 2887 return PETSC_SUCCESS; 2888 } 2889 2890 /* 2891 We wish to solve 2892 2893 min_y || y - x ||^2 subject to f(y) = 0 2894 2895 Let g(y) = grad(f). The minimization problem is equivalent to asking to satisfy 2896 f(y) = 0 and (y-x) is parallel to g(y). We do this by using Householder QR to obtain a basis for the 2897 tangent space and ask for both components in the tangent space to be zero. 2898 2899 Take g to be a column vector and compute the "full QR" factorization Q R = g, 2900 where Q = I - 2 n n^T is a symmetric orthogonal matrix. 2901 The first column of Q is parallel to g so the remaining two columns span the null space. 2902 Let Qn = Q[:,1:] be those remaining columns. Then Qn Qn^T is an orthogonal projector into the tangent space. 2903 Since Q is symmetric, this is equivalent to multiplying by Q and taking the last two entries. 2904 In total, we have a system of 3 equations in 3 unknowns: 2905 2906 f(y) = 0 1 equation 2907 Qn^T (y - x) = 0 2 equations 2908 2909 Here, we compute the residual and Jacobian of this system. 2910 */ 2911 static void TPSNearestPointResJac(TPSEvaluateFunc feval, const PetscScalar x[], const PetscScalar y[], PetscScalar res[], PetscScalar J[]) 2912 { 2913 PetscReal yreal[3] = {PetscRealPart(y[0]), PetscRealPart(y[1]), PetscRealPart(y[2])}; 2914 PetscReal d[3] = {PetscRealPart(y[0] - x[0]), PetscRealPart(y[1] - x[1]), PetscRealPart(y[2] - x[2])}; 2915 PetscReal f, grad[3], n[3], norm, norm_y[3], nd, nd_y[3], sign; 2916 PetscReal n_y[3][3] = { 2917 {0, 0, 0}, 2918 {0, 0, 0}, 2919 {0, 0, 0} 2920 }; 2921 2922 feval(yreal, &f, grad, n_y); 2923 2924 for (PetscInt i = 0; i < 3; i++) n[i] = grad[i]; 2925 norm = PetscSqrtReal(PetscSqr(n[0]) + PetscSqr(n[1]) + PetscSqr(n[2])); 2926 for (PetscInt i = 0; i < 3; i++) norm_y[i] = 1. / norm * n[i] * n_y[i][i]; 2927 2928 // Define the Householder reflector 2929 sign = n[0] >= 0 ? 1. : -1.; 2930 n[0] += norm * sign; 2931 for (PetscInt i = 0; i < 3; i++) n_y[0][i] += norm_y[i] * sign; 2932 2933 norm = PetscSqrtReal(PetscSqr(n[0]) + PetscSqr(n[1]) + PetscSqr(n[2])); 2934 norm_y[0] = 1. / norm * (n[0] * n_y[0][0]); 2935 norm_y[1] = 1. / norm * (n[0] * n_y[0][1] + n[1] * n_y[1][1]); 2936 norm_y[2] = 1. / norm * (n[0] * n_y[0][2] + n[2] * n_y[2][2]); 2937 2938 for (PetscInt i = 0; i < 3; i++) { 2939 n[i] /= norm; 2940 for (PetscInt j = 0; j < 3; j++) { 2941 // note that n[i] is n_old[i]/norm when executing the code below 2942 n_y[i][j] = n_y[i][j] / norm - n[i] / norm * norm_y[j]; 2943 } 2944 } 2945 2946 nd = n[0] * d[0] + n[1] * d[1] + n[2] * d[2]; 2947 for (PetscInt i = 0; i < 3; i++) nd_y[i] = n[i] + n_y[0][i] * d[0] + n_y[1][i] * d[1] + n_y[2][i] * d[2]; 2948 2949 res[0] = f; 2950 res[1] = d[1] - 2 * n[1] * nd; 2951 res[2] = d[2] - 2 * n[2] * nd; 2952 // J[j][i] is J_{ij} (column major) 2953 for (PetscInt j = 0; j < 3; j++) { 2954 J[0 + j * 3] = grad[j]; 2955 J[1 + j * 3] = (j == 1) * 1. - 2 * (n_y[1][j] * nd + n[1] * nd_y[j]); 2956 J[2 + j * 3] = (j == 2) * 1. - 2 * (n_y[2][j] * nd + n[2] * nd_y[j]); 2957 } 2958 } 2959 2960 /* 2961 Project x to the nearest point on the implicit surface using Newton's method. 2962 */ 2963 static PetscErrorCode TPSNearestPoint(TPSEvaluateFunc feval, PetscScalar x[]) 2964 { 2965 PetscScalar y[3] = {x[0], x[1], x[2]}; // Initial guess 2966 2967 PetscFunctionBegin; 2968 for (PetscInt iter = 0; iter < 10; iter++) { 2969 PetscScalar res[3], J[9]; 2970 PetscReal resnorm; 2971 TPSNearestPointResJac(feval, x, y, res, J); 2972 resnorm = PetscSqrtReal(PetscSqr(PetscRealPart(res[0])) + PetscSqr(PetscRealPart(res[1])) + PetscSqr(PetscRealPart(res[2]))); 2973 if (0) { // Turn on this monitor if you need to confirm quadratic convergence 2974 PetscCall(PetscPrintf(PETSC_COMM_SELF, "[%" PetscInt_FMT "] res [%g %g %g]\n", iter, (double)PetscRealPart(res[0]), (double)PetscRealPart(res[1]), (double)PetscRealPart(res[2]))); 2975 } 2976 if (resnorm < PETSC_SMALL) break; 2977 2978 // Take the Newton step 2979 PetscCall(PetscKernel_A_gets_inverse_A_3(J, 0., PETSC_FALSE, NULL)); 2980 PetscKernel_v_gets_v_minus_A_times_w_3(y, J, res); 2981 } 2982 for (PetscInt i = 0; i < 3; i++) x[i] = y[i]; 2983 PetscFunctionReturn(PETSC_SUCCESS); 2984 } 2985 2986 const char *const DMPlexTPSTypes[] = {"SCHWARZ_P", "GYROID", "DMPlexTPSType", "DMPLEX_TPS_", NULL}; 2987 2988 static PetscErrorCode DMPlexCreateTPSMesh_Internal(DM dm, DMPlexTPSType tpstype, const PetscInt extent[], const DMBoundaryType periodic[], PetscBool tps_distribute, PetscInt refinements, PetscInt layers, PetscReal thickness) 2989 { 2990 PetscMPIInt rank; 2991 PetscInt topoDim = 2, spaceDim = 3, numFaces = 0, numVertices = 0, numEdges = 0; 2992 PetscInt(*edges)[2] = NULL, *edgeSets = NULL; 2993 PetscInt *cells_flat = NULL; 2994 PetscReal *vtxCoords = NULL; 2995 TPSEvaluateFunc evalFunc = NULL; 2996 PetscSimplePointFunc normalFunc = NULL; 2997 DMLabel label; 2998 2999 PetscFunctionBegin; 3000 PetscCall(PetscLogEventBegin(DMPLEX_Generate, dm, 0, 0, 0)); 3001 PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank)); 3002 PetscCheck((layers != 0) ^ (thickness == 0.), PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_INCOMP, "Layers %" PetscInt_FMT " must be nonzero iff thickness %g is nonzero", layers, (double)thickness); 3003 switch (tpstype) { 3004 case DMPLEX_TPS_SCHWARZ_P: 3005 PetscCheck(!periodic || (periodic[0] == DM_BOUNDARY_NONE && periodic[1] == DM_BOUNDARY_NONE && periodic[2] == DM_BOUNDARY_NONE), PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Schwarz P does not support periodic meshes"); 3006 if (rank == 0) { 3007 PetscInt(*cells)[6][4][4] = NULL; // [junction, junction-face, cell, conn] 3008 PetscInt Njunctions = 0, Ncuts = 0, Npipes[3], vcount; 3009 PetscReal L = 1; 3010 3011 Npipes[0] = (extent[0] + 1) * extent[1] * extent[2]; 3012 Npipes[1] = extent[0] * (extent[1] + 1) * extent[2]; 3013 Npipes[2] = extent[0] * extent[1] * (extent[2] + 1); 3014 Njunctions = extent[0] * extent[1] * extent[2]; 3015 Ncuts = 2 * (extent[0] * extent[1] + extent[1] * extent[2] + extent[2] * extent[0]); 3016 numVertices = 4 * (Npipes[0] + Npipes[1] + Npipes[2]) + 8 * Njunctions; 3017 PetscCall(PetscMalloc1(3 * numVertices, &vtxCoords)); 3018 PetscCall(PetscMalloc1(Njunctions, &cells)); 3019 PetscCall(PetscMalloc1(Ncuts * 4, &edges)); 3020 PetscCall(PetscMalloc1(Ncuts * 4, &edgeSets)); 3021 // x-normal pipes 3022 vcount = 0; 3023 for (PetscInt i = 0; i < extent[0] + 1; i++) { 3024 for (PetscInt j = 0; j < extent[1]; j++) { 3025 for (PetscInt k = 0; k < extent[2]; k++) { 3026 for (PetscInt l = 0; l < 4; l++) { 3027 vtxCoords[vcount++] = (2 * i - 1) * L; 3028 vtxCoords[vcount++] = 2 * j * L + PetscCosReal((2 * l + 1) * PETSC_PI / 4) * L / 2; 3029 vtxCoords[vcount++] = 2 * k * L + PetscSinReal((2 * l + 1) * PETSC_PI / 4) * L / 2; 3030 } 3031 } 3032 } 3033 } 3034 // y-normal pipes 3035 for (PetscInt i = 0; i < extent[0]; i++) { 3036 for (PetscInt j = 0; j < extent[1] + 1; j++) { 3037 for (PetscInt k = 0; k < extent[2]; k++) { 3038 for (PetscInt l = 0; l < 4; l++) { 3039 vtxCoords[vcount++] = 2 * i * L + PetscSinReal((2 * l + 1) * PETSC_PI / 4) * L / 2; 3040 vtxCoords[vcount++] = (2 * j - 1) * L; 3041 vtxCoords[vcount++] = 2 * k * L + PetscCosReal((2 * l + 1) * PETSC_PI / 4) * L / 2; 3042 } 3043 } 3044 } 3045 } 3046 // z-normal pipes 3047 for (PetscInt i = 0; i < extent[0]; i++) { 3048 for (PetscInt j = 0; j < extent[1]; j++) { 3049 for (PetscInt k = 0; k < extent[2] + 1; k++) { 3050 for (PetscInt l = 0; l < 4; l++) { 3051 vtxCoords[vcount++] = 2 * i * L + PetscCosReal((2 * l + 1) * PETSC_PI / 4) * L / 2; 3052 vtxCoords[vcount++] = 2 * j * L + PetscSinReal((2 * l + 1) * PETSC_PI / 4) * L / 2; 3053 vtxCoords[vcount++] = (2 * k - 1) * L; 3054 } 3055 } 3056 } 3057 } 3058 // junctions 3059 for (PetscInt i = 0; i < extent[0]; i++) { 3060 for (PetscInt j = 0; j < extent[1]; j++) { 3061 for (PetscInt k = 0; k < extent[2]; k++) { 3062 const PetscInt J = (i * extent[1] + j) * extent[2] + k, Jvoff = (Npipes[0] + Npipes[1] + Npipes[2]) * 4 + J * 8; 3063 PetscCheck(vcount / 3 == Jvoff, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Unexpected vertex count"); 3064 for (PetscInt ii = 0; ii < 2; ii++) { 3065 for (PetscInt jj = 0; jj < 2; jj++) { 3066 for (PetscInt kk = 0; kk < 2; kk++) { 3067 double Ls = (1 - sqrt(2) / 4) * L; 3068 vtxCoords[vcount++] = 2 * i * L + (2 * ii - 1) * Ls; 3069 vtxCoords[vcount++] = 2 * j * L + (2 * jj - 1) * Ls; 3070 vtxCoords[vcount++] = 2 * k * L + (2 * kk - 1) * Ls; 3071 } 3072 } 3073 } 3074 const PetscInt jfaces[3][2][4] = { 3075 {{3, 1, 0, 2}, {7, 5, 4, 6}}, // x-aligned 3076 {{5, 4, 0, 1}, {7, 6, 2, 3}}, // y-aligned 3077 {{6, 2, 0, 4}, {7, 3, 1, 5}} // z-aligned 3078 }; 3079 const PetscInt pipe_lo[3] = {// vertex numbers of pipes 3080 ((i * extent[1] + j) * extent[2] + k) * 4, ((i * (extent[1] + 1) + j) * extent[2] + k + Npipes[0]) * 4, ((i * extent[1] + j) * (extent[2] + 1) + k + Npipes[0] + Npipes[1]) * 4}; 3081 const PetscInt pipe_hi[3] = {// vertex numbers of pipes 3082 (((i + 1) * extent[1] + j) * extent[2] + k) * 4, ((i * (extent[1] + 1) + j + 1) * extent[2] + k + Npipes[0]) * 4, ((i * extent[1] + j) * (extent[2] + 1) + k + 1 + Npipes[0] + Npipes[1]) * 4}; 3083 for (PetscInt dir = 0; dir < 3; dir++) { // x,y,z 3084 const PetscInt ijk[3] = {i, j, k}; 3085 for (PetscInt l = 0; l < 4; l++) { // rotations 3086 cells[J][dir * 2 + 0][l][0] = pipe_lo[dir] + l; 3087 cells[J][dir * 2 + 0][l][1] = Jvoff + jfaces[dir][0][l]; 3088 cells[J][dir * 2 + 0][l][2] = Jvoff + jfaces[dir][0][(l - 1 + 4) % 4]; 3089 cells[J][dir * 2 + 0][l][3] = pipe_lo[dir] + (l - 1 + 4) % 4; 3090 cells[J][dir * 2 + 1][l][0] = Jvoff + jfaces[dir][1][l]; 3091 cells[J][dir * 2 + 1][l][1] = pipe_hi[dir] + l; 3092 cells[J][dir * 2 + 1][l][2] = pipe_hi[dir] + (l - 1 + 4) % 4; 3093 cells[J][dir * 2 + 1][l][3] = Jvoff + jfaces[dir][1][(l - 1 + 4) % 4]; 3094 if (ijk[dir] == 0) { 3095 edges[numEdges][0] = pipe_lo[dir] + l; 3096 edges[numEdges][1] = pipe_lo[dir] + (l + 1) % 4; 3097 edgeSets[numEdges] = dir * 2 + 1; 3098 numEdges++; 3099 } 3100 if (ijk[dir] + 1 == extent[dir]) { 3101 edges[numEdges][0] = pipe_hi[dir] + l; 3102 edges[numEdges][1] = pipe_hi[dir] + (l + 1) % 4; 3103 edgeSets[numEdges] = dir * 2 + 2; 3104 numEdges++; 3105 } 3106 } 3107 } 3108 } 3109 } 3110 } 3111 PetscCheck(numEdges == Ncuts * 4, PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Edge count %" PetscInt_FMT " incompatible with number of cuts %" PetscInt_FMT, numEdges, Ncuts); 3112 numFaces = 24 * Njunctions; 3113 cells_flat = cells[0][0][0]; 3114 } 3115 evalFunc = TPSEvaluate_SchwarzP; 3116 normalFunc = TPSExtrudeNormalFunc_SchwarzP; 3117 break; 3118 case DMPLEX_TPS_GYROID: 3119 if (rank == 0) { 3120 // This is a coarse mesh approximation of the gyroid shifted to being the zero of the level set 3121 // 3122 // sin(pi*x)*cos(pi*(y+1/2)) + sin(pi*(y+1/2))*cos(pi*(z+1/4)) + sin(pi*(z+1/4))*cos(x) 3123 // 3124 // on the cell [0,2]^3. 3125 // 3126 // Think about dividing that cell into four columns, and focus on the column [0,1]x[0,1]x[0,2]. 3127 // If you looked at the gyroid in that column at different slices of z you would see that it kind of spins 3128 // like a boomerang: 3129 // 3130 // z = 0 z = 1/4 z = 1/2 z = 3/4 // 3131 // ----- ------- ------- ------- // 3132 // // 3133 // + + + + + + + \ + // 3134 // \ / \ // 3135 // \ `-_ _-' / } // 3136 // *-_ `-' _-' / // 3137 // + `-+ + + +-' + + / + // 3138 // // 3139 // // 3140 // z = 1 z = 5/4 z = 3/2 z = 7/4 // 3141 // ----- ------- ------- ------- // 3142 // // 3143 // +-_ + + + + _-+ + / + // 3144 // `-_ _-_ _-` / // 3145 // \ _-' `-_ / { // 3146 // \ / \ // 3147 // + + + + + + + \ + // 3148 // 3149 // 3150 // This course mesh approximates each of these slices by two line segments, 3151 // and then connects the segments in consecutive layers with quadrilateral faces. 3152 // All of the end points of the segments are multiples of 1/4 except for the 3153 // point * in the picture for z = 0 above and the similar points in other layers. 3154 // That point is at (gamma, gamma, 0), where gamma is calculated below. 3155 // 3156 // The column [1,2]x[1,2]x[0,2] looks the same as this column; 3157 // The columns [1,2]x[0,1]x[0,2] and [0,1]x[1,2]x[0,2] are mirror images. 3158 // 3159 // As for how this method turned into the names given to the vertices: 3160 // that was not systematic, it was just the way it worked out in my handwritten notes. 3161 3162 PetscInt facesPerBlock = 64; 3163 PetscInt vertsPerBlock = 56; 3164 PetscInt extentPlus[3]; 3165 PetscInt numBlocks, numBlocksPlus; 3166 const PetscInt A = 0, B = 1, C = 2, D = 3, E = 4, F = 5, G = 6, H = 7, II = 8, J = 9, K = 10, L = 11, M = 12, N = 13, O = 14, P = 15, Q = 16, R = 17, S = 18, T = 19, U = 20, V = 21, W = 22, X = 23, Y = 24, Z = 25, Ap = 26, Bp = 27, Cp = 28, Dp = 29, Ep = 30, Fp = 31, Gp = 32, Hp = 33, Ip = 34, Jp = 35, Kp = 36, Lp = 37, Mp = 38, Np = 39, Op = 40, Pp = 41, Qp = 42, Rp = 43, Sp = 44, Tp = 45, Up = 46, Vp = 47, Wp = 48, Xp = 49, Yp = 50, Zp = 51, Aq = 52, Bq = 53, Cq = 54, Dq = 55; 3167 const PetscInt pattern[64][4] = { 3168 /* face to vertex within the coarse discretization of a single gyroid block */ 3169 /* layer 0 */ 3170 {A, C, K, G }, 3171 {C, B, II, K }, 3172 {D, A, H, L }, 3173 {B + 56 * 1, D, L, J }, 3174 {E, B + 56 * 1, J, N }, 3175 {A + 56 * 2, E, N, H + 56 * 2 }, 3176 {F, A + 56 * 2, G + 56 * 2, M }, 3177 {B, F, M, II }, 3178 /* layer 1 */ 3179 {G, K, Q, O }, 3180 {K, II, P, Q }, 3181 {L, H, O + 56 * 1, R }, 3182 {J, L, R, P }, 3183 {N, J, P, S }, 3184 {H + 56 * 2, N, S, O + 56 * 3 }, 3185 {M, G + 56 * 2, O + 56 * 2, T }, 3186 {II, M, T, P }, 3187 /* layer 2 */ 3188 {O, Q, Y, U }, 3189 {Q, P, W, Y }, 3190 {R, O + 56 * 1, U + 56 * 1, Ap }, 3191 {P, R, Ap, W }, 3192 {S, P, X, Bp }, 3193 {O + 56 * 3, S, Bp, V + 56 * 1 }, 3194 {T, O + 56 * 2, V, Z }, 3195 {P, T, Z, X }, 3196 /* layer 3 */ 3197 {U, Y, Ep, Dp }, 3198 {Y, W, Cp, Ep }, 3199 {Ap, U + 56 * 1, Dp + 56 * 1, Gp }, 3200 {W, Ap, Gp, Cp }, 3201 {Bp, X, Cp + 56 * 2, Fp }, 3202 {V + 56 * 1, Bp, Fp, Dp + 56 * 1}, 3203 {Z, V, Dp, Hp }, 3204 {X, Z, Hp, Cp + 56 * 2}, 3205 /* layer 4 */ 3206 {Dp, Ep, Mp, Kp }, 3207 {Ep, Cp, Ip, Mp }, 3208 {Gp, Dp + 56 * 1, Lp, Np }, 3209 {Cp, Gp, Np, Jp }, 3210 {Fp, Cp + 56 * 2, Jp + 56 * 2, Pp }, 3211 {Dp + 56 * 1, Fp, Pp, Lp }, 3212 {Hp, Dp, Kp, Op }, 3213 {Cp + 56 * 2, Hp, Op, Ip + 56 * 2}, 3214 /* layer 5 */ 3215 {Kp, Mp, Sp, Rp }, 3216 {Mp, Ip, Qp, Sp }, 3217 {Np, Lp, Rp, Tp }, 3218 {Jp, Np, Tp, Qp + 56 * 1}, 3219 {Pp, Jp + 56 * 2, Qp + 56 * 3, Up }, 3220 {Lp, Pp, Up, Rp }, 3221 {Op, Kp, Rp, Vp }, 3222 {Ip + 56 * 2, Op, Vp, Qp + 56 * 2}, 3223 /* layer 6 */ 3224 {Rp, Sp, Aq, Yp }, 3225 {Sp, Qp, Wp, Aq }, 3226 {Tp, Rp, Yp, Cq }, 3227 {Qp + 56 * 1, Tp, Cq, Wp + 56 * 1}, 3228 {Up, Qp + 56 * 3, Xp + 56 * 1, Dq }, 3229 {Rp, Up, Dq, Zp }, 3230 {Vp, Rp, Zp, Bq }, 3231 {Qp + 56 * 2, Vp, Bq, Xp }, 3232 /* layer 7 (the top is the periodic image of the bottom of layer 0) */ 3233 {Yp, Aq, C + 56 * 4, A + 56 * 4 }, 3234 {Aq, Wp, B + 56 * 4, C + 56 * 4 }, 3235 {Cq, Yp, A + 56 * 4, D + 56 * 4 }, 3236 {Wp + 56 * 1, Cq, D + 56 * 4, B + 56 * 5 }, 3237 {Dq, Xp + 56 * 1, B + 56 * 5, E + 56 * 4 }, 3238 {Zp, Dq, E + 56 * 4, A + 56 * 6 }, 3239 {Bq, Zp, A + 56 * 6, F + 56 * 4 }, 3240 {Xp, Bq, F + 56 * 4, B + 56 * 4 } 3241 }; 3242 const PetscReal gamma = PetscAcosReal((PetscSqrtReal(3.) - 1.) / PetscSqrtReal(2.)) / PETSC_PI; 3243 const PetscReal patternCoords[56][3] = { 3244 {1., 0., 0. }, /* A */ 3245 {0., 1., 0. }, /* B */ 3246 {gamma, gamma, 0. }, /* C */ 3247 {1 + gamma, 1 - gamma, 0. }, /* D */ 3248 {2 - gamma, 2 - gamma, 0. }, /* E */ 3249 {1 - gamma, 1 + gamma, 0. }, /* F */ 3250 3251 {.5, 0, .25 }, /* G */ 3252 {1.5, 0., .25 }, /* H */ 3253 {.5, 1., .25 }, /* II */ 3254 {1.5, 1., .25 }, /* J */ 3255 {.25, .5, .25 }, /* K */ 3256 {1.25, .5, .25 }, /* L */ 3257 {.75, 1.5, .25 }, /* M */ 3258 {1.75, 1.5, .25 }, /* N */ 3259 3260 {0., 0., .5 }, /* O */ 3261 {1., 1., .5 }, /* P */ 3262 {gamma, 1 - gamma, .5 }, /* Q */ 3263 {1 + gamma, gamma, .5 }, /* R */ 3264 {2 - gamma, 1 + gamma, .5 }, /* S */ 3265 {1 - gamma, 2 - gamma, .5 }, /* T */ 3266 3267 {0., .5, .75 }, /* U */ 3268 {0., 1.5, .75 }, /* V */ 3269 {1., .5, .75 }, /* W */ 3270 {1., 1.5, .75 }, /* X */ 3271 {.5, .75, .75 }, /* Y */ 3272 {.5, 1.75, .75 }, /* Z */ 3273 {1.5, .25, .75 }, /* Ap */ 3274 {1.5, 1.25, .75 }, /* Bp */ 3275 3276 {1., 0., 1. }, /* Cp */ 3277 {0., 1., 1. }, /* Dp */ 3278 {1 - gamma, 1 - gamma, 1. }, /* Ep */ 3279 {1 + gamma, 1 + gamma, 1. }, /* Fp */ 3280 {2 - gamma, gamma, 1. }, /* Gp */ 3281 {gamma, 2 - gamma, 1. }, /* Hp */ 3282 3283 {.5, 0., 1.25}, /* Ip */ 3284 {1.5, 0., 1.25}, /* Jp */ 3285 {.5, 1., 1.25}, /* Kp */ 3286 {1.5, 1., 1.25}, /* Lp */ 3287 {.75, .5, 1.25}, /* Mp */ 3288 {1.75, .5, 1.25}, /* Np */ 3289 {.25, 1.5, 1.25}, /* Op */ 3290 {1.25, 1.5, 1.25}, /* Pp */ 3291 3292 {0., 0., 1.5 }, /* Qp */ 3293 {1., 1., 1.5 }, /* Rp */ 3294 {1 - gamma, gamma, 1.5 }, /* Sp */ 3295 {2 - gamma, 1 - gamma, 1.5 }, /* Tp */ 3296 {1 + gamma, 2 - gamma, 1.5 }, /* Up */ 3297 {gamma, 1 + gamma, 1.5 }, /* Vp */ 3298 3299 {0., .5, 1.75}, /* Wp */ 3300 {0., 1.5, 1.75}, /* Xp */ 3301 {1., .5, 1.75}, /* Yp */ 3302 {1., 1.5, 1.75}, /* Zp */ 3303 {.5, .25, 1.75}, /* Aq */ 3304 {.5, 1.25, 1.75}, /* Bq */ 3305 {1.5, .75, 1.75}, /* Cq */ 3306 {1.5, 1.75, 1.75}, /* Dq */ 3307 }; 3308 PetscInt(*cells)[64][4] = NULL; 3309 PetscBool *seen; 3310 PetscInt *vertToTrueVert; 3311 PetscInt count; 3312 3313 for (PetscInt i = 0; i < 3; i++) extentPlus[i] = extent[i] + 1; 3314 numBlocks = 1; 3315 for (PetscInt i = 0; i < 3; i++) numBlocks *= extent[i]; 3316 numBlocksPlus = 1; 3317 for (PetscInt i = 0; i < 3; i++) numBlocksPlus *= extentPlus[i]; 3318 numFaces = numBlocks * facesPerBlock; 3319 PetscCall(PetscMalloc1(numBlocks, &cells)); 3320 PetscCall(PetscCalloc1(numBlocksPlus * vertsPerBlock, &seen)); 3321 for (PetscInt k = 0; k < extent[2]; k++) { 3322 for (PetscInt j = 0; j < extent[1]; j++) { 3323 for (PetscInt i = 0; i < extent[0]; i++) { 3324 for (PetscInt f = 0; f < facesPerBlock; f++) { 3325 for (PetscInt v = 0; v < 4; v++) { 3326 PetscInt vertRaw = pattern[f][v]; 3327 PetscInt blockidx = vertRaw / 56; 3328 PetscInt patternvert = vertRaw % 56; 3329 PetscInt xplus = (blockidx & 1); 3330 PetscInt yplus = (blockidx & 2) >> 1; 3331 PetscInt zplus = (blockidx & 4) >> 2; 3332 PetscInt zcoord = (periodic && periodic[2] == DM_BOUNDARY_PERIODIC) ? ((k + zplus) % extent[2]) : (k + zplus); 3333 PetscInt ycoord = (periodic && periodic[1] == DM_BOUNDARY_PERIODIC) ? ((j + yplus) % extent[1]) : (j + yplus); 3334 PetscInt xcoord = (periodic && periodic[0] == DM_BOUNDARY_PERIODIC) ? ((i + xplus) % extent[0]) : (i + xplus); 3335 PetscInt vert = ((zcoord * extentPlus[1] + ycoord) * extentPlus[0] + xcoord) * 56 + patternvert; 3336 3337 cells[(k * extent[1] + j) * extent[0] + i][f][v] = vert; 3338 seen[vert] = PETSC_TRUE; 3339 } 3340 } 3341 } 3342 } 3343 } 3344 for (PetscInt i = 0; i < numBlocksPlus * vertsPerBlock; i++) 3345 if (seen[i]) numVertices++; 3346 count = 0; 3347 PetscCall(PetscMalloc1(numBlocksPlus * vertsPerBlock, &vertToTrueVert)); 3348 PetscCall(PetscMalloc1(numVertices * 3, &vtxCoords)); 3349 for (PetscInt i = 0; i < numBlocksPlus * vertsPerBlock; i++) vertToTrueVert[i] = -1; 3350 for (PetscInt k = 0; k < extentPlus[2]; k++) { 3351 for (PetscInt j = 0; j < extentPlus[1]; j++) { 3352 for (PetscInt i = 0; i < extentPlus[0]; i++) { 3353 for (PetscInt v = 0; v < vertsPerBlock; v++) { 3354 PetscInt vIdx = ((k * extentPlus[1] + j) * extentPlus[0] + i) * vertsPerBlock + v; 3355 3356 if (seen[vIdx]) { 3357 PetscInt thisVert; 3358 3359 vertToTrueVert[vIdx] = thisVert = count++; 3360 3361 for (PetscInt d = 0; d < 3; d++) vtxCoords[3 * thisVert + d] = patternCoords[v][d]; 3362 vtxCoords[3 * thisVert + 0] += i * 2; 3363 vtxCoords[3 * thisVert + 1] += j * 2; 3364 vtxCoords[3 * thisVert + 2] += k * 2; 3365 } 3366 } 3367 } 3368 } 3369 } 3370 for (PetscInt i = 0; i < numBlocks; i++) { 3371 for (PetscInt f = 0; f < facesPerBlock; f++) { 3372 for (PetscInt v = 0; v < 4; v++) cells[i][f][v] = vertToTrueVert[cells[i][f][v]]; 3373 } 3374 } 3375 PetscCall(PetscFree(vertToTrueVert)); 3376 PetscCall(PetscFree(seen)); 3377 cells_flat = cells[0][0]; 3378 numEdges = 0; 3379 for (PetscInt i = 0; i < numFaces; i++) { 3380 for (PetscInt e = 0; e < 4; e++) { 3381 PetscInt ev[] = {cells_flat[i * 4 + e], cells_flat[i * 4 + ((e + 1) % 4)]}; 3382 const PetscReal *evCoords[] = {&vtxCoords[3 * ev[0]], &vtxCoords[3 * ev[1]]}; 3383 3384 for (PetscInt d = 0; d < 3; d++) { 3385 if (!periodic || periodic[0] != DM_BOUNDARY_PERIODIC) { 3386 if (evCoords[0][d] == 0. && evCoords[1][d] == 0.) numEdges++; 3387 if (evCoords[0][d] == 2. * extent[d] && evCoords[1][d] == 2. * extent[d]) numEdges++; 3388 } 3389 } 3390 } 3391 } 3392 PetscCall(PetscMalloc1(numEdges, &edges)); 3393 PetscCall(PetscMalloc1(numEdges, &edgeSets)); 3394 for (PetscInt edge = 0, i = 0; i < numFaces; i++) { 3395 for (PetscInt e = 0; e < 4; e++) { 3396 PetscInt ev[] = {cells_flat[i * 4 + e], cells_flat[i * 4 + ((e + 1) % 4)]}; 3397 const PetscReal *evCoords[] = {&vtxCoords[3 * ev[0]], &vtxCoords[3 * ev[1]]}; 3398 3399 for (PetscInt d = 0; d < 3; d++) { 3400 if (!periodic || periodic[d] != DM_BOUNDARY_PERIODIC) { 3401 if (evCoords[0][d] == 0. && evCoords[1][d] == 0.) { 3402 edges[edge][0] = ev[0]; 3403 edges[edge][1] = ev[1]; 3404 edgeSets[edge++] = 2 * d; 3405 } 3406 if (evCoords[0][d] == 2. * extent[d] && evCoords[1][d] == 2. * extent[d]) { 3407 edges[edge][0] = ev[0]; 3408 edges[edge][1] = ev[1]; 3409 edgeSets[edge++] = 2 * d + 1; 3410 } 3411 } 3412 } 3413 } 3414 } 3415 } 3416 evalFunc = TPSEvaluate_Gyroid; 3417 normalFunc = TPSExtrudeNormalFunc_Gyroid; 3418 break; 3419 } 3420 3421 PetscCall(DMSetDimension(dm, topoDim)); 3422 if (rank == 0) PetscCall(DMPlexBuildFromCellList(dm, numFaces, numVertices, 4, cells_flat)); 3423 else PetscCall(DMPlexBuildFromCellList(dm, 0, 0, 0, NULL)); 3424 PetscCall(PetscFree(cells_flat)); 3425 { 3426 DM idm; 3427 PetscCall(DMPlexInterpolate(dm, &idm)); 3428 PetscCall(DMPlexReplace_Internal(dm, &idm)); 3429 } 3430 if (rank == 0) PetscCall(DMPlexBuildCoordinatesFromCellList(dm, spaceDim, vtxCoords)); 3431 else PetscCall(DMPlexBuildCoordinatesFromCellList(dm, spaceDim, NULL)); 3432 PetscCall(PetscFree(vtxCoords)); 3433 3434 PetscCall(DMCreateLabel(dm, "Face Sets")); 3435 PetscCall(DMGetLabel(dm, "Face Sets", &label)); 3436 for (PetscInt e = 0; e < numEdges; e++) { 3437 PetscInt njoin; 3438 const PetscInt *join, verts[] = {numFaces + edges[e][0], numFaces + edges[e][1]}; 3439 PetscCall(DMPlexGetJoin(dm, 2, verts, &njoin, &join)); 3440 PetscCheck(njoin == 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Expected unique join of vertices %" PetscInt_FMT " and %" PetscInt_FMT, edges[e][0], edges[e][1]); 3441 PetscCall(DMLabelSetValue(label, join[0], edgeSets[e])); 3442 PetscCall(DMPlexRestoreJoin(dm, 2, verts, &njoin, &join)); 3443 } 3444 PetscCall(PetscFree(edges)); 3445 PetscCall(PetscFree(edgeSets)); 3446 if (tps_distribute) { 3447 DM pdm = NULL; 3448 PetscPartitioner part; 3449 3450 PetscCall(DMPlexGetPartitioner(dm, &part)); 3451 PetscCall(PetscPartitionerSetFromOptions(part)); 3452 PetscCall(DMPlexDistribute(dm, 0, NULL, &pdm)); 3453 if (pdm) PetscCall(DMPlexReplace_Internal(dm, &pdm)); 3454 // Do not auto-distribute again 3455 PetscCall(DMPlexDistributeSetDefault(dm, PETSC_FALSE)); 3456 } 3457 3458 PetscCall(DMPlexSetRefinementUniform(dm, PETSC_TRUE)); 3459 for (PetscInt refine = 0; refine < refinements; refine++) { 3460 PetscInt m; 3461 DM dmf; 3462 Vec X; 3463 PetscScalar *x; 3464 PetscCall(DMRefine(dm, MPI_COMM_NULL, &dmf)); 3465 PetscCall(DMPlexReplace_Internal(dm, &dmf)); 3466 3467 PetscCall(DMGetCoordinatesLocal(dm, &X)); 3468 PetscCall(VecGetLocalSize(X, &m)); 3469 PetscCall(VecGetArray(X, &x)); 3470 for (PetscInt i = 0; i < m; i += 3) PetscCall(TPSNearestPoint(evalFunc, &x[i])); 3471 PetscCall(VecRestoreArray(X, &x)); 3472 } 3473 3474 // Face Sets has already been propagated to new vertices during refinement; this propagates to the initial vertices. 3475 PetscCall(DMGetLabel(dm, "Face Sets", &label)); 3476 PetscCall(DMPlexLabelComplete(dm, label)); 3477 3478 PetscCall(PetscLogEventEnd(DMPLEX_Generate, dm, 0, 0, 0)); 3479 3480 if (thickness > 0) { 3481 DM edm, cdm, ecdm; 3482 DMPlexTransform tr; 3483 const char *prefix; 3484 PetscOptions options; 3485 // Code from DMPlexExtrude 3486 PetscCall(DMPlexTransformCreate(PetscObjectComm((PetscObject)dm), &tr)); 3487 PetscCall(DMPlexTransformSetDM(tr, dm)); 3488 PetscCall(DMPlexTransformSetType(tr, DMPLEXEXTRUDE)); 3489 PetscCall(PetscObjectGetOptionsPrefix((PetscObject)dm, &prefix)); 3490 PetscCall(PetscObjectSetOptionsPrefix((PetscObject)tr, prefix)); 3491 PetscCall(PetscObjectGetOptions((PetscObject)dm, &options)); 3492 PetscCall(PetscObjectSetOptions((PetscObject)tr, options)); 3493 PetscCall(DMPlexTransformExtrudeSetLayers(tr, layers)); 3494 PetscCall(DMPlexTransformExtrudeSetThickness(tr, thickness)); 3495 PetscCall(DMPlexTransformExtrudeSetTensor(tr, PETSC_FALSE)); 3496 PetscCall(DMPlexTransformExtrudeSetSymmetric(tr, PETSC_TRUE)); 3497 PetscCall(DMPlexTransformExtrudeSetNormalFunction(tr, normalFunc)); 3498 PetscCall(DMPlexTransformSetFromOptions(tr)); 3499 PetscCall(PetscObjectSetOptions((PetscObject)tr, NULL)); 3500 PetscCall(DMPlexTransformSetUp(tr)); 3501 PetscCall(PetscObjectViewFromOptions((PetscObject)tr, NULL, "-dm_plex_tps_transform_view")); 3502 PetscCall(DMPlexTransformApply(tr, dm, &edm)); 3503 PetscCall(DMCopyDisc(dm, edm)); 3504 PetscCall(DMGetCoordinateDM(dm, &cdm)); 3505 PetscCall(DMGetCoordinateDM(edm, &ecdm)); 3506 PetscCall(DMCopyDisc(cdm, ecdm)); 3507 PetscCall(DMPlexTransformCreateDiscLabels(tr, edm)); 3508 PetscCall(DMPlexTransformDestroy(&tr)); 3509 if (edm) { 3510 ((DM_Plex *)edm->data)->printFEM = ((DM_Plex *)dm->data)->printFEM; 3511 ((DM_Plex *)edm->data)->printL2 = ((DM_Plex *)dm->data)->printL2; 3512 ((DM_Plex *)edm->data)->printLocate = ((DM_Plex *)dm->data)->printLocate; 3513 } 3514 PetscCall(DMPlexReplace_Internal(dm, &edm)); 3515 } 3516 PetscFunctionReturn(PETSC_SUCCESS); 3517 } 3518 3519 /*@ 3520 DMPlexCreateTPSMesh - Create a distributed, interpolated mesh of a triply-periodic surface 3521 3522 Collective 3523 3524 Input Parameters: 3525 + comm - The communicator for the `DM` object 3526 . tpstype - Type of triply-periodic surface 3527 . extent - Array of length 3 containing number of periods in each direction 3528 . periodic - array of length 3 with periodicity, or `NULL` for non-periodic 3529 . tps_distribute - Distribute 2D manifold mesh prior to refinement and extrusion (more scalable) 3530 . refinements - Number of factor-of-2 refinements of 2D manifold mesh 3531 . layers - Number of cell layers extruded in normal direction 3532 - thickness - Thickness in normal direction 3533 3534 Output Parameter: 3535 . dm - The `DM` object 3536 3537 Level: beginner 3538 3539 Notes: 3540 This meshes the surface of the Schwarz P or Gyroid surfaces. Schwarz P is is the simplest member of the triply-periodic minimal surfaces. 3541 https://en.wikipedia.org/wiki/Schwarz_minimal_surface#Schwarz_P_(%22Primitive%22) and can be cut with "clean" boundaries. 3542 The Gyroid (https://en.wikipedia.org/wiki/Gyroid) is another triply-periodic minimal surface with applications in additive manufacturing; it is much more difficult to "cut" since there are no planes of symmetry. 3543 Our implementation creates a very coarse mesh of the surface and refines (by 4-way splitting) as many times as requested. 3544 On each refinement, all vertices are projected to their nearest point on the surface. 3545 This projection could readily be extended to related surfaces. 3546 3547 The face (edge) sets for the Schwarz P surface are numbered 1(-x), 2(+x), 3(-y), 4(+y), 5(-z), 6(+z). 3548 When the mesh is refined, "Face Sets" contain the new vertices (created during refinement). Use `DMPlexLabelComplete()` to propagate to coarse-level vertices. 3549 3550 Developer Notes: 3551 The Gyroid mesh does not currently mark boundary sets. 3552 3553 References: 3554 . * - Maskery et al, Insights into the mechanical properties of several triply periodic minimal surface lattice structures made by polymer additive manufacturing, 2017. 3555 https://doi.org/10.1016/j.polymer.2017.11.049 3556 3557 .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexCreateSphereMesh()`, `DMSetType()`, `DMCreate()` 3558 @*/ 3559 PetscErrorCode DMPlexCreateTPSMesh(MPI_Comm comm, DMPlexTPSType tpstype, const PetscInt extent[], const DMBoundaryType periodic[], PetscBool tps_distribute, PetscInt refinements, PetscInt layers, PetscReal thickness, DM *dm) 3560 { 3561 PetscFunctionBegin; 3562 PetscCall(DMCreate(comm, dm)); 3563 PetscCall(DMSetType(*dm, DMPLEX)); 3564 PetscCall(DMPlexCreateTPSMesh_Internal(*dm, tpstype, extent, periodic, tps_distribute, refinements, layers, thickness)); 3565 PetscFunctionReturn(PETSC_SUCCESS); 3566 } 3567 3568 /*@ 3569 DMPlexCreateSphereMesh - Creates a mesh on the d-dimensional sphere, S^d. 3570 3571 Collective 3572 3573 Input Parameters: 3574 + comm - The communicator for the `DM` object 3575 . dim - The dimension 3576 . simplex - Use simplices, or tensor product cells 3577 - R - The radius 3578 3579 Output Parameter: 3580 . dm - The `DM` object 3581 3582 Level: beginner 3583 3584 .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexCreateBallMesh()`, `DMPlexCreateBoxMesh()`, `DMSetType()`, `DMCreate()` 3585 @*/ 3586 PetscErrorCode DMPlexCreateSphereMesh(MPI_Comm comm, PetscInt dim, PetscBool simplex, PetscReal R, DM *dm) 3587 { 3588 PetscFunctionBegin; 3589 PetscAssertPointer(dm, 5); 3590 PetscCall(DMCreate(comm, dm)); 3591 PetscCall(DMSetType(*dm, DMPLEX)); 3592 PetscCall(DMPlexCreateSphereMesh_Internal(*dm, dim, simplex, R)); 3593 PetscFunctionReturn(PETSC_SUCCESS); 3594 } 3595 3596 static PetscErrorCode DMPlexCreateBallMesh_Internal(DM dm, PetscInt dim, PetscReal R) 3597 { 3598 DM sdm, vol; 3599 DMLabel bdlabel; 3600 3601 PetscFunctionBegin; 3602 PetscCall(DMCreate(PetscObjectComm((PetscObject)dm), &sdm)); 3603 PetscCall(DMSetType(sdm, DMPLEX)); 3604 PetscCall(PetscObjectSetOptionsPrefix((PetscObject)sdm, "bd_")); 3605 PetscCall(DMPlexCreateSphereMesh_Internal(sdm, dim - 1, PETSC_TRUE, R)); 3606 PetscCall(DMSetFromOptions(sdm)); 3607 PetscCall(DMViewFromOptions(sdm, NULL, "-dm_view")); 3608 PetscCall(DMPlexGenerate(sdm, NULL, PETSC_TRUE, &vol)); 3609 PetscCall(DMDestroy(&sdm)); 3610 PetscCall(DMPlexReplace_Internal(dm, &vol)); 3611 PetscCall(DMCreateLabel(dm, "marker")); 3612 PetscCall(DMGetLabel(dm, "marker", &bdlabel)); 3613 PetscCall(DMPlexMarkBoundaryFaces(dm, PETSC_DETERMINE, bdlabel)); 3614 PetscCall(DMPlexLabelComplete(dm, bdlabel)); 3615 PetscFunctionReturn(PETSC_SUCCESS); 3616 } 3617 3618 /*@ 3619 DMPlexCreateBallMesh - Creates a simplex mesh on the d-dimensional ball, B^d. 3620 3621 Collective 3622 3623 Input Parameters: 3624 + comm - The communicator for the `DM` object 3625 . dim - The dimension 3626 - R - The radius 3627 3628 Output Parameter: 3629 . dm - The `DM` object 3630 3631 Options Database Key: 3632 . bd_dm_refine - This will refine the surface mesh preserving the sphere geometry 3633 3634 Level: beginner 3635 3636 .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexCreateSphereMesh()`, `DMPlexCreateBoxMesh()`, `DMSetType()`, `DMCreate()` 3637 @*/ 3638 PetscErrorCode DMPlexCreateBallMesh(MPI_Comm comm, PetscInt dim, PetscReal R, DM *dm) 3639 { 3640 PetscFunctionBegin; 3641 PetscCall(DMCreate(comm, dm)); 3642 PetscCall(DMSetType(*dm, DMPLEX)); 3643 PetscCall(DMPlexCreateBallMesh_Internal(*dm, dim, R)); 3644 PetscFunctionReturn(PETSC_SUCCESS); 3645 } 3646 3647 static PetscErrorCode DMPlexCreateReferenceCell_Internal(DM rdm, DMPolytopeType ct) 3648 { 3649 PetscFunctionBegin; 3650 switch (ct) { 3651 case DM_POLYTOPE_POINT: { 3652 PetscInt numPoints[1] = {1}; 3653 PetscInt coneSize[1] = {0}; 3654 PetscInt cones[1] = {0}; 3655 PetscInt coneOrientations[1] = {0}; 3656 PetscScalar vertexCoords[1] = {0.0}; 3657 3658 PetscCall(DMSetDimension(rdm, 0)); 3659 PetscCall(DMPlexCreateFromDAG(rdm, 0, numPoints, coneSize, cones, coneOrientations, vertexCoords)); 3660 } break; 3661 case DM_POLYTOPE_SEGMENT: { 3662 PetscInt numPoints[2] = {2, 1}; 3663 PetscInt coneSize[3] = {2, 0, 0}; 3664 PetscInt cones[2] = {1, 2}; 3665 PetscInt coneOrientations[2] = {0, 0}; 3666 PetscScalar vertexCoords[2] = {-1.0, 1.0}; 3667 3668 PetscCall(DMSetDimension(rdm, 1)); 3669 PetscCall(DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords)); 3670 } break; 3671 case DM_POLYTOPE_POINT_PRISM_TENSOR: { 3672 PetscInt numPoints[2] = {2, 1}; 3673 PetscInt coneSize[3] = {2, 0, 0}; 3674 PetscInt cones[2] = {1, 2}; 3675 PetscInt coneOrientations[2] = {0, 0}; 3676 PetscScalar vertexCoords[2] = {-1.0, 1.0}; 3677 3678 PetscCall(DMSetDimension(rdm, 1)); 3679 PetscCall(DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords)); 3680 } break; 3681 case DM_POLYTOPE_TRIANGLE: { 3682 PetscInt numPoints[2] = {3, 1}; 3683 PetscInt coneSize[4] = {3, 0, 0, 0}; 3684 PetscInt cones[3] = {1, 2, 3}; 3685 PetscInt coneOrientations[3] = {0, 0, 0}; 3686 PetscScalar vertexCoords[6] = {-1.0, -1.0, 1.0, -1.0, -1.0, 1.0}; 3687 3688 PetscCall(DMSetDimension(rdm, 2)); 3689 PetscCall(DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords)); 3690 } break; 3691 case DM_POLYTOPE_QUADRILATERAL: { 3692 PetscInt numPoints[2] = {4, 1}; 3693 PetscInt coneSize[5] = {4, 0, 0, 0, 0}; 3694 PetscInt cones[4] = {1, 2, 3, 4}; 3695 PetscInt coneOrientations[4] = {0, 0, 0, 0}; 3696 PetscScalar vertexCoords[8] = {-1.0, -1.0, 1.0, -1.0, 1.0, 1.0, -1.0, 1.0}; 3697 3698 PetscCall(DMSetDimension(rdm, 2)); 3699 PetscCall(DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords)); 3700 } break; 3701 case DM_POLYTOPE_SEG_PRISM_TENSOR: { 3702 PetscInt numPoints[2] = {4, 1}; 3703 PetscInt coneSize[5] = {4, 0, 0, 0, 0}; 3704 PetscInt cones[4] = {1, 2, 3, 4}; 3705 PetscInt coneOrientations[4] = {0, 0, 0, 0}; 3706 PetscScalar vertexCoords[8] = {-1.0, -1.0, 1.0, -1.0, -1.0, 1.0, 1.0, 1.0}; 3707 3708 PetscCall(DMSetDimension(rdm, 2)); 3709 PetscCall(DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords)); 3710 } break; 3711 case DM_POLYTOPE_TETRAHEDRON: { 3712 PetscInt numPoints[2] = {4, 1}; 3713 PetscInt coneSize[5] = {4, 0, 0, 0, 0}; 3714 PetscInt cones[4] = {1, 2, 3, 4}; 3715 PetscInt coneOrientations[4] = {0, 0, 0, 0}; 3716 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}; 3717 3718 PetscCall(DMSetDimension(rdm, 3)); 3719 PetscCall(DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords)); 3720 } break; 3721 case DM_POLYTOPE_HEXAHEDRON: { 3722 PetscInt numPoints[2] = {8, 1}; 3723 PetscInt coneSize[9] = {8, 0, 0, 0, 0, 0, 0, 0, 0}; 3724 PetscInt cones[8] = {1, 2, 3, 4, 5, 6, 7, 8}; 3725 PetscInt coneOrientations[8] = {0, 0, 0, 0, 0, 0, 0, 0}; 3726 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, -1.0, -1.0, 1.0, 1.0, -1.0, 1.0, 1.0, 1.0, 1.0, -1.0, 1.0, 1.0}; 3727 3728 PetscCall(DMSetDimension(rdm, 3)); 3729 PetscCall(DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords)); 3730 } break; 3731 case DM_POLYTOPE_TRI_PRISM: { 3732 PetscInt numPoints[2] = {6, 1}; 3733 PetscInt coneSize[7] = {6, 0, 0, 0, 0, 0, 0}; 3734 PetscInt cones[6] = {1, 2, 3, 4, 5, 6}; 3735 PetscInt coneOrientations[6] = {0, 0, 0, 0, 0, 0}; 3736 PetscScalar vertexCoords[18] = {-1.0, -1.0, -1.0, -1.0, 1.0, -1.0, 1.0, -1.0, -1.0, -1.0, -1.0, 1.0, 1.0, -1.0, 1.0, -1.0, 1.0, 1.0}; 3737 3738 PetscCall(DMSetDimension(rdm, 3)); 3739 PetscCall(DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords)); 3740 } break; 3741 case DM_POLYTOPE_TRI_PRISM_TENSOR: { 3742 PetscInt numPoints[2] = {6, 1}; 3743 PetscInt coneSize[7] = {6, 0, 0, 0, 0, 0, 0}; 3744 PetscInt cones[6] = {1, 2, 3, 4, 5, 6}; 3745 PetscInt coneOrientations[6] = {0, 0, 0, 0, 0, 0}; 3746 PetscScalar vertexCoords[18] = {-1.0, -1.0, -1.0, 1.0, -1.0, -1.0, -1.0, 1.0, -1.0, -1.0, -1.0, 1.0, 1.0, -1.0, 1.0, -1.0, 1.0, 1.0}; 3747 3748 PetscCall(DMSetDimension(rdm, 3)); 3749 PetscCall(DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords)); 3750 } break; 3751 case DM_POLYTOPE_QUAD_PRISM_TENSOR: { 3752 PetscInt numPoints[2] = {8, 1}; 3753 PetscInt coneSize[9] = {8, 0, 0, 0, 0, 0, 0, 0, 0}; 3754 PetscInt cones[8] = {1, 2, 3, 4, 5, 6, 7, 8}; 3755 PetscInt coneOrientations[8] = {0, 0, 0, 0, 0, 0, 0, 0}; 3756 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, -1.0, -1.0, 1.0, 1.0, -1.0, 1.0, 1.0, 1.0, 1.0, -1.0, 1.0, 1.0}; 3757 3758 PetscCall(DMSetDimension(rdm, 3)); 3759 PetscCall(DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords)); 3760 } break; 3761 case DM_POLYTOPE_PYRAMID: { 3762 PetscInt numPoints[2] = {5, 1}; 3763 PetscInt coneSize[6] = {5, 0, 0, 0, 0, 0}; 3764 PetscInt cones[5] = {1, 2, 3, 4, 5}; 3765 PetscInt coneOrientations[8] = {0, 0, 0, 0, 0, 0, 0, 0}; 3766 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, 0.0, 0.0, 1.0}; 3767 3768 PetscCall(DMSetDimension(rdm, 3)); 3769 PetscCall(DMPlexCreateFromDAG(rdm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords)); 3770 } break; 3771 default: 3772 SETERRQ(PetscObjectComm((PetscObject)rdm), PETSC_ERR_ARG_WRONG, "Cannot create reference cell for cell type %s", DMPolytopeTypes[ct]); 3773 } 3774 { 3775 PetscInt Nv, v; 3776 3777 /* Must create the celltype label here so that we do not automatically try to compute the types */ 3778 PetscCall(DMCreateLabel(rdm, "celltype")); 3779 PetscCall(DMPlexSetCellType(rdm, 0, ct)); 3780 PetscCall(DMPlexGetChart(rdm, NULL, &Nv)); 3781 for (v = 1; v < Nv; ++v) PetscCall(DMPlexSetCellType(rdm, v, DM_POLYTOPE_POINT)); 3782 } 3783 PetscCall(DMPlexInterpolateInPlace_Internal(rdm)); 3784 PetscCall(PetscObjectSetName((PetscObject)rdm, DMPolytopeTypes[ct])); 3785 PetscFunctionReturn(PETSC_SUCCESS); 3786 } 3787 3788 /*@ 3789 DMPlexCreateReferenceCell - Create a `DMPLEX` with the appropriate FEM reference cell 3790 3791 Collective 3792 3793 Input Parameters: 3794 + comm - The communicator 3795 - ct - The cell type of the reference cell 3796 3797 Output Parameter: 3798 . refdm - The reference cell 3799 3800 Level: intermediate 3801 3802 .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexCreateBoxMesh()` 3803 @*/ 3804 PetscErrorCode DMPlexCreateReferenceCell(MPI_Comm comm, DMPolytopeType ct, DM *refdm) 3805 { 3806 PetscFunctionBegin; 3807 PetscCall(DMCreate(comm, refdm)); 3808 PetscCall(DMSetType(*refdm, DMPLEX)); 3809 PetscCall(DMPlexCreateReferenceCell_Internal(*refdm, ct)); 3810 PetscFunctionReturn(PETSC_SUCCESS); 3811 } 3812 3813 static PetscErrorCode DMPlexCreateBoundaryLabel_Private(DM dm, const char name[]) 3814 { 3815 DM plex; 3816 DMLabel label; 3817 PetscBool hasLabel; 3818 3819 PetscFunctionBegin; 3820 PetscCall(DMHasLabel(dm, name, &hasLabel)); 3821 if (hasLabel) PetscFunctionReturn(PETSC_SUCCESS); 3822 PetscCall(DMCreateLabel(dm, name)); 3823 PetscCall(DMGetLabel(dm, name, &label)); 3824 PetscCall(DMConvert(dm, DMPLEX, &plex)); 3825 PetscCall(DMPlexMarkBoundaryFaces(plex, 1, label)); 3826 PetscCall(DMPlexLabelComplete(plex, label)); 3827 PetscCall(DMDestroy(&plex)); 3828 PetscFunctionReturn(PETSC_SUCCESS); 3829 } 3830 3831 /* 3832 We use the last coordinate as the radius, the inner radius is lower[dim-1] and the outer radius is upper[dim-1]. Then we map the first coordinate around the circle. 3833 3834 (x, y) -> (r, theta) = (x[1], (x[0] - lower[0]) * 2\pi/(upper[0] - lower[0])) 3835 */ 3836 static void boxToAnnulus(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 3837 { 3838 const PetscReal low = PetscRealPart(constants[0]); 3839 const PetscReal upp = PetscRealPart(constants[1]); 3840 const PetscReal r = PetscRealPart(u[1]); 3841 const PetscReal th = 2. * PETSC_PI * (PetscRealPart(u[0]) - low) / (upp - low); 3842 3843 f0[0] = r * PetscCosReal(th); 3844 f0[1] = r * PetscSinReal(th); 3845 } 3846 3847 const char *const DMPlexShapes[] = {"box", "box_surface", "ball", "sphere", "cylinder", "schwarz_p", "gyroid", "doublet", "annulus", "hypercubic", "zbox", "unknown", "DMPlexShape", "DM_SHAPE_", NULL}; 3848 3849 static PetscErrorCode DMPlexCreateFromOptions_Internal(PetscOptionItems *PetscOptionsObject, PetscBool *useCoordSpace, DM dm) 3850 { 3851 DMPlexShape shape = DM_SHAPE_BOX; 3852 DMPolytopeType cell = DM_POLYTOPE_TRIANGLE; 3853 PetscInt dim = 2; 3854 PetscBool simplex = PETSC_TRUE, interpolate = PETSC_TRUE, adjCone = PETSC_FALSE, adjClosure = PETSC_TRUE, refDomain = PETSC_FALSE; 3855 PetscBool flg, flg2, fflg, bdfflg, nameflg; 3856 MPI_Comm comm; 3857 char filename[PETSC_MAX_PATH_LEN] = "<unspecified>"; 3858 char bdFilename[PETSC_MAX_PATH_LEN] = "<unspecified>"; 3859 char plexname[PETSC_MAX_PATH_LEN] = ""; 3860 3861 PetscFunctionBegin; 3862 PetscCall(PetscLogEventBegin(DMPLEX_CreateFromOptions, dm, 0, 0, 0)); 3863 PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 3864 /* TODO Turn this into a registration interface */ 3865 PetscCall(PetscOptionsString("-dm_plex_filename", "File containing a mesh", "DMPlexCreateFromFile", filename, filename, sizeof(filename), &fflg)); 3866 PetscCall(PetscOptionsString("-dm_plex_boundary_filename", "File containing a mesh boundary", "DMPlexCreateFromFile", bdFilename, bdFilename, sizeof(bdFilename), &bdfflg)); 3867 PetscCall(PetscOptionsString("-dm_plex_name", "Name of the mesh in the file", "DMPlexCreateFromFile", plexname, plexname, sizeof(plexname), &nameflg)); 3868 PetscCall(PetscOptionsEnum("-dm_plex_cell", "Cell shape", "", DMPolytopeTypes, (PetscEnum)cell, (PetscEnum *)&cell, NULL)); 3869 PetscCall(PetscOptionsBool("-dm_plex_reference_cell_domain", "Use a reference cell domain", "", refDomain, &refDomain, NULL)); 3870 PetscCall(PetscOptionsEnum("-dm_plex_shape", "Shape for built-in mesh", "", DMPlexShapes, (PetscEnum)shape, (PetscEnum *)&shape, &flg)); 3871 PetscCall(PetscOptionsBoundedInt("-dm_plex_dim", "Topological dimension of the mesh", "DMGetDimension", dim, &dim, &flg, 0)); 3872 PetscCheck(dim >= 0, comm, PETSC_ERR_ARG_OUTOFRANGE, "Dimension %" PetscInt_FMT " should be in [0, infinity)", dim); 3873 PetscCall(PetscOptionsBool("-dm_plex_simplex", "Mesh cell shape", "", simplex, &simplex, &flg)); 3874 PetscCall(PetscOptionsBool("-dm_plex_interpolate", "Flag to create edges and faces automatically", "", interpolate, &interpolate, &flg)); 3875 PetscCall(PetscOptionsBool("-dm_plex_adj_cone", "Set adjacency direction", "DMSetBasicAdjacency", adjCone, &adjCone, &flg)); 3876 PetscCall(PetscOptionsBool("-dm_plex_adj_closure", "Set adjacency size", "DMSetBasicAdjacency", adjClosure, &adjClosure, &flg2)); 3877 if (flg || flg2) PetscCall(DMSetBasicAdjacency(dm, adjCone, adjClosure)); 3878 3879 switch (cell) { 3880 case DM_POLYTOPE_POINT: 3881 case DM_POLYTOPE_SEGMENT: 3882 case DM_POLYTOPE_POINT_PRISM_TENSOR: 3883 case DM_POLYTOPE_TRIANGLE: 3884 case DM_POLYTOPE_QUADRILATERAL: 3885 case DM_POLYTOPE_TETRAHEDRON: 3886 case DM_POLYTOPE_HEXAHEDRON: 3887 *useCoordSpace = PETSC_TRUE; 3888 break; 3889 default: 3890 *useCoordSpace = PETSC_FALSE; 3891 break; 3892 } 3893 3894 if (fflg) { 3895 DM dmnew; 3896 3897 PetscCall(DMPlexCreateFromFile(PetscObjectComm((PetscObject)dm), filename, plexname, interpolate, &dmnew)); 3898 PetscCall(DMPlexCopy_Internal(dm, PETSC_FALSE, PETSC_FALSE, dmnew)); 3899 PetscCall(DMPlexReplace_Internal(dm, &dmnew)); 3900 } else if (refDomain) { 3901 PetscCall(DMPlexCreateReferenceCell_Internal(dm, cell)); 3902 } else if (bdfflg) { 3903 DM bdm, dmnew; 3904 3905 PetscCall(DMPlexCreateFromFile(PetscObjectComm((PetscObject)dm), bdFilename, plexname, interpolate, &bdm)); 3906 PetscCall(PetscObjectSetOptionsPrefix((PetscObject)bdm, "bd_")); 3907 PetscCall(DMSetFromOptions(bdm)); 3908 PetscCall(DMPlexGenerate(bdm, NULL, interpolate, &dmnew)); 3909 PetscCall(DMDestroy(&bdm)); 3910 PetscCall(DMPlexCopy_Internal(dm, PETSC_FALSE, PETSC_FALSE, dmnew)); 3911 PetscCall(DMPlexReplace_Internal(dm, &dmnew)); 3912 } else { 3913 PetscCall(PetscObjectSetName((PetscObject)dm, DMPlexShapes[shape])); 3914 switch (shape) { 3915 case DM_SHAPE_BOX: 3916 case DM_SHAPE_ZBOX: 3917 case DM_SHAPE_ANNULUS: { 3918 PetscInt faces[3] = {0, 0, 0}; 3919 PetscReal lower[3] = {0, 0, 0}; 3920 PetscReal upper[3] = {1, 1, 1}; 3921 DMBoundaryType bdt[3] = {DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE}; 3922 PetscBool isAnnular = shape == DM_SHAPE_ANNULUS ? PETSC_TRUE : PETSC_FALSE; 3923 PetscInt i, n; 3924 3925 n = dim; 3926 for (i = 0; i < dim; ++i) faces[i] = (dim == 1 ? 1 : 4 - dim); 3927 PetscCall(PetscOptionsIntArray("-dm_plex_box_faces", "Number of faces along each dimension", "", faces, &n, &flg)); 3928 n = 3; 3929 PetscCall(PetscOptionsRealArray("-dm_plex_box_lower", "Lower left corner of box", "", lower, &n, &flg)); 3930 PetscCheck(!flg || !(n != dim), comm, PETSC_ERR_ARG_SIZ, "Lower box point had %" PetscInt_FMT " values, should have been %" PetscInt_FMT, n, dim); 3931 n = 3; 3932 PetscCall(PetscOptionsRealArray("-dm_plex_box_upper", "Upper right corner of box", "", upper, &n, &flg)); 3933 PetscCheck(!flg || !(n != dim), comm, PETSC_ERR_ARG_SIZ, "Upper box point had %" PetscInt_FMT " values, should have been %" PetscInt_FMT, n, dim); 3934 n = 3; 3935 PetscCall(PetscOptionsEnumArray("-dm_plex_box_bd", "Boundary type for each dimension", "", DMBoundaryTypes, (PetscEnum *)bdt, &n, &flg)); 3936 PetscCheck(!flg || !(n != dim), comm, PETSC_ERR_ARG_SIZ, "Box boundary types had %" PetscInt_FMT " values, should have been %" PetscInt_FMT, n, dim); 3937 3938 PetscCheck(!isAnnular || dim == 2, comm, PETSC_ERR_ARG_OUTOFRANGE, "Only two dimensional annuli have been implemented"); 3939 if (isAnnular) 3940 for (i = 0; i < dim - 1; ++i) bdt[i] = DM_BOUNDARY_PERIODIC; 3941 3942 switch (cell) { 3943 case DM_POLYTOPE_TRI_PRISM_TENSOR: 3944 PetscCall(DMPlexCreateWedgeBoxMesh_Internal(dm, faces, lower, upper, bdt)); 3945 if (!interpolate) { 3946 DM udm; 3947 3948 PetscCall(DMPlexUninterpolate(dm, &udm)); 3949 PetscCall(DMPlexReplace_Internal(dm, &udm)); 3950 } 3951 break; 3952 default: 3953 PetscCall(DMPlexCreateBoxMesh_Internal(dm, shape, dim, simplex, faces, lower, upper, bdt, interpolate)); 3954 break; 3955 } 3956 if (isAnnular) { 3957 DM cdm; 3958 PetscDS cds; 3959 PetscScalar bounds[2] = {lower[0], upper[0]}; 3960 3961 // Fix coordinates for annular region 3962 PetscCall(DMSetPeriodicity(dm, NULL, NULL, NULL)); 3963 PetscCall(DMSetCellCoordinatesLocal(dm, NULL)); 3964 PetscCall(DMSetCellCoordinates(dm, NULL)); 3965 PetscCall(DMPlexCreateCoordinateSpace(dm, 1, PETSC_TRUE, NULL)); 3966 PetscCall(DMGetCoordinateDM(dm, &cdm)); 3967 PetscCall(DMGetDS(cdm, &cds)); 3968 PetscCall(PetscDSSetConstants(cds, 2, bounds)); 3969 PetscCall(DMPlexRemapGeometry(dm, 0.0, boxToAnnulus)); 3970 } 3971 } break; 3972 case DM_SHAPE_BOX_SURFACE: { 3973 PetscInt faces[3] = {0, 0, 0}; 3974 PetscReal lower[3] = {0, 0, 0}; 3975 PetscReal upper[3] = {1, 1, 1}; 3976 PetscInt i, n; 3977 3978 n = dim + 1; 3979 for (i = 0; i < dim + 1; ++i) faces[i] = (dim + 1 == 1 ? 1 : 4 - (dim + 1)); 3980 PetscCall(PetscOptionsIntArray("-dm_plex_box_faces", "Number of faces along each dimension", "", faces, &n, &flg)); 3981 n = 3; 3982 PetscCall(PetscOptionsRealArray("-dm_plex_box_lower", "Lower left corner of box", "", lower, &n, &flg)); 3983 PetscCheck(!flg || !(n != dim + 1), comm, PETSC_ERR_ARG_SIZ, "Lower box point had %" PetscInt_FMT " values, should have been %" PetscInt_FMT, n, dim + 1); 3984 n = 3; 3985 PetscCall(PetscOptionsRealArray("-dm_plex_box_upper", "Upper right corner of box", "", upper, &n, &flg)); 3986 PetscCheck(!flg || !(n != dim + 1), comm, PETSC_ERR_ARG_SIZ, "Upper box point had %" PetscInt_FMT " values, should have been %" PetscInt_FMT, n, dim + 1); 3987 PetscCall(DMPlexCreateBoxSurfaceMesh_Internal(dm, dim + 1, faces, lower, upper, interpolate)); 3988 } break; 3989 case DM_SHAPE_SPHERE: { 3990 PetscReal R = 1.0; 3991 3992 PetscCall(PetscOptionsReal("-dm_plex_sphere_radius", "Radius of the sphere", "", R, &R, &flg)); 3993 PetscCall(DMPlexCreateSphereMesh_Internal(dm, dim, simplex, R)); 3994 } break; 3995 case DM_SHAPE_BALL: { 3996 PetscReal R = 1.0; 3997 3998 PetscCall(PetscOptionsReal("-dm_plex_ball_radius", "Radius of the ball", "", R, &R, &flg)); 3999 PetscCall(DMPlexCreateBallMesh_Internal(dm, dim, R)); 4000 } break; 4001 case DM_SHAPE_CYLINDER: { 4002 DMBoundaryType bdt = DM_BOUNDARY_NONE; 4003 PetscInt Nw = 6; 4004 4005 PetscCall(PetscOptionsEnum("-dm_plex_cylinder_bd", "Boundary type in the z direction", "", DMBoundaryTypes, (PetscEnum)bdt, (PetscEnum *)&bdt, NULL)); 4006 PetscCall(PetscOptionsInt("-dm_plex_cylinder_num_wedges", "Number of wedges around the cylinder", "", Nw, &Nw, NULL)); 4007 switch (cell) { 4008 case DM_POLYTOPE_TRI_PRISM_TENSOR: 4009 PetscCall(DMPlexCreateWedgeCylinderMesh_Internal(dm, Nw, interpolate)); 4010 break; 4011 default: 4012 PetscCall(DMPlexCreateHexCylinderMesh_Internal(dm, bdt)); 4013 break; 4014 } 4015 } break; 4016 case DM_SHAPE_SCHWARZ_P: // fallthrough 4017 case DM_SHAPE_GYROID: { 4018 PetscInt extent[3] = {1, 1, 1}, refine = 0, layers = 0, three; 4019 PetscReal thickness = 0.; 4020 DMBoundaryType periodic[3] = {DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE}; 4021 DMPlexTPSType tps_type = shape == DM_SHAPE_SCHWARZ_P ? DMPLEX_TPS_SCHWARZ_P : DMPLEX_TPS_GYROID; 4022 PetscBool tps_distribute; 4023 PetscCall(PetscOptionsIntArray("-dm_plex_tps_extent", "Number of replicas for each of three dimensions", NULL, extent, (three = 3, &three), NULL)); 4024 PetscCall(PetscOptionsInt("-dm_plex_tps_refine", "Number of refinements", NULL, refine, &refine, NULL)); 4025 PetscCall(PetscOptionsEnumArray("-dm_plex_tps_periodic", "Periodicity in each of three dimensions", NULL, DMBoundaryTypes, (PetscEnum *)periodic, (three = 3, &three), NULL)); 4026 PetscCall(PetscOptionsInt("-dm_plex_tps_layers", "Number of layers in volumetric extrusion (or zero to not extrude)", NULL, layers, &layers, NULL)); 4027 PetscCall(PetscOptionsReal("-dm_plex_tps_thickness", "Thickness of volumetric extrusion", NULL, thickness, &thickness, NULL)); 4028 PetscCall(DMPlexDistributeGetDefault(dm, &tps_distribute)); 4029 PetscCall(PetscOptionsBool("-dm_plex_tps_distribute", "Distribute the 2D mesh prior to refinement and extrusion", NULL, tps_distribute, &tps_distribute, NULL)); 4030 PetscCall(DMPlexCreateTPSMesh_Internal(dm, tps_type, extent, periodic, tps_distribute, refine, layers, thickness)); 4031 } break; 4032 case DM_SHAPE_DOUBLET: { 4033 DM dmnew; 4034 PetscReal rl = 0.0; 4035 4036 PetscCall(PetscOptionsReal("-dm_plex_doublet_refinementlimit", "Refinement limit", NULL, rl, &rl, NULL)); 4037 PetscCall(DMPlexCreateDoublet(PetscObjectComm((PetscObject)dm), dim, simplex, interpolate, rl, &dmnew)); 4038 PetscCall(DMPlexCopy_Internal(dm, PETSC_FALSE, PETSC_FALSE, dmnew)); 4039 PetscCall(DMPlexReplace_Internal(dm, &dmnew)); 4040 } break; 4041 case DM_SHAPE_HYPERCUBIC: { 4042 PetscInt *edges; 4043 PetscReal *lower, *upper; 4044 DMBoundaryType *bdt; 4045 PetscInt n, d; 4046 4047 *useCoordSpace = PETSC_FALSE; 4048 PetscCall(PetscMalloc4(dim, &edges, dim, &lower, dim, &upper, dim, &bdt)); 4049 for (d = 0; d < dim; ++d) { 4050 edges[d] = 1; 4051 lower[d] = 0.; 4052 upper[d] = 1.; 4053 bdt[d] = DM_BOUNDARY_PERIODIC; 4054 } 4055 n = dim; 4056 PetscCall(PetscOptionsIntArray("-dm_plex_box_faces", "Number of faces along each dimension", "", edges, &n, &flg)); 4057 n = dim; 4058 PetscCall(PetscOptionsRealArray("-dm_plex_box_lower", "Lower left corner of box", "", lower, &n, &flg)); 4059 PetscCheck(!flg || n == dim, comm, PETSC_ERR_ARG_SIZ, "Lower box point had %" PetscInt_FMT " values, should have been %" PetscInt_FMT, n, dim); 4060 n = dim; 4061 PetscCall(PetscOptionsRealArray("-dm_plex_box_upper", "Upper right corner of box", "", upper, &n, &flg)); 4062 PetscCheck(!flg || n == dim, comm, PETSC_ERR_ARG_SIZ, "Upper box point had %" PetscInt_FMT " values, should have been %" PetscInt_FMT, n, dim); 4063 n = dim; 4064 PetscCall(PetscOptionsEnumArray("-dm_plex_box_bd", "Boundary type for each dimension", "", DMBoundaryTypes, (PetscEnum *)bdt, &n, &flg)); 4065 PetscCheck(!flg || n == dim, comm, PETSC_ERR_ARG_SIZ, "Box boundary types had %" PetscInt_FMT " values, should have been %" PetscInt_FMT, n, dim); 4066 PetscCall(DMPlexCreateHypercubicMesh_Internal(dm, dim, lower, upper, edges, bdt)); 4067 PetscCall(PetscFree4(edges, lower, upper, bdt)); 4068 } break; 4069 default: 4070 SETERRQ(comm, PETSC_ERR_SUP, "Domain shape %s is unsupported", DMPlexShapes[shape]); 4071 } 4072 } 4073 PetscCall(DMPlexSetRefinementUniform(dm, PETSC_TRUE)); 4074 if (!((PetscObject)dm)->name && nameflg) PetscCall(PetscObjectSetName((PetscObject)dm, plexname)); 4075 PetscCall(PetscLogEventEnd(DMPLEX_CreateFromOptions, dm, 0, 0, 0)); 4076 PetscFunctionReturn(PETSC_SUCCESS); 4077 } 4078 4079 PetscErrorCode DMSetFromOptions_NonRefinement_Plex(DM dm, PetscOptionItems *PetscOptionsObject) 4080 { 4081 DM_Plex *mesh = (DM_Plex *)dm->data; 4082 PetscBool flg, flg2; 4083 char bdLabel[PETSC_MAX_PATH_LEN]; 4084 4085 PetscFunctionBegin; 4086 /* Handle viewing */ 4087 PetscCall(PetscOptionsBool("-dm_plex_print_set_values", "Output all set values info", "DMPlexMatSetClosure", PETSC_FALSE, &mesh->printSetValues, NULL)); 4088 PetscCall(PetscOptionsBoundedInt("-dm_plex_print_fem", "Debug output level for all fem computations", "DMPlexSNESComputeResidualFEM", 0, &mesh->printFEM, NULL, 0)); 4089 PetscCall(PetscOptionsBoundedInt("-dm_plex_print_fvm", "Debug output level for all fvm computations", "DMPlexSNESComputeResidualFVM", 0, &mesh->printFVM, NULL, 0)); 4090 PetscCall(PetscOptionsReal("-dm_plex_print_tol", "Tolerance for FEM output", "DMPlexSNESComputeResidualFEM", mesh->printTol, &mesh->printTol, NULL)); 4091 PetscCall(PetscOptionsBoundedInt("-dm_plex_print_l2", "Debug output level all L2 diff computations", "DMComputeL2Diff", 0, &mesh->printL2, NULL, 0)); 4092 PetscCall(PetscOptionsBoundedInt("-dm_plex_print_locate", "Debug output level all point location computations", "DMLocatePoints", 0, &mesh->printLocate, NULL, 0)); 4093 PetscCall(DMMonitorSetFromOptions(dm, "-dm_plex_monitor_throughput", "Monitor the simulation throughput", "DMPlexMonitorThroughput", DMPlexMonitorThroughput, NULL, &flg)); 4094 if (flg) PetscCall(PetscLogDefaultBegin()); 4095 /* Labeling */ 4096 PetscCall(PetscOptionsString("-dm_plex_boundary_label", "Label to mark the mesh boundary", "", bdLabel, bdLabel, sizeof(bdLabel), &flg)); 4097 if (flg) PetscCall(DMPlexCreateBoundaryLabel_Private(dm, bdLabel)); 4098 /* Point Location */ 4099 PetscCall(PetscOptionsBool("-dm_plex_hash_location", "Use grid hashing for point location", "DMInterpolate", PETSC_FALSE, &mesh->useHashLocation, NULL)); 4100 /* Partitioning and distribution */ 4101 PetscCall(PetscOptionsBool("-dm_plex_partition_balance", "Attempt to evenly divide points on partition boundary between processes", "DMPlexSetPartitionBalance", PETSC_FALSE, &mesh->partitionBalance, NULL)); 4102 /* Generation and remeshing */ 4103 PetscCall(PetscOptionsBool("-dm_plex_remesh_bd", "Allow changes to the boundary on remeshing", "DMAdapt", PETSC_FALSE, &mesh->remeshBd, NULL)); 4104 /* Projection behavior */ 4105 PetscCall(PetscOptionsBoundedInt("-dm_plex_max_projection_height", "Maximum mesh point height used to project locally", "DMPlexSetMaxProjectionHeight", 0, &mesh->maxProjectionHeight, NULL, 0)); 4106 PetscCall(PetscOptionsBool("-dm_plex_regular_refinement", "Use special nested projection algorithm for regular refinement", "DMPlexSetRegularRefinement", mesh->regularRefinement, &mesh->regularRefinement, NULL)); 4107 /* Checking structure */ 4108 { 4109 PetscBool all = PETSC_FALSE; 4110 4111 PetscCall(PetscOptionsBool("-dm_plex_check_all", "Perform all basic checks", "DMPlexCheck", PETSC_FALSE, &all, NULL)); 4112 if (all) { 4113 PetscCall(DMPlexCheck(dm)); 4114 } else { 4115 PetscCall(PetscOptionsBool("-dm_plex_check_symmetry", "Check that the adjacency information in the mesh is symmetric", "DMPlexCheckSymmetry", PETSC_FALSE, &flg, &flg2)); 4116 if (flg && flg2) PetscCall(DMPlexCheckSymmetry(dm)); 4117 PetscCall(PetscOptionsBool("-dm_plex_check_skeleton", "Check that each cell has the correct number of vertices (only for homogeneous simplex or tensor meshes)", "DMPlexCheckSkeleton", PETSC_FALSE, &flg, &flg2)); 4118 if (flg && flg2) PetscCall(DMPlexCheckSkeleton(dm, 0)); 4119 PetscCall(PetscOptionsBool("-dm_plex_check_faces", "Check that the faces of each cell give a vertex order this is consistent with what we expect from the cell type", "DMPlexCheckFaces", PETSC_FALSE, &flg, &flg2)); 4120 if (flg && flg2) PetscCall(DMPlexCheckFaces(dm, 0)); 4121 PetscCall(PetscOptionsBool("-dm_plex_check_geometry", "Check that cells have positive volume", "DMPlexCheckGeometry", PETSC_FALSE, &flg, &flg2)); 4122 if (flg && flg2) PetscCall(DMPlexCheckGeometry(dm)); 4123 PetscCall(PetscOptionsBool("-dm_plex_check_pointsf", "Check some necessary conditions for PointSF", "DMPlexCheckPointSF", PETSC_FALSE, &flg, &flg2)); 4124 if (flg && flg2) PetscCall(DMPlexCheckPointSF(dm, NULL, PETSC_FALSE)); 4125 PetscCall(PetscOptionsBool("-dm_plex_check_interface_cones", "Check points on inter-partition interfaces have conforming order of cone points", "DMPlexCheckInterfaceCones", PETSC_FALSE, &flg, &flg2)); 4126 if (flg && flg2) PetscCall(DMPlexCheckInterfaceCones(dm)); 4127 } 4128 PetscCall(PetscOptionsBool("-dm_plex_check_cell_shape", "Check cell shape", "DMPlexCheckCellShape", PETSC_FALSE, &flg, &flg2)); 4129 if (flg && flg2) PetscCall(DMPlexCheckCellShape(dm, PETSC_TRUE, PETSC_DETERMINE)); 4130 } 4131 { 4132 PetscReal scale = 1.0; 4133 4134 PetscCall(PetscOptionsReal("-dm_plex_scale", "Scale factor for mesh coordinates", "DMPlexScale", scale, &scale, &flg)); 4135 if (flg) { 4136 Vec coordinates, coordinatesLocal; 4137 4138 PetscCall(DMGetCoordinates(dm, &coordinates)); 4139 PetscCall(DMGetCoordinatesLocal(dm, &coordinatesLocal)); 4140 PetscCall(VecScale(coordinates, scale)); 4141 PetscCall(VecScale(coordinatesLocal, scale)); 4142 } 4143 } 4144 PetscCall(PetscPartitionerSetFromOptions(mesh->partitioner)); 4145 PetscFunctionReturn(PETSC_SUCCESS); 4146 } 4147 4148 PetscErrorCode DMSetFromOptions_Overlap_Plex(DM dm, PetscOptionItems *PetscOptionsObject, PetscInt *overlap) 4149 { 4150 PetscInt numOvLabels = 16, numOvExLabels = 16; 4151 char *ovLabelNames[16], *ovExLabelNames[16]; 4152 PetscInt numOvValues = 16, numOvExValues = 16, l; 4153 PetscBool flg; 4154 4155 PetscFunctionBegin; 4156 PetscCall(PetscOptionsBoundedInt("-dm_distribute_overlap", "The size of the overlap halo", "DMPlexDistribute", *overlap, overlap, NULL, 0)); 4157 PetscCall(PetscOptionsStringArray("-dm_distribute_overlap_labels", "List of overlap label names", "DMPlexDistribute", ovLabelNames, &numOvLabels, &flg)); 4158 if (!flg) numOvLabels = 0; 4159 if (numOvLabels) { 4160 ((DM_Plex *)dm->data)->numOvLabels = numOvLabels; 4161 for (l = 0; l < numOvLabels; ++l) { 4162 PetscCall(DMGetLabel(dm, ovLabelNames[l], &((DM_Plex *)dm->data)->ovLabels[l])); 4163 PetscCheck(((DM_Plex *)dm->data)->ovLabels[l], PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid label name %s", ovLabelNames[l]); 4164 PetscCall(PetscFree(ovLabelNames[l])); 4165 } 4166 PetscCall(PetscOptionsIntArray("-dm_distribute_overlap_values", "List of overlap label values", "DMPlexDistribute", ((DM_Plex *)dm->data)->ovValues, &numOvValues, &flg)); 4167 if (!flg) numOvValues = 0; 4168 PetscCheck(numOvLabels == numOvValues, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "The number of labels %" PetscInt_FMT " must match the number of values %" PetscInt_FMT, numOvLabels, numOvValues); 4169 4170 PetscCall(PetscOptionsStringArray("-dm_distribute_overlap_exclude_labels", "List of overlap exclude label names", "DMPlexDistribute", ovExLabelNames, &numOvExLabels, &flg)); 4171 if (!flg) numOvExLabels = 0; 4172 ((DM_Plex *)dm->data)->numOvExLabels = numOvExLabels; 4173 for (l = 0; l < numOvExLabels; ++l) { 4174 PetscCall(DMGetLabel(dm, ovExLabelNames[l], &((DM_Plex *)dm->data)->ovExLabels[l])); 4175 PetscCheck(((DM_Plex *)dm->data)->ovExLabels[l], PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid label name %s", ovExLabelNames[l]); 4176 PetscCall(PetscFree(ovExLabelNames[l])); 4177 } 4178 PetscCall(PetscOptionsIntArray("-dm_distribute_overlap_exclude_values", "List of overlap exclude label values", "DMPlexDistribute", ((DM_Plex *)dm->data)->ovExValues, &numOvExValues, &flg)); 4179 if (!flg) numOvExValues = 0; 4180 PetscCheck(numOvExLabels == numOvExValues, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "The number of exclude labels %" PetscInt_FMT " must match the number of values %" PetscInt_FMT, numOvExLabels, numOvExValues); 4181 } 4182 PetscFunctionReturn(PETSC_SUCCESS); 4183 } 4184 4185 static PetscErrorCode DMSetFromOptions_Plex(DM dm, PetscOptionItems *PetscOptionsObject) 4186 { 4187 PetscFunctionList ordlist; 4188 char oname[256]; 4189 DMPlexReorderDefaultFlag reorder; 4190 PetscReal volume = -1.0; 4191 PetscInt prerefine = 0, refine = 0, r, coarsen = 0, overlap = 0, extLayers = 0, dim; 4192 PetscBool uniformOrig = PETSC_FALSE, created = PETSC_FALSE, uniform = PETSC_TRUE, distribute, saveSF = PETSC_FALSE, interpolate = PETSC_TRUE, coordSpace = PETSC_TRUE, remap = PETSC_TRUE, ghostCells = PETSC_FALSE, isHierarchy, ignoreModel = PETSC_FALSE, flg; 4193 4194 PetscFunctionBegin; 4195 PetscOptionsHeadBegin(PetscOptionsObject, "DMPlex Options"); 4196 if (dm->cloneOpts) goto non_refine; 4197 /* Handle automatic creation */ 4198 PetscCall(DMGetDimension(dm, &dim)); 4199 if (dim < 0) { 4200 PetscCall(DMPlexCreateFromOptions_Internal(PetscOptionsObject, &coordSpace, dm)); 4201 created = PETSC_TRUE; 4202 } 4203 PetscCall(DMGetDimension(dm, &dim)); 4204 /* Handle interpolation before distribution */ 4205 PetscCall(PetscOptionsBool("-dm_plex_interpolate_pre", "Flag to interpolate mesh before distribution", "", interpolate, &interpolate, &flg)); 4206 if (flg) { 4207 DMPlexInterpolatedFlag interpolated; 4208 4209 PetscCall(DMPlexIsInterpolated(dm, &interpolated)); 4210 if (interpolated == DMPLEX_INTERPOLATED_FULL && !interpolate) { 4211 DM udm; 4212 4213 PetscCall(DMPlexUninterpolate(dm, &udm)); 4214 PetscCall(DMPlexReplace_Internal(dm, &udm)); 4215 } else if (interpolated != DMPLEX_INTERPOLATED_FULL && interpolate) { 4216 DM idm; 4217 4218 PetscCall(DMPlexInterpolate(dm, &idm)); 4219 PetscCall(DMPlexReplace_Internal(dm, &idm)); 4220 } 4221 } 4222 /* Handle DMPlex refinement before distribution */ 4223 PetscCall(PetscOptionsBool("-dm_refine_ignore_model", "Flag to ignore the geometry model when refining", "DMCreate", ignoreModel, &ignoreModel, &flg)); 4224 if (flg) ((DM_Plex *)dm->data)->ignoreModel = ignoreModel; 4225 PetscCall(DMPlexGetRefinementUniform(dm, &uniformOrig)); 4226 PetscCall(PetscOptionsBoundedInt("-dm_refine_pre", "The number of refinements before distribution", "DMCreate", prerefine, &prerefine, NULL, 0)); 4227 PetscCall(PetscOptionsBool("-dm_refine_remap_pre", "Flag to control coordinate remapping", "DMCreate", remap, &remap, NULL)); 4228 PetscCall(PetscOptionsBool("-dm_refine_uniform_pre", "Flag for uniform refinement before distribution", "DMCreate", uniform, &uniform, &flg)); 4229 if (flg) PetscCall(DMPlexSetRefinementUniform(dm, uniform)); 4230 PetscCall(PetscOptionsReal("-dm_refine_volume_limit_pre", "The maximum cell volume after refinement before distribution", "DMCreate", volume, &volume, &flg)); 4231 if (flg) { 4232 PetscCall(DMPlexSetRefinementUniform(dm, PETSC_FALSE)); 4233 PetscCall(DMPlexSetRefinementLimit(dm, volume)); 4234 prerefine = PetscMax(prerefine, 1); 4235 } 4236 for (r = 0; r < prerefine; ++r) { 4237 DM rdm; 4238 PetscPointFunc coordFunc = ((DM_Plex *)dm->data)->coordFunc; 4239 4240 PetscCall(DMSetFromOptions_NonRefinement_Plex(dm, PetscOptionsObject)); 4241 PetscCall(DMRefine(dm, PetscObjectComm((PetscObject)dm), &rdm)); 4242 PetscCall(DMPlexReplace_Internal(dm, &rdm)); 4243 PetscCall(DMSetFromOptions_NonRefinement_Plex(dm, PetscOptionsObject)); 4244 if (coordFunc && remap) { 4245 PetscCall(DMPlexRemapGeometry(dm, 0.0, coordFunc)); 4246 ((DM_Plex *)dm->data)->coordFunc = coordFunc; 4247 } 4248 } 4249 PetscCall(DMPlexSetRefinementUniform(dm, uniformOrig)); 4250 /* Handle DMPlex extrusion before distribution */ 4251 PetscCall(PetscOptionsBoundedInt("-dm_extrude", "The number of layers to extrude", "", extLayers, &extLayers, NULL, 0)); 4252 if (extLayers) { 4253 DM edm; 4254 4255 PetscCall(DMExtrude(dm, extLayers, &edm)); 4256 PetscCall(DMPlexReplace_Internal(dm, &edm)); 4257 ((DM_Plex *)dm->data)->coordFunc = NULL; 4258 PetscCall(DMSetFromOptions_NonRefinement_Plex(dm, PetscOptionsObject)); 4259 extLayers = 0; 4260 PetscCall(DMGetDimension(dm, &dim)); 4261 } 4262 /* Handle DMPlex reordering before distribution */ 4263 PetscCall(DMPlexReorderGetDefault(dm, &reorder)); 4264 PetscCall(MatGetOrderingList(&ordlist)); 4265 PetscCall(PetscStrncpy(oname, MATORDERINGNATURAL, sizeof(oname))); 4266 PetscCall(PetscOptionsFList("-dm_plex_reorder", "Set mesh reordering type", "DMPlexGetOrdering", ordlist, MATORDERINGNATURAL, oname, sizeof(oname), &flg)); 4267 if (reorder == DMPLEX_REORDER_DEFAULT_TRUE || flg) { 4268 DM pdm; 4269 IS perm; 4270 4271 PetscCall(DMPlexGetOrdering(dm, oname, NULL, &perm)); 4272 PetscCall(DMPlexPermute(dm, perm, &pdm)); 4273 PetscCall(ISDestroy(&perm)); 4274 PetscCall(DMPlexReplace_Internal(dm, &pdm)); 4275 PetscCall(DMSetFromOptions_NonRefinement_Plex(dm, PetscOptionsObject)); 4276 } 4277 /* Handle DMPlex distribution */ 4278 PetscCall(DMPlexDistributeGetDefault(dm, &distribute)); 4279 PetscCall(PetscOptionsBool("-dm_distribute", "Flag to redistribute a mesh among processes", "DMPlexDistribute", distribute, &distribute, NULL)); 4280 PetscCall(PetscOptionsBool("-dm_distribute_save_sf", "Flag to save the migration SF", "DMPlexSetMigrationSF", saveSF, &saveSF, NULL)); 4281 PetscCall(DMSetFromOptions_Overlap_Plex(dm, PetscOptionsObject, &overlap)); 4282 if (distribute) { 4283 DM pdm = NULL; 4284 PetscPartitioner part; 4285 PetscSF sfMigration; 4286 4287 PetscCall(DMPlexGetPartitioner(dm, &part)); 4288 PetscCall(PetscPartitionerSetFromOptions(part)); 4289 PetscCall(DMPlexDistribute(dm, overlap, &sfMigration, &pdm)); 4290 if (pdm) PetscCall(DMPlexReplace_Internal(dm, &pdm)); 4291 if (saveSF) PetscCall(DMPlexSetMigrationSF(dm, sfMigration)); 4292 PetscCall(PetscSFDestroy(&sfMigration)); 4293 } 4294 /* Must check CEED options before creating function space for coordinates */ 4295 { 4296 PetscBool useCeed = PETSC_FALSE, flg; 4297 4298 PetscCall(PetscOptionsBool("-dm_plex_use_ceed", "Use LibCEED as the FEM backend", "DMPlexSetUseCeed", useCeed, &useCeed, &flg)); 4299 if (flg) PetscCall(DMPlexSetUseCeed(dm, useCeed)); 4300 } 4301 /* Create coordinate space */ 4302 if (created) { 4303 DM_Plex *mesh = (DM_Plex *)dm->data; 4304 PetscInt degree = 1, deg; 4305 PetscInt height = 0; 4306 DM cdm; 4307 PetscBool flg; 4308 4309 PetscCall(PetscOptionsBool("-dm_coord_space", "Use an FEM space for coordinates", "", coordSpace, &coordSpace, &flg)); 4310 PetscCall(PetscOptionsInt("-dm_coord_petscspace_degree", "FEM degree for coordinate space", "", degree, °ree, NULL)); 4311 PetscCall(DMGetCoordinateDegree_Internal(dm, °)); 4312 if (coordSpace && deg <= 1) PetscCall(DMPlexCreateCoordinateSpace(dm, degree, PETSC_TRUE, mesh->coordFunc)); 4313 PetscCall(DMGetCoordinateDM(dm, &cdm)); 4314 if (flg && !coordSpace) { 4315 PetscDS cds; 4316 PetscObject obj; 4317 PetscClassId id; 4318 4319 PetscCall(DMGetDS(cdm, &cds)); 4320 PetscCall(PetscDSGetDiscretization(cds, 0, &obj)); 4321 PetscCall(PetscObjectGetClassId(obj, &id)); 4322 if (id == PETSCFE_CLASSID) { 4323 PetscContainer dummy; 4324 4325 PetscCall(PetscContainerCreate(PETSC_COMM_SELF, &dummy)); 4326 PetscCall(PetscObjectSetName((PetscObject)dummy, "coordinates")); 4327 PetscCall(DMSetField(cdm, 0, NULL, (PetscObject)dummy)); 4328 PetscCall(PetscContainerDestroy(&dummy)); 4329 PetscCall(DMClearDS(cdm)); 4330 } 4331 mesh->coordFunc = NULL; 4332 } 4333 PetscCall(PetscOptionsBool("-dm_sparse_localize", "Localize only necessary cells", "", dm->sparseLocalize, &dm->sparseLocalize, &flg)); 4334 PetscCall(PetscOptionsInt("-dm_localize_height", "Localize edges and faces in addition to cells", "", height, &height, &flg)); 4335 if (flg) PetscCall(DMPlexSetMaxProjectionHeight(cdm, height)); 4336 PetscCall(DMLocalizeCoordinates(dm)); 4337 } 4338 /* Handle DMPlex refinement */ 4339 remap = PETSC_TRUE; 4340 PetscCall(PetscOptionsBoundedInt("-dm_refine", "The number of uniform refinements", "DMCreate", refine, &refine, NULL, 0)); 4341 PetscCall(PetscOptionsBool("-dm_refine_remap", "Flag to control coordinate remapping", "DMCreate", remap, &remap, NULL)); 4342 PetscCall(PetscOptionsBoundedInt("-dm_refine_hierarchy", "The number of uniform refinements", "DMCreate", refine, &refine, &isHierarchy, 0)); 4343 if (refine) PetscCall(DMPlexSetRefinementUniform(dm, PETSC_TRUE)); 4344 if (refine && isHierarchy) { 4345 DM *dms, coarseDM; 4346 4347 PetscCall(DMGetCoarseDM(dm, &coarseDM)); 4348 PetscCall(PetscObjectReference((PetscObject)coarseDM)); 4349 PetscCall(PetscMalloc1(refine, &dms)); 4350 PetscCall(DMRefineHierarchy(dm, refine, dms)); 4351 /* Total hack since we do not pass in a pointer */ 4352 PetscCall(DMPlexSwap_Static(dm, dms[refine - 1])); 4353 if (refine == 1) { 4354 PetscCall(DMSetCoarseDM(dm, dms[0])); 4355 PetscCall(DMPlexSetRegularRefinement(dm, PETSC_TRUE)); 4356 } else { 4357 PetscCall(DMSetCoarseDM(dm, dms[refine - 2])); 4358 PetscCall(DMPlexSetRegularRefinement(dm, PETSC_TRUE)); 4359 PetscCall(DMSetCoarseDM(dms[0], dms[refine - 1])); 4360 PetscCall(DMPlexSetRegularRefinement(dms[0], PETSC_TRUE)); 4361 } 4362 PetscCall(DMSetCoarseDM(dms[refine - 1], coarseDM)); 4363 PetscCall(PetscObjectDereference((PetscObject)coarseDM)); 4364 /* Free DMs */ 4365 for (r = 0; r < refine; ++r) { 4366 PetscCall(DMSetFromOptions_NonRefinement_Plex(dms[r], PetscOptionsObject)); 4367 PetscCall(DMDestroy(&dms[r])); 4368 } 4369 PetscCall(PetscFree(dms)); 4370 } else { 4371 for (r = 0; r < refine; ++r) { 4372 DM rdm; 4373 PetscPointFunc coordFunc = ((DM_Plex *)dm->data)->coordFunc; 4374 4375 PetscCall(DMSetFromOptions_NonRefinement_Plex(dm, PetscOptionsObject)); 4376 PetscCall(DMRefine(dm, PetscObjectComm((PetscObject)dm), &rdm)); 4377 /* Total hack since we do not pass in a pointer */ 4378 PetscCall(DMPlexReplace_Internal(dm, &rdm)); 4379 PetscCall(DMSetFromOptions_NonRefinement_Plex(dm, PetscOptionsObject)); 4380 if (coordFunc && remap) { 4381 PetscCall(DMPlexRemapGeometry(dm, 0.0, coordFunc)); 4382 ((DM_Plex *)dm->data)->coordFunc = coordFunc; 4383 } 4384 } 4385 } 4386 /* Handle DMPlex coarsening */ 4387 PetscCall(PetscOptionsBoundedInt("-dm_coarsen", "Coarsen the mesh", "DMCreate", coarsen, &coarsen, NULL, 0)); 4388 PetscCall(PetscOptionsBoundedInt("-dm_coarsen_hierarchy", "The number of coarsenings", "DMCreate", coarsen, &coarsen, &isHierarchy, 0)); 4389 if (coarsen && isHierarchy) { 4390 DM *dms; 4391 4392 PetscCall(PetscMalloc1(coarsen, &dms)); 4393 PetscCall(DMCoarsenHierarchy(dm, coarsen, dms)); 4394 /* Free DMs */ 4395 for (r = 0; r < coarsen; ++r) { 4396 PetscCall(DMSetFromOptions_NonRefinement_Plex(dms[r], PetscOptionsObject)); 4397 PetscCall(DMDestroy(&dms[r])); 4398 } 4399 PetscCall(PetscFree(dms)); 4400 } else { 4401 for (r = 0; r < coarsen; ++r) { 4402 DM cdm; 4403 PetscPointFunc coordFunc = ((DM_Plex *)dm->data)->coordFunc; 4404 4405 PetscCall(DMSetFromOptions_NonRefinement_Plex(dm, PetscOptionsObject)); 4406 PetscCall(DMCoarsen(dm, PetscObjectComm((PetscObject)dm), &cdm)); 4407 /* Total hack since we do not pass in a pointer */ 4408 PetscCall(DMPlexReplace_Internal(dm, &cdm)); 4409 PetscCall(DMSetFromOptions_NonRefinement_Plex(dm, PetscOptionsObject)); 4410 if (coordFunc) { 4411 PetscCall(DMPlexRemapGeometry(dm, 0.0, coordFunc)); 4412 ((DM_Plex *)dm->data)->coordFunc = coordFunc; 4413 } 4414 } 4415 } 4416 /* Handle ghost cells */ 4417 PetscCall(PetscOptionsBool("-dm_plex_create_fv_ghost_cells", "Flag to create finite volume ghost cells on the boundary", "DMCreate", ghostCells, &ghostCells, NULL)); 4418 if (ghostCells) { 4419 DM gdm; 4420 char lname[PETSC_MAX_PATH_LEN]; 4421 4422 lname[0] = '\0'; 4423 PetscCall(PetscOptionsString("-dm_plex_fv_ghost_cells_label", "Label name for ghost cells boundary", "DMCreate", lname, lname, sizeof(lname), &flg)); 4424 PetscCall(DMPlexConstructGhostCells(dm, flg ? lname : NULL, NULL, &gdm)); 4425 PetscCall(DMPlexReplace_Internal(dm, &gdm)); 4426 } 4427 /* Handle 1D order */ 4428 if (reorder != DMPLEX_REORDER_DEFAULT_FALSE && dim == 1) { 4429 DM cdm, rdm; 4430 PetscDS cds; 4431 PetscObject obj; 4432 PetscClassId id = PETSC_OBJECT_CLASSID; 4433 IS perm; 4434 PetscInt Nf; 4435 PetscBool distributed; 4436 4437 PetscCall(DMPlexIsDistributed(dm, &distributed)); 4438 PetscCall(DMGetCoordinateDM(dm, &cdm)); 4439 PetscCall(DMGetDS(cdm, &cds)); 4440 PetscCall(PetscDSGetNumFields(cds, &Nf)); 4441 if (Nf) { 4442 PetscCall(PetscDSGetDiscretization(cds, 0, &obj)); 4443 PetscCall(PetscObjectGetClassId(obj, &id)); 4444 } 4445 if (!distributed && id != PETSCFE_CLASSID) { 4446 PetscCall(DMPlexGetOrdering1D(dm, &perm)); 4447 PetscCall(DMPlexPermute(dm, perm, &rdm)); 4448 PetscCall(DMPlexReplace_Internal(dm, &rdm)); 4449 PetscCall(ISDestroy(&perm)); 4450 } 4451 } 4452 /* Handle */ 4453 non_refine: 4454 PetscCall(DMSetFromOptions_NonRefinement_Plex(dm, PetscOptionsObject)); 4455 PetscOptionsHeadEnd(); 4456 PetscFunctionReturn(PETSC_SUCCESS); 4457 } 4458 4459 static PetscErrorCode DMCreateGlobalVector_Plex(DM dm, Vec *vec) 4460 { 4461 PetscFunctionBegin; 4462 PetscCall(DMCreateGlobalVector_Section_Private(dm, vec)); 4463 /* PetscCall(VecSetOperation(*vec, VECOP_DUPLICATE, (void(*)(void)) VecDuplicate_MPI_DM)); */ 4464 PetscCall(VecSetOperation(*vec, VECOP_VIEW, (void (*)(void))VecView_Plex)); 4465 PetscCall(VecSetOperation(*vec, VECOP_VIEWNATIVE, (void (*)(void))VecView_Plex_Native)); 4466 PetscCall(VecSetOperation(*vec, VECOP_LOAD, (void (*)(void))VecLoad_Plex)); 4467 PetscCall(VecSetOperation(*vec, VECOP_LOADNATIVE, (void (*)(void))VecLoad_Plex_Native)); 4468 PetscFunctionReturn(PETSC_SUCCESS); 4469 } 4470 4471 static PetscErrorCode DMCreateLocalVector_Plex(DM dm, Vec *vec) 4472 { 4473 PetscFunctionBegin; 4474 PetscCall(DMCreateLocalVector_Section_Private(dm, vec)); 4475 PetscCall(VecSetOperation(*vec, VECOP_VIEW, (void (*)(void))VecView_Plex_Local)); 4476 PetscCall(VecSetOperation(*vec, VECOP_LOAD, (void (*)(void))VecLoad_Plex_Local)); 4477 PetscFunctionReturn(PETSC_SUCCESS); 4478 } 4479 4480 static PetscErrorCode DMGetDimPoints_Plex(DM dm, PetscInt dim, PetscInt *pStart, PetscInt *pEnd) 4481 { 4482 PetscInt depth, d; 4483 4484 PetscFunctionBegin; 4485 PetscCall(DMPlexGetDepth(dm, &depth)); 4486 if (depth == 1) { 4487 PetscCall(DMGetDimension(dm, &d)); 4488 if (dim == 0) PetscCall(DMPlexGetDepthStratum(dm, dim, pStart, pEnd)); 4489 else if (dim == d) PetscCall(DMPlexGetDepthStratum(dm, 1, pStart, pEnd)); 4490 else { 4491 *pStart = 0; 4492 *pEnd = 0; 4493 } 4494 } else { 4495 PetscCall(DMPlexGetDepthStratum(dm, dim, pStart, pEnd)); 4496 } 4497 PetscFunctionReturn(PETSC_SUCCESS); 4498 } 4499 4500 static PetscErrorCode DMGetNeighbors_Plex(DM dm, PetscInt *nranks, const PetscMPIInt *ranks[]) 4501 { 4502 PetscSF sf; 4503 PetscInt niranks, njranks, n; 4504 const PetscMPIInt *iranks, *jranks; 4505 DM_Plex *data = (DM_Plex *)dm->data; 4506 4507 PetscFunctionBegin; 4508 PetscCall(DMGetPointSF(dm, &sf)); 4509 if (!data->neighbors) { 4510 PetscCall(PetscSFSetUp(sf)); 4511 PetscCall(PetscSFGetRootRanks(sf, &njranks, &jranks, NULL, NULL, NULL)); 4512 PetscCall(PetscSFGetLeafRanks(sf, &niranks, &iranks, NULL, NULL)); 4513 PetscCall(PetscMalloc1(njranks + niranks + 1, &data->neighbors)); 4514 PetscCall(PetscArraycpy(data->neighbors + 1, jranks, njranks)); 4515 PetscCall(PetscArraycpy(data->neighbors + njranks + 1, iranks, niranks)); 4516 n = njranks + niranks; 4517 PetscCall(PetscSortRemoveDupsMPIInt(&n, data->neighbors + 1)); 4518 /* The following cast should never fail: can't have more neighbors than PETSC_MPI_INT_MAX */ 4519 PetscCall(PetscMPIIntCast(n, data->neighbors)); 4520 } 4521 if (nranks) *nranks = data->neighbors[0]; 4522 if (ranks) { 4523 if (data->neighbors[0]) *ranks = data->neighbors + 1; 4524 else *ranks = NULL; 4525 } 4526 PetscFunctionReturn(PETSC_SUCCESS); 4527 } 4528 4529 PETSC_INTERN PetscErrorCode DMInterpolateSolution_Plex(DM, DM, Mat, Vec, Vec); 4530 4531 static PetscErrorCode DMInitialize_Plex(DM dm) 4532 { 4533 PetscFunctionBegin; 4534 dm->ops->view = DMView_Plex; 4535 dm->ops->load = DMLoad_Plex; 4536 dm->ops->setfromoptions = DMSetFromOptions_Plex; 4537 dm->ops->clone = DMClone_Plex; 4538 dm->ops->setup = DMSetUp_Plex; 4539 dm->ops->createlocalsection = DMCreateLocalSection_Plex; 4540 dm->ops->createdefaultconstraints = DMCreateDefaultConstraints_Plex; 4541 dm->ops->createglobalvector = DMCreateGlobalVector_Plex; 4542 dm->ops->createlocalvector = DMCreateLocalVector_Plex; 4543 dm->ops->getlocaltoglobalmapping = NULL; 4544 dm->ops->createfieldis = NULL; 4545 dm->ops->createcoordinatedm = DMCreateCoordinateDM_Plex; 4546 dm->ops->createcoordinatefield = DMCreateCoordinateField_Plex; 4547 dm->ops->getcoloring = NULL; 4548 dm->ops->creatematrix = DMCreateMatrix_Plex; 4549 dm->ops->createinterpolation = DMCreateInterpolation_Plex; 4550 dm->ops->createmassmatrix = DMCreateMassMatrix_Plex; 4551 dm->ops->createmassmatrixlumped = DMCreateMassMatrixLumped_Plex; 4552 dm->ops->createinjection = DMCreateInjection_Plex; 4553 dm->ops->refine = DMRefine_Plex; 4554 dm->ops->coarsen = DMCoarsen_Plex; 4555 dm->ops->refinehierarchy = DMRefineHierarchy_Plex; 4556 dm->ops->coarsenhierarchy = DMCoarsenHierarchy_Plex; 4557 dm->ops->extrude = DMExtrude_Plex; 4558 dm->ops->globaltolocalbegin = NULL; 4559 dm->ops->globaltolocalend = NULL; 4560 dm->ops->localtoglobalbegin = NULL; 4561 dm->ops->localtoglobalend = NULL; 4562 dm->ops->destroy = DMDestroy_Plex; 4563 dm->ops->createsubdm = DMCreateSubDM_Plex; 4564 dm->ops->createsuperdm = DMCreateSuperDM_Plex; 4565 dm->ops->getdimpoints = DMGetDimPoints_Plex; 4566 dm->ops->locatepoints = DMLocatePoints_Plex; 4567 dm->ops->projectfunctionlocal = DMProjectFunctionLocal_Plex; 4568 dm->ops->projectfunctionlabellocal = DMProjectFunctionLabelLocal_Plex; 4569 dm->ops->projectfieldlocal = DMProjectFieldLocal_Plex; 4570 dm->ops->projectfieldlabellocal = DMProjectFieldLabelLocal_Plex; 4571 dm->ops->projectbdfieldlabellocal = DMProjectBdFieldLabelLocal_Plex; 4572 dm->ops->computel2diff = DMComputeL2Diff_Plex; 4573 dm->ops->computel2gradientdiff = DMComputeL2GradientDiff_Plex; 4574 dm->ops->computel2fielddiff = DMComputeL2FieldDiff_Plex; 4575 dm->ops->getneighbors = DMGetNeighbors_Plex; 4576 PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMPlexInsertBoundaryValues_C", DMPlexInsertBoundaryValues_Plex)); 4577 PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMPlexInsertTimeDerviativeBoundaryValues_C", DMPlexInsertTimeDerivativeBoundaryValues_Plex)); 4578 PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMSetUpGLVisViewer_C", DMSetUpGLVisViewer_Plex)); 4579 PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMCreateNeumannOverlap_C", DMCreateNeumannOverlap_Plex)); 4580 PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMPlexGetOverlap_C", DMPlexGetOverlap_Plex)); 4581 PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMPlexDistributeGetDefault_C", DMPlexDistributeGetDefault_Plex)); 4582 PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMPlexDistributeSetDefault_C", DMPlexDistributeSetDefault_Plex)); 4583 PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMPlexReorderGetDefault_C", DMPlexReorderGetDefault_Plex)); 4584 PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMPlexReorderSetDefault_C", DMPlexReorderSetDefault_Plex)); 4585 PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMInterpolateSolution_C", DMInterpolateSolution_Plex)); 4586 PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMPlexGetOverlap_C", DMPlexGetOverlap_Plex)); 4587 PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMPlexSetOverlap_C", DMPlexSetOverlap_Plex)); 4588 PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMPlexGetUseCeed_C", DMPlexGetUseCeed_Plex)); 4589 PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMPlexSetUseCeed_C", DMPlexSetUseCeed_Plex)); 4590 PetscFunctionReturn(PETSC_SUCCESS); 4591 } 4592 4593 PETSC_INTERN PetscErrorCode DMClone_Plex(DM dm, DM *newdm) 4594 { 4595 DM_Plex *mesh = (DM_Plex *)dm->data; 4596 PetscSF face_sf; 4597 4598 PetscFunctionBegin; 4599 mesh->refct++; 4600 (*newdm)->data = mesh; 4601 PetscCall(DMPlexGetIsoperiodicFaceSF(dm, &face_sf)); 4602 PetscCall(DMPlexSetIsoperiodicFaceSF(*newdm, face_sf)); 4603 PetscCall(PetscObjectChangeTypeName((PetscObject)*newdm, DMPLEX)); 4604 PetscCall(DMInitialize_Plex(*newdm)); 4605 PetscFunctionReturn(PETSC_SUCCESS); 4606 } 4607 4608 /*MC 4609 DMPLEX = "plex" - A `DM` object that encapsulates an unstructured mesh, or CW Complex, which can be expressed using a Hasse Diagram. 4610 In the local representation, `Vec`s contain all unknowns in the interior and shared boundary. This is 4611 specified by a PetscSection object. Ownership in the global representation is determined by 4612 ownership of the underlying `DMPLEX` points. This is specified by another `PetscSection` object. 4613 4614 Options Database Keys: 4615 + -dm_refine_pre - Refine mesh before distribution 4616 + -dm_refine_uniform_pre - Choose uniform or generator-based refinement 4617 + -dm_refine_volume_limit_pre - Cell volume limit after pre-refinement using generator 4618 . -dm_distribute - Distribute mesh across processes 4619 . -dm_distribute_overlap - Number of cells to overlap for distribution 4620 . -dm_refine - Refine mesh after distribution 4621 . -dm_plex_hash_location - Use grid hashing for point location 4622 . -dm_plex_hash_box_faces <n,m,p> - The number of divisions in each direction of the grid hash 4623 . -dm_plex_partition_balance - Attempt to evenly divide points on partition boundary between processes 4624 . -dm_plex_remesh_bd - Allow changes to the boundary on remeshing 4625 . -dm_plex_max_projection_height - Maximum mesh point height used to project locally 4626 . -dm_plex_regular_refinement - Use special nested projection algorithm for regular refinement 4627 . -dm_plex_check_all - Perform all checks below 4628 . -dm_plex_check_symmetry - Check that the adjacency information in the mesh is symmetric 4629 . -dm_plex_check_skeleton <celltype> - Check that each cell has the correct number of vertices 4630 . -dm_plex_check_faces <celltype> - Check that the faces of each cell give a vertex order this is consistent with what we expect from the cell type 4631 . -dm_plex_check_geometry - Check that cells have positive volume 4632 . -dm_view :mesh.tex:ascii_latex - View the mesh in LaTeX/TikZ 4633 . -dm_plex_view_scale <num> - Scale the TikZ 4634 . -dm_plex_print_fem <num> - View FEM assembly information, such as element vectors and matrices 4635 - -dm_plex_print_fvm <num> - View FVM assembly information, such as flux updates 4636 4637 Level: intermediate 4638 4639 .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMType`, `DMPlexCreate()`, `DMCreate()`, `DMSetType()`, `PetscSection` 4640 M*/ 4641 4642 PETSC_EXTERN PetscErrorCode DMCreate_Plex(DM dm) 4643 { 4644 DM_Plex *mesh; 4645 PetscInt unit; 4646 4647 PetscFunctionBegin; 4648 PetscCall(PetscCitationsRegister(PlexCitation, &Plexcite)); 4649 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 4650 PetscCall(PetscNew(&mesh)); 4651 dm->data = mesh; 4652 4653 mesh->refct = 1; 4654 PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)dm), &mesh->coneSection)); 4655 PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)dm), &mesh->supportSection)); 4656 mesh->refinementUniform = PETSC_TRUE; 4657 mesh->refinementLimit = -1.0; 4658 mesh->distDefault = PETSC_TRUE; 4659 mesh->reorderDefault = DMPLEX_REORDER_DEFAULT_NOTSET; 4660 mesh->distributionName = NULL; 4661 mesh->interpolated = DMPLEX_INTERPOLATED_INVALID; 4662 mesh->interpolatedCollective = DMPLEX_INTERPOLATED_INVALID; 4663 4664 PetscCall(PetscPartitionerCreate(PetscObjectComm((PetscObject)dm), &mesh->partitioner)); 4665 mesh->remeshBd = PETSC_FALSE; 4666 4667 for (unit = 0; unit < NUM_PETSC_UNITS; ++unit) mesh->scale[unit] = 1.0; 4668 4669 mesh->depthState = -1; 4670 mesh->celltypeState = -1; 4671 mesh->printTol = 1.0e-10; 4672 4673 PetscCall(DMInitialize_Plex(dm)); 4674 PetscFunctionReturn(PETSC_SUCCESS); 4675 } 4676 4677 /*@ 4678 DMPlexCreate - Creates a `DMPLEX` object, which encapsulates an unstructured mesh, or CW complex, which can be expressed using a Hasse Diagram. 4679 4680 Collective 4681 4682 Input Parameter: 4683 . comm - The communicator for the `DMPLEX` object 4684 4685 Output Parameter: 4686 . mesh - The `DMPLEX` object 4687 4688 Level: beginner 4689 4690 .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMType`, `DMCreate()`, `DMSetType()` 4691 @*/ 4692 PetscErrorCode DMPlexCreate(MPI_Comm comm, DM *mesh) 4693 { 4694 PetscFunctionBegin; 4695 PetscAssertPointer(mesh, 2); 4696 PetscCall(DMCreate(comm, mesh)); 4697 PetscCall(DMSetType(*mesh, DMPLEX)); 4698 PetscFunctionReturn(PETSC_SUCCESS); 4699 } 4700 4701 /*@C 4702 DMPlexBuildFromCellListParallel - Build distributed `DMPLEX` topology from a list of vertices for each cell (common mesh generator output) 4703 4704 Collective; No Fortran Support 4705 4706 Input Parameters: 4707 + dm - The `DM` 4708 . numCells - The number of cells owned by this process 4709 . numVertices - The number of vertices to be owned by this process, or `PETSC_DECIDE` 4710 . NVertices - The global number of vertices, or `PETSC_DETERMINE` 4711 . numCorners - The number of vertices for each cell 4712 - cells - An array of numCells*numCorners numbers, the global vertex numbers for each cell 4713 4714 Output Parameters: 4715 + vertexSF - (Optional) `PetscSF` describing complete vertex ownership 4716 - verticesAdjSaved - (Optional) vertex adjacency array 4717 4718 Level: advanced 4719 4720 Notes: 4721 Two triangles sharing a face 4722 .vb 4723 4724 2 4725 / | \ 4726 / | \ 4727 / | \ 4728 0 0 | 1 3 4729 \ | / 4730 \ | / 4731 \ | / 4732 1 4733 .ve 4734 would have input 4735 .vb 4736 numCells = 2, numVertices = 4 4737 cells = [0 1 2 1 3 2] 4738 .ve 4739 which would result in the `DMPLEX` 4740 .vb 4741 4742 4 4743 / | \ 4744 / | \ 4745 / | \ 4746 2 0 | 1 5 4747 \ | / 4748 \ | / 4749 \ | / 4750 3 4751 .ve 4752 4753 Vertices are implicitly numbered consecutively 0,...,NVertices. 4754 Each rank owns a chunk of numVertices consecutive vertices. 4755 If numVertices is `PETSC_DECIDE`, PETSc will distribute them as evenly as possible using PetscLayout. 4756 If NVertices is `PETSC_DETERMINE` and numVertices is PETSC_DECIDE, NVertices is computed by PETSc as the maximum vertex index in cells + 1. 4757 If only NVertices is `PETSC_DETERMINE`, it is computed as the sum of numVertices over all ranks. 4758 4759 The cell distribution is arbitrary non-overlapping, independent of the vertex distribution. 4760 4761 .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexBuildFromCellList()`, `DMPlexCreateFromCellListParallelPetsc()`, `DMPlexBuildCoordinatesFromCellListParallel()`, 4762 `PetscSF` 4763 @*/ 4764 PetscErrorCode DMPlexBuildFromCellListParallel(DM dm, PetscInt numCells, PetscInt numVertices, PetscInt NVertices, PetscInt numCorners, const PetscInt cells[], PetscSF *vertexSF, PetscInt **verticesAdjSaved) 4765 { 4766 PetscSF sfPoint; 4767 PetscLayout layout; 4768 PetscInt numVerticesAdj, *verticesAdj, *cones, c, p; 4769 4770 PetscFunctionBegin; 4771 PetscValidLogicalCollectiveInt(dm, NVertices, 4); 4772 PetscCall(PetscLogEventBegin(DMPLEX_BuildFromCellList, dm, 0, 0, 0)); 4773 /* Get/check global number of vertices */ 4774 { 4775 PetscInt NVerticesInCells, i; 4776 const PetscInt len = numCells * numCorners; 4777 4778 /* NVerticesInCells = max(cells) + 1 */ 4779 NVerticesInCells = PETSC_MIN_INT; 4780 for (i = 0; i < len; i++) 4781 if (cells[i] > NVerticesInCells) NVerticesInCells = cells[i]; 4782 ++NVerticesInCells; 4783 PetscCall(MPIU_Allreduce(MPI_IN_PLACE, &NVerticesInCells, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)dm))); 4784 4785 if (numVertices == PETSC_DECIDE && NVertices == PETSC_DECIDE) NVertices = NVerticesInCells; 4786 else 4787 PetscCheck(NVertices == PETSC_DECIDE || NVertices >= NVerticesInCells, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Specified global number of vertices %" PetscInt_FMT " must be greater than or equal to the number of vertices in cells %" PetscInt_FMT, NVertices, NVerticesInCells); 4788 } 4789 /* Count locally unique vertices */ 4790 { 4791 PetscHSetI vhash; 4792 PetscInt off = 0; 4793 4794 PetscCall(PetscHSetICreate(&vhash)); 4795 for (c = 0; c < numCells; ++c) { 4796 for (p = 0; p < numCorners; ++p) PetscCall(PetscHSetIAdd(vhash, cells[c * numCorners + p])); 4797 } 4798 PetscCall(PetscHSetIGetSize(vhash, &numVerticesAdj)); 4799 if (!verticesAdjSaved) PetscCall(PetscMalloc1(numVerticesAdj, &verticesAdj)); 4800 else verticesAdj = *verticesAdjSaved; 4801 PetscCall(PetscHSetIGetElems(vhash, &off, verticesAdj)); 4802 PetscCall(PetscHSetIDestroy(&vhash)); 4803 PetscCheck(off == numVerticesAdj, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid number of local vertices %" PetscInt_FMT " should be %" PetscInt_FMT, off, numVerticesAdj); 4804 } 4805 PetscCall(PetscSortInt(numVerticesAdj, verticesAdj)); 4806 /* Create cones */ 4807 PetscCall(DMPlexSetChart(dm, 0, numCells + numVerticesAdj)); 4808 for (c = 0; c < numCells; ++c) PetscCall(DMPlexSetConeSize(dm, c, numCorners)); 4809 PetscCall(DMSetUp(dm)); 4810 PetscCall(DMPlexGetCones(dm, &cones)); 4811 for (c = 0; c < numCells; ++c) { 4812 for (p = 0; p < numCorners; ++p) { 4813 const PetscInt gv = cells[c * numCorners + p]; 4814 PetscInt lv; 4815 4816 /* Positions within verticesAdj form 0-based local vertex numbering; 4817 we need to shift it by numCells to get correct DAG points (cells go first) */ 4818 PetscCall(PetscFindInt(gv, numVerticesAdj, verticesAdj, &lv)); 4819 PetscCheck(lv >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Could not find global vertex %" PetscInt_FMT " in local connectivity", gv); 4820 cones[c * numCorners + p] = lv + numCells; 4821 } 4822 } 4823 /* Build point sf */ 4824 PetscCall(PetscLayoutCreate(PetscObjectComm((PetscObject)dm), &layout)); 4825 PetscCall(PetscLayoutSetSize(layout, NVertices)); 4826 PetscCall(PetscLayoutSetLocalSize(layout, numVertices)); 4827 PetscCall(PetscLayoutSetBlockSize(layout, 1)); 4828 PetscCall(PetscSFCreateByMatchingIndices(layout, numVerticesAdj, verticesAdj, NULL, numCells, numVerticesAdj, verticesAdj, NULL, numCells, vertexSF, &sfPoint)); 4829 PetscCall(PetscLayoutDestroy(&layout)); 4830 if (!verticesAdjSaved) PetscCall(PetscFree(verticesAdj)); 4831 PetscCall(PetscObjectSetName((PetscObject)sfPoint, "point SF")); 4832 if (dm->sf) { 4833 const char *prefix; 4834 4835 PetscCall(PetscObjectGetOptionsPrefix((PetscObject)dm->sf, &prefix)); 4836 PetscCall(PetscObjectSetOptionsPrefix((PetscObject)sfPoint, prefix)); 4837 } 4838 PetscCall(DMSetPointSF(dm, sfPoint)); 4839 PetscCall(PetscSFDestroy(&sfPoint)); 4840 if (vertexSF) PetscCall(PetscObjectSetName((PetscObject)(*vertexSF), "Vertex Ownership SF")); 4841 /* Fill in the rest of the topology structure */ 4842 PetscCall(DMPlexSymmetrize(dm)); 4843 PetscCall(DMPlexStratify(dm)); 4844 PetscCall(PetscLogEventEnd(DMPLEX_BuildFromCellList, dm, 0, 0, 0)); 4845 PetscFunctionReturn(PETSC_SUCCESS); 4846 } 4847 4848 /*@C 4849 DMPlexBuildCoordinatesFromCellListParallel - Build `DM` coordinates from a list of coordinates for each owned vertex (common mesh generator output) 4850 4851 Collective; No Fortran Support 4852 4853 Input Parameters: 4854 + dm - The `DM` 4855 . spaceDim - The spatial dimension used for coordinates 4856 . sfVert - `PetscSF` describing complete vertex ownership 4857 - vertexCoords - An array of numVertices*spaceDim numbers, the coordinates of each vertex 4858 4859 Level: advanced 4860 4861 .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexBuildCoordinatesFromCellList()`, `DMPlexCreateFromCellListParallelPetsc()`, `DMPlexBuildFromCellListParallel()` 4862 @*/ 4863 PetscErrorCode DMPlexBuildCoordinatesFromCellListParallel(DM dm, PetscInt spaceDim, PetscSF sfVert, const PetscReal vertexCoords[]) 4864 { 4865 PetscSection coordSection; 4866 Vec coordinates; 4867 PetscScalar *coords; 4868 PetscInt numVertices, numVerticesAdj, coordSize, v, vStart, vEnd; 4869 4870 PetscFunctionBegin; 4871 PetscCall(PetscLogEventBegin(DMPLEX_BuildCoordinatesFromCellList, dm, 0, 0, 0)); 4872 PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd)); 4873 PetscCheck(vStart >= 0 && vEnd >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "DM is not set up properly. DMPlexBuildFromCellList() should be called first."); 4874 PetscCall(DMSetCoordinateDim(dm, spaceDim)); 4875 PetscCall(PetscSFGetGraph(sfVert, &numVertices, &numVerticesAdj, NULL, NULL)); 4876 PetscCheck(vEnd - vStart == numVerticesAdj, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Supplied sfVert has wrong number of leaves = %" PetscInt_FMT " != %" PetscInt_FMT " = vEnd - vStart", numVerticesAdj, vEnd - vStart); 4877 PetscCall(DMGetCoordinateSection(dm, &coordSection)); 4878 PetscCall(PetscSectionSetNumFields(coordSection, 1)); 4879 PetscCall(PetscSectionSetFieldComponents(coordSection, 0, spaceDim)); 4880 PetscCall(PetscSectionSetChart(coordSection, vStart, vEnd)); 4881 for (v = vStart; v < vEnd; ++v) { 4882 PetscCall(PetscSectionSetDof(coordSection, v, spaceDim)); 4883 PetscCall(PetscSectionSetFieldDof(coordSection, v, 0, spaceDim)); 4884 } 4885 PetscCall(PetscSectionSetUp(coordSection)); 4886 PetscCall(PetscSectionGetStorageSize(coordSection, &coordSize)); 4887 PetscCall(VecCreate(PetscObjectComm((PetscObject)dm), &coordinates)); 4888 PetscCall(VecSetBlockSize(coordinates, spaceDim)); 4889 PetscCall(PetscObjectSetName((PetscObject)coordinates, "coordinates")); 4890 PetscCall(VecSetSizes(coordinates, coordSize, PETSC_DETERMINE)); 4891 PetscCall(VecSetType(coordinates, VECSTANDARD)); 4892 PetscCall(VecGetArray(coordinates, &coords)); 4893 { 4894 MPI_Datatype coordtype; 4895 4896 /* Need a temp buffer for coords if we have complex/single */ 4897 PetscCallMPI(MPI_Type_contiguous(spaceDim, MPIU_SCALAR, &coordtype)); 4898 PetscCallMPI(MPI_Type_commit(&coordtype)); 4899 #if defined(PETSC_USE_COMPLEX) 4900 { 4901 PetscScalar *svertexCoords; 4902 PetscInt i; 4903 PetscCall(PetscMalloc1(numVertices * spaceDim, &svertexCoords)); 4904 for (i = 0; i < numVertices * spaceDim; i++) svertexCoords[i] = vertexCoords[i]; 4905 PetscCall(PetscSFBcastBegin(sfVert, coordtype, svertexCoords, coords, MPI_REPLACE)); 4906 PetscCall(PetscSFBcastEnd(sfVert, coordtype, svertexCoords, coords, MPI_REPLACE)); 4907 PetscCall(PetscFree(svertexCoords)); 4908 } 4909 #else 4910 PetscCall(PetscSFBcastBegin(sfVert, coordtype, vertexCoords, coords, MPI_REPLACE)); 4911 PetscCall(PetscSFBcastEnd(sfVert, coordtype, vertexCoords, coords, MPI_REPLACE)); 4912 #endif 4913 PetscCallMPI(MPI_Type_free(&coordtype)); 4914 } 4915 PetscCall(VecRestoreArray(coordinates, &coords)); 4916 PetscCall(DMSetCoordinatesLocal(dm, coordinates)); 4917 PetscCall(VecDestroy(&coordinates)); 4918 PetscCall(PetscLogEventEnd(DMPLEX_BuildCoordinatesFromCellList, dm, 0, 0, 0)); 4919 PetscFunctionReturn(PETSC_SUCCESS); 4920 } 4921 4922 /*@ 4923 DMPlexCreateFromCellListParallelPetsc - Create distributed `DMPLEX` from a list of vertices for each cell (common mesh generator output) 4924 4925 Collective 4926 4927 Input Parameters: 4928 + comm - The communicator 4929 . dim - The topological dimension of the mesh 4930 . numCells - The number of cells owned by this process 4931 . numVertices - The number of vertices owned by this process, or `PETSC_DECIDE` 4932 . NVertices - The global number of vertices, or `PETSC_DECIDE` 4933 . numCorners - The number of vertices for each cell 4934 . interpolate - Flag indicating that intermediate mesh entities (faces, edges) should be created automatically 4935 . cells - An array of numCells*numCorners numbers, the global vertex numbers for each cell 4936 . spaceDim - The spatial dimension used for coordinates 4937 - vertexCoords - An array of numVertices*spaceDim numbers, the coordinates of each vertex 4938 4939 Output Parameters: 4940 + dm - The `DM` 4941 . vertexSF - (Optional) `PetscSF` describing complete vertex ownership 4942 - verticesAdj - (Optional) vertex adjacency array 4943 4944 Level: intermediate 4945 4946 Notes: 4947 This function is just a convenient sequence of `DMCreate()`, `DMSetType()`, `DMSetDimension()`, 4948 `DMPlexBuildFromCellListParallel()`, `DMPlexInterpolate()`, `DMPlexBuildCoordinatesFromCellListParallel()` 4949 4950 See `DMPlexBuildFromCellListParallel()` for an example and details about the topology-related parameters. 4951 4952 See `DMPlexBuildCoordinatesFromCellListParallel()` for details about the geometry-related parameters. 4953 4954 .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexCreateFromCellListPetsc()`, `DMPlexBuildFromCellListParallel()`, `DMPlexBuildCoordinatesFromCellListParallel()`, `DMPlexCreateFromDAG()`, `DMPlexCreate()` 4955 @*/ 4956 PetscErrorCode DMPlexCreateFromCellListParallelPetsc(MPI_Comm comm, PetscInt dim, PetscInt numCells, PetscInt numVertices, PetscInt NVertices, PetscInt numCorners, PetscBool interpolate, const PetscInt cells[], PetscInt spaceDim, const PetscReal vertexCoords[], PetscSF *vertexSF, PetscInt **verticesAdj, DM *dm) 4957 { 4958 PetscSF sfVert; 4959 4960 PetscFunctionBegin; 4961 PetscCall(DMCreate(comm, dm)); 4962 PetscCall(DMSetType(*dm, DMPLEX)); 4963 PetscValidLogicalCollectiveInt(*dm, dim, 2); 4964 PetscValidLogicalCollectiveInt(*dm, spaceDim, 9); 4965 PetscCall(DMSetDimension(*dm, dim)); 4966 PetscCall(DMPlexBuildFromCellListParallel(*dm, numCells, numVertices, NVertices, numCorners, cells, &sfVert, verticesAdj)); 4967 if (interpolate) { 4968 DM idm; 4969 4970 PetscCall(DMPlexInterpolate(*dm, &idm)); 4971 PetscCall(DMDestroy(dm)); 4972 *dm = idm; 4973 } 4974 PetscCall(DMPlexBuildCoordinatesFromCellListParallel(*dm, spaceDim, sfVert, vertexCoords)); 4975 if (vertexSF) *vertexSF = sfVert; 4976 else PetscCall(PetscSFDestroy(&sfVert)); 4977 PetscFunctionReturn(PETSC_SUCCESS); 4978 } 4979 4980 /*@C 4981 DMPlexBuildFromCellList - Build `DMPLEX` topology from a list of vertices for each cell (common mesh generator output) 4982 4983 Collective; No Fortran Support 4984 4985 Input Parameters: 4986 + dm - The `DM` 4987 . numCells - The number of cells owned by this process 4988 . numVertices - The number of vertices owned by this process, or `PETSC_DETERMINE` 4989 . numCorners - The number of vertices for each cell 4990 - cells - An array of numCells*numCorners numbers, the global vertex numbers for each cell 4991 4992 Level: advanced 4993 4994 Notes: 4995 Two triangles sharing a face 4996 .vb 4997 4998 2 4999 / | \ 5000 / | \ 5001 / | \ 5002 0 0 | 1 3 5003 \ | / 5004 \ | / 5005 \ | / 5006 1 5007 .ve 5008 would have input 5009 .vb 5010 numCells = 2, numVertices = 4 5011 cells = [0 1 2 1 3 2] 5012 .ve 5013 which would result in the `DMPLEX` 5014 .vb 5015 5016 4 5017 / | \ 5018 / | \ 5019 / | \ 5020 2 0 | 1 5 5021 \ | / 5022 \ | / 5023 \ | / 5024 3 5025 .ve 5026 5027 If numVertices is `PETSC_DETERMINE`, it is computed by PETSc as the maximum vertex index in cells + 1. 5028 5029 .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexBuildFromCellListParallel()`, `DMPlexBuildCoordinatesFromCellList()`, `DMPlexCreateFromCellListPetsc()` 5030 @*/ 5031 PetscErrorCode DMPlexBuildFromCellList(DM dm, PetscInt numCells, PetscInt numVertices, PetscInt numCorners, const PetscInt cells[]) 5032 { 5033 PetscInt *cones, c, p, dim; 5034 5035 PetscFunctionBegin; 5036 PetscCall(PetscLogEventBegin(DMPLEX_BuildFromCellList, dm, 0, 0, 0)); 5037 PetscCall(DMGetDimension(dm, &dim)); 5038 /* Get/check global number of vertices */ 5039 { 5040 PetscInt NVerticesInCells, i; 5041 const PetscInt len = numCells * numCorners; 5042 5043 /* NVerticesInCells = max(cells) + 1 */ 5044 NVerticesInCells = PETSC_MIN_INT; 5045 for (i = 0; i < len; i++) 5046 if (cells[i] > NVerticesInCells) NVerticesInCells = cells[i]; 5047 ++NVerticesInCells; 5048 5049 if (numVertices == PETSC_DECIDE) numVertices = NVerticesInCells; 5050 else 5051 PetscCheck(numVertices >= NVerticesInCells, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Specified number of vertices %" PetscInt_FMT " must be greater than or equal to the number of vertices in cells %" PetscInt_FMT, numVertices, NVerticesInCells); 5052 } 5053 PetscCall(DMPlexSetChart(dm, 0, numCells + numVertices)); 5054 for (c = 0; c < numCells; ++c) PetscCall(DMPlexSetConeSize(dm, c, numCorners)); 5055 PetscCall(DMSetUp(dm)); 5056 PetscCall(DMPlexGetCones(dm, &cones)); 5057 for (c = 0; c < numCells; ++c) { 5058 for (p = 0; p < numCorners; ++p) cones[c * numCorners + p] = cells[c * numCorners + p] + numCells; 5059 } 5060 PetscCall(DMPlexSymmetrize(dm)); 5061 PetscCall(DMPlexStratify(dm)); 5062 PetscCall(PetscLogEventEnd(DMPLEX_BuildFromCellList, dm, 0, 0, 0)); 5063 PetscFunctionReturn(PETSC_SUCCESS); 5064 } 5065 5066 /*@C 5067 DMPlexBuildCoordinatesFromCellList - Build `DM` coordinates from a list of coordinates for each owned vertex (common mesh generator output) 5068 5069 Collective; No Fortran Support 5070 5071 Input Parameters: 5072 + dm - The `DM` 5073 . spaceDim - The spatial dimension used for coordinates 5074 - vertexCoords - An array of numVertices*spaceDim numbers, the coordinates of each vertex 5075 5076 Level: advanced 5077 5078 .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexBuildCoordinatesFromCellListParallel()`, `DMPlexCreateFromCellListPetsc()`, `DMPlexBuildFromCellList()` 5079 @*/ 5080 PetscErrorCode DMPlexBuildCoordinatesFromCellList(DM dm, PetscInt spaceDim, const PetscReal vertexCoords[]) 5081 { 5082 PetscSection coordSection; 5083 Vec coordinates; 5084 DM cdm; 5085 PetscScalar *coords; 5086 PetscInt v, vStart, vEnd, d; 5087 5088 PetscFunctionBegin; 5089 PetscCall(PetscLogEventBegin(DMPLEX_BuildCoordinatesFromCellList, dm, 0, 0, 0)); 5090 PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd)); 5091 PetscCheck(vStart >= 0 && vEnd >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "DM is not set up properly. DMPlexBuildFromCellList() should be called first."); 5092 PetscCall(DMSetCoordinateDim(dm, spaceDim)); 5093 PetscCall(DMGetCoordinateSection(dm, &coordSection)); 5094 PetscCall(PetscSectionSetNumFields(coordSection, 1)); 5095 PetscCall(PetscSectionSetFieldComponents(coordSection, 0, spaceDim)); 5096 PetscCall(PetscSectionSetChart(coordSection, vStart, vEnd)); 5097 for (v = vStart; v < vEnd; ++v) { 5098 PetscCall(PetscSectionSetDof(coordSection, v, spaceDim)); 5099 PetscCall(PetscSectionSetFieldDof(coordSection, v, 0, spaceDim)); 5100 } 5101 PetscCall(PetscSectionSetUp(coordSection)); 5102 5103 PetscCall(DMGetCoordinateDM(dm, &cdm)); 5104 PetscCall(DMCreateLocalVector(cdm, &coordinates)); 5105 PetscCall(VecSetBlockSize(coordinates, spaceDim)); 5106 PetscCall(PetscObjectSetName((PetscObject)coordinates, "coordinates")); 5107 PetscCall(VecGetArrayWrite(coordinates, &coords)); 5108 for (v = 0; v < vEnd - vStart; ++v) { 5109 for (d = 0; d < spaceDim; ++d) coords[v * spaceDim + d] = vertexCoords[v * spaceDim + d]; 5110 } 5111 PetscCall(VecRestoreArrayWrite(coordinates, &coords)); 5112 PetscCall(DMSetCoordinatesLocal(dm, coordinates)); 5113 PetscCall(VecDestroy(&coordinates)); 5114 PetscCall(PetscLogEventEnd(DMPLEX_BuildCoordinatesFromCellList, dm, 0, 0, 0)); 5115 PetscFunctionReturn(PETSC_SUCCESS); 5116 } 5117 5118 /*@ 5119 DMPlexCreateFromCellListPetsc - Create `DMPLEX` from a list of vertices for each cell (common mesh generator output), but only process 0 takes in the input 5120 5121 Collective 5122 5123 Input Parameters: 5124 + comm - The communicator 5125 . dim - The topological dimension of the mesh 5126 . numCells - The number of cells, only on process 0 5127 . numVertices - The number of vertices owned by this process, or `PETSC_DECIDE`, only on process 0 5128 . numCorners - The number of vertices for each cell, only on process 0 5129 . interpolate - Flag indicating that intermediate mesh entities (faces, edges) should be created automatically 5130 . cells - An array of numCells*numCorners numbers, the vertices for each cell, only on process 0 5131 . spaceDim - The spatial dimension used for coordinates 5132 - vertexCoords - An array of numVertices*spaceDim numbers, the coordinates of each vertex, only on process 0 5133 5134 Output Parameter: 5135 . dm - The `DM`, which only has points on process 0 5136 5137 Level: intermediate 5138 5139 Notes: 5140 This function is just a convenient sequence of `DMCreate()`, `DMSetType()`, `DMSetDimension()`, `DMPlexBuildFromCellList()`, 5141 `DMPlexInterpolate()`, `DMPlexBuildCoordinatesFromCellList()` 5142 5143 See `DMPlexBuildFromCellList()` for an example and details about the topology-related parameters. 5144 See `DMPlexBuildCoordinatesFromCellList()` for details about the geometry-related parameters. 5145 See `DMPlexCreateFromCellListParallelPetsc()` for parallel input 5146 5147 .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexCreateFromCellListParallelPetsc()`, `DMPlexBuildFromCellList()`, `DMPlexBuildCoordinatesFromCellList()`, `DMPlexCreateFromDAG()`, `DMPlexCreate()` 5148 @*/ 5149 PetscErrorCode DMPlexCreateFromCellListPetsc(MPI_Comm comm, PetscInt dim, PetscInt numCells, PetscInt numVertices, PetscInt numCorners, PetscBool interpolate, const PetscInt cells[], PetscInt spaceDim, const PetscReal vertexCoords[], DM *dm) 5150 { 5151 PetscMPIInt rank; 5152 5153 PetscFunctionBegin; 5154 PetscCheck(dim, comm, PETSC_ERR_ARG_OUTOFRANGE, "This is not appropriate for 0-dimensional meshes. Consider either creating the DM using DMPlexCreateFromDAG(), by hand, or using DMSwarm."); 5155 PetscCallMPI(MPI_Comm_rank(comm, &rank)); 5156 PetscCall(DMCreate(comm, dm)); 5157 PetscCall(DMSetType(*dm, DMPLEX)); 5158 PetscCall(DMSetDimension(*dm, dim)); 5159 if (rank == 0) PetscCall(DMPlexBuildFromCellList(*dm, numCells, numVertices, numCorners, cells)); 5160 else PetscCall(DMPlexBuildFromCellList(*dm, 0, 0, 0, NULL)); 5161 if (interpolate) { 5162 DM idm; 5163 5164 PetscCall(DMPlexInterpolate(*dm, &idm)); 5165 PetscCall(DMDestroy(dm)); 5166 *dm = idm; 5167 } 5168 if (rank == 0) PetscCall(DMPlexBuildCoordinatesFromCellList(*dm, spaceDim, vertexCoords)); 5169 else PetscCall(DMPlexBuildCoordinatesFromCellList(*dm, spaceDim, NULL)); 5170 PetscFunctionReturn(PETSC_SUCCESS); 5171 } 5172 5173 /*@ 5174 DMPlexCreateFromDAG - This takes as input the adjacency-list representation of the Directed Acyclic Graph (Hasse Diagram) encoding a mesh, and produces a `DM` 5175 5176 Input Parameters: 5177 + dm - The empty `DM` object, usually from `DMCreate()` and `DMSetDimension()` 5178 . depth - The depth of the DAG 5179 . numPoints - Array of size depth + 1 containing the number of points at each `depth` 5180 . coneSize - The cone size of each point 5181 . cones - The concatenation of the cone points for each point, the cone list must be oriented correctly for each point 5182 . coneOrientations - The orientation of each cone point 5183 - vertexCoords - An array of `numPoints`[0]*spacedim numbers representing the coordinates of each vertex, with spacedim the value set via `DMSetCoordinateDim()` 5184 5185 Output Parameter: 5186 . dm - The `DM` 5187 5188 Level: advanced 5189 5190 Note: 5191 Two triangles sharing a face would have input 5192 .vb 5193 depth = 1, numPoints = [4 2], coneSize = [3 3 0 0 0 0] 5194 cones = [2 3 4 3 5 4], coneOrientations = [0 0 0 0 0 0] 5195 vertexCoords = [-1.0 0.0 0.0 -1.0 0.0 1.0 1.0 0.0] 5196 .ve 5197 which would result in the DMPlex 5198 .vb 5199 4 5200 / | \ 5201 / | \ 5202 / | \ 5203 2 0 | 1 5 5204 \ | / 5205 \ | / 5206 \ | / 5207 3 5208 .ve 5209 Notice that all points are numbered consecutively, unlike `DMPlexCreateFromCellListPetsc()` 5210 5211 .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexCreateFromCellListPetsc()`, `DMPlexCreate()` 5212 @*/ 5213 PetscErrorCode DMPlexCreateFromDAG(DM dm, PetscInt depth, const PetscInt numPoints[], const PetscInt coneSize[], const PetscInt cones[], const PetscInt coneOrientations[], const PetscScalar vertexCoords[]) 5214 { 5215 Vec coordinates; 5216 PetscSection coordSection; 5217 PetscScalar *coords; 5218 PetscInt coordSize, firstVertex = -1, pStart = 0, pEnd = 0, p, v, dim, dimEmbed, d, off; 5219 5220 PetscFunctionBegin; 5221 PetscCall(DMGetDimension(dm, &dim)); 5222 PetscCall(DMGetCoordinateDim(dm, &dimEmbed)); 5223 PetscCheck(dimEmbed >= dim, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Embedding dimension %" PetscInt_FMT " cannot be less than intrinsic dimension %" PetscInt_FMT, dimEmbed, dim); 5224 for (d = 0; d <= depth; ++d) pEnd += numPoints[d]; 5225 PetscCall(DMPlexSetChart(dm, pStart, pEnd)); 5226 for (p = pStart; p < pEnd; ++p) { 5227 PetscCall(DMPlexSetConeSize(dm, p, coneSize[p - pStart])); 5228 if (firstVertex < 0 && !coneSize[p - pStart]) firstVertex = p - pStart; 5229 } 5230 PetscCheck(firstVertex >= 0 || !numPoints[0], PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Expected %" PetscInt_FMT " vertices but could not find any", numPoints[0]); 5231 PetscCall(DMSetUp(dm)); /* Allocate space for cones */ 5232 for (p = pStart, off = 0; p < pEnd; off += coneSize[p - pStart], ++p) { 5233 PetscCall(DMPlexSetCone(dm, p, &cones[off])); 5234 PetscCall(DMPlexSetConeOrientation(dm, p, &coneOrientations[off])); 5235 } 5236 PetscCall(DMPlexSymmetrize(dm)); 5237 PetscCall(DMPlexStratify(dm)); 5238 /* Build coordinates */ 5239 PetscCall(DMGetCoordinateSection(dm, &coordSection)); 5240 PetscCall(PetscSectionSetNumFields(coordSection, 1)); 5241 PetscCall(PetscSectionSetFieldComponents(coordSection, 0, dimEmbed)); 5242 PetscCall(PetscSectionSetChart(coordSection, firstVertex, firstVertex + numPoints[0])); 5243 for (v = firstVertex; v < firstVertex + numPoints[0]; ++v) { 5244 PetscCall(PetscSectionSetDof(coordSection, v, dimEmbed)); 5245 PetscCall(PetscSectionSetFieldDof(coordSection, v, 0, dimEmbed)); 5246 } 5247 PetscCall(PetscSectionSetUp(coordSection)); 5248 PetscCall(PetscSectionGetStorageSize(coordSection, &coordSize)); 5249 PetscCall(VecCreate(PETSC_COMM_SELF, &coordinates)); 5250 PetscCall(PetscObjectSetName((PetscObject)coordinates, "coordinates")); 5251 PetscCall(VecSetSizes(coordinates, coordSize, PETSC_DETERMINE)); 5252 PetscCall(VecSetBlockSize(coordinates, dimEmbed)); 5253 PetscCall(VecSetType(coordinates, VECSTANDARD)); 5254 if (vertexCoords) { 5255 PetscCall(VecGetArray(coordinates, &coords)); 5256 for (v = 0; v < numPoints[0]; ++v) { 5257 PetscInt off; 5258 5259 PetscCall(PetscSectionGetOffset(coordSection, v + firstVertex, &off)); 5260 for (d = 0; d < dimEmbed; ++d) coords[off + d] = vertexCoords[v * dimEmbed + d]; 5261 } 5262 } 5263 PetscCall(VecRestoreArray(coordinates, &coords)); 5264 PetscCall(DMSetCoordinatesLocal(dm, coordinates)); 5265 PetscCall(VecDestroy(&coordinates)); 5266 PetscFunctionReturn(PETSC_SUCCESS); 5267 } 5268 5269 /* 5270 DMPlexCreateCellVertexFromFile - Create a `DMPLEX` mesh from a simple cell-vertex file. 5271 5272 Collective 5273 5274 + comm - The MPI communicator 5275 . filename - Name of the .dat file 5276 - interpolate - Create faces and edges in the mesh 5277 5278 Output Parameter: 5279 . dm - The `DM` object representing the mesh 5280 5281 Level: beginner 5282 5283 Note: 5284 The format is the simplest possible: 5285 .vb 5286 Ne 5287 v0 v1 ... vk 5288 Nv 5289 x y z marker 5290 .ve 5291 5292 Developer Note: 5293 Should use a `PetscViewer` not a filename 5294 5295 .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexCreateFromFile()`, `DMPlexCreateMedFromFile()`, `DMPlexCreateGmsh()`, `DMPlexCreate()` 5296 */ 5297 static PetscErrorCode DMPlexCreateCellVertexFromFile(MPI_Comm comm, const char filename[], PetscBool interpolate, DM *dm) 5298 { 5299 DMLabel marker; 5300 PetscViewer viewer; 5301 Vec coordinates; 5302 PetscSection coordSection; 5303 PetscScalar *coords; 5304 char line[PETSC_MAX_PATH_LEN]; 5305 PetscInt dim = 3, cdim = 3, coordSize, v, c, d; 5306 PetscMPIInt rank; 5307 int snum, Nv, Nc, Ncn, Nl; 5308 5309 PetscFunctionBegin; 5310 PetscCallMPI(MPI_Comm_rank(comm, &rank)); 5311 PetscCall(PetscViewerCreate(comm, &viewer)); 5312 PetscCall(PetscViewerSetType(viewer, PETSCVIEWERASCII)); 5313 PetscCall(PetscViewerFileSetMode(viewer, FILE_MODE_READ)); 5314 PetscCall(PetscViewerFileSetName(viewer, filename)); 5315 if (rank == 0) { 5316 PetscCall(PetscViewerRead(viewer, line, 4, NULL, PETSC_STRING)); 5317 snum = sscanf(line, "%d %d %d %d", &Nc, &Nv, &Ncn, &Nl); 5318 PetscCheck(snum == 4, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse cell-vertex file: %s", line); 5319 } else { 5320 Nc = Nv = Ncn = Nl = 0; 5321 } 5322 PetscCall(DMCreate(comm, dm)); 5323 PetscCall(DMSetType(*dm, DMPLEX)); 5324 PetscCall(DMPlexSetChart(*dm, 0, Nc + Nv)); 5325 PetscCall(DMSetDimension(*dm, dim)); 5326 PetscCall(DMSetCoordinateDim(*dm, cdim)); 5327 /* Read topology */ 5328 if (rank == 0) { 5329 char format[PETSC_MAX_PATH_LEN]; 5330 PetscInt cone[8]; 5331 int vbuf[8], v; 5332 5333 for (c = 0; c < Ncn; ++c) { 5334 format[c * 3 + 0] = '%'; 5335 format[c * 3 + 1] = 'd'; 5336 format[c * 3 + 2] = ' '; 5337 } 5338 format[Ncn * 3 - 1] = '\0'; 5339 for (c = 0; c < Nc; ++c) PetscCall(DMPlexSetConeSize(*dm, c, Ncn)); 5340 PetscCall(DMSetUp(*dm)); 5341 for (c = 0; c < Nc; ++c) { 5342 PetscCall(PetscViewerRead(viewer, line, Ncn, NULL, PETSC_STRING)); 5343 switch (Ncn) { 5344 case 2: 5345 snum = sscanf(line, format, &vbuf[0], &vbuf[1]); 5346 break; 5347 case 3: 5348 snum = sscanf(line, format, &vbuf[0], &vbuf[1], &vbuf[2]); 5349 break; 5350 case 4: 5351 snum = sscanf(line, format, &vbuf[0], &vbuf[1], &vbuf[2], &vbuf[3]); 5352 break; 5353 case 6: 5354 snum = sscanf(line, format, &vbuf[0], &vbuf[1], &vbuf[2], &vbuf[3], &vbuf[4], &vbuf[5]); 5355 break; 5356 case 8: 5357 snum = sscanf(line, format, &vbuf[0], &vbuf[1], &vbuf[2], &vbuf[3], &vbuf[4], &vbuf[5], &vbuf[6], &vbuf[7]); 5358 break; 5359 default: 5360 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "No cell shape with %d vertices", Ncn); 5361 } 5362 PetscCheck(snum == Ncn, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse cell-vertex file: %s", line); 5363 for (v = 0; v < Ncn; ++v) cone[v] = vbuf[v] + Nc; 5364 /* Hexahedra are inverted */ 5365 if (Ncn == 8) { 5366 PetscInt tmp = cone[1]; 5367 cone[1] = cone[3]; 5368 cone[3] = tmp; 5369 } 5370 PetscCall(DMPlexSetCone(*dm, c, cone)); 5371 } 5372 } 5373 PetscCall(DMPlexSymmetrize(*dm)); 5374 PetscCall(DMPlexStratify(*dm)); 5375 /* Read coordinates */ 5376 PetscCall(DMGetCoordinateSection(*dm, &coordSection)); 5377 PetscCall(PetscSectionSetNumFields(coordSection, 1)); 5378 PetscCall(PetscSectionSetFieldComponents(coordSection, 0, cdim)); 5379 PetscCall(PetscSectionSetChart(coordSection, Nc, Nc + Nv)); 5380 for (v = Nc; v < Nc + Nv; ++v) { 5381 PetscCall(PetscSectionSetDof(coordSection, v, cdim)); 5382 PetscCall(PetscSectionSetFieldDof(coordSection, v, 0, cdim)); 5383 } 5384 PetscCall(PetscSectionSetUp(coordSection)); 5385 PetscCall(PetscSectionGetStorageSize(coordSection, &coordSize)); 5386 PetscCall(VecCreate(PETSC_COMM_SELF, &coordinates)); 5387 PetscCall(PetscObjectSetName((PetscObject)coordinates, "coordinates")); 5388 PetscCall(VecSetSizes(coordinates, coordSize, PETSC_DETERMINE)); 5389 PetscCall(VecSetBlockSize(coordinates, cdim)); 5390 PetscCall(VecSetType(coordinates, VECSTANDARD)); 5391 PetscCall(VecGetArray(coordinates, &coords)); 5392 if (rank == 0) { 5393 char format[PETSC_MAX_PATH_LEN]; 5394 double x[3]; 5395 int l, val[3]; 5396 5397 if (Nl) { 5398 for (l = 0; l < Nl; ++l) { 5399 format[l * 3 + 0] = '%'; 5400 format[l * 3 + 1] = 'd'; 5401 format[l * 3 + 2] = ' '; 5402 } 5403 format[Nl * 3 - 1] = '\0'; 5404 PetscCall(DMCreateLabel(*dm, "marker")); 5405 PetscCall(DMGetLabel(*dm, "marker", &marker)); 5406 } 5407 for (v = 0; v < Nv; ++v) { 5408 PetscCall(PetscViewerRead(viewer, line, 3 + Nl, NULL, PETSC_STRING)); 5409 snum = sscanf(line, "%lg %lg %lg", &x[0], &x[1], &x[2]); 5410 PetscCheck(snum == 3, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse cell-vertex file: %s", line); 5411 switch (Nl) { 5412 case 0: 5413 snum = 0; 5414 break; 5415 case 1: 5416 snum = sscanf(line, format, &val[0]); 5417 break; 5418 case 2: 5419 snum = sscanf(line, format, &val[0], &val[1]); 5420 break; 5421 case 3: 5422 snum = sscanf(line, format, &val[0], &val[1], &val[2]); 5423 break; 5424 default: 5425 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Request support for %d labels", Nl); 5426 } 5427 PetscCheck(snum == Nl, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse cell-vertex file: %s", line); 5428 for (d = 0; d < cdim; ++d) coords[v * cdim + d] = x[d]; 5429 for (l = 0; l < Nl; ++l) PetscCall(DMLabelSetValue(marker, v + Nc, val[l])); 5430 } 5431 } 5432 PetscCall(VecRestoreArray(coordinates, &coords)); 5433 PetscCall(DMSetCoordinatesLocal(*dm, coordinates)); 5434 PetscCall(VecDestroy(&coordinates)); 5435 PetscCall(PetscViewerDestroy(&viewer)); 5436 if (interpolate) { 5437 DM idm; 5438 DMLabel bdlabel; 5439 5440 PetscCall(DMPlexInterpolate(*dm, &idm)); 5441 PetscCall(DMDestroy(dm)); 5442 *dm = idm; 5443 5444 if (!Nl) { 5445 PetscCall(DMCreateLabel(*dm, "marker")); 5446 PetscCall(DMGetLabel(*dm, "marker", &bdlabel)); 5447 PetscCall(DMPlexMarkBoundaryFaces(*dm, PETSC_DETERMINE, bdlabel)); 5448 PetscCall(DMPlexLabelComplete(*dm, bdlabel)); 5449 } 5450 } 5451 PetscFunctionReturn(PETSC_SUCCESS); 5452 } 5453 5454 /*@C 5455 DMPlexCreateFromFile - This takes a filename and produces a `DM` 5456 5457 Collective 5458 5459 Input Parameters: 5460 + comm - The communicator 5461 . filename - A file name 5462 . plexname - The object name of the resulting `DM`, also used for intra-datafile lookup by some formats 5463 - interpolate - Flag to create intermediate mesh pieces (edges, faces) 5464 5465 Output Parameter: 5466 . dm - The `DM` 5467 5468 Options Database Key: 5469 . -dm_plex_create_from_hdf5_xdmf - use the `PETSC_VIEWER_HDF5_XDMF` format for reading HDF5 5470 5471 Use `-dm_plex_create_ prefix` to pass options to the internal `PetscViewer`, e.g. 5472 $ -dm_plex_create_viewer_hdf5_collective 5473 5474 Level: beginner 5475 5476 Notes: 5477 Using `PETSCVIEWERHDF5` type with `PETSC_VIEWER_HDF5_PETSC` format, one can save multiple `DMPLEX` 5478 meshes in a single HDF5 file. This in turn requires one to name the `DMPLEX` object with `PetscObjectSetName()` 5479 before saving it with `DMView()` and before loading it with `DMLoad()` for identification of the mesh object. 5480 The input parameter name is thus used to name the `DMPLEX` object when `DMPlexCreateFromFile()` internally 5481 calls `DMLoad()`. Currently, name is ignored for other viewer types and/or formats. 5482 5483 .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexCreateFromDAG()`, `DMPlexCreateFromCellListPetsc()`, `DMPlexCreate()`, `PetscObjectSetName()`, `DMView()`, `DMLoad()` 5484 @*/ 5485 PetscErrorCode DMPlexCreateFromFile(MPI_Comm comm, const char filename[], const char plexname[], PetscBool interpolate, DM *dm) 5486 { 5487 const char extGmsh[] = ".msh"; 5488 const char extGmsh2[] = ".msh2"; 5489 const char extGmsh4[] = ".msh4"; 5490 const char extCGNS[] = ".cgns"; 5491 const char extExodus[] = ".exo"; 5492 const char extExodus_e[] = ".e"; 5493 const char extGenesis[] = ".gen"; 5494 const char extFluent[] = ".cas"; 5495 const char extHDF5[] = ".h5"; 5496 const char extXDMFHDF5[] = ".xdmf.h5"; 5497 const char extMed[] = ".med"; 5498 const char extPLY[] = ".ply"; 5499 const char extEGADSLite[] = ".egadslite"; 5500 const char extEGADS[] = ".egads"; 5501 const char extIGES[] = ".igs"; 5502 const char extSTEP[] = ".stp"; 5503 const char extCV[] = ".dat"; 5504 size_t len; 5505 PetscBool isGmsh, isGmsh2, isGmsh4, isCGNS, isExodus, isGenesis, isFluent, isHDF5, isMed, isPLY, isEGADSLite, isEGADS, isIGES, isSTEP, isCV, isXDMFHDF5; 5506 PetscMPIInt rank; 5507 5508 PetscFunctionBegin; 5509 PetscAssertPointer(filename, 2); 5510 if (plexname) PetscAssertPointer(plexname, 3); 5511 PetscAssertPointer(dm, 5); 5512 PetscCall(DMInitializePackage()); 5513 PetscCall(PetscLogEventBegin(DMPLEX_CreateFromFile, 0, 0, 0, 0)); 5514 PetscCallMPI(MPI_Comm_rank(comm, &rank)); 5515 PetscCall(PetscStrlen(filename, &len)); 5516 PetscCheck(len, comm, PETSC_ERR_ARG_WRONG, "Filename must be a valid path"); 5517 5518 #define CheckExtension(extension__, is_extension__) \ 5519 do { \ 5520 PetscAssert(sizeof(extension__), comm, PETSC_ERR_PLIB, "Zero-size extension: %s", extension__); \ 5521 /* don't count the null-terminator at the end */ \ 5522 const size_t ext_len = sizeof(extension__) - 1; \ 5523 if (len < ext_len) { \ 5524 is_extension__ = PETSC_FALSE; \ 5525 } else { \ 5526 PetscCall(PetscStrncmp(filename + len - ext_len, extension__, ext_len, &is_extension__)); \ 5527 } \ 5528 } while (0) 5529 5530 CheckExtension(extGmsh, isGmsh); 5531 CheckExtension(extGmsh2, isGmsh2); 5532 CheckExtension(extGmsh4, isGmsh4); 5533 CheckExtension(extCGNS, isCGNS); 5534 CheckExtension(extExodus, isExodus); 5535 if (!isExodus) CheckExtension(extExodus_e, isExodus); 5536 CheckExtension(extGenesis, isGenesis); 5537 CheckExtension(extFluent, isFluent); 5538 CheckExtension(extHDF5, isHDF5); 5539 CheckExtension(extMed, isMed); 5540 CheckExtension(extPLY, isPLY); 5541 CheckExtension(extEGADSLite, isEGADSLite); 5542 CheckExtension(extEGADS, isEGADS); 5543 CheckExtension(extIGES, isIGES); 5544 CheckExtension(extSTEP, isSTEP); 5545 CheckExtension(extCV, isCV); 5546 CheckExtension(extXDMFHDF5, isXDMFHDF5); 5547 5548 #undef CheckExtension 5549 5550 if (isGmsh || isGmsh2 || isGmsh4) { 5551 PetscCall(DMPlexCreateGmshFromFile(comm, filename, interpolate, dm)); 5552 } else if (isCGNS) { 5553 PetscCall(DMPlexCreateCGNSFromFile(comm, filename, interpolate, dm)); 5554 } else if (isExodus || isGenesis) { 5555 PetscCall(DMPlexCreateExodusFromFile(comm, filename, interpolate, dm)); 5556 } else if (isFluent) { 5557 PetscCall(DMPlexCreateFluentFromFile(comm, filename, interpolate, dm)); 5558 } else if (isHDF5) { 5559 PetscViewer viewer; 5560 5561 /* PETSC_VIEWER_HDF5_XDMF is used if the filename ends with .xdmf.h5, or if -dm_plex_create_from_hdf5_xdmf option is present */ 5562 PetscCall(PetscOptionsGetBool(NULL, NULL, "-dm_plex_create_from_hdf5_xdmf", &isXDMFHDF5, NULL)); 5563 PetscCall(PetscViewerCreate(comm, &viewer)); 5564 PetscCall(PetscViewerSetType(viewer, PETSCVIEWERHDF5)); 5565 PetscCall(PetscViewerSetOptionsPrefix(viewer, "dm_plex_create_")); 5566 PetscCall(PetscViewerSetFromOptions(viewer)); 5567 PetscCall(PetscViewerFileSetMode(viewer, FILE_MODE_READ)); 5568 PetscCall(PetscViewerFileSetName(viewer, filename)); 5569 5570 PetscCall(DMCreate(comm, dm)); 5571 PetscCall(PetscObjectSetName((PetscObject)(*dm), plexname)); 5572 PetscCall(DMSetType(*dm, DMPLEX)); 5573 if (isXDMFHDF5) PetscCall(PetscViewerPushFormat(viewer, PETSC_VIEWER_HDF5_XDMF)); 5574 PetscCall(DMLoad(*dm, viewer)); 5575 if (isXDMFHDF5) PetscCall(PetscViewerPopFormat(viewer)); 5576 PetscCall(PetscViewerDestroy(&viewer)); 5577 5578 if (interpolate) { 5579 DM idm; 5580 5581 PetscCall(DMPlexInterpolate(*dm, &idm)); 5582 PetscCall(DMDestroy(dm)); 5583 *dm = idm; 5584 } 5585 } else if (isMed) { 5586 PetscCall(DMPlexCreateMedFromFile(comm, filename, interpolate, dm)); 5587 } else if (isPLY) { 5588 PetscCall(DMPlexCreatePLYFromFile(comm, filename, interpolate, dm)); 5589 } else if (isEGADSLite || isEGADS || isIGES || isSTEP) { 5590 if (isEGADSLite) PetscCall(DMPlexCreateEGADSLiteFromFile(comm, filename, dm)); 5591 else PetscCall(DMPlexCreateEGADSFromFile(comm, filename, dm)); 5592 if (!interpolate) { 5593 DM udm; 5594 5595 PetscCall(DMPlexUninterpolate(*dm, &udm)); 5596 PetscCall(DMDestroy(dm)); 5597 *dm = udm; 5598 } 5599 } else if (isCV) { 5600 PetscCall(DMPlexCreateCellVertexFromFile(comm, filename, interpolate, dm)); 5601 } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cannot load file %s: unrecognized extension", filename); 5602 PetscCall(PetscStrlen(plexname, &len)); 5603 if (len) PetscCall(PetscObjectSetName((PetscObject)(*dm), plexname)); 5604 PetscCall(PetscLogEventEnd(DMPLEX_CreateFromFile, 0, 0, 0, 0)); 5605 PetscFunctionReturn(PETSC_SUCCESS); 5606 } 5607 5608 /*@C 5609 DMPlexCreateEphemeral - This takes a `DMPlexTransform` and a base `DMPlex` and produces an ephemeral `DM`, meaning one that is created on the fly in response to queries. 5610 5611 Input Parameters: 5612 + tr - The `DMPlexTransform` 5613 - prefix - An options prefix, or NULL 5614 5615 Output Parameter: 5616 . dm - The `DM` 5617 5618 Level: beginner 5619 5620 Notes: 5621 An emphemeral mesh is one that is not stored concretely, as in the default `DMPLEX` implementation, but rather is produced on the fly in response to queries, using information from the transform and the base mesh. 5622 5623 .seealso: `DMPlexCreateFromFile`, `DMPlexCreateFromDAG()`, `DMPlexCreateFromCellListPetsc()`, `DMPlexCreate()` 5624 @*/ 5625 PetscErrorCode DMPlexCreateEphemeral(DMPlexTransform tr, const char prefix[], DM *dm) 5626 { 5627 DM bdm, bcdm, cdm; 5628 Vec coordinates, coordinatesNew; 5629 PetscSection cs; 5630 PetscInt dim, cdim, Nl; 5631 5632 PetscFunctionBegin; 5633 PetscCall(DMCreate(PetscObjectComm((PetscObject)tr), dm)); 5634 PetscCall(DMSetType(*dm, DMPLEX)); 5635 ((DM_Plex *)(*dm)->data)->interpolated = DMPLEX_INTERPOLATED_FULL; 5636 // Handle coordinates 5637 PetscCall(DMPlexTransformGetDM(tr, &bdm)); 5638 PetscCall(DMGetCoordinateDim(bdm, &cdim)); 5639 PetscCall(DMSetCoordinateDim(*dm, cdim)); 5640 PetscCall(DMGetDimension(bdm, &dim)); 5641 PetscCall(DMSetDimension(*dm, dim)); 5642 PetscCall(DMGetCoordinateDM(bdm, &bcdm)); 5643 PetscCall(DMGetCoordinateDM(*dm, &cdm)); 5644 PetscCall(DMCopyDisc(bcdm, cdm)); 5645 PetscCall(DMGetLocalSection(cdm, &cs)); 5646 PetscCall(PetscSectionSetNumFields(cs, 1)); 5647 PetscCall(PetscSectionSetFieldComponents(cs, 0, cdim)); 5648 PetscCall(DMGetCoordinatesLocal(bdm, &coordinates)); 5649 PetscCall(VecDuplicate(coordinates, &coordinatesNew)); 5650 PetscCall(VecCopy(coordinates, coordinatesNew)); 5651 PetscCall(DMSetCoordinatesLocal(*dm, coordinatesNew)); 5652 PetscCall(VecDestroy(&coordinatesNew)); 5653 5654 PetscCall(PetscObjectReference((PetscObject)tr)); 5655 PetscCall(DMPlexTransformDestroy(&((DM_Plex *)(*dm)->data)->tr)); 5656 ((DM_Plex *)(*dm)->data)->tr = tr; 5657 PetscCall(DMPlexDistributeSetDefault(*dm, PETSC_FALSE)); 5658 PetscCall(PetscObjectSetOptionsPrefix((PetscObject)*dm, prefix)); 5659 PetscCall(DMSetFromOptions(*dm)); 5660 5661 PetscCall(DMGetNumLabels(bdm, &Nl)); 5662 for (PetscInt l = 0; l < Nl; ++l) { 5663 DMLabel label, labelNew; 5664 const char *lname; 5665 PetscBool isDepth, isCellType; 5666 5667 PetscCall(DMGetLabelName(bdm, l, &lname)); 5668 PetscCall(PetscStrcmp(lname, "depth", &isDepth)); 5669 if (isDepth) continue; 5670 PetscCall(PetscStrcmp(lname, "celltype", &isCellType)); 5671 if (isCellType) continue; 5672 PetscCall(DMCreateLabel(*dm, lname)); 5673 PetscCall(DMGetLabel(bdm, lname, &label)); 5674 PetscCall(DMGetLabel(*dm, lname, &labelNew)); 5675 PetscCall(DMLabelSetType(labelNew, DMLABELEPHEMERAL)); 5676 PetscCall(DMLabelEphemeralSetLabel(labelNew, label)); 5677 PetscCall(DMLabelEphemeralSetTransform(labelNew, tr)); 5678 PetscCall(DMLabelSetUp(labelNew)); 5679 } 5680 PetscFunctionReturn(PETSC_SUCCESS); 5681 } 5682