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 (rank == 0) {
135:     char args[32];

137:     /* Take away 'Q' for verbose output */
138:     PetscCall(PetscStrncpy(args, "pqezQ", sizeof(args)));
139:     if (createConvexHull) PetscCall(PetscStrlcat(args, "c", sizeof(args)));
140:     if (constrained) PetscCall(PetscStrncpy(args, "zepDQ", sizeof(args)));
141:     if (mesh->triangleOpts) {
142:       triangulate(mesh->triangleOpts, &in, &out, NULL);
143:     } else {
144:       triangulate(args, &in, &out, NULL);
145:     }
146:   }
147:   PetscCall(PetscFree(in.pointlist));
148:   PetscCall(PetscFree(in.pointmarkerlist));
149:   PetscCall(PetscFree(in.segmentlist));
150:   PetscCall(PetscFree(in.segmentmarkerlist));
151:   PetscCall(PetscFree(in.holelist));

153:   {
154:     DMLabel        glabel      = NULL;
155:     DMLabel        glabel2     = NULL;
156:     const PetscInt numCorners  = 3;
157:     const PetscInt numCells    = out.numberoftriangles;
158:     const PetscInt numVertices = out.numberofpoints;
159:     PetscInt      *cells;
160:     PetscReal     *meshCoords;

162:     if (sizeof(PetscReal) == sizeof(out.pointlist[0])) {
163:       meshCoords = (PetscReal *)out.pointlist;
164:     } else {
165:       PetscInt i;

167:       PetscCall(PetscMalloc1(dim * numVertices, &meshCoords));
168:       for (i = 0; i < dim * numVertices; i++) meshCoords[i] = (PetscReal)out.pointlist[i];
169:     }
170:     if (sizeof(PetscInt) == sizeof(out.trianglelist[0])) {
171:       cells = (PetscInt *)out.trianglelist;
172:     } else {
173:       PetscCall(PetscMalloc1(numCells * numCorners, &cells));
174:       for (PetscInt i = 0; i < numCells * numCorners; i++) cells[i] = (PetscInt)out.trianglelist[i];
175:     }
176:     PetscCall(DMPlexCreateFromCellListPetsc(comm, dim, numCells, numVertices, numCorners, interpolate, cells, dim, meshCoords, dm));
177:     if (sizeof(PetscReal) != sizeof(out.pointlist[0])) PetscCall(PetscFree(meshCoords));
178:     if (sizeof(PetscInt) != sizeof(out.trianglelist[0])) PetscCall(PetscFree(cells));
179:     if (label) {
180:       PetscCall(DMCreateLabel(*dm, labelName));
181:       PetscCall(DMGetLabel(*dm, labelName, &glabel));
182:     }
183:     if (label2) {
184:       PetscCall(DMCreateLabel(*dm, labelName2));
185:       PetscCall(DMGetLabel(*dm, labelName2, &glabel2));
186:     }
187:     /* Set labels */
188:     for (v = 0; v < numVertices; ++v) {
189:       if (out.pointmarkerlist[v] && glabel) PetscCall(DMLabelSetValue(glabel, v + numCells, out.pointmarkerlist[v]));
190:     }
191:     if (interpolate) {
192:       for (e = 0; e < out.numberofedges; e++) {
193:         if (out.edgemarkerlist[e]) {
194:           const PetscInt  vertices[2] = {out.edgelist[e * 2 + 0] + numCells, out.edgelist[e * 2 + 1] + numCells};
195:           const PetscInt *edges;
196:           PetscInt        numEdges;

198:           PetscCall(DMPlexGetJoin(*dm, 2, vertices, &numEdges, &edges));
199:           PetscCheck(numEdges == 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Two vertices must cover only one edge, not %" PetscInt_FMT, numEdges);
200:           if (glabel) PetscCall(DMLabelSetValue(glabel, edges[0], out.edgemarkerlist[e]));
201:           if (glabel2) PetscCall(DMLabelSetValue(glabel2, edges[0], out.edgemarkerlist[e]));
202:           PetscCall(DMPlexRestoreJoin(*dm, 2, vertices, &numEdges, &edges));
203:         }
204:       }
205:     }
206:     PetscCall(DMPlexSetRefinementUniform(*dm, PETSC_FALSE));
207:   }
208: #if 0 /* Do not currently support holes */
209:   PetscCall(DMPlexCopyHoles(*dm, boundary));
210: #endif
211:   PetscCall(FiniOutput_Triangle(&out));
212:   PetscFunctionReturn(PETSC_SUCCESS);
213: }

215: PETSC_EXTERN PetscErrorCode DMPlexRefine_Triangle(DM dm, PetscReal *inmaxVolumes, DM *dmRefined)
216: {
217:   MPI_Comm             comm;
218:   PetscInt             dim       = 2;
219:   const char          *labelName = "marker";
220:   struct triangulateio in;
221:   struct triangulateio out;
222:   DMLabel              label;
223:   PetscInt             vStart, vEnd, v, gcStart, cStart, cEnd, c, depth, depthGlobal;
224:   PetscMPIInt          rank;
225:   double              *maxVolumes;

227:   PetscFunctionBegin;
228:   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
229:   PetscCallMPI(MPI_Comm_rank(comm, &rank));
230:   PetscCall(InitInput_Triangle(&in));
231:   PetscCall(InitOutput_Triangle(&out));
232:   PetscCall(DMPlexGetDepth(dm, &depth));
233:   PetscCallMPI(MPIU_Allreduce(&depth, &depthGlobal, 1, MPIU_INT, MPI_MAX, comm));
234:   PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
235:   PetscCall(DMGetLabel(dm, labelName, &label));

237:   PetscCall(PetscCIntCast(vEnd - vStart, &in.numberofpoints));
238:   if (in.numberofpoints > 0) {
239:     PetscSection coordSection;
240:     Vec          coordinates;
241:     PetscScalar *array;

243:     PetscCall(PetscMalloc1(in.numberofpoints * dim, &in.pointlist));
244:     PetscCall(PetscMalloc1(in.numberofpoints, &in.pointmarkerlist));
245:     PetscCall(DMGetCoordinatesLocal(dm, &coordinates));
246:     PetscCall(DMGetCoordinateSection(dm, &coordSection));
247:     PetscCall(VecGetArray(coordinates, &array));
248:     for (v = vStart; v < vEnd; ++v) {
249:       const PetscInt idx = v - vStart;
250:       PetscInt       off, d, val;

252:       PetscCall(PetscSectionGetOffset(coordSection, v, &off));
253:       for (d = 0; d < dim; ++d) in.pointlist[idx * dim + d] = PetscRealPart(array[off + d]);
254:       if (label) {
255:         PetscCall(DMLabelGetValue(label, v, &val));
256:         PetscCall(PetscCIntCast(val, &in.pointmarkerlist[idx]));
257:       }
258:     }
259:     PetscCall(VecRestoreArray(coordinates, &array));
260:   }
261:   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
262:   PetscCall(DMPlexGetCellTypeStratum(dm, DM_POLYTOPE_FV_GHOST, &gcStart, NULL));
263:   if (gcStart >= 0) cEnd = gcStart;

265:   in.numberofcorners = 3;
266:   PetscCall(PetscCIntCast(cEnd - cStart, &in.numberoftriangles));

268: #if !defined(PETSC_USE_REAL_DOUBLE)
269:   PetscCall(PetscMalloc1(cEnd - cStart, &maxVolumes));
270:   for (c = 0; c < cEnd - cStart; ++c) maxVolumes[c] = (double)inmaxVolumes[c];
271: #else
272:   maxVolumes = inmaxVolumes;
273: #endif

275:   in.trianglearealist = (double *)maxVolumes;
276:   if (in.numberoftriangles > 0) {
277:     PetscCall(PetscMalloc1(in.numberoftriangles * in.numberofcorners, &in.trianglelist));
278:     for (c = cStart; c < cEnd; ++c) {
279:       const PetscInt idx     = c - cStart;
280:       PetscInt      *closure = NULL;
281:       PetscInt       closureSize;

283:       PetscCall(DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure));
284:       PetscCheck(!(closureSize != 4) || !(closureSize != 7), comm, PETSC_ERR_ARG_WRONG, "Mesh has cell which is not a triangle, %" PetscInt_FMT " vertices in closure", closureSize);
285:       for (v = 0; v < 3; ++v) PetscCall(PetscCIntCast(closure[(v + closureSize - 3) * 2] - vStart, &in.trianglelist[idx * in.numberofcorners + v]));
286:       PetscCall(DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure));
287:     }
288:   }
289:   /* TODO: Segment markers are missing on input */
290:   if (rank == 0) {
291:     char args[32];

293:     /* Take away 'Q' for verbose output */
294:     PetscCall(PetscStrncpy(args, "pqezQra", sizeof(args)));
295:     triangulate(args, &in, &out, NULL);
296:   }
297:   PetscCall(PetscFree(in.pointlist));
298:   PetscCall(PetscFree(in.pointmarkerlist));
299:   PetscCall(PetscFree(in.segmentlist));
300:   PetscCall(PetscFree(in.segmentmarkerlist));
301:   PetscCall(PetscFree(in.trianglelist));

303:   {
304:     DMLabel        rlabel      = NULL;
305:     const PetscInt numCorners  = 3;
306:     const PetscInt numCells    = out.numberoftriangles;
307:     const PetscInt numVertices = out.numberofpoints;
308:     PetscInt      *cells;
309:     PetscReal     *meshCoords;
310:     PetscBool      interpolate = depthGlobal > 1 ? PETSC_TRUE : PETSC_FALSE;

312:     if (sizeof(PetscReal) == sizeof(out.pointlist[0])) {
313:       meshCoords = (PetscReal *)out.pointlist;
314:     } else {
315:       PetscCall(PetscMalloc1(dim * numVertices, &meshCoords));
316:       for (PetscInt i = 0; i < dim * numVertices; i++) meshCoords[i] = (PetscReal)out.pointlist[i];
317:     }
318:     if (sizeof(PetscInt) == sizeof(out.trianglelist[0])) {
319:       cells = (PetscInt *)out.trianglelist;
320:     } else {
321:       PetscCall(PetscMalloc1(numCells * numCorners, &cells));
322:       for (PetscInt i = 0; i < numCells * numCorners; i++) cells[i] = (PetscInt)out.trianglelist[i];
323:     }

325:     PetscCall(DMPlexCreateFromCellListPetsc(comm, dim, numCells, numVertices, numCorners, interpolate, cells, dim, meshCoords, dmRefined));
326:     if (label) {
327:       PetscCall(DMCreateLabel(*dmRefined, labelName));
328:       PetscCall(DMGetLabel(*dmRefined, labelName, &rlabel));
329:     }
330:     if (sizeof(PetscReal) != sizeof(out.pointlist[0])) PetscCall(PetscFree(meshCoords));
331:     if (sizeof(PetscInt) != sizeof(out.trianglelist[0])) PetscCall(PetscFree(cells));
332:     /* Set labels */
333:     for (v = 0; v < numVertices; ++v) {
334:       if (out.pointmarkerlist[v] && rlabel) PetscCall(DMLabelSetValue(rlabel, v + numCells, out.pointmarkerlist[v]));
335:     }
336:     if (interpolate) {
337:       for (PetscInt e = 0; e < out.numberofedges; e++) {
338:         if (out.edgemarkerlist[e]) {
339:           const PetscInt  vertices[2] = {out.edgelist[e * 2 + 0] + numCells, out.edgelist[e * 2 + 1] + numCells};
340:           const PetscInt *edges;
341:           PetscInt        numEdges;

343:           PetscCall(DMPlexGetJoin(*dmRefined, 2, vertices, &numEdges, &edges));
344:           PetscCheck(numEdges == 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Two vertices must cover only one edge, not %" PetscInt_FMT, numEdges);
345:           if (rlabel) PetscCall(DMLabelSetValue(rlabel, edges[0], out.edgemarkerlist[e]));
346:           PetscCall(DMPlexRestoreJoin(*dmRefined, 2, vertices, &numEdges, &edges));
347:         }
348:       }
349:     }
350:     PetscCall(DMPlexSetRefinementUniform(*dmRefined, PETSC_FALSE));
351:   }
352: #if 0 /* Do not currently support holes */
353:   PetscCall(DMPlexCopyHoles(*dm, boundary));
354: #endif
355:   PetscCall(FiniOutput_Triangle(&out));
356: #if !defined(PETSC_USE_REAL_DOUBLE)
357:   PetscCall(PetscFree(maxVolumes));
358: #endif
359:   PetscFunctionReturn(PETSC_SUCCESS);
360: }