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 %D", 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 %D", 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 { 248 PetscCall(DMPlexComputeCellGeometryFVM(*dm, c, NULL, centroid, NULL)); 249 } 250 for (b = 0; b < Nb; ++b) { 251 if (EG_inTopology(bodies[b], centroid) == EGADS_SUCCESS) break; 252 } 253 if (b < Nb) { 254 PetscInt cval = b, eVal, fVal; 255 PetscInt *closure = NULL, Ncl, cl; 256 257 PetscCall(DMLabelSetValue(bodyLabel, c, cval)); 258 PetscCall(DMPlexGetTransitiveClosure(*dm, c, PETSC_TRUE, &Ncl, &closure)); 259 for (cl = 0; cl < Ncl; ++cl) { 260 const PetscInt p = closure[cl*2]; 261 262 if (p >= eStart && p < eEnd) { 263 PetscCall(DMLabelGetValue(bodyLabel, p, &eVal)); 264 if (eVal < 0) PetscCall(DMLabelSetValue(bodyLabel, p, cval)); 265 } 266 if (p >= fStart && p < fEnd) { 267 PetscCall(DMLabelGetValue(bodyLabel, p, &fVal)); 268 if (fVal < 0) PetscCall(DMLabelSetValue(bodyLabel, p, cval)); 269 } 270 } 271 PetscCall(DMPlexRestoreTransitiveClosure(*dm, c, PETSC_TRUE, &Ncl, &closure)); 272 } 273 } 274 } 275 } 276 #endif 277 PetscCall(DMPlexSetRefinementUniform(*dm, PETSC_FALSE)); 278 } 279 280 PetscCall(DMUniversalLabelDestroy(&universal)); 281 PetscCall(PLCDestroy(&in)); 282 PetscCall(PLCDestroy(&out)); 283 PetscFunctionReturn(0); 284 } 285 286 PETSC_EXTERN PetscErrorCode DMPlexRefine_CTetgen(DM dm, PetscReal *maxVolumes, DM *dmRefined) 287 { 288 MPI_Comm comm; 289 const PetscInt dim = 3; 290 PLC *in, *out; 291 DMUniversalLabel universal; 292 PetscInt vStart, vEnd, v, eStart, eEnd, e, fStart, fEnd, f, cStart, cEnd, c, verbose = 0; 293 DMPlexInterpolatedFlag isInterpolated; 294 PetscMPIInt rank; 295 296 PetscFunctionBegin; 297 PetscCall(PetscOptionsGetInt(NULL,((PetscObject) dm)->prefix, "-ctetgen_verbose", &verbose, NULL)); 298 PetscCall(PetscObjectGetComm((PetscObject)dm,&comm)); 299 PetscCallMPI(MPI_Comm_rank(comm, &rank)); 300 PetscCall(DMPlexIsInterpolatedCollective(dm, &isInterpolated)); 301 PetscCall(DMUniversalLabelCreate(dm, &universal)); 302 303 PetscCall(PLCCreate(&in)); 304 PetscCall(PLCCreate(&out)); 305 306 PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd)); 307 in->numberofpoints = vEnd - vStart; 308 if (in->numberofpoints > 0) { 309 PetscSection coordSection; 310 Vec coordinates; 311 PetscScalar *array; 312 313 PetscCall(PetscMalloc1(in->numberofpoints*dim, &in->pointlist)); 314 PetscCall(PetscMalloc1(in->numberofpoints, &in->pointmarkerlist)); 315 PetscCall(DMGetCoordinatesLocal(dm, &coordinates)); 316 PetscCall(DMGetCoordinateSection(dm, &coordSection)); 317 PetscCall(VecGetArray(coordinates, &array)); 318 for (v = vStart; v < vEnd; ++v) { 319 const PetscInt idx = v - vStart; 320 PetscInt off, d, m; 321 322 PetscCall(PetscSectionGetOffset(coordSection, v, &off)); 323 for (d = 0; d < dim; ++d) in->pointlist[idx*dim + d] = PetscRealPart(array[off+d]); 324 PetscCall(DMLabelGetValue(universal->label, v, &m)); 325 in->pointmarkerlist[idx] = (int) m; 326 } 327 PetscCall(VecRestoreArray(coordinates, &array)); 328 } 329 330 PetscCall(DMPlexGetDepthStratum(dm, 1, &eStart, &eEnd)); 331 in->numberofedges = eEnd - eStart; 332 if (isInterpolated == DMPLEX_INTERPOLATED_FULL && in->numberofedges > 0) { 333 PetscCall(PetscMalloc1(in->numberofedges * 2, &in->edgelist)); 334 PetscCall(PetscMalloc1(in->numberofedges, &in->edgemarkerlist)); 335 for (e = eStart; e < eEnd; ++e) { 336 const PetscInt idx = e - eStart; 337 const PetscInt *cone; 338 PetscInt coneSize, val; 339 340 PetscCall(DMPlexGetConeSize(dm, e, &coneSize)); 341 PetscCall(DMPlexGetCone(dm, e, &cone)); 342 in->edgelist[idx*2] = cone[0] - vStart; 343 in->edgelist[idx*2 + 1] = cone[1] - vStart; 344 345 PetscCall(DMLabelGetValue(universal->label, e, &val)); 346 in->edgemarkerlist[idx] = (int) val; 347 } 348 } 349 350 PetscCall(DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd)); 351 in->numberoftrifaces = 0; 352 for (f = fStart; f < fEnd; ++f) { 353 PetscInt supportSize; 354 355 PetscCall(DMPlexGetSupportSize(dm, f, &supportSize)); 356 if (supportSize == 1) ++in->numberoftrifaces; 357 } 358 if (isInterpolated == DMPLEX_INTERPOLATED_FULL && in->numberoftrifaces > 0) { 359 PetscInt tf = 0; 360 361 PetscCall(PetscMalloc1(in->numberoftrifaces*3, &in->trifacelist)); 362 PetscCall(PetscMalloc1(in->numberoftrifaces, &in->trifacemarkerlist)); 363 for (f = fStart; f < fEnd; ++f) { 364 PetscInt *points = NULL; 365 PetscInt supportSize, numPoints, p, Nv = 0, val; 366 367 PetscCall(DMPlexGetSupportSize(dm, f, &supportSize)); 368 if (supportSize != 1) continue; 369 PetscCall(DMPlexGetTransitiveClosure(dm, f, PETSC_TRUE, &numPoints, &points)); 370 for (p = 0; p < numPoints*2; p += 2) { 371 const PetscInt point = points[p]; 372 if ((point >= vStart) && (point < vEnd)) in->trifacelist[tf*3 + Nv++] = point - vStart; 373 } 374 PetscCall(DMPlexRestoreTransitiveClosure(dm, f, PETSC_TRUE, &numPoints, &points)); 375 PetscCheck(Nv == 3,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Face %D has %D vertices, not 3", f, Nv); 376 PetscCall(DMLabelGetValue(universal->label, f, &val)); 377 in->trifacemarkerlist[tf] = (int) val; 378 ++tf; 379 } 380 } 381 382 PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 383 in->numberofcorners = 4; 384 in->numberoftetrahedra = cEnd - cStart; 385 in->tetrahedronvolumelist = maxVolumes; 386 if (in->numberoftetrahedra > 0) { 387 PetscCall(PetscMalloc1(in->numberoftetrahedra*in->numberofcorners, &in->tetrahedronlist)); 388 for (c = cStart; c < cEnd; ++c) { 389 const PetscInt idx = c - cStart; 390 PetscInt *closure = NULL; 391 PetscInt closureSize; 392 393 PetscCall(DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure)); 394 PetscCheck((closureSize == 5) || (closureSize == 15),comm, PETSC_ERR_ARG_WRONG, "Mesh has cell which is not a tetrahedron, %" PetscInt_FMT " vertices in closure", closureSize); 395 for (v = 0; v < 4; ++v) in->tetrahedronlist[idx*in->numberofcorners + v] = closure[(v+closureSize-4)*2] - vStart; 396 PetscCall(DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure)); 397 } 398 } 399 400 if (rank == 0) { 401 TetGenOpts t; 402 403 PetscCall(TetGenOptsInitialize(&t)); 404 405 t.in = dm; /* Should go away */ 406 t.refine = 1; 407 t.varvolume = 1; 408 t.quality = 1; 409 t.edgesout = 1; 410 t.zeroindex = 1; 411 t.quiet = 1; 412 t.verbose = verbose; /* Change this */ 413 414 PetscCall(TetGenCheckOpts(&t)); 415 PetscCall(TetGenTetrahedralize(&t, in, out)); 416 } 417 418 in->tetrahedronvolumelist = NULL; 419 { 420 const PetscInt numCorners = 4; 421 const PetscInt numCells = out->numberoftetrahedra; 422 const PetscInt numVertices = out->numberofpoints; 423 PetscReal *meshCoords = NULL; 424 PetscInt *cells = NULL; 425 PetscBool interpolate = isInterpolated == DMPLEX_INTERPOLATED_FULL ? PETSC_TRUE : PETSC_FALSE; 426 427 if (sizeof (PetscReal) == sizeof (out->pointlist[0])) { 428 meshCoords = (PetscReal *) out->pointlist; 429 } else { 430 PetscInt i; 431 432 PetscCall(PetscMalloc1(dim * numVertices, &meshCoords)); 433 for (i = 0; i < dim * numVertices; ++i) meshCoords[i] = (PetscReal) out->pointlist[i]; 434 } 435 if (sizeof (PetscInt) == sizeof (out->tetrahedronlist[0])) { 436 cells = (PetscInt *) out->tetrahedronlist; 437 } else { 438 PetscInt i; 439 440 PetscCall(PetscMalloc1(numCells * numCorners, &cells)); 441 for (i = 0; i < numCells * numCorners; ++i) cells[i] = (PetscInt) out->tetrahedronlist[i]; 442 } 443 444 PetscCall(DMPlexInvertCells_CTetgen(numCells, numCorners, cells)); 445 PetscCall(DMPlexCreateFromCellListPetsc(comm, dim, numCells, numVertices, numCorners, interpolate, cells, dim, meshCoords, dmRefined)); 446 if (sizeof (PetscReal) != sizeof (out->pointlist[0])) PetscCall(PetscFree(meshCoords)); 447 if (sizeof (PetscInt) != sizeof (out->tetrahedronlist[0])) PetscCall(PetscFree(cells)); 448 449 /* Set labels */ 450 PetscCall(DMUniversalLabelCreateLabels(universal, PETSC_TRUE, *dmRefined)); 451 for (v = 0; v < numVertices; ++v) { 452 if (out->pointmarkerlist[v]) { 453 PetscCall(DMUniversalLabelSetLabelValue(universal, *dmRefined, PETSC_TRUE, v+numCells, out->pointmarkerlist[v])); 454 } 455 } 456 if (interpolate) { 457 PetscInt e, f; 458 459 for (e = 0; e < out->numberofedges; e++) { 460 if (out->edgemarkerlist[e]) { 461 const PetscInt vertices[2] = {out->edgelist[e*2+0]+numCells, out->edgelist[e*2+1]+numCells}; 462 const PetscInt *edges; 463 PetscInt numEdges; 464 465 PetscCall(DMPlexGetJoin(*dmRefined, 2, vertices, &numEdges, &edges)); 466 PetscCheck(numEdges == 1,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Two vertices must cover only one edge, not %D", numEdges); 467 PetscCall(DMUniversalLabelSetLabelValue(universal, *dmRefined, PETSC_TRUE, edges[0], out->edgemarkerlist[e])); 468 PetscCall(DMPlexRestoreJoin(*dmRefined, 2, vertices, &numEdges, &edges)); 469 } 470 } 471 for (f = 0; f < out->numberoftrifaces; f++) { 472 if (out->trifacemarkerlist[f]) { 473 const PetscInt vertices[3] = {out->trifacelist[f*3+0]+numCells, out->trifacelist[f*3+1]+numCells, out->trifacelist[f*3+2]+numCells}; 474 const PetscInt *faces; 475 PetscInt numFaces; 476 477 PetscCall(DMPlexGetFullJoin(*dmRefined, 3, vertices, &numFaces, &faces)); 478 PetscCheck(numFaces == 1,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Three vertices must cover only one face, not %D", numFaces); 479 PetscCall(DMUniversalLabelSetLabelValue(universal, *dmRefined, PETSC_TRUE, faces[0], out->trifacemarkerlist[f])); 480 PetscCall(DMPlexRestoreJoin(*dmRefined, 3, vertices, &numFaces, &faces)); 481 } 482 } 483 } 484 485 #ifdef PETSC_HAVE_EGADS 486 { 487 DMLabel bodyLabel; 488 PetscContainer modelObj; 489 PetscInt cStart, cEnd, c, eStart, eEnd, fStart, fEnd; 490 ego *bodies; 491 ego model, geom; 492 int Nb, oclass, mtype, *senses; 493 494 /* Get Attached EGADS Model from Original DMPlex */ 495 PetscCall(PetscObjectQuery((PetscObject) dm, "EGADS Model", (PetscObject *) &modelObj)); 496 if (modelObj) { 497 PetscCall(PetscContainerGetPointer(modelObj, (void **) &model)); 498 PetscCall(EG_getTopology(model, &geom, &oclass, &mtype, NULL, &Nb, &bodies, &senses)); 499 /* Transfer EGADS Model to Volumetric Mesh */ 500 PetscCall(PetscObjectCompose((PetscObject) *dmRefined, "EGADS Model", (PetscObject) modelObj)); 501 502 /* Set Cell Labels */ 503 PetscCall(DMGetLabel(*dmRefined, "EGADS Body ID", &bodyLabel)); 504 PetscCall(DMPlexGetHeightStratum(*dmRefined, 0, &cStart, &cEnd)); 505 PetscCall(DMPlexGetHeightStratum(*dmRefined, 1, &fStart, &fEnd)); 506 PetscCall(DMPlexGetDepthStratum(*dmRefined, 1, &eStart, &eEnd)); 507 508 for (c = cStart; c < cEnd; ++c) { 509 PetscReal centroid[3] = {0., 0., 0.}; 510 PetscInt b; 511 512 /* Deterimine what body the cell's centroid is located in */ 513 if (!interpolate) { 514 PetscSection coordSection; 515 Vec coordinates; 516 PetscScalar *coords = NULL; 517 PetscInt coordSize, s, d; 518 519 PetscCall(DMGetCoordinatesLocal(*dmRefined, &coordinates)); 520 PetscCall(DMGetCoordinateSection(*dmRefined, &coordSection)); 521 PetscCall(DMPlexVecGetClosure(*dmRefined, coordSection, coordinates, c, &coordSize, &coords)); 522 for (s = 0; s < coordSize; ++s) for (d = 0; d < dim; ++d) centroid[d] += coords[s*dim+d]; 523 PetscCall(DMPlexVecRestoreClosure(*dmRefined, coordSection, coordinates, c, &coordSize, &coords)); 524 } else { 525 PetscCall(DMPlexComputeCellGeometryFVM(*dmRefined, c, NULL, centroid, NULL)); 526 } 527 for (b = 0; b < Nb; ++b) { 528 if (EG_inTopology(bodies[b], centroid) == EGADS_SUCCESS) break; 529 } 530 if (b < Nb) { 531 PetscInt cval = b, eVal, fVal; 532 PetscInt *closure = NULL, Ncl, cl; 533 534 PetscCall(DMLabelSetValue(bodyLabel, c, cval)); 535 PetscCall(DMPlexGetTransitiveClosure(*dmRefined, c, PETSC_TRUE, &Ncl, &closure)); 536 for (cl = 0; cl < Ncl; cl += 2) { 537 const PetscInt p = closure[cl]; 538 539 if (p >= eStart && p < eEnd) { 540 PetscCall(DMLabelGetValue(bodyLabel, p, &eVal)); 541 if (eVal < 0) PetscCall(DMLabelSetValue(bodyLabel, p, cval)); 542 } 543 if (p >= fStart && p < fEnd) { 544 PetscCall(DMLabelGetValue(bodyLabel, p, &fVal)); 545 if (fVal < 0) PetscCall(DMLabelSetValue(bodyLabel, p, cval)); 546 } 547 } 548 PetscCall(DMPlexRestoreTransitiveClosure(*dmRefined, c, PETSC_TRUE, &Ncl, &closure)); 549 } 550 } 551 } 552 } 553 #endif 554 PetscCall(DMPlexSetRefinementUniform(*dmRefined, PETSC_FALSE)); 555 } 556 PetscCall(DMUniversalLabelDestroy(&universal)); 557 PetscCall(PLCDestroy(&in)); 558 PetscCall(PLCDestroy(&out)); 559 PetscFunctionReturn(0); 560 } 561