Actual source code: Generator.hh
1: #ifndef included_ALE_Generator_hh
2: #define included_ALE_Generator_hh
4: #ifndef included_ALE_Distribution_hh
5: #include <Distribution.hh>
6: #endif
8: #ifdef PETSC_HAVE_TRIANGLE
9: #include <triangle.h>
10: #endif
11: #ifdef PETSC_HAVE_TETGEN
12: #include <tetgen.h>
13: #endif
15: namespace ALE {
16: #ifdef PETSC_HAVE_TRIANGLE
17: namespace Triangle {
18: class Generator {
19: typedef ALE::Mesh Mesh;
20: public:
21: static void initInput(struct triangulateio *inputCtx) {
22: inputCtx->numberofpoints = 0;
23: inputCtx->numberofpointattributes = 0;
24: inputCtx->pointlist = NULL;
25: inputCtx->pointattributelist = NULL;
26: inputCtx->pointmarkerlist = NULL;
27: inputCtx->numberofsegments = 0;
28: inputCtx->segmentlist = NULL;
29: inputCtx->segmentmarkerlist = NULL;
30: inputCtx->numberoftriangleattributes = 0;
31: inputCtx->trianglelist = NULL;
32: inputCtx->numberofholes = 0;
33: inputCtx->holelist = NULL;
34: inputCtx->numberofregions = 0;
35: inputCtx->regionlist = NULL;
36: };
37: static void initOutput(struct triangulateio *outputCtx) {
38: outputCtx->pointlist = NULL;
39: outputCtx->pointattributelist = NULL;
40: outputCtx->pointmarkerlist = NULL;
41: outputCtx->trianglelist = NULL;
42: outputCtx->triangleattributelist = NULL;
43: outputCtx->neighborlist = NULL;
44: outputCtx->segmentlist = NULL;
45: outputCtx->segmentmarkerlist = NULL;
46: outputCtx->edgelist = NULL;
47: outputCtx->edgemarkerlist = NULL;
48: };
49: static void finiOutput(struct triangulateio *outputCtx) {
50: free(outputCtx->pointmarkerlist);
51: free(outputCtx->edgelist);
52: free(outputCtx->edgemarkerlist);
53: free(outputCtx->trianglelist);
54: free(outputCtx->neighborlist);
55: };
58: static Obj<Mesh> generateMesh(const Obj<Mesh>& boundary, const bool interpolate = false, const bool constrained = false) {
59: int dim = 2;
60: Obj<Mesh> mesh = new Mesh(boundary->comm(), dim, boundary->debug());
61: const Obj<Mesh::sieve_type>& sieve = boundary->getSieve();
62: const bool createConvexHull = false;
63: struct triangulateio in;
64: struct triangulateio out;
65: PetscErrorCode ierr;
67: initInput(&in);
68: initOutput(&out);
69: const Obj<Mesh::label_sequence>& vertices = boundary->depthStratum(0);
70: const Obj<Mesh::label_type>& markers = boundary->getLabel("marker");
71: const Obj<Mesh::real_section_type>& coordinates = boundary->getRealSection("coordinates");
72: const Obj<Mesh::numbering_type>& vNumbering = boundary->getFactory()->getLocalNumbering(boundary, 0);
74: in.numberofpoints = vertices->size();
75: if (in.numberofpoints > 0) {
76: PetscMalloc(in.numberofpoints * dim * sizeof(double), &in.pointlist);
77: PetscMalloc(in.numberofpoints * sizeof(int), &in.pointmarkerlist);
78: for(Mesh::label_sequence::iterator v_iter = vertices->begin(); v_iter != vertices->end(); ++v_iter) {
79: const Mesh::real_section_type::value_type *array = coordinates->restrictPoint(*v_iter);
80: const int idx = vNumbering->getIndex(*v_iter);
82: for(int d = 0; d < dim; d++) {
83: in.pointlist[idx*dim + d] = array[d];
84: }
85: in.pointmarkerlist[idx] = boundary->getValue(markers, *v_iter);
86: }
87: }
88: const Obj<Mesh::label_sequence>& edges = boundary->depthStratum(1);
89: const Obj<Mesh::numbering_type>& eNumbering = boundary->getFactory()->getLocalNumbering(boundary, 1);
91: in.numberofsegments = edges->size();
92: if (in.numberofsegments > 0) {
93: PetscMalloc(in.numberofsegments * 2 * sizeof(int), &in.segmentlist);
94: PetscMalloc(in.numberofsegments * sizeof(int), &in.segmentmarkerlist);
95: for(Mesh::label_sequence::iterator e_iter = edges->begin(); e_iter != edges->end(); ++e_iter) {
96: const Obj<Mesh::sieve_type::traits::coneSequence>& cone = sieve->cone(*e_iter);
97: const int idx = eNumbering->getIndex(*e_iter);
98: int v = 0;
99:
100: for(Mesh::sieve_type::traits::coneSequence::iterator c_iter = cone->begin(); c_iter != cone->end(); ++c_iter) {
101: in.segmentlist[idx*dim + (v++)] = vNumbering->getIndex(*c_iter);
102: }
103: in.segmentmarkerlist[idx] = boundary->getValue(markers, *e_iter);
104: }
105: }
106: const Mesh::holes_type& holes = boundary->getHoles();
108: in.numberofholes = holes.size();
109: if (in.numberofholes > 0) {
110: PetscMalloc(in.numberofholes*dim * sizeof(double), &in.holelist);
111: for(int h = 0; h < in.numberofholes; ++h) {
112: for(int d = 0; d < dim; ++d) {
113: in.holelist[h*dim+d] = holes[h][d];
114: }
115: }
116: }
117: if (mesh->commRank() == 0) {
118: std::string args("pqezQ");
120: if (createConvexHull) {
121: args += "c";
122: }
123: if (constrained) {
124: args = "zepQ";
125: }
126: triangulate((char *) args.c_str(), &in, &out, NULL);
127: }
129: if (in.pointlist) {PetscFree(in.pointlist);}
130: if (in.pointmarkerlist) {PetscFree(in.pointmarkerlist);}
131: if (in.segmentlist) {PetscFree(in.segmentlist);}
132: if (in.segmentmarkerlist) {PetscFree(in.segmentmarkerlist);}
133: if (in.holelist) {PetscFree(in.holelist);}
134: const Obj<Mesh::sieve_type> newSieve = new Mesh::sieve_type(mesh->comm(), mesh->debug());
135: int numCorners = 3;
136: int numCells = out.numberoftriangles;
137: int *cells = out.trianglelist;
138: int numVertices = out.numberofpoints;
139: double *coords = out.pointlist;
141: ALE::SieveBuilder<Mesh>::buildTopology(newSieve, dim, numCells, cells, numVertices, interpolate, numCorners, -1, mesh->getArrowSection("orientation"));
142: mesh->setSieve(newSieve);
143: mesh->stratify();
144: ALE::SieveBuilder<Mesh>::buildCoordinates(mesh, dim, coords);
145: const Obj<Mesh::label_type>& newMarkers = mesh->createLabel("marker");
147: if (mesh->commRank() == 0) {
148: for(int v = 0; v < out.numberofpoints; v++) {
149: if (out.pointmarkerlist[v]) {
150: mesh->setValue(newMarkers, v+out.numberoftriangles, out.pointmarkerlist[v]);
151: }
152: }
153: if (interpolate) {
154: for(int e = 0; e < out.numberofedges; e++) {
155: if (out.edgemarkerlist[e]) {
156: const Mesh::point_type vertexA(out.edgelist[e*2+0]+out.numberoftriangles);
157: const Mesh::point_type vertexB(out.edgelist[e*2+1]+out.numberoftriangles);
158: const Obj<Mesh::sieve_type::supportSet> edge = newSieve->nJoin(vertexA, vertexB, 1);
160: mesh->setValue(newMarkers, *(edge->begin()), out.edgemarkerlist[e]);
161: }
162: }
163: }
164: }
165: mesh->copyHoles(boundary);
166: finiOutput(&out);
167: return mesh;
168: };
169: };
170: class Refiner {
171: public:
172: static Obj<Mesh> refineMesh(const Obj<Mesh>& serialMesh, const double maxVolumes[], const bool interpolate = false) {
173: const int dim = serialMesh->getDimension();
174: const Obj<Mesh> refMesh = new Mesh(serialMesh->comm(), dim, serialMesh->debug());
175: const Obj<Mesh::sieve_type>& serialSieve = serialMesh->getSieve();
176: struct triangulateio in;
177: struct triangulateio out;
178: PetscErrorCode ierr;
180: Generator::initInput(&in);
181: Generator::initOutput(&out);
182: const Obj<Mesh::label_sequence>& vertices = serialMesh->depthStratum(0);
183: const Obj<Mesh::label_type>& markers = serialMesh->getLabel("marker");
184: const Obj<Mesh::real_section_type>& coordinates = serialMesh->getRealSection("coordinates");
185: const Obj<Mesh::numbering_type>& vNumbering = serialMesh->getFactory()->getLocalNumbering(serialMesh, 0);
187: in.numberofpoints = vertices->size();
188: if (in.numberofpoints > 0) {
189: PetscMalloc(in.numberofpoints * dim * sizeof(double), &in.pointlist);
190: PetscMalloc(in.numberofpoints * sizeof(int), &in.pointmarkerlist);
191: for(Mesh::label_sequence::iterator v_iter = vertices->begin(); v_iter != vertices->end(); ++v_iter) {
192: const Mesh::real_section_type::value_type *array = coordinates->restrictPoint(*v_iter);
193: const int idx = vNumbering->getIndex(*v_iter);
195: for(int d = 0; d < dim; d++) {
196: in.pointlist[idx*dim + d] = array[d];
197: }
198: in.pointmarkerlist[idx] = serialMesh->getValue(markers, *v_iter);
199: }
200: }
201: const Obj<Mesh::label_sequence>& faces = serialMesh->heightStratum(0);
202: const Obj<Mesh::numbering_type>& fNumbering = serialMesh->getFactory()->getLocalNumbering(serialMesh, serialMesh->depth());
204: in.numberofcorners = 3;
205: in.numberoftriangles = faces->size();
206: in.trianglearealist = (double *) maxVolumes;
207: if (in.numberoftriangles > 0) {
208: PetscMalloc(in.numberoftriangles*in.numberofcorners * sizeof(int), &in.trianglelist);
209: if (serialMesh->depth() == 1) {
210: for(Mesh::label_sequence::iterator f_iter = faces->begin(); f_iter != faces->end(); ++f_iter) {
211: const Obj<Mesh::sieve_type::traits::coneSequence>& cone = serialSieve->cone(*f_iter);
212: const int idx = fNumbering->getIndex(*f_iter);
213: int v = 0;
215: for(Mesh::sieve_type::traits::coneSequence::iterator c_iter = cone->begin(); c_iter != cone->end(); ++c_iter) {
216: in.trianglelist[idx*in.numberofcorners + v++] = vNumbering->getIndex(*c_iter);
217: }
218: }
219: } else if (serialMesh->depth() == 2) {
220: for(Mesh::label_sequence::iterator f_iter = faces->begin(); f_iter != faces->end(); ++f_iter) {
221: typedef ALE::SieveAlg<Mesh> sieve_alg_type;
222: const Obj<sieve_alg_type::coneArray>& cone = sieve_alg_type::nCone(serialMesh, *f_iter, 2);
223: const int idx = fNumbering->getIndex(*f_iter);
224: int v = 0;
226: for(Mesh::sieve_type::coneArray::iterator c_iter = cone->begin(); c_iter != cone->end(); ++c_iter) {
227: in.trianglelist[idx*in.numberofcorners + v++] = vNumbering->getIndex(*c_iter);
228: }
229: }
230: } else {
231: throw ALE::Exception("Invalid sieve: Cannot gives sieves of arbitrary depth to Triangle");
232: }
233: }
234: if (serialMesh->depth() == 2) {
235: const Obj<Mesh::label_sequence>& edges = serialMesh->depthStratum(1);
236: #define NEW_LABEL
237: #ifdef NEW_LABEL
238: for(Mesh::label_sequence::iterator e_iter = edges->begin(); e_iter != edges->end(); ++e_iter) {
239: if (serialMesh->getValue(markers, *e_iter)) {
240: in.numberofsegments++;
241: }
242: }
243: std::cout << "Number of segments: " << in.numberofsegments << std::endl;
244: if (in.numberofsegments > 0) {
245: int s = 0;
247: PetscMalloc(in.numberofsegments * 2 * sizeof(int), &in.segmentlist);
248: PetscMalloc(in.numberofsegments * sizeof(int), &in.segmentmarkerlist);
249: for(Mesh::label_sequence::iterator e_iter = edges->begin(); e_iter != edges->end(); ++e_iter) {
250: const int edgeMarker = serialMesh->getValue(markers, *e_iter);
252: if (edgeMarker) {
253: const Obj<Mesh::sieve_type::traits::coneSequence>& cone = serialSieve->cone(*e_iter);
254: int p = 0;
256: for(Mesh::sieve_type::traits::coneSequence::iterator v_iter = cone->begin(); v_iter != cone->end(); ++v_iter) {
257: in.segmentlist[s*2 + (p++)] = vNumbering->getIndex(*v_iter);
258: }
259: in.segmentmarkerlist[s++] = edgeMarker;
260: }
261: }
262: }
263: #else
264: const Obj<Mesh::label_type::baseSequence>& boundary = markers->base();
266: in.numberofsegments = 0;
267: for(Mesh::label_type::baseSequence::iterator b_iter = boundary->begin(); b_iter != boundary->end(); ++b_iter) {
268: for(Mesh::label_sequence::iterator e_iter = edges->begin(); e_iter != edges->end(); ++e_iter) {
269: if (*b_iter == *e_iter) {
270: in.numberofsegments++;
271: }
272: }
273: }
274: if (in.numberofsegments > 0) {
275: int s = 0;
277: PetscMalloc(in.numberofsegments * 2 * sizeof(int), &in.segmentlist);
278: PetscMalloc(in.numberofsegments * sizeof(int), &in.segmentmarkerlist);
279: for(Mesh::label_type::baseSequence::iterator b_iter = boundary->begin(); b_iter != boundary->end(); ++b_iter) {
280: for(Mesh::label_sequence::iterator e_iter = edges->begin(); e_iter != edges->end(); ++e_iter) {
281: if (*b_iter == *e_iter) {
282: const Obj<Mesh::sieve_type::traits::coneSequence>& cone = serialSieve->cone(*e_iter);
283: int p = 0;
285: for(Mesh::sieve_type::traits::coneSequence::iterator v_iter = cone->begin(); v_iter != cone->end(); ++v_iter) {
286: in.segmentlist[s*2 + (p++)] = vNumbering->getIndex(*v_iter);
287: }
288: in.segmentmarkerlist[s++] = serialMesh->getValue(markers, *e_iter);
289: }
290: }
291: }
292: }
293: #endif
294: }
296: in.numberofholes = 0;
297: if (in.numberofholes > 0) {
298: PetscMalloc(in.numberofholes * dim * sizeof(int), &in.holelist);
299: }
300: if (serialMesh->commRank() == 0) {
301: std::string args("pqezQra");
303: triangulate((char *) args.c_str(), &in, &out, NULL);
304: }
305: if (in.pointlist) {PetscFree(in.pointlist);}
306: if (in.pointmarkerlist) {PetscFree(in.pointmarkerlist);}
307: if (in.segmentlist) {PetscFree(in.segmentlist);}
308: if (in.segmentmarkerlist) {PetscFree(in.segmentmarkerlist);}
309: if (in.trianglelist) {PetscFree(in.trianglelist);}
310: const Obj<Mesh::sieve_type> newSieve = new Mesh::sieve_type(serialMesh->comm(), serialMesh->debug());
311: int numCorners = 3;
312: int numCells = out.numberoftriangles;
313: int *cells = out.trianglelist;
314: int numVertices = out.numberofpoints;
315: double *coords = out.pointlist;
317: ALE::SieveBuilder<Mesh>::buildTopology(newSieve, dim, numCells, cells, numVertices, interpolate, numCorners, -1, refMesh->getArrowSection("orientation"));
318: refMesh->setSieve(newSieve);
319: refMesh->stratify();
320: ALE::SieveBuilder<Mesh>::buildCoordinates(refMesh, dim, coords);
321: const Obj<Mesh::label_type>& newMarkers = refMesh->createLabel("marker");
323: if (refMesh->commRank() == 0) {
324: for(int v = 0; v < out.numberofpoints; v++) {
325: if (out.pointmarkerlist[v]) {
326: refMesh->setValue(newMarkers, v+out.numberoftriangles, out.pointmarkerlist[v]);
327: }
328: }
329: if (interpolate) {
330: for(int e = 0; e < out.numberofedges; e++) {
331: if (out.edgemarkerlist[e]) {
332: const Mesh::point_type vertexA(out.edgelist[e*2+0]+out.numberoftriangles);
333: const Mesh::point_type vertexB(out.edgelist[e*2+1]+out.numberoftriangles);
334: const Obj<Mesh::sieve_type::supportSet> edge = newSieve->nJoin(vertexA, vertexB, 1);
336: refMesh->setValue(newMarkers, *(edge->begin()), out.edgemarkerlist[e]);
337: }
338: }
339: }
340: }
342: Generator::finiOutput(&out);
343: if (refMesh->commSize() > 1) {
344: return ALE::Distribution<Mesh>::distributeMesh(refMesh);
345: }
346: return refMesh;
347: };
348: static Obj<Mesh> refineMesh(const Obj<Mesh>& mesh, const Obj<Mesh::real_section_type>& maxVolumes, const bool interpolate = false) {
349: Obj<Mesh> serialMesh = ALE::Distribution<Mesh>::unifyMesh(mesh);
350: const Obj<Mesh::real_section_type> serialMaxVolumes = ALE::Distribution<Mesh>::distributeSection(maxVolumes, serialMesh, serialMesh->getDistSendOverlap(), serialMesh->getDistRecvOverlap());
352: return refineMesh(serialMesh, serialMaxVolumes->restrict(), interpolate);
353: };
354: static Obj<Mesh> refineMesh(const Obj<Mesh>& mesh, const double maxVolume, const bool interpolate = false) {
355: Obj<Mesh> serialMesh;
356: if (mesh->commSize() > 1) {
357: serialMesh = ALE::Distribution<Mesh>::unifyMesh(mesh);
358: } else {
359: serialMesh = mesh;
360: }
361: const int numFaces = serialMesh->heightStratum(0)->size();
362: double *serialMaxVolumes = new double[numFaces];
364: for(int f = 0; f < numFaces; f++) {
365: serialMaxVolumes[f] = maxVolume;
366: }
367: const Obj<Mesh> refMesh = refineMesh(serialMesh, serialMaxVolumes, interpolate);
368: delete [] serialMaxVolumes;
369: return refMesh;
370: };
371: static Obj<Mesh> refineMeshLocal(const Obj<Mesh>& mesh, const double maxVolumes[], const bool interpolate = false) {
372: const int dim = mesh->getDimension();
373: const Obj<Mesh> refMesh = new Mesh(mesh->comm(), dim, mesh->debug());
374: const Obj<Mesh::sieve_type>& sieve = mesh->getSieve();
375: struct triangulateio in;
376: struct triangulateio out;
377: PetscErrorCode ierr;
379: Generator::initInput(&in);
380: Generator::initOutput(&out);
381: const Obj<Mesh::label_sequence>& vertices = mesh->depthStratum(0);
382: const Obj<Mesh::label_type>& markers = mesh->getLabel("marker");
383: const Obj<Mesh::real_section_type>& coordinates = mesh->getRealSection("coordinates");
384: const Obj<Mesh::numbering_type>& vNumbering = mesh->getFactory()->getLocalNumbering(mesh, 0);
386: in.numberofpoints = vertices->size();
387: if (in.numberofpoints > 0) {
388: PetscMalloc(in.numberofpoints * dim * sizeof(double), &in.pointlist);
389: PetscMalloc(in.numberofpoints * sizeof(int), &in.pointmarkerlist);
390: for(Mesh::label_sequence::iterator v_iter = vertices->begin(); v_iter != vertices->end(); ++v_iter) {
391: const Mesh::real_section_type::value_type *array = coordinates->restrictPoint(*v_iter);
392: const int idx = vNumbering->getIndex(*v_iter);
394: for(int d = 0; d < dim; d++) {
395: in.pointlist[idx*dim + d] = array[d];
396: }
397: in.pointmarkerlist[idx] = mesh->getValue(markers, *v_iter);
398: }
399: }
400: const Obj<Mesh::label_sequence>& faces = mesh->heightStratum(0);
401: const Obj<Mesh::numbering_type>& fNumbering = mesh->getFactory()->getLocalNumbering(mesh, mesh->depth());
403: in.numberofcorners = 3;
404: in.numberoftriangles = faces->size();
405: in.trianglearealist = (double *) maxVolumes;
406: if (in.numberoftriangles > 0) {
407: PetscMalloc(in.numberoftriangles*in.numberofcorners * sizeof(int), &in.trianglelist);
408: if (mesh->depth() == 1) {
409: for(Mesh::label_sequence::iterator f_iter = faces->begin(); f_iter != faces->end(); ++f_iter) {
410: const Obj<Mesh::sieve_type::traits::coneSequence>& cone = sieve->cone(*f_iter);
411: const int idx = fNumbering->getIndex(*f_iter);
412: int v = 0;
414: for(Mesh::sieve_type::traits::coneSequence::iterator c_iter = cone->begin(); c_iter != cone->end(); ++c_iter) {
415: in.trianglelist[idx*in.numberofcorners + v++] = vNumbering->getIndex(*c_iter);
416: }
417: }
418: } else if (mesh->depth() == 2) {
419: for(Mesh::label_sequence::iterator f_iter = faces->begin(); f_iter != faces->end(); ++f_iter) {
420: typedef ALE::SieveAlg<Mesh> sieve_alg_type;
421: const Obj<sieve_alg_type::coneArray>& cone = sieve_alg_type::nCone(mesh, *f_iter, 2);
422: const int idx = fNumbering->getIndex(*f_iter);
423: int v = 0;
425: for(Mesh::sieve_type::coneArray::iterator c_iter = cone->begin(); c_iter != cone->end(); ++c_iter) {
426: in.trianglelist[idx*in.numberofcorners + v++] = vNumbering->getIndex(*c_iter);
427: }
428: }
429: } else {
430: throw ALE::Exception("Invalid sieve: Cannot gives sieves of arbitrary depth to Triangle");
431: }
432: }
433: if (mesh->depth() == 2) {
434: const Obj<Mesh::label_sequence>& edges = mesh->depthStratum(1);
435: for(Mesh::label_sequence::iterator e_iter = edges->begin(); e_iter != edges->end(); ++e_iter) {
436: if (mesh->getValue(markers, *e_iter)) {
437: in.numberofsegments++;
438: }
439: }
440: std::cout << "Number of segments: " << in.numberofsegments << std::endl;
441: if (in.numberofsegments > 0) {
442: int s = 0;
444: PetscMalloc(in.numberofsegments * 2 * sizeof(int), &in.segmentlist);
445: PetscMalloc(in.numberofsegments * sizeof(int), &in.segmentmarkerlist);
446: for(Mesh::label_sequence::iterator e_iter = edges->begin(); e_iter != edges->end(); ++e_iter) {
447: const int edgeMarker = mesh->getValue(markers, *e_iter);
449: if (edgeMarker) {
450: const Obj<Mesh::sieve_type::traits::coneSequence>& cone = sieve->cone(*e_iter);
451: int p = 0;
453: for(Mesh::sieve_type::traits::coneSequence::iterator v_iter = cone->begin(); v_iter != cone->end(); ++v_iter) {
454: in.segmentlist[s*2 + (p++)] = vNumbering->getIndex(*v_iter);
455: }
456: in.segmentmarkerlist[s++] = edgeMarker;
457: }
458: }
459: }
460: }
462: in.numberofholes = 0;
463: if (in.numberofholes > 0) {
464: PetscMalloc(in.numberofholes * dim * sizeof(int), &in.holelist);
465: }
466: std::string args("pqezQra");
468: triangulate((char *) args.c_str(), &in, &out, NULL);
469: if (in.pointlist) {PetscFree(in.pointlist);}
470: if (in.pointmarkerlist) {PetscFree(in.pointmarkerlist);}
471: if (in.segmentlist) {PetscFree(in.segmentlist);}
472: if (in.segmentmarkerlist) {PetscFree(in.segmentmarkerlist);}
473: if (in.trianglelist) {PetscFree(in.trianglelist);}
474: const Obj<Mesh::sieve_type> newSieve = new Mesh::sieve_type(mesh->comm(), mesh->debug());
475: int numCorners = 3;
476: int numCells = out.numberoftriangles;
477: int *cells = out.trianglelist;
478: int numVertices = out.numberofpoints;
479: double *coords = out.pointlist;
481: ALE::SieveBuilder<Mesh>::buildTopologyMultiple(newSieve, dim, numCells, cells, numVertices, interpolate, numCorners, -1, refMesh->getArrowSection("orientation"));
482: refMesh->setSieve(newSieve);
483: refMesh->stratify();
484: ALE::SieveBuilder<Mesh>::buildCoordinatesMultiple(refMesh, dim, coords);
485: const Obj<Mesh::label_type>& newMarkers = refMesh->createLabel("marker");
487: for(int v = 0; v < out.numberofpoints; v++) {
488: if (out.pointmarkerlist[v]) {
489: refMesh->setValue(newMarkers, v+out.numberoftriangles, out.pointmarkerlist[v]);
490: }
491: }
492: if (interpolate) {
493: for(int e = 0; e < out.numberofedges; e++) {
494: if (out.edgemarkerlist[e]) {
495: const Mesh::point_type vertexA(out.edgelist[e*2+0]+out.numberoftriangles);
496: const Mesh::point_type vertexB(out.edgelist[e*2+1]+out.numberoftriangles);
497: const Obj<Mesh::sieve_type::supportSet> edge = newSieve->nJoin(vertexA, vertexB, 1);
499: refMesh->setValue(newMarkers, *(edge->begin()), out.edgemarkerlist[e]);
500: }
501: }
502: }
504: Generator::finiOutput(&out);
505: return refMesh;
506: };
507: static Obj<Mesh> refineMeshLocal(const Obj<Mesh>& mesh, const double maxVolume, const bool interpolate = false) {
508: const int numLocalFaces = mesh->heightStratum(0)->size();
509: double *localMaxVolumes = new double[numLocalFaces];
511: for(int f = 0; f < numLocalFaces; f++) {
512: localMaxVolumes[f] = maxVolume;
513: }
514: const Obj<Mesh> refMesh = refineMeshLocal(mesh, localMaxVolumes, interpolate);
515: const Obj<Mesh::sieve_type> refSieve = refMesh->getSieve();
516: delete [] localMaxVolumes;
517: #if 0
518: typedef typename ALE::New::Completion<Mesh, typename Mesh::sieve_type::point_type> sieveCompletion;
519: // This is where we enforce consistency over the overlap
520: // We need somehow to update the overlap to account for the new stuff
521: //
522: // 1) Since we are refining only, the vertices are invariant
523: // 2) We need to make another label for the interprocess boundaries so
524: // that Triangle will respect them
525: // 3) We then throw all that label into the new overlap
526: //
527: // Alternative: Figure out explicitly which segments were refined, and then
528: // communicated the refinement over the old overlap. Use this info to locally
529: // construct the new overlap and flip to get a decent mesh
530: sieveCompletion::scatterCones(refSieve, refSieve, reMesh->getDistSendOverlap(), refMesh->getDistRecvOverlap(), refMesh);
531: #endif
532: return refMesh;
533: };
534: };
535: };
536: #endif
537: #ifdef PETSC_HAVE_TETGEN
538: namespace TetGen {
539: class Generator {
540: public:
541: static Obj<Mesh> generateMesh(const Obj<Mesh>& boundary, const bool interpolate = false, const bool constrained = false) {
542: typedef ALE::SieveAlg<Mesh> sieve_alg_type;
543: const int dim = 3;
544: Obj<Mesh> mesh = new Mesh(boundary->comm(), dim, boundary->debug());
545: const PetscMPIInt rank = mesh->commRank();
546: bool createConvexHull = false;
547: ::tetgenio in;
548: ::tetgenio out;
550: const Obj<Mesh::label_sequence>& vertices = boundary->depthStratum(0);
551: const Obj<Mesh::numbering_type>& vNumbering = boundary->getFactory()->getLocalNumbering(boundary, 0);
552: const Obj<Mesh::real_section_type>& coordinates = boundary->getRealSection("coordinates");
553: const Obj<Mesh::label_type>& markers = boundary->getLabel("marker");
556: in.numberofpoints = vertices->size();
557: if (in.numberofpoints > 0) {
558: in.pointlist = new double[in.numberofpoints*dim];
559: in.pointmarkerlist = new int[in.numberofpoints];
560: for(Mesh::label_sequence::iterator v_iter = vertices->begin(); v_iter != vertices->end(); ++v_iter) {
561: const Mesh::real_section_type::value_type *array = coordinates->restrictPoint(*v_iter);
562: const int idx = vNumbering->getIndex(*v_iter);
564: for(int d = 0; d < dim; d++) {
565: in.pointlist[idx*dim + d] = array[d];
566: }
567: in.pointmarkerlist[idx] = boundary->getValue(markers, *v_iter);
568: }
569: }
570:
571: if (boundary->depth() != 0) { //our boundary mesh COULD be just a pointset; in which case depth = height = 0;
572: const Obj<Mesh::label_sequence>& facets = boundary->depthStratum(boundary->depth());
573: //PetscPrintf(boundary->comm(), "%d facets on the boundary\n", facets->size());
574: const Obj<Mesh::numbering_type>& fNumbering = boundary->getFactory()->getLocalNumbering(boundary, boundary->depth());
576: in.numberoffacets = facets->size();
577: if (in.numberoffacets > 0) {
578: in.facetlist = new tetgenio::facet[in.numberoffacets];
579: in.facetmarkerlist = new int[in.numberoffacets];
580: for(Mesh::label_sequence::iterator f_iter = facets->begin(); f_iter != facets->end(); ++f_iter) {
581: const Obj<sieve_alg_type::coneArray>& cone = sieve_alg_type::nCone(boundary, *f_iter, boundary->depth());
582: const int idx = fNumbering->getIndex(*f_iter);
583:
584: in.facetlist[idx].numberofpolygons = 1;
585: in.facetlist[idx].polygonlist = new tetgenio::polygon[in.facetlist[idx].numberofpolygons];
586: in.facetlist[idx].numberofholes = 0;
587: in.facetlist[idx].holelist = NULL;
588:
589: tetgenio::polygon *poly = in.facetlist[idx].polygonlist;
590: int c = 0;
592: poly->numberofvertices = cone->size();
593: poly->vertexlist = new int[poly->numberofvertices];
594: for(sieve_alg_type::coneArray::iterator c_iter = cone->begin(); c_iter != cone->end(); ++c_iter) {
595: const int vIdx = vNumbering->getIndex(*c_iter);
596:
597: poly->vertexlist[c++] = vIdx;
598: }
599: in.facetmarkerlist[idx] = boundary->getValue(markers, *f_iter);
600: }
601: }
602: }else {
603: createConvexHull = true;
604: }
606: in.numberofholes = 0;
607: if (rank == 0) {
608: // Normal operation
609: std::string args("pqezQ");
610: //constrained operation
611: if (constrained) {
612: args = "pezQ";
613: if (createConvexHull) {
614: args = "ezQ";
615: //PetscPrintf(boundary->comm(), "createConvexHull\n");
616: }
617: }
618: // Just make tetrahedrons
619: // std::string args("efzV");
620: // Adds a center point
621: // std::string args("pqezQi");
622: // in.numberofaddpoints = 1;
623: // in.addpointlist = new double[in.numberofaddpoints*dim];
624: // in.addpointlist[0] = 0.5;
625: // in.addpointlist[1] = 0.5;
626: // in.addpointlist[2] = 0.5;
628: //if (createConvexHull) args += "c"; NOT SURE, but this was basically unused before, and the convex hull should be filled in if "p" isn't included.
629: ::tetrahedralize((char *) args.c_str(), &in, &out);
630: }
631: const Obj<Mesh::sieve_type> newSieve = new Mesh::sieve_type(mesh->comm(), mesh->debug());
632: int numCorners = 4;
633: int numCells = out.numberoftetrahedra;
634: int *cells = out.tetrahedronlist;
635: int numVertices = out.numberofpoints;
636: double *coords = out.pointlist;
638: if (!interpolate) {
639: for(int c = 0; c < numCells; ++c) {
640: int tmp = cells[c*4+0];
641: cells[c*4+0] = cells[c*4+1];
642: cells[c*4+1] = tmp;
643: }
644: }
645: ALE::SieveBuilder<Mesh>::buildTopology(newSieve, dim, numCells, cells, numVertices, interpolate, numCorners, -1, mesh->getArrowSection("orientation"));
646: mesh->setSieve(newSieve);
647: mesh->stratify();
648: ALE::SieveBuilder<Mesh>::buildCoordinates(mesh, dim, coords);
649: const Obj<Mesh::label_type>& newMarkers = mesh->createLabel("marker");
650:
651: for(int v = 0; v < out.numberofpoints; v++) {
652: if (out.pointmarkerlist[v]) {
653: mesh->setValue(newMarkers, v+out.numberoftetrahedra, out.pointmarkerlist[v]);
654: }
655: }
656: if (interpolate) {
657: if (out.edgemarkerlist) {
658: for(int e = 0; e < out.numberofedges; e++) {
659: if (out.edgemarkerlist[e]) {
660: Mesh::point_type endpointA(out.edgelist[e*2+0]+out.numberoftetrahedra);
661: Mesh::point_type endpointB(out.edgelist[e*2+1]+out.numberoftetrahedra);
662: Obj<Mesh::sieve_type::supportSet> edge = newSieve->nJoin(endpointA, endpointB, 1);
664: mesh->setValue(newMarkers, *edge->begin(), out.edgemarkerlist[e]);
665: }
666: }
667: }
668: if (out.trifacemarkerlist) {
669: // Work around TetGen bug for raw tetrahedralization
670: // The boundary faces are 0,1,4,5,8,9,11,12,13,15,16,17
671: // for(int f = 0; f < out.numberoftrifaces; f++) {
672: // if (out.trifacemarkerlist[f]) {
673: // out.trifacemarkerlist[f] = 0;
674: // } else {
675: // out.trifacemarkerlist[f] = 1;
676: // }
677: // }
678: for(int f = 0; f < out.numberoftrifaces; f++) {
679: if (out.trifacemarkerlist[f]) {
680: Mesh::point_type cornerA(out.trifacelist[f*3+0]+out.numberoftetrahedra);
681: Mesh::point_type cornerB(out.trifacelist[f*3+1]+out.numberoftetrahedra);
682: Mesh::point_type cornerC(out.trifacelist[f*3+2]+out.numberoftetrahedra);
683: Obj<Mesh::sieve_type::supportSet> corners = Mesh::sieve_type::supportSet();
684: Obj<Mesh::sieve_type::supportSet> edges = Mesh::sieve_type::supportSet();
685: corners->clear();corners->insert(cornerA);corners->insert(cornerB);
686: edges->insert(*newSieve->nJoin1(corners)->begin());
687: corners->clear();corners->insert(cornerB);corners->insert(cornerC);
688: edges->insert(*newSieve->nJoin1(corners)->begin());
689: corners->clear();corners->insert(cornerC);corners->insert(cornerA);
690: edges->insert(*newSieve->nJoin1(corners)->begin());
691: const Mesh::point_type face = *newSieve->nJoin1(edges)->begin();
692: const int faceMarker = out.trifacemarkerlist[f];
693: const Obj<Mesh::coneArray> closure = sieve_alg_type::closure(mesh, face);
694: const Mesh::coneArray::iterator end = closure->end();
696: for(Mesh::coneArray::iterator cl_iter = closure->begin(); cl_iter != end; ++cl_iter) {
697: mesh->setValue(newMarkers, *cl_iter, faceMarker);
698: }
699: }
700: }
701: }
702: }
703: return mesh;
704: };
705: };
706: class Refiner {
707: public:
708: static Obj<Mesh> refineMesh(const Obj<Mesh>& serialMesh, const double maxVolumes[], const bool interpolate = false) {
709: typedef ALE::SieveAlg<Mesh> sieve_alg_type;
710: const int dim = serialMesh->getDimension();
711: const int depth = serialMesh->depth();
712: const Obj<Mesh> refMesh = new Mesh(serialMesh->comm(), dim, serialMesh->debug());
713: ::tetgenio in;
714: ::tetgenio out;
716: const Obj<Mesh::label_sequence>& vertices = serialMesh->depthStratum(0);
717: const Obj<Mesh::label_type>& markers = serialMesh->getLabel("marker");
718: const Obj<Mesh::real_section_type>& coordinates = serialMesh->getRealSection("coordinates");
719: const Obj<Mesh::numbering_type>& vNumbering = serialMesh->getFactory()->getLocalNumbering(serialMesh, 0);
721: in.numberofpoints = vertices->size();
722: if (in.numberofpoints > 0) {
723: in.pointlist = new double[in.numberofpoints*dim];
724: in.pointmarkerlist = new int[in.numberofpoints];
725: for(Mesh::label_sequence::iterator v_iter = vertices->begin(); v_iter != vertices->end(); ++v_iter) {
726: const Mesh::real_section_type::value_type *array = coordinates->restrictPoint(*v_iter);
727: const int idx = vNumbering->getIndex(*v_iter);
729: for(int d = 0; d < dim; d++) {
730: in.pointlist[idx*dim + d] = array[d];
731: }
732: in.pointmarkerlist[idx] = serialMesh->getValue(markers, *v_iter);
733: }
734: }
735: const Obj<Mesh::label_sequence>& cells = serialMesh->heightStratum(0);
736: const Obj<Mesh::numbering_type>& cNumbering = serialMesh->getFactory()->getLocalNumbering(serialMesh, depth);
738: in.numberofcorners = 4;
739: in.numberoftetrahedra = cells->size();
740: in.tetrahedronvolumelist = (double *) maxVolumes;
741: if (in.numberoftetrahedra > 0) {
742: in.tetrahedronlist = new int[in.numberoftetrahedra*in.numberofcorners];
743: for(Mesh::label_sequence::iterator c_iter = cells->begin(); c_iter != cells->end(); ++c_iter) {
744: typedef ALE::SieveAlg<Mesh> sieve_alg_type;
745: const Obj<sieve_alg_type::coneArray>& cone = sieve_alg_type::nCone(serialMesh, *c_iter, depth);
746: const int idx = cNumbering->getIndex(*c_iter);
747: int v = 0;
749: for(Mesh::sieve_type::coneArray::iterator v_iter = cone->begin(); v_iter != cone->end(); ++v_iter) {
750: in.tetrahedronlist[idx*in.numberofcorners + v++] = vNumbering->getIndex(*v_iter);
751: }
752: }
753: }
754: if (serialMesh->depth() == 3) {
755: const Obj<Mesh::label_sequence>& boundary = serialMesh->getLabelStratum("marker", 1);
757: in.numberoftrifaces = 0;
758: for(Mesh::label_sequence::iterator b_iter = boundary->begin(); b_iter != boundary->end(); ++b_iter) {
759: if (serialMesh->height(*b_iter) == 1) {
760: in.numberoftrifaces++;
761: }
762: }
763: if (in.numberoftrifaces > 0) {
764: int f = 0;
766: in.trifacelist = new int[in.numberoftrifaces*3];
767: in.trifacemarkerlist = new int[in.numberoftrifaces];
768: for(Mesh::label_sequence::iterator b_iter = boundary->begin(); b_iter != boundary->end(); ++b_iter) {
769: if (serialMesh->height(*b_iter) == 1) {
770: const Obj<Mesh::coneArray>& cone = sieve_alg_type::nCone(serialMesh, *b_iter, 2);
771: int p = 0;
773: for(Mesh::coneArray::iterator v_iter = cone->begin(); v_iter != cone->end(); ++v_iter) {
774: in.trifacelist[f*3 + (p++)] = vNumbering->getIndex(*v_iter);
775: }
776: in.trifacemarkerlist[f++] = serialMesh->getValue(markers, *b_iter);
777: }
778: }
779: }
780: }
782: in.numberofholes = 0;
783: if (serialMesh->commRank() == 0) {
784: std::string args("qezQra");
786: ::tetrahedralize((char *) args.c_str(), &in, &out);
787: }
788: in.tetrahedronvolumelist = NULL;
789: const Obj<Mesh::sieve_type> newSieve = new Mesh::sieve_type(refMesh->comm(), refMesh->debug());
790: int numCorners = 4;
791: int numCells = out.numberoftetrahedra;
792: int *newCells = out.tetrahedronlist;
793: int numVertices = out.numberofpoints;
794: double *coords = out.pointlist;
796: if (!interpolate) {
797: for(int c = 0; c < numCells; ++c) {
798: int tmp = newCells[c*4+0];
799: newCells[c*4+0] = newCells[c*4+1];
800: newCells[c*4+1] = tmp;
801: }
802: }
803: ALE::SieveBuilder<Mesh>::buildTopology(newSieve, dim, numCells, newCells, numVertices, interpolate, numCorners, -1, refMesh->getArrowSection("orientation"));
804: refMesh->setSieve(newSieve);
805: refMesh->stratify();
806: ALE::SieveBuilder<Mesh>::buildCoordinates(refMesh, dim, coords);
807: const Obj<Mesh::label_type>& newMarkers = refMesh->createLabel("marker");
810: for(int v = 0; v < out.numberofpoints; v++) {
811: if (out.pointmarkerlist[v]) {
812: refMesh->setValue(newMarkers, v+out.numberoftetrahedra, out.pointmarkerlist[v]);
813: }
814: }
815: if (interpolate) {
816: if (out.edgemarkerlist) {
817: for(int e = 0; e < out.numberofedges; e++) {
818: if (out.edgemarkerlist[e]) {
819: Mesh::point_type endpointA(out.edgelist[e*2+0]+out.numberoftetrahedra);
820: Mesh::point_type endpointB(out.edgelist[e*2+1]+out.numberoftetrahedra);
821: Obj<Mesh::sieve_type::supportSet> edge = newSieve->nJoin(endpointA, endpointB, 1);
823: refMesh->setValue(newMarkers, *edge->begin(), out.edgemarkerlist[e]);
824: }
825: }
826: }
827: if (out.trifacemarkerlist) {
828: for(int f = 0; f < out.numberoftrifaces; f++) {
829: if (out.trifacemarkerlist[f]) {
830: Mesh::point_type cornerA(out.trifacelist[f*3+0]+out.numberoftetrahedra);
831: Mesh::point_type cornerB(out.trifacelist[f*3+1]+out.numberoftetrahedra);
832: Mesh::point_type cornerC(out.trifacelist[f*3+2]+out.numberoftetrahedra);
833: Obj<Mesh::sieve_type::supportSet> corners = Mesh::sieve_type::supportSet();
834: Obj<Mesh::sieve_type::supportSet> edges = Mesh::sieve_type::supportSet();
835: corners->clear();corners->insert(cornerA);corners->insert(cornerB);
836: edges->insert(*newSieve->nJoin1(corners)->begin());
837: corners->clear();corners->insert(cornerB);corners->insert(cornerC);
838: edges->insert(*newSieve->nJoin1(corners)->begin());
839: corners->clear();corners->insert(cornerC);corners->insert(cornerA);
840: edges->insert(*newSieve->nJoin1(corners)->begin());
841: const Mesh::point_type face = *newSieve->nJoin1(edges)->begin();
842: const int faceMarker = out.trifacemarkerlist[f];
843: const Obj<Mesh::coneArray> closure = sieve_alg_type::closure(refMesh, face);
844: const Mesh::coneArray::iterator end = closure->end();
846: for(Mesh::coneArray::iterator cl_iter = closure->begin(); cl_iter != end; ++cl_iter) {
847: refMesh->setValue(newMarkers, *cl_iter, faceMarker);
848: }
849: }
850: }
851: }
852: }
853: if (refMesh->commSize() > 1) {
854: return ALE::Distribution<Mesh>::distributeMesh(refMesh);
855: }
856: return refMesh;
857: };
858: static Obj<Mesh> refineMesh(const Obj<Mesh>& mesh, const Obj<Mesh::real_section_type>& maxVolumes, const bool interpolate = false) {
859: Obj<Mesh> serialMesh = ALE::Distribution<Mesh>::unifyMesh(mesh);
860: const Obj<Mesh::real_section_type> serialMaxVolumes = ALE::Distribution<Mesh>::distributeSection(maxVolumes, serialMesh, serialMesh->getDistSendOverlap(), serialMesh->getDistRecvOverlap());
862: return refineMesh(serialMesh, serialMaxVolumes->restrict(), interpolate);
863: };
864: static Obj<Mesh> refineMesh(const Obj<Mesh>& mesh, const double maxVolume, const bool interpolate = false) {
865: Obj<Mesh> serialMesh;
866: if (mesh->commSize() > 1) {
867: serialMesh = ALE::Distribution<Mesh>::unifyMesh(mesh);
868: } else {
869: serialMesh = mesh;
870: }
871: const int numCells = serialMesh->heightStratum(0)->size();
872: double *serialMaxVolumes = new double[numCells];
874: for(int c = 0; c < numCells; c++) {
875: serialMaxVolumes[c] = maxVolume;
876: }
877: const Obj<Mesh> refMesh = refineMesh(serialMesh, serialMaxVolumes, interpolate);
878: delete [] serialMaxVolumes;
879: return refMesh;
880: };
881: };
882: };
883: #endif
884: class Generator {
885: typedef ALE::Mesh Mesh;
886: public:
887: static Obj<Mesh> generateMesh(const Obj<Mesh>& boundary, const bool interpolate = false, const bool constrained = false) {
888: int dim = boundary->getDimension();
890: if (dim == 1) {
891: #ifdef PETSC_HAVE_TRIANGLE
892: return ALE::Triangle::Generator::generateMesh(boundary, interpolate, constrained);
893: #else
894: throw ALE::Exception("Mesh generation currently requires Triangle to be installed. Use --download-triangle during configure.");
895: #endif
896: } else if (dim == 2) {
897: #ifdef PETSC_HAVE_TETGEN
898: return ALE::TetGen::Generator::generateMesh(boundary, interpolate, constrained);
899: #else
900: throw ALE::Exception("Mesh generation currently requires TetGen to be installed. Use --download-tetgen during configure.");
901: #endif
902: }
903: return NULL;
904: };
905: static Obj<Mesh> refineMesh(const Obj<Mesh>& mesh, const Obj<Mesh::real_section_type>& maxVolumes, const bool interpolate = false) {
906: int dim = mesh->getDimension();
908: if (dim == 2) {
909: #ifdef PETSC_HAVE_TRIANGLE
910: return ALE::Triangle::Refiner::refineMesh(mesh, maxVolumes, interpolate);
911: #else
912: throw ALE::Exception("Mesh refinement currently requires Triangle to be installed. Use --download-triangle during configure.");
913: #endif
914: } else if (dim == 3) {
915: #ifdef PETSC_HAVE_TETGEN
916: return ALE::TetGen::Refiner::refineMesh(mesh, maxVolumes, interpolate);
917: #else
918: throw ALE::Exception("Mesh refinement currently requires TetGen to be installed. Use --download-tetgen during configure.");
919: #endif
920: }
921: return NULL;
922: };
923: static Obj<Mesh> refineMesh(const Obj<Mesh>& mesh, const double maxVolume, const bool interpolate = false) {
924: int dim = mesh->getDimension();
926: if (dim == 2) {
927: #ifdef PETSC_HAVE_TRIANGLE
928: return ALE::Triangle::Refiner::refineMesh(mesh, maxVolume, interpolate);
929: #else
930: throw ALE::Exception("Mesh refinement currently requires Triangle to be installed. Use --download-triangle during configure.");
931: #endif
932: } else if (dim == 3) {
933: #ifdef PETSC_HAVE_TETGEN
934: return ALE::TetGen::Refiner::refineMesh(mesh, maxVolume, interpolate);
935: #else
936: throw ALE::Exception("Mesh refinement currently requires TetGen to be installed. Use --download-tetgen during configure.");
937: #endif
938: }
939: return NULL;
940: };
941: static Obj<Mesh> refineMeshLocal(const Obj<Mesh>& mesh, const double maxVolume, const bool interpolate = false) {
942: int dim = mesh->getDimension();
944: if (dim == 2) {
945: #ifdef PETSC_HAVE_TRIANGLE
946: return ALE::Triangle::Refiner::refineMeshLocal(mesh, maxVolume, interpolate);
947: #else
948: throw ALE::Exception("Mesh refinement currently requires Triangle to be installed. Use --download-triangle during configure.");
949: #endif
950: } else if (dim == 3) {
951: #ifdef PETSC_HAVE_TETGEN
952: return ALE::TetGen::Refiner::refineMesh(mesh, maxVolume, interpolate);
953: #else
954: throw ALE::Exception("Mesh refinement currently requires TetGen to be installed. Use --download-tetgen during configure.");
955: #endif
956: }
957: return NULL;
958: };
959: };
960: }
962: #endif