Actual source code: trigenerate.c
1: #include <petsc/private/dmpleximpl.h>
3: #if !defined(ANSI_DECLARATORS)
4: #define ANSI_DECLARATORS
5: #endif
6: #include <triangle.h>
8: static PetscErrorCode InitInput_Triangle(struct triangulateio *inputCtx)
9: {
10: PetscFunctionBegin;
11: inputCtx->numberofpoints = 0;
12: inputCtx->numberofpointattributes = 0;
13: inputCtx->pointlist = NULL;
14: inputCtx->pointattributelist = NULL;
15: inputCtx->pointmarkerlist = NULL;
16: inputCtx->numberofsegments = 0;
17: inputCtx->segmentlist = NULL;
18: inputCtx->segmentmarkerlist = NULL;
19: inputCtx->numberoftriangleattributes = 0;
20: inputCtx->trianglelist = NULL;
21: inputCtx->numberofholes = 0;
22: inputCtx->holelist = NULL;
23: inputCtx->numberofregions = 0;
24: inputCtx->regionlist = NULL;
25: PetscFunctionReturn(PETSC_SUCCESS);
26: }
28: static PetscErrorCode InitOutput_Triangle(struct triangulateio *outputCtx)
29: {
30: PetscFunctionBegin;
31: outputCtx->numberofpoints = 0;
32: outputCtx->pointlist = NULL;
33: outputCtx->pointattributelist = NULL;
34: outputCtx->pointmarkerlist = NULL;
35: outputCtx->numberoftriangles = 0;
36: outputCtx->trianglelist = NULL;
37: outputCtx->triangleattributelist = NULL;
38: outputCtx->neighborlist = NULL;
39: outputCtx->segmentlist = NULL;
40: outputCtx->segmentmarkerlist = NULL;
41: outputCtx->numberofedges = 0;
42: outputCtx->edgelist = NULL;
43: outputCtx->edgemarkerlist = NULL;
44: PetscFunctionReturn(PETSC_SUCCESS);
45: }
47: static PetscErrorCode FiniOutput_Triangle(struct triangulateio *outputCtx)
48: {
49: PetscFunctionBegin;
50: free(outputCtx->pointlist);
51: free(outputCtx->pointmarkerlist);
52: free(outputCtx->segmentlist);
53: free(outputCtx->segmentmarkerlist);
54: free(outputCtx->edgelist);
55: free(outputCtx->edgemarkerlist);
56: free(outputCtx->trianglelist);
57: free(outputCtx->neighborlist);
58: PetscFunctionReturn(PETSC_SUCCESS);
59: }
61: PETSC_EXTERN PetscErrorCode DMPlexGenerate_Triangle(DM boundary, PetscBool interpolate, DM *dm)
62: {
63: MPI_Comm comm;
64: DM_Plex *mesh = (DM_Plex *)boundary->data;
65: PetscInt dim = 2;
66: const PetscBool createConvexHull = PETSC_FALSE;
67: const PetscBool constrained = PETSC_FALSE;
68: const char *labelName = "marker";
69: const char *labelName2 = "Face Sets";
70: struct triangulateio in;
71: struct triangulateio out;
72: DMLabel label, label2;
73: PetscInt vStart, vEnd, v, eStart, eEnd, e;
74: PetscMPIInt rank;
75: PetscBool flg;
76: char opts[64];
78: PetscFunctionBegin;
79: PetscCall(PetscObjectGetComm((PetscObject)boundary, &comm));
80: PetscCallMPI(MPI_Comm_rank(comm, &rank));
81: PetscCall(InitInput_Triangle(&in));
82: PetscCall(InitOutput_Triangle(&out));
83: PetscCall(DMPlexGetDepthStratum(boundary, 0, &vStart, &vEnd));
84: PetscCall(DMGetLabel(boundary, labelName, &label));
85: PetscCall(DMGetLabel(boundary, labelName2, &label2));
86: PetscCall(PetscOptionsGetString(((PetscObject)boundary)->options, ((PetscObject)boundary)->prefix, "-dm_plex_generate_triangle_opts", opts, sizeof(opts), &flg));
87: if (flg) PetscCall(DMPlexTriangleSetOptions(boundary, opts));
89: PetscCall(PetscCIntCast(vEnd - vStart, &in.numberofpoints));
90: if (in.numberofpoints > 0) {
91: PetscSection coordSection;
92: Vec coordinates;
93: PetscScalar *array;
95: PetscCall(PetscMalloc1(in.numberofpoints * dim, &in.pointlist));
96: PetscCall(PetscMalloc1(in.numberofpoints, &in.pointmarkerlist));
97: PetscCall(DMGetCoordinatesLocal(boundary, &coordinates));
98: PetscCall(DMGetCoordinateSection(boundary, &coordSection));
99: PetscCall(VecGetArray(coordinates, &array));
100: for (v = vStart; v < vEnd; ++v) {
101: const PetscInt idx = v - vStart;
102: PetscInt val, off, d;
104: PetscCall(PetscSectionGetOffset(coordSection, v, &off));
105: for (d = 0; d < dim; ++d) in.pointlist[idx * dim + d] = PetscRealPart(array[off + d]);
106: if (label) {
107: PetscCall(DMLabelGetValue(label, v, &val));
108: PetscCall(PetscCIntCast(val, &in.pointmarkerlist[idx]));
109: }
110: }
111: PetscCall(VecRestoreArray(coordinates, &array));
112: }
113: PetscCall(DMPlexGetHeightStratum(boundary, 0, &eStart, &eEnd));
114: PetscCall(PetscCIntCast(eEnd - eStart, &in.numberofsegments));
115: if (in.numberofsegments > 0) {
116: PetscCall(PetscMalloc1(in.numberofsegments * 2, &in.segmentlist));
117: PetscCall(PetscMalloc1(in.numberofsegments, &in.segmentmarkerlist));
118: for (e = eStart; e < eEnd; ++e) {
119: const PetscInt idx = e - eStart;
120: const PetscInt *cone;
121: PetscInt val;
123: PetscCall(DMPlexGetCone(boundary, e, &cone));
125: PetscCall(PetscCIntCast(cone[0] - vStart, &in.segmentlist[idx * 2 + 0]));
126: PetscCall(PetscCIntCast(cone[1] - vStart, &in.segmentlist[idx * 2 + 1]));
128: if (label) {
129: PetscCall(DMLabelGetValue(label, e, &val));
130: PetscCall(PetscCIntCast(val, &in.segmentmarkerlist[idx]));
131: }
132: }
133: }
134: #if 0 /* Do not currently support holes */
135: PetscReal *holeCoords;
136: PetscInt h, d;
138: PetscCall(DMPlexGetHoles(boundary, &in.numberofholes, &holeCords));
139: if (in.numberofholes > 0) {
140: PetscCall(PetscMalloc1(in.numberofholes*dim, &in.holelist));
141: for (h = 0; h < in.numberofholes; ++h) {
142: for (d = 0; d < dim; ++d) in.holelist[h*dim+d] = holeCoords[h*dim+d];
143: }
144: }
145: #endif
146: if (rank == 0) {
147: char args[32];
149: /* Take away 'Q' for verbose output */
150: PetscCall(PetscStrncpy(args, "pqezQ", sizeof(args)));
151: if (createConvexHull) PetscCall(PetscStrlcat(args, "c", sizeof(args)));
152: if (constrained) PetscCall(PetscStrncpy(args, "zepDQ", sizeof(args)));
153: if (mesh->triangleOpts) {
154: triangulate(mesh->triangleOpts, &in, &out, NULL);
155: } else {
156: triangulate(args, &in, &out, NULL);
157: }
158: }
159: PetscCall(PetscFree(in.pointlist));
160: PetscCall(PetscFree(in.pointmarkerlist));
161: PetscCall(PetscFree(in.segmentlist));
162: PetscCall(PetscFree(in.segmentmarkerlist));
163: PetscCall(PetscFree(in.holelist));
165: {
166: DMLabel glabel = NULL;
167: DMLabel glabel2 = NULL;
168: const PetscInt numCorners = 3;
169: const PetscInt numCells = out.numberoftriangles;
170: const PetscInt numVertices = out.numberofpoints;
171: PetscInt *cells;
172: PetscReal *meshCoords;
174: if (sizeof(PetscReal) == sizeof(out.pointlist[0])) {
175: meshCoords = (PetscReal *)out.pointlist;
176: } else {
177: PetscInt i;
179: PetscCall(PetscMalloc1(dim * numVertices, &meshCoords));
180: for (i = 0; i < dim * numVertices; i++) meshCoords[i] = (PetscReal)out.pointlist[i];
181: }
182: if (sizeof(PetscInt) == sizeof(out.trianglelist[0])) {
183: cells = (PetscInt *)out.trianglelist;
184: } else {
185: PetscInt i;
187: PetscCall(PetscMalloc1(numCells * numCorners, &cells));
188: for (i = 0; i < numCells * numCorners; i++) cells[i] = (PetscInt)out.trianglelist[i];
189: }
190: PetscCall(DMPlexCreateFromCellListPetsc(comm, dim, numCells, numVertices, numCorners, interpolate, cells, dim, meshCoords, dm));
191: if (sizeof(PetscReal) != sizeof(out.pointlist[0])) PetscCall(PetscFree(meshCoords));
192: if (sizeof(PetscInt) != sizeof(out.trianglelist[0])) PetscCall(PetscFree(cells));
193: if (label) {
194: PetscCall(DMCreateLabel(*dm, labelName));
195: PetscCall(DMGetLabel(*dm, labelName, &glabel));
196: }
197: if (label2) {
198: PetscCall(DMCreateLabel(*dm, labelName2));
199: PetscCall(DMGetLabel(*dm, labelName2, &glabel2));
200: }
201: /* Set labels */
202: for (v = 0; v < numVertices; ++v) {
203: if (out.pointmarkerlist[v] && glabel) PetscCall(DMLabelSetValue(glabel, v + numCells, out.pointmarkerlist[v]));
204: }
205: if (interpolate) {
206: for (e = 0; e < out.numberofedges; e++) {
207: if (out.edgemarkerlist[e]) {
208: const PetscInt vertices[2] = {out.edgelist[e * 2 + 0] + numCells, out.edgelist[e * 2 + 1] + numCells};
209: const PetscInt *edges;
210: PetscInt numEdges;
212: PetscCall(DMPlexGetJoin(*dm, 2, vertices, &numEdges, &edges));
213: PetscCheck(numEdges == 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Two vertices must cover only one edge, not %" PetscInt_FMT, numEdges);
214: if (glabel) PetscCall(DMLabelSetValue(glabel, edges[0], out.edgemarkerlist[e]));
215: if (glabel2) PetscCall(DMLabelSetValue(glabel2, edges[0], out.edgemarkerlist[e]));
216: PetscCall(DMPlexRestoreJoin(*dm, 2, vertices, &numEdges, &edges));
217: }
218: }
219: }
220: PetscCall(DMPlexSetRefinementUniform(*dm, PETSC_FALSE));
221: }
222: #if 0 /* Do not currently support holes */
223: PetscCall(DMPlexCopyHoles(*dm, boundary));
224: #endif
225: PetscCall(FiniOutput_Triangle(&out));
226: PetscFunctionReturn(PETSC_SUCCESS);
227: }
229: PETSC_EXTERN PetscErrorCode DMPlexRefine_Triangle(DM dm, PetscReal *inmaxVolumes, DM *dmRefined)
230: {
231: MPI_Comm comm;
232: PetscInt dim = 2;
233: const char *labelName = "marker";
234: struct triangulateio in;
235: struct triangulateio out;
236: DMLabel label;
237: PetscInt vStart, vEnd, v, gcStart, cStart, cEnd, c, depth, depthGlobal;
238: PetscMPIInt rank;
239: double *maxVolumes;
241: PetscFunctionBegin;
242: PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
243: PetscCallMPI(MPI_Comm_rank(comm, &rank));
244: PetscCall(InitInput_Triangle(&in));
245: PetscCall(InitOutput_Triangle(&out));
246: PetscCall(DMPlexGetDepth(dm, &depth));
247: PetscCallMPI(MPIU_Allreduce(&depth, &depthGlobal, 1, MPIU_INT, MPI_MAX, comm));
248: PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
249: PetscCall(DMGetLabel(dm, labelName, &label));
251: PetscCall(PetscCIntCast(vEnd - vStart, &in.numberofpoints));
252: if (in.numberofpoints > 0) {
253: PetscSection coordSection;
254: Vec coordinates;
255: PetscScalar *array;
257: PetscCall(PetscMalloc1(in.numberofpoints * dim, &in.pointlist));
258: PetscCall(PetscMalloc1(in.numberofpoints, &in.pointmarkerlist));
259: PetscCall(DMGetCoordinatesLocal(dm, &coordinates));
260: PetscCall(DMGetCoordinateSection(dm, &coordSection));
261: PetscCall(VecGetArray(coordinates, &array));
262: for (v = vStart; v < vEnd; ++v) {
263: const PetscInt idx = v - vStart;
264: PetscInt off, d, val;
266: PetscCall(PetscSectionGetOffset(coordSection, v, &off));
267: for (d = 0; d < dim; ++d) in.pointlist[idx * dim + d] = PetscRealPart(array[off + d]);
268: if (label) {
269: PetscCall(DMLabelGetValue(label, v, &val));
270: PetscCall(PetscCIntCast(val, &in.pointmarkerlist[idx]));
271: }
272: }
273: PetscCall(VecRestoreArray(coordinates, &array));
274: }
275: PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
276: PetscCall(DMPlexGetCellTypeStratum(dm, DM_POLYTOPE_FV_GHOST, &gcStart, NULL));
277: if (gcStart >= 0) cEnd = gcStart;
279: in.numberofcorners = 3;
280: PetscCall(PetscCIntCast(cEnd - cStart, &in.numberoftriangles));
282: #if !defined(PETSC_USE_REAL_DOUBLE)
283: PetscCall(PetscMalloc1(cEnd - cStart, &maxVolumes));
284: for (c = 0; c < cEnd - cStart; ++c) maxVolumes[c] = (double)inmaxVolumes[c];
285: #else
286: maxVolumes = inmaxVolumes;
287: #endif
289: in.trianglearealist = (double *)maxVolumes;
290: if (in.numberoftriangles > 0) {
291: PetscCall(PetscMalloc1(in.numberoftriangles * in.numberofcorners, &in.trianglelist));
292: for (c = cStart; c < cEnd; ++c) {
293: const PetscInt idx = c - cStart;
294: PetscInt *closure = NULL;
295: PetscInt closureSize;
297: PetscCall(DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure));
298: PetscCheck(!(closureSize != 4) || !(closureSize != 7), comm, PETSC_ERR_ARG_WRONG, "Mesh has cell which is not a triangle, %" PetscInt_FMT " vertices in closure", closureSize);
299: for (v = 0; v < 3; ++v) PetscCall(PetscCIntCast(closure[(v + closureSize - 3) * 2] - vStart, &in.trianglelist[idx * in.numberofcorners + v]));
300: PetscCall(DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure));
301: }
302: }
303: /* TODO: Segment markers are missing on input */
304: #if 0 /* Do not currently support holes */
305: PetscReal *holeCoords;
306: PetscInt h, d;
308: PetscCall(DMPlexGetHoles(boundary, &in.numberofholes, &holeCords));
309: if (in.numberofholes > 0) {
310: PetscCall(PetscMalloc1(in.numberofholes*dim, &in.holelist));
311: for (h = 0; h < in.numberofholes; ++h) {
312: for (d = 0; d < dim; ++d) in.holelist[h*dim+d] = holeCoords[h*dim+d];
313: }
314: }
315: #endif
316: if (rank == 0) {
317: char args[32];
319: /* Take away 'Q' for verbose output */
320: PetscCall(PetscStrncpy(args, "pqezQra", sizeof(args)));
321: triangulate(args, &in, &out, NULL);
322: }
323: PetscCall(PetscFree(in.pointlist));
324: PetscCall(PetscFree(in.pointmarkerlist));
325: PetscCall(PetscFree(in.segmentlist));
326: PetscCall(PetscFree(in.segmentmarkerlist));
327: PetscCall(PetscFree(in.trianglelist));
329: {
330: DMLabel rlabel = NULL;
331: const PetscInt numCorners = 3;
332: const PetscInt numCells = out.numberoftriangles;
333: const PetscInt numVertices = out.numberofpoints;
334: PetscInt *cells;
335: PetscReal *meshCoords;
336: PetscBool interpolate = depthGlobal > 1 ? PETSC_TRUE : PETSC_FALSE;
338: if (sizeof(PetscReal) == sizeof(out.pointlist[0])) {
339: meshCoords = (PetscReal *)out.pointlist;
340: } else {
341: PetscInt i;
343: PetscCall(PetscMalloc1(dim * numVertices, &meshCoords));
344: for (i = 0; i < dim * numVertices; i++) meshCoords[i] = (PetscReal)out.pointlist[i];
345: }
346: if (sizeof(PetscInt) == sizeof(out.trianglelist[0])) {
347: cells = (PetscInt *)out.trianglelist;
348: } else {
349: PetscInt i;
351: PetscCall(PetscMalloc1(numCells * numCorners, &cells));
352: for (i = 0; i < numCells * numCorners; i++) cells[i] = (PetscInt)out.trianglelist[i];
353: }
355: PetscCall(DMPlexCreateFromCellListPetsc(comm, dim, numCells, numVertices, numCorners, interpolate, cells, dim, meshCoords, dmRefined));
356: if (label) {
357: PetscCall(DMCreateLabel(*dmRefined, labelName));
358: PetscCall(DMGetLabel(*dmRefined, labelName, &rlabel));
359: }
360: if (sizeof(PetscReal) != sizeof(out.pointlist[0])) PetscCall(PetscFree(meshCoords));
361: if (sizeof(PetscInt) != sizeof(out.trianglelist[0])) PetscCall(PetscFree(cells));
362: /* Set labels */
363: for (v = 0; v < numVertices; ++v) {
364: if (out.pointmarkerlist[v] && rlabel) PetscCall(DMLabelSetValue(rlabel, v + numCells, out.pointmarkerlist[v]));
365: }
366: if (interpolate) {
367: PetscInt e;
369: for (e = 0; e < out.numberofedges; e++) {
370: if (out.edgemarkerlist[e]) {
371: const PetscInt vertices[2] = {out.edgelist[e * 2 + 0] + numCells, out.edgelist[e * 2 + 1] + numCells};
372: const PetscInt *edges;
373: PetscInt numEdges;
375: PetscCall(DMPlexGetJoin(*dmRefined, 2, vertices, &numEdges, &edges));
376: PetscCheck(numEdges == 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Two vertices must cover only one edge, not %" PetscInt_FMT, numEdges);
377: if (rlabel) PetscCall(DMLabelSetValue(rlabel, edges[0], out.edgemarkerlist[e]));
378: PetscCall(DMPlexRestoreJoin(*dmRefined, 2, vertices, &numEdges, &edges));
379: }
380: }
381: }
382: PetscCall(DMPlexSetRefinementUniform(*dmRefined, PETSC_FALSE));
383: }
384: #if 0 /* Do not currently support holes */
385: PetscCall(DMPlexCopyHoles(*dm, boundary));
386: #endif
387: PetscCall(FiniOutput_Triangle(&out));
388: #if !defined(PETSC_USE_REAL_DOUBLE)
389: PetscCall(PetscFree(maxVolumes));
390: #endif
391: PetscFunctionReturn(PETSC_SUCCESS);
392: }