1 #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 2 3 #ifdef PETSC_HAVE_EGADS 4 #include <egads.h> 5 #endif 6 7 #include <ctetgen.h> 8 9 /* This is to fix the tetrahedron orientation from TetGen */ 10 static PetscErrorCode DMPlexInvertCells_CTetgen(PetscInt numCells, PetscInt numCorners, PetscInt cells[]) 11 { 12 PetscInt bound = numCells*numCorners, coff; 13 14 PetscFunctionBegin; 15 #define SWAP(a,b) do { PetscInt tmp = (a); (a) = (b); (b) = tmp; } while (0) 16 for (coff = 0; coff < bound; coff += numCorners) SWAP(cells[coff],cells[coff+1]); 17 #undef SWAP 18 PetscFunctionReturn(0); 19 } 20 21 PETSC_EXTERN PetscErrorCode DMPlexGenerate_CTetgen(DM boundary, PetscBool interpolate, DM *dm) 22 { 23 MPI_Comm comm; 24 const PetscInt dim = 3; 25 PLC *in, *out; 26 DMUniversalLabel universal; 27 PetscInt vStart, vEnd, v, eStart, eEnd, e, fStart, fEnd, f, verbose = 0; 28 DMPlexInterpolatedFlag isInterpolated; 29 PetscMPIInt rank; 30 31 PetscFunctionBegin; 32 PetscCall(PetscOptionsGetInt(NULL,((PetscObject) boundary)->prefix, "-ctetgen_verbose", &verbose, NULL)); 33 PetscCall(PetscObjectGetComm((PetscObject)boundary,&comm)); 34 PetscCallMPI(MPI_Comm_rank(comm, &rank)); 35 PetscCall(DMPlexIsInterpolatedCollective(boundary, &isInterpolated)); 36 PetscCall(DMUniversalLabelCreate(boundary, &universal)); 37 38 PetscCall(PLCCreate(&in)); 39 PetscCall(PLCCreate(&out)); 40 41 PetscCall(DMPlexGetDepthStratum(boundary, 0, &vStart, &vEnd)); 42 in->numberofpoints = vEnd - vStart; 43 if (in->numberofpoints > 0) { 44 PetscSection coordSection; 45 Vec coordinates; 46 const PetscScalar *array; 47 48 PetscCall(PetscMalloc1(in->numberofpoints*dim, &in->pointlist)); 49 PetscCall(PetscMalloc1(in->numberofpoints, &in->pointmarkerlist)); 50 PetscCall(DMGetCoordinatesLocal(boundary, &coordinates)); 51 PetscCall(DMGetCoordinateSection(boundary, &coordSection)); 52 PetscCall(VecGetArrayRead(coordinates, &array)); 53 for (v = vStart; v < vEnd; ++v) { 54 const PetscInt idx = v - vStart; 55 PetscInt off, d, m; 56 57 PetscCall(PetscSectionGetOffset(coordSection, v, &off)); 58 for (d = 0; d < dim; ++d) in->pointlist[idx*dim + d] = PetscRealPart(array[off+d]); 59 PetscCall(DMLabelGetValue(universal->label, v, &m)); 60 in->pointmarkerlist[idx] = (int) m; 61 } 62 PetscCall(VecRestoreArrayRead(coordinates, &array)); 63 } 64 65 PetscCall(DMPlexGetHeightStratum(boundary, 1, &eStart, &eEnd)); 66 in->numberofedges = eEnd - eStart; 67 if (isInterpolated == DMPLEX_INTERPOLATED_FULL && in->numberofedges > 0) { 68 PetscCall(PetscMalloc1(in->numberofedges*2, &in->edgelist)); 69 PetscCall(PetscMalloc1(in->numberofedges, &in->edgemarkerlist)); 70 for (e = eStart; e < eEnd; ++e) { 71 const PetscInt idx = e - eStart; 72 const PetscInt *cone; 73 PetscInt coneSize, val; 74 75 PetscCall(DMPlexGetConeSize(boundary, e, &coneSize)); 76 PetscCall(DMPlexGetCone(boundary, e, &cone)); 77 in->edgelist[idx*2] = cone[0] - vStart; 78 in->edgelist[idx*2 + 1] = cone[1] - vStart; 79 80 PetscCall(DMLabelGetValue(universal->label, e, &val)); 81 in->edgemarkerlist[idx] = (int) val; 82 } 83 } 84 85 PetscCall(DMPlexGetHeightStratum(boundary, 0, &fStart, &fEnd)); 86 in->numberoffacets = fEnd - fStart; 87 if (in->numberoffacets > 0) { 88 PetscCall(PetscMalloc1(in->numberoffacets, &in->facetlist)); 89 PetscCall(PetscMalloc1(in->numberoffacets, &in->facetmarkerlist)); 90 for (f = fStart; f < fEnd; ++f) { 91 const PetscInt idx = f - fStart; 92 PetscInt *points = NULL, numPoints, p, numVertices = 0, v, m = -1; 93 polygon *poly; 94 95 in->facetlist[idx].numberofpolygons = 1; 96 PetscCall(PetscMalloc1(in->facetlist[idx].numberofpolygons, &in->facetlist[idx].polygonlist)); 97 in->facetlist[idx].numberofholes = 0; 98 in->facetlist[idx].holelist = NULL; 99 100 PetscCall(DMPlexGetTransitiveClosure(boundary, f, PETSC_TRUE, &numPoints, &points)); 101 for (p = 0; p < numPoints*2; p += 2) { 102 const PetscInt point = points[p]; 103 if ((point >= vStart) && (point < vEnd)) points[numVertices++] = point; 104 } 105 106 poly = in->facetlist[idx].polygonlist; 107 poly->numberofvertices = numVertices; 108 PetscCall(PetscMalloc1(poly->numberofvertices, &poly->vertexlist)); 109 for (v = 0; v < numVertices; ++v) { 110 const PetscInt vIdx = points[v] - vStart; 111 poly->vertexlist[v] = vIdx; 112 } 113 PetscCall(DMLabelGetValue(universal->label, f, &m)); 114 in->facetmarkerlist[idx] = (int) m; 115 PetscCall(DMPlexRestoreTransitiveClosure(boundary, f, PETSC_TRUE, &numPoints, &points)); 116 } 117 } 118 if (rank == 0) { 119 TetGenOpts t; 120 121 PetscCall(TetGenOptsInitialize(&t)); 122 t.in = boundary; /* Should go away */ 123 t.plc = 1; 124 t.quality = 1; 125 t.edgesout = 1; 126 t.zeroindex = 1; 127 t.quiet = 1; 128 t.verbose = verbose; 129 #if 0 130 #ifdef PETSC_HAVE_EGADS 131 /* Need to add in more TetGen code */ 132 t.nobisect = 1; /* Add Y to preserve Surface Mesh for EGADS */ 133 #endif 134 #endif 135 136 PetscCall(TetGenCheckOpts(&t)); 137 PetscCall(TetGenTetrahedralize(&t, in, out)); 138 } 139 { 140 const PetscInt numCorners = 4; 141 const PetscInt numCells = out->numberoftetrahedra; 142 const PetscInt numVertices = out->numberofpoints; 143 PetscReal *meshCoords = NULL; 144 PetscInt *cells = NULL; 145 146 if (sizeof (PetscReal) == sizeof (out->pointlist[0])) { 147 meshCoords = (PetscReal *) out->pointlist; 148 } else { 149 PetscInt i; 150 151 PetscCall(PetscMalloc1(dim * numVertices, &meshCoords)); 152 for (i = 0; i < dim * numVertices; ++i) meshCoords[i] = (PetscReal) out->pointlist[i]; 153 } 154 if (sizeof (PetscInt) == sizeof (out->tetrahedronlist[0])) { 155 cells = (PetscInt *) out->tetrahedronlist; 156 } else { 157 PetscInt i; 158 159 PetscCall(PetscMalloc1(numCells * numCorners, &cells)); 160 for (i = 0; i < numCells * numCorners; i++) cells[i] = (PetscInt) out->tetrahedronlist[i]; 161 } 162 163 PetscCall(DMPlexInvertCells_CTetgen(numCells, numCorners, cells)); 164 PetscCall(DMPlexCreateFromCellListPetsc(comm, dim, numCells, numVertices, numCorners, interpolate, cells, dim, meshCoords, dm)); 165 if (sizeof (PetscReal) != sizeof (out->pointlist[0])) { 166 PetscCall(PetscFree(meshCoords)); 167 } 168 if (sizeof (PetscInt) != sizeof (out->tetrahedronlist[0])) { 169 PetscCall(PetscFree(cells)); 170 } 171 172 /* Set labels */ 173 PetscCall(DMUniversalLabelCreateLabels(universal, PETSC_TRUE, *dm)); 174 for (v = 0; v < numVertices; ++v) { 175 if (out->pointmarkerlist[v]) { 176 PetscCall(DMUniversalLabelSetLabelValue(universal, *dm, PETSC_TRUE, v+numCells, out->pointmarkerlist[v])); 177 } 178 } 179 if (interpolate) { 180 PetscInt e; 181 182 for (e = 0; e < out->numberofedges; e++) { 183 if (out->edgemarkerlist[e]) { 184 const PetscInt vertices[2] = {out->edgelist[e*2+0]+numCells, out->edgelist[e*2+1]+numCells}; 185 const PetscInt *edges; 186 PetscInt numEdges; 187 188 PetscCall(DMPlexGetJoin(*dm, 2, vertices, &numEdges, &edges)); 189 PetscCheck(numEdges == 1,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Two vertices must cover only one edge, not %" PetscInt_FMT, numEdges); 190 PetscCall(DMUniversalLabelSetLabelValue(universal, *dm, PETSC_TRUE, edges[0], out->edgemarkerlist[e])); 191 PetscCall(DMPlexRestoreJoin(*dm, 2, vertices, &numEdges, &edges)); 192 } 193 } 194 for (f = 0; f < out->numberoftrifaces; f++) { 195 if (out->trifacemarkerlist[f]) { 196 const PetscInt vertices[3] = {out->trifacelist[f*3+0]+numCells, out->trifacelist[f*3+1]+numCells, out->trifacelist[f*3+2]+numCells}; 197 const PetscInt *faces; 198 PetscInt numFaces; 199 200 PetscCall(DMPlexGetFullJoin(*dm, 3, vertices, &numFaces, &faces)); 201 PetscCheck(numFaces == 1,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Three vertices must cover only one face, not %" PetscInt_FMT, numFaces); 202 PetscCall(DMUniversalLabelSetLabelValue(universal, *dm, PETSC_TRUE, faces[0], out->trifacemarkerlist[f])); 203 PetscCall(DMPlexRestoreJoin(*dm, 3, vertices, &numFaces, &faces)); 204 } 205 } 206 } 207 208 #ifdef PETSC_HAVE_EGADS 209 { 210 DMLabel bodyLabel; 211 PetscContainer modelObj; 212 PetscInt cStart, cEnd, c, eStart, eEnd, fStart, fEnd; 213 ego *bodies; 214 ego model, geom; 215 int Nb, oclass, mtype, *senses; 216 217 /* Get Attached EGADS Model from Original DMPlex */ 218 PetscCall(PetscObjectQuery((PetscObject) boundary, "EGADS Model", (PetscObject *) &modelObj)); 219 if (modelObj) { 220 PetscCall(PetscContainerGetPointer(modelObj, (void **) &model)); 221 PetscCall(EG_getTopology(model, &geom, &oclass, &mtype, NULL, &Nb, &bodies, &senses)); 222 /* Transfer EGADS Model to Volumetric Mesh */ 223 PetscCall(PetscObjectCompose((PetscObject) *dm, "EGADS Model", (PetscObject) modelObj)); 224 225 /* Set Cell Labels */ 226 PetscCall(DMGetLabel(*dm, "EGADS Body ID", &bodyLabel)); 227 PetscCall(DMPlexGetHeightStratum(*dm, 0, &cStart, &cEnd)); 228 PetscCall(DMPlexGetHeightStratum(*dm, 1, &fStart, &fEnd)); 229 PetscCall(DMPlexGetDepthStratum(*dm, 1, &eStart, &eEnd)); 230 231 for (c = cStart; c < cEnd; ++c) { 232 PetscReal centroid[3] = {0., 0., 0.}; 233 PetscInt b; 234 235 /* Deterimine what body the cell's centroid is located in */ 236 if (!interpolate) { 237 PetscSection coordSection; 238 Vec coordinates; 239 PetscScalar *coords = NULL; 240 PetscInt coordSize, s, d; 241 242 PetscCall(DMGetCoordinatesLocal(*dm, &coordinates)); 243 PetscCall(DMGetCoordinateSection(*dm, &coordSection)); 244 PetscCall(DMPlexVecGetClosure(*dm, coordSection, coordinates, c, &coordSize, &coords)); 245 for (s = 0; s < coordSize; ++s) for (d = 0; d < dim; ++d) centroid[d] += coords[s*dim+d]; 246 PetscCall(DMPlexVecRestoreClosure(*dm, coordSection, coordinates, c, &coordSize, &coords)); 247 } else PetscCall(DMPlexComputeCellGeometryFVM(*dm, c, NULL, centroid, NULL)); 248 for (b = 0; b < Nb; ++b) { 249 if (EG_inTopology(bodies[b], centroid) == EGADS_SUCCESS) break; 250 } 251 if (b < Nb) { 252 PetscInt cval = b, eVal, fVal; 253 PetscInt *closure = NULL, Ncl, cl; 254 255 PetscCall(DMLabelSetValue(bodyLabel, c, cval)); 256 PetscCall(DMPlexGetTransitiveClosure(*dm, c, PETSC_TRUE, &Ncl, &closure)); 257 for (cl = 0; cl < Ncl; ++cl) { 258 const PetscInt p = closure[cl*2]; 259 260 if (p >= eStart && p < eEnd) { 261 PetscCall(DMLabelGetValue(bodyLabel, p, &eVal)); 262 if (eVal < 0) PetscCall(DMLabelSetValue(bodyLabel, p, cval)); 263 } 264 if (p >= fStart && p < fEnd) { 265 PetscCall(DMLabelGetValue(bodyLabel, p, &fVal)); 266 if (fVal < 0) PetscCall(DMLabelSetValue(bodyLabel, p, cval)); 267 } 268 } 269 PetscCall(DMPlexRestoreTransitiveClosure(*dm, c, PETSC_TRUE, &Ncl, &closure)); 270 } 271 } 272 } 273 } 274 #endif 275 PetscCall(DMPlexSetRefinementUniform(*dm, PETSC_FALSE)); 276 } 277 278 PetscCall(DMUniversalLabelDestroy(&universal)); 279 PetscCall(PLCDestroy(&in)); 280 PetscCall(PLCDestroy(&out)); 281 PetscFunctionReturn(0); 282 } 283 284 PETSC_EXTERN PetscErrorCode DMPlexRefine_CTetgen(DM dm, PetscReal *maxVolumes, DM *dmRefined) 285 { 286 MPI_Comm comm; 287 const PetscInt dim = 3; 288 PLC *in, *out; 289 DMUniversalLabel universal; 290 PetscInt vStart, vEnd, v, eStart, eEnd, e, fStart, fEnd, f, cStart, cEnd, c, verbose = 0; 291 DMPlexInterpolatedFlag isInterpolated; 292 PetscMPIInt rank; 293 294 PetscFunctionBegin; 295 PetscCall(PetscOptionsGetInt(NULL,((PetscObject) dm)->prefix, "-ctetgen_verbose", &verbose, NULL)); 296 PetscCall(PetscObjectGetComm((PetscObject)dm,&comm)); 297 PetscCallMPI(MPI_Comm_rank(comm, &rank)); 298 PetscCall(DMPlexIsInterpolatedCollective(dm, &isInterpolated)); 299 PetscCall(DMUniversalLabelCreate(dm, &universal)); 300 301 PetscCall(PLCCreate(&in)); 302 PetscCall(PLCCreate(&out)); 303 304 PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd)); 305 in->numberofpoints = vEnd - vStart; 306 if (in->numberofpoints > 0) { 307 PetscSection coordSection; 308 Vec coordinates; 309 PetscScalar *array; 310 311 PetscCall(PetscMalloc1(in->numberofpoints*dim, &in->pointlist)); 312 PetscCall(PetscMalloc1(in->numberofpoints, &in->pointmarkerlist)); 313 PetscCall(DMGetCoordinatesLocal(dm, &coordinates)); 314 PetscCall(DMGetCoordinateSection(dm, &coordSection)); 315 PetscCall(VecGetArray(coordinates, &array)); 316 for (v = vStart; v < vEnd; ++v) { 317 const PetscInt idx = v - vStart; 318 PetscInt off, d, m; 319 320 PetscCall(PetscSectionGetOffset(coordSection, v, &off)); 321 for (d = 0; d < dim; ++d) in->pointlist[idx*dim + d] = PetscRealPart(array[off+d]); 322 PetscCall(DMLabelGetValue(universal->label, v, &m)); 323 in->pointmarkerlist[idx] = (int) m; 324 } 325 PetscCall(VecRestoreArray(coordinates, &array)); 326 } 327 328 PetscCall(DMPlexGetDepthStratum(dm, 1, &eStart, &eEnd)); 329 in->numberofedges = eEnd - eStart; 330 if (isInterpolated == DMPLEX_INTERPOLATED_FULL && in->numberofedges > 0) { 331 PetscCall(PetscMalloc1(in->numberofedges * 2, &in->edgelist)); 332 PetscCall(PetscMalloc1(in->numberofedges, &in->edgemarkerlist)); 333 for (e = eStart; e < eEnd; ++e) { 334 const PetscInt idx = e - eStart; 335 const PetscInt *cone; 336 PetscInt coneSize, val; 337 338 PetscCall(DMPlexGetConeSize(dm, e, &coneSize)); 339 PetscCall(DMPlexGetCone(dm, e, &cone)); 340 in->edgelist[idx*2] = cone[0] - vStart; 341 in->edgelist[idx*2 + 1] = cone[1] - vStart; 342 343 PetscCall(DMLabelGetValue(universal->label, e, &val)); 344 in->edgemarkerlist[idx] = (int) val; 345 } 346 } 347 348 PetscCall(DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd)); 349 in->numberoftrifaces = 0; 350 for (f = fStart; f < fEnd; ++f) { 351 PetscInt supportSize; 352 353 PetscCall(DMPlexGetSupportSize(dm, f, &supportSize)); 354 if (supportSize == 1) ++in->numberoftrifaces; 355 } 356 if (isInterpolated == DMPLEX_INTERPOLATED_FULL && in->numberoftrifaces > 0) { 357 PetscInt tf = 0; 358 359 PetscCall(PetscMalloc1(in->numberoftrifaces*3, &in->trifacelist)); 360 PetscCall(PetscMalloc1(in->numberoftrifaces, &in->trifacemarkerlist)); 361 for (f = fStart; f < fEnd; ++f) { 362 PetscInt *points = NULL; 363 PetscInt supportSize, numPoints, p, Nv = 0, val; 364 365 PetscCall(DMPlexGetSupportSize(dm, f, &supportSize)); 366 if (supportSize != 1) continue; 367 PetscCall(DMPlexGetTransitiveClosure(dm, f, PETSC_TRUE, &numPoints, &points)); 368 for (p = 0; p < numPoints*2; p += 2) { 369 const PetscInt point = points[p]; 370 if ((point >= vStart) && (point < vEnd)) in->trifacelist[tf*3 + Nv++] = point - vStart; 371 } 372 PetscCall(DMPlexRestoreTransitiveClosure(dm, f, PETSC_TRUE, &numPoints, &points)); 373 PetscCheck(Nv == 3,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Face %" PetscInt_FMT " has %" PetscInt_FMT " vertices, not 3", f, Nv); 374 PetscCall(DMLabelGetValue(universal->label, f, &val)); 375 in->trifacemarkerlist[tf] = (int) val; 376 ++tf; 377 } 378 } 379 380 PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 381 in->numberofcorners = 4; 382 in->numberoftetrahedra = cEnd - cStart; 383 in->tetrahedronvolumelist = maxVolumes; 384 if (in->numberoftetrahedra > 0) { 385 PetscCall(PetscMalloc1(in->numberoftetrahedra*in->numberofcorners, &in->tetrahedronlist)); 386 for (c = cStart; c < cEnd; ++c) { 387 const PetscInt idx = c - cStart; 388 PetscInt *closure = NULL; 389 PetscInt closureSize; 390 391 PetscCall(DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure)); 392 PetscCheck((closureSize == 5) || (closureSize == 15),comm, PETSC_ERR_ARG_WRONG, "Mesh has cell which is not a tetrahedron, %" PetscInt_FMT " vertices in closure", closureSize); 393 for (v = 0; v < 4; ++v) in->tetrahedronlist[idx*in->numberofcorners + v] = closure[(v+closureSize-4)*2] - vStart; 394 PetscCall(DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure)); 395 } 396 } 397 398 if (rank == 0) { 399 TetGenOpts t; 400 401 PetscCall(TetGenOptsInitialize(&t)); 402 403 t.in = dm; /* Should go away */ 404 t.refine = 1; 405 t.varvolume = 1; 406 t.quality = 1; 407 t.edgesout = 1; 408 t.zeroindex = 1; 409 t.quiet = 1; 410 t.verbose = verbose; /* Change this */ 411 412 PetscCall(TetGenCheckOpts(&t)); 413 PetscCall(TetGenTetrahedralize(&t, in, out)); 414 } 415 416 in->tetrahedronvolumelist = NULL; 417 { 418 const PetscInt numCorners = 4; 419 const PetscInt numCells = out->numberoftetrahedra; 420 const PetscInt numVertices = out->numberofpoints; 421 PetscReal *meshCoords = NULL; 422 PetscInt *cells = NULL; 423 PetscBool interpolate = isInterpolated == DMPLEX_INTERPOLATED_FULL ? PETSC_TRUE : PETSC_FALSE; 424 425 if (sizeof (PetscReal) == sizeof (out->pointlist[0])) { 426 meshCoords = (PetscReal *) out->pointlist; 427 } else { 428 PetscInt i; 429 430 PetscCall(PetscMalloc1(dim * numVertices, &meshCoords)); 431 for (i = 0; i < dim * numVertices; ++i) meshCoords[i] = (PetscReal) out->pointlist[i]; 432 } 433 if (sizeof (PetscInt) == sizeof (out->tetrahedronlist[0])) { 434 cells = (PetscInt *) out->tetrahedronlist; 435 } else { 436 PetscInt i; 437 438 PetscCall(PetscMalloc1(numCells * numCorners, &cells)); 439 for (i = 0; i < numCells * numCorners; ++i) cells[i] = (PetscInt) out->tetrahedronlist[i]; 440 } 441 442 PetscCall(DMPlexInvertCells_CTetgen(numCells, numCorners, cells)); 443 PetscCall(DMPlexCreateFromCellListPetsc(comm, dim, numCells, numVertices, numCorners, interpolate, cells, dim, meshCoords, dmRefined)); 444 if (sizeof (PetscReal) != sizeof (out->pointlist[0])) PetscCall(PetscFree(meshCoords)); 445 if (sizeof (PetscInt) != sizeof (out->tetrahedronlist[0])) PetscCall(PetscFree(cells)); 446 447 /* Set labels */ 448 PetscCall(DMUniversalLabelCreateLabels(universal, PETSC_TRUE, *dmRefined)); 449 for (v = 0; v < numVertices; ++v) { 450 if (out->pointmarkerlist[v]) { 451 PetscCall(DMUniversalLabelSetLabelValue(universal, *dmRefined, PETSC_TRUE, v+numCells, out->pointmarkerlist[v])); 452 } 453 } 454 if (interpolate) { 455 PetscInt e, f; 456 457 for (e = 0; e < out->numberofedges; e++) { 458 if (out->edgemarkerlist[e]) { 459 const PetscInt vertices[2] = {out->edgelist[e*2+0]+numCells, out->edgelist[e*2+1]+numCells}; 460 const PetscInt *edges; 461 PetscInt numEdges; 462 463 PetscCall(DMPlexGetJoin(*dmRefined, 2, vertices, &numEdges, &edges)); 464 PetscCheck(numEdges == 1,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Two vertices must cover only one edge, not %" PetscInt_FMT, numEdges); 465 PetscCall(DMUniversalLabelSetLabelValue(universal, *dmRefined, PETSC_TRUE, edges[0], out->edgemarkerlist[e])); 466 PetscCall(DMPlexRestoreJoin(*dmRefined, 2, vertices, &numEdges, &edges)); 467 } 468 } 469 for (f = 0; f < out->numberoftrifaces; f++) { 470 if (out->trifacemarkerlist[f]) { 471 const PetscInt vertices[3] = {out->trifacelist[f*3+0]+numCells, out->trifacelist[f*3+1]+numCells, out->trifacelist[f*3+2]+numCells}; 472 const PetscInt *faces; 473 PetscInt numFaces; 474 475 PetscCall(DMPlexGetFullJoin(*dmRefined, 3, vertices, &numFaces, &faces)); 476 PetscCheck(numFaces == 1,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Three vertices must cover only one face, not %" PetscInt_FMT, numFaces); 477 PetscCall(DMUniversalLabelSetLabelValue(universal, *dmRefined, PETSC_TRUE, faces[0], out->trifacemarkerlist[f])); 478 PetscCall(DMPlexRestoreJoin(*dmRefined, 3, vertices, &numFaces, &faces)); 479 } 480 } 481 } 482 483 #ifdef PETSC_HAVE_EGADS 484 { 485 DMLabel bodyLabel; 486 PetscContainer modelObj; 487 PetscInt cStart, cEnd, c, eStart, eEnd, fStart, fEnd; 488 ego *bodies; 489 ego model, geom; 490 int Nb, oclass, mtype, *senses; 491 492 /* Get Attached EGADS Model from Original DMPlex */ 493 PetscCall(PetscObjectQuery((PetscObject) dm, "EGADS Model", (PetscObject *) &modelObj)); 494 if (modelObj) { 495 PetscCall(PetscContainerGetPointer(modelObj, (void **) &model)); 496 PetscCall(EG_getTopology(model, &geom, &oclass, &mtype, NULL, &Nb, &bodies, &senses)); 497 /* Transfer EGADS Model to Volumetric Mesh */ 498 PetscCall(PetscObjectCompose((PetscObject) *dmRefined, "EGADS Model", (PetscObject) modelObj)); 499 500 /* Set Cell Labels */ 501 PetscCall(DMGetLabel(*dmRefined, "EGADS Body ID", &bodyLabel)); 502 PetscCall(DMPlexGetHeightStratum(*dmRefined, 0, &cStart, &cEnd)); 503 PetscCall(DMPlexGetHeightStratum(*dmRefined, 1, &fStart, &fEnd)); 504 PetscCall(DMPlexGetDepthStratum(*dmRefined, 1, &eStart, &eEnd)); 505 506 for (c = cStart; c < cEnd; ++c) { 507 PetscReal centroid[3] = {0., 0., 0.}; 508 PetscInt b; 509 510 /* Deterimine what body the cell's centroid is located in */ 511 if (!interpolate) { 512 PetscSection coordSection; 513 Vec coordinates; 514 PetscScalar *coords = NULL; 515 PetscInt coordSize, s, d; 516 517 PetscCall(DMGetCoordinatesLocal(*dmRefined, &coordinates)); 518 PetscCall(DMGetCoordinateSection(*dmRefined, &coordSection)); 519 PetscCall(DMPlexVecGetClosure(*dmRefined, coordSection, coordinates, c, &coordSize, &coords)); 520 for (s = 0; s < coordSize; ++s) for (d = 0; d < dim; ++d) centroid[d] += coords[s*dim+d]; 521 PetscCall(DMPlexVecRestoreClosure(*dmRefined, coordSection, coordinates, c, &coordSize, &coords)); 522 } else PetscCall(DMPlexComputeCellGeometryFVM(*dmRefined, c, NULL, centroid, NULL)); 523 for (b = 0; b < Nb; ++b) { 524 if (EG_inTopology(bodies[b], centroid) == EGADS_SUCCESS) break; 525 } 526 if (b < Nb) { 527 PetscInt cval = b, eVal, fVal; 528 PetscInt *closure = NULL, Ncl, cl; 529 530 PetscCall(DMLabelSetValue(bodyLabel, c, cval)); 531 PetscCall(DMPlexGetTransitiveClosure(*dmRefined, c, PETSC_TRUE, &Ncl, &closure)); 532 for (cl = 0; cl < Ncl; cl += 2) { 533 const PetscInt p = closure[cl]; 534 535 if (p >= eStart && p < eEnd) { 536 PetscCall(DMLabelGetValue(bodyLabel, p, &eVal)); 537 if (eVal < 0) PetscCall(DMLabelSetValue(bodyLabel, p, cval)); 538 } 539 if (p >= fStart && p < fEnd) { 540 PetscCall(DMLabelGetValue(bodyLabel, p, &fVal)); 541 if (fVal < 0) PetscCall(DMLabelSetValue(bodyLabel, p, cval)); 542 } 543 } 544 PetscCall(DMPlexRestoreTransitiveClosure(*dmRefined, c, PETSC_TRUE, &Ncl, &closure)); 545 } 546 } 547 } 548 } 549 #endif 550 PetscCall(DMPlexSetRefinementUniform(*dmRefined, PETSC_FALSE)); 551 } 552 PetscCall(DMUniversalLabelDestroy(&universal)); 553 PetscCall(PLCDestroy(&in)); 554 PetscCall(PLCDestroy(&out)); 555 PetscFunctionReturn(0); 556 } 557