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