Actual source code: Sieve.cxx
1: #define ALE_Sieve_cxx
3: #ifndef included_ALE_Sieve_hh
4: #include <Sieve.hh>
5: #endif
7: #include <stack>
9: namespace ALE {
10:
11: #undef __FUNCT__
13: Sieve::Sieve() : PreSieve(),_additionPolicy(additionPolicyAcyclic), _stratificationPolicy(stratificationPolicyOnMutation) {};
15: #undef __FUNCT__
17: Sieve::Sieve(MPI_Comm comm) : PreSieve(comm), _additionPolicy(additionPolicyAcyclic), _stratificationPolicy(stratificationPolicyOnMutation) {};
19: #undef __FUNCT__
21: Sieve::~Sieve(){};
23: #undef __FUNCT__
25: void Sieve::setComm(MPI_Comm c) {
26: Coaster::setComm(c);
27: }// Coaster::setComm()
30: #undef __FUNCT__
32: Sieve& Sieve::clear(){
33: PreSieve::clear();
34: this->_depth.clear();
35: this->_height.clear();
36: return *this;
37: }// Sieve::clear()
38:
40: #undef __FUNCT__
42: Sieve& Sieve::getLock(){
43: CHKCOMM(*this);
44: this->_lock++;
45: PetscErrorCode MPI_Barrier(this->getComm()); CHKERROR(ierr, "Error in MPI_Barrier");
46: if(this->_stratificationPolicy == stratificationPolicyOnLocking) {
47: this->__computeClosureHeights(this->cone(this->_leaves));
48: this->__computeStarDepths(this->support(this->_roots));
49: }
50: return *this;
51: };
52: #undef __FUNCT__
54: Sieve& Sieve::setAdditionPolicy(Sieve::AdditionPolicy policy) {
55: this->__checkLock();
56: if( (policy != additionPolicyAcyclic) && (policy != additionPolicyStratified) ) {
57: throw ALE::Exception("Unknown AdditionPolicy value");
58: }
59: this->_additionPolicy = policy;
60: return *this;
61: }// Sieve::setAdditionPolicy()
63: #undef __FUNCT__
65: Sieve::AdditionPolicy Sieve::getAdditionPolicy() {
66: return this->_additionPolicy;
67: }// Sieve::getAdditionPolicy()
69: #undef __FUNCT__
71: Sieve& Sieve::setStratificationPolicy(Sieve::StratificationPolicy policy) {
72: this->__checkLock();
73: if( (policy != stratificationPolicyOnLocking) && (policy != stratificationPolicyOnMutation) ) {
74: throw ALE::Exception("Unknown StratificationPolicy value");
75: }
76: this->_stratificationPolicy = policy;
77: return *this;
78: }// Sieve::setStratificationPolicy()
80: #undef __FUNCT__
82: Sieve::StratificationPolicy Sieve::getStratificationPolicy() {
83: return this->_stratificationPolicy;
84: }// Sieve::getStratificationPolicy()
86: #undef __FUNCT__
88: Sieve& Sieve::addArrow(const Point& i, const Point& j) {
89: ALE_LOG_STAGE_BEGIN;
90: CHKCOMM(*this);
91: this->__checkLock();
92: // Check whether the arrow addition would violate the addition policy.
93: // This can be done only if the heights and depths are up-to-date, which means 'stratificationPolycyOnMutation'.
94: if(this->_stratificationPolicy == stratificationPolicyOnMutation) {
95: Point_set iSet, jSet, iClosure, jStar;
96: int32_t iDepth, jDepth;
97: ostringstream txt;
98: switch(this->_additionPolicy) {
99: case additionPolicyAcyclic:
100: if (i == j) {
101: ostringstream ex;
102: ex << "[" << this->getCommRank() << "]: ";
103: ex << "Attempted arrow insertion (" << i.prefix << ", " << i.index << ") --> (" << j.prefix << ", " << j.index << ") ";
104: ex << "would lead to a cycle";
105: //throw Exception("Attempted arrow insertion would lead to a cycle");
106: throw Exception(ex.str().c_str());
107: };
108: //jSet.insert(j);
109: //iSet.insert(i);
110:
111: jStar = this->star(j);
112: iClosure = this->closure(i);
113: if((jStar.find(i) != jStar.end()) || (iClosure.find(j) != iClosure.end())) {
114: printf("Adding (%d, %d)-->(%d, %d)\n", i.prefix, i.index, j.prefix, j.index);
115: if (jStar.find(i) != jStar.end()) {
116: printf("Head found in tail star\n");
117: }
118: if (iClosure.find(j) != iClosure.end()) {
119: printf("Tail found in head closure\n");
120: }
121: ostringstream ex;
122: ex << "[" << this->getCommRank() << "]: ";
123: ex << "Attempted arrow insertion (" << i.prefix << ", " << i.index << ") --> (" << j.prefix << ", " << j.index << ") ";
124: ex << "would lead to a cycle";
125: //throw Exception("Attempted arrow insertion would lead to a cycle");
126: throw Exception(ex.str().c_str());
127: }
128: break;
129: case additionPolicyStratified:
130: iDepth = this->depth(i);
131: jDepth = this->depth(j);
132: // Recall that iDepth < 0 means i is not in the Sieve; likewise for jDepth
133: if( (iDepth >= 0) && (jDepth >= 0) && (iDepth >= jDepth) ){
134: ostringstream ex;
135: ex << "[" << this->getCommRank() << "]: ";
136: ex << "Attempted arrow insertion would violate stratification";
137: //throw Exception("Attempted arrow insertion would violate stratification");
138: throw Exception(ex.str().c_str());
139: }
140: break;
141: default:
142: ostringstream ex;
143: ex << "[" << this->getCommRank() << "]: ";
144: ex << "Unknown addition policy";
145: throw Exception("Unknown addition policy");
146: //throw Exception(ex.str().c_str());
147: }// switch(this->_additionPolicy)
148: }// if(this->_stratificationPolicy == stratificationPolicyOnMutation)
150: // Now add the points
151: this->addCapPoint(i);
152: this->addBasePoint(j);
153:
154: // Finally, add the arrow
155: PreSieve::addArrow(i,j); // need to use PreSieve::addArrow to make sure that roots and leaves are consistent
157: // The heights and depths are computed now or later, depending on the stratificationPolicy
158: if(this->_stratificationPolicy == stratificationPolicyOnMutation) {
159: // Only points in the star of j may have gotten their depth changed by introducing a path up through i
160: if(this->depth(j) != (this->depth(i)+1)) {
161: this->__computeStarDepths(j);
162: }
163: // Only points in the closure of i may have gotten their height changed by introducing a path down through j
164: if(this->height(i) != (this->height(j) + 1)) {
165: this->__computeClosureHeights(i);
166: }
167: }
168: ALE_LOG_STAGE_END;
169: return *this;
170: }// Sieve::addArrow()
173: #undef __FUNCT__
175: Sieve& Sieve::removeArrow(const Point& i, const Point& j, bool removeSingleton) {
176: ALE_LOG_STAGE_BEGIN;
177: CHKCOMM(*this);
178: this->__checkLock();
179: // need to use PreSieve::removeArrow to make sure that roots and leaves are consistent
180: PreSieve::removeArrow(i,j, removeSingleton);
181: // now we may need to recompute the heights and depths
182: if(this->_stratificationPolicy == stratificationPolicyOnMutation) {
183: if (this->spaceContains(i)) {
184: // Only the points in the closure of i may have gotten their height changed, now that the path down through j is unavailable
185: __computeClosureHeights(i);
186: }
187: if (this->spaceContains(j)) {
188: // Only the points in the star of j may have gotten their depth changed, now that the path up through i is unavailable
189: __computeStarDepths(j);
190: }
191: }
192: ALE_LOG_STAGE_END;
193: return *this;
194: }// Sieve::removeArrow()
196: #undef __FUNCT__
198: Sieve& Sieve::removeBasePoint(const Point& p, bool removeSingleton) {
199: ALE_LOG_STAGE_BEGIN;
200: this->__checkLock();
201: ALE::PreSieve::removeBasePoint(p, removeSingleton);
202: if (!this->capContains(p)) {
203: this->__setDepth(p, -1);
204: this->__setHeight(p, -1);
205: }
206: ALE_LOG_STAGE_END;
207: return *this;
208: }
209:
210: #undef __FUNCT__
212: Sieve& Sieve::addBasePoint(const Point& p) {
213: ALE_LOG_STAGE_BEGIN;
214: CHKCOMM(*this);
215: this->__checkLock();
216: // If the point is absent from the Sieve, after insertion it will have zero height/depth
217: if(!(this->spaceContains(p))) {
218: this->__setDepth(p,0);
219: this->__setHeight(p,0);
220: }
221: // We must use PreSieve methods to make sure roots and leaves are maintained consistently
222: PreSieve::addBasePoint(p);
223: ALE_LOG_STAGE_END;
224: return *this;
225: }// Sieve::addBasePoint()
227: #undef __FUNCT__
229: Sieve& Sieve::removeCapPoint(const Point& q, bool removeSingleton) {
230: ALE_LOG_STAGE_BEGIN;
231: this->__checkLock();
232: ALE::PreSieve::removeCapPoint(q, removeSingleton);
233: if (!this->baseContains(q)) {
234: this->__setDepth(q, -1);
235: this->__setHeight(q, -1);
236: }
237: ALE_LOG_STAGE_END;
238: return *this;
239: }
241: #undef __FUNCT__
243: Sieve& Sieve::addCapPoint(const Point& p) {
244: ALE_LOG_STAGE_BEGIN;
245: CHKCOMM(*this);
246: this->__checkLock();
247: // If the point is absent from the Sieve, after insertion it will have zero height/depth
248: if(!(this->spaceContains(p))) {
249: this->__setDepth(p,0);
250: this->__setHeight(p,0);
251: }
252: // We must use PreSieve methods to make sure roots and leaves are maintained consistently
253: PreSieve::addCapPoint(p);
254: ALE_LOG_STAGE_END;
255: return *this;
256: }// Sieve::addCapPoint()
259: #undef __FUNCT__
261: void Sieve::__setHeight(Point p, int32_t h){
262: ALE_LOG_STAGE_BEGIN;
263: this->_height[p] = h;
264: ALE_LOG_STAGE_END;
265: }// Sieve::__setHeight()
267: #undef __FUNCT__
269: void Sieve::__setDepth(Point p, int32_t d){
270: ALE_LOG_STAGE_BEGIN;
271: this->_depth[p] = d;
272: ALE_LOG_STAGE_END;
273: }// Sieve::__setDepth()
275: #undef __FUNCT__
277: void Sieve::__computeClosureHeights(Obj<Point_set> points) {
278: ALE_LOG_STAGE_BEGIN;
279: ALE_LOG_EVENT_BEGIN
280: // points contains points for the current height computation;
281: // mpoints keeps track of 'modified' points identified at the current stage,
282: // and through which recursion propagates
283: Obj<Point_set> mpoints = Point_set();
285: for(Point_set::iterator p_itor = points->begin(); p_itor != points->end(); p_itor++) {
286: // retrieve the current height of p
287: int32_t h0 = this->height(*p_itor);
288: // compute the max height of the points in the support of p
289: int32_t maxH = this->maxHeight(this->support(*p_itor));
290: // the height is maxH + 1
291: int32_t h1 = maxH + 1;
292: // if h0 differs from h1, set the height of p to h1 and add p to mpoints -- its height has been modified
293: if(h1 != h0) {
294: this->__setHeight(*p_itor,h1);
295: mpoints->insert(*p_itor);
296: }
297: }
298: // if the mpoints set is not empty, we recursively call __computeClosureHeights on points = cone(mpoints)
299: if(mpoints->size() > 0) {
300: this->__computeClosureHeights(this->cone(mpoints));
301: }
302: ALE_LOG_EVENT_END
303: ALE_LOG_STAGE_END;
304: }//Sieve::__computeClosureHeights()
306: #undef __FUNCT__
308: void Sieve::__computeStarDepths(Obj<Point_set> points) {
309: ALE_LOG_STAGE_BEGIN;
310: ALE_LOG_EVENT_BEGIN
311: // points contains points for the current depth computation;
312: // mpoints keeps track of 'modified' points identified at the current stage,
313: // and through which recursion propagates
314: Obj<Point_set> mpoints = Point_set();
316: for(Point_set::iterator p_itor = points->begin(); p_itor != points->end(); p_itor++) {
317: // retrieve the current depth of p
318: int32_t d0 = this->depth(*p_itor);
319: // compute the max depth of the points in the cone over p
320: int32_t maxD = this->maxDepth(this->cone(*p_itor));
321: // the new depth is maxD + 1
322: int32_t d1 = maxD + 1;
323: // if d0 differs from d1, set the depth of p to d1 and add p to mpoints -- its depth has been modified
324: if(d1 != d0) {
325: this->__setDepth(*p_itor,d1);
326: mpoints->insert(*p_itor);
327: }
328: }
329: // if the mpoints set is not empty, we recursively call __computeStarDepths on points = support(mpoints)
330: if(mpoints->size() > 0) {
331: this->__computeStarDepths(this->support(mpoints));
332: }
333: ALE_LOG_EVENT_END
334: ALE_LOG_STAGE_END;
335: }//Sieve::__computeStarDepths()
337: #undef __FUNCT__
339: int32_t Sieve::depth(const Point& p) {
340: int32_t depth;
342: ALE_LOG_STAGE_BEGIN;
343: ALE_LOG_EVENT_BEGIN
344: CHKCOMM(*this);
345: if(this->_depth.find(p) != this->_depth.end()) {
346: depth = this->_depth[p];
347: } else {
348: /* This accomdates Stacks, since spaceContains() can return true before the point is added to the Stack itself */
349: depth = -1;
350: }
351: ALE_LOG_EVENT_END
352: ALE_LOG_STAGE_END;
353: return depth;
354: }// Sieve::depth()
356: #undef __FUNCT__
358: int32_t Sieve::maxDepth(Point_set& points) {
359: int32_t max = -1;
360: for(Point_set::iterator p_itor = points.begin(); p_itor != points.end(); p_itor++) {
361: int32_t m = this->depth(*p_itor);
362: if(m > max) {
363: max = m;
364: }
365: }
366: return max;
367: }// Sieve:maxDepth()
369: #undef __FUNCT__
371: int32_t Sieve::maxDepth(Obj<Point_set> points) {
372: int32_t max = -1;
373: for(Point_set::iterator p_itor = points->begin(); p_itor != points->end(); p_itor++) {
374: int32_t m = this->depth(*p_itor);
375: if(m > max) {
376: max = m;
377: }
378: }
379: return max;
380: }// Sieve:maxDepth()
382: #undef __FUNCT__
384: int32_t Sieve::maxHeight(Point_set& points) {
385: int32_t max = -1;
386: for(Point_set::iterator p_itor = points.begin(); p_itor != points.end(); p_itor++) {
387: int32_t m = this->height(*p_itor);
388: if(m > max) {
389: max = m;
390: }
391: }
392: return max;
393: }// Sieve:maxHeight()
395: #undef __FUNCT__
397: int32_t Sieve::maxHeight(Obj<Point_set> points) {
398: int32_t max = -1;
399: for(Point_set::iterator p_itor = points->begin(); p_itor != points->end(); p_itor++) {
400: int32_t m = this->height(*p_itor);
401: if(m > max) {
402: max = m;
403: }
404: }
405: return max;
406: }// Sieve:maxHeight()
408: #undef __FUNCT__
410: int32_t Sieve::height(Point p) {
411: int32_t height;
412: ALE_LOG_STAGE_BEGIN;
413: CHKCOMM(*this);
414: if(this->_height.find(p) != this->_height.end()) {
415: height = this->_height[p];
416: } else {
417: /* This accomdates Stacks, since spaceContains() can return true before the point is added to the Stack itself */
418: height = -1;
419: }
420: ALE_LOG_STAGE_END;
421: return height;
422: }// Sieve::height()
424: #undef __FUNCT__
426: int32_t Sieve::diameter(Point p) {
427: int32_t diameter;
428: ALE_LOG_STAGE_BEGIN;
429: CHKCOMM(*this);
430: if(this->spaceContains(p)) {
431: diameter = this->depth(p) + this->height(p);
432: }
433: else {
434: diameter = -1;
435: }
436: ALE_LOG_STAGE_END;
437: return diameter;
438: }// Sieve::diameter(Point)
440: #undef __FUNCT__
442: int32_t Sieve::diameter() {
443: int globalDiameter;
444: ALE_LOG_STAGE_BEGIN;
445: CHKCOMM(*this);
446: int32_t diameter = 0;
447: // IMPROVE: PreSieve::space() should return an iterator instead
448: Point_set space = this->space();
449: for(Point_set::iterator s_itor = space.begin(); s_itor != space.end(); s_itor++) {
450: Point p = *s_itor;
451: int32_t pDiameter = this->diameter(p);
452: if(pDiameter > diameter) {
453: diameter = pDiameter;
454: }
455: }
456: int MPI_Allreduce(&diameter, &globalDiameter, 1, MPI_INT, MPI_MAX, this->comm);
457: CHKMPIERROR(ierr, ERRORMSG("Error in MPI_Allreduce"));
458: ALE_LOG_STAGE_END;
459: return globalDiameter;
460: }// Sieve::diameter()
463: #undef __FUNCT__
465: Obj<Point_set> Sieve::closure(Obj<Point_set> chain) {
466: Obj<Point_set> closure = Point_set();
468: ALE_LOG_STAGE_BEGIN;
469: ALE_LOG_EVENT_BEGIN
470: CHKCOMM(*this);
471: int32_t depth = this->maxDepth(chain);
472: if(depth >= 0) {
473: closure = this->nClosure(chain, depth);
474: }
475: ALE_LOG_EVENT_END
476: ALE_LOG_STAGE_END;
477: return closure;
478: }// Sieve::closure()
480: #undef __FUNCT__
482: Obj<Sieve> Sieve::closureSieve(Obj<Point_set> chain, Obj<Sieve> closure) {
483: ALE_LOG_STAGE_BEGIN;
484: CHKCOMM(*this);
485: if(closure.isNull()) {
486: closure = Obj<Sieve>(new Sieve(this->getComm()));
487: }
488: int32_t depth = this->maxDepth(chain.object());
489: if(depth >= 0) {
490: this->nClosurePreSieve(chain,depth,closure);
491: }
492: ALE_LOG_STAGE_END;
493: return closure;
494: }// Sieve::closureSieve()
496: #undef __FUNCT__
498: Point_set Sieve::star(Point_set& chain) {
499: Point_set star;
500: ALE_LOG_STAGE_BEGIN;
501: CHKCOMM(*this);
502: int32_t height = this->maxHeight(chain);
503: if(height >= 0) {
504: star = this->nStar(chain,height);
505: }
506: ALE_LOG_STAGE_END;
507: return star;
508: }// Sieve::star()
510: #undef __FUNCT__
512: Obj<Sieve> Sieve::starSieve(Obj<Point_set> chain, Obj<Sieve> star) {
513: ALE_LOG_STAGE_BEGIN;
514: CHKCOMM(*this);
515: if(star.isNull()) {
516: star = Obj<Sieve>(new Sieve(this->getComm()));
517: }
518: int32_t height = this->maxHeight(chain.object());
519: if(height >= 0) {
520: this->nStarPreSieve(chain,height,star);
521: }
522: ALE_LOG_STAGE_END;
523: return star;
524: }// Sieve::starSieve()
526: #undef __FUNCT__
528: Point_set Sieve::meet(Point_set c0, Point_set c1) {
529: // The strategy is to compute the intersection of cones over the chains, remove the intersection
530: // and use the remaining two parts -- two disjoined components of the symmetric difference of cones -- as the new chains.
531: // The intersections at each stage are accumulated and their union is the meet.
532: // The iteration stops when at least one of the chains is empty.
533: Point_set meet;
534: ALE_LOG_STAGE_BEGIN;
535: // Check if any the initial chains may be empty, so that we don't perform spurious iterations
536: if((c0.size() == 0) || (c1.size() == 0)) {
537: // return meet;
538: }
539: else {
540: while(1) {
541: Point_set *c = &c0;
542: Point_set *cc = &c1;
543: // Traverse the smallest cone set
544: if(cc->size() < c->size()) {
545: Point_set *tmp = c; c = cc; cc = tmp;
546: }
547: // Compute the intersection of c & cc and put it in meet at the same time removing it from c and cc
548: for(Point_set::iterator c_itor = c->begin(); c_itor != c->end(); c_itor++) {
549: if(cc->find(*c_itor)!= cc->end()) {
550: meet.insert(*c_itor);
551: cc->erase(*c_itor);
552: c->erase(c_itor);
553: }
554: }// for(Point_set::iterator c_itor = c->begin(); c_itor != c->end(); c_itor++)
555: // Replace each of the cones with a cone over it, and check if either is empty; if so, return what's in meet at the moment.
556: c0 = this->cone(c0);
557: if(c0.size() == 0) {
558: break;
559: }
560: c1 = this->cone(c1);
561: if(c1.size() == 0) {
562: break;
563: }
564: }// while(1)
565: }
566: ALE_LOG_STAGE_END;
567: return meet;
568: }// Sieve::meet()
570: #undef __FUNCT__
572: Point_set Sieve::meetAll(Point_set& chain) {
573: // The strategy is the same as in meet, except it is performed on an array of chains/cones--one per point in chain--simultaneously.
574: // This may be faster than applying 'meet' recursively, since intersections may become empty faster.
575: Point_set meets;
576: // Populate the 'cones' map, while checking if any of the initial cones may be empty,
577: // so that we don't perform spurious iterations. At the same time determine the cone of the smallest size.
578: Point__Point_set cones;
579: Point minp = *(chain.begin());
580: for(Point_set::iterator chain_itor = chain.begin(); chain_itor != chain.end(); chain_itor++) {
581: Point p = *chain_itor;
582: cones[p] = this->cone(p);
583: if(cones[p].size() == 0) {
584: return meets;
585: }
586: if(cones[p].size() < cones[minp].size()) {
587: minp = p;
588: }
589: }// for(Point_set::iterator chain_itor = chain.begin(); chain_itor != chain.end(); chain_itor++)
590:
591: while(1) {
592: Point_set *c = &cones[minp];
593: // Traverse the smallest cone set pointed to by c and compute the intersection of c with the other cones,
594: // putting it in meet at the same time removing it from c and the other cones.
595: for(Point_set::iterator c_itor = c->begin(); c_itor != c->end(); c_itor++) {
596: Point q = *c_itor;
597: // Keep a flag indicating whether q belongs to the intersection of all cones.
598: int qIsIn = 1;
599: for(Point_set::iterator chain_itor = chain.begin(); chain_itor != chain.end(); chain_itor++) {
600: Point p = *chain_itor;
601: // skip p equal to minp
602: if(p != minp) {
603: if(cones[p].find(q) == cones[p].end()) {// q is absent from the cone over p
604: qIsIn = 0;
605: break;
606: }
607: }
608: }// for(Point_set::iterator chain_itor = chain.begin(); chain_itor != chain.end(); chain_itor++)
609: // if a is in, add it to the meet and remove from all the cones
610: if(qIsIn) {
611: meets.insert(q);
612: for(Point_set::iterator chain_itor = chain.begin(); chain_itor != chain.end(); chain_itor++) {
613: Point p = *chain_itor;
614: cones[p].erase(q);
615: }
616: }
617: // Recall that erase(q) should not invalidate c_itor
618: }// for(Point_set::iterator c_itor = c->begin(); c_itor != c->end(); c_itor++)
620: // Replace each of the cones with a cone over it, and check if either is empty; if so, we are done -- return meets.
621: // At the same time, locate the point with the smallest cone over it.
622: for(Point_set::iterator chain_itor = chain.begin(); chain_itor != chain.end(); chain_itor++) {
623: Point p = *chain_itor;
624: // We check before and after computing the cone
625: if(cones[p].size() == 0) {
626: return meets;
627: }
628: cones[p] = this->cone(cones[p]);
629: if(cones[p].size() == 0) {
630: return meets;
631: }
632: if(cones[p].size() < cones[minp].size()) {
633: minp = p;
634: }
635: }// for(Point_set::iterator chain_itor = chain.begin(); chain_itor != chain.end(); chain_itor++)
637: }// while(1)
638: return meets;
639: }// Sieve::meetAll()
642: #undef __FUNCT__
644: Point_set Sieve::join(Point_set c0, Point_set c1) {
645: // The strategy is to compute the intersection of the supports of the two chains, remove the intersection
646: // and use the remaining two parts -- two disjoined components of the symmetric difference of the supports -- as the new chains.
647: // The intersections at each stage are accumulated and their union is the join.
648: // The iteration stops when at least one of the chains is empty.
649: Point_set join;
650: ALE_LOG_STAGE_BEGIN;
651: // Check if any the initial chains may be empty, so that we don't perform spurious iterations
652: if((c0.size() == 0) || (c1.size() == 0)) {
653: //return join;
654: }
655: else{
656: while(1) {
657: Point_set *s = &c0;
658: Point_set *ss = &c1;
659: // Traverse the smallest supp set
660: if(ss->size() < s->size()) {
661: Point_set *tmp = s; s = ss; ss = tmp;
662: }
663: // Compute the intersection of s & ss and put it in join at the same time removing it from s and ss
664: for(Point_set::iterator s_itor = s->begin(); s_itor != s->end(); s_itor++) {
665: if(ss->find(*s_itor)!= ss->end()) {
666: join.insert(*s_itor);
667: ss->erase(*s_itor);
668: s->erase(s_itor);
669: }
670: }// for(Point_set::iterator s_itor = s->begin(); s_itor != s->end(); s_itor++)
671: // Replace each of the chains with its support, and check if either is empty; if so, stop
672: c0 = this->support(c0);
673: if(c0.size() == 0) {
674: break;
675: }
676: c1 = this->support(c1);
677: if(c1.size() == 0) {
678: break;
679: }
680: }// while(1)
681: }
682: ALE_LOG_STAGE_END;
683: return join;
684: }// Sieve::join()
687: #undef __FUNCT__
689: Point_set Sieve::joinAll(Point_set& chain) {
690: // The strategy is the same as in join, except it is performed on an array of chains/supps--one per point in chain--simultaneously.
691: // This may be faster than applying 'join' recursively, since intersections may become empty faster.
692: Point_set joins;
693: // Populate the 'supps' map, while checking if any of the initial supports may be empty,
694: // so that we don't perform spurious iterations. At the same time determine the supp of the smallest size.
695: Point__Point_set supps;
696: Point minp = *(chain.begin());
697: for(Point_set::iterator chain_itor = chain.begin(); chain_itor != chain.end(); chain_itor++) {
698: Point p = *chain_itor;
699: supps[p] = this->support(p);
700: if(supps[p].size() == 0) {
701: return joins;
702: }
703: if(supps[p].size() < supps[minp].size()) {
704: minp = p;
705: }
706: }// for(Point_set::iterator chain_itor = chain.begin(); chain_itor != chain.end(); chain_itor++)
707:
708: while(1) {
709: Point_set *s = &supps[minp];
710: // Traverse the smallest supp set pointed to by s and compute the intersection of s with the other supps,
711: // putting it in join at the same time removing it from s and the other supps.
712: for(Point_set::iterator s_itor = s->begin(); s_itor != s->end(); s_itor++) {
713: Point q = *s_itor;
714: // Keep a flag indicating whether q belongs to the intersection of all supps.
715: int qIsIn = 1;
716: for(Point_set::iterator chain_itor = chain.begin(); chain_itor != chain.end(); chain_itor++) {
717: Point p = *chain_itor;
718: // skip p equal to minp
719: if(p != minp) {
720: if(supps[p].find(q) == supps[p].end()) {// q is absent from the support of p
721: qIsIn = 0;
722: break;
723: }
724: }
725: }// for(Point_set::iterator chain_itor = chain.begin(); chain_itor != chain.end(); chain_itor++)
726: // if a is in, add it to the join and remove from all the supps
727: if(qIsIn) {
728: joins.insert(q);
729: for(Point_set::iterator chain_itor = chain.begin(); chain_itor != chain.end(); chain_itor++) {
730: Point p = *chain_itor;
731: supps[p].erase(q);
732: }
733: }
734: // Recall that erase(q) should not invalidate s_itor
735: }// for(Point_set::iterator s_itor = s->begin(); s_itor != s->end(); s_itor++)
737: // Replace each of the supps with its support, and check if any of the new supps is empty; if so, stop.
738: // At the same time, locate the point with the smallest supp.
739: for(Point_set::iterator chain_itor = chain.begin(); chain_itor != chain.end(); chain_itor++) {
740: Point p = *chain_itor;
741: // We check before and after computing the support
742: if(supps[p].size() == 0) {
743: return joins;
744: }
745: supps[p] = this->support(supps[p]);
746: if(supps[p].size() == 0) {
747: return joins;
748: }
749: if(supps[p].size() < supps[minp].size()) {
750: minp = p;
751: }
752: }// for(Point_set::iterator chain_itor = chain.begin(); chain_itor != chain.end(); chain_itor++)
754: }// while(1)
755: return joins;
756: }// Sieve::joinAll()
758: #undef __FUNCT__
760: Point_set Sieve::roots(Point_set chain) {
761: // Compute the roots (nodes with empty cones over them) from which chain is reacheable.
762: // The strategy is to examine each node in the current chain, replace it by its cone, if the latter is non-empty,
763: // or remove the node from the chain, if it is a root, while saving it in the roots set.
764: Point_set roots;
765: ALE_LOG_STAGE_BEGIN;
766: Point_set cone;
767: Point_set *c = &chain;
768: Point_set *cc = &cone;
769: int stop = 0;
770: while(!stop) {
771: stop = 1; // presume that everything in c is a root
772: // Traverse chain and check the size of the cone over each point.
773: for(Point_set::iterator c_itor = c->begin(); c_itor != c->end(); c_itor++) {
774: Point p = *c_itor;
775: if(this->coneSize(p) == 0){
776: roots.insert(p);
777: }
778: else {
779: // add the cone over p to cc
780: Point_set pCone = this->cone(p);
781: for(Point_set::iterator pCone_itor = pCone.begin(); pCone_itor != pCone.end(); pCone_itor++) {
782: cc->insert(*pCone_itor);
783: }
784: // turn off the stopping flag, as we are going to examine the just-added cone during the next iteration
785: stop = 0;
786: }// if(this->coneSize(p) == 0)
787: }// for(Point_set::iterator c_itor = c->begin(); c_itor != c->end(); c_itor++)
788: // swap c and cc
789: Point_set *tmp = c; c = cc; cc = tmp;
790: }// while(!stop)
791: ALE_LOG_STAGE_END;
792: return roots;
793: }// Sieve::roots()
795: #undef __FUNCT__
797: Point_set Sieve::leaves(Point_set chain) {
798: // Compute the leaves (nodes with empty supports) reacheable from chain.
799: // The strategy is to examine each node in the current chain, replace it by its support, if the latter is non-empty,
800: // or remove the node from the chain, if it is a leaf, while saving it in the leaves set.
801: Point_set leaves;
802: ALE_LOG_STAGE_BEGIN;
803: Point_set supp;
804: Point_set *c = &chain;
805: Point_set *cc = &supp;
806: int stop = 0;
807: while(!stop) {
808: stop = 1; // presume that everything in c is a leaf
809: // Traverse chain and check the size of each point's support.
810: for(Point_set::iterator c_itor = c->begin(); c_itor != c->end(); c_itor++) {
811: Point p = *c_itor;
812: if(this->supportSize(p) == 0){
813: leaves.insert(p);
814: }
815: else {
816: // add the support of p to cc
817: Point_set pSupp = this->support(p);
818: for(Point_set::iterator pSupp_itor = pSupp.begin(); pSupp_itor != pSupp.end(); pSupp_itor++) {
819: cc->insert(*pSupp_itor);
820: }
821: // turn off the stopping flag, as we are going to examine the just-added support during the next iteration
822: stop = 0;
823: }// if(this->supportSize(p) == 0)
824: }// for(Point_set::iterator c_itor = c->begin(); c_itor != c->end(); c_itor++)
825: // swap c and cc
826: Point_set *tmp = c; c = cc; cc = tmp;
827: }// while(!stop)
828: ALE_LOG_STAGE_END;
829: return leaves;
830: }// Sieve::leaves()
833: #undef __FUNCT__
835: Point_set Sieve::depthStratum(int32_t depth) {
836: Point_set stratum;
837: ALE_LOG_STAGE_BEGIN;
838: CHKCOMM(*this);
839: for(Point__int::iterator d_itor = this->_depth.begin(); d_itor != this->_depth.end(); d_itor++) {
840: if (d_itor->second == depth) {
841: stratum.insert(d_itor->first);
842: }
843: }
844: ALE_LOG_STAGE_END;
845: return stratum;
846: }// Sieve::depthStratum()
848: #undef __FUNCT__
850: Point_set Sieve::heightStratum(int32_t height) {
851: Point_set stratum;
852: ALE_LOG_STAGE_BEGIN;
853: CHKCOMM(*this);
854: for(Point__int::iterator h_itor = this->_height.begin(); h_itor != this->_height.end(); h_itor++) {
855: if (h_itor->second == height) {
856: stratum.insert(h_itor->first);
857: }
858: }
859: ALE_LOG_STAGE_END;
860: return stratum;
861: }// Sieve::heightStratum()
865: void Sieve::view(const char *name) {
866: CHKCOMM(*this);
867: int32_t rank = this->commRank;
869: ostringstream txt, hdr;
870: hdr << "Viewing ";
871: if(name == NULL) {
872: hdr << "a ";
873: }
874: if(!this->isLocked()) {
875: hdr << "(locked) ";
876: }
877: if(name != NULL) {
878: hdr << "sieve '" << name << "' (square brackets contain depth, height pairs)\n";
879: }
880: hdr << "\n";
881: // Print header
882: PetscSynchronizedPrintf(this->comm, hdr.str().c_str()); CHKERROR(ierr, "Error in PetscPrintf");
883: // Use a string stream to accumulate output that is then submitted to PetscSynchronizedPrintf
884: Point_set points = this->space();
885: txt << "[" << rank << "]: space of size " << points.size() << " : ";
886: for(Point_set::iterator p_itor = points.begin(); p_itor != points.end(); p_itor++)
887: {
888: Point p = (*p_itor);
889: if(p_itor != points.begin()) {
890: txt << ", ";
891: }
892: txt << "(" << p.prefix << ", " << p.index << ")["<<this->depth(p) << ", " << this->height(p) << "]";
893: }
894: txt << "\n";
895: //
896: points = this->cap();
897: txt << "[" << rank << "]: cap of size " << points.size() << " : ";
898: for(Point_set::iterator p_itor = points.begin(); p_itor != points.end(); p_itor++)
899: {
900: Point p = (*p_itor);
901: if(p_itor != points.begin()) {
902: txt << ", ";
903: }
904: txt << "(" << p.prefix << ", " << p.index << ")";
905: }
906: txt << "\n";
907: //
908: points = this->base();
909: txt << "[" << rank << "]: base of size " << points.size() << " : ";
910: for(Point_set::iterator p_itor = points.begin(); p_itor != points.end(); p_itor++)
911: {
912: Point p = (*p_itor);
913: if(p_itor != points.begin()) {
914: txt << ", ";
915: }
916: txt << "(" << p.prefix << ", " << p.index << ")";
917: }
918: txt << "\n";
919: //
920: for(Point__Point_set::iterator cone_itor = this->_cone.begin(); cone_itor != this->_cone.end(); cone_itor++)
921: {
922: Point p = (*cone_itor).first;
923: txt << "[" << rank << "]: cone over ("<<p.prefix<<", "<<p.index<<"): ";
924: // Traverse the local cone over p
925: for(Point_set::iterator pCone_itor = this->_cone[p].begin(); pCone_itor != this->_cone[p].end(); pCone_itor++) {
926: Point q = *pCone_itor;
927: if(pCone_itor != this->_cone[p].begin()) {
928: txt << ", ";
929: }
930: txt << "(" << q.prefix << ", " << q.index << ")";
931: }
932: txt << "\n";
933: }
934: for(Point__Point_set::iterator support_itor = this->_support.begin(); support_itor != this->_support.end(); support_itor++)
935: {
936: Point p = (*support_itor).first;
937: txt << "[" << rank << "]: support of ("<<p.prefix<<", "<<p.index<<"): ";
938: // Traverse the local support of p
939: for(Point_set::iterator pSupport_itor = this->_support[p].begin(); pSupport_itor != this->_support[p].end(); pSupport_itor++) {
940: Point q = *pSupport_itor;
941: if(pSupport_itor != this->_support[p].begin()) {
942: txt << ", ";
943: }
944: txt << "(" << q.prefix << ", " << q.index << ")";
945: }
946: txt << "\n";
947: }
948: PetscSynchronizedPrintf(this->comm, txt.str().c_str());
949: CHKERROR(ierr, "Error in PetscSynchronizedPrintf");
950: PetscSynchronizedFlush(this->comm);
951: CHKERROR(ierr, "Error in PetscSynchronizedFlush");
952: #if 0
953: ostringstream heighthdr;
954: heighthdr << "Height of ";
955: if (name == NULL) {
956: heighthdr << "the Sieve";
957: } else {
958: heighthdr << name;
959: }
960: heighthdr << std::endl;
961: this->_height.view(heighthdr.str().c_str());
962: #endif
963: }// Sieve::view()
964:
965: #undef __FUNCT__
967: Sieve* Sieve::baseRestriction(Point_set& base) {
968: Sieve *s;
969: ALE_LOG_STAGE_BEGIN;
970: CHKCOMM(*this);
971: s = new Sieve(this->getComm());
972: for(Point_set::iterator b_itor = base.begin(); b_itor != base.end(); b_itor++){
973: Point p = *b_itor;
974: // is point p present in the base of *this?
975: if(this->_cone.find(p) != this->_cone.end()){
976: s->addCone(this->_cone[p],p);
977: }
978: }// for(Point_set::iterator b_itor = base.begin(); b_itor != base.end(); b_itor++){
979: ALE_LOG_STAGE_END;
980: return s;
981: }// Sieve::baseRestriction()
983: #undef __FUNCT__
985: Sieve* Sieve::capRestriction(Point_set& cap) {
986: Sieve *s;
987: ALE_LOG_STAGE_BEGIN;
988: CHKCOMM(*this);
989: s = new Sieve(this->getComm());
990: for(Point_set::iterator c_itor = cap.begin(); c_itor != cap.end(); c_itor++){
991: Point q = *c_itor;
992: // is point q present in the cap of *this?
993: if(this->_support.find(q) == this->_support.end()){
994: s->addSupport(q,this->_support[q]);
995: }
996: }// for(Point_set::iterator c_itor = cap.begin(); c_itor != cap.end(); c_itor++){
997: ALE_LOG_STAGE_END;
998: return s;
999: }// Sieve::capRestriction()
1001: } // namespace ALE
1003: #undef ALE_Sieve_cxx