Actual source code: Mesh.hh
1: #ifndef included_ALE_Mesh_hh
2: #define included_ALE_Mesh_hh
4: #ifndef included_ALE_CoSifter_hh
5: #include <CoSifter.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 {
17: namespace Two {
18: class Mesh {
19: public:
20: typedef ALE::def::Point point_type;
21: typedef std::vector<point_type> PointArray;
22: typedef ALE::def::Sieve<point_type,int> sieve_type;
23: typedef point_type patch_type;
24: typedef CoSifter<sieve_type, patch_type, point_type, int> bundle_type;
25: typedef CoSifter<sieve_type, patch_type, point_type, double> field_type;
26: int debug;
27: private:
28: Obj<sieve_type> topology;
29: Obj<field_type> coordinates;
30: Obj<field_type> boundary;
31: std::map<int, Obj<bundle_type> > bundles;
32: std::map<std::string, Obj<field_type> > fields;
33: MPI_Comm comm;
34: int rank;
35: int dim;
36: public:
37: Mesh(MPI_Comm comm, int dimension, int debug = 0) : debug(debug), dim(dimension) {
38: this->setComm(comm);
39: this->topology = sieve_type(debug);
40: this->coordinates = field_type(debug);
41: this->boundary = field_type(debug);
42: };
44: MPI_Comm getComm() {return this->comm;};
45: void setComm(MPI_Comm comm) {this->comm = comm; MPI_Comm_rank(comm, &this->rank);};
46: int getRank() {return this->rank;};
47: Obj<sieve_type> getTopology() {return this->topology;};
48: void setTopology(const Obj<sieve_type>& topology) {this->topology = topology;};
49: int getDimension() {return this->dim;};
50: void setDimension(int dim) {this->dim = dim;};
51: Obj<field_type> getCoordinates() {return this->coordinates;};
52: void setCoordinates(const Obj<field_type>& coordinates) {this->coordinates = coordinates;};
53: Obj<field_type> getBoundary() {return this->boundary;};
54: void setBoundary(const Obj<field_type>& boundary) {this->boundary = boundary;};
55: Obj<bundle_type> getBundle(const int dim) {
56: if (this->bundles.find(dim) == this->bundles.end()) {
57: Obj<bundle_type> bundle = bundle_type(debug);
59: // Need to globalize indices (that is what we might use the value ints for)
60: std::cout << "Creating new bundle for dim " << dim << std::endl;
61: bundle->setTopology(this->topology);
62: bundle->setPatch(this->topology->leaves(), bundle_type::patch_type());
63: bundle->setFiberDimensionByDepth(bundle_type::patch_type(), dim, 1);
64: bundle->orderPatches();
65: // "element" reorder is in vertexBundle by default, and intermediate bundles could be handled by a cell tuple
66: this->bundles[dim] = bundle;
67: }
68: return this->bundles[dim];
69: };
70: Obj<field_type> getField(const std::string& name) {
71: if (this->fields.find(name) == this->fields.end()) {
72: Obj<field_type> field = field_type(debug);
74: std::cout << "Creating new field " << name << std::endl;
75: field->setTopology(this->topology);
76: this->fields[name] = field;
77: }
78: return this->fields[name];
79: }
81: void buildFaces(int dim, std::map<int, int*> *curSimplex, Obj<PointArray> boundary, point_type& simplex) {
82: Obj<PointArray> faces = PointArray();
84: if (debug > 1) {std::cout << " Building faces for boundary(size " << boundary->size() << "), dim " << dim << std::endl;}
85: if (dim > 1) {
86: // Use the cone construction
87: for(PointArray::iterator b_itor = boundary->begin(); b_itor != boundary->end(); ++b_itor) {
88: Obj<PointArray> faceBoundary = PointArray();
89: point_type face;
91: for(PointArray::iterator i_itor = boundary->begin(); i_itor != boundary->end(); ++i_itor) {
92: if (i_itor != b_itor) {
93: faceBoundary->push_back(*i_itor);
94: }
95: }
96: if (debug > 1) {std::cout << " boundary point " << *b_itor << std::endl;}
97: this->buildFaces(dim-1, curSimplex, faceBoundary, face);
98: faces->push_back(face);
99: }
100: } else {
101: if (debug > 1) {std::cout << " Just set faces to boundary in 1d" << std::endl;}
102: faces = boundary;
103: }
104: if (debug > 1) {
105: for(PointArray::iterator f_itor = faces->begin(); f_itor != faces->end(); ++f_itor) {
106: std::cout << " face point " << *f_itor << std::endl;
107: }
108: }
109: // We always create the toplevel, so we could shortcircuit somehow
110: // Should not have to loop here since the meet of just 2 boundary elements is an element
111: PointArray::iterator f_itor = faces->begin();
112: point_type start = *f_itor;
113: f_itor++;
114: point_type next = *f_itor;
115: Obj<ALE::def::PointSet> preElement = this->topology->nJoin(start, next, 1);
117: if (preElement->size() > 0) {
118: simplex = *preElement->begin();
119: if (debug > 1) {std::cout << " Found old simplex " << simplex << std::endl;}
120: } else {
121: int color = 0;
123: simplex = point_type(0, (*(*curSimplex)[dim])++);
124: for(PointArray::iterator f_itor = faces->begin(); f_itor != faces->end(); ++f_itor) {
125: this->topology->addArrow(*f_itor, simplex, color++);
126: }
127: if (debug > 1) {std::cout << " Added simplex " << simplex << " dim " << dim << std::endl;}
128: }
129: };
133: // Build a topology from a connectivity description
134: // (0, 0) ... (0, numSimplices-1): dim-dimensional simplices
135: // (0, numSimplices) ... (0, numVertices): vertices
136: // The other simplices are numbered as they are requested
137: void buildTopology(int numSimplices, int simplices[], int numVertices, bool interpolate = true) {
138: ALE_LOG_EVENT_BEGIN;
139: // Create a map from dimension to the current element number for that dimension
140: std::map<int,int*> curElement = std::map<int,int*>();
141: int curSimplex = 0;
142: int curVertex = numSimplices;
143: int newElement = numSimplices+numVertices;
144: Obj<PointArray> boundary = PointArray();
146: curElement[0] = &curVertex;
147: curElement[dim] = &curSimplex;
148: for(int d = 1; d < dim; d++) {
149: curElement[d] = &newElement;
150: }
151: for(int s = 0; s < numSimplices; s++) {
152: point_type simplex(0, s);
154: // Build the simplex
155: if (interpolate) {
156: boundary->clear();
157: for(int b = 0; b < dim+1; b++) {
158: point_type vertex(0, simplices[s*(dim+1)+b]+numSimplices);
160: if (debug > 1) {std::cout << "Adding boundary node " << vertex << std::endl;}
161: boundary->push_back(vertex);
162: }
163: if (debug) {std::cout << "simplex " << s << " boundary size " << boundary->size() << std::endl;}
165: this->buildFaces(this->dim, &curElement, boundary, simplex);
166: } else {
167: for(int b = 0; b < dim+1; b++) {
168: point_type p(0, simplices[s*(dim+1)+b]+numSimplices);
170: this->topology->addArrow(p, simplex, b);
171: }
172: }
173: }
174: ALE_LOG_EVENT_END;
175: };
179: void createVertexBundle(int numSimplices, int simplices[]) {
180: ALE_LOG_EVENT_BEGIN;
181: Obj<bundle_type> vertexBundle = this->getBundle(0);
182: Obj<sieve_type::heightSequence> elements = this->topology->heightStratum(0);
183: std::string orderName("element");
185: for(sieve_type::heightSequence::iterator e_iter = elements->begin(); e_iter != elements->end(); ++e_iter) {
186: // setFiberDimensionByDepth() does not work here since we only want it to apply to the patch cone
187: // What we really need is the depthStratum relative to the patch
188: Obj<PointArray> patch = PointArray();
190: for(int b = 0; b < dim+1; b++) {
191: patch->push_back(point_type(0, simplices[(*e_iter).index*(this->dim+1)+b]+numSimplices));
192: }
193: vertexBundle->setPatch(orderName, patch, *e_iter);
194: for(PointArray::iterator p_iter = patch->begin(); p_iter != patch->end(); ++p_iter) {
195: vertexBundle->setFiberDimension(orderName, *e_iter, *p_iter, 1);
196: }
197: }
198: vertexBundle->orderPatches(orderName);
199: ALE_LOG_EVENT_END;
200: };
204: void createSerialCoordinates(int embedDim, int numSimplices, double coords[]) {
205: ALE_LOG_EVENT_BEGIN;
206: patch_type patch;
208: this->coordinates->setTopology(this->topology);
209: this->coordinates->setPatch(this->topology->leaves(), patch);
210: this->coordinates->setFiberDimensionByDepth(patch, 0, embedDim);
211: this->coordinates->orderPatches();
212: Obj<sieve_type::depthSequence> vertices = this->topology->depthStratum(0);
213: for(sieve_type::depthSequence::iterator v_itor = vertices->begin(); v_itor != vertices->end(); v_itor++) {
214: this->coordinates->update(patch, *v_itor, &coords[((*v_itor).index - numSimplices)*embedDim]);
215: }
216: Obj<bundle_type> vertexBundle = this->getBundle(0);
217: Obj<sieve_type::heightSequence> elements = this->topology->heightStratum(0);
218: std::string orderName("element");
220: for(sieve_type::heightSequence::iterator e_iter = elements->begin(); e_iter != elements->end(); e_iter++) {
221: // setFiberDimensionByDepth() does not work here since we only want it to apply to the patch cone
222: // What we really need is the depthStratum relative to the patch
223: Obj<bundle_type::order_type::coneSequence> cone = vertexBundle->getPatch(*e_iter);
225: this->coordinates->setPatch(orderName, cone, *e_iter);
226: for(bundle_type::order_type::coneSequence::iterator c_iter = cone->begin(); c_iter != cone->end(); ++c_iter) {
227: this->coordinates->setFiberDimension(orderName, *e_iter, *c_iter, embedDim);
228: }
229: }
230: this->coordinates->orderPatches(orderName);
231: ALE_LOG_EVENT_END;
232: };
234: // Create a serial mesh
235: void populate(int numSimplices, int simplices[], int numVertices, double coords[], bool interpolate = true) {
236: this->topology->setStratification(false);
237: if (this->getRank() == 0) {
238: this->buildTopology(numSimplices, simplices, numVertices, interpolate);
239: }
240: this->topology->stratify();
241: this->topology->setStratification(true);
242: this->createVertexBundle(numSimplices, simplices);
243: this->createSerialCoordinates(this->dim, numSimplices, coords);
244: };
245: };
247: class Generator {
248: #ifdef PETSC_HAVE_TRIANGLE
249: static void initInput_Triangle(struct triangulateio *inputCtx) {
250: inputCtx->numberofpoints = 0;
251: inputCtx->numberofpointattributes = 0;
252: inputCtx->pointlist = NULL;
253: inputCtx->pointattributelist = NULL;
254: inputCtx->pointmarkerlist = NULL;
255: inputCtx->numberofsegments = 0;
256: inputCtx->segmentlist = NULL;
257: inputCtx->segmentmarkerlist = NULL;
258: inputCtx->numberoftriangleattributes = 0;
259: inputCtx->numberofholes = 0;
260: inputCtx->holelist = NULL;
261: inputCtx->numberofregions = 0;
262: inputCtx->regionlist = NULL;
263: };
264: static void initOutput_Triangle(struct triangulateio *outputCtx) {
265: outputCtx->pointlist = NULL;
266: outputCtx->pointattributelist = NULL;
267: outputCtx->pointmarkerlist = NULL;
268: outputCtx->trianglelist = NULL;
269: outputCtx->triangleattributelist = NULL;
270: outputCtx->neighborlist = NULL;
271: outputCtx->segmentlist = NULL;
272: outputCtx->segmentmarkerlist = NULL;
273: outputCtx->edgelist = NULL;
274: outputCtx->edgemarkerlist = NULL;
275: };
276: static void finiOutput_Triangle(struct triangulateio *outputCtx) {
277: free(outputCtx->pointmarkerlist);
278: free(outputCtx->edgelist);
279: free(outputCtx->edgemarkerlist);
280: free(outputCtx->trianglelist);
281: free(outputCtx->neighborlist);
282: };
285: static Obj<Mesh> generate_Triangle(Obj<Mesh> boundary, bool interpolate) {
286: struct triangulateio in;
287: struct triangulateio out;
288: int dim = 2;
289: Obj<Mesh> m = Mesh(boundary->getComm(), dim);
290: Obj<Mesh::sieve_type> bdTopology = boundary->getTopology();
291: PetscMPIInt rank;
292: PetscErrorCode ierr;
294: MPI_Comm_rank(boundary->getComm(), &rank);
295: initInput_Triangle(&in);
296: initOutput_Triangle(&out);
297: if (rank == 0) {
298: std::string args("pqenzQ");
299: bool createConvexHull = false;
300: Obj<Mesh::sieve_type::depthSequence> vertices = bdTopology->depthStratum(0);
301: Obj<Mesh::bundle_type> vertexBundle = boundary->getBundle(0);
302: Mesh::field_type::patch_type patch;
304: in.numberofpoints = vertices->size();
305: if (in.numberofpoints > 0) {
306: Obj<Mesh::field_type> coordinates = boundary->getCoordinates();
308: PetscMalloc(in.numberofpoints * dim * sizeof(double), &in.pointlist);
309: PetscMalloc(in.numberofpoints * sizeof(int), &in.pointmarkerlist);
310: for(Mesh::sieve_type::depthSequence::iterator v_itor = vertices->begin(); v_itor != vertices->end(); v_itor++) {
311: const Mesh::field_type::index_type& interval = coordinates->getIndex(patch, *v_itor);
312: const Mesh::field_type::value_type *array = coordinates->restrict(patch, *v_itor);
314: for(int d = 0; d < interval.index; d++) {
315: in.pointlist[interval.prefix + d] = array[d];
316: }
317: const Mesh::field_type::index_type& vInterval = vertexBundle->getIndex(patch, *v_itor);
318: in.pointmarkerlist[vInterval.prefix] = v_itor.getMarker();
319: }
320: }
322: Obj<Mesh::sieve_type::depthSequence> edges = bdTopology->depthStratum(1);
323: Obj<Mesh::bundle_type> edgeBundle = boundary->getBundle(1);
325: in.numberofsegments = edges->size();
326: if (in.numberofsegments > 0) {
327: PetscMalloc(in.numberofsegments * 2 * sizeof(int), &in.segmentlist);
328: PetscMalloc(in.numberofsegments * sizeof(int), &in.segmentmarkerlist);
329: for(Mesh::sieve_type::depthSequence::iterator e_itor = edges->begin(); e_itor != edges->end(); e_itor++) {
330: const Mesh::field_type::index_type& interval = edgeBundle->getIndex(patch, *e_itor);
331: Obj<Mesh::sieve_type::coneSequence> cone = bdTopology->cone(*e_itor);
332: int p = 0;
333:
334: for(Mesh::sieve_type::coneSequence::iterator c_itor = cone->begin(); c_itor != cone->end(); c_itor++) {
335: const Mesh::field_type::index_type& vInterval = vertexBundle->getIndex(patch, *c_itor);
337: in.segmentlist[interval.prefix * 2 + (p++)] = vInterval.prefix;
338: }
339: in.segmentmarkerlist[interval.prefix] = e_itor.getMarker();
340: }
341: }
343: in.numberofholes = 0;
344: if (in.numberofholes > 0) {
345: PetscMalloc(in.numberofholes * dim * sizeof(int), &in.holelist);
346: }
347: if (createConvexHull) {
348: args += "c";
349: }
350: triangulate((char *) args.c_str(), &in, &out, NULL);
352: PetscFree(in.pointlist);
353: PetscFree(in.pointmarkerlist);
354: PetscFree(in.segmentlist);
355: PetscFree(in.segmentmarkerlist);
356: }
357: m->populate(out.numberoftriangles, out.trianglelist, out.numberofpoints, out.pointlist, interpolate);
359: if (rank == 0) {
360: Obj<Mesh::sieve_type> topology = m->getTopology();
362: for(int v = 0; v < out.numberofpoints; v++) {
363: if (out.pointmarkerlist[v]) {
364: topology->setMarker(Mesh::point_type(0, v + out.numberoftriangles), out.pointmarkerlist[v]);
365: }
366: }
367: for(int e = 0; e < out.numberofedges; e++) {
368: if (out.edgemarkerlist[e]) {
369: Mesh::point_type endpointA(0, out.edgelist[e*2+0] + out.numberoftriangles);
370: Mesh::point_type endpointB(0, out.edgelist[e*2+1] + out.numberoftriangles);
371: Obj<ALE::def::PointSet> join = topology->nJoin(endpointA, endpointB, 1);
373: topology->setMarker(*join->begin(), out.edgemarkerlist[e]);
374: }
375: }
376: }
378: finiOutput_Triangle(&out);
379: return m;
380: };
381: #endif
382: #ifdef PETSC_HAVE_TETGEN
383: static Obj<Mesh> generate_TetGen(Obj<Mesh> boundary, bool interpolate) {
384: ::tetgenio in;
385: ::tetgenio out;
386: int dim = 3;
387: Obj<Mesh> m = Mesh(boundary->getComm(), dim);
388: Obj<Mesh::sieve_type> bdTopology = boundary->getTopology();
389: PetscMPIInt rank;
390: PetscErrorCode ierr;
392: MPI_Comm_rank(boundary->getComm(), &rank);
394: if (rank == 0) {
395: std::string args("pqenzQ");
396: bool createConvexHull = false;
397: Obj<Mesh::sieve_type::depthSequence> vertices = bdTopology->depthStratum(0);
398: Obj<Mesh::bundle_type> vertexBundle = boundary->getBundle(0);
399: Mesh::field_type::patch_type patch;
401: in.numberofpoints = vertices->size();
402: if (in.numberofpoints > 0) {
403: Obj<Mesh::field_type> coordinates = boundary->getCoordinates();
405: in.pointlist = new double[in.numberofpoints*dim];
406: in.pointmarkerlist = new int[in.numberofpoints];
407: for(Mesh::sieve_type::depthSequence::iterator v_itor = vertices->begin(); v_itor != vertices->end(); ++v_itor) {
408: const Mesh::field_type::index_type& interval = coordinates->getIndex(patch, *v_itor);
409: const Mesh::field_type::value_type *array = coordinates->restrict(patch, *v_itor);
411: for(int d = 0; d < interval.index; d++) {
412: in.pointlist[interval.prefix + d] = array[d];
413: }
414: const Mesh::field_type::index_type& vInterval = vertexBundle->getIndex(patch, *v_itor);
415: in.pointmarkerlist[vInterval.prefix] = v_itor.getMarker();
416: }
417: }
419: Obj<Mesh::sieve_type::heightSequence> facets = bdTopology->heightStratum(0);
420: Obj<Mesh::bundle_type> facetBundle = boundary->getBundle(bdTopology->depth());
422: in.numberoffacets = facets->size();
423: if (in.numberoffacets > 0) {
424: in.facetlist = new tetgenio::facet[in.numberoffacets];
425: in.facetmarkerlist = new int[in.numberoffacets];
426: for(Mesh::sieve_type::heightSequence::iterator f_itor = facets->begin(); f_itor != facets->end(); ++f_itor) {
427: const Mesh::field_type::index_type& interval = facetBundle->getIndex(patch, *f_itor);
428: Obj<Mesh::bundle_type::order_type::coneSequence> cone = vertexBundle->getPatch("element", *f_itor);
430: in.facetlist[interval.prefix].numberofpolygons = 1;
431: in.facetlist[interval.prefix].polygonlist = new tetgenio::polygon[in.facetlist[interval.prefix].numberofpolygons];
432: in.facetlist[interval.prefix].numberofholes = 0;
433: in.facetlist[interval.prefix].holelist = NULL;
435: tetgenio::polygon *poly = in.facetlist[interval.prefix].polygonlist;
436: int c = 0;
438: poly->numberofvertices = cone->size();
439: poly->vertexlist = new int[poly->numberofvertices];
440: // The "element" reorder should be fused with the structural order
441: for(Mesh::bundle_type::order_type::coneSequence::iterator c_itor = cone->begin(); c_itor != cone->end(); ++c_itor) {
442: const Mesh::field_type::index_type& vInterval = vertexBundle->getIndex(patch, *c_itor);
444: poly->vertexlist[c++] = vInterval.prefix;
445: }
446: in.facetmarkerlist[interval.prefix] = f_itor.getMarker();
447: }
448: }
450: in.numberofholes = 0;
451: if (createConvexHull) args += "c";
452: ::tetrahedralize((char *) args.c_str(), &in, &out);
453: }
454: m->populate(out.numberoftetrahedra, out.tetrahedronlist, out.numberofpoints, out.pointlist, interpolate);
455:
456: if (rank == 0) {
457: Obj<Mesh::sieve_type> topology = m->getTopology();
459: for(int v = 0; v < out.numberofpoints; v++) {
460: if (out.pointmarkerlist[v]) {
461: topology->setMarker(Mesh::point_type(0, v + out.numberoftetrahedra), out.pointmarkerlist[v]);
462: }
463: }
464: if (out.edgemarkerlist) {
465: for(int e = 0; e < out.numberofedges; e++) {
466: if (out.edgemarkerlist[e]) {
467: Mesh::point_type endpointA(0, out.edgelist[e*2+0] + out.numberoftetrahedra);
468: Mesh::point_type endpointB(0, out.edgelist[e*2+1] + out.numberoftetrahedra);
469: Obj<ALE::def::PointSet> join = topology->nJoin(endpointA, endpointB, 1);
471: topology->setMarker(*join->begin(), out.edgemarkerlist[e]);
472: }
473: }
474: }
475: if (out.trifacemarkerlist) {
476: for(int f = 0; f < out.numberoftrifaces; f++) {
477: if (out.trifacemarkerlist[f]) {
478: Obj<ALE::def::PointSet> point = ALE::def::PointSet();
479: Obj<ALE::def::PointSet> edge = ALE::def::PointSet();
480: Mesh::point_type cornerA(0, out.trifacelist[f*3+0] + out.numberoftetrahedra);
481: Mesh::point_type cornerB(0, out.trifacelist[f*3+1] + out.numberoftetrahedra);
482: Mesh::point_type cornerC(0, out.trifacelist[f*3+2] + out.numberoftetrahedra);
483: point->insert(cornerA);
484: edge->insert(cornerB);
485: edge->insert(cornerC);
486: Obj<ALE::def::PointSet> join = topology->nJoin(point, edge, 2);
488: topology->setMarker(*join->begin(), out.trifacemarkerlist[f]);
489: }
490: }
491: }
492: }
493: return m;
494: };
495: #endif
496: public:
497: static Obj<Mesh> generate(Obj<Mesh> boundary, bool interpolate = true) {
498: Obj<Mesh> mesh;
499: int dim = boundary->getDimension();
501: if (dim == 1) {
502: #ifdef PETSC_HAVE_TRIANGLE
503: mesh = generate_Triangle(boundary, interpolate);
504: #else
505: throw ALE::Exception("Mesh generation currently requires Triangle to be installed. Use --download-triangle during configure.");
506: #endif
507: } else if (dim == 2) {
508: #ifdef PETSC_HAVE_TETGEN
509: mesh = generate_TetGen(boundary, interpolate);
510: #else
511: throw ALE::Exception("Mesh generation currently requires TetGen to be installed. Use --download-tetgen during configure.");
512: #endif
513: }
514: return mesh;
515: };
516: private:
517: #ifdef PETSC_HAVE_TRIANGLE
518: static Obj<Mesh> refine_Triangle(Obj<Mesh> mesh, double maxAreas[]) {
519: struct triangulateio in;
520: struct triangulateio out;
521: int dim = 2;
522: Obj<Mesh> m = Mesh(mesh->getComm(), dim);
523: // FIX: Need to globalize
524: PetscInt numElements = mesh->getTopology()->heightStratum(0)->size();
525: PetscMPIInt rank;
526: PetscErrorCode ierr;
528: MPI_Comm_rank(mesh->getComm(), &rank);
529: initInput_Triangle(&in);
530: initOutput_Triangle(&out);
531: if (rank == 0) {
532: PetscMalloc(numElements * sizeof(double), &in.trianglearealist);
533: }
534: {
535: // Scatter in local area constraints
536: #ifdef PARALLEL
537: Vec locAreas;
538: VecScatter areaScatter;
540: VecCreateSeqWithArray(PETSC_COMM_SELF, numElements, areas, &locAreas);
541: MeshCreateMapping(oldMesh, elementBundle, partitionTypes, serialElementBundle, &areaScatter);
542: VecScatterBegin(maxAreas, locAreas, INSERT_VALUES, SCATTER_FORWARD, areaScatter);
543: VecScatterEnd(maxAreas, locAreas, INSERT_VALUES, SCATTER_FORWARD, areaScatter);
544: VecDestroy(locAreas);
545: VecScatterDestroy(areaScatter);
546: #else
547: for(int i = 0; i < numElements; i++) {
548: in.trianglearealist[i] = maxAreas[i];
549: }
550: #endif
551: }
553: #ifdef PARALLEL
554: Obj<Mesh> serialMesh = this->unify(mesh);
555: #else
556: Obj<Mesh> serialMesh = mesh;
557: #endif
558: Obj<Mesh::sieve_type> serialTopology = serialMesh->getTopology();
560: if (rank == 0) {
561: std::string args("pqenzQra");
562: Obj<Mesh::sieve_type::heightSequence> faces = serialTopology->heightStratum(0);
563: Obj<Mesh::sieve_type::depthSequence> vertices = serialTopology->depthStratum(0);
564: Obj<Mesh::bundle_type> vertexBundle = serialMesh->getBundle(0);
565: Obj<Mesh::field_type> coordinates = serialMesh->getCoordinates();
566: Mesh::field_type::patch_type patch;
567: int f = 0;
569: in.numberofpoints = vertices->size();
570: PetscMalloc(in.numberofpoints * dim * sizeof(double), &in.pointlist);
571: PetscMalloc(in.numberofpoints * sizeof(int), &in.pointmarkerlist);
572: for(Mesh::sieve_type::depthSequence::iterator v_itor = vertices->begin(); v_itor != vertices->end(); v_itor++) {
573: const Mesh::field_type::index_type& interval = coordinates->getIndex(patch, *v_itor);
574: const Mesh::field_type::value_type *array = coordinates->restrict(patch, *v_itor);
576: for(int d = 0; d < interval.index; d++) {
577: in.pointlist[interval.prefix + d] = array[d];
578: }
579: const Mesh::field_type::index_type& vInterval = vertexBundle->getIndex(patch, *v_itor);
580: in.pointmarkerlist[vInterval.prefix] = v_itor.getMarker();
581: }
583: in.numberofcorners = 3;
584: in.numberoftriangles = faces->size();
585: PetscMalloc(in.numberoftriangles * in.numberofcorners * sizeof(int), &in.trianglelist);
586: for(Mesh::sieve_type::heightSequence::iterator f_itor = faces->begin(); f_itor != faces->end(); f_itor++) {
587: Obj<Mesh::field_type::IndexArray> intervals = vertexBundle->getIndices("element", *f_itor);
588: int v = 0;
590: for(Mesh::field_type::IndexArray::iterator i_itor = intervals->begin(); i_itor != intervals->end(); i_itor++) {
591: if (i_itor->index) {
592: in.trianglelist[f * in.numberofcorners + v++] = i_itor->prefix;
593: }
594: }
595: f++;
596: }
598: Obj<Mesh::sieve_type::depthMarkerSequence> segments = serialTopology->depthStratum(1, 1);
599: Obj<Mesh::bundle_type> segmentBundle = Mesh::bundle_type();
601: segmentBundle->setTopology(serialTopology);
602: segmentBundle->setPatch(segments, patch);
603: segmentBundle->setFiberDimensionByDepth(patch, 1, 1);
604: segmentBundle->orderPatches();
605: in.numberofsegments = segments->size();
606: if (in.numberofsegments > 0) {
607: PetscMalloc(in.numberofsegments * 2 * sizeof(int), &in.segmentlist);
608: PetscMalloc(in.numberofsegments * sizeof(int), &in.segmentmarkerlist);
609: for(Mesh::sieve_type::depthMarkerSequence::iterator s_itor = segments->begin(); s_itor != segments->end(); s_itor++) {
610: const Mesh::field_type::index_type& interval = segmentBundle->getIndex(patch, *s_itor);
611: Obj<Mesh::sieve_type::coneSequence> cone = serialTopology->cone(*s_itor);
612: int p = 0;
613:
614: for(Mesh::sieve_type::coneSequence::iterator c_itor = cone->begin(); c_itor != cone->end(); c_itor++) {
615: const Mesh::field_type::index_type& vInterval = vertexBundle->getIndex(patch, *c_itor);
617: in.segmentlist[interval.prefix * 2 + (p++)] = vInterval.prefix;
618: }
619: in.segmentmarkerlist[interval.prefix] = s_itor.getMarker();
620: }
621: }
623: in.numberofholes = 0;
624: if (in.numberofholes > 0) {
625: PetscMalloc(in.numberofholes * dim * sizeof(int), &in.holelist);
626: }
627: triangulate((char *) args.c_str(), &in, &out, NULL);
628: PetscFree(in.trianglearealist);
629: PetscFree(in.pointlist);
630: PetscFree(in.pointmarkerlist);
631: PetscFree(in.segmentlist);
632: PetscFree(in.segmentmarkerlist);
633: }
634: m->populate(out.numberoftriangles, out.trianglelist, out.numberofpoints, out.pointlist);
635: //m->distribute(m);
637: // Need to make boundary
639: finiOutput_Triangle(&out);
640: return m;
641: };
642: #endif
643: #ifdef PETSC_HAVE_TETGEN
644: static Obj<Mesh> refine_TetGen(Obj<Mesh> mesh, double maxAreas[]) {
645: ::tetgenio in;
646: ::tetgenio out;
647: int dim = 3;
648: Obj<Mesh> m = Mesh(mesh->getComm(), dim);
649: // FIX: Need to globalize
650: PetscInt numElements = mesh->getTopology()->heightStratum(0)->size();
651: PetscMPIInt rank;
654: MPI_Comm_rank(mesh->getComm(), &rank);
656: if (rank == 0) {
657: in.tetrahedronvolumelist = new double[numElements];
658: }
659: {
660: // Scatter in local area constraints
661: #ifdef PARALLEL
662: Vec locAreas;
663: VecScatter areaScatter;
665: VecCreateSeqWithArray(PETSC_COMM_SELF, numElements, areas, &locAreas);
666: MeshCreateMapping(oldMesh, elementBundle, partitionTypes, serialElementBundle, &areaScatter);
667: VecScatterBegin(maxAreas, locAreas, INSERT_VALUES, SCATTER_FORWARD, areaScatter);
668: VecScatterEnd(maxAreas, locAreas, INSERT_VALUES, SCATTER_FORWARD, areaScatter);
669: VecDestroy(locAreas);
670: VecScatterDestroy(areaScatter);
671: #else
672: for(int i = 0; i < numElements; i++) {
673: in.tetrahedronvolumelist[i] = maxAreas[i];
674: }
675: #endif
676: }
678: #ifdef PARALLEL
679: Obj<Mesh> serialMesh = this->unify(mesh);
680: #else
681: Obj<Mesh> serialMesh = mesh;
682: #endif
683: Obj<Mesh::sieve_type> serialTopology = serialMesh->getTopology();
685: if (rank == 0) {
686: std::string args("qenzQra");
687: Obj<Mesh::sieve_type::heightSequence> cells = serialTopology->heightStratum(0);
688: Obj<Mesh::sieve_type::depthSequence> vertices = serialTopology->depthStratum(0);
689: Obj<Mesh::bundle_type> vertexBundle = serialMesh->getBundle(0);
690: Obj<Mesh::field_type> coordinates = serialMesh->getCoordinates();
691: Mesh::field_type::patch_type patch;
692: int c = 0;
694: in.numberofpoints = vertices->size();
695: in.pointlist = new double[in.numberofpoints*dim];
696: in.pointmarkerlist = new int[in.numberofpoints];
697: for(Mesh::sieve_type::depthSequence::iterator v_itor = vertices->begin(); v_itor != vertices->end(); ++v_itor) {
698: const Mesh::field_type::index_type& interval = coordinates->getIndex(patch, *v_itor);
699: const Mesh::field_type::value_type *array = coordinates->restrict(patch, *v_itor);
701: for(int d = 0; d < interval.index; d++) {
702: in.pointlist[interval.prefix + d] = array[d];
703: }
704: const Mesh::field_type::index_type& vInterval = vertexBundle->getIndex(patch, *v_itor);
705: in.pointmarkerlist[vInterval.prefix] = v_itor.getMarker();
706: }
708: in.numberofcorners = 4;
709: in.numberoftetrahedra = cells->size();
710: in.tetrahedronlist = new int[in.numberoftetrahedra*in.numberofcorners];
711: for(Mesh::sieve_type::heightSequence::iterator c_itor = cells->begin(); c_itor != cells->end(); ++c_itor) {
712: Obj<Mesh::field_type::IndexArray> intervals = vertexBundle->getIndices("element", *c_itor);
713: int v = 0;
715: for(Mesh::field_type::IndexArray::iterator i_itor = intervals->begin(); i_itor != intervals->end(); i_itor++) {
716: if (i_itor->index) {
717: in.tetrahedronlist[c * in.numberofcorners + v++] = i_itor->prefix;
718: }
719: }
720: c++;
721: }
723: in.numberofholes = 0;
724: ::tetrahedralize((char *) args.c_str(), &in, &out);
725: }
726: m->populate(out.numberoftetrahedra, out.tetrahedronlist, out.numberofpoints, out.pointlist);
727:
728: if (rank == 0) {
729: Obj<Mesh::sieve_type> topology = m->getTopology();
731: for(int v = 0; v < out.numberofpoints; v++) {
732: if (out.pointmarkerlist[v]) {
733: topology->setMarker(Mesh::point_type(0, v + out.numberoftetrahedra), out.pointmarkerlist[v]);
734: }
735: }
736: if (out.edgemarkerlist) {
737: for(int e = 0; e < out.numberofedges; e++) {
738: if (out.edgemarkerlist[e]) {
739: Mesh::point_type endpointA(0, out.edgelist[e*2+0] + out.numberoftetrahedra);
740: Mesh::point_type endpointB(0, out.edgelist[e*2+1] + out.numberoftetrahedra);
741: Obj<ALE::def::PointSet> join = topology->nJoin(endpointA, endpointB, 1);
743: topology->setMarker(*join->begin(), out.edgemarkerlist[e]);
744: }
745: }
746: }
747: if (out.trifacemarkerlist) {
748: for(int f = 0; f < out.numberoftrifaces; f++) {
749: if (out.trifacemarkerlist[f]) {
750: Obj<ALE::def::PointSet> point = ALE::def::PointSet();
751: Obj<ALE::def::PointSet> edge = ALE::def::PointSet();
752: Mesh::point_type cornerA(0, out.edgelist[f*3+0] + out.numberoftetrahedra);
753: Mesh::point_type cornerB(0, out.edgelist[f*3+1] + out.numberoftetrahedra);
754: Mesh::point_type cornerC(0, out.edgelist[f*3+2] + out.numberoftetrahedra);
755: point->insert(cornerA);
756: edge->insert(cornerB);
757: edge->insert(cornerC);
758: Obj<ALE::def::PointSet> join = topology->nJoin(point, edge, 2);
760: topology->setMarker(*join->begin(), out.trifacemarkerlist[f]);
761: }
762: }
763: }
764: }
765: return m;
766: };
767: #endif
768: public:
769: static Obj<Mesh> refine(Obj<Mesh> mesh, double maxArea) {
770: int numElements = mesh->getTopology()->heightStratum(0)->size();
771: double *maxAreas = new double[numElements];
772: for(int e = 0; e < numElements; e++) {
773: maxAreas[e] = maxArea;
774: }
775: Obj<Mesh> refinedMesh = refine(mesh, maxAreas);
777: delete [] maxAreas;
778: return refinedMesh;
779: };
780: static Obj<Mesh> refine(Obj<Mesh> mesh, double maxAreas[]) {
781: Obj<Mesh> refinedMesh;
782: int dim = mesh->getDimension();
784: if (dim == 2) {
785: #ifdef PETSC_HAVE_TRIANGLE
786: refinedMesh = refine_Triangle(mesh, maxAreas);
787: #else
788: throw ALE::Exception("Mesh refinement currently requires Triangle to be installed. Use --download-triangle during configure.");
789: #endif
790: } else if (dim == 3) {
791: #ifdef PETSC_HAVE_TETGEN
792: refinedMesh = refine_TetGen(mesh, maxAreas);
793: #else
794: throw ALE::Exception("Mesh generation currently requires TetGen to be installed. Use --download-tetgen during configure.");
795: #endif
796: }
797: return refinedMesh;
798: };
799: };
800: }
802: namespace def {
803: // A computational mesh
804: class Mesh {
805: public:
806: typedef Point point_type;
807: typedef Sieve<point_type,int> sieve_type;
808: typedef CoSieve<sieve_type, point_type, int, int> ordering_type;
809: typedef CoSieve<sieve_type, int, point_type, double> coordinate_type;
810: typedef CoSieve<sieve_type, int, point_type, int> bundle_type;
811: int debug;
812: private:
813: Obj<sieve_type> topology;
814: Obj<sieve_type> orientation;
815: Obj<ordering_type> ordering;
816: Obj<coordinate_type> coordinates;
817: Obj<coordinate_type> boundary;
818: std::map<std::string, Obj<coordinate_type> > fields;
819: std::map<int, Obj<bundle_type> > bundles;
820: MPI_Comm comm;
821: int dim;
822: public:
823: Mesh(MPI_Comm c, int dimension, int debug = 0) : debug(debug), comm(c), dim(dimension) {
824: this->topology = sieve_type(debug);
825: this->orientation = sieve_type(debug);
826: this->ordering = ordering_type(debug);
827: this->coordinates = coordinate_type(debug);
828: this->boundary = coordinate_type(debug);
829: };
831: MPI_Comm getComm() {return this->comm;};
832: Obj<sieve_type> getTopology() {return this->topology;};
833: Obj<sieve_type> getOrientation() {return this->orientation;};
834: Obj<ordering_type> getOrdering() {return this->ordering;};
835: int getDimension() {return this->dim;};
836: Obj<coordinate_type> getCoordinates() {return this->coordinates;};
837: Obj<coordinate_type> getBoundary() {return this->boundary;};
838: Obj<coordinate_type> getField() {return this->fields.begin()->second;};
839: Obj<coordinate_type> getField(const std::string& name) {return this->fields[name];};
840: void setField(const std::string& name, const Obj<coordinate_type>& field) {this->fields[name] = field;};
841: Obj<bundle_type> getBundle(const int dim) {
842: if (this->bundles.find(dim) == this->bundles.end()) {
843: Obj<ALE::def::Mesh::bundle_type> bundle = ALE::def::Mesh::bundle_type();
845: // Need to globalize indices (that is what we might use the value ints for)
846: std::cout << "Creating new bundle for dim " << dim << std::endl;
847: bundle->debug = this->debug;
848: bundle->setTopology(this->topology);
849: bundle->setPatch(this->topology->leaves(), 0);
850: bundle->setIndexDimensionByDepth(dim, 1);
851: bundle->orderPatches();
852: this->bundles[dim] = bundle;
853: }
854: return this->bundles[dim];
855: }
857: void buildFaces(int dim, std::map<int, int*> *curSimplex, Obj<PointSet> boundary, point_type& simplex) {
858: Obj<PointSet> faces = PointSet();
860: if (debug > 1) {std::cout << " Building faces for boundary(size " << boundary->size() << "), dim " << dim << std::endl;}
861: if (dim > 1) {
862: // Use the cone construction
863: for(PointSet::iterator b_itor = boundary->begin(); b_itor != boundary->end(); ++b_itor) {
864: Obj<PointSet> faceBoundary = PointSet();
865: point_type face;
867: faceBoundary.copy(boundary);
868: if (debug > 1) {std::cout << " boundary point " << *b_itor << std::endl;}
869: faceBoundary->erase(*b_itor);
870: this->buildFaces(dim-1, curSimplex, faceBoundary, face);
871: faces->insert(face);
872: }
873: } else {
874: if (debug > 1) {std::cout << " Just set faces to boundary in 1d" << std::endl;}
875: faces = boundary;
876: }
877: if (debug > 1) {
878: for(PointSet::iterator f_itor = faces->begin(); f_itor != faces->end(); ++f_itor) {
879: std::cout << " face point " << *f_itor << std::endl;
880: }
881: }
882: // We always create the toplevel, so we could shortcircuit somehow
883: // Should not have to loop here since the meet of just 2 boundary elements is an element
884: PointSet::iterator f_itor = faces->begin();
885: point_type start = *f_itor;
886: f_itor++;
887: point_type next = *f_itor;
888: Obj<PointSet> preElement = this->topology->nJoin(start, next, 1);
890: if (preElement->size() > 0) {
891: simplex = *preElement->begin();
892: if (debug > 1) {std::cout << " Found old simplex " << simplex << std::endl;}
893: } else {
894: simplex = point_type(0, (*(*curSimplex)[dim])++);
895: this->topology->addCone(faces, simplex);
896: if (debug > 1) {std::cout << " Added simplex " << simplex << " dim " << dim << std::endl;}
897: }
898: };
901: // Build a topology from a connectivity description
902: // (0, 0) ... (0, numSimplices-1): dim-dimensional simplices
903: // (0, numSimplices) ... (0, numVertices): vertices
904: // The other simplices are numbered as they are requested
905: void buildTopology(int numSimplices, int simplices[], int numVertices, bool interpolate = true) {
906: ALE_LOG_EVENT_BEGIN;
907: // Create a map from dimension to the current element number for that dimension
908: std::map<int,int*> curElement = std::map<int,int*>();
909: int curSimplex = 0;
910: int curVertex = numSimplices;
911: int newElement = numSimplices+numVertices;
912: Obj<PointSet> boundary = PointSet();
913: Obj<PointSet> cellTuple = PointSet();
915: curElement[0] = &curVertex;
916: curElement[dim] = &curSimplex;
917: for(int d = 1; d < dim; d++) {
918: curElement[d] = &newElement;
919: }
920: for(int s = 0; s < numSimplices; s++) {
921: point_type simplex(0, s);
923: // Build the simplex
924: if (interpolate) {
925: boundary->clear();
926: for(int b = 0; b < dim+1; b++) {
927: point_type vertex(0, simplices[s*(dim+1)+b]+numSimplices);
929: if (debug > 1) {std::cout << "Adding boundary node " << vertex << std::endl;}
930: boundary->insert(vertex);
931: }
932: if (debug) {std::cout << "simplex " << s << " boundary size " << boundary->size() << std::endl;}
934: this->buildFaces(this->dim, &curElement, boundary, simplex);
936: // Orient the simplex
937: point_type element(0, simplices[s*(dim+1)+0]+numSimplices);
938: cellTuple->clear();
939: cellTuple->insert(element);
940: for(int b = 1; b < dim+1; b++) {
941: point_type next(0, simplices[s*(dim+1)+b]+numSimplices);
942: Obj<PointSet> join = this->topology->nJoin(element, next, b);
944: if (join->size() == 0) {
945: if (debug) {std::cout << "element " << element << " next " << next << std::endl;}
946: throw ALE::Exception("Invalid join");
947: }
948: element = *join->begin();
949: cellTuple->insert(element);
950: }
951: this->orientation->addCone(cellTuple, simplex);
952: } else {
953: for(int b = 0; b < dim+1; b++) {
954: point_type p(0, simplices[s*(dim+1)+b]+numSimplices);
956: this->topology->addArrow(p, simplex);
957: this->orientation->addArrow(p, simplex, b);
958: }
959: }
960: }
961: ALE_LOG_EVENT_END;
962: };
965: void createSerialCoordinates(int embedDim, int numElements, double coords[]) {
966: ALE_LOG_EVENT_BEGIN;
967: this->topology->debug = this->debug;
968: this->coordinates->setTopology(this->topology);
969: this->coordinates->setPatch(this->topology->leaves(), 0);
970: this->coordinates->setIndexDimensionByDepth(0, embedDim);
971: this->coordinates->orderPatches();
972: Obj<sieve_type::depthSequence> vertices = this->topology->depthStratum(0);
973: for(sieve_type::depthSequence::iterator v_itor = vertices->begin(); v_itor != vertices->end(); v_itor++) {
974: this->coordinates->update(0, *v_itor, &coords[((*v_itor).index - numElements)*embedDim]);
975: }
976: ALE_LOG_EVENT_END;
977: };
978: // I think that the boundary should be marked in the Sieve
979: // It could be done with point markers like depth/height, but is that right?
980: // Should also put in boundary edges
981: void createBoundary(int numBoundaryVertices, int numBoundaryComponents, int boundaryVertices[], double boundaryValues[]) {
982: //FIX: Need to globalize
983: int numElements = this->topology->heightStratum(0)->size();
985: this->boundary->debug = this->debug;
986: this->boundary->setTopology(this->topology);
987: this->boundary->setPatch(this->topology->leaves(), 0);
988: // Reverse order allows newer conditions to override older, as required by PyLith
989: for(int v = numBoundaryVertices-1; v >= 0; v--) {
990: point_type vertex(0, boundaryVertices[v*(numBoundaryComponents+1)] + numElements);
992: if (this->boundary->getIndexDimension(0, vertex) == 0) {
993: for(int c = 0; c < numBoundaryComponents; c++) {
994: if (boundaryVertices[v*(numBoundaryComponents+1)+c+1]) {
995: this->boundary->setIndexDimension(0, vertex, c+1, 1);
996: }
997: }
998: }
999: }
1000: this->boundary->orderPatches();
1001: for(int v = 0; v < numBoundaryVertices; v++) {
1002: point_type vertex(0, boundaryVertices[v*(numBoundaryComponents+1)] + numElements);
1004: this->boundary->update(0, vertex, &boundaryValues[v*numBoundaryComponents]);
1005: }
1006: };
1007: void populate(int numSimplices, int simplices[], int numVertices, double coords[], bool interpolate = true) {
1008: PetscMPIInt rank;
1010: MPI_Comm_rank(this->comm, &rank);
1011: /* Create serial sieve */
1012: this->topology->debug = this->debug;
1013: this->topology->setStratification(false);
1014: if (rank == 0) {
1015: this->buildTopology(numSimplices, simplices, numVertices, interpolate);
1016: }
1017: this->topology->stratify();
1018: this->topology->setStratification(true);
1019: this->createSerialCoordinates(this->dim, numSimplices, coords);
1020: };
1021: };
1023: // Creation
1024: class PyLithBuilder {
1025: static inline void ignoreComments(char *buf, PetscInt bufSize, FILE *f) {
1026: while((fgets(buf, bufSize, f) != NULL) && (buf[0] == '#')) {}
1027: };
1029: static void readConnectivity(MPI_Comm comm, const std::string& filename, int dim, bool useZeroBase, int& numElements, int *vertices[]) {
1030: PetscViewer viewer;
1031: FILE *f;
1032: PetscInt maxCells = 1024, cellCount = 0;
1033: PetscInt *verts;
1034: char buf[2048];
1035: PetscInt c;
1036: PetscInt commRank;
1039: MPI_Comm_rank(comm, &commRank);
1040: if (dim != 3) {
1041: throw ALE::Exception("PyLith only works in 3D");
1042: }
1043: if (commRank != 0) return;
1044: PetscViewerCreate(PETSC_COMM_SELF, &viewer);
1045: PetscViewerSetType(viewer, PETSC_VIEWER_ASCII);
1046: PetscViewerFileSetMode(viewer, FILE_MODE_READ);
1047: PetscViewerFileSetName(viewer, filename.c_str());
1048: PetscViewerASCIIGetPointer(viewer, &f);
1049: /* Ignore comments */
1050: ignoreComments(buf, 2048, f);
1051: PetscMalloc(maxCells*(dim+1) * sizeof(PetscInt), &verts);
1052: do {
1053: const char *v = strtok(buf, " ");
1054: int elementType;
1056: if (cellCount == maxCells) {
1057: PetscInt *vtmp;
1059: vtmp = verts;
1060: PetscMalloc(maxCells*2*(dim+1) * sizeof(PetscInt), &verts);
1061: PetscMemcpy(verts, vtmp, maxCells*(dim+1) * sizeof(PetscInt));
1062: PetscFree(vtmp);
1063: maxCells *= 2;
1064: }
1065: /* Ignore cell number */
1066: v = strtok(NULL, " ");
1067: /* Verify element type is linear tetrahedron */
1068: elementType = atoi(v);
1069: if (elementType != 5) {
1070: throw ALE::Exception("We only accept linear tetrahedra right now");
1071: }
1072: v = strtok(NULL, " ");
1073: /* Ignore material type */
1074: v = strtok(NULL, " ");
1075: /* Ignore infinite domain element code */
1076: v = strtok(NULL, " ");
1077: for(c = 0; c <= dim; c++) {
1078: int vertex = atoi(v);
1079:
1080: if (!useZeroBase) vertex -= 1;
1081: verts[cellCount*(dim+1)+c] = vertex;
1082: v = strtok(NULL, " ");
1083: }
1084: cellCount++;
1085: } while(fgets(buf, 2048, f) != NULL);
1086: PetscViewerDestroy(viewer);
1087: numElements = cellCount;
1088: *vertices = verts;
1089: };
1090: static void readCoordinates(MPI_Comm comm, const std::string& filename, int dim, int& numVertices, double *coordinates[]) {
1091: PetscViewer viewer;
1092: FILE *f;
1093: PetscInt maxVerts = 1024, vertexCount = 0;
1094: PetscScalar *coords;
1095: char buf[2048];
1096: PetscInt c;
1097: PetscInt commRank;
1100: MPI_Comm_rank(comm, &commRank);
1101: if (commRank == 0) {
1102: PetscViewerCreate(PETSC_COMM_SELF, &viewer);
1103: PetscViewerSetType(viewer, PETSC_VIEWER_ASCII);
1104: PetscViewerFileSetMode(viewer, FILE_MODE_READ);
1105: PetscViewerFileSetName(viewer, filename.c_str());
1106: PetscViewerASCIIGetPointer(viewer, &f);
1107: /* Ignore comments and units line */
1108: ignoreComments(buf, 2048, f);
1109: PetscMalloc(maxVerts*dim * sizeof(PetscScalar), &coords);
1110: /* Ignore comments */
1111: ignoreComments(buf, 2048, f);
1112: do {
1113: const char *x = strtok(buf, " ");
1115: if (vertexCount == maxVerts) {
1116: PetscScalar *ctmp;
1118: ctmp = coords;
1119: PetscMalloc(maxVerts*2*dim * sizeof(PetscScalar), &coords);
1120: PetscMemcpy(coords, ctmp, maxVerts*dim * sizeof(PetscScalar));
1121: PetscFree(ctmp);
1122: maxVerts *= 2;
1123: }
1124: /* Ignore vertex number */
1125: x = strtok(NULL, " ");
1126: for(c = 0; c < dim; c++) {
1127: coords[vertexCount*dim+c] = atof(x);
1128: x = strtok(NULL, " ");
1129: }
1130: vertexCount++;
1131: } while(fgets(buf, 2048, f) != NULL);
1132: PetscViewerDestroy(viewer);
1133: numVertices = vertexCount;
1134: *coordinates = coords;
1135: }
1136: };
1137: public:
1138: PyLithBuilder() {};
1139: virtual ~PyLithBuilder() {};
1141: static Obj<Mesh> create(MPI_Comm comm, const std::string& baseFilename, bool interpolate = true, int debug = 0) {
1142: int dim = 3;
1143: bool useZeroBase = false;
1144: Obj<Mesh> mesh = Mesh(comm, dim);
1145: int *vertices;
1146: double *coordinates;
1147: int numElements, numVertices;
1149: mesh->debug = debug;
1150: readConnectivity(comm, baseFilename+".connect", dim, useZeroBase, numElements, &vertices);
1151: readCoordinates(comm, baseFilename+".coord", dim, numVertices, &coordinates);
1152: mesh->populate(numElements, vertices, numVertices, coordinates, interpolate);
1153: return mesh;
1154: };
1156: static Obj<ALE::Two::Mesh> createNew(MPI_Comm comm, const std::string& baseFilename, bool interpolate = true, int debug = 0) {
1157: int dim = 3;
1158: bool useZeroBase = false;
1159: Obj<ALE::Two::Mesh> mesh = ALE::Two::Mesh(comm, dim);
1160: int *vertices;
1161: double *coordinates;
1162: int numElements, numVertices;
1164: mesh->debug = debug;
1165: readConnectivity(comm, baseFilename+".connect", dim, useZeroBase, numElements, &vertices);
1166: readCoordinates(comm, baseFilename+".coord", dim, numVertices, &coordinates);
1167: mesh->populate(numElements, vertices, numVertices, coordinates, interpolate);
1168: return mesh;
1169: };
1170: };
1172: class PCICEBuilder {
1173: static void readConnectivity(MPI_Comm comm, const std::string& filename, int dim, bool useZeroBase, int& numElements, int *vertices[]) {
1174: PetscViewer viewer;
1175: FILE *f;
1176: PetscInt numCells, cellCount = 0;
1177: PetscInt *verts;
1178: char buf[2048];
1179: PetscInt c;
1180: PetscInt commRank;
1183: MPI_Comm_rank(comm, &commRank);
1185: if (commRank != 0) return;
1186: PetscViewerCreate(PETSC_COMM_SELF, &viewer);
1187: PetscViewerSetType(viewer, PETSC_VIEWER_ASCII);
1188: PetscViewerFileSetMode(viewer, FILE_MODE_READ);
1189: PetscViewerFileSetName(viewer, filename.c_str());
1190: PetscViewerASCIIGetPointer(viewer, &f);
1191: numCells = atoi(fgets(buf, 2048, f));
1192: PetscMalloc(numCells*(dim+1) * sizeof(PetscInt), &verts);
1193: while(fgets(buf, 2048, f) != NULL) {
1194: const char *v = strtok(buf, " ");
1195:
1196: /* Ignore cell number */
1197: v = strtok(NULL, " ");
1198: for(c = 0; c <= dim; c++) {
1199: int vertex = atoi(v);
1200:
1201: if (!useZeroBase) vertex -= 1;
1202: verts[cellCount*(dim+1)+c] = vertex;
1203: v = strtok(NULL, " ");
1204: }
1205: cellCount++;
1206: }
1207: PetscViewerDestroy(viewer);
1208: numElements = numCells;
1209: *vertices = verts;
1210: };
1211: static void readCoordinates(MPI_Comm comm, const std::string& filename, int dim, int& numVertices, double *coordinates[]) {
1212: PetscViewer viewer;
1213: FILE *f;
1214: PetscInt numVerts, vertexCount = 0;
1215: PetscScalar *coords;
1216: char buf[2048];
1217: PetscInt c;
1218: PetscInt commRank;
1221: MPI_Comm_rank(comm, &commRank);
1223: if (commRank != 0) return;
1224: PetscViewerCreate(PETSC_COMM_SELF, &viewer);
1225: PetscViewerSetType(viewer, PETSC_VIEWER_ASCII);
1226: PetscViewerFileSetMode(viewer, FILE_MODE_READ);
1227: PetscViewerFileSetName(viewer, filename.c_str());
1228: PetscViewerASCIIGetPointer(viewer, &f);
1229: numVerts = atoi(fgets(buf, 2048, f));
1230: PetscMalloc(numVerts*dim * sizeof(PetscScalar), &coords);
1231: while(fgets(buf, 2048, f) != NULL) {
1232: const char *x = strtok(buf, " ");
1233:
1234: /* Ignore vertex number */
1235: x = strtok(NULL, " ");
1236: for(c = 0; c < dim; c++) {
1237: coords[vertexCount*dim+c] = atof(x);
1238: x = strtok(NULL, " ");
1239: }
1240: vertexCount++;
1241: }
1242: PetscViewerDestroy(viewer);
1243: numVertices = numVerts;
1244: *coordinates = coords;
1245: };
1246: public:
1247: PCICEBuilder() {};
1248: virtual ~PCICEBuilder() {};
1250: static Obj<Mesh> create(MPI_Comm comm, const std::string& baseFilename, int dim, bool useZeroBase = false, int debug = 0) {
1251: Obj<Mesh> mesh = Mesh(comm, dim);
1252: int *vertices;
1253: double *coordinates;
1254: int numElements, numVertices;
1256: mesh->debug = debug;
1257: readConnectivity(comm, baseFilename+".lcon", dim, useZeroBase, numElements, &vertices);
1258: readCoordinates(comm, baseFilename+".nodes", dim, numVertices, &coordinates);
1259: mesh->populate(numElements, vertices, numVertices, coordinates);
1260: return mesh;
1261: };
1263: static Obj<ALE::Two::Mesh> createNew(MPI_Comm comm, const std::string& baseFilename, int dim, bool useZeroBase = false, int debug = 0) {
1264: Obj<ALE::Two::Mesh> mesh = ALE::Two::Mesh(comm, dim, debug);
1265: int *vertices;
1266: double *coordinates;
1267: int numElements, numVertices;
1269: readConnectivity(comm, baseFilename+".lcon", dim, useZeroBase, numElements, &vertices);
1270: readCoordinates(comm, baseFilename+".nodes", dim, numVertices, &coordinates);
1271: mesh->populate(numElements, vertices, numVertices, coordinates);
1272: return mesh;
1273: };
1274: };
1276: class Generator {
1277: #ifdef PETSC_HAVE_TRIANGLE
1278: static void initInput_Triangle(struct triangulateio *inputCtx) {
1279: inputCtx->numberofpoints = 0;
1280: inputCtx->numberofpointattributes = 0;
1281: inputCtx->pointlist = NULL;
1282: inputCtx->pointattributelist = NULL;
1283: inputCtx->pointmarkerlist = NULL;
1284: inputCtx->numberofsegments = 0;
1285: inputCtx->segmentlist = NULL;
1286: inputCtx->segmentmarkerlist = NULL;
1287: inputCtx->numberoftriangleattributes = 0;
1288: inputCtx->numberofholes = 0;
1289: inputCtx->holelist = NULL;
1290: inputCtx->numberofregions = 0;
1291: inputCtx->regionlist = NULL;
1292: };
1293: static void initOutput_Triangle(struct triangulateio *outputCtx) {
1294: outputCtx->pointlist = NULL;
1295: outputCtx->pointattributelist = NULL;
1296: outputCtx->pointmarkerlist = NULL;
1297: outputCtx->trianglelist = NULL;
1298: outputCtx->triangleattributelist = NULL;
1299: outputCtx->neighborlist = NULL;
1300: outputCtx->segmentlist = NULL;
1301: outputCtx->segmentmarkerlist = NULL;
1302: outputCtx->edgelist = NULL;
1303: outputCtx->edgemarkerlist = NULL;
1304: };
1305: static void finiOutput_Triangle(struct triangulateio *outputCtx) {
1306: free(outputCtx->pointmarkerlist);
1307: free(outputCtx->edgelist);
1308: free(outputCtx->edgemarkerlist);
1309: free(outputCtx->trianglelist);
1310: free(outputCtx->neighborlist);
1311: };
1314: static Obj<Mesh> generate_Triangle(Obj<Mesh> boundary) {
1315: struct triangulateio in;
1316: struct triangulateio out;
1317: int dim = 2;
1318: Obj<Mesh> m = Mesh(boundary->getComm(), dim);
1319: Obj<Mesh::sieve_type> bdTopology = boundary->getTopology();
1320: Obj<Mesh::sieve_type> bdOrientation = boundary->getOrientation();
1321: PetscMPIInt rank;
1322: PetscErrorCode ierr;
1324: MPI_Comm_rank(boundary->getComm(), &rank);
1325: initInput_Triangle(&in);
1326: initOutput_Triangle(&out);
1327: if (rank == 0) {
1328: std::string args("pqenzQ");
1329: bool createConvexHull = false;
1330: Obj<Mesh::sieve_type::depthSequence> vertices = bdTopology->depthStratum(0);
1331: Obj<Mesh::bundle_type> vertexBundle = boundary->getBundle(0);
1333: in.numberofpoints = vertices->size();
1334: if (in.numberofpoints > 0) {
1335: Obj<Mesh::coordinate_type> coordinates = boundary->getCoordinates();
1337: PetscMalloc(in.numberofpoints * dim * sizeof(double), &in.pointlist);
1338: PetscMalloc(in.numberofpoints * sizeof(int), &in.pointmarkerlist);
1339: for(Mesh::sieve_type::depthSequence::iterator v_itor = vertices->begin(); v_itor != vertices->end(); v_itor++) {
1340: const Mesh::coordinate_type::index_type& interval = coordinates->getIndex(0, *v_itor);
1341: const Mesh::coordinate_type::value_type *array = coordinates->restrict(0, *v_itor);
1343: for(int d = 0; d < interval.index; d++) {
1344: in.pointlist[interval.prefix + d] = array[d];
1345: }
1346: const Mesh::coordinate_type::index_type& vInterval = vertexBundle->getIndex(0, *v_itor);
1347: in.pointmarkerlist[vInterval.prefix] = v_itor.getMarker();
1348: }
1349: }
1351: Obj<Mesh::sieve_type::depthSequence> edges = bdTopology->depthStratum(1);
1352: Obj<Mesh::bundle_type> edgeBundle = boundary->getBundle(1);
1354: in.numberofsegments = edges->size();
1355: if (in.numberofsegments > 0) {
1356: PetscMalloc(in.numberofsegments * 2 * sizeof(int), &in.segmentlist);
1357: PetscMalloc(in.numberofsegments * sizeof(int), &in.segmentmarkerlist);
1358: for(Mesh::sieve_type::depthSequence::iterator e_itor = edges->begin(); e_itor != edges->end(); e_itor++) {
1359: const Mesh::coordinate_type::index_type& interval = edgeBundle->getIndex(0, *e_itor);
1360: Obj<Mesh::sieve_type::coneSequence> cone = bdTopology->cone(*e_itor);
1361: int p = 0;
1362:
1363: for(Mesh::sieve_type::coneSequence::iterator c_itor = cone->begin(); c_itor != cone->end(); c_itor++) {
1364: const Mesh::coordinate_type::index_type& vInterval = vertexBundle->getIndex(0, *c_itor);
1366: in.segmentlist[interval.prefix * 2 + (p++)] = vInterval.prefix;
1367: }
1368: in.segmentmarkerlist[interval.prefix] = e_itor.getMarker();
1369: }
1370: }
1372: in.numberofholes = 0;
1373: if (in.numberofholes > 0) {
1374: PetscMalloc(in.numberofholes * dim * sizeof(int), &in.holelist);
1375: }
1376: if (createConvexHull) {
1377: args += "c";
1378: }
1379: triangulate((char *) args.c_str(), &in, &out, NULL);
1381: PetscFree(in.pointlist);
1382: PetscFree(in.pointmarkerlist);
1383: PetscFree(in.segmentlist);
1384: PetscFree(in.segmentmarkerlist);
1385: }
1386: m->populate(out.numberoftriangles, out.trianglelist, out.numberofpoints, out.pointlist);
1388: if (rank == 0) {
1389: Obj<Mesh::sieve_type> topology = m->getTopology();
1391: for(int v = 0; v < out.numberofpoints; v++) {
1392: if (out.pointmarkerlist[v]) {
1393: topology->setMarker(Mesh::point_type(0, v + out.numberoftriangles), out.pointmarkerlist[v]);
1394: }
1395: }
1396: for(int e = 0; e < out.numberofedges; e++) {
1397: if (out.edgemarkerlist[e]) {
1398: Mesh::point_type endpointA(0, out.edgelist[e*2+0] + out.numberoftriangles);
1399: Mesh::point_type endpointB(0, out.edgelist[e*2+1] + out.numberoftriangles);
1400: Obj<PointSet> join = topology->nJoin(endpointA, endpointB, 1);
1402: topology->setMarker(*join->begin(), out.edgemarkerlist[e]);
1403: }
1404: }
1405: }
1407: finiOutput_Triangle(&out);
1408: return m;
1409: };
1410: #endif
1411: #ifdef PETSC_HAVE_TETGEN
1414: static Obj<Mesh> generate_TetGen(Obj<Mesh> boundary) {
1415: ::tetgenio in;
1416: ::tetgenio out;
1417: int dim = 3;
1418: Obj<Mesh> m = Mesh(boundary->getComm(), dim);
1419: Obj<Mesh::sieve_type> bdTopology = boundary->getTopology();
1420: Obj<Mesh::sieve_type> bdOrientation = boundary->getOrientation();
1421: Obj<Mesh::ordering_type> bdOrdering = boundary->getOrdering();
1422: PetscMPIInt rank;
1423: PetscErrorCode ierr;
1425: MPI_Comm_rank(boundary->getComm(), &rank);
1427: if (rank == 0) {
1428: std::string args("pqenzQ");
1429: bool createConvexHull = false;
1430: Obj<Mesh::sieve_type::depthSequence> vertices = bdTopology->depthStratum(0);
1431: Obj<Mesh::bundle_type> vertexBundle = boundary->getBundle(0);
1433: in.numberofpoints = vertices->size();
1434: if (in.numberofpoints > 0) {
1435: Obj<Mesh::coordinate_type> coordinates = boundary->getCoordinates();
1437: in.pointlist = new double[in.numberofpoints*dim];
1438: in.pointmarkerlist = new int[in.numberofpoints];
1439: for(Mesh::sieve_type::depthSequence::iterator v_itor = vertices->begin(); v_itor != vertices->end(); ++v_itor) {
1440: const Mesh::coordinate_type::index_type& interval = coordinates->getIndex(0, *v_itor);
1441: const Mesh::coordinate_type::value_type *array = coordinates->restrict(0, *v_itor);
1443: for(int d = 0; d < interval.index; d++) {
1444: in.pointlist[interval.prefix + d] = array[d];
1445: }
1446: const Mesh::coordinate_type::index_type& vInterval = vertexBundle->getIndex(0, *v_itor);
1447: in.pointmarkerlist[vInterval.prefix] = v_itor.getMarker();
1448: }
1449: }
1451: Obj<Mesh::sieve_type::heightSequence> facets = bdTopology->heightStratum(0);
1452: Obj<Mesh::bundle_type> facetBundle = boundary->getBundle(bdTopology->depth());
1454: in.numberoffacets = facets->size();
1455: if (in.numberoffacets > 0) {
1456: in.facetlist = new tetgenio::facet[in.numberoffacets];
1457: in.facetmarkerlist = new int[in.numberoffacets];
1458: for(Mesh::sieve_type::heightSequence::iterator f_itor = facets->begin(); f_itor != facets->end(); ++f_itor) {
1459: const Mesh::coordinate_type::index_type& interval = facetBundle->getIndex(0, *f_itor);
1460: Obj<Mesh::ordering_type::patches_type::coneSequence> cone = bdOrdering->getPatch(*f_itor);
1462: in.facetlist[interval.prefix].numberofpolygons = 1;
1463: in.facetlist[interval.prefix].polygonlist = new tetgenio::polygon[in.facetlist[interval.prefix].numberofpolygons];
1464: in.facetlist[interval.prefix].numberofholes = 0;
1465: in.facetlist[interval.prefix].holelist = NULL;
1467: tetgenio::polygon *poly = in.facetlist[interval.prefix].polygonlist;
1468: int c = 0;
1470: poly->numberofvertices = cone->size();
1471: poly->vertexlist = new int[poly->numberofvertices];
1472: for(Mesh::ordering_type::patches_type::coneSequence::iterator c_itor = cone->begin(); c_itor != cone->end(); ++c_itor) {
1473: const Mesh::coordinate_type::index_type& vInterval = vertexBundle->getIndex(0, *c_itor);
1475: poly->vertexlist[c++] = vInterval.prefix;
1476: }
1477: in.facetmarkerlist[interval.prefix] = f_itor.getMarker();
1478: }
1479: }
1481: in.numberofholes = 0;
1482: if (createConvexHull) args += "c";
1483: ::tetrahedralize((char *) args.c_str(), &in, &out);
1484: }
1485: m->populate(out.numberoftetrahedra, out.tetrahedronlist, out.numberofpoints, out.pointlist);
1486:
1487: if (rank == 0) {
1488: Obj<Mesh::sieve_type> topology = m->getTopology();
1490: for(int v = 0; v < out.numberofpoints; v++) {
1491: if (out.pointmarkerlist[v]) {
1492: topology->setMarker(Mesh::point_type(0, v + out.numberoftetrahedra), out.pointmarkerlist[v]);
1493: }
1494: }
1495: if (out.edgemarkerlist) {
1496: for(int e = 0; e < out.numberofedges; e++) {
1497: if (out.edgemarkerlist[e]) {
1498: Mesh::point_type endpointA(0, out.edgelist[e*2+0] + out.numberoftetrahedra);
1499: Mesh::point_type endpointB(0, out.edgelist[e*2+1] + out.numberoftetrahedra);
1500: Obj<PointSet> join = topology->nJoin(endpointA, endpointB, 1);
1502: topology->setMarker(*join->begin(), out.edgemarkerlist[e]);
1503: }
1504: }
1505: }
1506: if (out.trifacemarkerlist) {
1507: for(int f = 0; f < out.numberoftrifaces; f++) {
1508: if (out.trifacemarkerlist[f]) {
1509: Obj<PointSet> point = PointSet();
1510: Obj<PointSet> edge = PointSet();
1511: Mesh::point_type cornerA(0, out.trifacelist[f*3+0] + out.numberoftetrahedra);
1512: Mesh::point_type cornerB(0, out.trifacelist[f*3+1] + out.numberoftetrahedra);
1513: Mesh::point_type cornerC(0, out.trifacelist[f*3+2] + out.numberoftetrahedra);
1514: point->insert(cornerA);
1515: edge->insert(cornerB);
1516: edge->insert(cornerC);
1517: Obj<PointSet> join = topology->nJoin(point, edge, 2);
1519: topology->setMarker(*join->begin(), out.trifacemarkerlist[f]);
1520: }
1521: }
1522: }
1523: }
1524: return m;
1525: };
1526: #endif
1527: public:
1528: static Obj<Mesh> generate(Obj<Mesh> boundary) {
1529: Obj<Mesh> mesh;
1530: int dim = boundary->getDimension();
1532: if (dim == 1) {
1533: #ifdef PETSC_HAVE_TRIANGLE
1534: mesh = generate_Triangle(boundary);
1535: #else
1536: throw ALE::Exception("Mesh generation currently requires Triangle to be installed. Use --download-triangle during configure.");
1537: #endif
1538: } else if (dim == 2) {
1539: #ifdef PETSC_HAVE_TETGEN
1540: mesh = generate_TetGen(boundary);
1541: #else
1542: throw ALE::Exception("Mesh generation currently requires TetGen to be installed. Use --download-tetgen during configure.");
1543: #endif
1544: }
1545: return mesh;
1546: };
1547: private:
1548: #ifdef PETSC_HAVE_TRIANGLE
1549: static Obj<Mesh> refine_Triangle(Obj<Mesh> mesh, double maxAreas[]) {
1550: struct triangulateio in;
1551: struct triangulateio out;
1552: int dim = 2;
1553: Obj<Mesh> m = Mesh(mesh->getComm(), dim);
1554: // FIX: Need to globalize
1555: PetscInt numElements = mesh->getTopology()->heightStratum(0)->size();
1556: PetscMPIInt rank;
1557: PetscErrorCode ierr;
1559: MPI_Comm_rank(mesh->getComm(), &rank);
1560: initInput_Triangle(&in);
1561: initOutput_Triangle(&out);
1562: if (rank == 0) {
1563: PetscMalloc(numElements * sizeof(double), &in.trianglearealist);
1564: }
1565: {
1566: // Scatter in local area constraints
1567: #ifdef PARALLEL
1568: Vec locAreas;
1569: VecScatter areaScatter;
1571: VecCreateSeqWithArray(PETSC_COMM_SELF, numElements, areas, &locAreas);
1572: MeshCreateMapping(oldMesh, elementBundle, partitionTypes, serialElementBundle, &areaScatter);
1573: VecScatterBegin(maxAreas, locAreas, INSERT_VALUES, SCATTER_FORWARD, areaScatter);
1574: VecScatterEnd(maxAreas, locAreas, INSERT_VALUES, SCATTER_FORWARD, areaScatter);
1575: VecDestroy(locAreas);
1576: VecScatterDestroy(areaScatter);
1577: #else
1578: for(int i = 0; i < numElements; i++) {
1579: in.trianglearealist[i] = maxAreas[i];
1580: }
1581: #endif
1582: }
1584: #ifdef PARALLEL
1585: Obj<Mesh> serialMesh = this->unify(mesh);
1586: #else
1587: Obj<Mesh> serialMesh = mesh;
1588: #endif
1589: Obj<Mesh::sieve_type> serialTopology = serialMesh->getTopology();
1590: Obj<Mesh::sieve_type> serialOrientation = serialMesh->getOrientation();
1592: if (rank == 0) {
1593: std::string args("pqenzQra");
1594: Obj<Mesh::sieve_type::heightSequence> faces = serialTopology->heightStratum(0);
1595: Obj<Mesh::sieve_type::depthSequence> vertices = serialTopology->depthStratum(0);
1596: Obj<Mesh::bundle_type> vertexBundle = serialMesh->getBundle(0);
1597: Obj<Mesh::coordinate_type> coordinates = serialMesh->getCoordinates();
1598: int f = 0;
1600: in.numberofpoints = vertices->size();
1601: PetscMalloc(in.numberofpoints * dim * sizeof(double), &in.pointlist);
1602: PetscMalloc(in.numberofpoints * sizeof(int), &in.pointmarkerlist);
1603: for(Mesh::sieve_type::depthSequence::iterator v_itor = vertices->begin(); v_itor != vertices->end(); v_itor++) {
1604: const Mesh::coordinate_type::index_type& interval = coordinates->getIndex(0, *v_itor);
1605: const Mesh::coordinate_type::value_type *array = coordinates->restrict(0, *v_itor);
1607: for(int d = 0; d < interval.index; d++) {
1608: in.pointlist[interval.prefix + d] = array[d];
1609: }
1610: const Mesh::coordinate_type::index_type& vInterval = vertexBundle->getIndex(0, *v_itor);
1611: in.pointmarkerlist[vInterval.prefix] = v_itor.getMarker();
1612: }
1614: in.numberofcorners = 3;
1615: in.numberoftriangles = faces->size();
1616: PetscMalloc(in.numberoftriangles * in.numberofcorners * sizeof(int), &in.trianglelist);
1617: for(Mesh::sieve_type::heightSequence::iterator f_itor = faces->begin(); f_itor != faces->end(); f_itor++) {
1618: Obj<Mesh::coordinate_type::IndexArray> intervals = vertexBundle->getOrderedIndices(0, serialOrientation->cone(*f_itor));
1619: int v = 0;
1621: for(Mesh::coordinate_type::IndexArray::iterator i_itor = intervals->begin(); i_itor != intervals->end(); i_itor++) {
1622: in.trianglelist[f * in.numberofcorners + v++] = i_itor->prefix;
1623: }
1624: f++;
1625: }
1627: Obj<Mesh::sieve_type::depthMarkerSequence> segments = serialTopology->depthStratum(1, 1);
1628: Obj<Mesh::bundle_type> segmentBundle = Mesh::bundle_type();
1630: segmentBundle->setTopology(serialTopology);
1631: segmentBundle->setPatch(segments, 0);
1632: segmentBundle->setIndexDimensionByDepth(1, 1);
1633: segmentBundle->orderPatches();
1634: in.numberofsegments = segments->size();
1635: if (in.numberofsegments > 0) {
1636: PetscMalloc(in.numberofsegments * 2 * sizeof(int), &in.segmentlist);
1637: PetscMalloc(in.numberofsegments * sizeof(int), &in.segmentmarkerlist);
1638: for(Mesh::sieve_type::depthMarkerSequence::iterator s_itor = segments->begin(); s_itor != segments->end(); s_itor++) {
1639: const Mesh::coordinate_type::index_type& interval = segmentBundle->getIndex(0, *s_itor);
1640: Obj<Mesh::sieve_type::coneSequence> cone = serialTopology->cone(*s_itor);
1641: int p = 0;
1642:
1643: for(Mesh::sieve_type::coneSequence::iterator c_itor = cone->begin(); c_itor != cone->end(); c_itor++) {
1644: const Mesh::coordinate_type::index_type& vInterval = vertexBundle->getIndex(0, *c_itor);
1646: in.segmentlist[interval.prefix * 2 + (p++)] = vInterval.prefix;
1647: }
1648: in.segmentmarkerlist[interval.prefix] = s_itor.getMarker();
1649: }
1650: }
1652: in.numberofholes = 0;
1653: if (in.numberofholes > 0) {
1654: PetscMalloc(in.numberofholes * dim * sizeof(int), &in.holelist);
1655: }
1656: triangulate((char *) args.c_str(), &in, &out, NULL);
1657: PetscFree(in.trianglearealist);
1658: PetscFree(in.pointlist);
1659: PetscFree(in.pointmarkerlist);
1660: PetscFree(in.segmentlist);
1661: PetscFree(in.segmentmarkerlist);
1662: }
1663: m->populate(out.numberoftriangles, out.trianglelist, out.numberofpoints, out.pointlist);
1664: //m->distribute(m);
1666: // Need to make boundary
1668: finiOutput_Triangle(&out);
1669: return m;
1670: };
1671: #endif
1672: #ifdef PETSC_HAVE_TETGEN
1673: static Obj<Mesh> refine_TetGen(Obj<Mesh> mesh, double maxAreas[]) {
1674: ::tetgenio in;
1675: ::tetgenio out;
1676: int dim = 3;
1677: Obj<Mesh> m = Mesh(mesh->getComm(), dim);
1678: // FIX: Need to globalize
1679: PetscInt numElements = mesh->getTopology()->heightStratum(0)->size();
1680: PetscMPIInt rank;
1683: MPI_Comm_rank(mesh->getComm(), &rank);
1685: if (rank == 0) {
1686: in.tetrahedronvolumelist = new double[numElements];
1687: }
1688: {
1689: // Scatter in local area constraints
1690: #ifdef PARALLEL
1691: Vec locAreas;
1692: VecScatter areaScatter;
1694: VecCreateSeqWithArray(PETSC_COMM_SELF, numElements, areas, &locAreas);
1695: MeshCreateMapping(oldMesh, elementBundle, partitionTypes, serialElementBundle, &areaScatter);
1696: VecScatterBegin(maxAreas, locAreas, INSERT_VALUES, SCATTER_FORWARD, areaScatter);
1697: VecScatterEnd(maxAreas, locAreas, INSERT_VALUES, SCATTER_FORWARD, areaScatter);
1698: VecDestroy(locAreas);
1699: VecScatterDestroy(areaScatter);
1700: #else
1701: for(int i = 0; i < numElements; i++) {
1702: in.tetrahedronvolumelist[i] = maxAreas[i];
1703: }
1704: #endif
1705: }
1707: #ifdef PARALLEL
1708: Obj<Mesh> serialMesh = this->unify(mesh);
1709: #else
1710: Obj<Mesh> serialMesh = mesh;
1711: #endif
1712: Obj<Mesh::sieve_type> serialTopology = serialMesh->getTopology();
1713: Obj<Mesh::sieve_type> serialOrientation = serialMesh->getOrientation();
1714: Obj<Mesh::ordering_type> serialOrdering = serialMesh->getOrdering();
1716: if (rank == 0) {
1717: std::string args("qenzQra");
1718: Obj<Mesh::sieve_type::heightSequence> cells = serialTopology->heightStratum(0);
1719: Obj<Mesh::sieve_type::depthSequence> vertices = serialTopology->depthStratum(0);
1720: Obj<Mesh::bundle_type> vertexBundle = serialMesh->getBundle(0);
1721: Obj<Mesh::coordinate_type> coordinates = serialMesh->getCoordinates();
1722: int c = 0;
1724: in.numberofpoints = vertices->size();
1725: in.pointlist = new double[in.numberofpoints*dim];
1726: in.pointmarkerlist = new int[in.numberofpoints];
1727: for(Mesh::sieve_type::depthSequence::iterator v_itor = vertices->begin(); v_itor != vertices->end(); ++v_itor) {
1728: const Mesh::coordinate_type::index_type& interval = coordinates->getIndex(0, *v_itor);
1729: const Mesh::coordinate_type::value_type *array = coordinates->restrict(0, *v_itor);
1731: for(int d = 0; d < interval.index; d++) {
1732: in.pointlist[interval.prefix + d] = array[d];
1733: }
1734: const Mesh::coordinate_type::index_type& vInterval = vertexBundle->getIndex(0, *v_itor);
1735: in.pointmarkerlist[vInterval.prefix] = v_itor.getMarker();
1736: }
1738: in.numberofcorners = 4;
1739: in.numberoftetrahedra = cells->size();
1740: in.tetrahedronlist = new int[in.numberoftetrahedra*in.numberofcorners];
1741: for(Mesh::sieve_type::heightSequence::iterator c_itor = cells->begin(); c_itor != cells->end(); ++c_itor) {
1742: Obj<Mesh::coordinate_type::IndexArray> intervals = vertexBundle->getOrderedIndices(0, serialOrientation->cone(*c_itor));
1743: int v = 0;
1745: for(Mesh::coordinate_type::IndexArray::iterator i_itor = intervals->begin(); i_itor != intervals->end(); i_itor++) {
1746: in.tetrahedronlist[c * in.numberofcorners + v++] = i_itor->prefix;
1747: }
1748: c++;
1749: }
1751: in.numberofholes = 0;
1752: ::tetrahedralize((char *) args.c_str(), &in, &out);
1753: }
1754: m->populate(out.numberoftetrahedra, out.tetrahedronlist, out.numberofpoints, out.pointlist);
1755:
1756: if (rank == 0) {
1757: Obj<Mesh::sieve_type> topology = m->getTopology();
1759: for(int v = 0; v < out.numberofpoints; v++) {
1760: if (out.pointmarkerlist[v]) {
1761: topology->setMarker(Mesh::point_type(0, v + out.numberoftetrahedra), out.pointmarkerlist[v]);
1762: }
1763: }
1764: if (out.edgemarkerlist) {
1765: for(int e = 0; e < out.numberofedges; e++) {
1766: if (out.edgemarkerlist[e]) {
1767: Mesh::point_type endpointA(0, out.edgelist[e*2+0] + out.numberoftetrahedra);
1768: Mesh::point_type endpointB(0, out.edgelist[e*2+1] + out.numberoftetrahedra);
1769: Obj<PointSet> join = topology->nJoin(endpointA, endpointB, 1);
1771: topology->setMarker(*join->begin(), out.edgemarkerlist[e]);
1772: }
1773: }
1774: }
1775: if (out.trifacemarkerlist) {
1776: for(int f = 0; f < out.numberoftrifaces; f++) {
1777: if (out.trifacemarkerlist[f]) {
1778: Obj<PointSet> point = PointSet();
1779: Obj<PointSet> edge = PointSet();
1780: Mesh::point_type cornerA(0, out.edgelist[f*3+0] + out.numberoftetrahedra);
1781: Mesh::point_type cornerB(0, out.edgelist[f*3+1] + out.numberoftetrahedra);
1782: Mesh::point_type cornerC(0, out.edgelist[f*3+2] + out.numberoftetrahedra);
1783: point->insert(cornerA);
1784: edge->insert(cornerB);
1785: edge->insert(cornerC);
1786: Obj<PointSet> join = topology->nJoin(point, edge, 2);
1788: topology->setMarker(*join->begin(), out.trifacemarkerlist[f]);
1789: }
1790: }
1791: }
1792: }
1793: return m;
1794: };
1795: #endif
1796: public:
1797: static Obj<Mesh> refine(Obj<Mesh> mesh, double maxArea) {
1798: int numElements = mesh->getTopology()->heightStratum(0)->size();
1799: double *maxAreas = new double[numElements];
1800: for(int e = 0; e < numElements; e++) {
1801: maxAreas[e] = maxArea;
1802: }
1803: Obj<Mesh> refinedMesh = refine(mesh, maxAreas);
1805: delete [] maxAreas;
1806: return refinedMesh;
1807: };
1808: static Obj<Mesh> refine(Obj<Mesh> mesh, double maxAreas[]) {
1809: Obj<Mesh> refinedMesh;
1810: int dim = mesh->getDimension();
1812: if (dim == 2) {
1813: #ifdef PETSC_HAVE_TRIANGLE
1814: refinedMesh = refine_Triangle(mesh, maxAreas);
1815: #else
1816: throw ALE::Exception("Mesh refinement currently requires Triangle to be installed. Use --download-triangle during configure.");
1817: #endif
1818: } else if (dim == 3) {
1819: #ifdef PETSC_HAVE_TETGEN
1820: refinedMesh = refine_TetGen(mesh, maxAreas);
1821: #else
1822: throw ALE::Exception("Mesh generation currently requires TetGen to be installed. Use --download-tetgen during configure.");
1823: #endif
1824: }
1825: return refinedMesh;
1826: };
1827: };
1828: } // namespace def
1829: } // namespace ALE
1831: #endif