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