Actual source code: Completion.hh
1: #ifndef included_ALE_Completion_hh
2: #define included_ALE_Completion_hh
4: #ifndef included_ALE_Sections_hh
5: #include <Sections.hh>
6: #endif
8: #include <iostream>
9: #include <fstream>
11: namespace ALE {
12: namespace New {
13: template<typename Bundle_, typename Value_>
14: class Completion {
15: public:
16: typedef int point_type;
17: typedef Value_ value_type;
18: typedef Bundle_ bundle_type;
19: typedef typename bundle_type::sieve_type sieve_type;
20: typedef typename ALE::DiscreteSieve<point_type> dsieve_type;
21: typedef typename ALE::Topology<int, dsieve_type> topology_type;
22: typedef typename ALE::Sifter<int, point_type, point_type> send_overlap_type;
23: typedef typename ALE::Sifter<point_type, int, point_type> recv_overlap_type;
24: typedef typename ALE::Field<send_overlap_type, int, ALE::Section<point_type, int> > send_sizer_type;
25: typedef typename ALE::Field<recv_overlap_type, int, ALE::Section<point_type, int> > recv_sizer_type;
26: typedef typename ALE::New::ConeSizeSection<bundle_type, sieve_type> cone_size_section;
27: typedef typename ALE::New::ConeSection<sieve_type> cone_section;
28: typedef typename ALE::New::SectionCompletion<bundle_type, value_type> completion;
29: public:
30: template<typename PartitionType>
31: static void scatterSieve(const Obj<bundle_type>& bundle, const Obj<sieve_type>& sieve, const int dim, const Obj<sieve_type>& sieveNew, const Obj<send_overlap_type>& sendOverlap, const Obj<recv_overlap_type>& recvOverlap, const int height, const int numCells, const PartitionType assignment[]) {
32: typedef typename ALE::Field<send_overlap_type, int, ALE::Section<point_type, value_type> > send_section_type;
33: typedef typename ALE::Field<recv_overlap_type, int, ALE::Section<point_type, value_type> > recv_section_type;
34: int rank = sieve->commRank();
35: int debug = sieve->debug();
37: // Create local sieve
38: const Obj<typename bundle_type::label_sequence>& cells = bundle->heightStratum(height);
39: int e = 0;
41: if (sieve->debug()) {
42: int e2 = 0;
43: for(typename bundle_type::label_sequence::iterator e_iter = cells->begin(); e_iter != cells->end(); ++e_iter) {
44: std::cout << "assignment["<<*e_iter<<"]" << assignment[e2++] << std::endl;
45: }
46: }
47: PetscTruth flg;
48: PetscOptionsHasName(PETSC_NULL, "-output_partition", &flg);
49: if (flg) {
50: ostringstream fname;
51: fname << "part." << sieve->commSize() << ".dat";
52: std::ofstream f(fname.str().c_str());
53: int e2 = 0;
54: f << sieve->commSize() << std::endl;
55: for(typename bundle_type::label_sequence::iterator e_iter = cells->begin(); e_iter != cells->end(); ++e_iter) {
56: f << assignment[e2++] << std::endl;
57: }
58: }
59: for(typename bundle_type::label_sequence::iterator e_iter = cells->begin(); e_iter != cells->end(); ++e_iter) {
60: if (assignment[e] == rank) {
61: Obj<typename sieve_type::coneSet> current = new typename sieve_type::coneSet();
62: Obj<typename sieve_type::coneSet> next = new typename sieve_type::coneSet();
63: Obj<typename sieve_type::coneSet> tmp;
65: current->insert(*e_iter);
66: while(current->size()) {
67: for(typename sieve_type::coneSet::const_iterator p_iter = current->begin(); p_iter != current->end(); ++p_iter) {
68: const Obj<typename sieve_type::traits::coneSequence>& cone = sieve->cone(*p_iter);
69:
70: for(typename sieve_type::traits::coneSequence::iterator c_iter = cone->begin(); c_iter != cone->end(); ++c_iter) {
71: sieveNew->addArrow(*c_iter, *p_iter, c_iter.color());
72: next->insert(*c_iter);
73: }
74: }
75: tmp = current; current = next; next = tmp;
76: next->clear();
77: }
78: if (height) {
79: current->insert(*e_iter);
80: while(current->size()) {
81: for(typename sieve_type::coneSet::const_iterator p_iter = current->begin(); p_iter != current->end(); ++p_iter) {
82: const Obj<typename sieve_type::traits::supportSequence>& support = sieve->support(*p_iter);
83:
84: for(typename sieve_type::traits::supportSequence::iterator s_iter = support->begin(); s_iter != support->end(); ++s_iter) {
85: sieveNew->addArrow(*p_iter, *s_iter, s_iter.color());
86: next->insert(*s_iter);
87: }
88: }
89: tmp = current; current = next; next = tmp;
90: next->clear();
91: }
92: }
93: }
94: e++;
95: }
96: // Complete sizer section
97: typedef typename ALE::New::PartitionSizeSection<bundle_type, PartitionType> partition_size_section;
98: typedef typename ALE::New::PartitionSection<bundle_type, PartitionType> partition_section;
99: Obj<topology_type> secTopology = completion::createSendTopology(sendOverlap);
100: Obj<partition_size_section> partitionSizeSection = new partition_size_section(bundle, height, numCells, assignment);
101: Obj<partition_section> partitionSection = new partition_section(bundle, height, numCells, assignment);
102: Obj<send_section_type> sendSection = new send_section_type(sieve->comm(), sieve->debug());
103: Obj<recv_section_type> recvSection = new recv_section_type(sieve->comm(), sendSection->getTag(), sieve->debug());
105: completion::completeSection(sendOverlap, recvOverlap, partitionSizeSection, partitionSection, sendSection, recvSection);
106: // Unpack the section into the overlap
107: sendOverlap->clear();
108: recvOverlap->clear();
109: const typename send_section_type::sheaf_type& sendPatches = sendSection->getPatches();
111: for(typename send_section_type::sheaf_type::const_iterator p_iter = sendPatches.begin(); p_iter != sendPatches.end(); ++p_iter) {
112: const typename send_section_type::patch_type rank = p_iter->first;
113: const Obj<typename send_section_type::section_type>& section = p_iter->second;
114: const typename send_section_type::section_type::chart_type chart = section->getChart();
116: for(typename send_section_type::section_type::chart_type::iterator c_iter = chart.begin(); c_iter != chart.end(); ++c_iter) {
117: const typename send_section_type::value_type *points = section->restrictPoint(*c_iter);
118: int size = section->getFiberDimension(*c_iter);
120: for(int p = 0; p < size; p++) {
121: sendOverlap->addArrow(points[p], rank, points[p]);
122: }
123: }
124: }
125: const typename recv_section_type::sheaf_type& recvPatches = recvSection->getPatches();
127: for(typename recv_section_type::sheaf_type::const_iterator p_iter = recvPatches.begin(); p_iter != recvPatches.end(); ++p_iter) {
128: const typename send_section_type::patch_type rank = p_iter->first;
129: const Obj<typename send_section_type::section_type>& section = p_iter->second;
130: const typename send_section_type::section_type::chart_type chart = section->getChart();
132: for(typename recv_section_type::section_type::chart_type::iterator c_iter = chart.begin(); c_iter != chart.end(); ++c_iter) {
133: const typename recv_section_type::value_type *points = section->restrictPoint(*c_iter);
134: int size = section->getFiberDimension(*c_iter);
136: for(int p = 0; p < size; p++) {
137: recvOverlap->addArrow(rank, points[p], points[p]);
138: }
139: }
140: }
141: if (debug) {
142: sendOverlap->view(std::cout, "Send overlap for points");
143: recvOverlap->view(std::cout, "Receive overlap for points");
144: }
145: // Receive the point section
146: ALE::New::Completion<bundle_type, value_type>::scatterCones(sieve, sieveNew, sendOverlap, recvOverlap, bundle, height);
147: if (height) {
148: ALE::New::Completion<bundle_type, value_type>::scatterSupports(sieve, sieveNew, sendOverlap, recvOverlap, bundle, bundle->depth()-height);
149: }
150: };
151: template<typename SifterType>
152: static void scatterCones(const Obj<SifterType>& sifter, const Obj<SifterType>& sifterNew, const Obj<send_overlap_type>& sendOverlap, const Obj<recv_overlap_type>& recvOverlap, const Obj<bundle_type>& bundle = NULL, const int minimumHeight = 0) {
153: typedef typename ALE::New::ConeSizeSection<bundle_type, SifterType> cone_size_section;
154: typedef typename ALE::New::ConeSection<SifterType> cone_section;
155: typedef typename ALE::Field<send_overlap_type, int, ALE::Section<point_type, value_type> > send_section_type;
156: typedef typename ALE::Field<recv_overlap_type, int, ALE::Section<point_type, value_type> > recv_section_type;
157: Obj<topology_type> secTopology = completion::createSendTopology(sendOverlap);
158: Obj<cone_size_section> coneSizeSection = new cone_size_section(bundle, sifter, minimumHeight);
159: Obj<cone_section> coneSection = new cone_section(sifter);
160: Obj<send_section_type> sendSection = new send_section_type(sifter->comm(), sifter->debug());
161: Obj<recv_section_type> recvSection = new recv_section_type(sifter->comm(), sendSection->getTag(), sifter->debug());
163: completion::completeSection(sendOverlap, recvOverlap, coneSizeSection, coneSection, sendSection, recvSection);
164: // Unpack the section into the sieve
165: const typename recv_section_type::sheaf_type& patches = recvSection->getPatches();
167: for(typename recv_section_type::sheaf_type::const_iterator p_iter = patches.begin(); p_iter != patches.end(); ++p_iter) {
168: const Obj<typename recv_section_type::section_type>& section = p_iter->second;
169: const typename recv_section_type::section_type::chart_type& chart = section->getChart();
171: for(typename recv_section_type::section_type::chart_type::iterator c_iter = chart.begin(); c_iter != chart.end(); ++c_iter) {
172: const typename recv_section_type::value_type *points = section->restrictPoint(*c_iter);
173: int size = section->getFiberDimension(*c_iter);
174: int c = 0;
176: for(int p = 0; p < size; p++) {
177: sifterNew->addArrow(points[p], *c_iter, c++);
178: }
179: }
180: }
181: };
182: template<typename SifterType>
183: static void scatterSupports(const Obj<SifterType>& sifter, const Obj<SifterType>& sifterNew, const Obj<send_overlap_type>& sendOverlap, const Obj<recv_overlap_type>& recvOverlap, const Obj<bundle_type>& bundle = NULL, const int minimumDepth = 0) {
184: typedef typename ALE::New::SupportSizeSection<bundle_type, SifterType> support_size_section;
185: typedef typename ALE::New::SupportSection<SifterType> support_section;
186: typedef typename ALE::Field<send_overlap_type, int, ALE::Section<point_type, value_type> > send_section_type;
187: typedef typename ALE::Field<recv_overlap_type, int, ALE::Section<point_type, value_type> > recv_section_type;
188: Obj<topology_type> secTopology = completion::createSendTopology(sendOverlap);
189: Obj<support_size_section> supportSizeSection = new support_size_section(bundle, sifter, minimumDepth);
190: Obj<support_section> supportSection = new support_section(sifter);
191: Obj<send_section_type> sendSection = new send_section_type(sifter->comm(), sifter->debug());
192: Obj<recv_section_type> recvSection = new recv_section_type(sifter->comm(), sendSection->getTag(), sifter->debug());
194: completion::completeSection(sendOverlap, recvOverlap, supportSizeSection, supportSection, sendSection, recvSection);
195: // Unpack the section into the sieve
196: const typename recv_section_type::sheaf_type& recvPatches = recvSection->getPatches();
198: for(typename recv_section_type::sheaf_type::const_iterator p_iter = recvPatches.begin(); p_iter != recvPatches.end(); ++p_iter) {
199: const Obj<typename send_section_type::section_type>& section = p_iter->second;
200: const typename send_section_type::section_type::chart_type chart = section->getChart();
202: for(typename recv_section_type::section_type::chart_type::iterator c_iter = chart.begin(); c_iter != chart.end(); ++c_iter) {
203: const typename recv_section_type::value_type *points = section->restrictPoint(*c_iter);
204: int size = section->getFiberDimension(*c_iter);
205: int c = 0;
207: for(int p = 0; p < size; p++) {
208: sifterNew->addArrow(*c_iter, points[p], c++);
209: }
210: }
211: }
212: };
213: };
215: template<typename Value_>
216: class ParallelFactory {
217: public:
218: typedef Value_ value_type;
219: protected:
220: int _debug;
221: MPI_Datatype _mpiType;
222: protected:
223: MPI_Datatype constructMPIType() {
224: if (sizeof(value_type) == 4) {
225: return MPI_INT;
226: } else if (sizeof(value_type) == 8) {
227: return MPI_DOUBLE;
228: } else if (sizeof(value_type) == 28) {
229: int blen[2];
230: MPI_Aint indices[2];
231: MPI_Datatype oldtypes[2], newtype;
232: blen[0] = 1; indices[0] = 0; oldtypes[0] = MPI_INT;
233: blen[1] = 3; indices[1] = sizeof(int); oldtypes[1] = MPI_DOUBLE;
234: MPI_Type_struct(2, blen, indices, oldtypes, &newtype);
235: MPI_Type_commit(&newtype);
236: return newtype;
237: } else if (sizeof(value_type) == 32) {
238: int blen[2];
239: MPI_Aint indices[2];
240: MPI_Datatype oldtypes[2], newtype;
241: blen[0] = 1; indices[0] = 0; oldtypes[0] = MPI_DOUBLE;
242: blen[1] = 3; indices[1] = sizeof(int); oldtypes[1] = MPI_DOUBLE;
243: MPI_Type_struct(2, blen, indices, oldtypes, &newtype);
244: MPI_Type_commit(&newtype);
245: return newtype;
246: }
247: throw ALE::Exception("Cannot determine MPI type for value type");
248: };
249: ParallelFactory(const int debug) : _debug(debug) {
250: this->_mpiType = this->constructMPIType();
251: };
252: public:
253: ~ParallelFactory() {};
254: public:
255: static const Obj<ParallelFactory>& singleton(const int debug, bool cleanup = false) {
256: static Obj<ParallelFactory> *_singleton = NULL;
258: if (cleanup) {
259: if (debug) {std::cout << "Destroying ParallelFactory" << std::endl;}
260: if (_singleton) {delete _singleton;}
261: _singleton = NULL;
262: } else if (_singleton == NULL) {
263: if (debug) {std::cout << "Creating new ParallelFactory" << std::endl;}
264: _singleton = new Obj<ParallelFactory>();
265: *_singleton = new ParallelFactory(debug);
266: }
267: return *_singleton;
268: };
269: public: // Accessors
270: int debug() const {return this->_debug;};
271: MPI_Datatype getMPIType() const {return this->_mpiType;};
272: };
273: }
274: }
276: namespace ALECompat {
277: namespace New {
278: template<typename Section_>
279: class PatchlessSection : public ALE::ParallelObject {
280: public:
281: typedef Section_ section_type;
282: typedef typename section_type::patch_type patch_type;
283: typedef typename section_type::sieve_type sieve_type;
284: typedef typename section_type::point_type point_type;
285: typedef typename section_type::value_type value_type;
286: typedef typename section_type::chart_type chart_type;
287: protected:
288: Obj<section_type> _section;
289: const patch_type _patch;
290: public:
291: PatchlessSection(const Obj<section_type>& section, const patch_type& patch) : ParallelObject(MPI_COMM_SELF, section->debug()), _section(section), _patch(patch) {};
292: virtual ~PatchlessSection() {};
293: public:
294: const chart_type& getPatch(const patch_type& patch) {
295: return this->_section->getAtlas()->getPatch(this->_patch);
296: };
297: bool hasPoint(const patch_type& patch, const point_type& point) {
298: return this->_section->hasPoint(patch, point);
299: };
300: const value_type *restrict(const patch_type& patch) {
301: return this->_section->restrict(this->_patch);
302: };
303: const value_type *restrict(const patch_type& patch, const point_type& p) {
304: return this->_section->restrict(this->_patch, p);
305: };
306: const value_type *restrictPoint(const patch_type& patch, const point_type& p) {
307: return this->_section->restrictPoint(this->_patch, p);
308: };
309: void update(const patch_type& patch, const point_type& p, const value_type v[]) {
310: this->_section->update(this->_patch, p, v);
311: };
312: void updateAdd(const patch_type& patch, const point_type& p, const value_type v[]) {
313: this->_section->updateAdd(this->_patch, p, v);
314: };
315: void updatePoint(const patch_type& patch, const point_type& p, const value_type v[]) {
316: this->_section->updatePoint(this->_patch, p, v);
317: };
318: template<typename Input>
319: void update(const patch_type& patch, const point_type& p, const Obj<Input>& v) {
320: this->_section->update(this->_patch, p, v);
321: };
322: public:
323: void view(const std::string& name, MPI_Comm comm = MPI_COMM_NULL) const {
324: this->_section->view(name, comm);
325: };
326: };
328: template<typename Topology_, typename Value_>
329: class Completion {
330: public:
331: typedef int point_type;
332: typedef Value_ value_type;
333: typedef Topology_ mesh_topology_type;
334: typedef typename mesh_topology_type::sieve_type sieve_type;
335: typedef typename ALE::DiscreteSieve<point_type> dsieve_type;
336: typedef typename ALE::Topology<int, dsieve_type> topology_type;
337: typedef typename ALE::Sifter<int, point_type, point_type> send_overlap_type;
338: typedef typename ALECompat::New::OverlapValues<send_overlap_type, topology_type, int> send_sizer_type;
339: typedef typename ALE::Sifter<point_type, int, point_type> recv_overlap_type;
340: typedef typename ALECompat::New::OverlapValues<recv_overlap_type, topology_type, int> recv_sizer_type;
341: typedef typename ALECompat::New::OldConstantSection<topology_type, int> constant_sizer;
342: typedef typename ALECompat::New::OldConstantSection<topology_type, value_type> constant_section;
343: typedef typename ALECompat::New::ConeSizeSection<topology_type, mesh_topology_type, sieve_type> cone_size_section;
344: typedef typename ALECompat::New::ConeSection<topology_type, sieve_type> cone_section;
345: typedef typename ALECompat::New::SectionCompletion<mesh_topology_type,value_type> completion;
346: public:
347: template<typename PartitionType>
348: static void scatterSieve(const Obj<mesh_topology_type>& topology, const Obj<sieve_type>& sieve, const int dim, const Obj<sieve_type>& sieveNew, const Obj<send_overlap_type>& sendOverlap, const Obj<recv_overlap_type>& recvOverlap, const int numCells, const PartitionType assignment[]) {
349: typedef typename ALECompat::New::OverlapValues<send_overlap_type, topology_type, value_type> send_section_type;
350: typedef typename ALECompat::New::OverlapValues<recv_overlap_type, topology_type, value_type> recv_section_type;
351: int rank = sieve->commRank();
352: int debug = sieve->debug();
354: // Create local sieve
355: const Obj<topology_type::label_sequence>& cells = topology->heightStratum(0, 0);
356: int e = 0;
358: if (topology->debug()) {
359: int e2 = 0;
360: for(topology_type::label_sequence::iterator e_iter = cells->begin(); e_iter != cells->end(); ++e_iter) {
361: std::cout << "assignment["<<*e_iter<<"]" << assignment[e2++] << std::endl;
362: }
363: }
364: PetscTruth flg;
365: PetscOptionsHasName(PETSC_NULL, "-output_partition", &flg);
366: if (flg) {
367: ostringstream fname;
368: fname << "part." << sieve->commSize() << ".dat";
369: std::ofstream f(fname.str().c_str());
370: int e2 = 0;
371: f << sieve->commSize() << std::endl;
372: for(topology_type::label_sequence::iterator e_iter = cells->begin(); e_iter != cells->end(); ++e_iter) {
373: f << assignment[e2++] << std::endl;
374: }
375: }
376: for(topology_type::label_sequence::iterator e_iter = cells->begin(); e_iter != cells->end(); ++e_iter) {
377: if (assignment[e] == rank) {
378: Obj<typename sieve_type::coneSet> current = new typename sieve_type::coneSet();
379: Obj<typename sieve_type::coneSet> next = new typename sieve_type::coneSet();
380: Obj<typename sieve_type::coneSet> tmp;
382: current->insert(*e_iter);
383: while(current->size()) {
384: for(typename sieve_type::coneSet::const_iterator p_iter = current->begin(); p_iter != current->end(); ++p_iter) {
385: const Obj<typename sieve_type::traits::coneSequence>& cone = sieve->cone(*p_iter);
386:
387: for(typename sieve_type::traits::coneSequence::iterator c_iter = cone->begin(); c_iter != cone->end(); ++c_iter) {
388: sieveNew->addArrow(*c_iter, *p_iter, c_iter.color());
389: next->insert(*c_iter);
390: }
391: }
392: tmp = current; current = next; next = tmp;
393: next->clear();
394: }
395: }
396: e++;
397: }
398: sieveNew->stratify();
399: // Complete sizer section
400: typedef typename ALECompat::New::PartitionSizeSection<topology_type, mesh_topology_type, PartitionType> partition_size_section;
401: typedef typename ALECompat::New::PartitionSection<topology_type, mesh_topology_type, PartitionType> partition_section;
402: Obj<topology_type> secTopology = completion::createSendTopology(sendOverlap);
403: Obj<partition_size_section> partitionSizeSection = new partition_size_section(secTopology, topology, 0, numCells, assignment);
404: Obj<partition_section> partitionSection = new partition_section(secTopology, topology, 0, numCells, assignment);
405: Obj<send_section_type> sendSection = new send_section_type(sieve->comm(), sieve->debug());
406: Obj<recv_section_type> recvSection = new recv_section_type(sieve->comm(), sendSection->getTag(), sieve->debug());
408: completion::completeSection(sendOverlap, recvOverlap, partitionSizeSection, partitionSection, sendSection, recvSection);
409: // Unpack the section into the overlap
410: sendOverlap->clear();
411: recvOverlap->clear();
412: const topology_type::sheaf_type& sendPatches = sendSection->getTopology()->getPatches();
414: for(topology_type::sheaf_type::const_iterator p_iter = sendPatches.begin(); p_iter != sendPatches.end(); ++p_iter) {
415: const Obj<topology_type::sieve_type::baseSequence>& base = p_iter->second->base();
417: for(topology_type::sieve_type::baseSequence::iterator b_iter = base->begin(); b_iter != base->end(); ++b_iter) {
418: const typename send_section_type::value_type *points = sendSection->restrict(p_iter->first, *b_iter);
419: int size = sendSection->size(p_iter->first, *b_iter);
421: for(int p = 0; p < size; p++) {
422: sendOverlap->addArrow(points[p], p_iter->first, points[p]);
423: }
424: }
425: }
426: const topology_type::sheaf_type& recvPatches = recvSection->getTopology()->getPatches();
428: for(topology_type::sheaf_type::const_iterator p_iter = recvPatches.begin(); p_iter != recvPatches.end(); ++p_iter) {
429: const Obj<topology_type::sieve_type::baseSequence>& base = p_iter->second->base();
430: int rank = p_iter->first;
432: for(topology_type::sieve_type::baseSequence::iterator b_iter = base->begin(); b_iter != base->end(); ++b_iter) {
433: const typename recv_section_type::value_type *points = recvSection->restrict(rank, *b_iter);
434: int size = recvSection->getFiberDimension(rank, *b_iter);
436: for(int p = 0; p < size; p++) {
437: recvOverlap->addArrow(rank, points[p], points[p]);
438: }
439: }
440: }
441: if (debug) {
442: sendOverlap->view(std::cout, "Send overlap for points");
443: recvOverlap->view(std::cout, "Receive overlap for points");
444: }
445: // Receive the point section
446: ALECompat::New::Completion<mesh_topology_type,value_type>::scatterCones(sieve, sieveNew, sendOverlap, recvOverlap);
447: sieveNew->stratify();
448: };
449: template<typename PartitionType>
450: static void scatterSieveByFace(const Obj<mesh_topology_type>& topology, const Obj<sieve_type>& sieve, const int dim, const Obj<sieve_type>& sieveNew, const Obj<send_overlap_type>& sendOverlap, const Obj<recv_overlap_type>& recvOverlap, const int numFaces, const PartitionType assignment[]) {
451: typedef typename ALECompat::New::OverlapValues<send_overlap_type, topology_type, value_type> send_section_type;
452: typedef typename ALECompat::New::OverlapValues<recv_overlap_type, topology_type, value_type> recv_section_type;
453: const typename topology_type::patch_type patch = 0;
454: int rank = sieve->commRank();
455: int debug = sieve->debug();
457: // Create local sieve
458: const Obj<topology_type::label_sequence>& faces = topology->heightStratum(patch, 1);
459: int f = 0;
461: for(topology_type::label_sequence::iterator f_iter = faces->begin(); f_iter != faces->end(); ++f_iter) {
462: if (assignment[f] == rank) {
463: Obj<typename sieve_type::coneSet> current = new typename sieve_type::coneSet();
464: Obj<typename sieve_type::coneSet> next = new typename sieve_type::coneSet();
465: Obj<typename sieve_type::coneSet> tmp;
467: current->insert(*f_iter);
468: while(current->size()) {
469: for(typename sieve_type::coneSet::const_iterator p_iter = current->begin(); p_iter != current->end(); ++p_iter) {
470: const Obj<typename sieve_type::traits::coneSequence>& cone = sieve->cone(*p_iter);
471:
472: for(typename sieve_type::traits::coneSequence::iterator c_iter = cone->begin(); c_iter != cone->end(); ++c_iter) {
473: sieveNew->addArrow(*c_iter, *p_iter, c_iter.color());
474: next->insert(*c_iter);
475: }
476: }
477: tmp = current; current = next; next = tmp;
478: next->clear();
479: }
480: const Obj<typename sieve_type::traits::supportSequence>& support = sieve->support(*f_iter);
482: for(typename sieve_type::traits::supportSequence::iterator s_iter = support->begin(); s_iter != support->end(); ++s_iter) {
483: sieveNew->addArrow(*f_iter, *s_iter, s_iter.color());
484: }
485: }
486: f++;
487: }
488: sieveNew->stratify();
489: // Complete sizer section
490: typedef typename ALECompat::New::PartitionSizeSection<topology_type, mesh_topology_type, PartitionType> partition_size_section;
491: typedef typename ALECompat::New::PartitionSection<topology_type, mesh_topology_type, PartitionType> partition_section;
492: Obj<topology_type> secTopology = completion::createSendTopology(sendOverlap);
493: Obj<partition_size_section> partitionSizeSection = new partition_size_section(secTopology, topology, 1, numFaces, assignment);
494: Obj<partition_section> partitionSection = new partition_section(secTopology, topology, 1, numFaces, assignment);
495: Obj<send_section_type> sendSection = new send_section_type(sieve->comm(), sieve->debug());
496: Obj<recv_section_type> recvSection = new recv_section_type(sieve->comm(), sendSection->getTag(), sieve->debug());
498: completion::completeSection(sendOverlap, recvOverlap, partitionSizeSection, partitionSection, sendSection, recvSection);
499: // Unpack the section into the overlap
500: sendOverlap->clear();
501: recvOverlap->clear();
502: const topology_type::sheaf_type& sendPatches = sendSection->getTopology()->getPatches();
504: for(topology_type::sheaf_type::const_iterator p_iter = sendPatches.begin(); p_iter != sendPatches.end(); ++p_iter) {
505: const Obj<topology_type::sieve_type::baseSequence>& base = p_iter->second->base();
507: for(topology_type::sieve_type::baseSequence::iterator b_iter = base->begin(); b_iter != base->end(); ++b_iter) {
508: const typename send_section_type::value_type *points = sendSection->restrict(p_iter->first, *b_iter);
509: int size = sendSection->size(p_iter->first, *b_iter);
511: for(int p = 0; p < size; p++) {
512: sendOverlap->addArrow(points[p], p_iter->first, points[p]);
513: }
514: }
515: }
516: const topology_type::sheaf_type& recvPatches = recvSection->getTopology()->getPatches();
518: for(topology_type::sheaf_type::const_iterator p_iter = recvPatches.begin(); p_iter != recvPatches.end(); ++p_iter) {
519: const Obj<topology_type::sieve_type::baseSequence>& base = p_iter->second->base();
520: int rank = p_iter->first;
522: for(topology_type::sieve_type::baseSequence::iterator b_iter = base->begin(); b_iter != base->end(); ++b_iter) {
523: const typename recv_section_type::value_type *points = recvSection->restrict(rank, *b_iter);
524: int size = recvSection->getFiberDimension(rank, *b_iter);
526: for(int p = 0; p < size; p++) {
527: recvOverlap->addArrow(rank, points[p], points[p]);
528: }
529: }
530: }
531: if (debug) {
532: sendOverlap->view(std::cout, "Send overlap for points");
533: recvOverlap->view(std::cout, "Receive overlap for points");
534: }
535: // Receive the point section
536: ALECompat::New::Completion<mesh_topology_type,value_type>::scatterCones(sieve, sieveNew, sendOverlap, recvOverlap, topology, 1);
537: ALECompat::New::Completion<mesh_topology_type,value_type>::scatterSupports(sieve, sieveNew, sendOverlap, recvOverlap, topology, topology->depth()-1);
538: sieveNew->stratify();
539: };
540: template<typename SifterType>
541: static void scatterCones(const Obj<SifterType>& sifter, const Obj<SifterType>& sifterNew, const Obj<send_overlap_type>& sendOverlap, const Obj<recv_overlap_type>& recvOverlap, const Obj<mesh_topology_type>& topology = NULL, const int minimumHeight = 0) {
542: typedef typename ALECompat::New::ConeSizeSection<topology_type, mesh_topology_type, SifterType> cone_size_section;
543: typedef typename ALECompat::New::ConeSection<topology_type, SifterType> cone_section;
544: typedef typename ALECompat::New::OverlapValues<send_overlap_type, topology_type, value_type> send_section_type;
545: typedef typename ALECompat::New::OverlapValues<recv_overlap_type, topology_type, value_type> recv_section_type;
546: Obj<topology_type> secTopology = completion::createSendTopology(sendOverlap);
547: Obj<cone_size_section> coneSizeSection = new cone_size_section(secTopology, topology, sifter, minimumHeight);
548: Obj<cone_section> coneSection = new cone_section(secTopology, sifter);
549: Obj<send_section_type> sendSection = new send_section_type(sifter->comm(), sifter->debug());
550: Obj<recv_section_type> recvSection = new recv_section_type(sifter->comm(), sendSection->getTag(), sifter->debug());
552: completion::completeSection(sendOverlap, recvOverlap, coneSizeSection, coneSection, sendSection, recvSection);
553: // Unpack the section into the sieve
554: const topology_type::sheaf_type& patches = recvSection->getTopology()->getPatches();
556: for(topology_type::sheaf_type::const_iterator p_iter = patches.begin(); p_iter != patches.end(); ++p_iter) {
557: const Obj<topology_type::sieve_type::baseSequence>& base = p_iter->second->base();
558: int rank = p_iter->first;
560: for(topology_type::sieve_type::baseSequence::iterator b_iter = base->begin(); b_iter != base->end(); ++b_iter) {
561: const typename recv_section_type::value_type *points = recvSection->restrict(rank, *b_iter);
562: int size = recvSection->getFiberDimension(rank, *b_iter);
563: int c = 0;
565: for(int p = 0; p < size; p++) {
566: sifterNew->addArrow(points[p], *b_iter, c++);
567: }
568: }
569: }
570: };
571: template<typename SifterType>
572: static void scatterSupports(const Obj<SifterType>& sifter, const Obj<SifterType>& sifterNew, const Obj<send_overlap_type>& sendOverlap, const Obj<recv_overlap_type>& recvOverlap, const Obj<mesh_topology_type>& topology = NULL, const int minimumDepth = 0) {
573: typedef typename ALECompat::New::SupportSizeSection<topology_type, mesh_topology_type, SifterType> support_size_section;
574: typedef typename ALECompat::New::SupportSection<topology_type, SifterType> support_section;
575: typedef typename ALECompat::New::OverlapValues<send_overlap_type, topology_type, value_type> send_section_type;
576: typedef typename ALECompat::New::OverlapValues<recv_overlap_type, topology_type, value_type> recv_section_type;
577: Obj<topology_type> secTopology = completion::createSendTopology(sendOverlap);
578: Obj<support_size_section> supportSizeSection = new support_size_section(secTopology, topology, sifter, minimumDepth);
579: Obj<support_section> supportSection = new support_section(secTopology, sifter);
580: Obj<send_section_type> sendSection = new send_section_type(sifter->comm(), sifter->debug());
581: Obj<recv_section_type> recvSection = new recv_section_type(sifter->comm(), sendSection->getTag(), sifter->debug());
583: completion::completeSection(sendOverlap, recvOverlap, supportSizeSection, supportSection, sendSection, recvSection);
584: // Unpack the section into the sieve
585: const topology_type::sheaf_type& patches = recvSection->getTopology()->getPatches();
587: for(topology_type::sheaf_type::const_iterator p_iter = patches.begin(); p_iter != patches.end(); ++p_iter) {
588: const Obj<topology_type::sieve_type::baseSequence>& base = p_iter->second->base();
589: int rank = p_iter->first;
591: for(topology_type::sieve_type::baseSequence::iterator b_iter = base->begin(); b_iter != base->end(); ++b_iter) {
592: const typename recv_section_type::value_type *points = recvSection->restrict(rank, *b_iter);
593: int size = recvSection->getFiberDimension(rank, *b_iter);
594: int c = 0;
596: for(int p = 0; p < size; p++) {
597: sifterNew->addArrow(*b_iter, points[p], c++);
598: }
599: }
600: }
601: };
602: };
604: template<typename Value_>
605: class ParallelFactory {
606: public:
607: typedef Value_ value_type;
608: protected:
609: int _debug;
610: MPI_Datatype _mpiType;
611: protected:
612: MPI_Datatype constructMPIType() {
613: if (sizeof(value_type) == 4) {
614: return MPI_INT;
615: } else if (sizeof(value_type) == 8) {
616: return MPI_DOUBLE;
617: } else if (sizeof(value_type) == 28) {
618: int blen[2];
619: MPI_Aint indices[2];
620: MPI_Datatype oldtypes[2], newtype;
621: blen[0] = 1; indices[0] = 0; oldtypes[0] = MPI_INT;
622: blen[1] = 3; indices[1] = sizeof(int); oldtypes[1] = MPI_DOUBLE;
623: MPI_Type_struct(2, blen, indices, oldtypes, &newtype);
624: MPI_Type_commit(&newtype);
625: return newtype;
626: } else if (sizeof(value_type) == 32) {
627: int blen[2];
628: MPI_Aint indices[2];
629: MPI_Datatype oldtypes[2], newtype;
630: blen[0] = 1; indices[0] = 0; oldtypes[0] = MPI_DOUBLE;
631: blen[1] = 3; indices[1] = sizeof(int); oldtypes[1] = MPI_DOUBLE;
632: MPI_Type_struct(2, blen, indices, oldtypes, &newtype);
633: MPI_Type_commit(&newtype);
634: return newtype;
635: }
636: throw ALE::Exception("Cannot determine MPI type for value type");
637: };
638: ParallelFactory(const int debug) : _debug(debug) {
639: this->_mpiType = this->constructMPIType();
640: };
641: public:
642: ~ParallelFactory() {};
643: public:
644: static const Obj<ParallelFactory>& singleton(const int debug, bool cleanup = false) {
645: static Obj<ParallelFactory> *_singleton = NULL;
647: if (cleanup) {
648: if (debug) {std::cout << "Destroying ParallelFactory" << std::endl;}
649: if (_singleton) {delete _singleton;}
650: _singleton = NULL;
651: } else if (_singleton == NULL) {
652: if (debug) {std::cout << "Creating new ParallelFactory" << std::endl;}
653: _singleton = new Obj<ParallelFactory>();
654: *_singleton = new ParallelFactory(debug);
655: }
656: return *_singleton;
657: };
658: public: // Accessors
659: int debug() const {return this->_debug;};
660: MPI_Datatype getMPIType() const {return this->_mpiType;};
661: };
662: }
663: }
664: #endif