Actual source code: Field.hh
1: #ifndef included_ALE_Field_hh
2: #define included_ALE_Field_hh
4: #ifndef included_ALE_SieveAlgorithms_hh
5: #include <SieveAlgorithms.hh>
6: #endif
10: // Sieve need point_type
11: // Section need point_type and value_type
12: // size(), restrict(), update() need orientation which is a Bundle (circular)
13: // Bundle is Sieve+Section
15: // An AbstractSection is a mapping from Sieve points to sets of values
16: // This is our presentation of a section of a fibre bundle,
17: // in which the Topology is the base space, and
18: // the value sets are vectors in the fiber spaces
19: // The main interface to values is through restrict() and update()
20: // This retrieval uses Sieve closure()
21: // We should call these rawRestrict/rawUpdate
22: // The Section must also be able to set/report the dimension of each fiber
23: // for which we use another section we call an \emph{Atlas}
24: // For some storage schemes, we also want offsets to go with these dimensions
25: // We must have some way of limiting the points associated with values
26: // so each section must support a getPatch() call returning the points with associated fibers
27: // I was using the Chart for this
28: // The Section must be able to participate in \emph{completion}
29: // This means restrict to a provided overlap, and exchange in the restricted sections
30: // Completion does not use hierarchy, so we see the Topology as a DiscreteTopology
31: namespace ALE {
32: template<typename Point_>
33: class DiscreteSieve {
34: public:
35: typedef Point_ point_type;
36: typedef std::vector<point_type> coneSequence;
37: typedef std::vector<point_type> coneSet;
38: typedef std::vector<point_type> coneArray;
39: typedef std::vector<point_type> supportSequence;
40: typedef std::vector<point_type> supportSet;
41: typedef std::vector<point_type> supportArray;
42: typedef std::set<point_type> points_type;
43: typedef points_type baseSequence;
44: typedef points_type capSequence;
45: protected:
46: Obj<points_type> _points;
47: Obj<coneSequence> _empty;
48: Obj<coneSequence> _return;
49: void _init() {
50: this->_points = new points_type();
51: this->_empty = new coneSequence();
52: this->_return = new coneSequence();
53: };
54: public:
55: DiscreteSieve() {
56: this->_init();
57: };
58: template<typename Input>
59: DiscreteSieve(const Obj<Input>& points) {
60: this->_init();
61: this->_points->insert(points->begin(), points->end());
62: };
63: virtual ~DiscreteSieve() {};
64: public:
65: void addPoint(const point_type& point) {
66: this->_points->insert(point);
67: };
68: template<typename Input>
69: void addPoints(const Obj<Input>& points) {
70: this->_points->insert(points->begin(), points->end());
71: };
72: const Obj<coneSequence>& cone(const point_type& p) {return this->_empty;};
73: template<typename Input>
74: const Obj<coneSequence>& cone(const Input& p) {return this->_empty;};
75: const Obj<coneSet>& nCone(const point_type& p, const int level) {
76: if (level == 0) {
77: return this->closure(p);
78: } else {
79: return this->_empty;
80: }
81: };
82: const Obj<coneArray>& closure(const point_type& p) {
83: this->_return->clear();
84: this->_return->push_back(p);
85: return this->_return;
86: };
87: const Obj<supportSequence>& support(const point_type& p) {return this->_empty;};
88: template<typename Input>
89: const Obj<supportSequence>& support(const Input& p) {return this->_empty;};
90: const Obj<supportSet>& nSupport(const point_type& p, const int level) {
91: if (level == 0) {
92: return this->star(p);
93: } else {
94: return this->_empty;
95: }
96: };
97: const Obj<supportArray>& star(const point_type& p) {
98: this->_return->clear();
99: this->_return->push_back(p);
100: return this->_return;
101: };
102: const Obj<capSequence>& roots() {return this->_points;};
103: const Obj<capSequence>& cap() {return this->_points;};
104: const Obj<baseSequence>& leaves() {return this->_points;};
105: const Obj<baseSequence>& base() {return this->_points;};
106: template<typename Color>
107: void addArrow(const point_type& p, const point_type& q, const Color& color) {
108: throw ALE::Exception("Cannot add an arrow to a DiscreteSieve");
109: };
110: void stratify() {};
111: void view(const std::string& name, MPI_Comm comm = MPI_COMM_NULL) const {
112: ostringstream txt;
113: int rank;
115: if (comm == MPI_COMM_NULL) {
116: comm = MPI_COMM_SELF;
117: rank = 0;
118: } else {
119: MPI_Comm_rank(comm, &rank);
120: }
121: if (name == "") {
122: if(rank == 0) {
123: txt << "viewing a DiscreteSieve" << std::endl;
124: }
125: } else {
126: if(rank == 0) {
127: txt << "viewing DiscreteSieve '" << name << "'" << std::endl;
128: }
129: }
130: for(typename points_type::const_iterator p_iter = this->_points->begin(); p_iter != this->_points->end(); ++p_iter) {
131: txt << " Point " << *p_iter << std::endl;
132: }
133: PetscSynchronizedPrintf(comm, txt.str().c_str());
134: PetscSynchronizedFlush(comm);
135: };
136: };
137: // A ConstantSection is the simplest Section
138: // All fibers are dimension 1
139: // All values are equal to a constant
140: // We need no value storage and no communication for completion
141: template<typename Point_, typename Value_>
142: class ConstantSection : public ALE::ParallelObject {
143: public:
144: typedef Point_ point_type;
145: typedef std::set<point_type> chart_type;
146: typedef Value_ value_type;
147: protected:
148: chart_type _chart;
149: value_type _value;
150: value_type _defaultValue;
151: public:
152: ConstantSection(MPI_Comm comm, const int debug = 0) : ParallelObject(comm, debug), _defaultValue(0) {};
153: ConstantSection(MPI_Comm comm, const value_type& value, const int debug) : ParallelObject(comm, debug), _value(value), _defaultValue(value) {};
154: ConstantSection(MPI_Comm comm, const value_type& value, const value_type& defaultValue, const int debug) : ParallelObject(comm, debug), _value(value), _defaultValue(defaultValue) {};
155: public: // Verifiers
156: void checkPoint(const point_type& point) const {
157: if (this->_chart.find(point) == this->_chart.end()) {
158: ostringstream msg;
159: msg << "Invalid section point " << point << std::endl;
160: throw ALE::Exception(msg.str().c_str());
161: }
162: };
163: void checkDimension(const int& dim) {
164: if (dim != 1) {
165: ostringstream msg;
166: msg << "Invalid fiber dimension " << dim << " must be 1" << std::endl;
167: throw ALE::Exception(msg.str().c_str());
168: }
169: };
170: bool hasPoint(const point_type& point) const {
171: return this->_chart.count(point) > 0;
172: };
173: public: // Accessors
174: const chart_type& getChart() {return this->_chart;};
175: void addPoint(const point_type& point) {
176: this->_chart.insert(point);
177: };
178: template<typename Points>
179: void addPoint(const Obj<Points>& points) {
180: this->_chart.insert(points->begin(), points->end());
181: };
182: void addPoint(const std::set<point_type>& points) {
183: this->_chart.insert(points.begin(), points.end());
184: };
185: value_type getDefaultValue() {return this->_defaultValue;};
186: void setDefaultValue(const value_type value) {this->_defaultValue = value;};
187: void copy(const Obj<ConstantSection>& section) {
188: const chart_type& chart = section->getChart();
190: this->addPoint(chart);
191: this->_value = section->restrict(*chart.begin())[0];
192: };
193: public: // Sizes
194: void clear() {
195: this->_chart.clear();
196: };
197: int getFiberDimension(const point_type& p) const {
198: if (this->hasPoint(p)) return 1;
199: return 0;
200: };
201: void setFiberDimension(const point_type& p, int dim) {
202: this->checkDimension(dim);
203: this->addPoint(p);
204: };
205: template<typename Sequence>
206: void setFiberDimension(const Obj<Sequence>& points, int dim) {
207: for(typename Sequence::iterator p_iter = points->begin(); p_iter != points->end(); ++p_iter) {
208: this->setFiberDimension(*p_iter, dim);
209: }
210: };
211: void addFiberDimension(const point_type& p, int dim) {
212: if (this->_chart.find(p) != this->_chart.end()) {
213: ostringstream msg;
214: msg << "Invalid addition to fiber dimension " << dim << " cannot exceed 1" << std::endl;
215: throw ALE::Exception(msg.str().c_str());
216: } else {
217: this->setFiberDimension(p, dim);
218: }
219: };
220: int size() {return this->_sheaf.size();};
221: int size(const point_type& p) {return this->getFiberDimension(p);};
222: public: // Restriction
223: const value_type *restrict(const point_type& p) const {
224: if (this->hasPoint(p)) {
225: return &this->_value;
226: }
227: return &this->_defaultValue;
228: };
229: const value_type *restrictPoint(const point_type& p) const {return this->restrict(p);};
230: void update(const point_type& p, const value_type v[]) {
231: this->_value = v[0];
232: };
233: void updatePoint(const point_type& p, const value_type v[]) {return this->update(p, v);};
234: void updateAdd(const point_type& p, const value_type v[]) {
235: this->_value += v[0];
236: };
237: void updateAddPoint(const point_type& p, const value_type v[]) {return this->updateAdd(p, v);};
238: public:
239: void view(const std::string& name, MPI_Comm comm = MPI_COMM_NULL) const {
240: ostringstream txt;
241: int rank;
243: if (comm == MPI_COMM_NULL) {
244: comm = this->comm();
245: rank = this->commRank();
246: } else {
247: MPI_Comm_rank(comm, &rank);
248: }
249: if (name == "") {
250: if(rank == 0) {
251: txt << "viewing a ConstantSection" << std::endl;
252: }
253: } else {
254: if(rank == 0) {
255: txt << "viewing ConstantSection '" << name << "'" << std::endl;
256: }
257: }
258: txt <<"["<<this->commRank()<<"]: Value " << this->_value << std::endl;
259: PetscSynchronizedPrintf(comm, txt.str().c_str());
260: PetscSynchronizedFlush(comm);
261: };
262: };
264: // A UniformSection often acts as an Atlas
265: // All fibers are the same dimension
266: // Note we can use a ConstantSection for this Atlas
267: // Each point may have a different vector
268: // Thus we need storage for values, and hence must implement completion
269: template<typename Point_, typename Value_, int fiberDim = 1>
270: class UniformSection : public ALE::ParallelObject {
271: public:
272: typedef Point_ point_type;
273: typedef ConstantSection<point_type, int> atlas_type;
274: typedef typename atlas_type::chart_type chart_type;
275: typedef Value_ value_type;
276: typedef struct {value_type v[fiberDim];} fiber_type;
277: typedef std::map<point_type, fiber_type> values_type;
278: protected:
279: Obj<atlas_type> _atlas;
280: values_type _array;
281: public:
282: UniformSection(MPI_Comm comm, const int debug = 0) : ParallelObject(comm, debug) {
283: this->_atlas = new atlas_type(comm, fiberDim, 0, debug);
284: };
285: UniformSection(const Obj<atlas_type>& atlas) : ParallelObject(atlas->comm(), atlas->debug()), _atlas(atlas) {
286: int dim = fiberDim;
287: this->_atlas->update(*this->_atlas->getChart().begin(), &dim);
288: };
289: public:
290: value_type *getRawArray(const int size) {
291: static value_type *array = NULL;
292: static int maxSize = 0;
294: if (size > maxSize) {
295: maxSize = size;
296: if (array) delete [] array;
297: array = new value_type[maxSize];
298: };
299: return array;
300: };
301: public: // Verifiers
302: bool hasPoint(const point_type& point) {
303: return this->_atlas->hasPoint(point);
304: };
305: void checkDimension(const int& dim) {
306: if (dim != fiberDim) {
307: ostringstream msg;
308: msg << "Invalid fiber dimension " << dim << " must be " << fiberDim << std::endl;
309: throw ALE::Exception(msg.str().c_str());
310: }
311: };
312: public: // Accessors
313: const chart_type& getChart() {return this->_atlas->getChart();};
314: const Obj<atlas_type>& getAtlas() {return this->_atlas;};
315: void setAtlas(const Obj<atlas_type>& atlas) {this->_atlas = atlas;};
316: void addPoint(const point_type& point) {
317: this->setFiberDimension(point, fiberDim);
318: };
319: template<typename Points>
320: void addPoint(const Obj<Points>& points) {
321: for(typename Points::iterator p_iter = points->begin(); p_iter != points->end(); ++p_iter) {
322: this->setFiberDimension(*p_iter, fiberDim);
323: }
324: };
325: void copy(const Obj<UniformSection>& section) {
326: this->getAtlas()->copy(section->getAtlas());
327: const chart_type& chart = section->getChart();
329: for(typename chart_type::const_iterator c_iter = chart.begin(); c_iter != chart.end(); ++c_iter) {
330: this->updatePoint(*c_iter, section->restrictPoint(*c_iter));
331: }
332: };
333: public: // Sizes
334: void clear() {
335: this->_atlas->clear();
336: this->_array.clear();
337: };
338: int getFiberDimension(const point_type& p) const {
339: return this->_atlas->restrictPoint(p)[0];
340: };
341: void setFiberDimension(const point_type& p, int dim) {
342: this->checkDimension(dim);
343: this->_atlas->addPoint(p);
344: this->_atlas->updatePoint(p, &dim);
345: };
346: template<typename Sequence>
347: void setFiberDimension(const Obj<Sequence>& points, int dim) {
348: for(typename Sequence::iterator p_iter = points->begin(); p_iter != points->end(); ++p_iter) {
349: this->setFiberDimension(*p_iter, dim);
350: }
351: };
352: void setFiberDimension(const std::set<point_type>& points, int dim) {
353: for(typename std::set<point_type>::iterator p_iter = points.begin(); p_iter != points.end(); ++p_iter) {
354: this->setFiberDimension(*p_iter, dim);
355: }
356: };
357: void addFiberDimension(const point_type& p, int dim) {
358: if (this->hasPoint(p)) {
359: ostringstream msg;
360: msg << "Invalid addition to fiber dimension " << dim << " cannot exceed " << fiberDim << std::endl;
361: throw ALE::Exception(msg.str().c_str());
362: } else {
363: this->setFiberDimension(p, dim);
364: }
365: };
366: int size() const {
367: const chart_type& points = this->_atlas->getChart();
368: int size = 0;
370: for(typename chart_type::iterator p_iter = points.begin(); p_iter != points.end(); ++p_iter) {
371: size += this->getFiberDimension(*p_iter);
372: }
373: return size;
374: };
375: int sizeWithBC() const {
376: return this->size();
377: };
378: public: // Restriction
379: // Return a pointer to the entire contiguous storage array
380: const values_type& restrict() {
381: return this->_array;
382: };
383: // Return only the values associated to this point, not its closure
384: const value_type *restrictPoint(const point_type& p) {
385: return this->_array[p].v;
386: };
387: // Update only the values associated to this point, not its closure
388: void updatePoint(const point_type& p, const value_type v[]) {
389: for(int i = 0; i < fiberDim; ++i) {
390: this->_array[p].v[i] = v[i];
391: }
392: };
393: // Update only the values associated to this point, not its closure
394: void updateAddPoint(const point_type& p, const value_type v[]) {
395: for(int i = 0; i < fiberDim; ++i) {
396: this->_array[p].v[i] += v[i];
397: }
398: };
399: void updatePointAll(const point_type& p, const value_type v[]) {
400: this->updatePoint(p, v);
401: };
402: public:
403: void view(const std::string& name, MPI_Comm comm = MPI_COMM_NULL) {
404: ostringstream txt;
405: int rank;
407: if (comm == MPI_COMM_NULL) {
408: comm = this->comm();
409: rank = this->commRank();
410: } else {
411: MPI_Comm_rank(comm, &rank);
412: }
413: if (name == "") {
414: if(rank == 0) {
415: txt << "viewing a UniformSection" << std::endl;
416: }
417: } else {
418: if(rank == 0) {
419: txt << "viewing UniformSection '" << name << "'" << std::endl;
420: }
421: }
422: const typename atlas_type::chart_type& chart = this->_atlas->getChart();
423: values_type& array = this->_array;
425: for(typename atlas_type::chart_type::const_iterator p_iter = chart.begin(); p_iter != chart.end(); ++p_iter) {
426: const point_type& point = *p_iter;
427: const typename atlas_type::value_type dim = this->_atlas->restrict(point)[0];
429: if (dim != 0) {
430: txt << "[" << this->commRank() << "]: " << point << " dim " << dim << " ";
431: for(int i = 0; i < dim; i++) {
432: txt << " " << array[point].v[i];
433: }
434: txt << std::endl;
435: }
436: }
437: if (chart.size() == 0) {
438: txt << "[" << this->commRank() << "]: empty" << std::endl;
439: }
440: PetscSynchronizedPrintf(comm, txt.str().c_str());
441: PetscSynchronizedFlush(comm);
442: };
443: };
444: // A Section is our most general construct (more general ones could be envisioned)
445: // The Atlas is a UniformSection of dimension 1 and value type Point
446: // to hold each fiber dimension and offsets into a contiguous patch array
447: template<typename Point_, typename Value_>
448: class Section : public ALE::ParallelObject {
449: public:
450: typedef Point_ point_type;
451: typedef ALE::Point index_type;
452: typedef UniformSection<point_type, index_type> atlas_type;
453: typedef typename atlas_type::chart_type chart_type;
454: typedef Value_ value_type;
455: typedef value_type * values_type;
456: protected:
457: Obj<atlas_type> _atlas;
458: Obj<atlas_type> _atlasNew;
459: values_type _array;
460: public:
461: Section(MPI_Comm comm, const int debug = 0) : ParallelObject(comm, debug) {
462: this->_atlas = new atlas_type(comm, debug);
463: this->_atlasNew = NULL;
464: this->_array = NULL;
465: };
466: Section(const Obj<atlas_type>& atlas) : ParallelObject(atlas->comm(), atlas->debug()), _atlas(atlas), _atlasNew(NULL), _array(NULL) {};
467: virtual ~Section() {
468: if (this->_array) {
469: delete [] this->_array;
470: this->_array = NULL;
471: }
472: };
473: public:
474: value_type *getRawArray(const int size) {
475: static value_type *array = NULL;
476: static int maxSize = 0;
478: if (size > maxSize) {
479: maxSize = size;
480: if (array) delete [] array;
481: array = new value_type[maxSize];
482: };
483: return array;
484: };
485: public: // Verifiers
486: bool hasPoint(const point_type& point) {
487: return this->_atlas->hasPoint(point);
488: };
489: public: // Accessors
490: const Obj<atlas_type>& getAtlas() {return this->_atlas;};
491: void setAtlas(const Obj<atlas_type>& atlas) {this->_atlas = atlas;};
492: const Obj<atlas_type>& getNewAtlas() {return this->_atlasNew;};
493: void setNewAtlas(const Obj<atlas_type>& atlas) {this->_atlasNew = atlas;};
494: const chart_type& getChart() const {return this->_atlas->getChart();};
495: public: // BC
496: // Returns the number of constraints on a point
497: const int getConstraintDimension(const point_type& p) const {
498: return std::max(0, -this->_atlas->restrictPoint(p)->prefix);
499: };
500: // Set the number of constraints on a point
501: // We only allow the entire point to be constrained, so these will be the
502: // only dofs on the point
503: void setConstraintDimension(const point_type& p, const int numConstraints) {
504: this->setFiberDimension(p, -numConstraints);
505: };
506: void copyBC(const Obj<Section>& section) {
507: const chart_type& chart = this->getChart();
509: for(typename chart_type::const_iterator p_iter = chart.begin(); p_iter != chart.end(); ++p_iter) {
510: if (this->getConstraintDimension(*p_iter)) {
511: this->updatePointBC(*p_iter, section->restrictPoint(*p_iter));
512: }
513: }
514: };
515: void defaultConstraintDof() {};
516: public: // Sizes
517: void clear() {
518: this->_atlas->clear();
519: delete [] this->_array;
520: this->_array = NULL;
521: };
522: // Return the total number of dofs on the point (free and constrained)
523: int getFiberDimension(const point_type& p) const {
524: return std::abs(this->_atlas->restrictPoint(p)->prefix);
525: };
526: int getFiberDimension(const Obj<atlas_type>& atlas, const point_type& p) const {
527: return std::abs(atlas->restrictPoint(p)->prefix);
528: };
529: // Set the total number of dofs on the point (free and constrained)
530: void setFiberDimension(const point_type& p, int dim) {
531: const index_type idx(dim, -1);
532: this->_atlas->addPoint(p);
533: this->_atlas->updatePoint(p, &idx);
534: };
535: template<typename Sequence>
536: void setFiberDimension(const Obj<Sequence>& points, int dim) {
537: for(typename Sequence::iterator p_iter = points->begin(); p_iter != points->end(); ++p_iter) {
538: this->setFiberDimension(*p_iter, dim);
539: }
540: };
541: void addFiberDimension(const point_type& p, int dim) {
542: if (this->_atlas->hasPoint(p)) {
543: const index_type values(dim, 0);
544: this->_atlas->updateAddPoint(p, &values);
545: } else {
546: this->setFiberDimension(p, dim);
547: }
548: };
549: // Return the number of constrained dofs on this point
550: // If constrained, this is equal to the fiber dimension
551: // Otherwise, 0
552: int getConstrainedFiberDimension(const point_type& p) const {
553: return std::max(0, this->_atlas->restrictPoint(p)->prefix);
554: };
555: // Return the total number of free dofs
556: int size() const {
557: const chart_type& points = this->getChart();
558: int size = 0;
560: for(typename chart_type::iterator p_iter = points.begin(); p_iter != points.end(); ++p_iter) {
561: size += this->getConstrainedFiberDimension(*p_iter);
562: }
563: return size;
564: };
565: // Return the total number of dofs (free and constrained)
566: int sizeWithBC() const {
567: const chart_type& points = this->getChart();
568: int size = 0;
570: for(typename chart_type::iterator p_iter = points.begin(); p_iter != points.end(); ++p_iter) {
571: size += this->getFiberDimension(*p_iter);
572: }
573: return size;
574: };
575: public: // Index retrieval
576: const int& getIndex(const point_type& p) {
577: return this->_atlas->restrictPoint(p)->index;
578: };
579: void setIndex(const point_type& p, const int& index) {
580: ((typename atlas_type::value_type *) this->_atlas->restrictPoint(p))->index = index;
581: };
582: void setIndexBC(const point_type& p, const int& index) {
583: this->setIndex(p, index);
584: };
585: void getIndices(const point_type& p, PetscInt indices[], PetscInt *indx, const int orientation = 1, const bool freeOnly = false) {
586: this->getIndices(p, this->getIndex(p), indices, indx, orientation, freeOnly);
587: };
588: template<typename Order_>
589: void getIndices(const point_type& p, const Obj<Order_>& order, PetscInt indices[], PetscInt *indx, const int orientation = 1, const bool freeOnly = false) {
590: this->getIndices(p, order->getIndex(p), indices, indx, orientation, freeOnly);
591: };
592: void getIndices(const point_type& p, const int start, PetscInt indices[], PetscInt *indx, const int orientation = 1, const bool freeOnly = false) {
593: const int& dim = this->getFiberDimension(p);
594: const int& cDim = this->getConstraintDimension(p);
595: const int end = start + dim;
597: if (!cDim) {
598: if (orientation >= 0) {
599: for(int i = start; i < end; ++i) {
600: indices[(*indx)++] = i;
601: }
602: } else {
603: for(int i = end-1; i >= start; --i) {
604: indices[(*indx)++] = i;
605: }
606: }
607: } else {
608: if (!freeOnly) {
609: for(int i = start; i < end; ++i) {
610: indices[(*indx)++] = -1;
611: }
612: }
613: }
614: };
615: public: // Allocation
616: void allocateStorage() {
617: const int totalSize = this->sizeWithBC();
619: this->_array = new value_type[totalSize];
620: PetscMemzero(this->_array, totalSize * sizeof(value_type));
621: };
622: void replaceStorage(value_type *newArray) {
623: delete [] this->_array;
624: this->_array = newArray;
625: this->_atlasNew = NULL;
626: };
627: void addPoint(const point_type& point, const int dim) {
628: if (dim == 0) return;
629: if (this->_atlasNew.isNull()) {
630: this->_atlasNew = new atlas_type(this->comm(), this->debug());
631: this->_atlasNew->copy(this->_atlas);
632: }
633: const index_type idx(dim, -1);
634: this->_atlasNew->addPoint(point);
635: this->_atlasNew->updatePoint(point, &idx);
636: };
637: template<typename Sequence>
638: void addPoints(const Obj<Sequence>& points, const int dim) {
639: for(typename Sequence::iterator p_iter = points->begin(); p_iter != points->end(); ++p_iter) {
640: this->addPoint(*p_iter, dim);
641: }
642: };
643: void orderPoints(const Obj<atlas_type>& atlas){
644: const chart_type& chart = this->getChart();
645: int offset = 0;
646: int bcOffset = this->size();
648: for(typename chart_type::iterator c_iter = chart.begin(); c_iter != chart.end(); ++c_iter) {
649: typename atlas_type::value_type idx = atlas->restrictPoint(*c_iter)[0];
650: const int& dim = idx.prefix;
652: if (dim >= 0) {
653: idx.index = offset;
654: atlas->updatePoint(*c_iter, &idx);
655: offset += dim;
656: } else {
657: idx.index = bcOffset;
658: atlas->updatePoint(*c_iter, &idx);
659: bcOffset -= dim;
660: }
661: }
662: };
663: void allocatePoint() {
664: this->orderPoints(this->_atlas);
665: this->allocateStorage();
666: };
667: public: // Restriction and Update
668: // Zero entries
669: void zero() {
670: memset(this->_array, 0, this->size()* sizeof(value_type));
671: };
672: // Return a pointer to the entire contiguous storage array
673: const value_type *restrict() {
674: return this->_array;
675: };
676: // Update the entire contiguous storage array
677: void update(const value_type v[]) {
678: const value_type *array = this->_array;
679: const int size = this->size();
681: for(int i = 0; i < size; i++) {
682: array[i] = v[i];
683: }
684: };
685: // Return the free values on a point
686: const value_type *restrictPoint(const point_type& p) {
687: return &(this->_array[this->_atlas->restrictPoint(p)[0].index]);
688: };
689: // Update the free values on a point
690: void updatePoint(const point_type& p, const value_type v[], const int orientation = 1) {
691: const index_type& idx = this->_atlas->restrictPoint(p)[0];
692: value_type *a = &(this->_array[idx.index]);
694: if (orientation >= 0) {
695: for(int i = 0; i < idx.prefix; ++i) {
696: a[i] = v[i];
697: }
698: } else {
699: const int last = idx.prefix-1;
701: for(int i = 0; i < idx.prefix; ++i) {
702: a[i] = v[last-i];
703: }
704: }
705: };
706: // Update the free values on a point
707: void updateAddPoint(const point_type& p, const value_type v[], const int orientation = 1) {
708: const index_type& idx = this->_atlas->restrictPoint(p)[0];
709: value_type *a = &(this->_array[idx.index]);
711: if (orientation >= 0) {
712: for(int i = 0; i < idx.prefix; ++i) {
713: a[i] += v[i];
714: }
715: } else {
716: const int last = idx.prefix-1;
718: for(int i = 0; i < idx.prefix; ++i) {
719: a[i] += v[last-i];
720: }
721: }
722: };
723: // Update only the constrained dofs on a point
724: void updatePointBC(const point_type& p, const value_type v[], const int orientation = 1) {
725: const index_type& idx = this->_atlas->restrictPoint(p)[0];
726: const int dim = -idx.prefix;
727: value_type *a = &(this->_array[idx.index]);
729: if (orientation >= 0) {
730: for(int i = 0; i < dim; ++i) {
731: a[i] = v[i];
732: }
733: } else {
734: const int last = dim-1;
736: for(int i = 0; i < dim; ++i) {
737: a[i] = v[last-i];
738: }
739: }
740: };
741: // Update all dofs on a point (free and constrained)
742: void updatePointAll(const point_type& p, const value_type v[], const int orientation = 1) {
743: const index_type& idx = this->_atlas->restrictPoint(p)[0];
744: const int dim = std::abs(idx.prefix);
745: value_type *a = &(this->_array[idx.index]);
747: if (orientation >= 0) {
748: for(int i = 0; i < dim; ++i) {
749: a[i] = v[i];
750: }
751: } else {
752: const int last = dim-1;
754: for(int i = 0; i < dim; ++i) {
755: a[i] = v[last-i];
756: }
757: }
758: };
759: public:
760: void view(const std::string& name, MPI_Comm comm = MPI_COMM_NULL) const {
761: ostringstream txt;
762: int rank;
764: if (comm == MPI_COMM_NULL) {
765: comm = this->comm();
766: rank = this->commRank();
767: } else {
768: MPI_Comm_rank(comm, &rank);
769: }
770: if(rank == 0) {
771: txt << "viewing Section " << this->getName() << std::endl;
772: if (name != "") {
773: txt << ": " << name << "'";
774: }
775: txt << std::endl;
776: }
777: const typename atlas_type::chart_type& chart = this->_atlas->getChart();
778: const value_type *array = this->_array;
780: for(typename atlas_type::chart_type::const_iterator p_iter = chart.begin(); p_iter != chart.end(); ++p_iter) {
781: const point_type& point = *p_iter;
782: const index_type& idx = this->_atlas->restrictPoint(point)[0];
783: const int dim = this->getFiberDimension(point);
785: if (idx.prefix != 0) {
786: txt << "[" << this->commRank() << "]: " << point << " dim " << idx.prefix << " offset " << idx.index << " ";
787: for(int i = 0; i < dim; i++) {
788: txt << " " << array[idx.index+i];
789: }
790: txt << std::endl;
791: }
792: }
793: PetscSynchronizedPrintf(comm, txt.str().c_str());
794: PetscSynchronizedFlush(comm);
795: };
796: };
797: // GeneralSection will support BC on a subset of unknowns on a point
798: // We make a separate constraint Atlas to mark constrained dofs on a point
799: // Storage will be contiguous by node, just as in Section
800: // This allows fast restrict(p)
801: // Then update() is accomplished by skipping constrained unknowns
802: // We must eliminate restrict() since it does not correspond to the constrained system
803: // Numbering will have to be rewritten to calculate correct mappings
804: // I think we can just generate multiple tuples per point
805: template<typename Point_, typename Value_>
806: class GeneralSection : public ALE::ParallelObject {
807: public:
808: typedef Point_ point_type;
809: typedef ALE::Point index_type;
810: typedef UniformSection<point_type, index_type> atlas_type;
811: typedef Section<point_type, int> bc_type;
812: typedef typename atlas_type::chart_type chart_type;
813: typedef Value_ value_type;
814: typedef value_type * values_type;
815: protected:
816: Obj<atlas_type> _atlas;
817: Obj<atlas_type> _atlasNew;
818: values_type _array;
819: Obj<bc_type> _bc;
820: public:
821: GeneralSection(MPI_Comm comm, const int debug = 0) : ParallelObject(comm, debug) {
822: this->_atlas = new atlas_type(comm, debug);
823: this->_atlasNew = NULL;
824: this->_array = NULL;
825: this->_bc = new bc_type(comm, debug);
826: };
827: GeneralSection(const Obj<atlas_type>& atlas) : ParallelObject(atlas->comm(), atlas->debug()), _atlas(atlas), _atlasNew(NULL), _array(NULL) {
828: this->_bc = new bc_type(comm, debug);
829: };
830: virtual ~GeneralSection() {
831: if (!this->_array) {
832: delete [] this->_array;
833: this->_array = NULL;
834: }
835: };
836: public:
837: value_type *getRawArray(const int size) {
838: static value_type *array = NULL;
839: static int maxSize = 0;
841: if (size > maxSize) {
842: maxSize = size;
843: if (array) delete [] array;
844: array = new value_type[maxSize];
845: };
846: return array;
847: };
848: public: // Verifiers
849: bool hasPoint(const point_type& point) {
850: return this->_atlas->hasPoint(point);
851: };
852: public: // Accessors
853: const Obj<atlas_type>& getAtlas() const {return this->_atlas;};
854: void setAtlas(const Obj<atlas_type>& atlas) {this->_atlas = atlas;};
855: const Obj<atlas_type>& getNewAtlas() const {return this->_atlasNew;};
856: void setNewAtlas(const Obj<atlas_type>& atlas) {this->_atlasNew = atlas;};
857: const Obj<bc_type>& getBC() const {return this->_bc;};
858: void setBC(const Obj<bc_type>& bc) {this->_bc = bc;};
859: const chart_type& getChart() const {return this->_atlas->getChart();};
860: public: // BC
861: // Returns the number of constraints on a point
862: const int getConstraintDimension(const point_type& p) const {
863: if (!this->_bc->hasPoint(p)) return 0;
864: return this->_bc->getFiberDimension(p);
865: };
866: // Set the number of constraints on a point
867: void setConstraintDimension(const point_type& p, const int numConstraints) {
868: this->_bc->setFiberDimension(p, numConstraints);
869: };
870: // Return the local dofs which are constrained on a point
871: const int *getConstraintDof(const point_type& p) {
872: return this->_bc->restrictPoint(p);
873: };
874: // Set the local dofs which are constrained on a point
875: void setConstraintDof(const point_type& p, const int dofs[]) {
876: this->_bc->updatePoint(p, dofs);
877: };
878: void copyBC(const Obj<GeneralSection>& section) {
879: this->setBC(section->getBC());
880: const chart_type& chart = this->getChart();
882: for(typename chart_type::const_iterator p_iter = chart.begin(); p_iter != chart.end(); ++p_iter) {
883: if (this->getConstraintDimension(*p_iter)) {
884: this->updatePointBC(*p_iter, section->restrictPoint(*p_iter));
885: }
886: }
887: };
888: void defaultConstraintDof() {
889: const chart_type& chart = this->getChart();
890: int size = 0;
892: for(typename chart_type::const_iterator p_iter = chart.begin(); p_iter != chart.end(); ++p_iter) {
893: size = std::max(size, this->getConstraintDimension(*p_iter));
894: }
895: int *dofs = new int[size];
896: for(typename chart_type::const_iterator p_iter = chart.begin(); p_iter != chart.end(); ++p_iter) {
897: const int cDim = this->getConstraintDimension(*p_iter);
899: if (cDim) {
900: for(int d = 0; d < cDim; ++d) {
901: dofs[d] = d;
902: }
903: this->_bc->updatePoint(*p_iter, dofs);
904: }
905: }
906: };
907: public: // Sizes
908: void clear() {
909: this->_atlas->clear();
910: delete [] this->_array;
911: this->_array = NULL;
912: this->_bc->clear();
913: };
914: // Return the total number of dofs on the point (free and constrained)
915: int getFiberDimension(const point_type& p) const {
916: return this->_atlas->restrictPoint(p)->prefix;
917: };
918: int getFiberDimension(const Obj<atlas_type>& atlas, const point_type& p) const {
919: return atlas->restrictPoint(p)->prefix;
920: };
921: // Set the total number of dofs on the point (free and constrained)
922: void setFiberDimension(const point_type& p, int dim) {
923: const index_type idx(dim, -1);
924: this->_atlas->addPoint(p);
925: this->_atlas->updatePoint(p, &idx);
926: };
927: template<typename Sequence>
928: void setFiberDimension(const Obj<Sequence>& points, int dim) {
929: for(typename Sequence::iterator p_iter = points->begin(); p_iter != points->end(); ++p_iter) {
930: this->setFiberDimension(*p_iter, dim);
931: }
932: };
933: void addFiberDimension(const point_type& p, int dim) {
934: if (this->_atlas->hasPoint(p)) {
935: const index_type values(dim, 0);
936: this->_atlas->updateAddPoint(p, &values);
937: } else {
938: this->setFiberDimension(p, dim);
939: }
940: };
941: // Return the number of constrained dofs on this point
942: int getConstrainedFiberDimension(const point_type& p) const {
943: return this->getFiberDimension(p) - this->getConstraintDimension(p);
944: };
945: // Return the total number of free dofs
946: int size() const {
947: const chart_type& points = this->getChart();
948: int size = 0;
950: for(typename chart_type::iterator p_iter = points.begin(); p_iter != points.end(); ++p_iter) {
951: size += this->getConstrainedFiberDimension(*p_iter);
952: }
953: return size;
954: };
955: // Return the total number of dofs (free and constrained)
956: int sizeWithBC() const {
957: const chart_type& points = this->getChart();
958: int size = 0;
960: for(typename chart_type::iterator p_iter = points.begin(); p_iter != points.end(); ++p_iter) {
961: size += this->getFiberDimension(*p_iter);
962: }
963: return size;
964: };
965: public: // Index retrieval
966: const int& getIndex(const point_type& p) {
967: return this->_atlas->restrictPoint(p)->index;
968: };
969: void setIndex(const point_type& p, const int& index) {
970: ((typename atlas_type::value_type *) this->_atlas->restrictPoint(p))->index = index;
971: };
972: void setIndexBC(const point_type& p, const int& index) {};
973: void getIndices(const point_type& p, PetscInt indices[], PetscInt *indx, const int orientation = 1, const bool freeOnly = false) {
974: this->getIndices(p, this->getIndex(p), indices, indx, orientation, freeOnly);
975: };
976: template<typename Order_>
977: void getIndices(const point_type& p, const Obj<Order_>& order, PetscInt indices[], PetscInt *indx, const int orientation = 1, const bool freeOnly = false) {
978: this->getIndices(p, order->getIndex(p), indices, indx, orientation, freeOnly);
979: };
980: void getIndices(const point_type& p, const int start, PetscInt indices[], PetscInt *indx, const int orientation = 1, const bool freeOnly = false) {
981: const int& dim = this->getFiberDimension(p);
982: const int& cDim = this->getConstraintDimension(p);
983: const int end = start + dim;
985: if (!cDim) {
986: if (orientation >= 0) {
987: for(int i = start; i < end; ++i) {
988: indices[(*indx)++] = i;
989: }
990: } else {
991: for(int i = end-1; i >= start; --i) {
992: indices[(*indx)++] = i;
993: }
994: }
995: } else {
996: const typename bc_type::value_type *cDof = this->getConstraintDof(p);
997: int cInd = 0;
999: if (orientation >= 0) {
1000: for(int i = start, k = 0; k < dim; ++k) {
1001: if ((cInd < cDim) && (k == cDof[cInd])) {
1002: if (!freeOnly) indices[(*indx)++] = -(k+1);
1003: ++cInd;
1004: } else {
1005: indices[(*indx)++] = i++;
1006: }
1007: }
1008: } else {
1009: const int tEnd = start + this->getConstrainedFiberDimension(p);
1011: for(int i = tEnd-1, k = 0; k < dim; ++k) {
1012: if ((cInd < cDim) && (k == cDof[cInd])) {
1013: if (!freeOnly) indices[(*indx)++] = -(dim-k+1);
1014: ++cInd;
1015: } else {
1016: indices[(*indx)++] = i--;
1017: }
1018: }
1019: }
1020: }
1021: };
1022: public: // Allocation
1023: void allocateStorage() {
1024: const int totalSize = this->sizeWithBC();
1026: this->_array = new value_type[totalSize];
1027: PetscMemzero(this->_array, totalSize * sizeof(value_type));
1028: this->_bc->allocatePoint();
1029: };
1030: void replaceStorage(value_type *newArray) {
1031: delete [] this->_array;
1032: this->_array = newArray;
1033: this->_atlas = this->_atlasNew;
1034: this->_atlasNew = NULL;
1035: };
1036: void addPoint(const point_type& point, const int dim) {
1037: if (dim == 0) return;
1038: if (this->_atlasNew.isNull()) {
1039: this->_atlasNew = new atlas_type(this->comm(), this->debug());
1040: this->_atlasNew->copy(this->_atlas);
1041: }
1042: const index_type idx(dim, -1);
1043: this->_atlasNew->addPoint(point);
1044: this->_atlasNew->updatePoint(point, &idx);
1045: };
1046: void orderPoints(const Obj<atlas_type>& atlas){
1047: const chart_type& chart = this->getChart();
1048: int offset = 0;
1050: for(typename chart_type::iterator c_iter = chart.begin(); c_iter != chart.end(); ++c_iter) {
1051: typename atlas_type::value_type idx = atlas->restrictPoint(*c_iter)[0];
1052: const int& dim = idx.prefix;
1054: idx.index = offset;
1055: atlas->updatePoint(*c_iter, &idx);
1056: offset += dim;
1057: }
1058: };
1059: void allocatePoint() {
1060: this->orderPoints(this->_atlas);
1061: this->allocateStorage();
1062: };
1063: public: // Restriction and Update
1064: // Zero entries
1065: void zero() {
1066: //memset(this->_array, 0, this->sizeWithBC()* sizeof(value_type));
1067: const chart_type& chart = this->getChart();
1069: for(typename chart_type::iterator c_iter = chart.begin(); c_iter != chart.end(); ++c_iter) {
1070: value_type *array = (value_type *) this->restrictPoint(*c_iter);
1071: const int& dim = this->getFiberDimension(*c_iter);
1072: const int& cDim = this->getConstraintDimension(*c_iter);
1074: if (!cDim) {
1075: memset(array, 0, dim * sizeof(value_type));
1076: } else {
1077: const typename bc_type::value_type *cDof = this->getConstraintDof(*c_iter);
1078: int cInd = 0;
1080: for(int i = 0; i < dim; ++i) {
1081: if ((cInd < cDim) && (i == cDof[cInd])) {++cInd; continue;}
1082: array[i] = 0.0;
1083: }
1084: }
1085: }
1086: };
1087: // Return the free values on a point
1088: const value_type *restrict() const {
1089: return this->_array;
1090: };
1091: // Return the free values on a point
1092: const value_type *restrictPoint(const point_type& p) const {
1093: return &(this->_array[this->_atlas->restrictPoint(p)[0].index]);
1094: };
1095: // Update the free values on a point
1096: void updatePoint(const point_type& p, const value_type v[], const int orientation = 1) {
1097: value_type *array = (value_type *) this->restrictPoint(p);
1098: const int& dim = this->getFiberDimension(p);
1099: const int& cDim = this->getConstraintDimension(p);
1101: if (!cDim) {
1102: if (orientation >= 0) {
1103: for(int i = 0; i < dim; ++i) {
1104: array[i] = v[i];
1105: }
1106: } else {
1107: const int last = dim-1;
1109: for(int i = 0; i < dim; ++i) {
1110: array[i] = v[last-i];
1111: }
1112: }
1113: } else {
1114: const typename bc_type::value_type *cDof = this->getConstraintDof(p);
1115: int cInd = 0;
1117: if (orientation >= 0) {
1118: for(int i = 0, k = -1; i < dim; ++i) {
1119: if ((cInd < cDim) && (i == cDof[cInd])) {++cInd; continue;}
1120: array[i] = v[++k];
1121: }
1122: } else {
1123: const int tDim = this->getConstrainedFiberDimension(p);
1125: for(int i = 0, k = 0; i < dim; ++i) {
1126: if ((cInd < cDim) && (i == cDof[cInd])) {++cInd; continue;}
1127: array[i] = v[tDim-k];
1128: ++k;
1129: }
1130: }
1131: }
1132: };
1133: // Update the free values on a point
1134: void updateAddPoint(const point_type& p, const value_type v[], const int orientation = 1) {
1135: value_type *array = (value_type *) this->restrictPoint(p);
1136: const int& dim = this->getFiberDimension(p);
1137: const int& cDim = this->getConstraintDimension(p);
1139: if (!cDim) {
1140: if (orientation >= 0) {
1141: for(int i = 0; i < dim; ++i) {
1142: array[i] += v[i];
1143: }
1144: } else {
1145: const int last = dim-1;
1147: for(int i = 0; i < dim; ++i) {
1148: array[i] += v[last-i];
1149: }
1150: }
1151: } else {
1152: const typename bc_type::value_type *cDof = this->getConstraintDof(p);
1153: int cInd = 0;
1155: if (orientation >= 0) {
1156: for(int i = 0, k = -1; i < dim; ++i) {
1157: if ((cInd < cDim) && (i == cDof[cInd])) {++cInd; continue;}
1158: array[i] += v[++k];
1159: }
1160: } else {
1161: const int tDim = this->getConstrainedFiberDimension(p);
1163: for(int i = 0, k = 0; i < dim; ++i) {
1164: if ((cInd < cDim) && (i == cDof[cInd])) {++cInd; continue;}
1165: array[i] += v[tDim-k];
1166: ++k;
1167: }
1168: }
1169: }
1170: };
1171: // Update only the constrained dofs on a point
1172: void updatePointBC(const point_type& p, const value_type v[], const int orientation = 1) {
1173: value_type *array = (value_type *) this->restrictPoint(p);
1174: const int& dim = this->getFiberDimension(p);
1175: const int& cDim = this->getConstraintDimension(p);
1177: if (cDim) {
1178: const typename bc_type::value_type *cDof = this->getConstraintDof(p);
1179: int cInd = 0;
1181: for(int i = 0, k = 0; i < dim; ++i) {
1182: if (cInd == cDim) break;
1183: if (i == cDof[cInd]) {
1184: array[i] = v[k];
1185: ++cInd;
1186: ++k;
1187: }
1188: }
1189: }
1190: };
1191: // Update all dofs on a point (free and constrained)
1192: void updatePointAll(const point_type& p, const value_type v[], const int orientation = 1) {
1193: value_type *array = (value_type *) this->restrictPoint(p);
1194: const int& dim = this->getFiberDimension(p);
1196: if (orientation >= 0) {
1197: for(int i = 0; i < dim; ++i) {
1198: array[i] = v[i];
1199: }
1200: } else {
1201: const int last = dim-1;
1203: for(int i = 0; i < dim; ++i) {
1204: array[i] = v[last-i];
1205: }
1206: }
1207: };
1208: public:
1209: void view(const std::string& name, MPI_Comm comm = MPI_COMM_NULL) const {
1210: ostringstream txt;
1211: int rank;
1213: if (comm == MPI_COMM_NULL) {
1214: comm = this->comm();
1215: rank = this->commRank();
1216: } else {
1217: MPI_Comm_rank(comm, &rank);
1218: }
1219: if (name == "") {
1220: if(rank == 0) {
1221: txt << "viewing a GeneralSection" << std::endl;
1222: }
1223: } else {
1224: if(rank == 0) {
1225: txt << "viewing GeneralSection '" << name << "'" << std::endl;
1226: }
1227: }
1228: const chart_type& chart = this->getChart();
1230: for(typename chart_type::const_iterator p_iter = chart.begin(); p_iter != chart.end(); ++p_iter) {
1231: const value_type *array = this->restrictPoint(*p_iter);
1232: const int& dim = this->getFiberDimension(*p_iter);
1234: if (dim != 0) {
1235: txt << "[" << this->commRank() << "]: " << *p_iter << " dim " << dim << " offset " << this->_atlas->restrictPoint(*p_iter)->index << " ";
1236: for(int i = 0; i < dim; i++) {
1237: txt << " " << array[i];
1238: }
1239: const int& dim = this->getConstraintDimension(*p_iter);
1241: if (dim) {
1242: const typename bc_type::value_type *bcArray = this->_bc->restrictPoint(*p_iter);
1244: txt << " constrained";
1245: for(int i = 0; i < dim; ++i) {
1246: txt << " " << bcArray[i];
1247: }
1248: }
1249: txt << std::endl;
1250: }
1251: }
1252: PetscSynchronizedPrintf(comm, txt.str().c_str());
1253: PetscSynchronizedFlush(comm);
1254: };
1255: };
1256: // A Field combines several sections
1257: template<typename Overlap_, typename Patch_, typename Section_>
1258: class Field : public ALE::ParallelObject {
1259: public:
1260: typedef Overlap_ overlap_type;
1261: typedef Patch_ patch_type;
1262: typedef Section_ section_type;
1263: typedef typename section_type::point_type point_type;
1264: typedef typename section_type::chart_type chart_type;
1265: typedef typename section_type::value_type value_type;
1266: typedef std::map<patch_type, Obj<section_type> > sheaf_type;
1267: typedef enum {SEND, RECEIVE} request_type;
1268: typedef std::map<patch_type, MPI_Request> requests_type;
1269: protected:
1270: sheaf_type _sheaf;
1271: int _tag;
1272: MPI_Datatype _datatype;
1273: requests_type _requests;
1274: public:
1275: Field(MPI_Comm comm, const int debug = 0) : ParallelObject(comm, debug) {
1276: this->_tag = this->getNewTag();
1277: this->_datatype = this->getMPIDatatype();
1278: };
1279: Field(MPI_Comm comm, const int tag, const int debug) : ParallelObject(comm, debug), _tag(tag) {
1280: this->_datatype = this->getMPIDatatype();
1281: };
1282: virtual ~Field() {};
1283: protected:
1284: MPI_Datatype getMPIDatatype() {
1285: if (sizeof(value_type) == 4) {
1286: return MPI_INT;
1287: } else if (sizeof(value_type) == 8) {
1288: return MPI_DOUBLE;
1289: } else if (sizeof(value_type) == 28) {
1290: int blen[2];
1291: MPI_Aint indices[2];
1292: MPI_Datatype oldtypes[2], newtype;
1293: blen[0] = 1; indices[0] = 0; oldtypes[0] = MPI_INT;
1294: blen[1] = 3; indices[1] = sizeof(int); oldtypes[1] = MPI_DOUBLE;
1295: MPI_Type_struct(2, blen, indices, oldtypes, &newtype);
1296: MPI_Type_commit(&newtype);
1297: return newtype;
1298: } else if (sizeof(value_type) == 32) {
1299: int blen[2];
1300: MPI_Aint indices[2];
1301: MPI_Datatype oldtypes[2], newtype;
1302: blen[0] = 1; indices[0] = 0; oldtypes[0] = MPI_DOUBLE;
1303: blen[1] = 3; indices[1] = sizeof(int); oldtypes[1] = MPI_DOUBLE;
1304: MPI_Type_struct(2, blen, indices, oldtypes, &newtype);
1305: MPI_Type_commit(&newtype);
1306: return newtype;
1307: }
1308: throw ALE::Exception("Cannot determine MPI type for value type");
1309: };
1310: int getNewTag() {
1311: static int tagKeyval = MPI_KEYVAL_INVALID;
1312: int *tagvalp = NULL, *maxval, flg;
1314: if (tagKeyval == MPI_KEYVAL_INVALID) {
1315: tagvalp = (int *) malloc(sizeof(int));
1316: MPI_Keyval_create(MPI_NULL_COPY_FN, Mesh_DelTag, &tagKeyval, (void *) NULL);
1317: MPI_Attr_put(this->_comm, tagKeyval, tagvalp);
1318: tagvalp[0] = 0;
1319: }
1320: MPI_Attr_get(this->_comm, tagKeyval, (void **) &tagvalp, &flg);
1321: if (tagvalp[0] < 1) {
1322: MPI_Attr_get(MPI_COMM_WORLD, MPI_TAG_UB, (void **) &maxval, &flg);
1323: tagvalp[0] = *maxval - 128; // hope that any still active tags were issued right at the beginning of the run
1324: }
1325: if (this->debug()) {
1326: std::cout << "[" << this->commRank() << "]Got new tag " << tagvalp[0] << std::endl;
1327: }
1328: return tagvalp[0]--;
1329: };
1330: public: // Verifiers
1331: void checkPatch(const patch_type& patch) const {
1332: if (this->_sheaf.find(patch) == this->_sheaf.end()) {
1333: ostringstream msg;
1334: msg << "Invalid field patch " << patch << std::endl;
1335: throw ALE::Exception(msg.str().c_str());
1336: }
1337: };
1338: bool hasPatch(const patch_type& patch) {
1339: if (this->_sheaf.find(patch) == this->_sheaf.end()) {
1340: return false;
1341: }
1342: return true;
1343: };
1344: public: // Accessors
1345: int getTag() const {return this->_tag;};
1346: void setTag(const int tag) {this->_tag = tag;};
1347: Obj<section_type>& getSection(const patch_type& patch) {
1348: if (this->_sheaf.find(patch) == this->_sheaf.end()) {
1349: this->_sheaf[patch] = new section_type(this->comm(), this->debug());
1350: }
1351: return this->_sheaf[patch];
1352: };
1353: void setSection(const patch_type& patch, const Obj<section_type>& section) {this->_sheaf[patch] = section;};
1354: const sheaf_type& getPatches() {
1355: return this->_sheaf;
1356: };
1357: void clear() {
1358: for(typename sheaf_type::const_iterator p_iter = this->_sheaf.begin(); p_iter != this->_sheaf.end(); ++p_iter) {
1359: p_iter->second->clear();
1360: }
1361: };
1362: public: // Adapter
1363: template<typename Topology_>
1364: void setTopology(const Obj<Topology_>& topology) {
1365: const typename Topology_::sheaf_type& patches = topology->getPatches();
1367: for(typename Topology_::sheaf_type::const_iterator p_iter = patches.begin(); p_iter != patches.end(); ++p_iter) {
1368: int rank = p_iter->first;
1369: const Obj<section_type>& section = this->getSection(rank);
1370: const Obj<typename Topology_::sieve_type::baseSequence>& base = p_iter->second->base();
1372: for(typename Topology_::sieve_type::baseSequence::iterator b_iter = base->begin(); b_iter != base->end(); ++b_iter) {
1373: section->setFiberDimension(*b_iter, 1);
1374: }
1375: }
1376: };
1377: void allocate() {
1378: for(typename sheaf_type::const_iterator p_iter = this->_sheaf.begin(); p_iter != this->_sheaf.end(); ++p_iter) {
1379: p_iter->second->allocatePoint();
1380: }
1381: };
1382: public: // Communication
1383: void construct(const int size) {
1384: const sheaf_type& patches = this->getPatches();
1386: for(typename sheaf_type::const_iterator p_iter = patches.begin(); p_iter != patches.end(); ++p_iter) {
1387: const patch_type rank = p_iter->first;
1388: const Obj<section_type>& section = this->getSection(rank);
1389: const chart_type& chart = section->getChart();
1391: for(typename chart_type::iterator c_iter = chart.begin(); c_iter != chart.end(); ++c_iter) {
1392: section->setFiberDimension(*c_iter, size);
1393: }
1394: }
1395: };
1396: template<typename Sizer>
1397: void construct(const Obj<Sizer>& sizer) {
1398: const sheaf_type& patches = this->getPatches();
1400: for(typename sheaf_type::const_iterator p_iter = patches.begin(); p_iter != patches.end(); ++p_iter) {
1401: const patch_type rank = p_iter->first;
1402: const Obj<section_type>& section = this->getSection(rank);
1403: const chart_type& chart = section->getChart();
1404:
1405: for(typename chart_type::iterator c_iter = chart.begin(); c_iter != chart.end(); ++c_iter) {
1406: section->setFiberDimension(*c_iter, *(sizer->getSection(rank)->restrictPoint(*c_iter)));
1407: }
1408: }
1409: };
1410: void constructCommunication(const request_type& requestType) {
1411: const sheaf_type& patches = this->getPatches();
1413: for(typename sheaf_type::const_iterator p_iter = patches.begin(); p_iter != patches.end(); ++p_iter) {
1414: const patch_type patch = p_iter->first;
1415: const Obj<section_type>& section = this->getSection(patch);
1416: MPI_Request request;
1418: if (requestType == RECEIVE) {
1419: if (this->_debug) {std::cout <<"["<<this->commRank()<<"] Receiving data(" << section->size() << ") from " << patch << " tag " << this->_tag << std::endl;}
1420: MPI_Recv_init((void *) section->restrict(), section->size(), this->_datatype, patch, this->_tag, this->comm(), &request);
1421: } else {
1422: if (this->_debug) {std::cout <<"["<<this->commRank()<<"] Sending data (" << section->size() << ") to " << patch << " tag " << this->_tag << std::endl;}
1423: MPI_Send_init((void *) section->restrict(), section->size(), this->_datatype, patch, this->_tag, this->comm(), &request);
1424: }
1425: this->_requests[patch] = request;
1426: }
1427: };
1428: void startCommunication() {
1429: const sheaf_type& patches = this->getPatches();
1431: for(typename sheaf_type::const_iterator p_iter = patches.begin(); p_iter != patches.end(); ++p_iter) {
1432: MPI_Request request = this->_requests[p_iter->first];
1434: MPI_Start(&request);
1435: }
1436: };
1437: void endCommunication() {
1438: const sheaf_type& patches = this->getPatches();
1439: MPI_Status status;
1441: for(typename sheaf_type::const_iterator p_iter = patches.begin(); p_iter != patches.end(); ++p_iter) {
1442: MPI_Request request = this->_requests[p_iter->first];
1444: MPI_Wait(&request, &status);
1445: }
1446: };
1447: public:
1448: void view(const std::string& name, MPI_Comm comm = MPI_COMM_NULL) const {
1449: ostringstream txt;
1450: int rank;
1452: if (comm == MPI_COMM_NULL) {
1453: comm = this->comm();
1454: rank = this->commRank();
1455: } else {
1456: MPI_Comm_rank(comm, &rank);
1457: }
1458: if (name == "") {
1459: if(rank == 0) {
1460: txt << "viewing a Field" << std::endl;
1461: }
1462: } else {
1463: if(rank == 0) {
1464: txt << "viewing Field '" << name << "'" << std::endl;
1465: }
1466: }
1467: PetscSynchronizedPrintf(comm, txt.str().c_str());
1468: PetscSynchronizedFlush(comm);
1469: for(typename sheaf_type::const_iterator p_iter = this->_sheaf.begin(); p_iter != this->_sheaf.end(); ++p_iter) {
1470: ostringstream txt1;
1472: txt1 << "[" << this->commRank() << "]: Patch " << p_iter->first << std::endl;
1473: PetscSynchronizedPrintf(comm, txt1.str().c_str());
1474: PetscSynchronizedFlush(comm);
1475: p_iter->second->view("field section", comm);
1476: }
1477: };
1478: };
1479: }
1481: namespace ALECompat {
1482: namespace New {
1483: using ALE::Obj;
1484: // A ConstantSection is the simplest Section
1485: // All fibers are dimension 1
1486: // All values are equal to a constant
1487: // We need no value storage and no communication for completion
1488: template<typename Topology_, typename Value_>
1489: class NewConstantSection : public ALE::ParallelObject {
1490: public:
1491: typedef Topology_ topology_type;
1492: typedef typename topology_type::patch_type patch_type;
1493: typedef typename topology_type::sieve_type sieve_type;
1494: typedef typename topology_type::point_type point_type;
1495: typedef std::set<point_type> chart_type;
1496: typedef std::map<patch_type, chart_type> atlas_type;
1497: typedef Value_ value_type;
1498: protected:
1499: Obj<topology_type> _topology;
1500: atlas_type _atlas;
1501: chart_type _emptyChart;
1502: value_type _value;
1503: value_type _defaultValue;
1504: public:
1505: NewConstantSection(MPI_Comm comm, const int debug = 0) : ParallelObject(comm, debug), _defaultValue(0) {
1506: this->_topology = new topology_type(comm, debug);
1507: };
1508: NewConstantSection(const Obj<topology_type>& topology) : ParallelObject(topology->comm(), topology->debug()), _topology(topology) {};
1509: NewConstantSection(const Obj<topology_type>& topology, const value_type& value) : ParallelObject(topology->comm(), topology->debug()), _topology(topology), _value(value), _defaultValue(value) {};
1510: NewConstantSection(const Obj<topology_type>& topology, const value_type& value, const value_type& defaultValue) : ParallelObject(topology->comm(), topology->debug()), _topology(topology), _value(value), _defaultValue(defaultValue) {};
1511: public: // Verifiers
1512: void checkPatch(const patch_type& patch) const {
1513: this->_topology->checkPatch(patch);
1514: if (this->_atlas.find(patch) == this->_atlas.end()) {
1515: ostringstream msg;
1516: msg << "Invalid atlas patch " << patch << std::endl;
1517: throw ALE::Exception(msg.str().c_str());
1518: }
1519: };
1520: void checkPoint(const patch_type& patch, const point_type& point) const {
1521: this->checkPatch(patch);
1522: if (this->_atlas.find(patch)->second.find(point) == this->_atlas.find(patch)->second.end()) {
1523: ostringstream msg;
1524: msg << "Invalid section point " << point << std::endl;
1525: throw ALE::Exception(msg.str().c_str());
1526: }
1527: };
1528: void checkDimension(const int& dim) {
1529: if (dim != 1) {
1530: ostringstream msg;
1531: msg << "Invalid fiber dimension " << dim << " must be 1" << std::endl;
1532: throw ALE::Exception(msg.str().c_str());
1533: }
1534: };
1535: bool hasPatch(const patch_type& patch) {
1536: if (this->_atlas.find(patch) == this->_atlas.end()) {
1537: return false;
1538: }
1539: return true;
1540: };
1541: bool hasPoint(const patch_type& patch, const point_type& point) const {
1542: this->checkPatch(patch);
1543: return this->_atlas.find(patch)->second.count(point) > 0;
1544: };
1545: bool hasPoint(const point_type& point) const {
1546: this->checkPatch(0);
1547: return this->_atlas.find(0)->second.count(point) > 0;
1548: };
1549: public: // Accessors
1550: const Obj<topology_type>& getTopology() const {return this->_topology;};
1551: void setTopology(const Obj<topology_type>& topology) {this->_topology = topology;};
1552: const chart_type& getPatch(const patch_type& patch) {
1553: if (this->hasPatch(patch)) {
1554: return this->_atlas[patch];
1555: }
1556: return this->_emptyChart;
1557: };
1558: void updatePatch(const patch_type& patch, const point_type& point) {
1559: this->_atlas[patch].insert(point);
1560: };
1561: template<typename Points>
1562: void updatePatch(const patch_type& patch, const Obj<Points>& points) {
1563: this->_atlas[patch].insert(points->begin(), points->end());
1564: };
1565: value_type getDefaultValue() {return this->_defaultValue;};
1566: void setDefaultValue(const value_type value) {this->_defaultValue = value;};
1567: public: // Sizes
1568: void clear() {
1569: this->_atlas.clear();
1570: };
1571: int getFiberDimension(const patch_type& patch, const point_type& p) const {
1572: if (this->hasPoint(patch, p)) return 1;
1573: return 0;
1574: };
1575: void setFiberDimension(const patch_type& patch, const point_type& p, int dim) {
1576: this->checkDimension(dim);
1577: this->updatePatch(patch, p);
1578: };
1579: template<typename Sequence>
1580: void setFiberDimension(const patch_type& patch, const Obj<Sequence>& points, int dim) {
1581: for(typename topology_type::label_sequence::iterator p_iter = points->begin(); p_iter != points->end(); ++p_iter) {
1582: this->setFiberDimension(patch, *p_iter, dim);
1583: }
1584: };
1585: void addFiberDimension(const patch_type& patch, const point_type& p, int dim) {
1586: if (this->hasPatch(patch) && (this->_atlas[patch].find(p) != this->_atlas[patch].end())) {
1587: ostringstream msg;
1588: msg << "Invalid addition to fiber dimension " << dim << " cannot exceed 1" << std::endl;
1589: throw ALE::Exception(msg.str().c_str());
1590: } else {
1591: this->setFiberDimension(patch, p, dim);
1592: }
1593: };
1594: void setFiberDimensionByDepth(const patch_type& patch, int depth, int dim) {
1595: this->setFiberDimension(patch, this->_topology->getLabelStratum(patch, "depth", depth), dim);
1596: };
1597: void setFiberDimensionByHeight(const patch_type& patch, int height, int dim) {
1598: this->setFiberDimension(patch, this->_topology->getLabelStratum(patch, "height", height), dim);
1599: };
1600: int size(const patch_type& patch) {return this->_atlas[patch].size();};
1601: int size(const patch_type& patch, const point_type& p) {return this->getFiberDimension(patch, p);};
1602: public: // Restriction
1603: const value_type *restrict(const patch_type& patch, const point_type& p) const {
1604: //std::cout <<"["<<this->commRank()<<"]: Constant restrict ("<<patch<<","<<p<<") from " << std::endl;
1605: //for(typename chart_type::iterator c_iter = this->_atlas.find(patch)->second.begin(); c_iter != this->_atlas.find(patch)->second.end(); ++c_iter) {
1606: // std::cout <<"["<<this->commRank()<<"]: point " << *c_iter << std::endl;
1607: //}
1608: if (this->hasPoint(patch, p)) {
1609: return &this->_value;
1610: }
1611: return &this->_defaultValue;
1612: };
1613: const value_type *restrictPoint(const patch_type& patch, const point_type& p) const {return this->restrict(patch, p);};
1614: const value_type *restrictPoint(const point_type& p) const {return this->restrict(0, p);};
1615: void update(const patch_type& patch, const point_type& p, const value_type v[]) {
1616: this->checkPatch(patch);
1617: this->_value = v[0];
1618: };
1619: void updatePoint(const patch_type& patch, const point_type& p, const value_type v[]) {return this->update(patch, p, v);};
1620: void updateAdd(const patch_type& patch, const point_type& p, const value_type v[]) {
1621: this->checkPatch(patch);
1622: this->_value += v[0];
1623: };
1624: void updateAddPoint(const patch_type& patch, const point_type& p, const value_type v[]) {return this->updateAdd(patch, p, v);};
1625: public:
1626: void copy(const Obj<NewConstantSection>& section) {
1627: const typename topology_type::sheaf_type& patches = this->_topology->getPatches();
1629: for(typename topology_type::sheaf_type::const_iterator p_iter = patches.begin(); p_iter != patches.end(); ++p_iter) {
1630: const patch_type& patch = p_iter->first;
1631: if (!section->hasPatch(patch)) continue;
1632: const chart_type& chart = section->getPatch(patch);
1634: for(typename chart_type::iterator c_iter = chart.begin(); c_iter != chart.end(); ++c_iter) {
1635: this->updatePatch(patch, *c_iter);
1636: }
1637: this->_value = section->restrict(patch, *chart.begin())[0];
1638: }
1639: };
1640: void view(const std::string& name, MPI_Comm comm = MPI_COMM_NULL) const {
1641: ostringstream txt;
1642: int rank;
1644: if (comm == MPI_COMM_NULL) {
1645: comm = this->comm();
1646: rank = this->commRank();
1647: } else {
1648: MPI_Comm_rank(comm, &rank);
1649: }
1650: if (name == "") {
1651: if(rank == 0) {
1652: txt << "viewing a NewConstantSection" << std::endl;
1653: }
1654: } else {
1655: if(rank == 0) {
1656: txt << "viewing NewConstantSection '" << name << "'" << std::endl;
1657: }
1658: }
1659: const typename topology_type::sheaf_type& sheaf = this->_topology->getPatches();
1661: for(typename topology_type::sheaf_type::const_iterator p_iter = sheaf.begin(); p_iter != sheaf.end(); ++p_iter) {
1662: txt <<"["<<this->commRank()<<"]: Patch " << p_iter->first << std::endl;
1663: txt <<"["<<this->commRank()<<"]: Value " << this->_value << std::endl;
1664: }
1665: PetscSynchronizedPrintf(comm, txt.str().c_str());
1666: PetscSynchronizedFlush(comm);
1667: };
1668: };
1670: // A UniformSection often acts as an Atlas
1671: // All fibers are the same dimension
1672: // Note we can use a ConstantSection for this Atlas
1673: // Each point may have a different vector
1674: // Thus we need storage for values, and hence must implement completion
1675: template<typename Topology_, typename Value_, int fiberDim = 1>
1676: class UniformSection : public ALE::ParallelObject {
1677: public:
1678: typedef Topology_ topology_type;
1679: typedef typename topology_type::patch_type patch_type;
1680: typedef typename topology_type::sieve_type sieve_type;
1681: typedef typename topology_type::point_type point_type;
1682: typedef NewConstantSection<topology_type, int> atlas_type;
1683: typedef typename atlas_type::chart_type chart_type;
1684: typedef Value_ value_type;
1685: typedef struct {value_type v[fiberDim];} fiber_type;
1686: typedef std::map<point_type, fiber_type> array_type;
1687: typedef std::map<patch_type, array_type> values_type;
1688: protected:
1689: Obj<atlas_type> _atlas;
1690: values_type _arrays;
1691: public:
1692: UniformSection(MPI_Comm comm, const int debug = 0) : ParallelObject(comm, debug) {
1693: this->_atlas = new atlas_type(comm, debug);
1694: };
1695: UniformSection(const Obj<topology_type>& topology) : ParallelObject(topology->comm(), topology->debug()) {
1696: this->_atlas = new atlas_type(topology, fiberDim, 0);
1697: };
1698: UniformSection(const Obj<atlas_type>& atlas) : ParallelObject(atlas->comm(), atlas->debug()), _atlas(atlas) {};
1699: protected:
1700: value_type *getRawArray(const int size) {
1701: static value_type *array = NULL;
1702: static int maxSize = 0;
1704: if (size > maxSize) {
1705: maxSize = size;
1706: if (array) delete [] array;
1707: array = new value_type[maxSize];
1708: };
1709: return array;
1710: };
1711: public: // Verifiers
1712: void checkPatch(const patch_type& patch) {
1713: this->_atlas->checkPatch(patch);
1714: if (this->_arrays.find(patch) == this->_arrays.end()) {
1715: ostringstream msg;
1716: msg << "Invalid section patch: " << patch << std::endl;
1717: throw ALE::Exception(msg.str().c_str());
1718: }
1719: };
1720: bool hasPatch(const patch_type& patch) {
1721: return this->_atlas->hasPatch(patch);
1722: };
1723: bool hasPoint(const patch_type& patch, const point_type& point) {
1724: return this->_atlas->hasPoint(patch, point);
1725: };
1726: bool hasPoint(const point_type& point) {
1727: return this->_atlas->hasPoint(0, point);
1728: };
1729: void checkDimension(const int& dim) {
1730: if (dim != fiberDim) {
1731: ostringstream msg;
1732: msg << "Invalid fiber dimension " << dim << " must be " << fiberDim << std::endl;
1733: throw ALE::Exception(msg.str().c_str());
1734: }
1735: };
1736: public: // Accessors
1737: const Obj<atlas_type>& getAtlas() {return this->_atlas;};
1738: void setAtlas(const Obj<atlas_type>& atlas) {this->_atlas = atlas;};
1739: const Obj<topology_type>& getTopology() {return this->_atlas->getTopology();};
1740: void setTopology(const Obj<topology_type>& topology) {this->_atlas->setTopology(topology);};
1741: const chart_type& getPatch(const patch_type& patch) {
1742: return this->_atlas->getPatch(patch);
1743: };
1744: void updatePatch(const patch_type& patch, const point_type& point) {
1745: this->setFiberDimension(patch, point, 1);
1746: };
1747: template<typename Points>
1748: void updatePatch(const patch_type& patch, const Obj<Points>& points) {
1749: for(typename Points::iterator p_iter = points->begin(); p_iter != points->end(); ++p_iter) {
1750: this->setFiberDimension(patch, *p_iter, 1);
1751: }
1752: };
1753: void copy(const Obj<UniformSection<Topology_, Value_, fiberDim> >& section) {
1754: this->getAtlas()->copy(section->getAtlas());
1755: const typename topology_type::sheaf_type& sheaf = section->getTopology()->getPatches();
1757: for(typename topology_type::sheaf_type::const_iterator s_iter = sheaf.begin(); s_iter != sheaf.end(); ++s_iter) {
1758: const patch_type& patch = s_iter->first;
1759: if (!section->hasPatch(patch)) continue;
1760: const chart_type& chart = section->getPatch(patch);
1762: for(typename chart_type::const_iterator c_iter = chart.begin(); c_iter != chart.end(); ++c_iter) {
1763: this->updatePoint(s_iter->first, *c_iter, section->restrictPoint(s_iter->first, *c_iter));
1764: }
1765: }
1766: };
1767: public: // Sizes
1768: void clear() {
1769: this->_atlas->clear();
1770: this->_arrays.clear();
1771: };
1772: int getFiberDimension(const patch_type& patch, const point_type& p) const {
1773: // Could check for non-existence here
1774: return this->_atlas->restrictPoint(patch, p)[0];
1775: };
1776: void setFiberDimension(const patch_type& patch, const point_type& p, int dim) {
1777: this->checkDimension(dim);
1778: this->_atlas->updatePatch(patch, p);
1779: this->_atlas->updatePoint(patch, p, &dim);
1780: };
1781: template<typename Sequence>
1782: void setFiberDimension(const patch_type& patch, const Obj<Sequence>& points, int dim) {
1783: for(typename Sequence::iterator p_iter = points->begin(); p_iter != points->end(); ++p_iter) {
1784: this->setFiberDimension(patch, *p_iter, dim);
1785: }
1786: };
1787: void setFiberDimension(const patch_type& patch, const std::set<point_type>& points, int dim) {
1788: for(typename std::set<point_type>::iterator p_iter = points.begin(); p_iter != points.end(); ++p_iter) {
1789: this->setFiberDimension(patch, *p_iter, dim);
1790: }
1791: };
1792: void addFiberDimension(const patch_type& patch, const point_type& p, int dim) {
1793: if (this->hasPatch(patch) && (this->_atlas[patch].find(p) != this->_atlas[patch].end())) {
1794: ostringstream msg;
1795: msg << "Invalid addition to fiber dimension " << dim << " cannot exceed " << fiberDim << std::endl;
1796: throw ALE::Exception(msg.str().c_str());
1797: } else {
1798: this->setFiberDimension(patch, p, dim);
1799: }
1800: };
1801: void setFiberDimensionByDepth(const patch_type& patch, int depth, int dim) {
1802: this->setFiberDimension(patch, this->getTopology()->getLabelStratum(patch, "depth", depth), dim);
1803: };
1804: void setFiberDimensionByHeight(const patch_type& patch, int height, int dim) {
1805: this->setFiberDimension(patch, this->getTopology()->getLabelStratum(patch, "height", height), dim);
1806: };
1807: int size(const patch_type& patch) {
1808: const typename atlas_type::chart_type& points = this->_atlas->getPatch(patch);
1809: int size = 0;
1811: for(typename atlas_type::chart_type::iterator p_iter = points.begin(); p_iter != points.end(); ++p_iter) {
1812: size += this->getFiberDimension(patch, *p_iter);
1813: }
1814: return size;
1815: };
1816: int size(const patch_type& patch, const point_type& p) {
1817: const typename atlas_type::chart_type& points = this->_atlas->getPatch(patch);
1818: const Obj<typename sieve_type::coneSet> closure = this->getTopology()->getPatch(patch)->closure(p);
1819: typename sieve_type::coneSet::iterator end = closure->end();
1820: int size = 0;
1822: for(typename sieve_type::coneSet::iterator c_iter = closure->begin(); c_iter != end; ++c_iter) {
1823: if (points.count(*c_iter)) {
1824: size += this->getFiberDimension(patch, *c_iter);
1825: }
1826: }
1827: return size;
1828: };
1829: void orderPatches() {};
1830: public: // Restriction
1831: // Return a pointer to the entire contiguous storage array
1832: const array_type& restrict(const patch_type& patch) {
1833: this->checkPatch(patch);
1834: return this->_arrays[patch];
1835: };
1836: // Return the values for the closure of this point
1837: // use a smart pointer?
1838: const value_type *restrict(const patch_type& patch, const point_type& p) {
1839: this->checkPatch(patch);
1840: const chart_type& chart = this->getPatch(patch);
1841: array_type& array = this->_arrays[patch];
1842: const int size = this->size(patch, p);
1843: value_type *values = this->getRawArray(size);
1844: int j = -1;
1846: // We could actually ask for the height of the individual point
1847: if (this->getTopology()->height(patch) < 2) {
1848: // Avoids only the copy of closure()
1849: const int& dim = this->_atlas->restrictPoint(patch, p)[0];
1851: if (chart.count(p)) {
1852: for(int i = 0; i < dim; ++i) {
1853: values[++j] = array[p].v[i];
1854: }
1855: }
1856: // Need only the cone
1857: const Obj<typename sieve_type::coneSequence>& cone = this->getTopology()->getPatch(patch)->cone(p);
1858: typename sieve_type::coneSequence::iterator end = cone->end();
1860: for(typename sieve_type::coneSequence::iterator p_iter = cone->begin(); p_iter != end; ++p_iter) {
1861: if (chart.count(*p_iter)) {
1862: const int& dim = this->_atlas->restrictPoint(patch, *p_iter)[0];
1864: for(int i = 0; i < dim; ++i) {
1865: values[++j] = array[*p_iter].v[i];
1866: }
1867: }
1868: }
1869: } else {
1870: // Right now, we have no way of consistently ordering the closure()
1871: const Obj<typename sieve_type::coneSet>& closure = this->getTopology()->getPatch(patch)->closure(p);
1872: typename sieve_type::coneSet::iterator end = closure->end();
1874: for(typename sieve_type::coneSet::iterator p_iter = closure->begin(); p_iter != end; ++p_iter) {
1875: if (chart.count(*p_iter)) {
1876: const int& dim = this->_atlas->restrictPoint(patch, *p_iter)[0];
1878: for(int i = 0; i < dim; ++i) {
1879: values[++j] = array[*p_iter].v[i];
1880: }
1881: }
1882: }
1883: }
1884: if (j != size-1) {
1885: ostringstream txt;
1887: txt << "Invalid restrict to point " << p << std::endl;
1888: txt << " j " << j << " should be " << (size-1) << std::endl;
1889: std::cout << txt.str();
1890: throw ALE::Exception(txt.str().c_str());
1891: }
1892: return values;
1893: };
1894: void update(const patch_type& patch, const point_type& p, const value_type v[]) {
1895: this->_atlas->checkPatch(patch);
1896: const chart_type& chart = this->getPatch(patch);
1897: array_type& array = this->_arrays[patch];
1898: int j = -1;
1900: if (this->getTopology()->height(patch) < 2) {
1901: // Only avoids the copy of closure()
1902: const int& dim = this->_atlas->restrict(patch, p)[0];
1904: if (chart.count(p)) {
1905: for(int i = 0; i < dim; ++i) {
1906: array[p].v[i] = v[++j];
1907: }
1908: }
1909: // Should be closure()
1910: const Obj<typename sieve_type::coneSequence>& cone = this->getTopology()->getPatch(patch)->cone(p);
1912: for(typename sieve_type::coneSequence::iterator p_iter = cone->begin(); p_iter != cone->end(); ++p_iter) {
1913: if (chart.count(*p_iter)) {
1914: const int& dim = this->_atlas->restrict(patch, *p_iter)[0];
1916: for(int i = 0; i < dim; ++i) {
1917: array[*p_iter].v[i] = v[++j];
1918: }
1919: }
1920: }
1921: } else {
1922: throw ALE::Exception("Update is not yet implemented for interpolated sieves");
1923: }
1924: };
1925: void updateAdd(const patch_type& patch, const point_type& p, const value_type v[]) {
1926: this->_atlas->checkPatch(patch);
1927: const chart_type& chart = this->getPatch(patch);
1928: array_type& array = this->_arrays[patch];
1929: int j = -1;
1931: if (this->getTopology()->height(patch) < 2) {
1932: // Only avoids the copy of closure()
1933: const int& dim = this->_atlas->restrict(patch, p)[0];
1935: if (chart.count(p)) {
1936: for(int i = 0; i < dim; ++i) {
1937: array[p].v[i] += v[++j];
1938: }
1939: }
1940: // Should be closure()
1941: const Obj<typename sieve_type::coneSequence>& cone = this->getTopology()->getPatch(patch)->cone(p);
1943: for(typename sieve_type::coneSequence::iterator p_iter = cone->begin(); p_iter != cone->end(); ++p_iter) {
1944: if (chart.count(*p_iter)) {
1945: const int& dim = this->_atlas->restrict(patch, *p_iter)[0];
1947: for(int i = 0; i < dim; ++i) {
1948: array[*p_iter].v[i] += v[++j];
1949: }
1950: }
1951: }
1952: } else {
1953: throw ALE::Exception("Not yet implemented for interpolated sieves");
1954: }
1955: };
1956: // Return only the values associated to this point, not its closure
1957: const value_type *restrictPoint(const patch_type& patch, const point_type& p) {
1958: this->checkPatch(patch);
1959: return this->_arrays[patch][p].v;
1960: };
1961: const value_type *restrictPoint(const point_type& p) {
1962: this->checkPatch(0);
1963: return this->_arrays[0][p].v;
1964: };
1965: // Update only the values associated to this point, not its closure
1966: void updatePoint(const patch_type& patch, const point_type& p, const value_type v[]) {
1967: this->_atlas->checkPatch(patch);
1968: for(int i = 0; i < fiberDim; ++i) {
1969: this->_arrays[patch][p].v[i] = v[i];
1970: }
1971: };
1972: // Update only the values associated to this point, not its closure
1973: void updateAddPoint(const patch_type& patch, const point_type& p, const value_type v[]) {
1974: this->_atlas->checkPatch(patch);
1975: for(int i = 0; i < fiberDim; ++i) {
1976: this->_arrays[patch][p].v[i] += v[i];
1977: }
1978: };
1979: public:
1980: void view(const std::string& name, MPI_Comm comm = MPI_COMM_NULL) {
1981: ostringstream txt;
1982: int rank;
1984: if (comm == MPI_COMM_NULL) {
1985: comm = this->comm();
1986: rank = this->commRank();
1987: } else {
1988: MPI_Comm_rank(comm, &rank);
1989: }
1990: if (name == "") {
1991: if(rank == 0) {
1992: txt << "viewing a UniformSection" << std::endl;
1993: }
1994: } else {
1995: if(rank == 0) {
1996: txt << "viewing UniformSection '" << name << "'" << std::endl;
1997: }
1998: }
1999: for(typename values_type::const_iterator a_iter = this->_arrays.begin(); a_iter != this->_arrays.end(); ++a_iter) {
2000: const patch_type& patch = a_iter->first;
2001: array_type& array = this->_arrays[patch];
2003: txt << "[" << this->commRank() << "]: Patch " << patch << std::endl;
2004: const typename atlas_type::chart_type& chart = this->_atlas->getPatch(patch);
2006: for(typename atlas_type::chart_type::const_iterator p_iter = chart.begin(); p_iter != chart.end(); ++p_iter) {
2007: const point_type& point = *p_iter;
2008: const typename atlas_type::value_type dim = this->_atlas->restrict(patch, point)[0];
2010: if (dim != 0) {
2011: txt << "[" << this->commRank() << "]: " << point << " dim " << dim << " ";
2012: for(int i = 0; i < dim; i++) {
2013: txt << " " << array[point].v[i];
2014: }
2015: txt << std::endl;
2016: }
2017: }
2018: }
2019: PetscSynchronizedPrintf(comm, txt.str().c_str());
2020: PetscSynchronizedFlush(comm);
2021: };
2022: };
2024: // A Section is our most general construct (more general ones could be envisioned)
2025: // The Atlas is a UniformSection of dimension 1 and value type Point
2026: // to hold each fiber dimension and offsets into a contiguous patch array
2027: template<typename Topology_, typename Value_>
2028: class Section : public ALE::ParallelObject {
2029: public:
2030: typedef Topology_ topology_type;
2031: typedef typename topology_type::patch_type patch_type;
2032: typedef typename topology_type::sieve_type sieve_type;
2033: typedef typename topology_type::point_type point_type;
2034: typedef ALE::Point index_type;
2035: typedef UniformSection<topology_type, index_type> atlas_type;
2036: typedef typename atlas_type::chart_type chart_type;
2037: typedef Value_ value_type;
2038: typedef value_type * array_type;
2039: typedef std::map<patch_type, array_type> values_type;
2040: typedef std::vector<index_type> IndexArray;
2041: protected:
2042: Obj<atlas_type> _atlas;
2043: Obj<atlas_type> _atlasNew;
2044: values_type _arrays;
2045: Obj<IndexArray> _indexArray;
2046: public:
2047: Section(MPI_Comm comm, const int debug = 0) : ParallelObject(comm, debug) {
2048: this->_atlas = new atlas_type(comm, debug);
2049: this->_atlasNew = NULL;
2050: this->_indexArray = new IndexArray();
2051: };
2052: Section(const Obj<topology_type>& topology) : ParallelObject(topology->comm(), topology->debug()), _atlasNew(NULL) {
2053: this->_atlas = new atlas_type(topology);
2054: this->_indexArray = new IndexArray();
2055: };
2056: Section(const Obj<atlas_type>& atlas) : ParallelObject(atlas->comm(), atlas->debug()), _atlas(atlas), _atlasNew(NULL) {
2057: this->_indexArray = new IndexArray();
2058: };
2059: virtual ~Section() {
2060: for(typename values_type::iterator a_iter = this->_arrays.begin(); a_iter != this->_arrays.end(); ++a_iter) {
2061: delete [] a_iter->second;
2062: a_iter->second = NULL;
2063: }
2064: };
2065: protected:
2066: value_type *getRawArray(const int size) {
2067: static value_type *array = NULL;
2068: static int maxSize = 0;
2070: if (size > maxSize) {
2071: maxSize = size;
2072: if (array) delete [] array;
2073: array = new value_type[maxSize];
2074: };
2075: return array;
2076: };
2077: public: // Verifiers
2078: void checkPatch(const patch_type& patch) {
2079: this->_atlas->checkPatch(patch);
2080: if (this->_arrays.find(patch) == this->_arrays.end()) {
2081: ostringstream msg;
2082: msg << "Invalid section patch: " << patch << std::endl;
2083: throw ALE::Exception(msg.str().c_str());
2084: }
2085: };
2086: bool hasPatch(const patch_type& patch) {
2087: return this->_atlas->hasPatch(patch);
2088: };
2089: public: // Accessors
2090: const Obj<atlas_type>& getAtlas() {return this->_atlas;};
2091: void setAtlas(const Obj<atlas_type>& atlas) {this->_atlas = atlas;};
2092: const Obj<topology_type>& getTopology() {return this->_atlas->getTopology();};
2093: void setTopology(const Obj<topology_type>& topology) {this->_atlas->setTopology(topology);};
2094: const chart_type& getPatch(const patch_type& patch) {
2095: return this->_atlas->getPatch(patch);
2096: };
2097: bool hasPoint(const patch_type& patch, const point_type& point) {
2098: return this->_atlas->hasPoint(patch, point);
2099: };
2100: public: // Sizes
2101: void clear() {
2102: this->_atlas->clear();
2103: this->_arrays.clear();
2104: };
2105: int getFiberDimension(const patch_type& patch, const point_type& p) const {
2106: // Could check for non-existence here
2107: return this->_atlas->restrictPoint(patch, p)->prefix;
2108: };
2109: int getFiberDimension(const Obj<atlas_type>& atlas, const patch_type& patch, const point_type& p) const {
2110: // Could check for non-existence here
2111: return atlas->restrictPoint(patch, p)->prefix;
2112: };
2113: void setFiberDimension(const patch_type& patch, const point_type& p, int dim) {
2114: const index_type idx(dim, -1);
2115: this->_atlas->updatePatch(patch, p);
2116: this->_atlas->updatePoint(patch, p, &idx);
2117: };
2118: template<typename Sequence>
2119: void setFiberDimension(const patch_type& patch, const Obj<Sequence>& points, int dim) {
2120: for(typename topology_type::label_sequence::iterator p_iter = points->begin(); p_iter != points->end(); ++p_iter) {
2121: this->setFiberDimension(patch, *p_iter, dim);
2122: }
2123: };
2124: void addFiberDimension(const patch_type& patch, const point_type& p, int dim) {
2125: if (this->_atlas->hasPatch(patch) && this->_atlas->hasPoint(patch, p)) {
2126: const index_type values(dim, 0);
2127: this->_atlas->updateAddPoint(patch, p, &values);
2128: } else {
2129: this->setFiberDimension(patch, p, dim);
2130: }
2131: };
2132: void setFiberDimensionByDepth(const patch_type& patch, int depth, int dim) {
2133: this->setFiberDimension(patch, this->getTopology()->getLabelStratum(patch, "depth", depth), dim);
2134: };
2135: void setFiberDimensionByHeight(const patch_type& patch, int height, int dim) {
2136: this->setFiberDimension(patch, this->getTopology()->getLabelStratum(patch, "height", height), dim);
2137: };
2138: int size(const patch_type& patch) {
2139: const typename atlas_type::chart_type& points = this->_atlas->getPatch(patch);
2140: int size = 0;
2142: for(typename atlas_type::chart_type::iterator p_iter = points.begin(); p_iter != points.end(); ++p_iter) {
2143: size += std::max(0, this->getFiberDimension(patch, *p_iter));
2144: }
2145: return size;
2146: };
2147: int sizeWithBC(const patch_type& patch) {
2148: const typename atlas_type::chart_type& points = this->_atlas->getPatch(patch);
2149: int size = 0;
2151: for(typename atlas_type::chart_type::iterator p_iter = points.begin(); p_iter != points.end(); ++p_iter) {
2152: size += std::abs(this->getFiberDimension(patch, *p_iter));
2153: }
2154: return size;
2155: };
2156: int size(const patch_type& patch, const point_type& p) {
2157: if (this->getTopology()->depth() > 1) throw ALE::Exception("Compatibility layer is not for interpolated meshes");
2158: const typename atlas_type::chart_type& points = this->_atlas->getPatch(patch);
2159: const Obj<typename sieve_type::coneSequence> closure = this->getTopology()->getPatch(patch)->cone(p);
2160: typename sieve_type::coneSequence::iterator end = closure->end();
2161: int size = 0;
2163: size += std::max(0, this->getFiberDimension(patch, p));
2164: for(typename sieve_type::coneSequence::iterator c_iter = closure->begin(); c_iter != end; ++c_iter) {
2165: if (points.count(*c_iter)) {
2166: size += std::max(0, this->getFiberDimension(patch, *c_iter));
2167: }
2168: }
2169: return size;
2170: };
2171: int sizeWithBC(const patch_type& patch, const point_type& p) {
2172: if (this->getTopology()->depth() > 1) throw ALE::Exception("Compatibility layer is not for interpolated meshes");
2173: const typename atlas_type::chart_type& points = this->_atlas->getPatch(patch);
2174: const Obj<typename sieve_type::coneSequence> closure = this->getTopology()->getPatch(patch)->cone(p);
2175: typename sieve_type::coneSequence::iterator end = closure->end();
2176: int size = 0;
2178: size += std::abs(this->getFiberDimension(patch, p));
2179: for(typename sieve_type::coneSequence::iterator c_iter = closure->begin(); c_iter != end; ++c_iter) {
2180: if (points.count(*c_iter)) {
2181: size += std::abs(this->getFiberDimension(patch, *c_iter));
2182: }
2183: }
2184: return size;
2185: };
2186: int size(const Obj<atlas_type>& atlas, const patch_type& patch) {
2187: const typename atlas_type::chart_type& points = atlas->getPatch(patch);
2188: int size = 0;
2190: for(typename atlas_type::chart_type::iterator p_iter = points.begin(); p_iter != points.end(); ++p_iter) {
2191: size += std::max(0, this->getFiberDimension(atlas, patch, *p_iter));
2192: }
2193: return size;
2194: };
2195: public: // Index retrieval
2196: const index_type& getIndex(const patch_type& patch, const point_type& p) {
2197: this->checkPatch(patch);
2198: return this->_atlas->restrictPoint(patch, p)[0];
2199: };
2200: template<typename Numbering>
2201: const index_type getIndex(const patch_type& patch, const point_type& p, const Obj<Numbering>& numbering) {
2202: this->checkPatch(patch);
2203: return index_type(this->getFiberDimension(patch, p), numbering->getIndex(p));
2204: };
2205: const Obj<IndexArray>& getIndices(const patch_type& patch, const point_type& p, const int level = -1) {
2206: this->_indexArray->clear();
2208: if (level == 0) {
2209: this->_indexArray->push_back(this->getIndex(patch, p));
2210: } else if ((level == 1) || (this->getTopology()->height(patch) == 1)) {
2211: const Obj<typename sieve_type::coneSequence>& cone = this->getTopology()->getPatch(patch)->cone(p);
2212: typename sieve_type::coneSequence::iterator end = cone->end();
2214: this->_indexArray->push_back(this->getIndex(patch, p));
2215: for(typename sieve_type::coneSequence::iterator p_iter = cone->begin(); p_iter != end; ++p_iter) {
2216: this->_indexArray->push_back(this->getIndex(patch, *p_iter));
2217: }
2218: } else if (level == -1) {
2219: #if 1
2220: throw ALE::Exception("Call should be moved to Bundle");
2221: #else
2222: const Obj<typename sieve_type::coneSet> closure = this->getTopology()->getPatch(patch)->closure(p);
2223: typename sieve_type::coneSet::iterator end = closure->end();
2225: for(typename sieve_type::coneSet::iterator p_iter = closure->begin(); p_iter != end; ++p_iter) {
2226: this->_indexArray->push_back(this->getIndex(patch, *p_iter));
2227: }
2228: #endif
2229: } else {
2230: const Obj<typename sieve_type::coneArray> cone = this->getTopology()->getPatch(patch)->nCone(p, level);
2231: typename sieve_type::coneArray::iterator end = cone->end();
2233: for(typename sieve_type::coneArray::iterator p_iter = cone->begin(); p_iter != end; ++p_iter) {
2234: this->_indexArray->push_back(this->getIndex(patch, *p_iter));
2235: }
2236: }
2237: return this->_indexArray;
2238: };
2239: template<typename Numbering>
2240: const Obj<IndexArray>& getIndices(const patch_type& patch, const point_type& p, const Obj<Numbering>& numbering, const int level = -1) {
2241: this->_indexArray->clear();
2243: if (level == 0) {
2244: this->_indexArray->push_back(this->getIndex(patch, p, numbering));
2245: } else if ((level == 1) || (this->getTopology()->height(patch) == 1)) {
2246: const Obj<typename sieve_type::coneSequence>& cone = this->getTopology()->getPatch(patch)->cone(p);
2247: typename sieve_type::coneSequence::iterator end = cone->end();
2249: this->_indexArray->push_back(this->getIndex(patch, p, numbering));
2250: for(typename sieve_type::coneSequence::iterator p_iter = cone->begin(); p_iter != end; ++p_iter) {
2251: this->_indexArray->push_back(this->getIndex(patch, *p_iter, numbering));
2252: }
2253: } else if (level == -1) {
2254: #if 1
2255: throw ALE::Exception("Call should be moved to Bundle");
2256: #else
2257: const Obj<typename sieve_type::coneSet> closure = this->getTopology()->getPatch(patch)->closure(p);
2258: typename sieve_type::coneSet::iterator end = closure->end();
2260: for(typename sieve_type::coneSet::iterator p_iter = closure->begin(); p_iter != end; ++p_iter) {
2261: this->_indexArray->push_back(this->getIndex(patch, *p_iter, numbering));
2262: }
2263: #endif
2264: } else {
2265: const Obj<typename sieve_type::coneArray> cone = this->getTopology()->getPatch(patch)->nCone(p, level);
2266: typename sieve_type::coneArray::iterator end = cone->end();
2268: for(typename sieve_type::coneArray::iterator p_iter = cone->begin(); p_iter != end; ++p_iter) {
2269: this->_indexArray->push_back(this->getIndex(patch, *p_iter, numbering));
2270: }
2271: }
2272: return this->_indexArray;
2273: };
2274: public: // Allocation
2275: void orderPoint(const Obj<atlas_type>& atlas, const Obj<sieve_type>& sieve, const patch_type& patch, const point_type& point, int& offset, int& bcOffset, const bool postponeGhosts = false) {
2276: const Obj<typename sieve_type::coneSequence>& cone = sieve->cone(point);
2277: typename sieve_type::coneSequence::iterator end = cone->end();
2278: index_type idx = atlas->restrictPoint(patch, point)[0];
2279: const int& dim = idx.prefix;
2280: const index_type defaultIdx(0, -1);
2282: if (atlas->getPatch(patch).count(point) == 0) {
2283: idx = defaultIdx;
2284: }
2285: if (idx.index == -1) {
2286: for(typename sieve_type::coneSequence::iterator c_iter = cone->begin(); c_iter != end; ++c_iter) {
2287: if (this->_debug > 1) {std::cout << " Recursing to " << *c_iter << std::endl;}
2288: this->orderPoint(atlas, sieve, patch, *c_iter, offset, bcOffset);
2289: }
2290: if (dim > 0) {
2291: bool number = true;
2293: // Maybe use template specialization here
2294: if (postponeGhosts && this->getTopology()->getSendOverlap()->capContains(point)) {
2295: const Obj<typename topology_type::send_overlap_type::supportSequence>& ranks = this->getTopology()->getSendOverlap()->support(point);
2297: for(typename topology_type::send_overlap_type::supportSequence::iterator r_iter = ranks->begin(); r_iter != ranks->end(); ++r_iter) {
2298: if (this->commRank() > *r_iter) {
2299: number = false;
2300: break;
2301: }
2302: }
2303: }
2304: if (number) {
2305: if (this->_debug > 1) {std::cout << " Ordering point " << point << " at " << offset << std::endl;}
2306: idx.index = offset;
2307: atlas->updatePoint(patch, point, &idx);
2308: offset += dim;
2309: } else {
2310: if (this->_debug > 1) {std::cout << " Ignoring ghost point " << point << std::endl;}
2311: }
2312: } else if (dim < 0) {
2313: if (this->_debug > 1) {std::cout << " Ordering boundary point " << point << " at " << bcOffset << std::endl;}
2314: idx.index = bcOffset;
2315: atlas->updatePoint(patch, point, &idx);
2316: bcOffset += dim;
2317: }
2318: }
2319: }
2320: void orderPatch(const Obj<atlas_type>& atlas, const patch_type& patch, int& offset, int& bcOffset, const bool postponeGhosts = false) {
2321: const typename atlas_type::chart_type& chart = atlas->getPatch(patch);
2323: for(typename atlas_type::chart_type::const_iterator p_iter = chart.begin(); p_iter != chart.end(); ++p_iter) {
2324: if (this->_debug > 1) {std::cout << "Ordering closure of point " << *p_iter << std::endl;}
2325: this->orderPoint(atlas, this->getTopology()->getPatch(patch), patch, *p_iter, offset, bcOffset, postponeGhosts);
2326: }
2327: for(typename atlas_type::chart_type::const_iterator p_iter = chart.begin(); p_iter != chart.end(); ++p_iter) {
2328: index_type idx = atlas->restrictPoint(patch, *p_iter)[0];
2329: const int& dim = idx.prefix;
2331: if (dim < 0) {
2332: if (this->_debug > 1) {std::cout << "Correcting boundary offset of point " << *p_iter << std::endl;}
2333: idx.index = offset - (idx.index+2);
2334: atlas->updatePoint(patch, *p_iter, &idx);
2335: }
2336: }
2337: };
2338: void orderPatches(const Obj<atlas_type>& atlas, const bool postponeGhosts = false) {
2339: const typename topology_type::sheaf_type& patches = this->getTopology()->getPatches();
2341: for(typename topology_type::sheaf_type::const_iterator p_iter = patches.begin(); p_iter != patches.end(); ++p_iter) {
2342: if (this->_debug > 1) {std::cout << "Ordering patch " << p_iter->first << std::endl;}
2343: int offset = 0, bcOffset = -2;
2345: if (!atlas->hasPatch(p_iter->first)) continue;
2346: this->orderPatch(atlas, p_iter->first, offset, bcOffset, postponeGhosts);
2347: }
2348: };
2349: void orderPatches(const bool postponeGhosts = false) {
2350: this->orderPatches(this->_atlas, postponeGhosts);
2351: };
2352: void allocateStorage() {
2353: const typename topology_type::sheaf_type& patches = this->getTopology()->getPatches();
2355: for(typename topology_type::sheaf_type::const_iterator p_iter = patches.begin(); p_iter != patches.end(); ++p_iter) {
2356: if (!this->_atlas->hasPatch(p_iter->first)) continue;
2357: this->_arrays[p_iter->first] = new value_type[this->sizeWithBC(p_iter->first)];
2358: PetscMemzero(this->_arrays[p_iter->first], this->sizeWithBC(p_iter->first) * sizeof(value_type));
2359: }
2360: };
2361: void allocate(const bool postponeGhosts = false) {
2362: bool doGhosts = false;
2364: if (postponeGhosts && !this->getTopology()->getSendOverlap().isNull()) {
2365: doGhosts = true;
2366: }
2367: this->orderPatches(doGhosts);
2368: if (doGhosts) {
2369: const typename topology_type::sheaf_type& patches = this->getTopology()->getPatches();
2371: for(typename topology_type::sheaf_type::const_iterator p_iter = patches.begin(); p_iter != patches.end(); ++p_iter) {
2372: if (this->_debug > 1) {std::cout << "Ordering patch " << p_iter->first << " for ghosts" << std::endl;}
2373: const typename atlas_type::chart_type& points = this->_atlas->getPatch(p_iter->first);
2374: int offset = 0, bcOffset = -2;
2376: for(typename atlas_type::chart_type::iterator point = points.begin(); point != points.end(); ++point) {
2377: const index_type& idx = this->_atlas->restrictPoint(p_iter->first, *point)[0];
2379: offset = std::max(offset, idx.index + std::abs(idx.prefix));
2380: }
2381: if (!this->_atlas->hasPatch(p_iter->first)) continue;
2382: this->orderPatch(this->_atlas, p_iter->first, offset, bcOffset);
2383: if (offset != this->sizeWithBC(p_iter->first)) throw ALE::Exception("Inconsistent array sizes in section");
2384: }
2385: }
2386: this->allocateStorage();
2387: };
2388: void addPoint(const patch_type& patch, const point_type& point, const int dim) {
2389: if (dim == 0) return;
2390: //const typename atlas_type::chart_type& chart = this->_atlas->getPatch(patch);
2392: //if (chart.find(point) == chart.end()) {
2393: if (this->_atlasNew.isNull()) {
2394: this->_atlasNew = new atlas_type(this->getTopology());
2395: this->_atlasNew->copy(this->_atlas);
2396: }
2397: const index_type idx(dim, -1);
2398: this->_atlasNew->updatePatch(patch, point);
2399: this->_atlasNew->updatePoint(patch, point, &idx);
2400: };
2401: void reallocate() {
2402: if (this->_atlasNew.isNull()) return;
2403: const typename topology_type::sheaf_type& patches = this->getTopology()->getPatches();
2405: // Since copy() preserves offsets, we must reinitialize them before ordering
2406: for(typename topology_type::sheaf_type::const_iterator p_iter = patches.begin(); p_iter != patches.end(); ++p_iter) {
2407: const patch_type& patch = p_iter->first;
2408: const typename atlas_type::chart_type& chart = this->_atlasNew->getPatch(patch);
2409: index_type defaultIdx(0, -1);
2411: for(typename atlas_type::chart_type::const_iterator c_iter = chart.begin(); c_iter != chart.end(); ++c_iter) {
2412: defaultIdx.prefix = this->_atlasNew->restrictPoint(patch, *c_iter)[0].prefix;
2413: this->_atlasNew->updatePoint(patch, *c_iter, &defaultIdx);
2414: }
2415: }
2416: this->orderPatches(this->_atlasNew);
2417: // Copy over existing values
2418: for(typename topology_type::sheaf_type::const_iterator p_iter = patches.begin(); p_iter != patches.end(); ++p_iter) {
2419: const patch_type& patch = p_iter->first;
2420: value_type *newArray = new value_type[this->size(this->_atlasNew, patch)];
2422: if (!this->_atlas->hasPatch(patch)) {
2423: this->_arrays[patch] = newArray;
2424: continue;
2425: }
2426: const typename atlas_type::chart_type& chart = this->_atlas->getPatch(patch);
2427: const value_type *array = this->_arrays[patch];
2429: for(typename atlas_type::chart_type::const_iterator c_iter = chart.begin(); c_iter != chart.end(); ++c_iter) {
2430: const index_type& idx = this->_atlas->restrictPoint(patch, *c_iter)[0];
2431: const int size = idx.prefix;
2432: const int offset = idx.index;
2433: const int& newOffset = this->_atlasNew->restrictPoint(patch, *c_iter)[0].index;
2435: for(int i = 0; i < size; ++i) {
2436: newArray[newOffset+i] = array[offset+i];
2437: }
2438: }
2439: delete [] this->_arrays[patch];
2440: this->_arrays[patch] = newArray;
2441: }
2442: this->_atlas = this->_atlasNew;
2443: this->_atlasNew = NULL;
2444: };
2445: public: // Restriction and Update
2446: // Zero entries
2447: void zero(const patch_type& patch) {
2448: this->checkPatch(patch);
2449: memset(this->_arrays[patch], 0, this->size(patch)* sizeof(value_type));
2450: };
2451: // Return a pointer to the entire contiguous storage array
2452: const value_type *restrict(const patch_type& patch) {
2453: this->checkPatch(patch);
2454: return this->_arrays[patch];
2455: };
2456: // Update the entire contiguous storage array
2457: void update(const patch_type& patch, const value_type v[]) {
2458: const value_type *array = this->_arrays[patch];
2459: const int size = this->size(patch);
2461: for(int i = 0; i < size; i++) {
2462: array[i] = v[i];
2463: }
2464: };
2465: // Return the values for the closure of this point
2466: // use a smart pointer?
2467: const value_type *restrict(const patch_type& patch, const point_type& p) {
2468: this->checkPatch(patch);
2469: const value_type *a = this->_arrays[patch];
2470: const int size = this->sizeWithBC(patch, p);
2471: value_type *values = this->getRawArray(size);
2472: int j = -1;
2474: if (this->getTopology()->height(patch) < 2) {
2475: // Avoids the copy of both
2476: // points in topology->closure()
2477: // indices in _atlas->restrict()
2478: const index_type& pInd = this->_atlas->restrictPoint(patch, p)[0];
2480: for(int i = pInd.index; i < std::abs(pInd.prefix) + pInd.index; ++i) {
2481: values[++j] = a[i];
2482: }
2483: const Obj<typename sieve_type::coneSequence>& cone = this->getTopology()->getPatch(patch)->cone(p);
2484: typename sieve_type::coneSequence::iterator end = cone->end();
2486: for(typename sieve_type::coneSequence::iterator p_iter = cone->begin(); p_iter != end; ++p_iter) {
2487: const index_type& ind = this->_atlas->restrictPoint(patch, *p_iter)[0];
2488: const int& start = ind.index;
2489: const int& length = std::abs(ind.prefix);
2491: for(int i = start; i < start + length; ++i) {
2492: values[++j] = a[i];
2493: }
2494: }
2495: } else {
2496: const Obj<IndexArray>& ind = this->getIndices(patch, p);
2498: for(typename IndexArray::iterator i_iter = ind->begin(); i_iter != ind->end(); ++i_iter) {
2499: const int& start = i_iter->index;
2500: const int& length = std::abs(i_iter->prefix);
2502: for(int i = start; i < start + length; ++i) {
2503: values[++j] = a[i];
2504: }
2505: }
2506: }
2507: if (j != size-1) {
2508: ostringstream txt;
2510: txt << "Invalid restrict to point " << p << std::endl;
2511: txt << " j " << j << " should be " << (size-1) << std::endl;
2512: std::cout << txt.str();
2513: throw ALE::Exception(txt.str().c_str());
2514: }
2515: return values;
2516: };
2517: // Update the values for the closure of this point
2518: void update(const patch_type& patch, const point_type& p, const value_type v[]) {
2519: this->checkPatch(patch);
2520: value_type *a = this->_arrays[patch];
2521: int j = -1;
2523: if (this->getTopology()->height(patch) < 2) {
2524: // Avoids the copy of both
2525: // points in topology->closure()
2526: // indices in _atlas->restrict()
2527: const index_type& pInd = this->_atlas->restrictPoint(patch, p)[0];
2529: for(int i = pInd.index; i < pInd.prefix + pInd.index; ++i) {
2530: a[i] = v[++j];
2531: }
2532: j += std::max(0, -pInd.prefix);
2533: const Obj<typename sieve_type::coneSequence>& cone = this->getTopology()->getPatch(patch)->cone(p);
2534: typename sieve_type::coneSequence::iterator end = cone->end();
2536: for(typename sieve_type::coneSequence::iterator p_iter = cone->begin(); p_iter != end; ++p_iter) {
2537: const index_type& ind = this->_atlas->restrictPoint(patch, *p_iter)[0];
2538: const int& start = ind.index;
2539: const int& length = ind.prefix;
2541: for(int i = start; i < start + length; ++i) {
2542: a[i] = v[++j];
2543: }
2544: j += std::max(0, -length);
2545: }
2546: } else {
2547: const Obj<IndexArray>& ind = this->getIndices(patch, p);
2549: for(typename IndexArray::iterator i_iter = ind->begin(); i_iter != ind->end(); ++i_iter) {
2550: const int& start = i_iter->index;
2551: const int& length = i_iter->prefix;
2553: for(int i = start; i < start + length; ++i) {
2554: a[i] = v[++j];
2555: }
2556: j += std::max(0, -length);
2557: }
2558: }
2559: };
2560: // Update the values for the closure of this point
2561: void updateAdd(const patch_type& patch, const point_type& p, const value_type v[]) {
2562: this->checkPatch(patch);
2563: value_type *a = this->_arrays[patch];
2564: int j = -1;
2566: if (this->getTopology()->height(patch) < 2) {
2567: // Avoids the copy of both
2568: // points in topology->closure()
2569: // indices in _atlas->restrict()
2570: const index_type& pInd = this->_atlas->restrictPoint(patch, p)[0];
2572: for(int i = pInd.index; i < pInd.prefix + pInd.index; ++i) {
2573: a[i] += v[++j];
2574: }
2575: j += std::max(0, -pInd.prefix);
2576: const Obj<typename sieve_type::coneSequence>& cone = this->getTopology()->getPatch(patch)->cone(p);
2577: typename sieve_type::coneSequence::iterator end = cone->end();
2579: for(typename sieve_type::coneSequence::iterator p_iter = cone->begin(); p_iter != end; ++p_iter) {
2580: const index_type& ind = this->_atlas->restrictPoint(patch, *p_iter)[0];
2581: const int& start = ind.index;
2582: const int& length = ind.prefix;
2584: for(int i = start; i < start + length; ++i) {
2585: a[i] += v[++j];
2586: }
2587: j += std::max(0, -length);
2588: }
2589: } else {
2590: const Obj<IndexArray>& ind = this->getIndices(patch, p);
2592: for(typename IndexArray::iterator i_iter = ind->begin(); i_iter != ind->end(); ++i_iter) {
2593: const int& start = i_iter->index;
2594: const int& length = i_iter->prefix;
2596: for(int i = start; i < start + length; ++i) {
2597: a[i] += v[++j];
2598: }
2599: j += std::max(0, -length);
2600: }
2601: }
2602: };
2603: // Update the values for the closure of this point
2604: template<typename Input>
2605: void update(const patch_type& patch, const point_type& p, const Obj<Input>& v) {
2606: this->checkPatch(patch);
2607: value_type *a = this->_arrays[patch];
2609: if (this->getTopology()->height(patch) == 1) {
2610: // Only avoids the copy
2611: const index_type& pInd = this->_atlas->restrictPoint(patch, p)[0];
2612: typename Input::iterator v_iter = v->begin();
2613: typename Input::iterator v_end = v->end();
2615: for(int i = pInd.index; i < pInd.prefix + pInd.index; ++i) {
2616: a[i] = *v_iter;
2617: ++v_iter;
2618: }
2619: const Obj<typename sieve_type::coneSequence>& cone = this->getTopology()->getPatch(patch)->cone(p);
2620: typename sieve_type::coneSequence::iterator end = cone->end();
2622: for(typename sieve_type::coneSequence::iterator p_iter = cone->begin(); p_iter != end; ++p_iter) {
2623: const index_type& ind = this->_atlas->restrictPoint(patch, *p_iter)[0];
2624: const int& start = ind.index;
2625: const int& length = ind.prefix;
2627: for(int i = start; i < start + length; ++i) {
2628: a[i] = *v_iter;
2629: ++v_iter;
2630: }
2631: }
2632: } else {
2633: const Obj<IndexArray>& ind = this->getIndices(patch, p);
2634: typename Input::iterator v_iter = v->begin();
2635: typename Input::iterator v_end = v->end();
2637: for(typename IndexArray::iterator i_iter = ind->begin(); i_iter != ind->end(); ++i_iter) {
2638: const int& start = i_iter->index;
2639: const int& length = i_iter->prefix;
2641: for(int i = start; i < start + length; ++i) {
2642: a[i] = *v_iter;
2643: ++v_iter;
2644: }
2645: }
2646: }
2647: };
2648: void updateBC(const patch_type& patch, const point_type& p, const value_type v[]) {
2649: this->checkPatch(patch);
2650: value_type *a = this->_arrays[patch];
2651: int j = -1;
2653: if (this->getTopology()->height(patch) < 2) {
2654: // Avoids the copy of both
2655: // points in topology->closure()
2656: // indices in _atlas->restrict()
2657: const index_type& pInd = this->_atlas->restrictPoint(patch, p)[0];
2659: for(int i = pInd.index; i < std::abs(pInd.prefix) + pInd.index; ++i) {
2660: a[i] = v[++j];
2661: }
2662: const Obj<typename sieve_type::coneSequence>& cone = this->getTopology()->getPatch(patch)->cone(p);
2663: typename sieve_type::coneSequence::iterator end = cone->end();
2665: for(typename sieve_type::coneSequence::iterator p_iter = cone->begin(); p_iter != end; ++p_iter) {
2666: const index_type& ind = this->_atlas->restrictPoint(patch, *p_iter)[0];
2667: const int& start = ind.index;
2668: const int& length = std::abs(ind.prefix);
2670: for(int i = start; i < start + length; ++i) {
2671: a[i] = v[++j];
2672: }
2673: }
2674: } else {
2675: const Obj<IndexArray>& ind = this->getIndices(patch, p);
2677: for(typename IndexArray::iterator i_iter = ind->begin(); i_iter != ind->end(); ++i_iter) {
2678: const int& start = i_iter->index;
2679: const int& length = std::abs(i_iter->prefix);
2681: for(int i = start; i < start + length; ++i) {
2682: a[i] = v[++j];
2683: }
2684: }
2685: }
2686: };
2687: // Return only the values associated to this point, not its closure
2688: const value_type *restrictPoint(const patch_type& patch, const point_type& p) {
2689: this->checkPatch(patch);
2690: return &(this->_arrays[patch][this->_atlas->restrictPoint(patch, p)[0].index]);
2691: };
2692: // Update only the values associated to this point, not its closure
2693: void updatePoint(const patch_type& patch, const point_type& p, const value_type v[]) {
2694: this->checkPatch(patch);
2695: const index_type& idx = this->_atlas->restrictPoint(patch, p)[0];
2696: value_type *a = &(this->_arrays[patch][idx.index]);
2698: for(int i = 0; i < idx.prefix; ++i) {
2699: a[i] = v[i];
2700: }
2701: };
2702: // Update only the values associated to this point, not its closure
2703: void updateAddPoint(const patch_type& patch, const point_type& p, const value_type v[]) {
2704: this->checkPatch(patch);
2705: const index_type& idx = this->_atlas->restrictPoint(patch, p)[0];
2706: value_type *a = &(this->_arrays[patch][idx.index]);
2708: for(int i = 0; i < idx.prefix; ++i) {
2709: a[i] += v[i];
2710: }
2711: };
2712: void updatePointBC(const patch_type& patch, const point_type& p, const value_type v[]) {
2713: this->checkPatch(patch);
2714: const index_type& idx = this->_atlas->restrictPoint(patch, p)[0];
2715: value_type *a = &(this->_arrays[patch][idx.index]);
2717: for(int i = 0; i < std::abs(idx.prefix); ++i) {
2718: a[i] = v[i];
2719: }
2720: };
2721: public: // BC
2722: void copyBC(const Obj<Section>& section) {
2723: const typename topology_type::sheaf_type& patches = this->getTopology()->getPatches();
2725: for(typename topology_type::sheaf_type::const_iterator patch_iter = patches.begin(); patch_iter != patches.end(); ++patch_iter) {
2726: const typename atlas_type::chart_type& chart = this->_atlas->getPatch(patch_iter->first);
2728: for(typename atlas_type::chart_type::const_iterator p_iter = chart.begin(); p_iter != chart.end(); ++p_iter) {
2729: const index_type& idx = this->_atlas->restrictPoint(patch_iter->first, *p_iter)[0];
2731: if (idx.prefix < 0) {
2732: this->updatePointBC(patch_iter->first, *p_iter, section->restrictPoint(patch_iter->first, *p_iter));
2733: }
2734: }
2735: }
2736: };
2737: public:
2738: void view(const std::string& name, MPI_Comm comm = MPI_COMM_NULL) const {
2739: ostringstream txt;
2740: int rank;
2742: if (comm == MPI_COMM_NULL) {
2743: comm = this->comm();
2744: rank = this->commRank();
2745: } else {
2746: MPI_Comm_rank(comm, &rank);
2747: }
2748: if (name == "") {
2749: if(rank == 0) {
2750: txt << "viewing a Section" << std::endl;
2751: }
2752: } else {
2753: if(rank == 0) {
2754: txt << "viewing Section '" << name << "'" << std::endl;
2755: }
2756: }
2757: for(typename values_type::const_iterator a_iter = this->_arrays.begin(); a_iter != this->_arrays.end(); ++a_iter) {
2758: const patch_type& patch = a_iter->first;
2759: const value_type *array = a_iter->second;
2761: txt << "[" << this->commRank() << "]: Patch " << patch << std::endl;
2762: const typename atlas_type::chart_type& chart = this->_atlas->getPatch(patch);
2764: for(typename atlas_type::chart_type::const_iterator p_iter = chart.begin(); p_iter != chart.end(); ++p_iter) {
2765: const point_type& point = *p_iter;
2766: const index_type& idx = this->_atlas->restrictPoint(patch, point)[0];
2768: if (idx.prefix != 0) {
2769: txt << "[" << this->commRank() << "]: " << point << " dim " << idx.prefix << " offset " << idx.index << " ";
2770: for(int i = 0; i < std::abs(idx.prefix); i++) {
2771: txt << " " << array[idx.index+i];
2772: }
2773: txt << std::endl;
2774: }
2775: }
2776: }
2777: if (this->_arrays.empty()) {
2778: txt << "[" << this->commRank() << "]: empty" << std::endl;
2779: }
2780: PetscSynchronizedPrintf(comm, txt.str().c_str());
2781: PetscSynchronizedFlush(comm);
2782: };
2783: };
2785: // An Overlap is a Sifter describing the overlap of two Sieves
2786: // Each arrow is local point ---(remote point)---> remote rank right now
2787: // For XSifter, this should change to (local patch, local point) ---> (remote rank, remote patch, remote point)
2789: template<typename Topology_, typename Value_>
2790: class OldConstantSection : public ALE::ParallelObject {
2791: public:
2792: typedef OldConstantSection<Topology_, Value_> section_type;
2793: typedef Topology_ topology_type;
2794: typedef typename topology_type::patch_type patch_type;
2795: typedef typename topology_type::sieve_type sieve_type;
2796: typedef typename topology_type::point_type point_type;
2797: typedef Value_ value_type;
2798: protected:
2799: Obj<topology_type> _topology;
2800: const value_type _value;
2801: Obj<section_type> _section;
2802: public:
2803: OldConstantSection(MPI_Comm comm, const value_type value, const int debug = 0) : ParallelObject(comm, debug), _value(value) {
2804: this->_topology = new topology_type(comm, debug);
2805: this->_section = this;
2806: this->_section.addRef();
2807: };
2808: OldConstantSection(const Obj<topology_type>& topology, const value_type value) : ParallelObject(topology->comm(), topology->debug()), _topology(topology), _value(value) {
2809: this->_section = this;
2810: this->_section.addRef();
2811: };
2812: virtual ~OldConstantSection() {};
2813: public: // Verifiers
2814: bool hasPoint(const patch_type& patch, const point_type& point) const {return true;};
2815: public: // Restriction
2816: const value_type *restrict(const patch_type& patch) {return &this->_value;};
2817: const value_type *restrictPoint(const patch_type& patch, const point_type& p) {return &this->_value;};
2818: public: // Adapter
2819: const Obj<section_type>& getSection(const patch_type& patch) {return this->_section;};
2820: const value_type *restrictPoint(const point_type& p) {return &this->_value;};
2821: public:
2822: void view(const std::string& name, MPI_Comm comm = MPI_COMM_NULL) const {
2823: ostringstream txt;
2824: int rank;
2826: if (comm == MPI_COMM_NULL) {
2827: comm = this->comm();
2828: rank = this->commRank();
2829: } else {
2830: MPI_Comm_rank(comm, &rank);
2831: }
2832: if (name == "") {
2833: if(rank == 0) {
2834: txt << "viewing a OldConstantSection with value " << this->_value << std::endl;
2835: }
2836: } else {
2837: if(rank == 0) {
2838: txt << "viewing OldConstantSection '" << name << "' with value " << this->_value << std::endl;
2839: }
2840: }
2841: PetscSynchronizedPrintf(comm, txt.str().c_str());
2842: PetscSynchronizedFlush(comm);
2843: };
2844: };
2846: template<typename Overlap_, typename Topology_, typename Value_>
2847: class OverlapValues : public Section<Topology_, Value_> {
2848: public:
2849: typedef Overlap_ overlap_type;
2850: typedef Section<Topology_, Value_> base_type;
2851: typedef typename base_type::topology_type topology_type;
2852: typedef typename base_type::patch_type patch_type;
2853: typedef typename base_type::chart_type chart_type;
2854: typedef typename base_type::atlas_type atlas_type;
2855: typedef typename base_type::value_type value_type;
2856: typedef enum {SEND, RECEIVE} request_type;
2857: typedef std::map<patch_type, MPI_Request> requests_type;
2858: protected:
2859: int _tag;
2860: MPI_Datatype _datatype;
2861: requests_type _requests;
2862: public:
2863: OverlapValues(MPI_Comm comm, const int debug = 0) : Section<Topology_, Value_>(comm, debug) {
2864: this->_tag = this->getNewTag();
2865: this->_datatype = this->getMPIDatatype();
2866: };
2867: OverlapValues(MPI_Comm comm, const int tag, const int debug) : Section<Topology_, Value_>(comm, debug), _tag(tag) {
2868: this->_datatype = this->getMPIDatatype();
2869: };
2870: virtual ~OverlapValues() {};
2871: protected:
2872: MPI_Datatype getMPIDatatype() {
2873: if (sizeof(value_type) == 4) {
2874: return MPI_INT;
2875: } else if (sizeof(value_type) == 8) {
2876: return MPI_DOUBLE;
2877: } else if (sizeof(value_type) == 28) {
2878: int blen[2];
2879: MPI_Aint indices[2];
2880: MPI_Datatype oldtypes[2], newtype;
2881: blen[0] = 1; indices[0] = 0; oldtypes[0] = MPI_INT;
2882: blen[1] = 3; indices[1] = sizeof(int); oldtypes[1] = MPI_DOUBLE;
2883: MPI_Type_struct(2, blen, indices, oldtypes, &newtype);
2884: MPI_Type_commit(&newtype);
2885: return newtype;
2886: } else if (sizeof(value_type) == 32) {
2887: int blen[2];
2888: MPI_Aint indices[2];
2889: MPI_Datatype oldtypes[2], newtype;
2890: blen[0] = 1; indices[0] = 0; oldtypes[0] = MPI_DOUBLE;
2891: blen[1] = 3; indices[1] = sizeof(int); oldtypes[1] = MPI_DOUBLE;
2892: MPI_Type_struct(2, blen, indices, oldtypes, &newtype);
2893: MPI_Type_commit(&newtype);
2894: return newtype;
2895: }
2896: throw ALE::Exception("Cannot determine MPI type for value type");
2897: };
2898: int getNewTag() {
2899: static int tagKeyval = MPI_KEYVAL_INVALID;
2900: int *tagvalp = NULL, *maxval, flg;
2902: if (tagKeyval == MPI_KEYVAL_INVALID) {
2903: tagvalp = (int *) malloc(sizeof(int));
2904: MPI_Keyval_create(MPI_NULL_COPY_FN, Mesh_DelTag, &tagKeyval, (void *) NULL);
2905: MPI_Attr_put(this->_comm, tagKeyval, tagvalp);
2906: tagvalp[0] = 0;
2907: }
2908: MPI_Attr_get(this->_comm, tagKeyval, (void **) &tagvalp, &flg);
2909: if (tagvalp[0] < 1) {
2910: MPI_Attr_get(MPI_COMM_WORLD, MPI_TAG_UB, (void **) &maxval, &flg);
2911: tagvalp[0] = *maxval - 128; // hope that any still active tags were issued right at the beginning of the run
2912: }
2913: if (this->debug()) {
2914: std::cout << "[" << this->commRank() << "]Got new tag " << tagvalp[0] << std::endl;
2915: }
2916: return tagvalp[0]--;
2917: };
2918: public: // Accessors
2919: int getTag() const {return this->_tag;};
2920: void setTag(const int tag) {this->_tag = tag;};
2921: public:
2922: void construct(const int size) {
2923: const typename topology_type::sheaf_type& patches = this->getTopology()->getPatches();
2925: for(typename topology_type::sheaf_type::const_iterator p_iter = patches.begin(); p_iter != patches.end(); ++p_iter) {
2926: const Obj<typename topology_type::sieve_type::baseSequence>& base = p_iter->second->base();
2927: int rank = p_iter->first;
2929: for(typename topology_type::sieve_type::baseSequence::iterator b_iter = base->begin(); b_iter != base->end(); ++b_iter) {
2930: this->setFiberDimension(rank, *b_iter, size);
2931: }
2932: }
2933: };
2934: template<typename Sizer>
2935: void construct(const Obj<Sizer>& sizer) {
2936: const typename topology_type::sheaf_type& patches = this->getTopology()->getPatches();
2938: for(typename topology_type::sheaf_type::const_iterator p_iter = patches.begin(); p_iter != patches.end(); ++p_iter) {
2939: const Obj<typename topology_type::sieve_type::baseSequence>& base = p_iter->second->base();
2940: int rank = p_iter->first;
2942: for(typename topology_type::sieve_type::baseSequence::iterator b_iter = base->begin(); b_iter != base->end(); ++b_iter) {
2943: this->setFiberDimension(rank, *b_iter, *(sizer->restrictPoint(rank, *b_iter)));
2944: }
2945: }
2946: };
2947: void constructCommunication(const request_type& requestType) {
2948: const typename topology_type::sheaf_type& patches = this->getAtlas()->getTopology()->getPatches();
2950: for(typename topology_type::sheaf_type::const_iterator p_iter = patches.begin(); p_iter != patches.end(); ++p_iter) {
2951: const patch_type patch = p_iter->first;
2952: MPI_Request request;
2954: if (requestType == RECEIVE) {
2955: if (this->_debug) {std::cout <<"["<<this->commRank()<<"] Receiving data(" << this->size(patch) << ") from " << patch << " tag " << this->_tag << std::endl;}
2956: MPI_Recv_init(this->_arrays[patch], this->size(patch), this->_datatype, patch, this->_tag, this->_comm, &request);
2957: } else {
2958: if (this->_debug) {std::cout <<"["<<this->commRank()<<"] Sending data (" << this->size(patch) << ") to " << patch << " tag " << this->_tag << std::endl;}
2959: MPI_Send_init(this->_arrays[patch], this->size(patch), this->_datatype, patch, this->_tag, this->_comm, &request);
2960: }
2961: this->_requests[patch] = request;
2962: }
2963: };
2964: void startCommunication() {
2965: const typename topology_type::sheaf_type& patches = this->getAtlas()->getTopology()->getPatches();
2967: for(typename topology_type::sheaf_type::const_iterator p_iter = patches.begin(); p_iter != patches.end(); ++p_iter) {
2968: MPI_Request request = this->_requests[p_iter->first];
2970: MPI_Start(&request);
2971: }
2972: };
2973: void endCommunication() {
2974: const typename topology_type::sheaf_type& patches = this->getAtlas()->getTopology()->getPatches();
2975: MPI_Status status;
2977: for(typename topology_type::sheaf_type::const_iterator p_iter = patches.begin(); p_iter != patches.end(); ++p_iter) {
2978: MPI_Request request = this->_requests[p_iter->first];
2980: MPI_Wait(&request, &status);
2981: }
2982: };
2983: };
2984: }
2985: }
2986: #endif