Actual source code: CoSieve.hh
1: #ifndef included_ALE_CoSieve_hh
2: #define included_ALE_CoSieve_hh
4: #ifndef included_ALE_Sieve_hh
5: #include <Sieve.hh>
6: #endif
8: #ifndef included_ALE_Field_hh
9: #include <Field.hh>
10: #endif
14: namespace ALE {
15: // A Topology is a collection of Sieves
16: // Each Sieve has a label, which we call a \emph{patch}
17: // The collection itself we call a \emph{sheaf}
18: // The main operation we provide in Topology is the creation of a \emph{label}
19: // A label is a bidirectional mapping of Sieve points to integers, implemented with a Sifter
20: template<typename Patch_, typename Sieve_>
21: class Topology : public ALE::ParallelObject {
22: public:
23: typedef Patch_ patch_type;
24: typedef Sieve_ sieve_type;
25: typedef typename sieve_type::point_type point_type;
26: typedef typename std::map<patch_type, Obj<sieve_type> > sheaf_type;
27: typedef typename ALE::Sifter<int, point_type, int> patch_label_type;
28: typedef typename std::map<patch_type, Obj<patch_label_type> > label_type;
29: typedef typename std::map<patch_type, int> max_label_type;
30: typedef typename std::map<const std::string, label_type> labels_type;
31: typedef typename patch_label_type::supportSequence label_sequence;
32: typedef typename std::set<point_type> point_set_type;
33: typedef typename ALE::Sifter<int,point_type,point_type> send_overlap_type;
34: typedef typename ALE::Sifter<point_type,int,point_type> recv_overlap_type;
35: protected:
36: sheaf_type _sheaf;
37: labels_type _labels;
38: int _maxHeight;
39: max_label_type _maxHeights;
40: int _maxDepth;
41: max_label_type _maxDepths;
42: bool _calculatedOverlap;
43: Obj<send_overlap_type> _sendOverlap;
44: Obj<recv_overlap_type> _recvOverlap;
45: Obj<send_overlap_type> _distSendOverlap;
46: Obj<recv_overlap_type> _distRecvOverlap;
47: // Work space
48: Obj<point_set_type> _modifiedPoints;
49: public:
50: Topology(MPI_Comm comm, const int debug = 0) : ParallelObject(comm, debug), _maxHeight(-1), _maxDepth(-1), _calculatedOverlap(false) {
51: this->_sendOverlap = new send_overlap_type(this->comm(), this->debug());
52: this->_recvOverlap = new recv_overlap_type(this->comm(), this->debug());
53: this->_modifiedPoints = new point_set_type();
54: };
55: virtual ~Topology() {};
56: public: // Verifiers
57: void checkPatch(const patch_type& patch) {
58: if (this->_sheaf.find(patch) == this->_sheaf.end()) {
59: ostringstream msg;
60: msg << "Invalid topology patch: " << patch << std::endl;
61: throw ALE::Exception(msg.str().c_str());
62: }
63: };
64: void checkLabel(const std::string& name, const patch_type& patch) {
65: this->checkPatch(patch);
66: if ((this->_labels.find(name) == this->_labels.end()) || (this->_labels[name].find(patch) == this->_labels[name].end())) {
67: ostringstream msg;
68: msg << "Invalid label name: " << name << " for patch " << patch << std::endl;
69: throw ALE::Exception(msg.str().c_str());
70: }
71: };
72: bool hasPatch(const patch_type& patch) {
73: if (this->_sheaf.find(patch) != this->_sheaf.end()) {
74: return true;
75: }
76: return false;
77: };
78: bool hasLabel(const std::string& name, const patch_type& patch) {
79: if ((this->_labels.find(name) != this->_labels.end()) && (this->_labels[name].find(patch) != this->_labels[name].end())) {
80: return true;
81: }
82: return false;
83: };
84: public: // Accessors
85: const Obj<sieve_type>& getPatch(const patch_type& patch) {
86: this->checkPatch(patch);
87: return this->_sheaf[patch];
88: };
89: void setPatch(const patch_type& patch, const Obj<sieve_type>& sieve) {
90: this->_sheaf[patch] = sieve;
91: };
92: int getValue (const Obj<patch_label_type>& label, const point_type& point, const int defValue = 0) {
93: const Obj<typename patch_label_type::coneSequence>& cone = label->cone(point);
95: if (cone->size() == 0) return defValue;
96: return *cone->begin();
97: };
98: template<typename InputPoints>
99: int getMaxValue (const Obj<patch_label_type>& label, const Obj<InputPoints>& points, const int defValue = 0) {
100: int maxValue = defValue;
102: for(typename InputPoints::iterator p_iter = points->begin(); p_iter != points->end(); ++p_iter) {
103: maxValue = std::max(maxValue, this->getValue(label, *p_iter, defValue));
104: }
105: return maxValue;
106: };
107: void setValue(const Obj<patch_label_type>& label, const point_type& point, const int value) {
108: label->setCone(value, point);
109: };
110: const Obj<patch_label_type>& createLabel(const patch_type& patch, const std::string& name) {
111: this->checkPatch(patch);
112: this->_labels[name][patch] = new patch_label_type(this->comm(), this->debug());
113: return this->_labels[name][patch];
114: };
115: const Obj<patch_label_type>& getLabel(const patch_type& patch, const std::string& name) {
116: this->checkLabel(name, patch);
117: return this->_labels[name][patch];
118: };
119: const Obj<label_sequence>& getLabelStratum(const patch_type& patch, const std::string& name, int value) {
120: this->checkLabel(name, patch);
121: return this->_labels[name][patch]->support(value);
122: };
123: const sheaf_type& getPatches() {
124: return this->_sheaf;
125: };
126: const labels_type& getLabels() {
127: return this->_labels;
128: };
129: void clear() {
130: this->_sheaf.clear();
131: this->_labels.clear();
132: this->_maxHeight = -1;
133: this->_maxHeights.clear();
134: this->_maxDepth = -1;
135: this->_maxDepths.clear();
136: };
137: const Obj<send_overlap_type>& getSendOverlap() const {return this->_sendOverlap;};
138: void setSendOverlap(const Obj<send_overlap_type>& overlap) {this->_sendOverlap = overlap;};
139: const Obj<recv_overlap_type>& getRecvOverlap() const {return this->_recvOverlap;};
140: void setRecvOverlap(const Obj<recv_overlap_type>& overlap) {this->_recvOverlap = overlap;};
141: const Obj<send_overlap_type>& getDistSendOverlap() const {return this->_distSendOverlap;};
142: void setDistSendOverlap(const Obj<send_overlap_type>& overlap) {this->_distSendOverlap = overlap;};
143: const Obj<recv_overlap_type>& getDistRecvOverlap() const {return this->_distRecvOverlap;};
144: void setDistRecvOverlap(const Obj<recv_overlap_type>& overlap) {this->_distRecvOverlap = overlap;};
145: public: // Stratification
146: template<class InputPoints>
147: void computeHeight(const Obj<patch_label_type>& height, const Obj<sieve_type>& sieve, const Obj<InputPoints>& points, int& maxHeight) {
148: this->_modifiedPoints->clear();
150: for(typename InputPoints::iterator p_iter = points->begin(); p_iter != points->end(); ++p_iter) {
151: // Compute the max height of the points in the support of p, and add 1
152: int h0 = this->getValue(height, *p_iter, -1);
153: int h1 = this->getMaxValue(height, sieve->support(*p_iter), -1) + 1;
155: if(h1 != h0) {
156: this->setValue(height, *p_iter, h1);
157: if (h1 > maxHeight) maxHeight = h1;
158: this->_modifiedPoints->insert(*p_iter);
159: }
160: }
161: // FIX: We would like to avoid the copy here with cone()
162: if(this->_modifiedPoints->size() > 0) {
163: this->computeHeight(height, sieve, sieve->cone(this->_modifiedPoints), maxHeight);
164: }
165: };
166: void computeHeights() {
167: const std::string name("height");
169: this->_maxHeight = -1;
170: for(typename sheaf_type::iterator s_iter = this->_sheaf.begin(); s_iter != this->_sheaf.end(); ++s_iter) {
171: const Obj<patch_label_type>& label = this->createLabel(s_iter->first, name);
173: this->_maxHeights[s_iter->first] = -1;
174: this->computeHeight(label, s_iter->second, s_iter->second->leaves(), this->_maxHeights[s_iter->first]);
175: if (this->_maxHeights[s_iter->first] > this->_maxHeight) this->_maxHeight = this->_maxHeights[s_iter->first];
176: }
177: };
178: int height() const {return this->_maxHeight;};
179: int height(const patch_type& patch) {
180: this->checkPatch(patch);
181: return this->_maxHeights[patch];
182: };
183: int height(const patch_type& patch, const point_type& point) {
184: return this->getValue(this->_labels["height"][patch], point, -1);
185: };
186: const Obj<label_sequence>& heightStratum(const patch_type& patch, int height) {
187: return this->getLabelStratum(patch, "height", height);
188: };
189: template<class InputPoints>
190: void computeDepth(const Obj<patch_label_type>& depth, const Obj<sieve_type>& sieve, const Obj<InputPoints>& points, int& maxDepth) {
191: this->_modifiedPoints->clear();
193: for(typename InputPoints::iterator p_iter = points->begin(); p_iter != points->end(); ++p_iter) {
194: // Compute the max depth of the points in the cone of p, and add 1
195: int d0 = this->getValue(depth, *p_iter, -1);
196: int d1 = this->getMaxValue(depth, sieve->cone(*p_iter), -1) + 1;
198: if(d1 != d0) {
199: this->setValue(depth, *p_iter, d1);
200: if (d1 > maxDepth) maxDepth = d1;
201: this->_modifiedPoints->insert(*p_iter);
202: }
203: }
204: // FIX: We would like to avoid the copy here with support()
205: if(this->_modifiedPoints->size() > 0) {
206: this->computeDepth(depth, sieve, sieve->support(this->_modifiedPoints), maxDepth);
207: }
208: };
209: void computeDepths() {
210: const std::string name("depth");
212: this->_maxDepth = -1;
213: for(typename sheaf_type::iterator s_iter = this->_sheaf.begin(); s_iter != this->_sheaf.end(); ++s_iter) {
214: const Obj<patch_label_type>& label = this->createLabel(s_iter->first, name);
216: this->_maxDepths[s_iter->first] = -1;
217: this->computeDepth(label, s_iter->second, s_iter->second->roots(), this->_maxDepths[s_iter->first]);
218: if (this->_maxDepths[s_iter->first] > this->_maxDepth) this->_maxDepth = this->_maxDepths[s_iter->first];
219: }
220: };
221: int depth() const {return this->_maxDepth;};
222: int depth(const patch_type& patch) {
223: this->checkPatch(patch);
224: return this->_maxDepths[patch];
225: };
226: int depth(const patch_type& patch, const point_type& point) {
227: return this->getValue(this->_labels["depth"][patch], point, -1);
228: };
229: const Obj<label_sequence>& depthStratum(const patch_type& patch, int depth) {
230: return this->getLabelStratum(patch, "depth", depth);
231: };
234: void stratify() {
235: ALE_LOG_EVENT_BEGIN;
236: this->computeHeights();
237: this->computeDepths();
238: ALE_LOG_EVENT_END;
239: };
240: public: // Viewers
241: void view(const std::string& name, MPI_Comm comm = MPI_COMM_NULL) {
242: if (comm == MPI_COMM_NULL) {
243: comm = this->comm();
244: }
245: if (name == "") {
246: PetscPrintf(comm, "viewing a Topology\n");
247: } else {
248: PetscPrintf(comm, "viewing Topology '%s'\n", name.c_str());
249: }
250: PetscPrintf(comm, " maximum height %d maximum depth %d\n", this->height(), this->depth());
251: for(typename sheaf_type::const_iterator s_iter = this->_sheaf.begin(); s_iter != this->_sheaf.end(); ++s_iter) {
252: ostringstream txt;
254: txt << "Patch " << s_iter->first;
255: s_iter->second->view(txt.str().c_str(), comm);
256: PetscPrintf(comm, " maximum height %d maximum depth %d\n", this->height(s_iter->first), this->depth(s_iter->first));
257: }
258: for(typename labels_type::const_iterator l_iter = this->_labels.begin(); l_iter != this->_labels.end(); ++l_iter) {
259: PetscPrintf(comm, " label %s constructed\n", l_iter->first.c_str());
260: }
261: };
262: public:
263: void constructOverlap(const patch_type& patch) {
264: if (this->_calculatedOverlap) return;
265: if (this->hasPatch(patch)) {
266: this->constructOverlap(this->getPatch(patch)->base(), this->_sendOverlap, this->_recvOverlap);
267: this->constructOverlap(this->getPatch(patch)->cap(), this->_sendOverlap, this->_recvOverlap);
268: }
269: if (this->debug()) {
270: this->_sendOverlap->view("Send overlap");
271: this->_recvOverlap->view("Receive overlap");
272: }
273: this->_calculatedOverlap = true;
274: };
275: template<typename Sequence>
276: void constructOverlap(const Obj<Sequence>& points, const Obj<send_overlap_type>& sendOverlap, const Obj<recv_overlap_type>& recvOverlap) {
277: point_type *sendBuf = new point_type[points->size()];
278: int size = 0;
279: for(typename Sequence::iterator l_iter = points->begin(); l_iter != points->end(); ++l_iter) {
280: sendBuf[size++] = *l_iter;
281: }
282: int *sizes = new int[this->commSize()]; // The number of points coming from each process
283: int *offsets = new int[this->commSize()+1]; // Prefix sums for sizes
284: int *oldOffs = new int[this->commSize()+1]; // Temporary storage
285: point_type *remotePoints = NULL; // The points from each process
286: int *remoteRanks = NULL; // The rank and number of overlap points of each process that overlaps another
288: // Change to Allgather() for the correct binning algorithm
289: MPI_Gather(&size, 1, MPI_INT, sizes, 1, MPI_INT, 0, this->comm());
290: if (this->commRank() == 0) {
291: offsets[0] = 0;
292: for(int p = 1; p <= this->commSize(); p++) {
293: offsets[p] = offsets[p-1] + sizes[p-1];
294: }
295: remotePoints = new point_type[offsets[this->commSize()]];
296: }
297: MPI_Gatherv(sendBuf, size, MPI_INT, remotePoints, sizes, offsets, MPI_INT, 0, this->comm());
298: std::map<int, std::map<int, std::set<point_type> > > overlapInfo; // Maps (p,q) to their set of overlap points
300: if (this->commRank() == 0) {
301: for(int p = 0; p < this->commSize(); p++) {
302: std::sort(&remotePoints[offsets[p]], &remotePoints[offsets[p+1]]);
303: }
304: for(int p = 0; p <= this->commSize(); p++) {
305: oldOffs[p] = offsets[p];
306: }
307: for(int p = 0; p < this->commSize(); p++) {
308: for(int q = p+1; q < this->commSize(); q++) {
309: std::set_intersection(&remotePoints[oldOffs[p]], &remotePoints[oldOffs[p+1]],
310: &remotePoints[oldOffs[q]], &remotePoints[oldOffs[q+1]],
311: std::insert_iterator<std::set<point_type> >(overlapInfo[p][q], overlapInfo[p][q].begin()));
312: overlapInfo[q][p] = overlapInfo[p][q];
313: }
314: sizes[p] = overlapInfo[p].size()*2;
315: offsets[p+1] = offsets[p] + sizes[p];
316: }
317: remoteRanks = new int[offsets[this->commSize()]];
318: int k = 0;
319: for(int p = 0; p < this->commSize(); p++) {
320: for(typename std::map<int, std::set<point_type> >::iterator r_iter = overlapInfo[p].begin(); r_iter != overlapInfo[p].end(); ++r_iter) {
321: remoteRanks[k*2] = r_iter->first;
322: remoteRanks[k*2+1] = r_iter->second.size();
323: k++;
324: }
325: }
326: }
327: int numOverlaps; // The number of processes overlapping this process
328: MPI_Scatter(sizes, 1, MPI_INT, &numOverlaps, 1, MPI_INT, 0, this->comm());
329: int *overlapRanks = new int[numOverlaps]; // The rank and overlap size for each overlapping process
330: MPI_Scatterv(remoteRanks, sizes, offsets, MPI_INT, overlapRanks, numOverlaps, MPI_INT, 0, this->comm());
331: point_type *sendPoints = NULL; // The points to send to each process
332: if (this->commRank() == 0) {
333: for(int p = 0, k = 0; p < this->commSize(); p++) {
334: sizes[p] = 0;
335: for(int r = 0; r < (int) overlapInfo[p].size(); r++) {
336: sizes[p] += remoteRanks[k*2+1];
337: k++;
338: }
339: offsets[p+1] = offsets[p] + sizes[p];
340: }
341: sendPoints = new point_type[offsets[this->commSize()]];
342: for(int p = 0, k = 0; p < this->commSize(); p++) {
343: for(typename std::map<int, std::set<point_type> >::iterator r_iter = overlapInfo[p].begin(); r_iter != overlapInfo[p].end(); ++r_iter) {
344: int rank = r_iter->first;
345: for(typename std::set<point_type>::iterator p_iter = (overlapInfo[p][rank]).begin(); p_iter != (overlapInfo[p][rank]).end(); ++p_iter) {
346: sendPoints[k++] = *p_iter;
347: }
348: }
349: }
350: }
351: int numOverlapPoints = 0;
352: for(int r = 0; r < numOverlaps/2; r++) {
353: numOverlapPoints += overlapRanks[r*2+1];
354: }
355: point_type *overlapPoints = new point_type[numOverlapPoints];
356: MPI_Scatterv(sendPoints, sizes, offsets, MPI_INT, overlapPoints, numOverlapPoints, MPI_INT, 0, this->comm());
358: for(int r = 0, k = 0; r < numOverlaps/2; r++) {
359: int rank = overlapRanks[r*2];
361: for(int p = 0; p < overlapRanks[r*2+1]; p++) {
362: point_type point = overlapPoints[k++];
364: sendOverlap->addArrow(point, rank, point);
365: recvOverlap->addArrow(rank, point, point);
366: }
367: }
369: delete [] overlapPoints;
370: delete [] overlapRanks;
371: delete [] sizes;
372: delete [] offsets;
373: delete [] oldOffs;
374: if (this->commRank() == 0) {
375: delete [] remoteRanks;
376: delete [] remotePoints;
377: delete [] sendPoints;
378: }
379: };
380: };
382: template<typename Bundle_>
383: class SieveBuilder {
384: public:
385: typedef Bundle_ bundle_type;
386: typedef typename bundle_type::sieve_type sieve_type;
387: typedef typename bundle_type::arrow_section_type arrow_section_type;
388: typedef std::vector<typename sieve_type::point_type> PointArray;
389: public:
390: static void buildHexFaces(Obj<sieve_type> sieve, int dim, std::map<int, int*>& curElement, std::map<int,PointArray>& bdVertices, std::map<int,PointArray>& faces, typename sieve_type::point_type& cell) {
391: int debug = sieve->debug();
393: if (debug > 1) {std::cout << " Building hex faces for boundary of " << cell << " (size " << bdVertices[dim].size() << "), dim " << dim << std::endl;}
394: faces[dim].clear();
395: if (dim > 3) {
396: throw ALE::Exception("Cannot do hexes of dimension greater than three");
397: } else if (dim > 2) {
398: int nodes[24] = {0, 1, 2, 3, 4, 5, 6, 7, 0, 1, 5, 4,
399: 1, 2, 6, 5, 2, 3, 7, 6, 3, 0, 4, 7};
401: for(int b = 0; b < 6; b++) {
402: typename sieve_type::point_type face;
404: bdVertices[dim-1].clear();
405: for(int c = 0; c < 4; c++) {
406: bdVertices[dim-1].push_back(bdVertices[dim][nodes[b*4+c]]);
407: }
408: if (debug > 1) {std::cout << " boundary hex face " << b << std::endl;}
409: buildHexFaces(sieve, dim-1, curElement, bdVertices, faces, face);
410: if (debug > 1) {std::cout << " added face " << face << std::endl;}
411: faces[dim].push_back(face);
412: }
413: } else if (dim > 1) {
414: int boundarySize = bdVertices[dim].size();
416: for(int b = 0; b < boundarySize; b++) {
417: typename sieve_type::point_type face;
419: bdVertices[dim-1].clear();
420: for(int c = 0; c < 2; c++) {
421: bdVertices[dim-1].push_back(bdVertices[dim][(b+c)%boundarySize]);
422: }
423: if (debug > 1) {
424: std::cout << " boundary point " << bdVertices[dim][b] << std::endl;
425: std::cout << " boundary vertices";
426: for(int c = 0; c < (int) bdVertices[dim-1].size(); c++) {
427: std::cout << " " << bdVertices[dim-1][c];
428: }
429: std::cout << std::endl;
430: }
431: buildHexFaces(sieve, dim-1, curElement, bdVertices, faces, face);
432: if (debug > 1) {std::cout << " added face " << face << std::endl;}
433: faces[dim].push_back(face);
434: }
435: } else {
436: if (debug > 1) {std::cout << " Just set faces to boundary in 1d" << std::endl;}
437: faces[dim].insert(faces[dim].end(), bdVertices[dim].begin(), bdVertices[dim].end());
438: }
439: if (debug > 1) {
440: for(typename PointArray::iterator f_iter = faces[dim].begin(); f_iter != faces[dim].end(); ++f_iter) {
441: std::cout << " face point " << *f_iter << std::endl;
442: }
443: }
444: // We always create the toplevel, so we could short circuit somehow
445: // Should not have to loop here since the meet of just 2 boundary elements is an element
446: typename PointArray::iterator f_itor = faces[dim].begin();
447: const typename sieve_type::point_type& start = *f_itor;
448: const typename sieve_type::point_type& next = *(++f_itor);
449: Obj<typename sieve_type::supportSet> preElement = sieve->nJoin(start, next, 1);
451: if (preElement->size() > 0) {
452: cell = *preElement->begin();
453: if (debug > 1) {std::cout << " Found old cell " << cell << std::endl;}
454: } else {
455: int color = 0;
457: cell = typename sieve_type::point_type((*curElement[dim])++);
458: for(typename PointArray::iterator f_itor = faces[dim].begin(); f_itor != faces[dim].end(); ++f_itor) {
459: sieve->addArrow(*f_itor, cell, color++);
460: }
461: if (debug > 1) {std::cout << " Added cell " << cell << " dim " << dim << std::endl;}
462: }
463: };
464: static void buildFaces(Obj<sieve_type> sieve, int dim, std::map<int, int*>& curElement, std::map<int,PointArray>& bdVertices, std::map<int,PointArray>& faces, typename sieve_type::point_type& cell) {
465: int debug = sieve->debug();
467: if (debug > 1) {
468: if (cell >= 0) {
469: std::cout << " Building faces for boundary of " << cell << " (size " << bdVertices[dim].size() << "), dim " << dim << std::endl;
470: } else {
471: std::cout << " Building faces for boundary of undetermined cell (size " << bdVertices[dim].size() << "), dim " << dim << std::endl;
472: }
473: }
474: faces[dim].clear();
475: if (dim > 1) {
476: // Use the cone construction
477: for(typename PointArray::iterator b_itor = bdVertices[dim].begin(); b_itor != bdVertices[dim].end(); ++b_itor) {
478: typename sieve_type::point_type face = -1;
480: bdVertices[dim-1].clear();
481: for(typename PointArray::iterator i_itor = bdVertices[dim].begin(); i_itor != bdVertices[dim].end(); ++i_itor) {
482: if (i_itor != b_itor) {
483: bdVertices[dim-1].push_back(*i_itor);
484: }
485: }
486: if (debug > 1) {std::cout << " boundary point " << *b_itor << std::endl;}
487: buildFaces(sieve, dim-1, curElement, bdVertices, faces, face);
488: if (debug > 1) {std::cout << " added face " << face << std::endl;}
489: faces[dim].push_back(face);
490: }
491: } else {
492: if (debug > 1) {std::cout << " Just set faces to boundary in 1d" << std::endl;}
493: faces[dim].insert(faces[dim].end(), bdVertices[dim].begin(), bdVertices[dim].end());
494: }
495: if (debug > 1) {
496: for(typename PointArray::iterator f_iter = faces[dim].begin(); f_iter != faces[dim].end(); ++f_iter) {
497: std::cout << " face point " << *f_iter << std::endl;
498: }
499: }
500: // We always create the toplevel, so we could short circuit somehow
501: // Should not have to loop here since the meet of just 2 boundary elements is an element
502: typename PointArray::iterator f_itor = faces[dim].begin();
503: const typename sieve_type::point_type& start = *f_itor;
504: const typename sieve_type::point_type& next = *(++f_itor);
505: Obj<typename sieve_type::supportSet> preElement = sieve->nJoin(start, next, 1);
507: if (preElement->size() > 0) {
508: cell = *preElement->begin();
509: if (debug > 1) {std::cout << " Found old cell " << cell << std::endl;}
510: } else {
511: int color = 0;
513: cell = typename sieve_type::point_type((*curElement[dim])++);
514: for(typename PointArray::iterator f_itor = faces[dim].begin(); f_itor != faces[dim].end(); ++f_itor) {
515: sieve->addArrow(*f_itor, cell, color++);
516: }
517: if (debug > 1) {std::cout << " Added cell " << cell << " dim " << dim << std::endl;}
518: }
519: };
523: // Build a topology from a connectivity description
524: // (0, 0) ... (0, numCells-1): dim-dimensional cells
525: // (0, numCells) ... (0, numVertices): vertices
526: // The other cells are numbered as they are requested
527: static void buildTopology(Obj<sieve_type> sieve, int dim, int numCells, int cells[], int numVertices, bool interpolate = true, int corners = -1, int firstVertex = -1, Obj<arrow_section_type> orientation = NULL) {
528: int debug = sieve->debug();
530: ALE_LOG_EVENT_BEGIN;
531: if (sieve->commRank() != 0) {
532: ALE_LOG_EVENT_END;
533: return;
534: }
535: if (firstVertex < 0) firstVertex = numCells;
536: // Create a map from dimension to the current element number for that dimension
537: std::map<int,int*> curElement;
538: std::map<int,PointArray> bdVertices;
539: std::map<int,PointArray> faces;
540: int curCell = 0;
541: int curVertex = firstVertex;
542: int newElement = firstVertex+numVertices;
544: if (corners < 0) corners = dim+1;
545: curElement[0] = &curVertex;
546: curElement[dim] = &curCell;
547: for(int d = 1; d < dim; d++) {
548: curElement[d] = &newElement;
549: }
550: for(int c = 0; c < numCells; c++) {
551: typename sieve_type::point_type cell(c);
553: // Build the cell
554: if (interpolate) {
555: bdVertices[dim].clear();
556: for(int b = 0; b < corners; b++) {
557: // This ordering produces the same vertex order as the uninterpolated mesh
558: //typename sieve_type::point_type vertex(cells[c*corners+(b+corners-1)%corners]+firstVertex);
559: typename sieve_type::point_type vertex(cells[c*corners+b]+firstVertex);
561: if (debug > 1) {std::cout << "Adding boundary vertex " << vertex << std::endl;}
562: bdVertices[dim].push_back(vertex);
563: }
564: if (debug) {std::cout << "cell " << cell << " num boundary vertices " << bdVertices[dim].size() << std::endl;}
566: if (corners != dim+1) {
567: buildHexFaces(sieve, dim, curElement, bdVertices, faces, cell);
568: } else {
569: buildFaces(sieve, dim, curElement, bdVertices, faces, cell);
570: }
571: if ((dim == 2) && (!orientation.isNull())) {
572: if (debug > 1) {std::cout << "Orienting cell " << cell << std::endl;}
573: const Obj<typename sieve_type::traits::coneSequence>& cone = sieve->cone(cell);
574: const typename sieve_type::traits::coneSequence::iterator end = cone->end();
576: for(typename sieve_type::traits::coneSequence::iterator p_iter = cone->begin(); p_iter != end; ++p_iter) {
577: if (debug > 1) {std::cout << " edge " << *p_iter << std::endl;}
578: const Obj<typename sieve_type::traits::coneSequence>& vertices = sieve->cone(*p_iter);
579: typename sieve_type::traits::coneSequence::iterator vertex = vertices->begin();
580: MinimalArrow<typename sieve_type::point_type,typename sieve_type::point_type> arrow(*p_iter, cell);
581: int indA, indB, value;
583: orientation->addPoint(arrow);
584: for(indA = 0; indA < corners; indA++) {if (*vertex == cells[c*corners+indA] + numCells) break;}
585: if (debug > 1) {std::cout << " vertexA " << *vertex << " indA " << indA <<std::endl;}
586: ++vertex;
587: for(indB = 0; indB < corners; indB++) {if (*vertex == cells[c*corners+indB] + numCells) break;}
588: if (debug > 1) {std::cout << " vertexB " << *vertex << " indB " << indB <<std::endl;}
589: if ((indA == corners) || (indB == corners) || (indA == indB)) {throw ALE::Exception("Invalid edge endpoints");}
590: if ((indB - indA == 1) || (indA - indB == 2)) {
591: value = 1;
592: } else {
593: value = -1;
594: }
595: if (debug > 1) {std::cout << " value " << value <<std::endl;}
596: orientation->updatePoint(arrow, &value);
597: }
598: }
599: } else {
600: for(int b = 0; b < corners; b++) {
601: sieve->addArrow(typename sieve_type::point_type(cells[c*corners+b]+firstVertex), cell, b);
602: }
603: if (debug) {
604: if (debug > 1) {
605: for(int b = 0; b < corners; b++) {
606: std::cout << " Adding vertex " << typename sieve_type::point_type(cells[c*corners+b]+firstVertex) << std::endl;
607: }
608: }
609: std::cout << "Adding cell " << cell << " dim " << dim << std::endl;
610: }
611: }
612: }
613: ALE_LOG_EVENT_END;
614: };
615: static void buildCoordinates(const Obj<Bundle_>& bundle, const int embedDim, const double coords[]) {
616: const Obj<typename Bundle_::real_section_type>& coordinates = bundle->getRealSection("coordinates");
617: const Obj<typename Bundle_::label_sequence>& vertices = bundle->depthStratum(0);
618: const int numCells = bundle->heightStratum(0)->size();
620: coordinates->setFiberDimension(vertices, embedDim);
621: bundle->allocate(coordinates);
622: for(typename Bundle_::label_sequence::iterator v_iter = vertices->begin(); v_iter != vertices->end(); ++v_iter) {
623: coordinates->updatePoint(*v_iter, &(coords[(*v_iter - numCells)*embedDim]));
624: }
625: };
626: };
628: // An Overlap is a Sifter describing the overlap of two Sieves
629: // Each arrow is local point ---(remote point)---> remote rank right now
630: // For XSifter, this should change to (local patch, local point) ---> (remote rank, remote patch, remote point)
631: }
633: #endif