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