Actual source code: IndexBundle.cxx
1: #define ALE_IndexBundle_cxx
3: #ifndef included_ALE_IndexBundle_hh
4: #include <IndexBundle.hh>
5: #endif
7: #include <stack>
8: #include <queue>
10: namespace ALE {
12: #undef __FUNCT__
14: IndexBundle& IndexBundle::setAssemblyPolicy(BundleAssemblyPolicy policy){
15: if((policy != ADDITION) && (policy != INSERTION)) {
16: throw(Exception("Invalid BundleAssemblyPolicy value"));
17: }
18: this->_assemblyPolicy = policy;
19: return *this;
20: }//IndexBundle::setAssemblyPolicy()
22: #undef __FUNCT__
24: IndexBundle& IndexBundle::setTopology(Obj<Sieve> topology){
25: CHKCOMMS(*this, *topology.pointer());
26: // we have to reset all other PreSieves/Stacks
27: __reset(topology);
28: return *this;
29: }// IndexBundle::setTopology()
31: #undef __FUNCT__
33: IndexBundle& IndexBundle::setFiberDimension(Point element, int32_t d) {
34: ALE_LOG_STAGE_BEGIN;
35: // Setting a fiber dimension makes the bundle dirty
36: if(d < 0) {
37: throw Exception("Invalid dimension");
38: }
39: if(this->getTopology()->spaceContains(element)) {
40: Point dim(-1, d);
41: // add a point with the dimension in the index to the dimensions PreSieve
42: this->_dimensionsToElements->top()->addPoint(dim);
43: // add a single arrow from dim to the element in the dimensionsToElements Stack
44: this->_dimensionsToElements->setCone(dim, element);
45: this->__markDirty();
46: //return *this;
47: }
48: else {
49: throw ALE::Exception("Non-existent element");
50: }
51: ALE_LOG_STAGE_END;
52: return *this;
53: }//IndexBundle::setFiberDimension()
54:
55: #undef __FUNCT__
57: IndexBundle& IndexBundle::__setFiberDimensionByStratum(int stratumType, int32_t stratumIndex, int32_t dim) {
58: ALE_LOG_STAGE_BEGIN;
59: Obj<Point_set> stratum;
60: Obj<Sieve> topology = this->getTopology();
61: switch(stratumType) {
62: case stratumTypeDepth:
63: stratum = topology->depthStratum(stratumIndex);
64: break;
65: case stratumTypeHeight:
66: stratum = topology->heightStratum(stratumIndex);
67: break;
68: default:
69: throw ALE::Exception("Unknown stratum type");
70: }
71: for(Point_set::iterator s_itor = stratum->begin(); s_itor != stratum->end(); s_itor++) {
72: this->setFiberDimension(*s_itor,dim);
73: }
74: ALE_LOG_STAGE_END;
75: return *this;
76: }//IndexBundle::__setFiberDimensionByStratum()
77:
78: #undef __FUNCT__
80: int32_t IndexBundle::getFiberDimension(Obj<Point_set> base) {
81: int32_t dim = 0;
82: ALE_LOG_STAGE_BEGIN;
83: base = this->__validateChain(base);
84: for(Point_set::iterator base_itor = base->begin(); base_itor != base->end(); base_itor++) {
85: Point e = *base_itor;
86: Point_set dims = this->_dimensionsToElements->cone(e);
87: // make sure there is only one dim element
88: if(dims.size() > 1) {
89: ostringstream ex;
90: ex << "Multiple dimension designators found for element (" << e.prefix << ", " << e.index << ")";
91: throw(Exception(ex.str().c_str()));
92: }
93: // no dimension designators means 'dimension zero'
94: if(dims.size() != 0) {
95: // The actual dimension is stored in the point's index
96: Point dimPoint = *dims.begin();
97: dim += dimPoint.index;
98: }
99: }// for(Point_set::iterator base_itor = base->begin(); base_itor != base->end(); base++) {
100: ALE_LOG_STAGE_END;
101: return dim;
103: }//IndexBundle::getFiberDimension()
105: #undef __FUNCT__
107: int32_t IndexBundle::getBundleDimension(Obj<Point_set> base) {
108: int32_t dim = 0;
109: ALE_LOG_STAGE_BEGIN;
110: base = this->__validateChain(base);
111: // Traverse the closure of base and add up fiber dimensions
112: while(base->size() > 0){
113: for(Point_set::iterator base_itor = base->begin(); base_itor != base->end(); base_itor++) {
114: Point p = *base_itor;
115: dim += this->getFiberDimension(p);
116: }
117: base = this->getTopology()->cone(base);
118: }
119: ALE_LOG_STAGE_END;
120: return dim;
122: }//IndexBundle::getBundleDimension()
124: #undef __FUNCT__
126: Obj<PreSieve> IndexBundle::__computeIndices(Obj<Point_set> supports, Obj<Point_set> base, bool includeBoundary, Obj<Point_set> exclusion) {
127: Obj<PreSieve> indices = PreSieve(MPI_COMM_SELF);
128: ALE_LOG_STAGE_BEGIN;
129: base = this->__validateChain(base);
130: supports = this->__validateChain(supports);
131: // IMPROVE: we could make this subroutine consult cache, if base is singleton
132: // Traverse the boundary of s -- the closure of s except s itself -- in a depth-first search and compute the indices for the
133: // boundary fibers. For an already seen point ss, its offset is in seen[ss].prefix, or in 'off' for a newly seen point.
134: // 'off' is updated during the calculation by accumulating the dimensions of the fibers over all newly encountered elements.
135: // An element ss that has been seen AND indexed has a nonzero seen[ss].index equal to its dimension.
136: // If requested, we will cache the results in arrows to e.
138: Point__Point seen;
139: // Traverse the closure of base in a depth-first search storing the offsets of each element ss we see during the search
140: // in seen[ss].prefix.
141: // Continue searching until we encounter s -- one of supports' points, lookup or compute s's offset and store its indices
142: // in 'indices' as well as seen[ss].
143: // If none of supports are not found in the closure of base, return an empty index set.
145: // To enable the depth-first order traversal without recursion, we employ a stack.
146: std::stack<Point> stk;
147: int32_t off = 0; // keeps track of the offset within the bundle
149: // We traverse the (sub)bundle over base until one of s is encountered
150: while(1) { // use 'break' to exit the loop
152: // First we stack base elements up on stk in the reverse order to ensure offsets
153: // are computed in the bundle order.
154: for(Point_set::reverse_iterator base_ritor = base->rbegin(); base_ritor != base->rend(); base_ritor++) {
155: stk.push(*base_ritor);
156: }
158: // If the stack is empty, we are done, albeit with a negative result.
159: if(stk.empty()) { // if stk is empty
160: break;
161: }
162: else { // if stk is not empty
163: Point ss = stk.top(); stk.pop();
164: int32_t ss_dim = this->getFiberDimension(ss);
165: int32_t ss_off;
166: // If ss has already been seen, we use the stored offset.
167: if(seen.find(ss) != seen.end()){ // seen ss already
168: ss_off = seen[ss].prefix;
169: }
170: else { // not seen ss yet
171: // Compute ss_off
172: ss_off = off;
173: // Store ss_off, but ss_dim in 'seen'
174: seen[ss] = Point(ss_off, 0);
175: // off is only incremented when we encounter a not yet seen ss.
176: off += ss_dim;
177: }
178: // ss in supports ?
179: Point_set::iterator s_itor = supports->find(ss);
180: if(s_itor != supports->end()) {
181: Point s = *s_itor;
182: // If s (aka ss) has a nonzero dimension but has not been indexed yet
183: if((ss_dim > 0) && (seen[ss].index == 0)) {
184: // store ss_off, ss_dim both in indices and seen[ss]
185: Point ssPoint(ss_off, ss_dim);
186: indices->addCone(ssPoint, ss);
187: seen[ss] = Point(ss_off, ss_dim);
188: }
189: // If the boundary of s should be included in the calculation, add in the boundary indices
190: if(includeBoundary) {
191: Obj<PreSieve> boundary_indices = this->__computeBoundaryIndices(s, seen, off);
192: indices->add(boundary_indices.object());
193: }
194: } // ss in supports
195: // Regardless of whether ss is in supports or not, we need to index the elements in its closure to obtain a correct off.
196: // Compute the cone over ss, which will replace ss on stk at the beginning of the next iteration.
197: // IMPROVE: can probably reduce runtime by putting this line inside the 'if not seen' clause'.
198: // Also, this can break the numbering if we have given a very specific base and od not want other points numbered
199: base = this->getTopology()->cone(ss);
200: if (!exclusion.isNull()) {
201: base->subtract(exclusion);
202: }
203: }// if stk is not empty
204: }// while(1)
205: ALE_LOG_STAGE_END;
206: return indices;
207: }//IndexBundle::__computeIndices()
210: #undef __FUNCT__
212: Obj<PreSieve> IndexBundle::__computeBoundaryIndices(Point s, Point__Point& seen, int32_t& off) {
214: Obj<PreSieve> indices(new PreSieve(MPI_COMM_SELF));
215: ALE_LOG_STAGE_BEGIN;
216: // Traverse the boundary of s -- the closure of s except s itself -- in a depth-first search and compute the indices for the
217: // boundary fibers. For an already seen point ss, its offset is in seen[ss].prefix, or in 'off' for a newly seen point.
218: // 'off' is updated during the calculation by accumulating the dimensions of the fibers over all newly encountered elements.
219: // An element ss that has been seen AND indexed has a nonzero seen[ss].index equal to its dimension.
220: // If requested, we will cache the results in arrows to e.
222: // To enable the depth-first order traversal without recursion, we employ a stack.
223: std::stack<Point> stk;
224: // Initialize base with the cone over s
225: Obj<Point_set> base = this->getTopology()->cone(s);
226: while(1) { // use 'break' to exit the loop
227: // First we stack base elements up on stk in the reverse order to ensure that offsets
228: // are computed in the bundle order.
229: for(Point_set::reverse_iterator base_ritor = base->rbegin(); base_ritor != base->rend(); base_ritor++) {
230: stk.push(*base_ritor);
231: }
232: // If the stack is empty, we are done.
233: if(stk.empty()) { // if stk is empty
234: break;
235: }
236: else { // if stk is not empty
237: Point ss = stk.top(); stk.pop();
238: int32_t ss_dim = this->getFiberDimension(ss);
239: int32_t ss_off;
240: if(seen.find(ss) != seen.end()) { // already seen ss
241: // use the stored offset
242: ss_off = seen[ss].prefix;
243: if(seen[ss].index == 0) { // but not yet indexed ss
244: seen[ss] = Point(ss_off, ss_dim);
245: if(ss_dim > 0) {
246: Point ssPoint = Point(ss_off, ss_dim);
247: indices->addCone(ssPoint, ss);
248: }
249: } // not yet indexed ss
250: } // already seen ss
251: else // not seen ss yet
252: {
253: // use the computed offset
254: ss_off = off;
255: // Mark ss as already seen and store its offset and dimension
256: seen[ss] = Point(ss_off, ss_dim);
257: if(ss_dim > 0) {
258: Point ssPoint(ss_off, ss_dim);
259: indices->addCone(ssPoint, ss);
260: }
261: off += ss_dim;
262: } // not seen ss yet
263: base = this->getTopology()->cone(ss);
264: }// if stk is not empty
265: }// while(1)
266: ALE_LOG_STAGE_END;
267: return indices;
268: }//IndexBundle::__computeBoundaryIndices()
271: #undef __FUNCT__
273: Point IndexBundle::getFiberInterval(Point support, Obj<Point_set> base) {
274: Point interval;
275: ALE_LOG_STAGE_BEGIN;
276: Obj<PreSieve> indices = this->__computeIndices(Point_set(support), base);
277: if(indices->cap().size() == 0) {
278: interval = Point(-1,0);
279: } else {
280: interval = *(indices->cap().begin());
281: }
282: ALE_LOG_STAGE_END;
283: return interval;
284: }//IndexBundle::getFiberInterval()
287: #undef __FUNCT__
289: Obj<PreSieve> IndexBundle::getFiberIndices(Obj<Point_set> supports, Obj<Point_set> base, Obj<Point_set> exclusion) {
290: Obj<PreSieve> indices;
291: ALE_LOG_STAGE_BEGIN;
292: indices = this->__computeIndices(supports, base, false, exclusion);
293: ALE_LOG_STAGE_END;
294: return indices;
295: }//IndexBundle::getFiberIndices()
298: #undef __FUNCT__
300: Obj<PreSieve> IndexBundle::getClosureIndices(Obj<Point_set> supports, Obj<Point_set> base) {
301: bool includingBoundary = true;
302: Obj<PreSieve> indices = this->__computeIndices(supports, base, includingBoundary);
303: return indices;
304: }//IndexBundle::getClosureIndices()
306: #undef __FUNCT__
308: ALE::Point IndexBundle::__getMaxDepthElement(Obj<Point_set> points) {
309: ALE::Point max(-1, 0);
310: int32_t maxDepth = -1;
312: for(ALE::Point_set::iterator e_itor = points->begin(); e_itor != points->end(); e_itor++) {
313: int32_t depth = this->getTopology()->depth(*e_itor);
315: if (depth > maxDepth) {
316: maxDepth = depth;
317: max = *e_itor;
318: }
319: }
320: return max;
321: }
323: #undef __FUNCT__
325: ALE::int__Point IndexBundle::__checkOrderChain(Obj<Point_set> order, int& minDepth, int& maxDepth) {
326: Obj<Sieve> topology = this->getTopology();
327: int__Point dElement;
328: minDepth = 0;
329: maxDepth = 0;
331: // A topology cell-tuple contains one element per dimension, so we order the points by depth.
332: for(Point_set::iterator ord_itor = order->begin(); ord_itor != order->end(); ord_itor++) {
333: Point e = *ord_itor;
334: int32_t depth = topology->depth(e);
336: if (depth < 0) {
337: throw Exception("Invalid element: negative depth returned");
338: }
339: if (depth > maxDepth) {
340: maxDepth = depth;
341: }
342: if (depth < minDepth) {
343: minDepth = depth;
344: }
345: dElement[depth] = e;
346: }
347: // Verify that the chain is a "baricentric chain", i.e. it starts at depth 0
348: if(minDepth != 0) {
349: throw Exception("Invalid order chain: minimal depth is nonzero");
350: }
351: // Verify the chain has an element of every depth between 0 and maxDepth
352: // and that each element at depth d is in the cone of the element at depth d+1
353: for(int32_t d = 0; d <= maxDepth; d++) {
354: int__Point::iterator d_itor = dElement.find(d);
356: if(d_itor == dElement.end()){
357: ostringstream ex;
358: ex << "[" << this->getCommRank() << "]: ";
359: ex << "Missing Point at depth " << d;
360: throw Exception(ex.str().c_str());
361: }
362: if(d > 0) {
363: if(!topology->coneContains(dElement[d], dElement[d-1])){
364: ostringstream ex;
365: ex << "[" << this->getCommRank() << "]: ";
366: ex << "point (" << dElement[d-1].prefix << ", " << dElement[d-1].index << ") at depth " << d-1 << " not in the cone of ";
367: ex << "point (" << dElement[d].prefix << ", " << dElement[d].index << ") at depth " << d;
368: throw Exception(ex.str().c_str());
369: }
370: }
371: }
372: return dElement;
373: }
375: #undef __FUNCT__
377: void IndexBundle::__orderElement(int dim, ALE::Point element, std::map<int, std::queue<Point> > *ordered, ALE::Obj<ALE::Point_set> elementsOrdered) {
378: if (elementsOrdered->find(element) != elementsOrdered->end()) return;
379: (*ordered)[dim].push(element);
380: elementsOrdered->insert(element);
381: if (this->verbosity > 10) {printf(" ordered element (%d, %d) dim %d\n", element.prefix, element.index, dim);}
382: }
384: #undef __FUNCT__
386: ALE::Point IndexBundle::__orderCell(int dim, int__Point *orderChain, std::map<int, std::queue<Point> > *ordered, ALE::Obj<ALE::Point_set> elementsOrdered) {
387: Obj<Sieve> closure = this->getTopology()->closureSieve(Point_set((*orderChain)[dim]), Obj<Sieve>(new Sieve(MPI_COMM_SELF)));
388: ALE::Point last;
390: if (this->verbosity > 10) {
391: printf("Ordering cell (%d, %d) dim %d\n", (*orderChain)[dim].prefix, (*orderChain)[dim].index, dim);
392: for(int d = 0; d < dim; d++) {
393: printf(" orderChain[%d] (%d, %d)\n", d, (*orderChain)[d].prefix, (*orderChain)[d].index);
394: }
395: }
396: if (dim == 0) {
397: last = (*orderChain)[0];
398: this->__orderElement(0, last, ordered, elementsOrdered);
399: return last;
400: } else if (dim == 1) {
401: Obj<Point_set> flip = closure->cone((*orderChain)[1]);
403: this->__orderElement(0, (*orderChain)[0], ordered, elementsOrdered);
404: flip->erase((*orderChain)[0]);
405: last = *flip->begin();
406: this->__orderElement(0, last, ordered, elementsOrdered);
407: this->__orderElement(1, (*orderChain)[1], ordered, elementsOrdered);
408: (*orderChain)[dim-1] = last;
409: return last;
410: }
411: do {
412: last = this->__orderCell(dim-1, orderChain, ordered, elementsOrdered);
413: if (this->verbosity > 10) {printf(" last (%d, %d)\n", last.prefix, last.index);}
414: Obj<Point_set> faces = closure->support(last);
415: faces->erase((*orderChain)[dim-1]);
416: if (faces->size() != 1) {
417: throw Exception("Last edge did not separate two faces");
418: }
419: last = (*orderChain)[dim-1];
420: (*orderChain)[dim-1] = *faces->begin();
421: } while(elementsOrdered->find((*orderChain)[dim-1]) == elementsOrdered->end());
422: if (this->verbosity > 10) {
423: printf("Finish ordering for cell (%d, %d)\n", (*orderChain)[dim].prefix, (*orderChain)[dim].index);
424: printf(" with last (%d, %d)\n", last.prefix, last.index);
425: }
426: (*orderChain)[dim-1] = last;
427: this->__orderElement(dim, (*orderChain)[dim], ordered, elementsOrdered);
428: return last;
429: }
431: #undef __FUNCT__
433: Obj<Point_array> IndexBundle::getOrderedIndices(Obj<Point_set> order, Obj<PreSieve> indices) {
434: // IMPROVE: the ordering algorithm employed here works only for 'MeshSieves' --
435: // Sieves with the property that any pair of elements of depth d > 0 share at most
436: // one element of depth d-1. 'MeshSieve' would check that this property is preserved
437: // under arrow additions.
438: // We should define class 'MeshBundle' that would take 'MeshSieves' as topology and
439: // move this method there.
440: int minDepth, maxDepth;
441: ALE::int__Point dElement = this->__checkOrderChain(order, minDepth, maxDepth);
442: // We store the elements ordered in each dimension
443: std::map<int, std::queue<Point> > ordered;
444: // Set of the elements already ordered
445: ALE::Point_set elementsOrdered;
447: // Order elements in the closure
448: if (this->verbosity > 10) {printf("Ordering (%d, %d)\n", dElement[maxDepth].prefix, dElement[maxDepth].index);}
449: ALE::Point last = this->__orderCell(maxDepth, &dElement, &ordered, elementsOrdered);
451: // Generate indices from ordered elements
452: Obj<Point_array> indexArray(new Point_array);
453: for(int d = 0; d <= maxDepth; d++) {
454: while(!ordered[d].empty()) {
455: Obj<Point_set> indCone = indices->cone(ordered[d].front());
457: ordered[d].pop();
458: if (this->verbosity > 10) {printf(" indices (%d, %d)\n", indCone->begin()->prefix, indCone->begin()->index);}
459: if (indCone->begin()->index > 0) {
460: indexArray->push_back(*indCone->begin());
461: }
462: }
463: }
464: return indexArray;
465: }
467: #undef __FUNCT__
469: void IndexBundle::computeOverlapIndices() {
470: ALE_LOG_EVENT_BEGIN
471: this->_overlapOwnership = this->getTopology()->baseFootprint(PreSieve::completionTypePoint, PreSieve::footprintTypeCone, NULL)->left();
472: if (this->verbosity > 10) {this->_overlapOwnership->view("Overlap ownership");}
473: this->_localOverlapIndices->setBottom(this->_overlapOwnership);
474: this->_localOverlapIndices->setTop(new PreSieve(this->getComm()));
475: // Traverse the points in the overlapOwnership base, compute the local indices over it and attach them using _localOverlapIndices
476: ALE::Obj<ALE::Point_set> base = this->_overlapOwnership->base();
477: for(Point_set::iterator o_itor = base->begin(); o_itor != base->end(); o_itor++) {
478: Point e = *o_itor;
479: Point interval = this->getFiberInterval(e);
480: this->_localOverlapIndices->top()->addCapPoint(interval);
481: this->_localOverlapIndices->addCone(interval, e);
482: }
483: if (this->verbosity > 10) {this->_localOverlapIndices->view("Local overlap indices");}
484: // Now we do the completion
485: this->_localOverlapIndices->coneCompletion(PreSieve::completionTypeArrow, PreSieve::footprintTypeCone, this->_remoteOverlapIndices);
486: if (this->verbosity > 10) {this->_remoteOverlapIndices->view("Remote overlap indices");}
487: ALE_LOG_EVENT_END
488: }//IndexBundle::computeOverlapIndices()
490: #undef __FUNCT__
492: int32_t IndexBundle::getOverlapSize() {
493: return getFiberDimension(this->_overlapOwnership->base());
494: }
495:
496: #undef __FUNCT__
498: Obj<Point_set> IndexBundle::getOverlapOwners(Point e) {
499: return this->_overlapOwnership->cone(e);
500: }//IndexBundle::getOverlapOwners()
502: #undef __FUNCT__
504: Obj<PreSieve> IndexBundle::getOverlapFiberIndices(Obj<Point_set> supports, int32_t rank) {
505: supports = this->__validateChain(supports);
506: Obj<PreSieve> indices(new PreSieve(this->getComm()));
507: if (rank == this->commRank) {
508: for(ALE::Point_set::iterator e_itor = supports->begin(); e_itor != supports->end(); e_itor++) {
509: ALE::Point e = *e_itor;
510: //ALE::Obj<ALE::Point_set> ind = this->_localOverlapIndices->cone(e);
511: ALE::Point_set ind = this->_localOverlapIndices->cone(e);
513: indices->addCone(ind, e);
514: }
515: } else {
516: throw Exception("Not supported: Remote overlap indices not done");
517: }
518: return indices;
519: }//IndexBundle::getOverlapFiberIndices()
521: #undef __FUNCT__
523: Obj<PreSieve> IndexBundle::getOverlapClosureIndices(Obj<Point_set> supports, int32_t rank) {
524: supports = this->__validateChain(supports);
525: Obj<PreSieve> indices(new PreSieve(this->getComm()));
527: if (rank == this->commRank) {
528: Obj<Point_set> points = this->getTopology()->closure(supports);
530: for(ALE::Point_set::iterator e_itor = points->begin(); e_itor != points->end(); e_itor++) {
531: ALE::Point e = *e_itor;
532: //ALE::Obj<ALE::Point_set> ind = this->_localOverlapIndices->cone(e);
533: ALE::Point_set ind = this->_localOverlapIndices->cone(e);
535: indices->addCone(ind, e);
536: }
537: } else {
538: throw Exception("Not supported: Remote overlap indices not done");
539: }
540: return indices;
541: }
543: #undef __FUNCT__
545: ALE::PointType IndexBundle::__getPointType(ALE::Point point) {
546: ALE::Obj<ALE::Point_set> owners = getOverlapOwners(point);
547: ALE::PointType pointType = localPoint;
549: if (owners->size()) {
550: for(ALE::Point_set::iterator o_itor = owners->begin(); o_itor != owners->end(); o_itor++) {
551: if ((*o_itor).index < this->commRank) {
552: pointType = rentedPoint;
553: break;
554: }
555: }
556: if (pointType == localPoint) {
557: pointType = leasedPoint;
558: }
559: }
560: return pointType;
561: }
563: #undef __FUNCT__
565: /*
566: This PreSieve classifies points in the topology. Marker points are in the base and topological points are in the cap.
568: (rank, PointType) is covered by topological points of that type
569: ( p, p) is covered by topological points leased from process p
570: (-p-1, p) is covered by topological points rented to process p
571: */
572: ALE::Obj<ALE::PreSieve> IndexBundle::__computePointTypes() {
573: ALE::Obj<ALE::Point_set> space = this->getTopology()->space();
574: ALE::Obj<ALE::PreSieve> pointTypes(new PreSieve(this->comm));
576: for(ALE::Point_set::iterator e_itor = space->begin(); e_itor != space->end(); e_itor++) {
577: ALE::Point point = *e_itor;
578: ALE::Obj<ALE::Point_set> owners = getOverlapOwners(point);
579: ALE::PointType pointType = localPoint;
580: ALE::Point typePoint;
582: if (owners->size()) {
583: for(ALE::Point_set::iterator o_itor = owners->begin(); o_itor != owners->end(); o_itor++) {
584: if ((*o_itor).index < this->commRank) {
585: typePoint = ALE::Point(-(*o_itor).prefix-1, (*o_itor).index);
587: pointType = rentedPoint;
588: pointTypes->addSupport(point, typePoint);
589: break;
590: }
591: }
592: if (pointType == localPoint) {
593: pointType = leasedPoint;
594: pointTypes->addSupport(point, owners);
595: }
596: }
597: typePoint = ALE::Point(this->commRank, pointType);
598: pointTypes->addCone(point, typePoint);
599: }
600: this->_pointTypes = pointTypes;
601: if (this->verbosity > 10) {pointTypes->view("Point types");}
602: return pointTypes;
603: }
605: #undef __FUNCT__
607: ALE::Obj<ALE::PreSieve> IndexBundle::getPointTypes() {
608: return this->_pointTypes;
609: }
611: #undef __FUNCT__
613: void IndexBundle::computeGlobalIndices() {
614: ALE_LOG_EVENT_BEGIN
615: // Make local indices
616: ALE::Point_set localTypes;
617: localTypes.insert(ALE::Point(this->commRank, ALE::localPoint));
618: localTypes.insert(ALE::Point(this->commRank, ALE::leasedPoint));
619: ALE::Obj<ALE::PreSieve> pointTypes = this->__computePointTypes();
620: ALE::Obj<ALE::Point_set> localPoints = pointTypes->cone(localTypes);
621: ALE::Obj<ALE::Point_set> rentedPoints = pointTypes->cone(ALE::Point(this->commRank, ALE::rentedPoint));
622: if (this->verbosity > 10) {localPoints->view("Global local points");}
623: ALE::Obj<ALE::PreSieve> localIndices = this->getFiberIndices(localPoints, localPoints, rentedPoints);
624: if (this->verbosity > 10) {localIndices->view("Global local indices");}
625: int localSize = this->getFiberDimension(localPoints);
626: if (this->verbosity > 10) {
627: PetscSynchronizedPrintf(this->comm, "[%d]Local size %d\n", this->commRank, localSize);
628: PetscSynchronizedFlush(this->comm);
629: }
631: // Calculate global size
632: ALE::Obj<ALE::PreSieve> globalIndices(new PreSieve(this->comm));
633: this->_firstGlobalIndex = new int[this->commSize+1];
634: int MPI_Allgather(&localSize, 1, MPI_INT, &(this->_firstGlobalIndex[1]), 1, MPI_INT, this->comm);
635: CHKMPIERROR(ierr, ERRORMSG("Error in MPI_Allgather"));
636: this->_firstGlobalIndex[0] = 0;
637: for(int p = 0; p < this->commSize; p++) {
638: this->_firstGlobalIndex[p+1] = this->_firstGlobalIndex[p+1] + this->_firstGlobalIndex[p];
639: }
641: // Add rented points
642: for(ALE::Point_set::iterator r_itor = rentedPoints->begin(); r_itor != rentedPoints->end(); r_itor++) {
643: ALE::Point p = *r_itor;
644: int32_t dim = this->getFiberDimension(p);
645: ALE::Point interval(localSize, dim);
647: localIndices->addCone(interval, p);
648: localSize += dim;
649: }
651: // Make global indices
652: for(ALE::Point_set::iterator e_itor = localPoints->begin(); e_itor != localPoints->end(); e_itor++) {
653: ALE::Obj<ALE::Point_set> cone = localIndices->cone(*e_itor);
654: ALE::Point globalIndex(-1, 0);
655: ALE::Point point = *e_itor;
657: if (cone->size()) {
658: globalIndex = *cone->begin();
659: if (this->verbosity > 10) {
660: PetscSynchronizedPrintf(this->comm, "[%d] local interval (%d, %d) for point (%d, %d)\n", this->commRank, globalIndex.prefix, globalIndex.index, (*e_itor).prefix, (*e_itor).index);
661: }
662: globalIndex.prefix += this->_firstGlobalIndex[this->commRank];
663: if (this->verbosity > 10) {
664: PetscSynchronizedPrintf(this->comm, "[%d] global interval (%d, %d) for point (%d, %d)\n", this->commRank, globalIndex.prefix, globalIndex.index, (*e_itor).prefix, (*e_itor).index);
665: }
666: }
667: globalIndices->addCone(globalIndex, point);
668: }
669: if (this->verbosity > 10) {
670: PetscSynchronizedPrintf(this->comm, "[%d]Global size %d\n", this->commRank, this->_firstGlobalIndex[this->commSize]);
671: PetscSynchronizedFlush(this->comm);
672: }
674: // FIX: Communicate remote indices
675: MPI_Request *requests = new MPI_Request[this->commSize];
676: MPI_Status *statuses = new MPI_Status[this->commSize];
677: int **recvIntervals = new int *[this->commSize];
678: for(int p = 0; p < this->commSize; p++) {
679: int size;
680: if (p == this->commRank) {
681: size = 0;
682: } else {
683: ALE::Obj<ALE::Point_set> rentedPoints = pointTypes->cone(ALE::Point(-p-1, p));
684: size = rentedPoints->size()*2;
685: }
686: recvIntervals[p] = new int[size+1];
688: MPI_Irecv(recvIntervals[p], size, MPI_INT, p, 1, this->comm, &(requests[p]));
689: CHKMPIERROR(ierr, ERRORMSG("Error in MPI_Irecv"));
690: if (this->verbosity > 10) {PetscSynchronizedPrintf(this->comm, "[%d] rented size %d for proc %d\n", this->commRank, size, p);}
691: }
692: for(int p = 0; p < this->commSize; p++) {
693: if (p == this->commRank) {
694: MPI_Send(&p, 0, MPI_INT, p, 1, this->comm);
695: CHKMPIERROR(ierr, ERRORMSG("Error in MPI_Send"));
696: continue;
697: }
698: ALE::Obj<ALE::Point_set> leasedPoints = pointTypes->cone(ALE::Point(p, p));
699: int size = leasedPoints->size()*2;
700: int *intervals = new int[size+1];
701: int i = 0;
703: for(ALE::Point_set::iterator e_itor = leasedPoints->begin(); e_itor != leasedPoints->end(); e_itor++) {
704: ALE::Obj<ALE::Point_set> cone = globalIndices->cone(*e_itor);
705: ALE::Point interval(-1, 0);
707: if (cone->size()) {
708: interval = *cone->begin();
709: }
711: intervals[i++] = interval.prefix;
712: intervals[i++] = interval.index;
713: }
714: MPI_Send(intervals, size, MPI_INT, p, 1, this->comm);
715: CHKMPIERROR(ierr, ERRORMSG("Error in MPI_Send"));
716: delete intervals;
717: if (this->verbosity > 10) {PetscSynchronizedPrintf(this->comm, "[%d] leased size %d for proc %d\n", this->commRank, size, p);}
718: }
719: MPI_Waitall(this->commSize, requests, statuses); CHKMPIERROR(ierr, ERRORMSG("Error in MPI_Waitall"));
720: for(int p = 0; p < this->commSize; p++) {
721: if (p == this->commRank) {
722: delete [] recvIntervals[p];
723: continue;
724: }
725: ALE::Obj<ALE::Point_set> rentedPoints = pointTypes->cone(ALE::Point(-p-1, p));
726: int i = 0;
728: for(ALE::Point_set::iterator e_itor = rentedPoints->begin(); e_itor != rentedPoints->end(); e_itor++) {
729: ALE::Point interval(recvIntervals[p][i], recvIntervals[p][i+1]);
730: ALE::Point point = *e_itor;
732: globalIndices->addCone(interval, point);
733: if (this->verbosity > 10) {PetscSynchronizedPrintf(this->comm, "[%d]Set global indices of (%d, %d) to (%d, %d)\n", this->commRank, point.prefix, point.index, interval.prefix, interval.index);}
734: i += 2;
735: }
736: delete recvIntervals[p];
737: }
738: if (this->verbosity > 10) {PetscSynchronizedFlush(this->comm);}
739: delete [] requests;
740: delete [] statuses;
741: delete [] recvIntervals;
742: this->_localIndices = localIndices;
743: this->_globalIndices = globalIndices;
744: ALE_LOG_EVENT_END
745: }
747: #undef __FUNCT__
749: int32_t IndexBundle::getLocalSize() {
750: ALE::Point_set localTypes;
751: localTypes.insert(ALE::Point(this->commRank, ALE::localPoint));
752: localTypes.insert(ALE::Point(this->commRank, ALE::leasedPoint));
753: return getFiberDimension(this->_pointTypes->cone(localTypes));
754: }
756: #undef __FUNCT__
758: int* IndexBundle::getLocalSizes() {
759: return this->_firstGlobalIndex;
760: }
762: #undef __FUNCT__
764: int32_t IndexBundle::getGlobalSize() {
765: int localSize = getLocalSize(), globalSize;
766: int MPI_Allreduce(&localSize, &globalSize, 1, MPI_INT, MPI_SUM, this->comm);
767: CHKMPIERROR(ierr, ERRORMSG("Error in MPI_Allreduce"));
768: return globalSize;
769: }
771: #undef __FUNCT__
773: int32_t IndexBundle::getRemoteSize() {
774: return getFiberDimension(this->_pointTypes->cone(ALE::Point(this->commRank, ALE::rentedPoint)));
775: }
777: #undef __FUNCT__
779: ALE::Obj<ALE::PreSieve> IndexBundle::getGlobalIndices() {
780: return this->_globalIndices;
781: }
783: #undef __FUNCT__
785: ALE::Obj<ALE::PreSieve> IndexBundle::getLocalIndices() {
786: return this->_localIndices;
787: }
789: #undef __FUNCT__
791: Point IndexBundle::getGlobalFiberInterval(Point support) {
792: ALE::Obj<ALE::Point_set> indices = this->_globalIndices->cone(support);
793: Point interval;
795: if(indices->size() == 0) {
796: interval = Point(-1,0);
797: } else {
798: interval = *(indices->begin());
799: }
800: return interval;
801: }
803: #undef __FUNCT__
805: ALE::Obj<ALE::PreSieve> IndexBundle::getGlobalFiberIndices(Obj<Point_set> supports) {
806: supports = this->__validateChain(supports);
807: Obj<PreSieve> indices(new PreSieve(MPI_COMM_SELF));
809: for(ALE::Point_set::iterator e_itor = supports->begin(); e_itor != supports->end(); e_itor++) {
810: ALE::Point point = *e_itor;
811: ALE::Point_set cone = this->_globalIndices->cone(point);
813: indices->addCone(cone, point);
814: }
815: return indices;
816: }
818: #undef __FUNCT__
820: ALE::Obj<ALE::PreSieve> IndexBundle::getLocalFiberIndices(Obj<Point_set> supports) {
821: supports = this->__validateChain(supports);
822: Obj<PreSieve> indices(new PreSieve(MPI_COMM_SELF));
824: for(ALE::Point_set::iterator e_itor = supports->begin(); e_itor != supports->end(); e_itor++) {
825: ALE::Point point = *e_itor;
826: ALE::Point_set cone = this->_localIndices->cone(point);
828: indices->addCone(cone, point);
829: }
830: return indices;
831: }
833: #undef __FUNCT__
835: ALE::Obj<ALE::PreSieve> IndexBundle::getGlobalClosureIndices(Obj<Point_set> supports) {
836: supports = this->__validateChain(supports);
837: Obj<PreSieve> indices(new PreSieve(MPI_COMM_SELF));
838: Obj<Point_set> points = this->getTopology()->closure(supports);
840: for(ALE::Point_set::iterator e_itor = points->begin(); e_itor != points->end(); e_itor++) {
841: ALE::Point point = *e_itor;
842: ALE::Point_set cone = this->_globalIndices->cone(point);
844: indices->addCone(cone, point);
845: }
846: return indices;
847: }
849: #undef __FUNCT__
851: ALE::Obj<ALE::PreSieve> IndexBundle::getLocalClosureIndices(Obj<Point_set> supports) {
852: supports = this->__validateChain(supports);
853: Obj<PreSieve> indices(new PreSieve(MPI_COMM_SELF));
854: Obj<Point_set> points = this->getTopology()->closure(supports);
856: for(ALE::Point_set::iterator e_itor = points->begin(); e_itor != points->end(); e_itor++) {
857: ALE::Point point = *e_itor;
858: ALE::Point_set cone = this->_localIndices->cone(point);
860: indices->addCone(cone, point);
861: }
862: return indices;
863: }
865: #undef __FUNCT__
867: void IndexBundle::__postIntervalRequests(ALE::Obj<ALE::PreSieve> pointTypes, int__Point rentMarkers, MPI_Request *intervalRequests[], int **receivedIntervals[]) {
868: MPI_Request *requests = new MPI_Request[this->commSize];
869: int **recvIntervals = new int *[this->commSize];
872: for(int p = 0; p < this->commSize; p++) {
873: int size;
874: if (p == this->commRank) {
875: size = 0;
876: } else {
877: ALE::Obj<ALE::Point_set> rentedPoints = pointTypes->cone(rentMarkers[p]);
878: size = rentedPoints->size()*2;
879: }
880: recvIntervals[p] = new int[size+1];
882: MPI_Irecv(recvIntervals[p], size, MPI_INT, p, 1, this->comm, &(requests[p]));
883: CHKMPIERROR(ierr, ERRORMSG("Error in MPI_Irecv"));
884: if (this->verbosity > 10) {PetscSynchronizedPrintf(this->comm, "[%d] rented size %d for proc %d\n", this->commRank, size, p);}
885: }
886: PetscSynchronizedFlush(this->comm);
887: *intervalRequests = requests;
888: *receivedIntervals = recvIntervals;
889: }
891: #undef __FUNCT__
893: void IndexBundle::__sendIntervals(ALE::Obj<ALE::PreSieve> pointTypes, int__Point leaseMarkers, ALE::Obj<ALE::PreSieve> indices) {
896: for(int p = 0; p < this->commSize; p++) {
897: if (p == this->commRank) {
898: MPI_Send(&p, 0, MPI_INT, p, 1, this->comm);
899: CHKMPIERROR(ierr, ERRORMSG("Error in MPI_Send"));
900: continue;
901: }
902: ALE::Obj<ALE::Point_set> leasedPoints = pointTypes->cone(leaseMarkers[p]);
903: int size = leasedPoints->size()*2;
904: int *intervals = new int[size+1];
905: int i = 0;
907: for(ALE::Point_set::iterator e_itor = leasedPoints->begin(); e_itor != leasedPoints->end(); e_itor++) {
908: ALE::Obj<ALE::Point_set> cone = indices->cone(*e_itor);
909: ALE::Point interval(-1, 0);
911: if (cone->size()) {
912: interval = *cone->begin();
913: }
914: intervals[i++] = interval.prefix;
915: intervals[i++] = interval.index;
916: }
917: MPI_Send(intervals, size, MPI_INT, p, 1, this->comm);
918: CHKMPIERROR(ierr, ERRORMSG("Error in MPI_Send"));
919: delete intervals;
920: if (this->verbosity > 10) {PetscSynchronizedPrintf(this->comm, "[%d] leased size %d for proc %d\n", this->commRank, size, p);}
921: }
922: PetscSynchronizedFlush(this->comm);
923: }
925: #undef __FUNCT__
927: void IndexBundle::__receiveIntervals(ALE::Obj<ALE::PreSieve> pointTypes, int__Point rentMarkers, MPI_Request *requests, int *recvIntervals[], ALE::Obj<ALE::PreSieve> indices) {
928: MPI_Status *statuses = new MPI_Status[this->commSize];
931: MPI_Waitall(this->commSize, requests, statuses); CHKMPIERROR(ierr, ERRORMSG("Error in MPI_Waitall"));
932: for(int p = 0; p < this->commSize; p++) {
933: if (p == this->commRank) {
934: delete [] recvIntervals[p];
935: continue;
936: }
937: ALE::Obj<ALE::Point_set> rentedPoints = pointTypes->cone(rentMarkers[p]);
938: int i = 0;
940: for(ALE::Point_set::iterator e_itor = rentedPoints->begin(); e_itor != rentedPoints->end(); e_itor++) {
941: ALE::Point interval(recvIntervals[p][i], recvIntervals[p][i+1]);
942: ALE::Point point = *e_itor;
944: indices->addCone(interval, point);
945: if (this->verbosity > 10) {
946: PetscSynchronizedPrintf(this->comm, "[%d]Set indices of (%d, %d) to (%d, %d)\n", this->commRank, point.prefix, point.index, interval.prefix, interval.index);
947: }
948: i += 2;
949: }
950: delete recvIntervals[p];
951: }
952: PetscSynchronizedFlush(this->comm);
953: delete [] requests;
954: delete [] statuses;
955: delete [] recvIntervals;
956: }
958: #undef __FUNCT__
960: /*
961: Creates a mapping (scatter) from this bundle to the target bundle
963: This currently uses global indices to interface with VecScatter, but could
964: use another communication structure with purely local indexing.
965: */
966: ALE::Obj<ALE::Stack> IndexBundle::computeMappingIndices(ALE::Obj<ALE::PreSieve> pointTypes, ALE::Obj<ALE::IndexBundle> target) {
967: if (this->verbosity > 10) {pointTypes->view("Mapping point types");}
968: ALE::Obj<ALE::Stack> stack(new ALE::Stack(this->comm));
970: // Make global source indices (local + leased)
971: ALE::Point_set sourceTypes;
972: sourceTypes.insert(ALE::Point(this->commRank, ALE::localPoint));
973: sourceTypes.insert(ALE::Point(this->commRank, ALE::leasedPoint));
974: ALE::Obj<ALE::Point_set> sourcePoints = pointTypes->cone(sourceTypes);
975: // Need to implement copy(), then use restrictBase()
976: ALE::Obj<ALE::PreSieve> sourceIndices(new ALE::PreSieve(this->comm));
977: for(ALE::Point_set::iterator e_itor = sourcePoints->begin(); e_itor != sourcePoints->end(); e_itor++) {
978: ALE::Point point = *e_itor;
979: ALE::Point_set cone = this->_globalIndices->cone(*e_itor);
980: sourceIndices->addCone(cone, point);
981: }
982: int sourceSize = this->getFiberDimension(sourcePoints);
983: if (this->verbosity > 10) {
984: PetscSynchronizedPrintf(this->comm, "[%d]Source size %d\n", this->commRank, sourceSize);
985: PetscSynchronizedFlush(this->comm);
986: }
988: // Make initial global target indices (local)
989: ALE::Point_set targetTypes;
990: targetTypes.insert(ALE::Point(this->commRank, ALE::localPoint));
991: ALE::Obj<ALE::Point_set> targetPoints = pointTypes->cone(targetTypes);
992: // Need to implement copy(), then use restrictBase()
993: ALE::Obj<ALE::PreSieve> targetGlobalIndices = target->getGlobalIndices();
994: ALE::Obj<ALE::PreSieve> targetIndices(new ALE::PreSieve(this->comm));
995: for(ALE::Point_set::iterator e_itor = targetPoints->begin(); e_itor != targetPoints->end(); e_itor++) {
996: ALE::Point point = *e_itor;
997: ALE::Point_set cone = targetGlobalIndices->cone(*e_itor);
998: targetIndices->addCone(cone, point);
999: }
1001: int__Point rentMarkers, leaseMarkers;
1002: MPI_Request *requests;
1003: int **recvIntervals;
1005: for(int p = 0; p < this->commSize; p++) {
1006: rentMarkers[p] = ALE::Point(-p-1, p);
1007: leaseMarkers[p] = ALE::Point(p, p);
1008: }
1010: // Send leased indices from source to target, Accept rented indices at target from source
1011: this->__postIntervalRequests(pointTypes, leaseMarkers, &requests, &recvIntervals);
1012: this->__sendIntervals(pointTypes, rentMarkers, targetGlobalIndices);
1013: this->__receiveIntervals(pointTypes, leaseMarkers, requests, recvIntervals, targetIndices);
1015: if (this->verbosity > 10) {
1016: sourceIndices->view("Source indices");
1017: targetIndices->view("Target indices");
1018: }
1019: stack->setTop(sourceIndices);
1020: stack->setBottom(targetIndices);
1021: return stack;
1022: }
1024: #undef __FUNCT__
1026: Point IndexBundle::__getArrow(Point e1, Point e) {
1027: Point_set arrows = this->__getArrows()->nMeet(this->_arrowsToStarts->cone(e1),this->_arrowsToEnds->cone(e),0);
1028: if(arrows.size() > 1) {
1029: throw(Exception("Multiple arrows attached to an element pair"));
1030: }
1031: Point arrow;
1032: if(arrows.size() == 0) {
1033: // We must add a new arrow. How do we ensure it is indeed new?
1034: // We insist of prefix == this->commRank; then all __getArrows()->cap() points are sorted by the index,
1035: // so we take the last + 1.
1036: arrow.prefix = this->commRank;
1037: arrow.index = (--(this->__getArrows()->cap().end()))->index + 1;
1038: this->__getArrows()->addPoint(arrow);
1039: this->_arrowsToStarts->addCone(arrow,e1);
1040: this->_arrowsToEnds->addCone(arrow,e);
1041: }
1042: else {// There already exists a unique arrow, so return it.
1043: arrow = *(arrows.begin());
1044: }
1045: return arrow;
1046:
1047: }// IndexBundle::__getArrow()
1050: #undef __FUNCT__
1052: Point IndexBundle::__getArrowInterval(Point e1, Point e) {
1053: // If the bundle is dirty, we reset all the indices first
1054: if(!this->__isClean()) {
1055: this->__resetArrowIndices();
1056: this->__markClean();
1057: return Point(-1,0);
1058: }
1059: // Retrieve the arrow between e1 and e
1060: Point arrow = this->__getArrow(e1,e);
1061: // Now retrieve the index set arrached to arrow
1062: Obj<Point_set> indices = this->_indicesToArrows->cone(arrow);
1063: if(indices->size() != 1) { // either nothing or too many things cached
1064: return Point(-1,0);
1065: }
1066: return *(indices->begin());
1067: }// IndexBundle::__getArrowInterval()
1069: #undef __FUNCT__
1071: void IndexBundle::__setArrowInterval(Point e1, Point e, Point interval) {
1072: // If the bundle is dirty, we reset all the indices first
1073: if(!this->__isClean()) {
1074: this->__resetArrowIndices();
1075: this->__markClean();
1076: }
1077: // First we retrieve the arrow
1078: Point arrow = this->__getArrow(e1, e);
1079: // Now set the arrow index interval
1080: this->_indicesToArrows->setCone(interval, arrow);
1081: }// IndexBundle::__setArrowInterval()
1085: void IndexBundle::view(const char *name) {
1086: CHKCOMM(*this);
1088: ostringstream txt, hdr;
1089: hdr << "Viewing ";
1090: if(name != NULL) {
1091: hdr << name;
1092: }
1093: hdr << " IndexBundle\n";
1094: // Print header
1095: PetscPrintf(this->comm, hdr.str().c_str()); CHKERROR(ierr, "Error in PetscPrintf");
1096:
1097: this->_dimensionsToElements->view("Dimension Assignment Stack");
1099: }// IndexBundle::view()
1102: } // namespace ALE
1104: #undef ALE_IndexBundle_cxx