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) {
 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:         }

107:         in.numberofholes = 0;
108:         if (in.numberofholes > 0) {
109:           PetscMalloc(in.numberofholes*dim * sizeof(int), &in.holelist);
110:         }
111:         if (mesh->commRank() == 0) {
112:           std::string args("pqezQ");

114:           if (createConvexHull) {
115:             args += "c";
116:           }
117:           triangulate((char *) args.c_str(), &in, &out, NULL);
118:         }

120:         if (in.pointlist)         {PetscFree(in.pointlist);}
121:         if (in.pointmarkerlist)   {PetscFree(in.pointmarkerlist);}
122:         if (in.segmentlist)       {PetscFree(in.segmentlist);}
123:         if (in.segmentmarkerlist) {PetscFree(in.segmentmarkerlist);}
124:         const Obj<Mesh::sieve_type> newSieve = new Mesh::sieve_type(mesh->comm(), mesh->debug());
125:         int     numCorners  = 3;
126:         int     numCells    = out.numberoftriangles;
127:         int    *cells       = out.trianglelist;
128:         int     numVertices = out.numberofpoints;
129:         double *coords      = out.pointlist;

131:         ALE::SieveBuilder<Mesh>::buildTopology(newSieve, dim, numCells, cells, numVertices, interpolate, numCorners, -1, mesh->getArrowSection("orientation"));
132:         mesh->setSieve(newSieve);
133:         mesh->stratify();
134:         ALE::SieveBuilder<Mesh>::buildCoordinates(mesh, dim, coords);
135:         const Obj<Mesh::label_type>& newMarkers = mesh->createLabel("marker");

137:         if (mesh->commRank() == 0) {
138:           for(int v = 0; v < out.numberofpoints; v++) {
139:             if (out.pointmarkerlist[v]) {
140:               mesh->setValue(newMarkers, v+out.numberoftriangles, out.pointmarkerlist[v]);
141:             }
142:           }
143:           if (interpolate) {
144:             for(int e = 0; e < out.numberofedges; e++) {
145:               if (out.edgemarkerlist[e]) {
146:                 const Mesh::point_type vertexA(out.edgelist[e*2+0]+out.numberoftriangles);
147:                 const Mesh::point_type vertexB(out.edgelist[e*2+1]+out.numberoftriangles);
148:                 const Obj<Mesh::sieve_type::supportSet> edge = newSieve->nJoin(vertexA, vertexB, 1);

150:                 mesh->setValue(newMarkers, *(edge->begin()), out.edgemarkerlist[e]);
151:               }
152:             }
153:           }
154:         }
155:         finiOutput(&out);
156:         return mesh;
157:       };
158:     };
159:     class Refiner {
160:     public:
161:       static Obj<Mesh> refineMesh(const Obj<Mesh>& serialMesh, const double maxVolumes[], const bool interpolate = false) {
162:         const int                    dim         = serialMesh->getDimension();
163:         const Obj<Mesh>              refMesh     = new Mesh(serialMesh->comm(), dim, serialMesh->debug());
164:         const Obj<Mesh::sieve_type>& serialSieve = serialMesh->getSieve();
165:         struct triangulateio in;
166:         struct triangulateio out;
167:         PetscErrorCode       ierr;

169:         Generator::initInput(&in);
170:         Generator::initOutput(&out);
171:         const Obj<Mesh::label_sequence>&    vertices    = serialMesh->depthStratum(0);
172:         const Obj<Mesh::label_type>&        markers     = serialMesh->getLabel("marker");
173:         const Obj<Mesh::real_section_type>& coordinates = serialMesh->getRealSection("coordinates");
174:         const Obj<Mesh::numbering_type>&    vNumbering  = serialMesh->getFactory()->getLocalNumbering(serialMesh, 0);

176:         in.numberofpoints = vertices->size();
177:         if (in.numberofpoints > 0) {
178:           PetscMalloc(in.numberofpoints * dim * sizeof(double), &in.pointlist);
179:           PetscMalloc(in.numberofpoints * sizeof(int), &in.pointmarkerlist);
180:           for(Mesh::label_sequence::iterator v_iter = vertices->begin(); v_iter != vertices->end(); ++v_iter) {
181:             const Mesh::real_section_type::value_type *array = coordinates->restrictPoint(*v_iter);
182:             const int                                  idx   = vNumbering->getIndex(*v_iter);

184:             for(int d = 0; d < dim; d++) {
185:               in.pointlist[idx*dim + d] = array[d];
186:             }
187:             in.pointmarkerlist[idx] = serialMesh->getValue(markers, *v_iter);
188:           }
189:         }
190:         const Obj<Mesh::label_sequence>& faces      = serialMesh->heightStratum(0);
191:         const Obj<Mesh::numbering_type>& fNumbering = serialMesh->getFactory()->getLocalNumbering(serialMesh, serialMesh->depth());

193:         in.numberofcorners   = 3;
194:         in.numberoftriangles = faces->size();
195:         in.trianglearealist  = (double *) maxVolumes;
196:         if (in.numberoftriangles > 0) {
197:           PetscMalloc(in.numberoftriangles*in.numberofcorners * sizeof(int), &in.trianglelist);
198:           if (serialMesh->depth() == 1) {
199:             for(Mesh::label_sequence::iterator f_iter = faces->begin(); f_iter != faces->end(); ++f_iter) {
200:               const Obj<Mesh::sieve_type::traits::coneSequence>& cone = serialSieve->cone(*f_iter);
201:               const int                                          idx  = fNumbering->getIndex(*f_iter);
202:               int                                                v    = 0;

204:               for(Mesh::sieve_type::traits::coneSequence::iterator c_iter = cone->begin(); c_iter != cone->end(); ++c_iter) {
205:                 in.trianglelist[idx*in.numberofcorners + v++] = vNumbering->getIndex(*c_iter);
206:               }
207:             }
208:           } else if (serialMesh->depth() == 2) {
209:             for(Mesh::label_sequence::iterator f_iter = faces->begin(); f_iter != faces->end(); ++f_iter) {
210:               typedef ALE::SieveAlg<Mesh> sieve_alg_type;
211:               const Obj<sieve_alg_type::coneArray>& cone = sieve_alg_type::nCone(serialMesh, *f_iter, 2);
212:               const int                             idx  = fNumbering->getIndex(*f_iter);
213:               int                                   v    = 0;

215:               for(Mesh::sieve_type::coneArray::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 {
220:             throw ALE::Exception("Invalid sieve: Cannot gives sieves of arbitrary depth to Triangle");
221:           }
222:         }
223:         if (serialMesh->depth() == 2) {
224:           const Obj<Mesh::label_sequence>&           edges    = serialMesh->depthStratum(1);
225:           const Obj<Mesh::label_type::baseSequence>& boundary = markers->base();

227:           in.numberofsegments = 0;
228:           for(Mesh::label_type::baseSequence::iterator b_iter = boundary->begin(); b_iter != boundary->end(); ++b_iter) {
229:             for(Mesh::label_sequence::iterator e_iter = edges->begin(); e_iter != edges->end(); ++e_iter) {
230:               if (*b_iter == *e_iter) {
231:                 in.numberofsegments++;
232:               }
233:             }
234:           }
235:           if (in.numberofsegments > 0) {
236:             int s = 0;

238:             PetscMalloc(in.numberofsegments * 2 * sizeof(int), &in.segmentlist);
239:             PetscMalloc(in.numberofsegments * sizeof(int), &in.segmentmarkerlist);
240:             for(Mesh::label_type::baseSequence::iterator b_iter = boundary->begin(); b_iter != boundary->end(); ++b_iter) {
241:               for(Mesh::label_sequence::iterator e_iter = edges->begin(); e_iter != edges->end(); ++e_iter) {
242:                 if (*b_iter == *e_iter) {
243:                   const Obj<Mesh::sieve_type::traits::coneSequence>& cone = serialSieve->cone(*e_iter);
244:                   int                                                            p    = 0;

246:                   for(Mesh::sieve_type::traits::coneSequence::iterator v_iter = cone->begin(); v_iter != cone->end(); ++v_iter) {
247:                     in.segmentlist[s*2 + (p++)] = vNumbering->getIndex(*v_iter);
248:                   }
249:                   in.segmentmarkerlist[s++] = serialMesh->getValue(markers, *e_iter);
250:                 }
251:               }
252:             }
253:           }
254:         }

256:         in.numberofholes = 0;
257:         if (in.numberofholes > 0) {
258:           PetscMalloc(in.numberofholes * dim * sizeof(int), &in.holelist);
259:         }
260:         if (serialMesh->commRank() == 0) {
261:           std::string args("pqezQra");

263:           triangulate((char *) args.c_str(), &in, &out, NULL);
264:         }
265:         if (in.pointlist)         {PetscFree(in.pointlist);}
266:         if (in.pointmarkerlist)   {PetscFree(in.pointmarkerlist);}
267:         if (in.segmentlist)       {PetscFree(in.segmentlist);}
268:         if (in.segmentmarkerlist) {PetscFree(in.segmentmarkerlist);}
269:         if (in.trianglelist)      {PetscFree(in.trianglelist);}
270:         const Obj<Mesh::sieve_type> newSieve = new Mesh::sieve_type(serialMesh->comm(), serialMesh->debug());
271:         int     numCorners  = 3;
272:         int     numCells    = out.numberoftriangles;
273:         int    *cells       = out.trianglelist;
274:         int     numVertices = out.numberofpoints;
275:         double *coords      = out.pointlist;

277:         ALE::SieveBuilder<Mesh>::buildTopology(newSieve, dim, numCells, cells, numVertices, interpolate, numCorners, -1, refMesh->getArrowSection("orientation"));
278:         refMesh->setSieve(newSieve);
279:         refMesh->stratify();
280:         ALE::SieveBuilder<Mesh>::buildCoordinates(refMesh, dim, coords);
281:         const Obj<Mesh::label_type>& newMarkers = refMesh->createLabel("marker");

283:         if (refMesh->commRank() == 0) {
284:           for(int v = 0; v < out.numberofpoints; v++) {
285:             if (out.pointmarkerlist[v]) {
286:               refMesh->setValue(newMarkers, v+out.numberoftriangles, out.pointmarkerlist[v]);
287:             }
288:           }
289:           if (interpolate) {
290:             for(int e = 0; e < out.numberofedges; e++) {
291:               if (out.edgemarkerlist[e]) {
292:                 const Mesh::point_type vertexA(out.edgelist[e*2+0]+out.numberoftriangles);
293:                 const Mesh::point_type vertexB(out.edgelist[e*2+1]+out.numberoftriangles);
294:                 const Obj<Mesh::sieve_type::supportSet> edge = newSieve->nJoin(vertexA, vertexB, 1);

296:                 refMesh->setValue(newMarkers, *(edge->begin()), out.edgemarkerlist[e]);
297:               }
298:             }
299:           }
300:         }

302:         Generator::finiOutput(&out);
303:         if (refMesh->commSize() > 1) {
304:           return ALE::Distribution<Mesh>::distributeMesh(refMesh);
305:         }
306:         return refMesh;
307:       };
308:       static Obj<Mesh> refineMesh(const Obj<Mesh>& mesh, const Obj<Mesh::real_section_type>& maxVolumes, const bool interpolate = false) {
309:         Obj<Mesh>                          serialMesh       = ALE::Distribution<Mesh>::unifyMesh(mesh);
310:         const Obj<Mesh::real_section_type> serialMaxVolumes = ALE::Distribution<Mesh>::distributeSection(maxVolumes, serialMesh, serialMesh->getDistSendOverlap(), serialMesh->getDistRecvOverlap());

312:         return refineMesh(serialMesh, serialMaxVolumes->restrict(), interpolate);
313:       };
314:       static Obj<Mesh> refineMesh(const Obj<Mesh>& mesh, const double maxVolume, const bool interpolate = false) {
315:         Obj<Mesh> serialMesh;
316:         if (mesh->commSize() > 1) {
317:           serialMesh = ALE::Distribution<Mesh>::unifyMesh(mesh);
318:         } else {
319:           serialMesh = mesh;
320:         }
321:         const int numFaces         = serialMesh->heightStratum(0)->size();
322:         double   *serialMaxVolumes = new double[numFaces];

324:         for(int f = 0; f < numFaces; f++) {
325:           serialMaxVolumes[f] = maxVolume;
326:         }
327:         const Obj<Mesh> refMesh = refineMesh(serialMesh, serialMaxVolumes, interpolate);
328:         delete [] serialMaxVolumes;
329:         return refMesh;
330:       };
331:     };
332:   };
333: #endif
334: #ifdef PETSC_HAVE_TETGEN
335:   namespace TetGen {
336:     class Generator {
337:     public:
338:       static Obj<Mesh> generateMesh(const Obj<Mesh>& boundary, const bool interpolate = false) {
339:         typedef ALE::SieveAlg<Mesh> sieve_alg_type;
340:         const int         dim   = 3;
341:         Obj<Mesh>         mesh  = new Mesh(boundary->comm(), dim, boundary->debug());
342:         const PetscMPIInt rank  = mesh->commRank();
343:         const bool        createConvexHull = false;
344:         ::tetgenio        in;
345:         ::tetgenio        out;

347:         const Obj<Mesh::label_sequence>&    vertices    = boundary->depthStratum(0);
348:         const Obj<Mesh::numbering_type>&    vNumbering  = boundary->getFactory()->getLocalNumbering(boundary, 0);
349:         const Obj<Mesh::real_section_type>& coordinates = boundary->getRealSection("coordinates");
350:         const Obj<Mesh::label_type>&        markers     = boundary->getLabel("marker");


353:         in.numberofpoints = vertices->size();
354:         if (in.numberofpoints > 0) {
355:           in.pointlist       = new double[in.numberofpoints*dim];
356:           in.pointmarkerlist = new int[in.numberofpoints];
357:           for(Mesh::label_sequence::iterator v_iter = vertices->begin(); v_iter != vertices->end(); ++v_iter) {
358:             const Mesh::real_section_type::value_type *array = coordinates->restrictPoint(*v_iter);
359:             const int                                  idx   = vNumbering->getIndex(*v_iter);

361:             for(int d = 0; d < dim; d++) {
362:               in.pointlist[idx*dim + d] = array[d];
363:             }
364:             in.pointmarkerlist[idx] = boundary->getValue(markers, *v_iter);
365:           }
366:         }

368:         const Obj<Mesh::label_sequence>& facets     = boundary->depthStratum(boundary->depth());
369:         const Obj<Mesh::numbering_type>& fNumbering = boundary->getFactory()->getLocalNumbering(boundary, boundary->depth());

371:         in.numberoffacets = facets->size();
372:         if (in.numberoffacets > 0) {
373:           in.facetlist       = new tetgenio::facet[in.numberoffacets];
374:           in.facetmarkerlist = new int[in.numberoffacets];
375:           for(Mesh::label_sequence::iterator f_iter = facets->begin(); f_iter != facets->end(); ++f_iter) {
376:             const Obj<sieve_alg_type::coneArray>& cone = sieve_alg_type::nCone(boundary, *f_iter, boundary->depth());
377:             const int                             idx  = fNumbering->getIndex(*f_iter);

379:             in.facetlist[idx].numberofpolygons = 1;
380:             in.facetlist[idx].polygonlist      = new tetgenio::polygon[in.facetlist[idx].numberofpolygons];
381:             in.facetlist[idx].numberofholes    = 0;
382:             in.facetlist[idx].holelist         = NULL;

384:             tetgenio::polygon *poly = in.facetlist[idx].polygonlist;
385:             int                c    = 0;

387:             poly->numberofvertices = cone->size();
388:             poly->vertexlist       = new int[poly->numberofvertices];
389:             for(sieve_alg_type::coneArray::iterator c_iter = cone->begin(); c_iter != cone->end(); ++c_iter) {
390:               const int vIdx = vNumbering->getIndex(*c_iter);

392:               poly->vertexlist[c++] = vIdx;
393:             }
394:             in.facetmarkerlist[idx] = boundary->getValue(markers, *f_iter);
395:           }
396:         }

398:         in.numberofholes = 0;
399:         if (rank == 0) {
400:           // Normal operation
401:           std::string args("pqezQ");
402:           // Just make tetrahedrons
403: //           std::string args("efzV");
404:           // Adds a center point
405: //           std::string args("pqezQi");
406: //           in.numberofaddpoints = 1;
407: //           in.addpointlist      = new double[in.numberofaddpoints*dim];
408: //           in.addpointlist[0]   = 0.5;
409: //           in.addpointlist[1]   = 0.5;
410: //           in.addpointlist[2]   = 0.5;

412:           if (createConvexHull) args += "c";
413:           ::tetrahedralize((char *) args.c_str(), &in, &out);
414:         }
415:         const Obj<Mesh::sieve_type> newSieve = new Mesh::sieve_type(mesh->comm(), mesh->debug());
416:         int     numCorners  = 4;
417:         int     numCells    = out.numberoftetrahedra;
418:         int    *cells       = out.tetrahedronlist;
419:         int     numVertices = out.numberofpoints;
420:         double *coords      = out.pointlist;

422:         ALE::SieveBuilder<Mesh>::buildTopology(newSieve, dim, numCells, cells, numVertices, interpolate, numCorners, -1, mesh->getArrowSection("orientation"));
423:         mesh->setSieve(newSieve);
424:         mesh->stratify();
425:         ALE::SieveBuilder<Mesh>::buildCoordinates(mesh, dim, coords);
426:         const Obj<Mesh::label_type>& newMarkers = mesh->createLabel("marker");
427: 
428:         for(int v = 0; v < out.numberofpoints; v++) {
429:           if (out.pointmarkerlist[v]) {
430:             mesh->setValue(newMarkers, v+out.numberoftetrahedra, out.pointmarkerlist[v]);
431:           }
432:         }
433:         if (interpolate) {
434:           if (out.edgemarkerlist) {
435:             for(int e = 0; e < out.numberofedges; e++) {
436:               if (out.edgemarkerlist[e]) {
437:                 Mesh::point_type endpointA(out.edgelist[e*2+0]+out.numberoftetrahedra);
438:                 Mesh::point_type endpointB(out.edgelist[e*2+1]+out.numberoftetrahedra);
439:                 Obj<Mesh::sieve_type::supportSet> edge = newSieve->nJoin(endpointA, endpointB, 1);

441:                 mesh->setValue(newMarkers, *edge->begin(), out.edgemarkerlist[e]);
442:               }
443:             }
444:           }
445:           if (out.trifacemarkerlist) {
446:             // Work around TetGen bug for raw tetrahedralization
447:             //   The boundary faces are 0,1,4,5,8,9,11,12,13,15,16,17
448: //             for(int f = 0; f < out.numberoftrifaces; f++) {
449: //               if (out.trifacemarkerlist[f]) {
450: //                 out.trifacemarkerlist[f] = 0;
451: //               } else {
452: //                 out.trifacemarkerlist[f] = 1;
453: //               }
454: //             }
455:             for(int f = 0; f < out.numberoftrifaces; f++) {
456:               if (out.trifacemarkerlist[f]) {
457:                 Mesh::point_type cornerA(out.trifacelist[f*3+0]+out.numberoftetrahedra);
458:                 Mesh::point_type cornerB(out.trifacelist[f*3+1]+out.numberoftetrahedra);
459:                 Mesh::point_type cornerC(out.trifacelist[f*3+2]+out.numberoftetrahedra);
460:                 Obj<Mesh::sieve_type::supportSet> corners = Mesh::sieve_type::supportSet();
461:                 Obj<Mesh::sieve_type::supportSet> edges   = Mesh::sieve_type::supportSet();
462:                 corners->clear();corners->insert(cornerA);corners->insert(cornerB);
463:                 edges->insert(*newSieve->nJoin1(corners)->begin());
464:                 corners->clear();corners->insert(cornerB);corners->insert(cornerC);
465:                 edges->insert(*newSieve->nJoin1(corners)->begin());
466:                 corners->clear();corners->insert(cornerC);corners->insert(cornerA);
467:                 edges->insert(*newSieve->nJoin1(corners)->begin());
468:                 const Mesh::point_type          face       = *newSieve->nJoin1(edges)->begin();
469:                 const int                       faceMarker = out.trifacemarkerlist[f];
470:                 const Obj<Mesh::coneArray>      closure    = sieve_alg_type::closure(mesh, face);
471:                 const Mesh::coneArray::iterator end        = closure->end();

473:                 for(Mesh::coneArray::iterator cl_iter = closure->begin(); cl_iter != end; ++cl_iter) {
474:                   mesh->setValue(newMarkers, *cl_iter, faceMarker);
475:                 }
476:               }
477:             }
478:           }
479:         }
480:         return mesh;
481:       };
482:     };
483:     class Refiner {
484:     public:
485:       static Obj<Mesh> refineMesh(const Obj<Mesh>& serialMesh, const double maxVolumes[], const bool interpolate = false) {
486:         typedef ALE::SieveAlg<Mesh> sieve_alg_type;
487:         const int       dim     = serialMesh->getDimension();
488:         const int       depth   = serialMesh->depth();
489:         const Obj<Mesh> refMesh = new Mesh(serialMesh->comm(), dim, serialMesh->debug());
490:         ::tetgenio      in;
491:         ::tetgenio      out;

493:         const Obj<Mesh::label_sequence>&    vertices    = serialMesh->depthStratum(0);
494:         const Obj<Mesh::label_type>&        markers     = serialMesh->getLabel("marker");
495:         const Obj<Mesh::real_section_type>& coordinates = serialMesh->getRealSection("coordinates");
496:         const Obj<Mesh::numbering_type>&    vNumbering  = serialMesh->getFactory()->getLocalNumbering(serialMesh, 0);

498:         in.numberofpoints = vertices->size();
499:         if (in.numberofpoints > 0) {
500:           in.pointlist       = new double[in.numberofpoints*dim];
501:           in.pointmarkerlist = new int[in.numberofpoints];
502:           for(Mesh::label_sequence::iterator v_iter = vertices->begin(); v_iter != vertices->end(); ++v_iter) {
503:             const Mesh::real_section_type::value_type *array = coordinates->restrictPoint(*v_iter);
504:             const int                                  idx   = vNumbering->getIndex(*v_iter);

506:             for(int d = 0; d < dim; d++) {
507:               in.pointlist[idx*dim + d] = array[d];
508:             }
509:             in.pointmarkerlist[idx] = serialMesh->getValue(markers, *v_iter);
510:           }
511:         }
512:         const Obj<Mesh::label_sequence>& cells      = serialMesh->heightStratum(0);
513:         const Obj<Mesh::numbering_type>& cNumbering = serialMesh->getFactory()->getLocalNumbering(serialMesh, depth);

515:         in.numberofcorners       = 4;
516:         in.numberoftetrahedra    = cells->size();
517:         in.tetrahedronvolumelist = (double *) maxVolumes;
518:         if (in.numberoftetrahedra > 0) {
519:           in.tetrahedronlist     = new int[in.numberoftetrahedra*in.numberofcorners];
520:           for(Mesh::label_sequence::iterator c_iter = cells->begin(); c_iter != cells->end(); ++c_iter) {
521:             typedef ALE::SieveAlg<Mesh> sieve_alg_type;
522:             const Obj<sieve_alg_type::coneArray>& cone = sieve_alg_type::nCone(serialMesh, *c_iter, depth);
523:             const int                             idx  = cNumbering->getIndex(*c_iter);
524:             int                                   v    = 0;

526:             for(Mesh::sieve_type::coneArray::iterator v_iter = cone->begin(); v_iter != cone->end(); ++v_iter) {
527:               in.tetrahedronlist[idx*in.numberofcorners + v++] = vNumbering->getIndex(*v_iter);
528:             }
529:           }
530:         }
531:         if (serialMesh->depth() == 3) {
532:           const Obj<Mesh::label_sequence>& boundary = serialMesh->getLabelStratum("marker", 1);

534:           in.numberoftrifaces = 0;
535:           for(Mesh::label_sequence::iterator b_iter = boundary->begin(); b_iter != boundary->end(); ++b_iter) {
536:             if (serialMesh->height(*b_iter) == 1) {
537:               in.numberoftrifaces++;
538:             }
539:           }
540:           if (in.numberoftrifaces > 0) {
541:             int f = 0;

543:             in.trifacelist       = new int[in.numberoftrifaces*3];
544:             in.trifacemarkerlist = new int[in.numberoftrifaces];
545:             for(Mesh::label_sequence::iterator b_iter = boundary->begin(); b_iter != boundary->end(); ++b_iter) {
546:               if (serialMesh->height(*b_iter) == 1) {
547:                 const Obj<Mesh::coneArray>& cone = sieve_alg_type::nCone(serialMesh, *b_iter, 2);
548:                 int                         p    = 0;

550:                 for(Mesh::coneArray::iterator v_iter = cone->begin(); v_iter != cone->end(); ++v_iter) {
551:                   in.trifacelist[f*3 + (p++)] = vNumbering->getIndex(*v_iter);
552:                 }
553:                 in.trifacemarkerlist[f++] = serialMesh->getValue(markers, *b_iter);
554:               }
555:             }
556:           }
557:         }

559:         in.numberofholes = 0;
560:         if (serialMesh->commRank() == 0) {
561:           std::string args("qezQra");

563:           ::tetrahedralize((char *) args.c_str(), &in, &out);
564:         }
565:         in.tetrahedronvolumelist = NULL;
566:         const Obj<Mesh::sieve_type> newSieve = new Mesh::sieve_type(refMesh->comm(), refMesh->debug());
567:         int     numCorners  = 4;
568:         int     numCells    = out.numberoftetrahedra;
569:         int    *newCells       = out.tetrahedronlist;
570:         int     numVertices = out.numberofpoints;
571:         double *coords      = out.pointlist;

573:         ALE::SieveBuilder<Mesh>::buildTopology(newSieve, dim, numCells, newCells, numVertices, interpolate, numCorners, -1, refMesh->getArrowSection("orientation"));
574:         refMesh->setSieve(newSieve);
575:         refMesh->stratify();
576:         ALE::SieveBuilder<Mesh>::buildCoordinates(refMesh, dim, coords);
577:         const Obj<Mesh::label_type>& newMarkers = refMesh->createLabel("marker");


580:         for(int v = 0; v < out.numberofpoints; v++) {
581:           if (out.pointmarkerlist[v]) {
582:             refMesh->setValue(newMarkers, v+out.numberoftetrahedra, out.pointmarkerlist[v]);
583:           }
584:         }
585:         if (interpolate) {
586:           if (out.edgemarkerlist) {
587:             for(int e = 0; e < out.numberofedges; e++) {
588:               if (out.edgemarkerlist[e]) {
589:                 Mesh::point_type endpointA(out.edgelist[e*2+0]+out.numberoftetrahedra);
590:                 Mesh::point_type endpointB(out.edgelist[e*2+1]+out.numberoftetrahedra);
591:                 Obj<Mesh::sieve_type::supportSet> edge = newSieve->nJoin(endpointA, endpointB, 1);

593:                 refMesh->setValue(newMarkers, *edge->begin(), out.edgemarkerlist[e]);
594:               }
595:             }
596:           }
597:           if (out.trifacemarkerlist) {
598:             for(int f = 0; f < out.numberoftrifaces; f++) {
599:               if (out.trifacemarkerlist[f]) {
600:                 Mesh::point_type cornerA(out.trifacelist[f*3+0]+out.numberoftetrahedra);
601:                 Mesh::point_type cornerB(out.trifacelist[f*3+1]+out.numberoftetrahedra);
602:                 Mesh::point_type cornerC(out.trifacelist[f*3+2]+out.numberoftetrahedra);
603:                 Obj<Mesh::sieve_type::supportSet> corners = Mesh::sieve_type::supportSet();
604:                 Obj<Mesh::sieve_type::supportSet> edges   = Mesh::sieve_type::supportSet();
605:                 corners->clear();corners->insert(cornerA);corners->insert(cornerB);
606:                 edges->insert(*newSieve->nJoin1(corners)->begin());
607:                 corners->clear();corners->insert(cornerB);corners->insert(cornerC);
608:                 edges->insert(*newSieve->nJoin1(corners)->begin());
609:                 corners->clear();corners->insert(cornerC);corners->insert(cornerA);
610:                 edges->insert(*newSieve->nJoin1(corners)->begin());
611:                 const Mesh::point_type          face       = *newSieve->nJoin1(edges)->begin();
612:                 const int                       faceMarker = out.trifacemarkerlist[f];
613:                 const Obj<Mesh::coneArray>      closure    = sieve_alg_type::closure(refMesh, face);
614:                 const Mesh::coneArray::iterator end        = closure->end();

616:                 for(Mesh::coneArray::iterator cl_iter = closure->begin(); cl_iter != end; ++cl_iter) {
617:                   refMesh->setValue(newMarkers, *cl_iter, faceMarker);
618:                 }
619:               }
620:             }
621:           }
622:         }
623:         if (refMesh->commSize() > 1) {
624:           return ALE::Distribution<Mesh>::distributeMesh(refMesh);
625:         }
626:         return refMesh;
627:       };
628:       static Obj<Mesh> refineMesh(const Obj<Mesh>& mesh, const Obj<Mesh::real_section_type>& maxVolumes, const bool interpolate = false) {
629:         Obj<Mesh>                          serialMesh       = ALE::Distribution<Mesh>::unifyMesh(mesh);
630:         const Obj<Mesh::real_section_type> serialMaxVolumes = ALE::Distribution<Mesh>::distributeSection(maxVolumes, serialMesh, serialMesh->getDistSendOverlap(), serialMesh->getDistRecvOverlap());

632:         return refineMesh(serialMesh, serialMaxVolumes->restrict(), interpolate);
633:       };
634:       static Obj<Mesh> refineMesh(const Obj<Mesh>& mesh, const double maxVolume, const bool interpolate = false) {
635:         Obj<Mesh> serialMesh;
636:         if (mesh->commSize() > 1) {
637:           serialMesh = ALE::Distribution<Mesh>::unifyMesh(mesh);
638:         } else {
639:           serialMesh = mesh;
640:         }
641:         const int numCells         = serialMesh->heightStratum(0)->size();
642:         double   *serialMaxVolumes = new double[numCells];

644:         for(int c = 0; c < numCells; c++) {
645:           serialMaxVolumes[c] = maxVolume;
646:         }
647:         const Obj<Mesh> refMesh = refineMesh(serialMesh, serialMaxVolumes, interpolate);
648:         delete [] serialMaxVolumes;
649:         return refMesh;
650:       };
651:     };
652:   };
653: #endif
654:   class Generator {
655:     typedef ALE::Mesh Mesh;
656:   public:
657:     static Obj<Mesh> generateMesh(const Obj<Mesh>& boundary, const bool interpolate = false) {
658:       int dim = boundary->getDimension();

660:       if (dim == 1) {
661: #ifdef PETSC_HAVE_TRIANGLE
662:         return ALE::Triangle::Generator::generateMesh(boundary, interpolate);
663: #else
664:         throw ALE::Exception("Mesh generation currently requires Triangle to be installed. Use --download-triangle during configure.");
665: #endif
666:       } else if (dim == 2) {
667: #ifdef PETSC_HAVE_TETGEN
668:         return ALE::TetGen::Generator::generateMesh(boundary, interpolate);
669: #else
670:         throw ALE::Exception("Mesh generation currently requires TetGen to be installed. Use --download-tetgen during configure.");
671: #endif
672:       }
673:       return NULL;
674:     };
675:     static Obj<Mesh> refineMesh(const Obj<Mesh>& mesh, const Obj<Mesh::real_section_type>& maxVolumes, const bool interpolate = false) {
676:       int dim = mesh->getDimension();

678:       if (dim == 2) {
679: #ifdef PETSC_HAVE_TRIANGLE
680:         return ALE::Triangle::Refiner::refineMesh(mesh, maxVolumes, interpolate);
681: #else
682:         throw ALE::Exception("Mesh refinement currently requires Triangle to be installed. Use --download-triangle during configure.");
683: #endif
684:       } else if (dim == 3) {
685: #ifdef PETSC_HAVE_TETGEN
686:         return ALE::TetGen::Refiner::refineMesh(mesh, maxVolumes, interpolate);
687: #else
688:         throw ALE::Exception("Mesh refinement currently requires TetGen to be installed. Use --download-tetgen during configure.");
689: #endif
690:       }
691:       return NULL;
692:     };
693:     static Obj<Mesh> refineMesh(const Obj<Mesh>& mesh, const double maxVolume, const bool interpolate = false) {
694:       int dim = mesh->getDimension();

696:       if (dim == 2) {
697: #ifdef PETSC_HAVE_TRIANGLE
698:         return ALE::Triangle::Refiner::refineMesh(mesh, maxVolume, interpolate);
699: #else
700:         throw ALE::Exception("Mesh refinement currently requires Triangle to be installed. Use --download-triangle during configure.");
701: #endif
702:       } else if (dim == 3) {
703: #ifdef PETSC_HAVE_TETGEN
704:         return ALE::TetGen::Refiner::refineMesh(mesh, maxVolume, interpolate);
705: #else
706:         throw ALE::Exception("Mesh refinement currently requires TetGen to be installed. Use --download-tetgen during configure.");
707: #endif
708:       }
709:       return NULL;
710:     };
711:   };
712: }

714: #endif