Actual source code: PreSieve.cxx
1: #define ALE_PreSieve_cxx
3: #ifndef included_ALE_PreSieve_hh
4: #include <PreSieve.hh>
5: #endif
7: #ifndef included_ALE_Stack_hh
8: #include <Stack.hh>
9: #endif
11: #include <stack>
12: #include <map>
14: namespace ALE {
16:
17: #undef __FUNCT__
19: PreSieve::PreSieve() : Coaster(), _cone(), _support() {};
21: #undef __FUNCT__
23: PreSieve::PreSieve(MPI_Comm comm) : Coaster(comm), _cone(), _support() {};
25: #undef __FUNCT__
27: PreSieve::~PreSieve(){};
30: #undef __FUNCT__
32: PreSieve& PreSieve::clear(){
33: /* This is not what we mean by clear: Coaster::clear(); */
34: this->_cone.clear();
35: this->_support.clear();
36: return *this;
37: }// PreSieve::clear()
38:
39: #undef __FUNCT__
41: PreSieve& PreSieve::addBasePoint(const Point& p) {
42: this->__checkLock();
43: if(!this->baseContains(p)) {
44: this->_cone[p] = Point_set();
45: // This is a little counterintuitive, but after the initial addition of a base point,
46: // it is initially a root -- no incoming arrows.
47: this->_roots.insert(p);
48: // If the point is not in the cap, it is also a leaf.
49: // It could be a leaf even while being in the cap, but then addCapPoint would take care of that.
50: if(!capContains(p)) {
51: this->_leaves.insert(p);
52: }
53: }
54: return *this;
55: }// PreSieve::addBasePoint()
57: #undef __FUNCT__
59: PreSieve& PreSieve::removeBasePoint(const Point& p, bool removeSingleton) {
60: this->__checkLock();
61: if (this->_cone.find(p) != this->_cone.end()) {
62: // IMPROVE: use 'coneView' and iterate over it to avoid copying the set
63: ALE::Obj<ALE::Point_set> cone = this->cone(p);
64: for(ALE::Point_set::iterator c_itor = cone->begin(); c_itor != cone->end(); c_itor++) {
65: ALE::Point cover = *c_itor;
67: this->removeArrow(cover, p, removeSingleton);
68: }
69: // Remove the cone over p
70: this->_cone.erase(p);
71: // After the removal of p from the base, if p is still in the cap, it must necessarily be a root,
72: // since there are no longer any arrows terminating at p.
73: if(this->capContains(p)) {
74: this->_roots.insert(p);
75: } else {
76: // otherwise the point is gone from the PreSieve, and is not a root, hence it better not be a leaf either
77: this->_leaves.erase(p);
78: }
79: }
80: return *this;
81: }// PreSieve::removeBasePoint()
83: #undef __FUNCT__
85: PreSieve& PreSieve::addCapPoint(const Point& q) {
86: CHKCOMM(*this);
87: this->__checkLock();
88: if(!this->capContains(q)) {
89: this->_support[q] = Point_set();
90: // This may appear counterintuitive, but upon initial addition to the cap the point is a leaf
91: this->_leaves.insert(q);
92: // If the point is not in the base, it is also a root. It may be a root while being in the base as well,
93: // but that is handled by addBasePoint
94: if(!this->baseContains(q)) {
95: this->_roots.insert(q);
96: }
97: }
98: return *this;
99: }// PreSieve::addCapPoint()
101: #undef __FUNCT__
103: PreSieve& PreSieve::removeCapPoint(const Point& q, bool removeSingleton) {
104: CHKCOMM(*this);
105: this->__checkLock();
106: if(this->capContains(q)) {
107: // IMPROVE: use 'supportView' and iterate over that, instead of copying a set
108: ALE::Obj<ALE::Point_set> support = this->support(q);
109: for(ALE::Point_set::iterator s_itor = support->begin(); s_itor != support->end(); s_itor++) {
110: ALE::Point s = *s_itor;
112: this->removeArrow(q, s, removeSingleton);
113: // We do not erase points from 'base' unless explicitly removed, even if no points terminate there
114: // WARNING: the following code looks suspicious
115: // if (this->cone(s).size() == 0) {
116: // if (!this->capContains(s)) {
117: // this->removeBasePoint(s);
118: // } else {
119: // this->_cone.erase(s);
120: // }
121: // }
122: }
123: // Remove the point from the cap and delete it the support under q
124: this->_support.erase(q);
125: // After removal from the cap the point that is still in the base becomes a leaf -- no outgoing arrows.
126: if(this->baseContains(q)) {
127: this->_leaves.insert(q);
128: } else {
129: // otherwise the point is gone from the PreSieve, and is not a leaf, hence it must not be a root either
130: this->_roots.erase(q);
131: }
132: }
133: return *this;
134: }// PreSieve::removeCapPoint()
136: #undef __FUNCT__
138: PreSieve& PreSieve::addArrow(const Point& i, const Point& j) {
139: ALE_LOG_STAGE_BEGIN;
140: CHKCOMM(*this);
141: this->__checkLock();
142: this->addBasePoint(j);
143: this->addCapPoint(i);
144: this->_cone[j].insert(i);
145: this->_support[i].insert(j);
146: this->_leaves.erase(i);
147: this->_roots.erase(j);
148: ALE_LOG_STAGE_END;
149: return *this;
150: }// PreSieve::addArrow()
152: #undef __FUNCT__
154: PreSieve& PreSieve::removeArrow(const Point& i, const Point& j, bool removeSingleton) {
155: ALE_LOG_STAGE_BEGIN;
156: CHKCOMM(*this);
157: this->__checkLock();
158: this->_cone[j].erase(i);
159: this->_support[i].erase(j);
160: // If this was the last arrow terminating at j, it becomes a root
161: if(this->coneSize(j) == 0) {
162: if (removeSingleton && (this->_leaves.find(j) != this->_leaves.end())) {
163: this->removePoint(j);
164: } else {
165: this->_roots.insert(j);
166: }
167: }
168: // If this was the last arrow emanating from i, it becomes a root
169: if(this->supportSize(i) == 0) {
170: if (removeSingleton && (this->_roots.find(i) != this->_roots.end())) {
171: this->removePoint(i);
172: } else {
173: this->_leaves.insert(i);
174: }
175: }
176: ALE_LOG_STAGE_END;
177: return *this;
178: }// PreSieve::removeArrow()
180: #undef __FUNCT__
182: PreSieve& PreSieve::invert() {
183: CHKCOMM(*this);
184: this->__checkLock();
185: // We keep track of the inversion in an new PreSieve object
186: PreSieve inv;
188: // Traverse the base of 'this'
189: for(Point__Point_set::iterator base_itor = this->_cone.begin(); base_itor != this->_cone.end(); base_itor++) {
190: Point p = base_itor->first;
191: // Traverse the cone over p in 'this'
192: for(Point_set::iterator pCone_itor = this->_cone[p].begin(); pCone_itor != this->_cone[p].end(); pCone_itor++) {
193: Point q = *pCone_itor;
194: inv.addArrow(p,q); // inversion takes place here
195: }
196: }
197: // Replace *this with inv
198: *this = inv;
199: return *this;
200: }// PreSieve::invert()
201:
203: #undef __FUNCT__
205: PreSieve& PreSieve::stackBelow(PreSieve& s) {
206: CHKCOMMS(*this,s);
207: this->__checkLock();
208: // Composition is by 'stacking': s is placed above 'this' so that the base of s is on par with the cap of 'this'.
209: // Then the arrows from s' cap are composed with 'this's arrows below.
210: // This can also be thought of as the 'right action by s': *this = *this o s
211: // so that the domain of the composition is the domain (cap) of s, and the range of the composition is the range (base) of 'this'.
212: // This is done so we can keep the base and hence the number of cones in this->_cone, although some cones may become
213: // empty upon composition.
215: // We keep the result of the composition in a new PreSieve C
216: PreSieve C(this->getComm());
217: // To construct a cone over the base points in the new PreSieve, we traverse the cones of *this.
218: for(Point__Point_set::iterator base_itor = this->_cone.begin(); base_itor != this->_cone.end(); base_itor++) {
219: Point p = (*base_itor).first;
220: // Add p to the base of C ahead of time in case the cone ends up being empty
221: C.addBasePoint(p);
222: // Traverse the cone over p in this
223: for(Point_set::iterator pCone_itor=this->_cone[p].begin(); pCone_itor != this->_cone[p].end(); pCone_itor++) {
224: Point q = *pCone_itor;
225: // For each q in the cone over p add it as the cone over p in C
226: Point_set qCone = s.cone(q);
227: C.addCone(qCone, p);
228: }// for(Point_set::iterator pCone_itor = this->_cone[p].begin(); pCone_itor != this->_cone[p].end(); pCone_itor++)
229: }// for(Point__Point_set::iterator base_itor = this->_cone.begin(); base_itor != this->_cone.end(); base_itor++)
230: // We use the cap of s for C
231: for(Point__Point_set::iterator cap_itor = s._support.begin(); cap_itor != s._support.end(); cap_itor++) {
232: Point q = (*cap_itor).first;
233: C.addCapPoint(q);
234: }
235: // Now C to *this
236: *this = C;
237: return *this;
238: }// PreSieve::stackBelow()
241: #undef __FUNCT__
243: PreSieve& PreSieve::stackAbove(PreSieve& s) {
244: CHKCOMMS(*this,s);
245: this->__checkLock();
246: // Composition is by 'stacking': s is placed below 'this' so that the cap of s is on par with the base of 'this'.
247: // Then the arrows from this's cap are composed with s's arrows below.
248: // This can also be thought of as the 'left action by s': *this = s o *this
249: // so that the domain of the composition is the domain (cap) of 'this', and the range of the composition
250: // is the range (base) of s.
252: // We keep the result of the composition in a new PreSieve ss
253: PreSieve ss;
254: // Now traverse s's base.
255: for(Point__Point_set::iterator base_itor = s._cone.begin(); base_itor != s._cone.end(); base_itor++) {
256: Point p = (*base_itor).first;
257: // Add p to the base of ss ahead of time in case the cone ends up being empty
258: ss.addBasePoint(p);
259: // Traverse the cone over p in s
260: for(Point_set::iterator pCone_itor=s._cone[p].begin(); pCone_itor != s._cone[p].end(); pCone_itor++) {
261: Point q = *pCone_itor;
262: // For each q in the cone over p in s take the cone over q in *this and add it as the cone over p in ss
263: Point_set qCone = this->cone(q);
264: ss.addCone(qCone,p);
265: }// for(Point_set::iterator pCone_itor = s._cone[p].begin(); pCone_itor != s._cone[p].end(); pCone_itor++)
266: }// for(Point__Point_set::iterator base_itor = s._cone.begin(); base_itor != s._cone.end(); base_itor++)
267: // We use the cap of *this for C
268: for(Point__Point_set::iterator cap_itor = this->_support.begin(); cap_itor != this->_support.end(); cap_itor++) {
269: Point q = (*cap_itor).first;
270: ss.addCapPoint(q);
271: }
272: // Now replace the guts of *this with those of ss
273: *this = ss;
274: return *this;
275: }// PreSieve::stackAbove()
277: #undef __FUNCT__
279: PreSieve& PreSieve::restrictBase(Point_set& base) {
280: this->__checkLock();
281: // IMPROVE: this should be improved to do the removal 'in-place', instead of creating a potentially huge 'removal' set
282: Point_set removal;
283: while(1) {
284: for(Point__Point_set::iterator b_itor = this->_cone.begin(); b_itor != this->_cone.end(); b_itor++){
285: Point p = b_itor->first;
286: // is point p absent from base?
287: if(base.find(p) == base.end()){
288: removal.insert(p);
289: }
290: }
291: if (!removal.size()) break;
292: for(Point_set::iterator r_itor = removal.begin(); r_itor != removal.end(); r_itor++){
293: Point p = *r_itor;
294: this->removeBasePoint(p, true);
295: }
296: removal.clear();
297: }
298: return *this;
299: }// PreSieve::restrictBase()
301: #undef __FUNCT__
303: PreSieve& PreSieve::excludeBase(Point_set& base) {
304: this->__checkLock();
305: // IMPROVE: this should be improved to do the removal 'in-place', instead of creating a potentially huge 'removal' set
306: Point_set removal;
307: for(Point__Point_set::iterator b_itor = this->_cone.begin(); b_itor != this->_cone.end(); b_itor++){
308: Point p = b_itor->first;
309: // is point p present in base?
310: if(base.find(p) != base.end()){
311: removal.insert(p);
312: }
313: }// for(Point_set::iterator b_itor = this->_cone.begin(); b_itor != this->_cone.end(); b_itor++){
314: for(Point_set::iterator r_itor = removal.begin(); r_itor != removal.end(); r_itor++){
315: Point p = *r_itor;
316: this->removeBasePoint(p);
317: }
318: return *this;
319: }// PreSieve::excludeBase()
321: #undef __FUNCT__
323: PreSieve* PreSieve::baseRestriction(Point_set& base) {
324: CHKCOMM(*this);
325: PreSieve *s = new PreSieve(this->getComm());
326: for(Point_set::iterator b_itor = base.begin(); b_itor != base.end(); b_itor++){
327: Point p = *b_itor;
328: // is point p present in the base of *this?
329: if(this->_cone.find(p) != this->_cone.end()){
330: s->addCone(this->_cone[p],p);
331: }
332: }// for(Point_set::iterator b_itor = base.begin(); b_itor != base.end(); b_itor++){
333: return s;
334: }// PreSieve::baseRestriction()
336: #undef __FUNCT__
338: PreSieve* PreSieve::baseExclusion(Point_set& base) {
339: CHKCOMM(*this);
340: PreSieve *s = new PreSieve(this->getComm());
341: this->__computeBaseExclusion(base, s);
342: return s;
343: }// PreSieve::baseExclusion()
345: #undef __FUNCT__
347: void PreSieve::__computeBaseExclusion(Point_set& base, PreSieve *s) {
348: // This function is used by Stack as well
349: for(Point__Point_set::iterator cone_itor = this->_cone.begin(); cone_itor != this->_cone.end(); cone_itor++){
350: Point p = cone_itor->first;
351: Point_set pCone = cone_itor->second;
352: // is point p absent from base?
353: if(base.find(p) == base.end()){
354: s->addCone(pCone,p);
355: }
356: }
357: }// PreSieve::__computeBaseExclusion()
359: #undef __FUNCT__
361: PreSieve& PreSieve::restrictCap(Point_set& cap) {
362: this->__checkLock();
363: // IMPROVE: this should be improved to do the removal 'in-place', instead of creating a potentially huge 'removal' set
364: Point_set removal;
365: for(Point__Point_set::iterator c_itor = this->_support.begin(); c_itor != this->_support.end(); c_itor++){
366: Point q = c_itor->first;
367: // is point q absent from cap?
368: if(cap.find(q) == cap.end()){
369: removal.insert(q);
370: }
371: }// for(Point_set::iterator c_itor = this->_support.begin(); c_itor != this->_support.end(); c_itor++){
372: for(Point_set::iterator r_itor = removal.begin(); r_itor != removal.end(); r_itor++){
373: Point r = *r_itor;
374: this->removeCapPoint(r, true);
375: }
376: return *this;
377: }// PreSieve::restrictCap()
379: #undef __FUNCT__
381: PreSieve& PreSieve::excludeCap(Point_set& cap) {
382: this->__checkLock();
383: // IMPROVE: this should be improved to do the removal 'in-place', instead of creating a potentially huge 'removal' set
384: Point_set removal;
385: for(Point__Point_set::iterator c_itor = this->_support.begin(); c_itor != this->_support.end(); c_itor++){
386: Point q = c_itor->first;
387: // is point q present in cap?
388: if(cap.find(q) != cap.end()){
389: removal.insert(q);
390: }
391: }// for(Point_set::iterator c_itor = this->_support.begin(); c_itor != this->_support.end(); c_itor++){
392: for(Point_set::iterator r_itor = removal.begin(); r_itor != removal.end(); r_itor++){
393: Point r = *r_itor;
394: this->removeCapPoint(r);
395: }
396: return *this;
397: }// PreSieve::excludeCap()
400: #undef __FUNCT__
402: PreSieve* PreSieve::capRestriction(Point_set& cap) {
403: CHKCOMM(*this);
404: PreSieve *s = new PreSieve(this->getComm());
405: for(Point_set::iterator c_itor = cap.begin(); c_itor != cap.end(); c_itor++){
406: Point q = *c_itor;
407: // is point q present in the cap of *this?
408: if(this->_support.find(q) != this->_support.end()){
409: s->addSupport(q,this->_support[q]);
410: }
411: }// for(Point_set::iterator c_itor = cap.begin(); c_itor != cap.end(); c_itor++){
412: return s;
413: }// PreSieve::capRestriction()
416: #undef __FUNCT__
418: PreSieve* PreSieve::capExclusion(Point_set& cap) {
419: CHKCOMM(*this);
420: PreSieve *s = new PreSieve(this->getComm());
421: this->__computeCapExclusion(cap, s);
422: return s;
423: }// PreSieve::capExclusion()
425: #undef __FUNCT__
427: void PreSieve::__computeCapExclusion(Point_set& cap, PreSieve *s) {
428: // This function is used by Stack as well
429: for(Point__Point_set::iterator support_itor = this->_support.begin(); support_itor != this->_support.end(); support_itor++){
430: Point q = support_itor->first;
431: Point_set qSupport = support_itor->second;
432: // is point q absent from cap?
433: if(cap.find(q) == cap.end()){
434: s->addSupport(q,qSupport);
435: }
436: }
437: }// PreSieve::__computeCapExclusion()
440: #undef __FUNCT__
442: Obj<Point_set> PreSieve::space() {
443: Obj<Point_set> space = Point_set();
444: ALE_LOG_STAGE_BEGIN;
445: CHKCOMM(*this);
446: // Adding both cap and base
448: for(Point__Point_set::iterator cap_itor = this->_support.begin(); cap_itor != this->_support.end(); cap_itor++){
449: space->insert(cap_itor->first);
450: }
451: for(Point__Point_set::iterator base_itor = this->_cone.begin(); base_itor != this->_cone.end(); base_itor++){
452: space->insert(base_itor->first);
453: }
454: ALE_LOG_STAGE_END;
455: return space;
456: }// PreSieve::space()
459: #undef __FUNCT__
461: Point_set PreSieve::base() {
462: Point_set base;
463: ALE_LOG_STAGE_BEGIN;
464: CHKCOMM(*this);
465: for(Point__Point_set::iterator cone_itor = this->_cone.begin(); cone_itor != this->_cone.end(); cone_itor++) {
466: base.insert((*cone_itor).first);
467: }
468: ALE_LOG_STAGE_END;
469: return base;
470: }// PreSieve::base()
472: #undef __FUNCT__
474: Point_set PreSieve::cap() {
475: Point_set cap;
476: ALE_LOG_STAGE_BEGIN;
477: CHKCOMM(*this);
478: for(Point__Point_set::iterator support_itor = this->_support.begin(); support_itor != this->_support.end(); support_itor++) {
479: cap.insert((*support_itor).first);
480: }
481: ALE_LOG_STAGE_END;
482: return cap;
483: }// PreSieve::cap()
486: #undef __FUNCT__
488: int32_t PreSieve::spaceSize() {
489: CHKCOMM(*this);
490: return this->space()->size();
491: }// PreSieve::spaceSize()
494: #undef __FUNCT__
496: int32_t *PreSieve::spaceSizes() {
497: CHKCOMM(*this);
498: // Allocate array of size commSize
499: int32_t spaceSize = this->space()->size();
500: int32_t *spaceSizes;
502: spaceSizes = (int32_t *) malloc(sizeof(int32_t)*this->commSize);
503: MPI_Allgather(&spaceSize, 1, MPIU_INT, spaceSizes, 1, MPIU_INT, comm); CHKERROR(ierr, "Error in MPI_Allgather");
504: return spaceSizes;
505: }// PreSieve::spaceSizes()
508: #undef __FUNCT__
510: int32_t PreSieve::baseSize() {
511: CHKCOMM(*this);
512: return this->_cone.size();
513: }// PreSieve::baseSize()
515: #undef __FUNCT__
517: int32_t *PreSieve::baseSizes() {
518: CHKCOMM(*this);
519: // Allocate array of size commSize
520: int32_t baseSize = this->baseSize();
521: int32_t *baseSizes;
523: baseSizes = (int32_t*)malloc(sizeof(int32_t)*this->commSize);
524: MPI_Allgather(&baseSize, 1, MPIU_INT, baseSizes, 1, MPIU_INT, comm); CHKERROR(ierr, "Error in MPI_Allgather");
525: return baseSizes;
526: }// PreSieve::baseSizes()
528: #undef __FUNCT__
530: int32_t PreSieve::capSize() {
531: CHKCOMM(*this);
532: return this->_support.size();
533: }// PreSieve::capSize()
535: #undef __FUNCT__
537: int32_t *PreSieve::capSizes() {
538: CHKCOMM(*this);
539: // Allocate array of size commSize
540: int32_t capSize = this->capSize();
541: int32_t *capSizes;
543: capSizes = (int32_t *)malloc(sizeof(int32_t)*this->commSize);
544: MPI_Allgather(&capSize, 1, MPIU_INT, capSizes, 1, MPIU_INT, comm); CHKERROR(ierr, "Error in MPI_Allgather");
545: return capSizes;
546: }// PreSieve::capSizes()
549: #undef __FUNCT__
551: int32_t PreSieve::coneSize(Point_set& c) {
552: CHKCOMM(*this);
553: int32_t coneSize = 0;
554: for(Point_set::iterator c_itor = c.begin(); c_itor != c.end(); c_itor++) {
555: Point p = *c_itor;
557: if (this->_cone.find(p) != this->_cone.end()) {
558: coneSize += this->_cone[p].size();
559: }
560: }
561: return coneSize;
562: }// PreSieve::coneSize()
565: #undef __FUNCT__
567: int32_t PreSieve::supportSize(const Point_set& c) {
568: CHKCOMM(*this);
569: int32_t supportSize = 0;
570: for(Point_set::iterator c_itor = c.begin(); c_itor != c.end(); c_itor++) {
571: Point q = *c_itor;
572: supportSize += this->_support[q].size();
573: }
574: return supportSize;
575: }// PreSieve::supportSize()
578: #undef __FUNCT__
580: int PreSieve::spaceContains(Point point) {
581: CHKCOMM(*this);
582: int flag;
583: flag = (this->capContains(point) || this->baseContains(point));
584: return flag;
585: }// PreSieve::spaceContains()
586:
587: #undef __FUNCT__
589: int PreSieve::baseContains(Point point) {
590: CHKCOMM(*this);
591: int flag;
592: //flag = ((this->_cone.find(point) != this->_cone.end()) && (this->_cone[point].size() > 0));
593: // SEMANTICS: we assume that the point is in the base regardless of whether it actually has arrows terminating at it
594: flag = ((this->_cone.find(point) != this->_cone.end()));
595: return flag;
596: }// PreSieve::baseContains()
598: #undef __FUNCT__
600: bool PreSieve::capContains(const Point& point) {
601: CHKCOMM(*this);
602: bool flag;
603: // SEMANTICS: we assume that the point is in the cap regardless of whether it actually has arrows emanating from it
604: //flag = ((this->_support.find(point) != this->_support.end()) && (this->_support[point].size() > 0));
605: flag = ((this->_support.find(point) != this->_support.end()));
606: return flag;
607: }// PreSieve::capContains()
608:
610: #undef __FUNCT__
612: int PreSieve::coneContains(Point& p, Point point) {
613: CHKCOMM(*this);
614: if (this->_cone.find(p) == this->_cone.end()) {
615: return 0;
616: }
617: Point_set& pCone = this->_cone[p];
618: int flag = 0;
619: if(pCone.find(point) != pCone.end()) {
620: flag = 1;
621: }
622: return flag;
623: }// PreSieve::coneContains()
625: #undef __FUNCT__
627: int PreSieve::supportContains(Point& q, Point point) {
628: CHKCOMM(*this);
629: int flag = 0;
630: Point_set& qSupport = this->_support[q];
631: if(qSupport.find(point) != qSupport.end()) {
632: flag = 1;
633: }
634: return flag;
635: }// PreSieve::supportContains()
637: #undef __FUNCT__
639: void PreSieve::addCone(const Point& cone, const Point& j) {
640: ALE_LOG_EVENT_BEGIN
641: CHKCOMM(*this);
642: this->__checkLock();
643: // Add j to the base in case coneSet is empty and no addArrow are executed
644: this->addBasePoint(j);
645: this->addArrow(cone, j);
646: ALE_LOG_EVENT_END
647: }
649: #undef __FUNCT__
651: void PreSieve::addCone(Obj<Point_set> coneSet, const Point& j) {
652: ALE_LOG_EVENT_BEGIN
653: CHKCOMM(*this);
654: this->__checkLock();
655: // Add j to the base in case coneSet is empty and no addArrow are executed
656: this->addBasePoint(j);
657: for(Point_set::iterator cone_itor = coneSet->begin(); cone_itor != coneSet->end(); cone_itor++) {
658: this->addArrow(*cone_itor, j);
659: }
660: ALE_LOG_EVENT_END
661: }
663: #undef __FUNCT__
665: Obj<Point_set> PreSieve::cone(const Point& point) {
666: static Obj<Point_set> c;
668: ALE_LOG_STAGE_BEGIN;
669: ALE_LOG_EVENT_BEGIN
670: if (this->_cone.find(point) != this->_cone.end()) {
671: c = this->_cone[point];
672: } else {
673: c = Point_set();
674: }
675: ALE_LOG_EVENT_END
676: ALE_LOG_STAGE_END;
677: return c;
678: }
680: #undef __FUNCT__
682: Obj<Point_set> PreSieve::nCone(Obj<Point_set> chain, int32_t n) {
683: Obj<Point_set> top(new Point_set);
684: ALE_LOG_STAGE_BEGIN;
685: CHKCOMM(*this);
686: // Compute the point set obtained by taking the cone recursively on a set of points in the base
687: // (i.e., the set of cap points resulting after each iteration is used again as the base for the next cone computation).
688: // Note: a 0-cone is the chain itself.
690: // We use two Point_set pointers and swap them at the beginning of each iteration
691: Obj<Point_set> bottom(new Point_set);
692: if(n == 0) {
693: top.copy(chain);
694: }
695: else {
696: top = chain;
697: }
698: // If no iterations are executed, chain is returned
699: for(int32_t i = 0; i < n; i++) {
700: // Swap pointers and clear top
701: Obj<Point_set> tmp = top; top = bottom; bottom = tmp;
702: // If top == chain, replace top->pointer() with another &Point_set; this avoids having to copy *chain.
703: if(top == chain) {
704: top = new Point_set;
705: }
706: else {
707: top->clear();
708: }
709: // Traverse the points in bottom
710: for(Point_set::iterator b_itor = bottom->begin(); b_itor != bottom->end(); b_itor++) {
711: Point p = *b_itor;
712: if (this->_cone.find(p) != this->_cone.end()) {
713: // Traverse the points in the cone over p
714: for(Point_set::iterator pCone_itor = this->_cone[p].begin(); pCone_itor != this->_cone[p].end(); pCone_itor++) {
715: Point q = *pCone_itor;
716: top->insert(q);
717: }
718: }
719: }
720: }
721: // IMPROVE: memory use can be imporoved (number of copies and alloc/dealloc reduced)
722: // if pointers to Point_set are used
723: ALE_LOG_STAGE_END;
724: return top;
725: }// PreSieve::nCone()
728: #undef __FUNCT__
730: Obj<Point_set> PreSieve::nClosure(Obj<Point_set> chain, int32_t n) {
731: Obj<Point_set> closure = Point_set();
732: ALE_LOG_STAGE_BEGIN;
733: CHKCOMM(*this);
734: // Compute the point set obtained by recursively accumulating the cone over all of points of a set in the base
735: // (i.e., the set of cap points resulting after each iteration is both stored in the resulting set and
736: // used again as the base of a cone computation).
737: // Note: a 0-closure is the chain itself.
739: // If no iterations are executed, chain is returned
740: closure.copy(chain); // copy the initial set
741: // We use two Point_set pointers and swap them at the beginning of each iteration
742: Obj<Point_set> top(chain);
743: Obj<Point_set> bottom = Point_set();
744: for(int32_t i = 0; i < n; i++) {
745: // Swap pointers and clear top
746: Obj<Point_set> tmp = top; top = bottom; bottom = tmp;
747: // If top == chain, replace top->pointer() with another &Point_set; this avoids having to copy *chain.
748: if(top == chain) {
749: top = new Point_set;
750: }
751: else {
752: top->clear();
753: }
754: // Traverse the points in bottom
755: for(Point_set::iterator b_itor = bottom->begin(); b_itor != bottom->end(); b_itor++) {
756: if (this->_cone.find(*b_itor) != this->_cone.end()) {
757: // Traverse the points in the cone over p
758: for(Point_set::iterator pCone_itor = this->_cone[*b_itor].begin(); pCone_itor != this->_cone[*b_itor].end(); pCone_itor++) {
759: top->insert(*pCone_itor);
760: closure->insert(*pCone_itor);
761: }
762: }
763: }
764: }
765: // IMPROVE: memory use can be imporoved (number of copies and alloc/dealloc reduced)
766: // if pointers to Point_sets are used
767: ALE_LOG_STAGE_END;
768: return closure;
769: }// PreSieve::nClosure()
772: #undef __FUNCT__
774: Obj<PreSieve> PreSieve::nClosurePreSieve(Obj<Point_set> chain, int32_t n, Obj<PreSieve> closure) {
775: ALE_LOG_STAGE_BEGIN;
776: CHKCOMM(*this);
777: if(closure.isNull()) {
778: closure = Obj<PreSieve>(new PreSieve(this->comm));
779: }
780: // Compute the PreSieve obtained by carving out the PreSieve 'over' a given chain up to height n.
781: // Note: a 0-closure is a PreSieve with the chain itself in the base, an empty cap and no arrows.
783: Obj<Point_set> base(new Point_set);
784: Obj<Point_set> top = chain;;
785: // Initially closure contains the chain only in the base, and top contains the chain itself
786: for(Point_set::iterator c_itor = chain->begin(); c_itor != chain->end(); c_itor++) {
787: Point p = *c_itor;
788: closure->addBasePoint(p);
789: }
790: // If no iterations are executed, chain is returned
791: for(int32_t i = 0; i < n; i++) {
792: // Swap base and top
793: Obj<Point_set> tmp = top; top = base; base = tmp;
794: // Traverse the points in the base
795: for(Point_set::iterator b_itor = base->begin(); b_itor != base->end(); b_itor++) {
796: Point p = *b_itor;
797: // Compute the cone over p and add it to closure
798: // IMPROVE: memory use can be improve if nCone returned a pointer to Point_set
799: Point_set pCone = this->nCone(p,1);
800: // Add the cone to top
801: top->join(pCone);
802: // Set the cone over p in the closure
803: closure->addCone(pCone,p);
804: }
805: }
806: ALE_LOG_STAGE_END;
807: return closure;
808: }// PreSieve::nClosurePreSieve()
810: #undef __FUNCT__
812: Obj<Point_set> PreSieve::support(const Point& point) {
813: static Obj<Point_set> s;
815: ALE_LOG_STAGE_BEGIN;
816: ALE_LOG_EVENT_BEGIN
817: if (this->_support.find(point) != this->_support.end()) {
818: s = this->_support[point];
819: } else {
820: s = Point_set();
821: }
822: ALE_LOG_EVENT_END
823: ALE_LOG_STAGE_END;
824: return s;
825: }
827: #undef __FUNCT__
829: Obj<Point_set> PreSieve::nSupport(const Obj<Point_set>& chain, int32_t n) {
830: Obj<Point_set> top;
831: Obj<Point_set> bottom;
833: ALE_LOG_STAGE_BEGIN;
834: CHKCOMM(*this);
835: // Compute the point set obtained by taking support recursively on a set of points in the cap
836: // (i.e., the set of base points resulting after each iteration is used again as the cap for the next support computation).
837: // Note: a 0-support is the chain itself.
839: // We use two Point_set pointers and swap them at the beginning of each iteration
840: top.create(Point_set());
841: bottom.create(Point_set());
842: if(n == 0) {
843: bottom.copy(chain);
844: } else {
845: bottom = chain;
846: }
847: // If no iterations are executed, chain is returned
848: for(int32_t i = 0; i < n; i++) {
849: // Swap pointers and clear bottom
850: Obj<Point_set> tmp = top; top = bottom; bottom = tmp;
851: // If bottom == chain, replace bottom->pointer() with another &Point_set; this avoids having to copy *chain.
852: if(bottom == chain) {
853: bottom.create(Point_set());
854: } else {
855: bottom->clear();
856: }
857: // Traverse the points in top
858: for(Point_set::iterator t_itor = top->begin(); t_itor != top->end(); t_itor++) {
859: Point q = *t_itor;
860: // Must check since map automatically inserts missing keys
861: if (this->_support.find(q) != this->_support.end()) {
862: // Traverse the points in the support under q
863: for(Point_set::iterator qSupport_itor = this->_support[q].begin(); qSupport_itor != this->_support[q].end(); qSupport_itor++) {
864: Point p = *qSupport_itor;
865: bottom->insert(p);
866: }
867: }
868: }
869: }
870: // IMPROVE: memory use can be improved (number of copies and alloc/dealloc reduced)
871: // if pointers to Point_set are used
872: ALE_LOG_STAGE_END;
873: return bottom;
874: }// PreSieve::nSupport()
876: #undef __FUNCT__
878: Obj<Point_set> PreSieve::nStar(const Obj<Point_set>& chain, int32_t n) {
879: Obj<Point_set> star(new Point_set);
880: ALE_LOG_STAGE_BEGIN;
881: CHKCOMM(*this);
882: // Compute the point set obtained by recursively accumulating the support of all of points of a set in the cap
883: // (i.e., the set of base points resulting after each iteration is both stored in the resulting set and
884: // used again as the cap of a support computation).
885: // Note: a 0-star is the chain itself.
887: // If no iterations are executed, chain is returned
888: star.copy(chain); // copy the initial set
889: // We use two Point_set pointers and swap them at the beginning of each iteration
890: Obj<Point_set> bottom(chain);
891: Obj<Point_set> top(new Point_set);
892: for(int32_t i = 0; i < n; i++) {
893: // Swap pointers and clear bottom
894: Obj<Point_set> tmp = top; top = bottom; bottom = tmp;
895: // If bottom == chain, replace bottom->pointer() with another &Point_set; this avoids having to copy *chain.
896: if(bottom == chain) {
897: bottom = new Point_set;
898: }
899: else {
900: bottom->clear();
901: }
902: // Traverse the points in top
903: for(Point_set::iterator t_itor = top->begin(); t_itor != top->end(); t_itor++) {
904: Point q = *t_itor;
905: // Traverse the points in the support of q
906: for(Point_set::iterator qSupport_itor = this->_support[q].begin(); qSupport_itor != this->_support[q].end(); qSupport_itor++) {
907: Point p = *qSupport_itor;
908: bottom->insert(p);
909: star->insert(p);
910: }
911: }
912: }
913: // IMPROVE: memory use can be imporoved (number of copies and alloc/dealloc reduced)
914: // if pointers to Point_sets are used
915: ALE_LOG_STAGE_END;
916: return star;
917: }// PreSieve::nStar()
920: #undef __FUNCT__
922: Point_set PreSieve::nStar(const Point_set& chain, int32_t n) {
923: return (Point_set) nStar(Obj<Point_set>(chain), n);
924: }// PreSieve::nStar()
927: #undef __FUNCT__
929: Obj<PreSieve> PreSieve::nStarPreSieve(Obj<Point_set> chain, int32_t n, Obj<PreSieve> star) {
930: ALE_LOG_STAGE_BEGIN;
931: CHKCOMM(*this);
932: // Compute the PreSieve obtained by accumulating intermediate kSupports (1<=k<=n) for a set of points in the cap.
933: // Note: a 0-star is the PreSieve containing the chain itself in the cap, an empty base and no arrows.
934: if(star.isNull()){
935: star = Obj<PreSieve>(new PreSieve(this->comm));
936: }
937: // We use a 'PreSieve-light' -- a map from points in the chain to the bottom-most supports computed so far.
938: Point__Point_set bottom;
939: // Initially star contains only the chain in the cap, and bottom has a singleton {q} set under each q.
940: for(Point_set::iterator c_itor = chain->begin(); c_itor != chain->end(); c_itor++) {
941: Point q = *c_itor;
942: Point_set qSet; qSet.insert(q);
943: bottom[q] = qSet;
944: Point_set emptySet;
945: star->addSupport(q,emptySet);
946: }
947: for(int32_t i = 0; i < n; i++) {
948: // Traverse the chain
949: for(Point_set::iterator c_itor = chain->begin(); c_itor != chain->end(); c_itor++) {
950: Point q = *c_itor;
951: // Compute the new bottom set as the 1-support of the current bottom set in 'this'
952: bottom[q] = this->nSupport(bottom[q],1);
953: // Add the new bottom set to the support of q in the support PreSieve
954: star->addSupport(q,bottom[q]);
955: }// for(Point_set::iterator c_itor = chain.begin(); c_itor != chain.end(); c_itor++)
956: }// for(int32_t i = 0; i < n; i++)
957: ALE_LOG_STAGE_END;
958: return star;
959: }// PreSieve::nStarPreSieve()
962: #undef __FUNCT__
964: Point_set PreSieve::nMeet(Point_set c0, Point_set c1, int32_t n) {
965: Point_set meet;
966: ALE_LOG_STAGE_BEGIN;
967: // The strategy is to compute the intersection of cones over the chains, remove the intersection
968: // and use the remaining two parts -- two disjoined components of the symmetric difference of cones -- as the new chains.
969: // The intersections at each stage are accumulated and their union is the meet.
970: // The iteration stops after n steps in addition to the meet of the initial chains or sooner if at least one of the chains is empty.
971: ALE::Obj<ALE::Point_set> cone;
972: // Check if any the initial chains may be empty, so that we don't perform spurious iterations
973: if((c0.size() == 0) || (c1.size() == 0)) {
974: //return meet;
975: }
976: else { // nonzero sizes
977: for(int32_t i = 0; i <= n; i++) {
978: Point_set *c = &c0;
979: Point_set *cc = &c1;
980: // Traverse the smallest cone set
981: if(cc->size() < c->size()) {
982: Point_set *tmp = c; c = cc; cc = tmp;
983: }
984: // Compute the intersection of c & cc and put it in meet at the same time removing it from c and cc
985: for(Point_set::iterator c_itor = c->begin(); c_itor != c->end(); c_itor++) {
986: if(cc->find(*c_itor)!= cc->end()) {
987: meet.insert(*c_itor);
988: cc->erase(*c_itor);
989: c->erase(c_itor);
990: }
991: }// for(Point_set::iterator c_itor = c->begin(); c_itor != c->end(); c_itor++)
992: // 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.
993: cone = this->cone(c0);
994: c0.insert(cone->begin(), cone->end());
995: if(c0.size() == 0) {
996: //return meet;
997: break;
998: }
999: cone = this->cone(c1);
1000: c1.insert(cone->begin(), cone->end());
1001: if(c1.size() == 0) {
1002: //return meet;
1003: break;
1004: }
1005: }// for(int32_t i = 0; i <= n; i++)
1006: }// nonzero sizes
1007: ALE_LOG_STAGE_END;
1008: return meet;
1009: }// PreSieve::nMeet()
1012: #undef __FUNCT__
1014: Point_set PreSieve::nMeetAll(Point_set& chain, int32_t n) {
1015: // The strategy is the same as in nMeet, except it is performed on an array of chains/cones--one per point in chain--simultaneously.
1016: // This may be faster than applying 'meet' recursively, since intersections may become empty faster.
1017: Point_set meets;
1018: // Populate the 'cones' map, while checking if any of the initial cones may be empty,
1019: // so that we don't perform spurious iterations. At the same time determine the cone of the smallest size.
1020: Point__Point_set cones;
1021: Point minp = *(chain.begin());
1022: for(Point_set::iterator chain_itor = chain.begin(); chain_itor != chain.end(); chain_itor++) {
1023: Point p = *chain_itor;
1024: cones[p] = this->cone(p);
1025: if(cones[p].size() == 0) {
1026: return meets;
1027: }
1028: if(cones[p].size() < cones[minp].size()) {
1029: minp = p;
1030: }
1031: }// for(Point_set::iterator chain_itor = chain.begin(); chain_itor != chain.end(); chain_itor++)
1032:
1033: for(int32_t i = 0; i <= n; i++) {
1034: Point_set *c = &cones[minp];
1035: // Traverse the smallest cone set pointed to by c and compute the intersection of c with the other cones,
1036: // putting it in meet at the same time removing it from c and the other cones.
1037: for(Point_set::iterator c_itor = c->begin(); c_itor != c->end(); c_itor++) {
1038: Point q = *c_itor;
1039: // Keep a flag indicating whether q belongs to the intersection of all cones.
1040: int qIsIn = 1;
1041: for(Point_set::iterator chain_itor = chain.begin(); chain_itor != chain.end(); chain_itor++) {
1042: Point p = *chain_itor;
1043: // skip p equal to minp
1044: if(p != minp) {
1045: if(cones[p].find(q) == cones[p].end()) {// q is absent from the cone over p
1046: qIsIn = 0;
1047: break;
1048: }
1049: }
1050: }// for(Point_set::iterator chain_itor = chain.begin(); chain_itor != chain.end(); chain_itor++)
1051: // if a is in, add it to the meet and remove from all the cones
1052: if(qIsIn) {
1053: meets.insert(q);
1054: for(Point_set::iterator chain_itor = chain.begin(); chain_itor != chain.end(); chain_itor++) {
1055: Point p = *chain_itor;
1056: cones[p].erase(q);
1057: }
1058: }
1059: // Recall that erase(q) should not invalidate c_itor
1060: }// for(Point_set::iterator c_itor = c->begin(); c_itor != c->end(); c_itor++)
1062: // Replace each of the cones with a cone over it, and check if either is empty; if so, we are done -- return meets.
1063: // At the same time, locate the point with the smallest cone over it.
1064: for(Point_set::iterator chain_itor = chain.begin(); chain_itor != chain.end(); chain_itor++) {
1065: Point p = *chain_itor;
1066: // We check before and after computing the cone
1067: if(cones[p].size() == 0) {
1068: return meets;
1069: }
1070: cones[p] = this->cone(cones[p]);
1071: if(cones[p].size() == 0) {
1072: return meets;
1073: }
1074: if(cones[p].size() < cones[minp].size()) {
1075: minp = p;
1076: }
1077: }// for(Point_set::iterator chain_itor = chain.begin(); chain_itor != chain.end(); chain_itor++)
1079: }
1080: return meets;
1081: }// PreSieve::nMeetAll()
1083: #undef __FUNCT__
1085: Point_set PreSieve::nJoin(Point_set c0, Point_set c1, int32_t n) {
1086: Point_set join;
1087: ALE_LOG_STAGE_BEGIN;
1088: // The strategy is to compute the intersection of the supports of the two chains, remove the intersection
1089: // and use the remaining two parts -- two disjoined components of the symmetric difference of the supports -- as the new chains.
1090: // The intersections at each stage are accumulated and their union is the join.
1091: // The iteration stops after n steps apart from the join of the initial chains or sooner if at least one of the chains is empty.
1092: ALE::Obj<ALE::Point_set> support;
1094: // Check if any the initial chains may be empty, so that we don't perform spurious iterations
1095: if((c0.size() == 0) || (c1.size() == 0)) {
1096: // return join;
1097: }
1098: else { // nonzero sizes
1099: for(int32_t i = 0; i <= n; i++) {
1100: Point_set *s = &c0;
1101: Point_set *ss = &c1;
1102: // Traverse the smallest supp set
1103: if(ss->size() < s->size()) {
1104: Point_set *tmp = s; s = ss; ss = tmp;
1105: }
1106: // Compute the intersection of s & ss and put it in join at the same time removing it from s and ss
1107: for(Point_set::iterator s_itor = s->begin(); s_itor != s->end(); s_itor++) {
1108: if(ss->find(*s_itor)!= ss->end()) {
1109: join.insert(*s_itor);
1110: ss->erase(*s_itor);
1111: s->erase(*s_itor);
1112: }
1113: }// for(Point_set::iterator s_itor = s->begin(); s_itor != s->end(); s_itor++)
1114: // Replace each of the chains with its support, and check if either is empty; if so, stop
1115: support = this->support(c0);
1116: c0.insert(support->begin(), support->end());
1117: if(c0.size() == 0) {
1118: // return join;
1119: break;
1120: }
1121: support = this->support(c1);
1122: c1.insert(support->begin(), support->end());
1123: if(c1.size() == 0) {
1124: // return join;
1125: break;
1126: }
1127: }// for(int32_t i = 0; i <= n; i++)
1128: }// nonzero sizes
1129: ALE_LOG_STAGE_END;
1130: return join;
1131: }// PreSieve::nJoin()
1134: #undef __FUNCT__
1136: Point_set PreSieve::nJoinAll(Point_set& chain, int32_t n) {
1137: // The strategy is the same as in nJoin, except it is performed on an array of chains/supps--one per point in chain--simultaneously.
1138: // This may be faster than applying 'join' recursively, since intersections may become empty faster.
1139: Point_set joins;
1141: // Populate the 'supps' map, while checking if any of the initial supports may be empty,
1142: // so that we don't perform spurious iterations. At the same time determine the supp of the smallest size.
1143: Point__Point_set supps;
1144: Point minp = *(chain.begin());
1145: for(Point_set::iterator chain_itor = chain.begin(); chain_itor != chain.end(); chain_itor++) {
1146: Point p = *chain_itor;
1147: supps[p] = this->support(p);
1148: if(supps[p].size() == 0) {
1149: return joins;
1150: }
1151: if(supps[p].size() < supps[minp].size()) {
1152: minp = p;
1153: }
1154: }// for(Point_set::iterator chain_itor = chain.begin(); chain_itor != chain.end(); chain_itor++)
1155:
1156: for(int32_t i = 0; i <= n; i++) {
1157: Point_set *s = &supps[minp];
1158: // Traverse the smallest supp set pointed to by s and compute the intersection of s with the other supps,
1159: // putting it in join at the same time removing it from s and the other supps.
1160: for(Point_set::iterator s_itor = s->begin(); s_itor != s->end(); s_itor++) {
1161: Point q = *s_itor;
1162: // Keep a flag indicating whether q belongs to the intersection of all supps.
1163: int qIsIn = 1;
1164: for(Point_set::iterator chain_itor = chain.begin(); chain_itor != chain.end(); chain_itor++) {
1165: Point p = *chain_itor;
1166: // skip p equal to minp
1167: if(p != minp) {
1168: if(supps[p].find(q) == supps[p].end()) {// q is absent from the support of p
1169: qIsIn = 0;
1170: break;
1171: }
1172: }
1173: }// for(Point_set::iterator chain_itor = chain.begin(); chain_itor != chain.end(); chain_itor++)
1174: // if a is in, add it to the join and remove from all the supps
1175: if(qIsIn) {
1176: joins.insert(q);
1177: for(Point_set::iterator chain_itor = chain.begin(); chain_itor != chain.end(); chain_itor++) {
1178: Point p = *chain_itor;
1179: supps[p].erase(q);
1180: }
1181: }
1182: // Recall that erase(q) should not invalidate s_itor
1183: }// for(Point_set::iterator s_itor = s->begin(); s_itor != s->end(); s_itor++)
1185: // Replace each of the supps with its support, and check if any of the new supps is empty; if so, stop.
1186: // At the same time, locate the point with the smallest supp.
1187: for(Point_set::iterator chain_itor = chain.begin(); chain_itor != chain.end(); chain_itor++) {
1188: Point p = *chain_itor;
1189: // We check before and after computing the support
1190: if(supps[p].size() == 0) {
1191: return joins;
1192: }
1193: supps[p] = this->support(supps[p]);
1194: if(supps[p].size() == 0) {
1195: return joins;
1196: }
1197: if(supps[p].size() < supps[minp].size()) {
1198: minp = p;
1199: }
1200: }// for(Point_set::iterator chain_itor = chain.begin(); chain_itor != chain.end(); chain_itor++)
1202: }
1203: return joins;
1204: }// PreSieve::nJoinAll()
1206: #undef __FUNCT__
1208: PreSieve *PreSieve::add(PreSieve& s) {
1209: ALE_LOG_STAGE_BEGIN;
1210: CHKCOMMS(*this,s);
1211: this->__checkLock();
1212: // Add the points and arrows of s to those of 'this'
1213: Point_set emptySet;
1214: for(Point__Point_set::iterator sBase_itor = s._cone.begin(); sBase_itor != s._cone.end(); sBase_itor++) {
1215: Point p = sBase_itor->first;
1216: Point_set pCone = s._cone[p];
1217: this->addCone(pCone, p);
1218: }
1219: // Make sure all of the cap of s is added (some cap points have no arrows).
1220: for(Point__Point_set::iterator cap_itor = s._support.begin(); cap_itor != s._support.end(); cap_itor++) {
1221: Point q = cap_itor->first;
1222: this->addCapPoint(q);
1223: }
1224: ALE_LOG_STAGE_END;
1225: return this;
1226: }// PreSieve::add()
1229: #undef __FUNCT__
1231: PreSieve& PreSieve::add(Obj<PreSieve> s) {
1232: ALE_LOG_STAGE_BEGIN;
1233: // IMPROVE: do something about CHKCOMMS so that I can use Obj<PreSieve> to do the checking.
1234: //CHKCOMMS(*this,s.object());
1235: this->__checkLock();
1236: // Add the points and arrows of s to those of 'this'
1237: Point_set emptySet;
1238: for(Point__Point_set::iterator sBase_itor = s->_cone.begin(); sBase_itor != s->_cone.end(); sBase_itor++) {
1239: Point p = sBase_itor->first;
1240: Point_set pCone = s->_cone[p];
1241: this->addCone(pCone, p);
1242: }
1243: // Make sure all of the cap of s is added (some cap points have no arrows).
1244: for(Point__Point_set::iterator cap_itor = s->_support.begin(); cap_itor != s->_support.end(); cap_itor++) {
1245: Point q = cap_itor->first;
1246: this->addCapPoint(q);
1247: }
1248: ALE_LOG_STAGE_END;
1249: return *this;
1250: }// PreSieve::add()
1255: void PreSieve::view(const char *name) {
1256: CHKCOMM(*this);
1257: int32_t rank = this->commRank;
1259: ostringstream txt1, txt2, txt3, txt4, hdr;
1260: hdr << "[" << rank << "]:Viewing ";
1261: if(name == NULL) {
1262: hdr << "a ";
1263: }
1264: if(!this->isLocked()) {
1265: hdr << "(locked) ";
1266: }
1267: if(name != NULL) {
1268: hdr << "presieve '" << name << "'\n";
1269: }
1270: // Print header
1271: PetscSynchronizedPrintf(this->comm, hdr.str().c_str()); CHKERROR(ierr, "Error in PetscPrintf");
1272: // Use a string stream to accumulate output that is then submitted to PetscSynchronizedPrintf
1273: Point_set points = this->space();
1274: txt1 << "[" << rank << "]: space of size " << points.size() << " : ";
1275: for(Point_set::iterator p_itor = points.begin(); p_itor != points.end(); p_itor++)
1276: {
1277: Point p = (*p_itor);
1278: if(p_itor != points.begin()) {
1279: txt1 << ", ";
1280: }
1281: txt1 << "(" << p.prefix << ", " << p.index << ")";
1282: }
1283: txt1 << "\n";
1284: PetscSynchronizedPrintf(this->comm, txt1.str().c_str());
1285: //
1286: points = this->cap();
1287: txt2 << "[" << rank << "]: cap of size " << points.size() << " : ";
1288: for(Point_set::iterator p_itor = points.begin(); p_itor != points.end(); p_itor++)
1289: {
1290: Point p = (*p_itor);
1291: if(p_itor != points.begin()) {
1292: txt2 << ", ";
1293: }
1294: txt2 << "(" << p.prefix << ", " << p.index << ")";
1295: }
1296: txt2 << "\n";
1297: PetscSynchronizedPrintf(this->comm, txt2.str().c_str());
1298: //
1299: points = this->base();
1300: txt3 << "[" << rank << "]: base of size " << points.size() << " : ";
1301: for(Point_set::iterator p_itor = points.begin(); p_itor != points.end(); p_itor++)
1302: {
1303: Point p = (*p_itor);
1304: if(p_itor != points.begin()) {
1305: txt3 << ", ";
1306: }
1307: txt3 << "(" << p.prefix << ", " << p.index << ")";
1308: }
1309: txt3 << "\n";
1310: PetscSynchronizedPrintf(this->comm, txt3.str().c_str());
1311: //
1312: for(Point__Point_set::iterator cone_itor = this->_cone.begin(); cone_itor != this->_cone.end(); cone_itor++)
1313: {
1314: Point p = (*cone_itor).first;
1315: txt4 << "[" << rank << "]: cone over (" << p.prefix << ", " << p.index << "): ";
1316: // Traverse the local cone over p
1317: for(Point_set::iterator pCone_itor = this->_cone[p].begin(); pCone_itor != this->_cone[p].end(); pCone_itor++) {
1318: Point q = *pCone_itor;
1319: if(pCone_itor != this->_cone[p].begin()) {
1320: txt4 << ", ";
1321: }
1322: txt4 << "(" << q.prefix << ", " << q.index << ")";
1323: }
1324: txt4 << "\n";
1325: }
1326: #if 0
1327: for(Point__Point_set::iterator support_itor = this->_support.begin(); support_itor != this->_support.end(); support_itor++)
1328: {
1329: Point p = (*support_itor).first;
1330: txt4 << "[" << rank << "]: support of (" << p.prefix << ", " << p.index << "): ";
1331: // Traverse the local support of p
1332: for(Point_set::iterator pSupport_itor = this->_support[p].begin(); pSupport_itor != this->_support[p].end(); pSupport_itor++) {
1333: Point q = *pSupport_itor;
1334: if(pSupport_itor != this->_support[p].begin()) {
1335: txt4 << ", ";
1336: }
1337: txt4 << "(" << q.prefix << ", " << q.index << ")";
1338: }
1339: txt4 << "\n";
1340: }
1341: #endif
1342: PetscSynchronizedPrintf(this->comm, txt4.str().c_str());
1343: CHKERROR(ierr, "Error in PetscSynchronizedPrintf");
1344: PetscSynchronizedFlush(this->comm);
1345: CHKERROR(ierr, "Error in PetscSynchronizedFlush");
1347: }// PreSieve::view()
1348:
1349:
1350: // -------------------------------------------------------------------------------------------------------------------------------- //
1351: #undef __FUNCT__
1353: Stack* PreSieve::coneCompletion(int completionType, int footprintType, PreSieve *c) {
1354: Stack *s;
1355: ALE_LOG_STAGE_BEGIN;
1356: s = this->__computeCompletion(completedSetCap, completionType, footprintType, c);
1357: ALE_LOG_STAGE_END;
1358: return s;
1359: }// PreSieve::coneCompletion()
1361: #undef __FUNCT__
1363: Stack* PreSieve::supportCompletion(int completionType, int footprintType, PreSieve *c) {
1364: Stack *s;
1365: ALE_LOG_STAGE_BEGIN;
1366: s = this->__computeCompletion(completedSetBase, completionType, footprintType, c);
1367: ALE_LOG_STAGE_END;
1368: return s;
1369: }// PreSieve::supportCompletion()
1373: Stack* PreSieve::baseFootprint(int completionType, int footprintType, PreSieve *f) {
1374: PreSieve ownership(this->comm);
1375: Point self(this->commRank, this->commRank);
1376: if(footprintType == footprintTypeNone) {
1377: return NULL;
1378: }
1379: else if(footprintType == footprintTypeCone) {
1380: // Create a temporary 'ownership' PreSieve with this->base as the base and the process ID in the cone of all this->base points
1381: Point_set base = this->base();
1382: ownership.addSupport(self, base);
1383: // Now compute the cone completion of ownership and place it in footprint.
1384: // No footprint computation is required during the completion, since the completion itself carries that information.
1385: return ownership.coneCompletion(completionType, footprintTypeNone, f);
1386: }
1387: else if(footprintType == footprintTypeSupport) {
1388: // Create a temporary 'ownership' PreSieve with this->base as the cap and the process ID in the support of all this->base points
1389: Point_set base = this->base();
1390: ownership.addCone(base, self);
1391: // Now compute the support completion of ownership and place it in footprint.
1392: // No footprint computation is required during the completion, since the completion itself carries that information.
1393: return ownership.coneCompletion(completionType, footprintTypeNone, f);
1394: }
1395: else {
1396: throw Exception("Unknown footprintType");
1397: }
1398:
1399: }// PreSieve::baseFootprint()
1403: Stack* PreSieve::capFootprint(int completionType, int footprintType, PreSieve *f) {
1404: PreSieve ownership(this->comm);
1405: Point self(this->commRank, this->commRank);
1406: if(footprintType == footprintTypeNone) {
1407: return NULL;
1408: }
1409: else if(footprintType == footprintTypeCone) {
1410: // Create a temporary 'ownership' PreSieve with this->cap as the base and the process ID in the cone of all this->cap points.
1411: Point_set cap = this->cap();
1412: ownership.addSupport(self, cap);
1413: // Now compute the cone completion of ownership and place it in footprint.
1414: // No footprint computation is required during the completion, since the completion itself carries that information.
1415: return ownership.coneCompletion(completionType, footprintTypeNone, f);
1416: }
1417: else if(footprintType == footprintTypeSupport) {
1418: // Create a temporary 'ownership' PreSieve with this->cap as the cap and the process ID in the support of all this->cap points.
1419: Point_set cap = this->cap();
1420: ownership.addCone(cap, self);
1421: // Now compute the support completion of ownership and place it in footprint.
1422: // No footprint computation is required during the completion, since the completion itself carries that information.
1423: return ownership.coneCompletion(completionType, footprintTypeNone, f);
1424: }
1425: else {
1426: throw Exception("Unknown footprintType");
1427: }
1428:
1429: }// PreSieve::capFootprint()
1433: Stack* PreSieve::spaceFootprint(int completionType, int footprintType, PreSieve *f) {
1434: PreSieve ownership(this->comm);
1435: Point self(this->commRank, this->commRank);
1436: if(footprintType == footprintTypeNone) {
1437: return NULL;
1438: }
1439: else if(footprintType == footprintTypeCone) {
1440: // Create a temporary 'ownership' PreSieve with this->space as the base and the process ID in the cone of all this->space points.
1441: Point_set space = this->space();
1442: ownership.addSupport(self, space);
1443: // Now compute the cone completion of ownership and place it in footprint.
1444: // No footprint computation is required during the completion, since the completion itself carries that information.
1445: return ownership.coneCompletion(completionType, footprintTypeNone, f);
1446: }
1447: else if(footprintType == footprintTypeSupport) {
1448: //Create a temporary 'ownership' PreSieve with this->space as the cap and the process ID in the support of all this->space points.
1449: Point_set space = this->space();
1450: ownership.addCone(space, self);
1451: // Now compute the support completion of ownership and place it in footprint.
1452: // No footprint computation is required during the completion, since the completion itself carries that information.
1453: return ownership.coneCompletion(completionType, footprintTypeNone, f);
1454: }
1455: else {
1456: throw Exception("Unknown footprintType");
1457: }
1458:
1459: }// PreSieve::spaceFootprint()
1461: #undef __FUNCT__
1463: void PreSieve::__determinePointOwners(Point_set& points, int32_t *LeaseData, std::map<Point, int32_t, Point::Cmp>& owner) {
1464: ALE_LOG_STAGE_BEGIN;
1465: CHKCOMM(*this);
1467: int32_t size = this->commSize;
1468: int32_t rank = this->commRank;
1470: // We need to partition global nodes among lessors, which we do by global prefix
1471: // First we determine the extent of global prefices and the bounds on the indices with each global prefix.
1472: int32_t minGlobalPrefix = 0;
1473: // Determine the local extent of global domains
1474: for(Point_set::iterator point_itor = points.begin(); point_itor != points.end(); point_itor++) {
1475: Point p = (*point_itor);
1476: if((p.prefix < 0) && (p.prefix < minGlobalPrefix)) {
1477: minGlobalPrefix = p.prefix;
1478: }
1479: }
1480: int32_t MinGlobalPrefix;
1481: MPI_Allreduce(&minGlobalPrefix, &MinGlobalPrefix, 1, MPIU_INT, MPI_MIN, this->comm); CHKERROR(ierr, "Error in MPI_Allreduce");
1483: int__int BaseLowerBound, BaseUpperBound; // global quantities computed from the local quantities below
1484: int__int BaseMaxSize; // the maximum size of the global base index space by global prefix
1485: int__int BaseSliceScale, BaseSliceSize, BaseSliceOffset;
1487: if(MinGlobalPrefix < 0) { // if we actually do have global base points
1488: // Determine the upper and lower bounds on the indices of base points with each global prefix.
1489: // We use maps to keep track of these quantities with different global prefices.
1490: int__int baseLowerBound, baseUpperBound; // local quantities
1491: // Initialize local bound maps with the upper below lower so we can later recognize omitted prefices.
1492: for(int d = -1; d >= MinGlobalPrefix; d--) {
1493: baseLowerBound[d] = 0; baseUpperBound[d] = -1;
1494: }
1495: // Compute local bounds
1496: for(Point_set::iterator point_itor = points.begin(); point_itor != points.end(); point_itor++) {
1497: Point p = (*point_itor);
1498: int32_t d = p.prefix;
1499: int32_t i = p.index;
1500: if(d < 0) { // it is indeed a global prefix
1501: if (i < baseLowerBound[d]) {
1502: baseLowerBound[d] = i;
1503: }
1504: if (i > baseUpperBound[d]) {
1505: baseUpperBound[d] = i;
1506: }
1507: }
1508: }
1509: // Compute global bounds
1510: for(int32_t d = -1; d >= MinGlobalPrefix; d--){
1511: int32_t lowerBound, upperBound, maxSize;
1512: MPI_Allreduce(&baseLowerBound[d],&lowerBound,1,MPIU_INT,MPI_MIN,this->comm);
1513: CHKERROR(ierr, "Error in MPI_Allreduce");
1514: MPI_Allreduce(&baseUpperBound[d],&upperBound,1,MPIU_INT,MPI_MAX,this->comm);
1515: CHKERROR(ierr, "Error in MPI_Allreduce");
1516: maxSize = upperBound - lowerBound + 1;
1517: if(maxSize > 0) { // there are actually some indices in this global prefix
1518: BaseLowerBound[d] = lowerBound;
1519: BaseUpperBound[d] = upperBound;
1520: BaseMaxSize[d] = maxSize;
1522: // Each processor (at least potentially) owns a slice of the base indices with each global indices.
1523: // The size of the slice with global prefix d is BaseMaxSize[d]/size + 1 (except if rank == size-1,
1524: // where the slice size can be smaller; +1 is for safety).
1525:
1526: // For a non-empty domain d we compute and store the slice size in BaseSliceScale[d] (the 'typical' slice size) and
1527: // BaseSliceSize[d] (the 'actual' slice size, which only differs from 'typical' for processor with rank == size -1 ).
1528: // Likewise, each processor has to keep track of the index offset for each slice it owns and stores it in BaseSliceOffset[d].
1529: BaseSliceScale[d] = BaseMaxSize[d]/size + 1;
1530: BaseSliceSize[d] = BaseSliceScale[d];
1531: if (rank == size-1) {
1532: BaseSliceSize[d] = BaseMaxSize[d] - BaseSliceScale[d]*(size-1);
1533: }
1534: BaseSliceSize[d] = PetscMax(1,BaseSliceSize[d]);
1535: BaseSliceOffset[d] = BaseLowerBound[d] + BaseSliceScale[d]*rank;
1536: }// for(int32_t d = -1; d >= MinGlobalPrefix; d--){
1537: }//
1538: }// if(MinGlobalDomain < 0)
1540: for (Point_set::iterator point_itor = points.begin(); point_itor != points.end(); point_itor++) {
1541: Point p = (*point_itor);
1542: // Determine which slice p falls into
1543: // ASSUMPTION on Point type
1544: int32_t d = p.prefix;
1545: int32_t i = p.index;
1546: int32_t proc;
1547: if(d < 0) { // global domain -- determine the owner by which slice p falls into
1548: proc = (i-BaseLowerBound[d])/BaseSliceScale[d];
1549: }
1550: else { // local domain -- must refer to a rank within the comm
1551: if(d >= size) {
1552: throw ALE::Exception("Local domain outside of comm size");
1553: }
1554: proc = d;
1555: }
1556: owner[p] = proc;
1557: LeaseData[2*proc+1] = 1; // processor owns at least one of ours (i.e., the number of leases from proc is 1)
1558: LeaseData[2*proc]++; // count of how many we lease from proc
1559: }
1560: ALE_LOG_STAGE_END;
1561: }// PreSieve::__determinePointOwners()
1563: #undef __FUNCT__
1565: Stack *PreSieve::__computeCompletion(int completedSet, int completionType, int footprintType, PreSieve *completion) {
1566: // Overall completion stack C
1567: Stack *C = new Stack(this->comm);
1569: // First, we set up the structure of the stack in which we return the completion information
1570: // This structure is controled by the flags in the calling sequence (completedSet, completionType, footprintType)
1571: // and is best described by a picture (see the Paper).
1573:
1574: // "Structure" presieve -- either the point presieve or the arrow stack (both defined below).
1575: Obj<PreSieve> S;
1577: // Point presieve P
1578: Obj<PreSieve> P(completion);
1579: if(completion != NULL) {
1580: // The completion sieve must have the same communicator
1581: CHKCOMMS(*this, *completion);
1582: }
1583: else {
1584: P = Obj<PreSieve>(new PreSieve(this->getComm()));
1585: }
1587: // Arrow stacks/presieve
1588: Obj<PreSieve> A;
1589: Obj<Stack> AAA, A0, A1;
1590: switch(completionType) {
1591: case completionTypePoint:
1592: S = P;
1593: break;
1594: case completionTypeArrow:
1595: A = new PreSieve(this->getComm());
1596: A0 = new Stack(this->getComm());
1597: A0->setTop(P); A0->setBottom(A);
1598: A1 = new Stack(this->getComm());
1599: A1->setTop(A); A1->setBottom(P);
1600: AAA = new Stack(this->getComm());
1601: AAA->setTop(A0); AAA->setBottom(A1);
1602: S = AAA;
1603: break;
1604: default:
1605: throw(ALE::Exception("Unknown completionType"));
1606: }
1607:
1608: // Footprint presieve
1609: Obj<PreSieve> F;
1610: // Depending on the footprintType flag, the structure presieve S sits either at the bottom or the top of C,
1611: // while F is null or not.
1612: switch(footprintType){
1613: case footprintTypeCone:
1614: F = Obj<PreSieve>(new PreSieve(this->getComm()));
1615: // S is the top of C
1616: C->setTop(S);
1617: // Footprint
1618: C->setBottom(F);
1619: break;
1620: case footprintTypeSupport:
1621: F = Obj<PreSieve>(new PreSieve(this->getComm()));
1622: // S is the bottom of C
1623: C->setBottom(S);
1624: // Footprint
1625: C->setTop(F);
1626: break;
1627: case footprintTypeNone:
1628: C->setTop(S);
1629: C->setBottom(S);
1630: break;
1631: default:
1632: throw(ALE::Exception("Unknown footprint type"));
1633: }
1634:
1635: // If this is a local sieve, there is no completion to be computed
1636: // Otherwise we have to compute the completion of the cone over each i from the cones over i
1637: // held by all the processes in the communicator, and at the same time contributing to whatever cones
1638: // those processors may be computing.
1639:
1640: if( (this->comm == MPI_COMM_SELF) || (this->commSize == 1) ) {
1641: return C;
1642: }
1644: int32_t size = this->commSize;
1645: int32_t rank = this->commRank;
1646:
1647: bool debug = this->verbosity > 10;
1648: bool debug2 = this->verbosity > 11;
1649:
1650: /* Glossary:
1651: -- Node -- a Point either in the base or in the cone; terms 'Node' and 'Point' (capitalized or not)
1652: will be used interchangeably (sp?).
1653:
1654: -- Owned node -- a that have been assigned to this processor just for this routine.
1655: It is either a node with this processor's rank as the prefix, or a global node that has been
1656: assigned to this processor during the 'preprocessing' stage.
1657:
1658: -- Domain -- the set of all points in the base with a given prefix
1660: -- Global domain/node -- a domain corresponding to a global (i.e., negative) prefix;
1661: a global node is a node from a global domain, equivalently, a node with global prefix.
1663: -- Rented node -- an owned node assigned to this processor that belongs to the base of another processor.
1665: -- Leased node -- a node belonging to this processor's base that is owned by another processor.
1667: -- Lessor -- a processor owning a node rented by this processor (hence leasing it to this processor).
1668:
1669: -- Renter -- a processor that has leased a node owned by this processor.
1670:
1671: -- Rental -- a renter--rented-node combination
1673: -- Sharer -- (in relation to a given rented node) a renter sharing the node with another renter
1674:
1675: -- Neighbor -- (in relation to a given leased node) a processor leasing the same node -- analogous to a sharer
1676: from a renter's point of view
1678: -- Total -- always refers to the cumulative number of entities across all processors
1679: (e.g., total number of nodes -- the sum of all nodes owned by each processor).
1680: */
1681:
1682: // ASSUMPTION: here we assume that Point == {int32_t prefix, index}; should the definition of Point change, so must this code.
1684: // Determine the owners of base nodes and collect the lease data for each processor:
1685: // the number of nodes leased and the number of leases (0 or 1).
1686: int32_t *LeaseData; // 2 ints per processor: number of leased nodes and number of leases (0 or 1).
1687: PetscMalloc(2*size*sizeof(PetscInt),&LeaseData);CHKERROR(ierr, "Error in PetscMalloc");
1688: PetscMemzero(LeaseData,2*size*sizeof(PetscInt));CHKERROR(ierr, "Error in PetscMemzero");
1689:
1690: // We also need to keep the set of points (base or cap) we use to determine the neighbors,
1691: // as well as a list of "cones" that are either cones or supports, depending on what the points are.
1692: Point_set *points;
1693: Point__Point_set *cones;
1694: if(completedSet == completedSetCap) {
1695: points = new Point_set();
1696: *points = this->base();
1697: cones = &this->_cone;
1698: }
1699: else if(completedSet == completedSetBase){
1700: points = new Point_set();
1701: *points = this->cap();
1702: cones = &this->_support;
1703: }
1704: else {
1705: throw ALE::Exception("Unknown completedSet");
1706: }
1707: // determine owners of each base node and save it in a map
1708: std::map<Point, int32_t, Point::Cmp> owner;
1709: this->__determinePointOwners(*points, LeaseData, owner);
1710:
1711: // Now we accumulate the max lease size and the total number of renters
1712: int32_t MaxLeaseSize, RenterCount;
1713: PetscMaxSum(this->comm,LeaseData,&MaxLeaseSize,&RenterCount);
1714: CHKERROR(ierr,"Error in PetscMaxSum");
1715: PetscInfo1(0,"PreSieve::__computeCompletion: Number of renters %d\n", RenterCount);
1716: CHKERROR(ierr,"Error in PetscInfo");
1718: if(debug) { /* -------------------------------------------------------------- */
1719: PetscSynchronizedPrintf(this->comm, "[%d]: __computeCompletion: RenterCount = %d, MaxLeaseSize = %d\n", rank, RenterCount, MaxLeaseSize);
1720: CHKERROR(ierr, "Error in PetscSynchronizedPrintf");
1721: PetscSynchronizedFlush(this->comm);
1722: CHKERROR(ierr, "Error in PetscSynchronizedFlush");
1723: } /* ----------------------------------------------------------------------- */
1724:
1725: // post receives for all Rented nodes; we will be receiving 3 data items per rented node,
1726: // and at most MaxLeaseSize of nodes per renter
1727: PetscMPIInt tag1;
1728: PetscObjectGetNewTag(this->petscObj, &tag1); CHKERROR(ierr, "Failded on PetscObjectGetNewTag");
1729: int32_t *RentedNodes;
1730: MPI_Request *Renter_waits;
1731: if(RenterCount){
1732: PetscMalloc((RenterCount)*(3*MaxLeaseSize+1)*sizeof(int32_t),&RentedNodes); CHKERROR(ierr,"Error in PetscMalloc");
1733: PetscMemzero(RentedNodes,(RenterCount)*(3*MaxLeaseSize+1)*sizeof(int32_t)); CHKERROR(ierr,"Error in PetscMemzero");
1734: PetscMalloc((RenterCount)*sizeof(MPI_Request),&Renter_waits); CHKERROR(ierr,"Error in PetscMalloc");
1735: }
1736: for (int32_t i=0; i<RenterCount; i++) {
1737: MPI_Irecv(RentedNodes+3*MaxLeaseSize*i,3*MaxLeaseSize,MPIU_INT,MPI_ANY_SOURCE,tag1,this->comm,Renter_waits+i);
1738: CHKERROR(ierr,"Error in MPI_Irecv");
1739: }
1740:
1741: int32_t LessorCount;
1742: LessorCount = 0; for (int32_t i=0; i<size; i++) LessorCount += LeaseData[2*i+1];
1743: PetscInfo1(0,"PreSieve::__computeCompletion: Number of lessors %d\n",LessorCount);
1744: CHKERROR(ierr,"Error in PetscInfo");
1745: if(debug) { /* -------------------------------------------------------------- */
1746: PetscSynchronizedPrintf(this->comm, "[%d]: __computeCompletion: LessorCount = %d\n", rank, LessorCount);
1747: CHKERROR(ierr, "Error in PetscSynchronizedPrintf");
1748: PetscSynchronizedFlush(this->comm);
1749: CHKERROR(ierr, "Error in PetscSynchronizedFlush");
1750: } /* ----------------------------------------------------------------------- */
1751:
1752: // We keep only the data about the real lessors -- those that own the nodes we lease
1753: int32_t *LeaseSizes, *Lessors;
1754: if(LessorCount) {
1755: PetscMalloc(sizeof(int32_t)*(LessorCount), &LeaseSizes); CHKERROR(ierr, "Error in PetscMalloc");
1756: PetscMalloc(sizeof(int32_t)*(LessorCount), &Lessors); CHKERROR(ierr, "Error in PetscMalloc");
1757: }
1758: // We also need to compute the inverse to the Lessors array, since we need to be able to convert i into cntr
1759: // after using the owner array. We use a map LessorIndex; it is likely to be small -- ASSUMPTION
1760: int__int LessorIndex;
1761: // Traverse all processes in ascending order
1762: int32_t cntr = 0; // keep track of entered records
1763: for(int32_t i = 0; i < size; i++) {
1764: if(LeaseData[2*i]) { // if there are nodes leased from process i, record it
1765: LeaseSizes[cntr] = LeaseData[2*i];
1766: Lessors[cntr] = i;
1767: LessorIndex[i] = cntr;
1768: cntr++;
1769: }
1770: }
1771: PetscFree(LeaseData); CHKERROR(ierr, "Error in PetscFree");
1772: if(debug2) { /* ----------------------------------- */
1773: ostringstream txt;
1774: txt << "[" << rank << "]: __computeCompletion: lessor data [index, rank, lease size]: ";
1775: for(int32_t i = 0; i < LessorCount; i++) {
1776: txt << "[" << i << ", " << Lessors[i] << ", " << LeaseSizes[i] << "] ";
1777: }
1778: txt << "\n";
1779: PetscSynchronizedPrintf(this->comm, txt.str().c_str()); CHKERROR(ierr, "Error in PetscSynchronizedPrintf");
1780: PetscSynchronizedFlush(this->comm);CHKERROR(ierr,"Error in PetscSynchronizedFlush");
1781: }/* ----------------------------------- */
1782: if(debug2) { /* ----------------------------------- */
1783: ostringstream txt;
1784: txt << "[" << rank << "]: __computeCompletion: LessorIndex: ";
1785: for(int__int::iterator li_itor = LessorIndex.begin(); li_itor!= LessorIndex.end(); li_itor++) {
1786: int32_t i = (*li_itor).first;
1787: int32_t j = (*li_itor).second;
1788: txt << i << "-->" << j << "; ";
1789: }
1790: txt << "\n";
1791: PetscSynchronizedPrintf(this->comm, txt.str().c_str()); CHKERROR(ierr, "Error in PetscSynchronizedPrintf");
1792: PetscSynchronizedFlush(this->comm);CHKERROR(ierr,"Error in PetscSynchronizedFlush");
1793: }/* ----------------------------------- */
1795:
1796: // pack messages containing lists of leased base nodes and their cone sizes to the lessors
1797: int32_t LeasedNodeCount = points->size(); // all points are considered leased from someone
1798: int32_t *LeasedNodes;
1799: int32_t *LessorOffsets;
1800: // We need 3 ints per leased node -- 2 per Point and 1 for the cone size
1801: if(LeasedNodeCount) {
1802: PetscMalloc((3*LeasedNodeCount)*sizeof(PetscInt),&LeasedNodes); CHKERROR(ierr,"Error in PetscMalloc");
1803: }
1804: if(LessorCount) {
1805: PetscMalloc((LessorCount)*sizeof(PetscInt),&LessorOffsets); CHKERROR(ierr,"Error in PetscMalloc");
1806: LessorOffsets[0] = 0;
1807: }
1808: for (int32_t i=1; i<LessorCount; i++) { LessorOffsets[i] = LessorOffsets[i-1] + 3*LeaseSizes[i-1];}
1809: for (Point_set::iterator point_itor = points->begin(); point_itor != points->end(); point_itor++) {
1810: Point p = (*point_itor);
1811: int32_t ow = owner[p];
1812: int32_t ind = LessorIndex[ow];
1813: LeasedNodes[LessorOffsets[ind]++] = p.prefix;
1814: LeasedNodes[LessorOffsets[ind]++] = p.index;
1815: LeasedNodes[LessorOffsets[ind]++] = (*cones)[p].size();
1816: }
1817: if(LessorCount) {
1818: LessorOffsets[0] = 0;
1819: }
1820: for (int32_t i=1; i<LessorCount; i++) { LessorOffsets[i] = LessorOffsets[i-1] + 3*LeaseSizes[i-1];}
1821:
1822: // send the messages to the lessors
1823: MPI_Request *Lessor_waits;
1824: if(LessorCount) {
1825: PetscMalloc((LessorCount)*sizeof(MPI_Request),&Lessor_waits);CHKERROR(ierr,"Error in PetscMalloc");
1826: }
1827: for (int32_t i=0; i<LessorCount; i++) {
1828: MPI_Isend(LeasedNodes+LessorOffsets[i],3*LeaseSizes[i],MPIU_INT,Lessors[i],tag1,this->comm,&Lessor_waits[i]);
1829: CHKERROR(ierr,"Error in MPI_Isend");
1830: }
1831:
1832: // wait on receive request and prepare to record the identities of the renters responding to the request and their lease sizes
1833: std::map<int32_t, int32_t> Renters, RenterLeaseSizes;
1834: // Prepare to compute the set of renters of each owned node along with the cone sizes held by those renters over the node.
1835: // Since we don't have a unique ordering on the owned nodes a priori, we will utilize a map.
1836: std::map<Point, int_pair_set, Point::Cmp> NodeRenters;
1837: cntr = RenterCount;
1838: while (cntr) {
1839: int32_t arrivalNumber;
1840: MPI_Status Renter_status;
1841: MPI_Waitany(RenterCount,Renter_waits,&arrivalNumber,&Renter_status);
1842: CHKMPIERROR(ierr,ERRORMSG("Error in MPI_Waitany"));
1843: int32_t renter = Renter_status.MPI_SOURCE;
1844: Renters[arrivalNumber] = renter;
1845: MPI_Get_count(&Renter_status,MPIU_INT,&RenterLeaseSizes[arrivalNumber]); CHKERROR(ierr,"Error in MPI_Get_count");
1846: // Since there are 3 ints per leased node, the lease size is computed by dividing the received count by 3;
1847: RenterLeaseSizes[arrivalNumber] = RenterLeaseSizes[arrivalNumber]/3;
1848: // Record the renters for each node
1849: for (int32_t i=0; i<RenterLeaseSizes[arrivalNumber]; i++) {
1850: // Compute the offset into the RentedNodes array for the arrived lease.
1851: int32_t LeaseOffset = arrivalNumber*3*MaxLeaseSize;
1852: // ASSUMPTION on Point type
1853: Point node = Point(RentedNodes[LeaseOffset + 3*i], RentedNodes[LeaseOffset + 3*i+1]);
1854: int32_t coneSize = RentedNodes[LeaseOffset + 3*i + 2];
1855: NodeRenters[node].insert(int_pair(renter,coneSize));
1856: }
1857: cntr--;
1858: }
1859:
1860: if (debug) { /* ----------------------------------- */
1861: // We need to collect all the data to be submitted to PetscSynchronizedPrintf
1862: // We use a C++ string streams for that
1863: ostringstream txt;
1864: for (std::map<Point, int_pair_set, Point::Cmp>::iterator nodeRenters_itor=NodeRenters.begin();nodeRenters_itor!= NodeRenters.end();nodeRenters_itor++) {
1865: Point node = (*nodeRenters_itor).first;
1866: int_pair_set renterSet = (*nodeRenters_itor).second;
1867: // ASSUMPTION on point type
1868: txt << "[" << rank << "]: __computeCompletion: node (" << node.prefix << "," << node.index << ") is rented by " << renterSet.size() << " renters (renter, cone size): ";
1869: for (int_pair_set::iterator renterSet_itor = renterSet.begin(); renterSet_itor != renterSet.end(); renterSet_itor++)
1870: {
1871: txt << "(" << (*renterSet_itor).first << "," << (*renterSet_itor).second << ") ";
1872: }
1873: txt << "\n";
1874: }
1875: // Now send the C-string behind txt to PetscSynchronizedPrintf
1876: PetscSynchronizedPrintf(this->comm, txt.str().c_str()); CHKERROR(ierr, "Error in PetscSynchronizedPrintf");
1877: PetscSynchronizedFlush(this->comm);CHKERROR(ierr,"Error in PetscSynchronizedFlush");
1878: }/* ----------------------------------- */
1879:
1880: // wait on the original sends to the lessors
1881: MPI_Status *Lessor_status;
1882: if (LessorCount) {
1883: PetscMalloc((LessorCount)*sizeof(MPI_Status),&Lessor_status); CHKERROR(ierr,"Error in PetscMalloc");
1884: MPI_Waitall(LessorCount,Lessor_waits,Lessor_status); CHKERROR(ierr,"Error in MPI_Waitall");
1885: }
1886:
1888: // Neighbor counts: here the renters receive from the lessors the number of other renters sharing each leased node.
1889: // Prepare to receive three integers per leased node: two for the node itself and one for the number of neighbors over that node.
1890: // The buffer has the same structure as LeasedNodes, hence LessorOffsets can be reused.
1891: // IMPROVE: can probably reduce the message size by a factor of 3 if we assume an ordering on the nodes received from each lessor.
1892: // ASSUMPTION on Point type
1893: int32_t *NeighborCounts;
1894: if(LeasedNodeCount) {
1895: PetscMalloc(3*(LeasedNodeCount)*sizeof(PetscInt),&NeighborCounts); CHKERROR(ierr,"Error in PetscMalloc");
1896: }
1897: // Post receives for NeighbornCounts
1898: PetscMPIInt tag2;
1899: PetscObjectGetNewTag(this->petscObj, &tag2); CHKERROR(ierr, "Failded on PetscObjectGetNewTag");
1900: for (int32_t i=0; i<LessorCount; i++) {
1901: MPI_Irecv(NeighborCounts+LessorOffsets[i],3*LeaseSizes[i],MPIU_INT,Lessors[i],tag2,this->comm,&Lessor_waits[i]);
1902: CHKERROR(ierr,"Error in MPI_Irecv");
1903: }
1904: // pack and send messages back to renters; we need to send 3 integers per rental (2 for Point, 1 for sharer count)
1905: // grouped by the renter
1906: // ASSUMPTION on Point type
1907: // first we compute the total number of rentals
1908: int32_t TotalRentalCount = 0;
1909: for(std::map<Point, int_pair_set, Point::Cmp>::iterator nodeRenters_itor=NodeRenters.begin();nodeRenters_itor!=NodeRenters.end();nodeRenters_itor++){
1910: TotalRentalCount += (*nodeRenters_itor).second.size();
1911: }
1912: if(debug2) {
1913: PetscSynchronizedPrintf(this->comm, "[%d]: TotalRentalCount %d\n", rank, TotalRentalCount);
1914: CHKERROR(ierr, "Error in PetscSynchronizedPrintf");
1915: PetscSynchronizedFlush(this->comm); CHKERROR(ierr, "PetscSynchronizedFlush");
1916: }/* ----------------------------------- */
1918: // Allocate sharer counts array for all rentals
1919: int32_t *SharerCounts;
1920: if(TotalRentalCount) {
1921: PetscMalloc(3*(TotalRentalCount)*sizeof(int32_t),&SharerCounts); CHKERROR(ierr,"Error in PetscMalloc");
1922: }
1923: // Renters are traversed in the order of their original arrival index by arrival number a
1924: int32_t RenterOffset = 0;
1925: cntr = 0;
1926: for(int32_t a = 0; a < RenterCount; a++) {
1927: // traverse the nodes leased by the renter
1928: int32_t RenterLeaseOffset = a*3*MaxLeaseSize;
1929: for(int32_t i = 0; i < RenterLeaseSizes[a]; i++) {
1930: // ASSUMPTION on Point type
1931: Point node;
1932: node.prefix = RentedNodes[RenterLeaseOffset + 3*i];
1933: node.index = RentedNodes[RenterLeaseOffset + 3*i + 1];
1934: SharerCounts[cntr++] = node.prefix;
1935: SharerCounts[cntr++] = node.index;
1936: // Decrement the sharer count by one not to count the current renter itself (with arrival number a).
1937: SharerCounts[cntr++] = NodeRenters[node].size()-1;
1938: }
1939: // Send message to renter
1940: MPI_Isend(SharerCounts+RenterOffset,3*RenterLeaseSizes[a],MPIU_INT,Renters[a],tag2,this->comm,Renter_waits+a);
1941: CHKERROR(ierr, "Error in MPI_Isend");
1942: // Offset is advanced by thrice the number of leased nodes, since we store 3 integers per leased node: Point and cone size
1943: RenterOffset = cntr;
1944: }
1945: // Wait on receives from lessors with the neighbor counts
1946: if (LessorCount) {
1947: MPI_Waitall(LessorCount,Lessor_waits,Lessor_status); CHKERROR(ierr,"Error in MPI_Waitall");
1948: }
1949: // Wait on the original sends to the renters
1950: MPI_Status *Renter_status;
1951: PetscMalloc((RenterCount)*sizeof(MPI_Status),&Renter_status);CHKERROR(ierr,"Error in PetscMalloc");
1952: if(RenterCount) {
1953: MPI_Waitall(RenterCount, Renter_waits, Renter_status);CHKERROR(ierr,"Error in MPI_Waitall");
1954: }
1955:
1956: if (debug) { /* ----------------------------------- */
1957: // Use a C++ string stream to report the numbers of shared nodes leased from each lessor
1958: ostringstream txt;
1959: cntr = 0;
1960: txt << "[" << rank << "]: __computeCompletion: neighbor counts by lessor-node [lessor rank, (node), neighbor count]: ";
1961: for(int32_t i = 0; i < LessorCount; i++) {
1962: // ASSUMPTION on point type
1963: for(int32_t j = 0; j < LeaseSizes[i]; j++)
1964: {
1965: int32_t prefix, index, sharerCount;
1966: prefix = NeighborCounts[cntr++];
1967: index = NeighborCounts[cntr++];
1968: sharerCount = NeighborCounts[cntr++];
1969: txt << "[" << Lessors[i] <<", (" << prefix << "," << index << "), " << sharerCount << "] ";
1970: }
1971: }
1972: txt << "\n";
1973: PetscSynchronizedPrintf(this->comm, txt.str().c_str()); CHKERROR(ierr, "Error in PetscSynchronizedPrintf");
1974: PetscSynchronizedFlush(this->comm); CHKERROR(ierr, "PetscSynchronizedFlush");
1975: }/* ----------------------------------- */
1978: // Now we allocate an array to receive the neighbor ranks and the remote cone sizes for each leased node,
1979: // hence, the total array size is 2*TotalNeighborCount.
1980: // Note that the lessor offsets must be recalculated, since they are no longer based on the number of nodes
1981: // leased from that lessor, but on the number of neighbor over the nodes leased from that lessor.
1983: // First we compute the numbers of neighbors over the nodes leased from a given lessor.
1984: int32_t TotalNeighborCount = 0;
1985: int32_t *NeighborCountsByLessor;
1986: if(LessorCount) {
1987: PetscMalloc((LessorCount)*sizeof(int32_t), &NeighborCountsByLessor); CHKERROR(ierr, "Error in PetscMalloc");
1988: }
1989: cntr = 0;
1990: for(int32_t i = 0; i < LessorCount; i++) {
1991: for(int32_t j = 0; j < LeaseSizes[i]; j++) {
1992: //ASSUMPTION on Point type affects NeighborCountsOffset size
1993: cntr += 2;
1994: TotalNeighborCount += NeighborCounts[cntr++];
1995: }
1996: if(i == 0) {
1997: NeighborCountsByLessor[i] = TotalNeighborCount;
1998: }
1999: else {
2000: NeighborCountsByLessor[i] = TotalNeighborCount - NeighborCountsByLessor[i-1];
2001: }
2002: }
2003: if (debug2) { /* ----------------------------------- */
2004: // Use a C++ string stream to report the numbers of shared nodes leased from each lessor
2005: ostringstream txt;
2006: cntr = 0;
2007: txt << "[" << rank << "]: __computeCompletion: NeighborCountsByLessor [rank, count]: ";
2008: for(int32_t i = 0; i < LessorCount; i++) {
2009: txt << "[" << Lessors[i] <<"," << NeighborCountsByLessor[i] << "]; ";
2010: }
2011: txt << "\n";
2012: PetscSynchronizedPrintf(this->comm, txt.str().c_str()); CHKERROR(ierr, "Error in PetscSynchronizedPrintf");
2013: PetscSynchronizedFlush(this->comm); CHKERROR(ierr, "PetscSynchronizedFlush");
2014: }/* ----------------------------------- */
2015: int32_t *Neighbors = 0;
2016: if(TotalNeighborCount) {
2017: PetscMalloc((2*TotalNeighborCount)*sizeof(int32_t),&Neighbors); CHKERROR(ierr,"Error in PetscMalloc");
2018: }
2019:
2020: // Post receives for Neighbors
2021: PetscMPIInt tag3;
2022: PetscObjectGetNewTag(this->petscObj, &tag3); CHKERROR(ierr, "Failded on PetscObjectGetNewTag");
2023: int32_t lessorOffset = 0;
2024: for(int32_t i=0; i<LessorCount; i++) {
2025: if(NeighborCountsByLessor[i]) { // We expect messages from lessors with a non-zero NeighborCountsByLessor entry only
2026: MPI_Irecv(Neighbors+lessorOffset,2*NeighborCountsByLessor[i],MPIU_INT,Lessors[i],tag3,this->comm,&Lessor_waits[i]);
2027: CHKERROR(ierr,"Error in MPI_Irecv");
2028: lessorOffset += 2*NeighborCountsByLessor[i];
2029: }
2030: }
2031: // Pack and send messages back to renters.
2032: // For each node p and each renter r (hence for each rental (p,r)) we must send to r a segment consisting of the list of all
2033: // (rr,cc) such that (p,rr) is a share and cc is the cone size over p at rr.
2034: // ALTERNATIVE, SCALABILITY:
2035: // 1. allocate an array capable of holding all messages to all renters and send one message per renter (more memory)
2036: // 2. allocate an array capable of holding all rentals for all nodes and send one message per share (more messages).
2037: // Here we choose 1 since we assume that the memory requirement is modest and the communication time is much more expensive,
2038: // however, this is likely to be application-dependent, and a switch should be introduced to change this behavior at will.
2039: // The rental segments are grouped by the renter recepient and within the renter by the node in the same order as SharerCounts.
2041: // We need to compute the send buffer size using the SharerCounts array.
2042: // Traverse the renters in order of their original arrival, indexed by the arrival number a, and then by the nodes leased by a.
2043: // Add up all entries equal to 2 mod 3 in SharerCounts (0 & 1 mod 3 are node IDs, ASSUMPTION on Point type) and double that number
2044: // to account for sharer ranks AND the cone sizes we are sending.
2045: int32_t SharersSize = 0; // 'Sharers' buffer size
2046: cntr = 0;
2047: for(int32_t a = 0; a < RenterCount; a++) {
2048: // traverse the number of nodes leased by the renter
2049: for(int32_t i = 0; i < RenterLeaseSizes[a]; i++) {
2050: SharersSize += SharerCounts[3*cntr+2];
2051: cntr++;
2052: }
2053: }
2054: SharersSize *= 2;
2055: // Allocate the Sharers array
2056: int32_t *Sharers;
2057: if(SharersSize) {
2058: PetscMalloc(SharersSize*sizeof(int32_t),&Sharers); CHKERROR(ierr,"Error in PetscMalloc");
2059: }
2060: // Now pack the messages and send them off.
2061: // Renters are traversed in the order of their original arrival index by arrival number a
2062: ostringstream txt; // DEBUG
2063: if(debug2) {
2064: txt << "[" << rank << "]: __computeCompletion: RenterCount = " << RenterCount << "\n";
2065: }
2066: RenterOffset = 0; // this is the current offset into Sharers needed for the send statement
2067: for(int32_t a = 0; a < RenterCount; a++) {//
2068: int32_t r = Renters[a];
2069: int32_t RenterLeaseOffset = a*3*MaxLeaseSize;
2070: int32_t SegmentSize = 0;
2071: // traverse the nodes leased by the renter
2072: for(int32_t i = 0; i < RenterLeaseSizes[a]; i++) {
2073: // Get a node p rented to r
2074: // ASSUMPTION on Point type
2075: Point p;
2076: p.prefix = RentedNodes[RenterLeaseOffset + 3*i];
2077: p.index = RentedNodes[RenterLeaseOffset + 3*i + 1];
2078: if(debug) {
2079: txt << "[" << rank << "]: __computeCompletion: renters sharing with " << r << " of node (" << p.prefix << "," << p.index << ") [rank, cone size]: ";
2080: }
2081: // now traverse the set of all the renters of p
2082: for(int_pair_set::iterator pRenters_itor=NodeRenters[p].begin(); pRenters_itor!=NodeRenters[p].end(); pRenters_itor++) {
2083: int32_t rr = (*pRenters_itor).first; // rank of a pRenter
2084: int32_t cc = (*pRenters_itor).second; // cone size over p at rr
2085: // skip r itself
2086: if(rr != r){
2087: Sharers[RenterOffset+SegmentSize++] = rr;
2088: Sharers[RenterOffset+SegmentSize++] = cc;
2089: if(debug) {
2090: txt << "[" << rr << "," << cc << "]; ";
2091: }
2092: }
2093: }// for(int_pair_set::iterator pRenters_itor=NodeRenters[p].begin(); pRenters_itor!=NodeRenters[p].end(); pRenters_itor++) {
2094: if(debug) {
2095: txt << "\n";
2096: }
2097: }// for(int32_t i = 0; i < RenterLeaseSizes[a]; i++) {
2098: // Send message to renter only if the segment size is positive
2099: if(SegmentSize > 0) {
2100: MPI_Isend(Sharers+RenterOffset,SegmentSize,MPIU_INT,Renters[a],tag3,this->comm,Renter_waits+a);
2101: CHKERROR(ierr, "Error in MPI_Isend");
2102: }
2103: // Offset is advanced by the segmentSize
2104: RenterOffset += SegmentSize;
2105: }// for(int32_t a = 0; a < RenterCount; a++) {
2106: if(debug) {
2107: PetscSynchronizedPrintf(this->comm, txt.str().c_str()); CHKERROR(ierr, "Error in PetscSynchronizedPrintf");
2108: PetscSynchronizedFlush(this->comm); CHKERROR(ierr, "PetscSynchronizedFlush");
2109: }
2110:
2111: // Wait on receives from lessors with the neighbor counts
2112: if (LessorCount) {
2113: MPI_Waitall(LessorCount,Lessor_waits,Lessor_status);CHKERROR(ierr,"Error in MPI_Waitall");
2114: }
2115: if (debug) { /* ----------------------------------- */
2116: // To report the neighbors at each lessor we use C++ a string stream
2117: ostringstream txt;
2118: int32_t cntr1 = 0;
2119: int32_t cntr2 = 0;
2120: for(int32_t i = 0; i < LessorCount; i++) {
2121: // ASSUMPTION on point type
2122: txt << "[" <<rank<< "]: __computeCompletion: neighbors over nodes leased from " <<Lessors[i]<< ":\n";
2123: int32_t activeLessor = 0;
2124: for(int32_t j = 0; j < LeaseSizes[i]; j++)
2125: {
2126: int32_t prefix, index, sharerCount;
2127: prefix = NeighborCounts[cntr1++];
2128: index = NeighborCounts[cntr1++];
2129: sharerCount = NeighborCounts[cntr1++];
2130: if(sharerCount > 0) {
2131: txt <<"[" << rank << "]:\t(" << prefix <<","<<index<<"): [rank, coneSize]: ";
2132: activeLessor++;
2133: }
2134: for(int32_t k = 0; k < sharerCount; k++) {
2135: int32_t sharer = Neighbors[cntr2++];
2136: int32_t coneSize = Neighbors[cntr2++];
2137: txt << "[" <<sharer <<", "<< coneSize << "] ";
2138: }
2139: }// for(int32_t j = 0; j < LeaseSizes[i]; j++)
2140: if(!activeLessor) {
2141: txt <<"[" << rank << "]:\tnone";
2142: }
2143: txt << "\n";
2144: }// for(int32_t i = 0; i < LessorCount; i++)
2145: PetscSynchronizedPrintf(this->comm,txt.str().c_str());CHKERROR(ierr,"Error in PetscSynchronizedPrintf");
2146: PetscSynchronizedFlush(this->comm);CHKERROR(ierr,"Error in PetscSynchronizedFlush");
2147: }/* ----------------------------------- */
2149: // This concludes the interaction of lessors and renters, and the exchange is completed by a peer-to-peer neighbor cone swap
2150: // (except we still have to wait on our last sends to the renters -- see below).
2151: // However, we don't free all of the arrays associated with the lessor-renter exchanges, since some of the data
2152: // still use those structures. Here are the arrays we can get rid of:
2153: if(RenterCount) {
2154: PetscFree(RentedNodes); CHKERROR(ierr, "Error in PetscFree");
2155: }
2156: if(SharersSize) {PetscFree(Sharers); CHKERROR(ierr, "Error in PetscFree");}
2157: if(LessorCount) {
2158: PetscFree(NeighborCountsByLessor); CHKERROR(ierr, "Error in PetscFree");
2159: PetscFree(Lessor_status); CHKERROR(ierr,"Error in PetscFree");
2160: PetscFree(Lessor_waits); CHKERROR(ierr,"Error in PetscFree");
2161: PetscFree(LessorOffsets); CHKERROR(ierr,"Error in PetscFree");
2162: PetscFree(LeaseSizes); CHKERROR(ierr,"Error in PetscFree");
2163: PetscFree(Lessors); CHKERROR(ierr,"Error in PetscFree");
2164: }
2165: if(LeasedNodeCount) {
2166: PetscFree(LeasedNodes); CHKERROR(ierr,"Error in PetscFree");
2167: }
2169: // First we determine the neighbors and the cones over each node to be received from or sent to each neigbor.
2170: std::map<int32_t, Point__int > NeighborPointConeSizeIn, NeighborPointConeSizeOut;
2171: // Traverse Neighbors and separate the data by neighbor over each point: NeighborPointConeSizeIn/Out.
2172: // cntr keeps track of the current position within the Neighbors array, node boundaries are delineated using NeighborCounts.
2173: // ASSUMPTION: 'Neighbors' stores node renter segments in the same order as NeighborCounts stores the node data.
2174: cntr = 0;
2175: for(int32_t i = 0; i < LeasedNodeCount; i++) {
2176: // ASSUMPTION on Point type
2177: Point p;
2178: p.prefix = NeighborCounts[3*i];
2179: p.index = NeighborCounts[3*i+1];
2180: int32_t pNeighborsCount = NeighborCounts[3*i+2]; // recall that NeighborCounts lists the number of neighbors after each node
2181: // extract the renters of p from Neighbors
2182: for(int32_t j = 0; j < pNeighborsCount; j++) {
2183: int32_t neighbor = Neighbors[cntr++];
2184: int32_t coneSize = Neighbors[cntr++];
2185: // Record the size of the cone over p coming in from neighbor
2186: NeighborPointConeSizeIn[neighbor][p]= coneSize;
2187: // Record the size of the cone over p going out to neighbor
2188: NeighborPointConeSizeOut[neighbor][p] = (*cones)[p].size();
2189: }
2190: }// for(int32_t i = 0; i < LeasedNodeCount; i++)
2193: // Compute total incoming cone sizes by neighbor and the total incomping cone size.
2194: // Also count the total number of neighbors we will be communicating with
2195: int32_t NeighborCount = 0;
2196: int__int NeighborConeSizeIn;
2197: int32_t ConeSizeIn = 0;
2198: ostringstream txt3;
2199: // Traverse all of the neighbors from whom we will be receiving cones.
2200: for(std::map<int32_t, Point__int>::iterator np_itor = NeighborPointConeSizeIn.begin();np_itor!=NeighborPointConeSizeIn.end(); np_itor++){
2201: int32_t neighbor = (*np_itor).first;
2202: NeighborConeSizeIn[neighbor] = 0;
2203: if(debug2) {
2204: txt3 << "[" << rank << "]: " << "__computeCompletion: NeighborPointConeSizeIn[" << neighbor << "]: ";
2205: }
2206: // Traverse all the points the cones over which we are receiving and add the cone sizes
2207: for(Point__int::iterator pConeSize_itor = (*np_itor).second.begin(); pConeSize_itor != (*np_itor).second.end(); pConeSize_itor++){
2208: NeighborConeSizeIn[neighbor] = NeighborConeSizeIn[neighbor] + (*pConeSize_itor).second;
2209: txt3 << "(" << (*pConeSize_itor).first.prefix << "," << (*pConeSize_itor).first.index << ")" << "->" << (*pConeSize_itor).second << "; ";
2210: }
2211: // Accumulate the total cone size
2212: ConeSizeIn += NeighborConeSizeIn[neighbor];
2213: NeighborCount++;
2214: txt3 << "NeighborConeSizeIn[" << neighbor << "]: " << NeighborConeSizeIn[neighbor] << "\n";
2215: }//for(std::map<int32_t, Point__int>::iterator np_itor=NeighborPointConeSizeIn.begin();np_itor!=NeighborPointConeSizeIn.end(); np_itor++)
2216: if(debug2) {
2217: PetscSynchronizedPrintf(this->comm,txt3.str().c_str());CHKERROR(ierr,"Error in PetscSynchronizedPrintf");
2218: PetscSynchronizedFlush(this->comm);CHKERROR(ierr,"Error in PetscSynchronizedFlush");
2219: }
2220: if(debug) {/* --------------------------------------------------------------------------------------------- */
2221: ostringstream txt;
2222: txt << "[" << rank << "]: __computeCompletion: total size of incoming cone: " << ConeSizeIn << "\n";
2223: for(int__int::iterator np_itor = NeighborConeSizeIn.begin();np_itor!=NeighborConeSizeIn.end();np_itor++)
2224: {
2225: int32_t neighbor = (*np_itor).first;
2226: int32_t coneSize = (*np_itor).second;
2227: txt << "[" << rank << "]: __computeCompletion: size of cone from " << neighbor << ": " << coneSize << "\n";
2229: }//int__int::iterator np_itor=NeighborConeSizeIn.begin();np_itor!=NeighborConeSizeIn.end();np_itor++)
2230: PetscSynchronizedPrintf(this->comm, txt.str().c_str());
2231: CHKERROR(ierr, "Error in PetscSynchronizedPrintf");
2232: PetscSynchronizedFlush(this->comm);
2233: CHKERROR(ierr, "Error in PetscSynchronizedFlush");
2234: }/* --------------------------------------------------------------------------------------------- */
2235: // Now we can allocate a receive buffer to receive all of the remote cones from neighbors
2236: int32_t *ConesIn;
2237: // ASSUMPTION on Point type affects ConesIn size -- double ConeSizeIn if Point uses two ints
2238: if(ConeSizeIn) {
2239: PetscMalloc(2*ConeSizeIn*sizeof(PetscInt),&ConesIn); CHKERROR(ierr,"Error in PetscMalloc");
2240: }
2241:
2242: // Wait on the original sends to the renters (the last vestige of the lessor-renter exchange epoch; we delayed it to afford the
2243: // greatest opportunity for a communication-computation overlap).
2244: if(RenterCount) {
2245: MPI_Waitall(RenterCount, Renter_waits, Renter_status); CHKERROR(ierr,"Error in MPI_Waitall");
2246: }
2247: if(RenterCount) {
2248: PetscFree(Renter_waits); CHKERROR(ierr, "Error in PetscFree");
2249: PetscFree(Renter_status); CHKERROR(ierr, "Error in PetscFree");
2250: }
2252: // Allocate receive requests
2253: MPI_Request *NeighborsIn_waits;
2254: if(NeighborCount) {
2255: PetscMalloc((NeighborCount)*sizeof(MPI_Request),&NeighborsIn_waits);CHKERROR(ierr,"Error in PetscMalloc");
2256: }
2257: // Post receives for ConesIn
2258: PetscMPIInt tag4;
2259: PetscObjectGetNewTag(this->petscObj, &tag4); CHKERROR(ierr, "Failded on PetscObjectGetNewTag");
2260: // Traverse all neighbors from whom we are receiving cones
2261: int32_t NeighborOffset = 0;
2262: int32_t n = 0;
2263: if(debug2) {
2264: PetscSynchronizedPrintf(this->comm, "[%d]: __computeCompletion: NeighborConeSizeIn.size() = %d\n",rank, NeighborConeSizeIn.size());
2265: CHKERROR(ierr, "Error in PetscSynchornizedPrintf");
2266: PetscSynchronizedFlush(this->comm);
2267: CHKERROR(ierr, "Error in PetscSynchornizedFlush");
2268: if(NeighborConeSizeIn.size()) {
2269: ierr=PetscSynchronizedPrintf(this->comm, "[%d]: __computeCompletion: *NeighborConeSizeIn.begin() = (%d,%d)\n",
2270: rank, (*NeighborConeSizeIn.begin()).first, (*NeighborConeSizeIn.begin()).second);
2271: CHKERROR(ierr, "Error in PetscSynchornizedPrintf");
2272: PetscSynchronizedFlush(this->comm);
2273: CHKERROR(ierr, "Error in PetscSynchornizedFlush");
2274:
2275: }
2276: }
2277: for(std::map<int32_t, int32_t>::iterator n_itor = NeighborConeSizeIn.begin(); n_itor!=NeighborConeSizeIn.end(); n_itor++) {
2278: int32_t neighbor = (*n_itor).first;
2279: int32_t coneSize = (*n_itor).second;
2280: // ASSUMPTION on Point type affects NeighborOffset and coneSize
2281: MPI_Irecv(ConesIn+NeighborOffset,2*coneSize,MPIU_INT,neighbor,tag4,this->comm, NeighborsIn_waits+n);
2282: CHKERROR(ierr, "Error in MPI_Irecv");
2283: // ASSUMPTION on Point type affects NeighborOffset
2284: NeighborOffset += 2*coneSize;
2285: n++;
2286: }
2288: // Compute the total outgoing cone sizes by neighbor and the total outgoing cone size.
2289: int__int NeighborConeSizeOut;
2290: int32_t ConeSizeOut = 0;
2291: // Traverse all the neighbors to whom we will be sending cones
2292: for(std::map<int32_t, Point__int>::iterator np_itor=NeighborPointConeSizeOut.begin();np_itor!=NeighborPointConeSizeOut.end();np_itor++){
2293: int32_t neighbor = (*np_itor).first;
2294: NeighborConeSizeOut[neighbor] = 0;
2295: // Traverse all the points cones over which we are sending and add up the cone sizes.
2296: for(Point__int::iterator pConeSize_itor = (*np_itor).second.begin(); pConeSize_itor != (*np_itor).second.end(); pConeSize_itor++){
2297: NeighborConeSizeOut[neighbor] = NeighborConeSizeOut[neighbor] + (*pConeSize_itor).second;
2298: }
2299: // Accumulate the total cone size
2300: ConeSizeOut += NeighborConeSizeOut[neighbor];
2301: }//for(std::map<int32_t,Point__int>::iterator np_itor=NeighborPointConeSizeOut.begin();np_itor!=NeighborPointConeSizeOut.end();np_itor++)
2304: // Now we can allocate a send buffer to send all of the remote cones to neighbors
2305: int32_t *ConesOut;
2306: // ASSUMPTION on Point type affects ConesOut size -- double ConeSizeOut if Point takes up 2 ints
2307: PetscMalloc(2*ConeSizeOut*sizeof(PetscInt),&ConesOut); CHKERROR(ierr,"Error in PetscMalloc");
2308: // Allocate send requests
2309: MPI_Request *NeighborsOut_waits;
2310: if(NeighborCount) {
2311: PetscMalloc((NeighborCount)*sizeof(MPI_Request),&NeighborsOut_waits);CHKERROR(ierr,"Error in PetscMalloc");
2312: }
2314: // Pack and send messages
2315: NeighborOffset = 0;
2316: cntr = 0;
2317: n = 0;
2318: ostringstream txt2;
2319: if(debug) {/* --------------------------------------------------------------------------------------------- */
2320: txt2 << "[" << rank << "]: __computeCompletion: total outgoing cone size: " << ConeSizeOut << "\n";
2321: }/* --------------------------------------------------------------------------------------------- */
2322: // Traverse all neighbors to whom we are sending cones
2323: for(std::map<int32_t, Point__int>::iterator np_itor=NeighborPointConeSizeOut.begin();np_itor!=NeighborPointConeSizeOut.end(); np_itor++){
2324: int32_t neighbor = (*np_itor).first;
2325: if(debug) { /* ------------------------------------------------------------ */
2326: txt2 << "[" << rank << "]: __computeCompletion: outgoing cones destined for " << neighbor << "\n";
2327: }/* ----------------------------------------------------------------------- */
2328: // ASSUMPTION: all Point__int maps are sorted the same way, so we safely can assume that
2329: // the receiver will be expecting points in the same order
2330: // Traverse all the points within this neighbor
2331: for(Point__int::iterator pConeSize_itor = (*np_itor).second.begin(); pConeSize_itor != (*np_itor).second.end(); pConeSize_itor++)
2332: {
2333: Point p = (*pConeSize_itor).first;
2334: if(debug) { /* ------------------------------------------------------------ */
2335: txt2 << "[" << rank << "]: \t cone over (" << p.prefix << ", " << p.index << "): ";
2336: }/* ----------------------------------------------------------------------- */
2337: // Traverse the local cone over p and store it in ConesOut
2338: for(Point_set::iterator cone_itor = (*cones)[p].begin(); cone_itor != (*cones)[p].end(); cone_itor++) {
2339: // ASSUMPTION on Point type affects storage of points in ConesOut
2340: Point q = *cone_itor;
2341: ConesOut[cntr++] = q.prefix;
2342: ConesOut[cntr++] = q.index;
2343: if(debug) { /* ------------------------------------------------------------ */
2344: txt2 << "(" << q.prefix << ", " << q.index << ") ";
2345: }/* ----------------------------------------------------------------------- */
2346: }
2347: if(debug) { /* ------------------------------------------------------------ */
2348: txt2 << "\n";
2349: }/* ----------------------------------------------------------------------- */
2350: }//for(Point__int::iterator pConeSize_itor = (*np_itor).second.begin(); pConeSize_itor!=(*np_itor).second.end();pConeSize_itor++)
2351: int32_t coneSize = NeighborConeSizeOut[neighbor];
2352: // ASSUMPTION on Point type affects the usage of coneSize and the calculation of Neighbor offset
2353: MPI_Isend(ConesOut+NeighborOffset,2*coneSize,MPIU_INT,neighbor,tag4,this->comm, NeighborsOut_waits+n);
2354: CHKERROR(ierr, "Error in MPI_Isend");
2355: // ASSUMPTION on Point type affects the computation of NeighborOffset -- double coneSize if Point uses up 2 ints
2356: NeighborOffset += 2*coneSize; // keep track of offset
2357: n++; // count neighbors
2358: }//for(std::map<int32_t,Point__int>::iterator np_itor=NeighborPointConeSizeIn.begin();np_itor!=NeighborPointConeSizeIn.end();np_itor++){
2359: if(debug) {/* --------------------------------------------------------------------------------------------- */
2360: PetscSynchronizedPrintf(this->comm, txt2.str().c_str());
2361: CHKERROR(ierr, "Error in PetscSynchronizedPrintf");
2362: PetscSynchronizedFlush(this->comm);
2363: CHKERROR(ierr, "Error in PetscSynchronizedFlush");
2364: }/* --------------------------------------------------------------------------------------------- */
2366: // Allocate a status array
2367: MPI_Status *Neighbor_status;
2368: if(NeighborCount) {
2369: PetscMalloc((NeighborCount)*sizeof(MPI_Status),&Neighbor_status);CHKERROR(ierr,"Error in PetscMalloc");
2370: }
2372: // Wait on the receives
2373: if(NeighborCount) {
2374: MPI_Waitall(NeighborCount, NeighborsIn_waits, Neighbor_status); CHKERROR(ierr,"Error in MPI_Waitall");
2375: }
2376:
2377: // Now we unpack the received cones and store them in the completion stack C.
2378: // Traverse all neighbors from whom we are expecting cones
2379: cntr = 0; // cone point counter
2380: for(std::map<int32_t,Point__int>::iterator np_itor=NeighborPointConeSizeIn.begin();np_itor!=NeighborPointConeSizeIn.end();np_itor++)
2381: {
2382: int32_t neighbor = (*np_itor).first;
2383: // Traverse all the points over which have received cones from neighbors
2384: // ASSUMPTION: points are sorted within each neighbor, so we are expecting points in the same order as they arrived in ConesIn
2385: for(Point__int::iterator pConeSize_itor = (*np_itor).second.begin(); pConeSize_itor != (*np_itor).second.end(); pConeSize_itor++){
2386: Point p = (*pConeSize_itor).first;
2387: Point s; // this is the structure point whose footprint will be recorded
2388: int32_t coneSize = (*pConeSize_itor).second;
2389: //Traverse the cone points received in ConesIn
2390: for(int32_t i = 0; i < coneSize; i++) {
2391: // Insert into C an arrow from/to each arrived point q to/from p, depending on completedSet: cap/base.
2392: // ASSUMPTION on Point type affects the storage of nodes in ConesIn and the usage of cntr
2393: Point q = Point(ConesIn[2*cntr], ConesIn[2*cntr+1]);
2394: Point p0; // initial point of the completion arrow
2395: Point p1; // final point of the completion arrow
2396: cntr++;
2397: if(completedSet == completedSetCap) {
2398: p0 = q;
2399: p1 = p;
2400: }
2401: else if (completedSet == completedSetBase) {
2402: p0 = p;
2403: p1 = q;
2404: }
2405: else {
2406: throw ALE::Exception("Unknown completedSet");
2407: }
2408: //
2409: P->addArrow(p0,p1);
2410: if(completionType == completionTypeArrow) {
2411: // create an arrow node and insert it into A
2412: Point a(this->commRank, cntr);
2413: A->addPoint(a);
2414: // now add the arrows from p0 to a and from a to p1
2415: A0->addArrow(p0,a); A1->addArrow(a,p1);
2416: // save the footprint point
2417: s = a;
2418: }
2419: else {
2420: // save the footprint point
2421: s = q;
2422: }
2423: if(footprintType != footprintTypeNone) {
2424: // Record the footprint of s
2425: Point f;
2426: //f.prefix = -(rank+1);
2427: f.prefix = rank;
2428: f.index = neighbor;
2429: // Insert the footprint point into F
2430: F->addPoint(f);
2431: // Add the arrow to/from f from/to s into S
2432: if(footprintType == footprintTypeSupport) {
2433: C->addArrow(f,s);
2434: }
2435: else if(footprintType == footprintTypeCone) {
2436: C->addArrow(s,f);
2437: }
2438: else {
2439: throw Exception("Unknown footprintType");
2440: }
2441: }
2442: }// for(int32_t i = 0; i < coneSize; i++)
2443: }
2444: }//for(std::map<int32_t, Point__int>::iterator np_itor=NeighborPointConeSizeIn.begin();np_itor!=NeighborPointConeSizeIn.end(); np_itor++)
2445: if (debug) { /* ----------------------------------- */
2446: ostringstream txt;
2447: cntr = 0;
2448: txt << "[" << rank << "]: debug begin\n";
2449: txt << "[" << rank << "]: __computeCompletion: completion cone: base of size " << C->_cone.size() << "\n";
2450: for(std::map<Point,Point_set,Point::Cmp>::iterator Cbase_itor=C->_cone.begin(); Cbase_itor!=C->_cone.end(); Cbase_itor++) {
2451: Point Cb = (*Cbase_itor).first;
2452: txt << "[" << rank << "]: \t(" << Cb.prefix << ", " << Cb.index << "): ";
2453: Point_set CbCone = (*Cbase_itor).second;
2454: txt << "cone of size " << CbCone.size() << ": ";
2455: for(Point_set::iterator CbCone_itor = CbCone.begin(); CbCone_itor != CbCone.end(); CbCone_itor++)
2456: {
2457: Point q = *CbCone_itor;
2458: txt << "(" << q.prefix <<"," << q.index << ") ";
2459: }
2460: txt << "\n";
2461: }
2462: PetscSynchronizedPrintf(this->comm, txt.str().c_str()); CHKERROR(ierr, "Error in PetscSynchronizedPrintf");
2463: PetscSynchronizedFlush(this->comm); CHKERROR(ierr, "PetscSynchronizedFlush");
2464: txt << "[" << rank << "]: debug end\n";
2465: }/* ----------------------------------- */
2467: // Wait on the original sends
2468: if(NeighborCount) {
2469: MPI_Waitall(NeighborCount, NeighborsOut_waits, Neighbor_status); CHKERROR(ierr,"Error in MPI_Waitall");
2470: }
2472: // Computation complete; freeing memory.
2473: // Some of these can probably be freed earlier, if memory is needed.
2474: // However, be careful while freeing memory that may be in use implicitly.
2475: // For instance, ConesOut is a send buffer and should probably be retained until all send requests have been waited on.
2476: if(NeighborCount){
2477: PetscFree(NeighborsOut_waits); CHKERROR(ierr, "Error in PetscFree");
2478: PetscFree(NeighborsIn_waits); CHKERROR(ierr, "Error in PetscFree");
2479: PetscFree(Neighbor_status); CHKERROR(ierr, "Error in PetscFree");
2480: }
2481: if(LeasedNodeCount) {
2482: PetscFree(NeighborCounts); CHKERROR(ierr,"Error in PetscFree");
2483: }
2485: if(TotalNeighborCount) {PetscFree(Neighbors); CHKERROR(ierr, "Error in PetscFree");}
2486: if(ConeSizeIn) {PetscFree(ConesIn); CHKERROR(ierr, "Error in PetscFree");}
2487: if(ConeSizeOut){PetscFree(ConesOut); CHKERROR(ierr, "Error in PetscFree");}
2488: if(TotalRentalCount){PetscFree(SharerCounts); CHKERROR(ierr, "Error in PetscFree");}
2490: if(completedSet == completedSetCap) {
2491: delete points;
2492: }
2493: else if(completedSet == completedSetBase) {
2494: delete cones;
2495: }
2496: else {
2497: throw ALE::Exception("Unknown completionType");
2498: }
2499: return C;
2500: // Done!
2501:
2502: }// PreSieve::__computeCompletion()
2505: } // namespace ALE
2507: #undef ALE_PreSieve_cxx