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: #define NEW_LABEL
226: #ifdef NEW_LABEL
227:           for(Mesh::label_sequence::iterator e_iter = edges->begin(); e_iter != edges->end(); ++e_iter) {
228:             if (serialMesh->getValue(markers, *e_iter)) {
229:               in.numberofsegments++;
230:             }
231:           }
232:           std::cout << "Number of segments: " << in.numberofsegments << std::endl;
233:           if (in.numberofsegments > 0) {
234:             int s = 0;

236:             PetscMalloc(in.numberofsegments * 2 * sizeof(int), &in.segmentlist);
237:             PetscMalloc(in.numberofsegments * sizeof(int), &in.segmentmarkerlist);
238:             for(Mesh::label_sequence::iterator e_iter = edges->begin(); e_iter != edges->end(); ++e_iter) {
239:               const int edgeMarker = serialMesh->getValue(markers, *e_iter);

241:               if (edgeMarker) {
242:                 const Obj<Mesh::sieve_type::traits::coneSequence>& cone = serialSieve->cone(*e_iter);
243:                 int                                                p    = 0;

245:                 for(Mesh::sieve_type::traits::coneSequence::iterator v_iter = cone->begin(); v_iter != cone->end(); ++v_iter) {
246:                   in.segmentlist[s*2 + (p++)] = vNumbering->getIndex(*v_iter);
247:                 }
248:                 in.segmentmarkerlist[s++] = edgeMarker;
249:               }
250:             }
251:           }
252: #else
253:           const Obj<Mesh::label_type::baseSequence>& boundary = markers->base();

255:           in.numberofsegments = 0;
256:           for(Mesh::label_type::baseSequence::iterator b_iter = boundary->begin(); b_iter != boundary->end(); ++b_iter) {
257:             for(Mesh::label_sequence::iterator e_iter = edges->begin(); e_iter != edges->end(); ++e_iter) {
258:               if (*b_iter == *e_iter) {
259:                 in.numberofsegments++;
260:               }
261:             }
262:           }
263:           if (in.numberofsegments > 0) {
264:             int s = 0;

266:             PetscMalloc(in.numberofsegments * 2 * sizeof(int), &in.segmentlist);
267:             PetscMalloc(in.numberofsegments * sizeof(int), &in.segmentmarkerlist);
268:             for(Mesh::label_type::baseSequence::iterator b_iter = boundary->begin(); b_iter != boundary->end(); ++b_iter) {
269:               for(Mesh::label_sequence::iterator e_iter = edges->begin(); e_iter != edges->end(); ++e_iter) {
270:                 if (*b_iter == *e_iter) {
271:                   const Obj<Mesh::sieve_type::traits::coneSequence>& cone = serialSieve->cone(*e_iter);
272:                   int                                                p    = 0;

274:                   for(Mesh::sieve_type::traits::coneSequence::iterator v_iter = cone->begin(); v_iter != cone->end(); ++v_iter) {
275:                     in.segmentlist[s*2 + (p++)] = vNumbering->getIndex(*v_iter);
276:                   }
277:                   in.segmentmarkerlist[s++] = serialMesh->getValue(markers, *e_iter);
278:                 }
279:               }
280:             }
281:           }
282: #endif
283:         }

285:         in.numberofholes = 0;
286:         if (in.numberofholes > 0) {
287:           PetscMalloc(in.numberofholes * dim * sizeof(int), &in.holelist);
288:         }
289:         if (serialMesh->commRank() == 0) {
290:           std::string args("pqezQra");

292:           triangulate((char *) args.c_str(), &in, &out, NULL);
293:         }
294:         if (in.pointlist)         {PetscFree(in.pointlist);}
295:         if (in.pointmarkerlist)   {PetscFree(in.pointmarkerlist);}
296:         if (in.segmentlist)       {PetscFree(in.segmentlist);}
297:         if (in.segmentmarkerlist) {PetscFree(in.segmentmarkerlist);}
298:         if (in.trianglelist)      {PetscFree(in.trianglelist);}
299:         const Obj<Mesh::sieve_type> newSieve = new Mesh::sieve_type(serialMesh->comm(), serialMesh->debug());
300:         int     numCorners  = 3;
301:         int     numCells    = out.numberoftriangles;
302:         int    *cells       = out.trianglelist;
303:         int     numVertices = out.numberofpoints;
304:         double *coords      = out.pointlist;

306:         ALE::SieveBuilder<Mesh>::buildTopology(newSieve, dim, numCells, cells, numVertices, interpolate, numCorners, -1, refMesh->getArrowSection("orientation"));
307:         refMesh->setSieve(newSieve);
308:         refMesh->stratify();
309:         ALE::SieveBuilder<Mesh>::buildCoordinates(refMesh, dim, coords);
310:         const Obj<Mesh::label_type>& newMarkers = refMesh->createLabel("marker");

312:         if (refMesh->commRank() == 0) {
313:           for(int v = 0; v < out.numberofpoints; v++) {
314:             if (out.pointmarkerlist[v]) {
315:               refMesh->setValue(newMarkers, v+out.numberoftriangles, out.pointmarkerlist[v]);
316:             }
317:           }
318:           if (interpolate) {
319:             for(int e = 0; e < out.numberofedges; e++) {
320:               if (out.edgemarkerlist[e]) {
321:                 const Mesh::point_type vertexA(out.edgelist[e*2+0]+out.numberoftriangles);
322:                 const Mesh::point_type vertexB(out.edgelist[e*2+1]+out.numberoftriangles);
323:                 const Obj<Mesh::sieve_type::supportSet> edge = newSieve->nJoin(vertexA, vertexB, 1);

325:                 refMesh->setValue(newMarkers, *(edge->begin()), out.edgemarkerlist[e]);
326:               }
327:             }
328:           }
329:         }

331:         Generator::finiOutput(&out);
332:         if (refMesh->commSize() > 1) {
333:           return ALE::Distribution<Mesh>::distributeMesh(refMesh);
334:         }
335:         return refMesh;
336:       };
337:       static Obj<Mesh> refineMesh(const Obj<Mesh>& mesh, const Obj<Mesh::real_section_type>& maxVolumes, const bool interpolate = false) {
338:         Obj<Mesh>                          serialMesh       = ALE::Distribution<Mesh>::unifyMesh(mesh);
339:         const Obj<Mesh::real_section_type> serialMaxVolumes = ALE::Distribution<Mesh>::distributeSection(maxVolumes, serialMesh, serialMesh->getDistSendOverlap(), serialMesh->getDistRecvOverlap());

341:         return refineMesh(serialMesh, serialMaxVolumes->restrict(), interpolate);
342:       };
343:       static Obj<Mesh> refineMesh(const Obj<Mesh>& mesh, const double maxVolume, const bool interpolate = false) {
344:         Obj<Mesh> serialMesh;
345:         if (mesh->commSize() > 1) {
346:           serialMesh = ALE::Distribution<Mesh>::unifyMesh(mesh);
347:         } else {
348:           serialMesh = mesh;
349:         }
350:         const int numFaces         = serialMesh->heightStratum(0)->size();
351:         double   *serialMaxVolumes = new double[numFaces];

353:         for(int f = 0; f < numFaces; f++) {
354:           serialMaxVolumes[f] = maxVolume;
355:         }
356:         const Obj<Mesh> refMesh = refineMesh(serialMesh, serialMaxVolumes, interpolate);
357:         delete [] serialMaxVolumes;
358:         return refMesh;
359:       };
360:       static Obj<Mesh> refineMeshLocal(const Obj<Mesh>& mesh, const double maxVolumes[], const bool interpolate = false) {
361:         const int                    dim     = mesh->getDimension();
362:         const Obj<Mesh>              refMesh = new Mesh(mesh->comm(), dim, mesh->debug());
363:         const Obj<Mesh::sieve_type>& sieve   = mesh->getSieve();
364:         struct triangulateio in;
365:         struct triangulateio out;
366:         PetscErrorCode       ierr;

368:         Generator::initInput(&in);
369:         Generator::initOutput(&out);
370:         const Obj<Mesh::label_sequence>&    vertices    = mesh->depthStratum(0);
371:         const Obj<Mesh::label_type>&        markers     = mesh->getLabel("marker");
372:         const Obj<Mesh::real_section_type>& coordinates = mesh->getRealSection("coordinates");
373:         const Obj<Mesh::numbering_type>&    vNumbering  = mesh->getFactory()->getLocalNumbering(mesh, 0);

375:         in.numberofpoints = vertices->size();
376:         if (in.numberofpoints > 0) {
377:           PetscMalloc(in.numberofpoints * dim * sizeof(double), &in.pointlist);
378:           PetscMalloc(in.numberofpoints * sizeof(int), &in.pointmarkerlist);
379:           for(Mesh::label_sequence::iterator v_iter = vertices->begin(); v_iter != vertices->end(); ++v_iter) {
380:             const Mesh::real_section_type::value_type *array = coordinates->restrictPoint(*v_iter);
381:             const int                                  idx   = vNumbering->getIndex(*v_iter);

383:             for(int d = 0; d < dim; d++) {
384:               in.pointlist[idx*dim + d] = array[d];
385:             }
386:             in.pointmarkerlist[idx] = mesh->getValue(markers, *v_iter);
387:           }
388:         }
389:         const Obj<Mesh::label_sequence>& faces      = mesh->heightStratum(0);
390:         const Obj<Mesh::numbering_type>& fNumbering = mesh->getFactory()->getLocalNumbering(mesh, mesh->depth());

392:         in.numberofcorners   = 3;
393:         in.numberoftriangles = faces->size();
394:         in.trianglearealist  = (double *) maxVolumes;
395:         if (in.numberoftriangles > 0) {
396:           PetscMalloc(in.numberoftriangles*in.numberofcorners * sizeof(int), &in.trianglelist);
397:           if (mesh->depth() == 1) {
398:             for(Mesh::label_sequence::iterator f_iter = faces->begin(); f_iter != faces->end(); ++f_iter) {
399:               const Obj<Mesh::sieve_type::traits::coneSequence>& cone = sieve->cone(*f_iter);
400:               const int                                          idx  = fNumbering->getIndex(*f_iter);
401:               int                                                v    = 0;

403:               for(Mesh::sieve_type::traits::coneSequence::iterator c_iter = cone->begin(); c_iter != cone->end(); ++c_iter) {
404:                 in.trianglelist[idx*in.numberofcorners + v++] = vNumbering->getIndex(*c_iter);
405:               }
406:             }
407:           } else if (mesh->depth() == 2) {
408:             for(Mesh::label_sequence::iterator f_iter = faces->begin(); f_iter != faces->end(); ++f_iter) {
409:               typedef ALE::SieveAlg<Mesh> sieve_alg_type;
410:               const Obj<sieve_alg_type::coneArray>& cone = sieve_alg_type::nCone(mesh, *f_iter, 2);
411:               const int                             idx  = fNumbering->getIndex(*f_iter);
412:               int                                   v    = 0;

414:               for(Mesh::sieve_type::coneArray::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 {
419:             throw ALE::Exception("Invalid sieve: Cannot gives sieves of arbitrary depth to Triangle");
420:           }
421:         }
422:         if (mesh->depth() == 2) {
423:           const Obj<Mesh::label_sequence>& edges = mesh->depthStratum(1);
424:           for(Mesh::label_sequence::iterator e_iter = edges->begin(); e_iter != edges->end(); ++e_iter) {
425:             if (mesh->getValue(markers, *e_iter)) {
426:               in.numberofsegments++;
427:             }
428:           }
429:           std::cout << "Number of segments: " << in.numberofsegments << std::endl;
430:           if (in.numberofsegments > 0) {
431:             int s = 0;

433:             PetscMalloc(in.numberofsegments * 2 * sizeof(int), &in.segmentlist);
434:             PetscMalloc(in.numberofsegments * sizeof(int), &in.segmentmarkerlist);
435:             for(Mesh::label_sequence::iterator e_iter = edges->begin(); e_iter != edges->end(); ++e_iter) {
436:               const int edgeMarker = mesh->getValue(markers, *e_iter);

438:               if (edgeMarker) {
439:                 const Obj<Mesh::sieve_type::traits::coneSequence>& cone = sieve->cone(*e_iter);
440:                 int                                                p    = 0;

442:                 for(Mesh::sieve_type::traits::coneSequence::iterator v_iter = cone->begin(); v_iter != cone->end(); ++v_iter) {
443:                   in.segmentlist[s*2 + (p++)] = vNumbering->getIndex(*v_iter);
444:                 }
445:                 in.segmentmarkerlist[s++] = edgeMarker;
446:               }
447:             }
448:           }
449:         }

451:         in.numberofholes = 0;
452:         if (in.numberofholes > 0) {
453:           PetscMalloc(in.numberofholes * dim * sizeof(int), &in.holelist);
454:         }
455:         std::string args("pqezQra");

457:         triangulate((char *) args.c_str(), &in, &out, NULL);
458:         if (in.pointlist)         {PetscFree(in.pointlist);}
459:         if (in.pointmarkerlist)   {PetscFree(in.pointmarkerlist);}
460:         if (in.segmentlist)       {PetscFree(in.segmentlist);}
461:         if (in.segmentmarkerlist) {PetscFree(in.segmentmarkerlist);}
462:         if (in.trianglelist)      {PetscFree(in.trianglelist);}
463:         const Obj<Mesh::sieve_type> newSieve = new Mesh::sieve_type(mesh->comm(), mesh->debug());
464:         int     numCorners  = 3;
465:         int     numCells    = out.numberoftriangles;
466:         int    *cells       = out.trianglelist;
467:         int     numVertices = out.numberofpoints;
468:         double *coords      = out.pointlist;

470:         ALE::SieveBuilder<Mesh>::buildTopologyMultiple(newSieve, dim, numCells, cells, numVertices, interpolate, numCorners, -1, refMesh->getArrowSection("orientation"));
471:         refMesh->setSieve(newSieve);
472:         refMesh->stratify();
473:         ALE::SieveBuilder<Mesh>::buildCoordinatesMultiple(refMesh, dim, coords);
474:         const Obj<Mesh::label_type>& newMarkers = refMesh->createLabel("marker");

476:         for(int v = 0; v < out.numberofpoints; v++) {
477:           if (out.pointmarkerlist[v]) {
478:             refMesh->setValue(newMarkers, v+out.numberoftriangles, out.pointmarkerlist[v]);
479:           }
480:         }
481:         if (interpolate) {
482:           for(int e = 0; e < out.numberofedges; e++) {
483:             if (out.edgemarkerlist[e]) {
484:               const Mesh::point_type vertexA(out.edgelist[e*2+0]+out.numberoftriangles);
485:               const Mesh::point_type vertexB(out.edgelist[e*2+1]+out.numberoftriangles);
486:               const Obj<Mesh::sieve_type::supportSet> edge = newSieve->nJoin(vertexA, vertexB, 1);

488:               refMesh->setValue(newMarkers, *(edge->begin()), out.edgemarkerlist[e]);
489:             }
490:           }
491:         }

493:         Generator::finiOutput(&out);
494:         return refMesh;
495:       };
496:       static Obj<Mesh> refineMeshLocal(const Obj<Mesh>& mesh, const double maxVolume, const bool interpolate = false) {
497:         const int numLocalFaces   = mesh->heightStratum(0)->size();
498:         double   *localMaxVolumes = new double[numLocalFaces];

500:         for(int f = 0; f < numLocalFaces; f++) {
501:           localMaxVolumes[f] = maxVolume;
502:         }
503:         const Obj<Mesh> refMesh = refineMeshLocal(mesh, localMaxVolumes, interpolate);
504:         const Obj<Mesh::sieve_type> refSieve = refMesh->getSieve();
505:         delete [] localMaxVolumes;
506: #if 0
507:         typedef typename ALE::New::Completion<Mesh, typename Mesh::sieve_type::point_type> sieveCompletion;
508:         // This is where we enforce consistency over the overlap
509:         //   We need somehow to update the overlap to account for the new stuff
510:         //
511:         //   1) Since we are refining only, the vertices are invariant
512:         //   2) We need to make another label for the interprocess boundaries so
513:         //      that Triangle will respect them
514:         //   3) We then throw all that label into the new overlap
515:         //
516:         // Alternative: Figure out explicitly which segments were refined, and then
517:         //   communicated the refinement over the old overlap. Use this info to locally
518:         //   construct the new overlap and flip to get a decent mesh
519:         sieveCompletion::scatterCones(refSieve, refSieve, reMesh->getDistSendOverlap(), refMesh->getDistRecvOverlap(), refMesh);
520: #endif
521:         return refMesh;
522:       };
523:     };
524:   };
525: #endif
526: #ifdef PETSC_HAVE_TETGEN
527:   namespace TetGen {
528:     class Generator {
529:     public:
530:       static Obj<Mesh> generateMesh(const Obj<Mesh>& boundary, const bool interpolate = false) {
531:         typedef ALE::SieveAlg<Mesh> sieve_alg_type;
532:         const int         dim   = 3;
533:         Obj<Mesh>         mesh  = new Mesh(boundary->comm(), dim, boundary->debug());
534:         const PetscMPIInt rank  = mesh->commRank();
535:         const bool        createConvexHull = false;
536:         ::tetgenio        in;
537:         ::tetgenio        out;

539:         const Obj<Mesh::label_sequence>&    vertices    = boundary->depthStratum(0);
540:         const Obj<Mesh::numbering_type>&    vNumbering  = boundary->getFactory()->getLocalNumbering(boundary, 0);
541:         const Obj<Mesh::real_section_type>& coordinates = boundary->getRealSection("coordinates");
542:         const Obj<Mesh::label_type>&        markers     = boundary->getLabel("marker");


545:         in.numberofpoints = vertices->size();
546:         if (in.numberofpoints > 0) {
547:           in.pointlist       = new double[in.numberofpoints*dim];
548:           in.pointmarkerlist = new int[in.numberofpoints];
549:           for(Mesh::label_sequence::iterator v_iter = vertices->begin(); v_iter != vertices->end(); ++v_iter) {
550:             const Mesh::real_section_type::value_type *array = coordinates->restrictPoint(*v_iter);
551:             const int                                  idx   = vNumbering->getIndex(*v_iter);

553:             for(int d = 0; d < dim; d++) {
554:               in.pointlist[idx*dim + d] = array[d];
555:             }
556:             in.pointmarkerlist[idx] = boundary->getValue(markers, *v_iter);
557:           }
558:         }

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

563:         in.numberoffacets = facets->size();
564:         if (in.numberoffacets > 0) {
565:           in.facetlist       = new tetgenio::facet[in.numberoffacets];
566:           in.facetmarkerlist = new int[in.numberoffacets];
567:           for(Mesh::label_sequence::iterator f_iter = facets->begin(); f_iter != facets->end(); ++f_iter) {
568:             const Obj<sieve_alg_type::coneArray>& cone = sieve_alg_type::nCone(boundary, *f_iter, boundary->depth());
569:             const int                             idx  = fNumbering->getIndex(*f_iter);

571:             in.facetlist[idx].numberofpolygons = 1;
572:             in.facetlist[idx].polygonlist      = new tetgenio::polygon[in.facetlist[idx].numberofpolygons];
573:             in.facetlist[idx].numberofholes    = 0;
574:             in.facetlist[idx].holelist         = NULL;

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

579:             poly->numberofvertices = cone->size();
580:             poly->vertexlist       = new int[poly->numberofvertices];
581:             for(sieve_alg_type::coneArray::iterator c_iter = cone->begin(); c_iter != cone->end(); ++c_iter) {
582:               const int vIdx = vNumbering->getIndex(*c_iter);

584:               poly->vertexlist[c++] = vIdx;
585:             }
586:             in.facetmarkerlist[idx] = boundary->getValue(markers, *f_iter);
587:           }
588:         }

590:         in.numberofholes = 0;
591:         if (rank == 0) {
592:           // Normal operation
593:           std::string args("pqezQ");
594:           // Just make tetrahedrons
595: //           std::string args("efzV");
596:           // Adds a center point
597: //           std::string args("pqezQi");
598: //           in.numberofaddpoints = 1;
599: //           in.addpointlist      = new double[in.numberofaddpoints*dim];
600: //           in.addpointlist[0]   = 0.5;
601: //           in.addpointlist[1]   = 0.5;
602: //           in.addpointlist[2]   = 0.5;

604:           if (createConvexHull) args += "c";
605:           ::tetrahedralize((char *) args.c_str(), &in, &out);
606:         }
607:         const Obj<Mesh::sieve_type> newSieve = new Mesh::sieve_type(mesh->comm(), mesh->debug());
608:         int     numCorners  = 4;
609:         int     numCells    = out.numberoftetrahedra;
610:         int    *cells       = out.tetrahedronlist;
611:         int     numVertices = out.numberofpoints;
612:         double *coords      = out.pointlist;

614:         if (!interpolate) {
615:           for(int c = 0; c < numCells; ++c) {
616:             int tmp = cells[c*4+0];
617:             cells[c*4+0] = cells[c*4+1];
618:             cells[c*4+1] = tmp;
619:           }
620:         }
621:         ALE::SieveBuilder<Mesh>::buildTopology(newSieve, dim, numCells, cells, numVertices, interpolate, numCorners, -1, mesh->getArrowSection("orientation"));
622:         mesh->setSieve(newSieve);
623:         mesh->stratify();
624:         ALE::SieveBuilder<Mesh>::buildCoordinates(mesh, dim, coords);
625:         const Obj<Mesh::label_type>& newMarkers = mesh->createLabel("marker");
626: 
627:         for(int v = 0; v < out.numberofpoints; v++) {
628:           if (out.pointmarkerlist[v]) {
629:             mesh->setValue(newMarkers, v+out.numberoftetrahedra, out.pointmarkerlist[v]);
630:           }
631:         }
632:         if (interpolate) {
633:           if (out.edgemarkerlist) {
634:             for(int e = 0; e < out.numberofedges; e++) {
635:               if (out.edgemarkerlist[e]) {
636:                 Mesh::point_type endpointA(out.edgelist[e*2+0]+out.numberoftetrahedra);
637:                 Mesh::point_type endpointB(out.edgelist[e*2+1]+out.numberoftetrahedra);
638:                 Obj<Mesh::sieve_type::supportSet> edge = newSieve->nJoin(endpointA, endpointB, 1);

640:                 mesh->setValue(newMarkers, *edge->begin(), out.edgemarkerlist[e]);
641:               }
642:             }
643:           }
644:           if (out.trifacemarkerlist) {
645:             // Work around TetGen bug for raw tetrahedralization
646:             //   The boundary faces are 0,1,4,5,8,9,11,12,13,15,16,17
647: //             for(int f = 0; f < out.numberoftrifaces; f++) {
648: //               if (out.trifacemarkerlist[f]) {
649: //                 out.trifacemarkerlist[f] = 0;
650: //               } else {
651: //                 out.trifacemarkerlist[f] = 1;
652: //               }
653: //             }
654:             for(int f = 0; f < out.numberoftrifaces; f++) {
655:               if (out.trifacemarkerlist[f]) {
656:                 Mesh::point_type cornerA(out.trifacelist[f*3+0]+out.numberoftetrahedra);
657:                 Mesh::point_type cornerB(out.trifacelist[f*3+1]+out.numberoftetrahedra);
658:                 Mesh::point_type cornerC(out.trifacelist[f*3+2]+out.numberoftetrahedra);
659:                 Obj<Mesh::sieve_type::supportSet> corners = Mesh::sieve_type::supportSet();
660:                 Obj<Mesh::sieve_type::supportSet> edges   = Mesh::sieve_type::supportSet();
661:                 corners->clear();corners->insert(cornerA);corners->insert(cornerB);
662:                 edges->insert(*newSieve->nJoin1(corners)->begin());
663:                 corners->clear();corners->insert(cornerB);corners->insert(cornerC);
664:                 edges->insert(*newSieve->nJoin1(corners)->begin());
665:                 corners->clear();corners->insert(cornerC);corners->insert(cornerA);
666:                 edges->insert(*newSieve->nJoin1(corners)->begin());
667:                 const Mesh::point_type          face       = *newSieve->nJoin1(edges)->begin();
668:                 const int                       faceMarker = out.trifacemarkerlist[f];
669:                 const Obj<Mesh::coneArray>      closure    = sieve_alg_type::closure(mesh, face);
670:                 const Mesh::coneArray::iterator end        = closure->end();

672:                 for(Mesh::coneArray::iterator cl_iter = closure->begin(); cl_iter != end; ++cl_iter) {
673:                   mesh->setValue(newMarkers, *cl_iter, faceMarker);
674:                 }
675:               }
676:             }
677:           }
678:         }
679:         return mesh;
680:       };
681:     };
682:     class Refiner {
683:     public:
684:       static Obj<Mesh> refineMesh(const Obj<Mesh>& serialMesh, const double maxVolumes[], const bool interpolate = false) {
685:         typedef ALE::SieveAlg<Mesh> sieve_alg_type;
686:         const int       dim     = serialMesh->getDimension();
687:         const int       depth   = serialMesh->depth();
688:         const Obj<Mesh> refMesh = new Mesh(serialMesh->comm(), dim, serialMesh->debug());
689:         ::tetgenio      in;
690:         ::tetgenio      out;

692:         const Obj<Mesh::label_sequence>&    vertices    = serialMesh->depthStratum(0);
693:         const Obj<Mesh::label_type>&        markers     = serialMesh->getLabel("marker");
694:         const Obj<Mesh::real_section_type>& coordinates = serialMesh->getRealSection("coordinates");
695:         const Obj<Mesh::numbering_type>&    vNumbering  = serialMesh->getFactory()->getLocalNumbering(serialMesh, 0);

697:         in.numberofpoints = vertices->size();
698:         if (in.numberofpoints > 0) {
699:           in.pointlist       = new double[in.numberofpoints*dim];
700:           in.pointmarkerlist = new int[in.numberofpoints];
701:           for(Mesh::label_sequence::iterator v_iter = vertices->begin(); v_iter != vertices->end(); ++v_iter) {
702:             const Mesh::real_section_type::value_type *array = coordinates->restrictPoint(*v_iter);
703:             const int                                  idx   = vNumbering->getIndex(*v_iter);

705:             for(int d = 0; d < dim; d++) {
706:               in.pointlist[idx*dim + d] = array[d];
707:             }
708:             in.pointmarkerlist[idx] = serialMesh->getValue(markers, *v_iter);
709:           }
710:         }
711:         const Obj<Mesh::label_sequence>& cells      = serialMesh->heightStratum(0);
712:         const Obj<Mesh::numbering_type>& cNumbering = serialMesh->getFactory()->getLocalNumbering(serialMesh, depth);

714:         in.numberofcorners       = 4;
715:         in.numberoftetrahedra    = cells->size();
716:         in.tetrahedronvolumelist = (double *) maxVolumes;
717:         if (in.numberoftetrahedra > 0) {
718:           in.tetrahedronlist     = new int[in.numberoftetrahedra*in.numberofcorners];
719:           for(Mesh::label_sequence::iterator c_iter = cells->begin(); c_iter != cells->end(); ++c_iter) {
720:             typedef ALE::SieveAlg<Mesh> sieve_alg_type;
721:             const Obj<sieve_alg_type::coneArray>& cone = sieve_alg_type::nCone(serialMesh, *c_iter, depth);
722:             const int                             idx  = cNumbering->getIndex(*c_iter);
723:             int                                   v    = 0;

725:             for(Mesh::sieve_type::coneArray::iterator v_iter = cone->begin(); v_iter != cone->end(); ++v_iter) {
726:               in.tetrahedronlist[idx*in.numberofcorners + v++] = vNumbering->getIndex(*v_iter);
727:             }
728:           }
729:         }
730:         if (serialMesh->depth() == 3) {
731:           const Obj<Mesh::label_sequence>& boundary = serialMesh->getLabelStratum("marker", 1);

733:           in.numberoftrifaces = 0;
734:           for(Mesh::label_sequence::iterator b_iter = boundary->begin(); b_iter != boundary->end(); ++b_iter) {
735:             if (serialMesh->height(*b_iter) == 1) {
736:               in.numberoftrifaces++;
737:             }
738:           }
739:           if (in.numberoftrifaces > 0) {
740:             int f = 0;

742:             in.trifacelist       = new int[in.numberoftrifaces*3];
743:             in.trifacemarkerlist = new int[in.numberoftrifaces];
744:             for(Mesh::label_sequence::iterator b_iter = boundary->begin(); b_iter != boundary->end(); ++b_iter) {
745:               if (serialMesh->height(*b_iter) == 1) {
746:                 const Obj<Mesh::coneArray>& cone = sieve_alg_type::nCone(serialMesh, *b_iter, 2);
747:                 int                         p    = 0;

749:                 for(Mesh::coneArray::iterator v_iter = cone->begin(); v_iter != cone->end(); ++v_iter) {
750:                   in.trifacelist[f*3 + (p++)] = vNumbering->getIndex(*v_iter);
751:                 }
752:                 in.trifacemarkerlist[f++] = serialMesh->getValue(markers, *b_iter);
753:               }
754:             }
755:           }
756:         }

758:         in.numberofholes = 0;
759:         if (serialMesh->commRank() == 0) {
760:           std::string args("qezQra");

762:           ::tetrahedralize((char *) args.c_str(), &in, &out);
763:         }
764:         in.tetrahedronvolumelist = NULL;
765:         const Obj<Mesh::sieve_type> newSieve = new Mesh::sieve_type(refMesh->comm(), refMesh->debug());
766:         int     numCorners  = 4;
767:         int     numCells    = out.numberoftetrahedra;
768:         int    *newCells       = out.tetrahedronlist;
769:         int     numVertices = out.numberofpoints;
770:         double *coords      = out.pointlist;

772:         if (!interpolate) {
773:           for(int c = 0; c < numCells; ++c) {
774:             int tmp = newCells[c*4+0];
775:             newCells[c*4+0] = newCells[c*4+1];
776:             newCells[c*4+1] = tmp;
777:           }
778:         }
779:         ALE::SieveBuilder<Mesh>::buildTopology(newSieve, dim, numCells, newCells, numVertices, interpolate, numCorners, -1, refMesh->getArrowSection("orientation"));
780:         refMesh->setSieve(newSieve);
781:         refMesh->stratify();
782:         ALE::SieveBuilder<Mesh>::buildCoordinates(refMesh, dim, coords);
783:         const Obj<Mesh::label_type>& newMarkers = refMesh->createLabel("marker");


786:         for(int v = 0; v < out.numberofpoints; v++) {
787:           if (out.pointmarkerlist[v]) {
788:             refMesh->setValue(newMarkers, v+out.numberoftetrahedra, out.pointmarkerlist[v]);
789:           }
790:         }
791:         if (interpolate) {
792:           if (out.edgemarkerlist) {
793:             for(int e = 0; e < out.numberofedges; e++) {
794:               if (out.edgemarkerlist[e]) {
795:                 Mesh::point_type endpointA(out.edgelist[e*2+0]+out.numberoftetrahedra);
796:                 Mesh::point_type endpointB(out.edgelist[e*2+1]+out.numberoftetrahedra);
797:                 Obj<Mesh::sieve_type::supportSet> edge = newSieve->nJoin(endpointA, endpointB, 1);

799:                 refMesh->setValue(newMarkers, *edge->begin(), out.edgemarkerlist[e]);
800:               }
801:             }
802:           }
803:           if (out.trifacemarkerlist) {
804:             for(int f = 0; f < out.numberoftrifaces; f++) {
805:               if (out.trifacemarkerlist[f]) {
806:                 Mesh::point_type cornerA(out.trifacelist[f*3+0]+out.numberoftetrahedra);
807:                 Mesh::point_type cornerB(out.trifacelist[f*3+1]+out.numberoftetrahedra);
808:                 Mesh::point_type cornerC(out.trifacelist[f*3+2]+out.numberoftetrahedra);
809:                 Obj<Mesh::sieve_type::supportSet> corners = Mesh::sieve_type::supportSet();
810:                 Obj<Mesh::sieve_type::supportSet> edges   = Mesh::sieve_type::supportSet();
811:                 corners->clear();corners->insert(cornerA);corners->insert(cornerB);
812:                 edges->insert(*newSieve->nJoin1(corners)->begin());
813:                 corners->clear();corners->insert(cornerB);corners->insert(cornerC);
814:                 edges->insert(*newSieve->nJoin1(corners)->begin());
815:                 corners->clear();corners->insert(cornerC);corners->insert(cornerA);
816:                 edges->insert(*newSieve->nJoin1(corners)->begin());
817:                 const Mesh::point_type          face       = *newSieve->nJoin1(edges)->begin();
818:                 const int                       faceMarker = out.trifacemarkerlist[f];
819:                 const Obj<Mesh::coneArray>      closure    = sieve_alg_type::closure(refMesh, face);
820:                 const Mesh::coneArray::iterator end        = closure->end();

822:                 for(Mesh::coneArray::iterator cl_iter = closure->begin(); cl_iter != end; ++cl_iter) {
823:                   refMesh->setValue(newMarkers, *cl_iter, faceMarker);
824:                 }
825:               }
826:             }
827:           }
828:         }
829:         if (refMesh->commSize() > 1) {
830:           return ALE::Distribution<Mesh>::distributeMesh(refMesh);
831:         }
832:         return refMesh;
833:       };
834:       static Obj<Mesh> refineMesh(const Obj<Mesh>& mesh, const Obj<Mesh::real_section_type>& maxVolumes, const bool interpolate = false) {
835:         Obj<Mesh>                          serialMesh       = ALE::Distribution<Mesh>::unifyMesh(mesh);
836:         const Obj<Mesh::real_section_type> serialMaxVolumes = ALE::Distribution<Mesh>::distributeSection(maxVolumes, serialMesh, serialMesh->getDistSendOverlap(), serialMesh->getDistRecvOverlap());

838:         return refineMesh(serialMesh, serialMaxVolumes->restrict(), interpolate);
839:       };
840:       static Obj<Mesh> refineMesh(const Obj<Mesh>& mesh, const double maxVolume, const bool interpolate = false) {
841:         Obj<Mesh> serialMesh;
842:         if (mesh->commSize() > 1) {
843:           serialMesh = ALE::Distribution<Mesh>::unifyMesh(mesh);
844:         } else {
845:           serialMesh = mesh;
846:         }
847:         const int numCells         = serialMesh->heightStratum(0)->size();
848:         double   *serialMaxVolumes = new double[numCells];

850:         for(int c = 0; c < numCells; c++) {
851:           serialMaxVolumes[c] = maxVolume;
852:         }
853:         const Obj<Mesh> refMesh = refineMesh(serialMesh, serialMaxVolumes, interpolate);
854:         delete [] serialMaxVolumes;
855:         return refMesh;
856:       };
857:     };
858:   };
859: #endif
860:   class Generator {
861:     typedef ALE::Mesh Mesh;
862:   public:
863:     static Obj<Mesh> generateMesh(const Obj<Mesh>& boundary, const bool interpolate = false) {
864:       int dim = boundary->getDimension();

866:       if (dim == 1) {
867: #ifdef PETSC_HAVE_TRIANGLE
868:         return ALE::Triangle::Generator::generateMesh(boundary, interpolate);
869: #else
870:         throw ALE::Exception("Mesh generation currently requires Triangle to be installed. Use --download-triangle during configure.");
871: #endif
872:       } else if (dim == 2) {
873: #ifdef PETSC_HAVE_TETGEN
874:         return ALE::TetGen::Generator::generateMesh(boundary, interpolate);
875: #else
876:         throw ALE::Exception("Mesh generation currently requires TetGen to be installed. Use --download-tetgen during configure.");
877: #endif
878:       }
879:       return NULL;
880:     };
881:     static Obj<Mesh> refineMesh(const Obj<Mesh>& mesh, const Obj<Mesh::real_section_type>& maxVolumes, const bool interpolate = false) {
882:       int dim = mesh->getDimension();

884:       if (dim == 2) {
885: #ifdef PETSC_HAVE_TRIANGLE
886:         return ALE::Triangle::Refiner::refineMesh(mesh, maxVolumes, interpolate);
887: #else
888:         throw ALE::Exception("Mesh refinement currently requires Triangle to be installed. Use --download-triangle during configure.");
889: #endif
890:       } else if (dim == 3) {
891: #ifdef PETSC_HAVE_TETGEN
892:         return ALE::TetGen::Refiner::refineMesh(mesh, maxVolumes, interpolate);
893: #else
894:         throw ALE::Exception("Mesh refinement currently requires TetGen to be installed. Use --download-tetgen during configure.");
895: #endif
896:       }
897:       return NULL;
898:     };
899:     static Obj<Mesh> refineMesh(const Obj<Mesh>& mesh, const double maxVolume, const bool interpolate = false) {
900:       int dim = mesh->getDimension();

902:       if (dim == 2) {
903: #ifdef PETSC_HAVE_TRIANGLE
904:         return ALE::Triangle::Refiner::refineMesh(mesh, maxVolume, interpolate);
905: #else
906:         throw ALE::Exception("Mesh refinement currently requires Triangle to be installed. Use --download-triangle during configure.");
907: #endif
908:       } else if (dim == 3) {
909: #ifdef PETSC_HAVE_TETGEN
910:         return ALE::TetGen::Refiner::refineMesh(mesh, maxVolume, interpolate);
911: #else
912:         throw ALE::Exception("Mesh refinement currently requires TetGen to be installed. Use --download-tetgen during configure.");
913: #endif
914:       }
915:       return NULL;
916:     };
917:     static Obj<Mesh> refineMeshLocal(const Obj<Mesh>& mesh, const double maxVolume, const bool interpolate = false) {
918:       int dim = mesh->getDimension();

920:       if (dim == 2) {
921: #ifdef PETSC_HAVE_TRIANGLE
922:         return ALE::Triangle::Refiner::refineMeshLocal(mesh, maxVolume, interpolate);
923: #else
924:         throw ALE::Exception("Mesh refinement currently requires Triangle to be installed. Use --download-triangle during configure.");
925: #endif
926:       } else if (dim == 3) {
927: #ifdef PETSC_HAVE_TETGEN
928:         return ALE::TetGen::Refiner::refineMesh(mesh, maxVolume, interpolate);
929: #else
930:         throw ALE::Exception("Mesh refinement currently requires TetGen to be installed. Use --download-tetgen during configure.");
931: #endif
932:       }
933:       return NULL;
934:     };
935:   };
936: }

938: #endif