Actual source code: Distribution.hh
1: #ifndef included_ALE_Distribution_hh
2: #define included_ALE_Distribution_hh
4: #ifndef included_ALE_Mesh_hh
5: #include <Mesh.hh>
6: #endif
8: #ifndef included_ALE_Partitioner_hh
9: #include <Partitioner.hh>
10: #endif
12: #ifndef included_ALE_Completion_hh
13: #include <Completion.hh>
14: #endif
16: // Attempt to unify all of the distribution mechanisms:
17: // one to many (distributeMesh)
18: // many to one (unifyMesh)
19: // many to many (Numbering)
20: // as well as things being distributed
21: // Section
22: // Sieve (This sends two sections, the points and cones)
23: // Numbering (Should be an integer section)
24: // Global Order (should be an integer section with extra methods)
25: //
26: // 0) Create the new object to hold the communicated data
27: //
28: // 1) Create Overlap
29: // There may be special ways to do this based upon what we know at the time
30: //
31: // 2) Create send and receive sections over the interface
32: // These have a flat topology now, consisting only of the overlap nodes
33: // We could make a full topology on the overlap (maybe it is necessary for higher order)
34: //
35: // 3) Communication section
36: // Create sizer sections on interface (uses constant sizer)
37: // Communicate sizes on interface (uses custom filler)
38: // Fill send section
39: // sendSection->startCommunication();
40: // recvSection->startCommunication();
41: // sendSection->endCommunication();
42: // recvSection->endCommunication();
43: //
44: // Create section on interface (uses previous sizer)
45: // Communicate values on interface (uses custom filler)
46: // Same stuff as above
47: //
48: // 4) Update new section with old local values (can be done in between the communication?)
49: // Loop over patches in new topology
50: // Loop over chart from patch in old atlas
51: // If this point is in the new sieve from patch
52: // Set to old fiber dimension
53: // Order and allocate new section
54: // Repeat loop, but update values
55: //
56: // 5) Update new section with old received values
57: // Loop over patches in discrete topology of receive section (these are ranks)
58: // Loop over base of discrete sieve (we should transform this to a chart to match above)
59: // Get new patch from overlap, or should the receive patches be <rank, patch>?
60: // Guaranteed to be in the new sieve from patch (but we could check anyway)
61: // Set to recevied fiber dimension
62: // Order and allocate new section
63: // Repeat loop, but update values
64: //
65: // 6) Synchronize PETSc tags (can I get around this?)
66: namespace ALE {
67: template<typename Bundle_>
68: class Distribution {
69: public:
70: typedef Bundle_ bundle_type;
71: typedef typename bundle_type::sieve_type sieve_type;
72: typedef typename bundle_type::point_type point_type;
73: typedef typename ALE::New::Completion<bundle_type, typename sieve_type::point_type> sieveCompletion;
74: typedef typename ALE::New::SectionCompletion<bundle_type, typename bundle_type::real_section_type::value_type> sectionCompletion;
75: typedef typename sectionCompletion::send_overlap_type send_overlap_type;
76: typedef typename sectionCompletion::recv_overlap_type recv_overlap_type;
77: public:
80: static void createPartitionOverlap(const Obj<bundle_type>& bundle, const Obj<send_overlap_type>& sendOverlap, const Obj<recv_overlap_type>& recvOverlap) {
81: const Obj<send_overlap_type>& topSendOverlap = bundle->getSendOverlap();
82: const Obj<recv_overlap_type>& topRecvOverlap = bundle->getRecvOverlap();
83: const Obj<typename send_overlap_type::traits::baseSequence> base = topSendOverlap->base();
84: const Obj<typename recv_overlap_type::traits::capSequence> cap = topRecvOverlap->cap();
85: const int rank = bundle->commRank();
87: if (base->empty()) {
88: if (rank == 0) {
89: for(int p = 1; p < bundle->commSize(); p++) {
90: // The arrow is from local partition point p (source) to remote partition point p (color) on rank p (target)
91: sendOverlap->addCone(p, p, p);
92: }
93: }
94: } else {
95: for(typename send_overlap_type::traits::baseSequence::iterator b_iter = base->begin(); b_iter != base->end(); ++b_iter) {
96: const int& p = *b_iter;
97: // The arrow is from local partition point p (source) to remote partition point p (color) on rank p (target)
98: sendOverlap->addCone(p, p, p);
99: }
100: }
101: if (cap->empty()) {
102: if (rank != 0) {
103: // The arrow is from local partition point rank (color) on rank 0 (source) to remote partition point rank (target)
104: recvOverlap->addCone(0, rank, rank);
105: }
106: } else {
107: for(typename recv_overlap_type::traits::capSequence::iterator c_iter = cap->begin(); c_iter != cap->end(); ++c_iter) {
108: const int& p = *c_iter;
109: // The arrow is from local partition point rank (color) on rank p (source) to remote partition point rank (target)
110: recvOverlap->addCone(p, rank, rank);
111: }
112: }
113: };
116: template<typename Partitioner>
117: static typename Partitioner::part_type *createAssignment(const Obj<bundle_type>& bundle, const int dim, const Obj<send_overlap_type>& sendOverlap, const Obj<recv_overlap_type>& recvOverlap, const int height = 0) {
118: // 1) Form partition point overlap a priori
119: createPartitionOverlap(bundle, sendOverlap, recvOverlap);
120: if (bundle->debug()) {
121: sendOverlap->view("Send overlap for partition");
122: recvOverlap->view("Receive overlap for partition");
123: }
124: // 2) Partition the mesh
125: if (height == 0) {
126: return Partitioner::partitionSieve(bundle, dim);
127: } else if (height == 1) {
128: return Partitioner::partitionSieveByFace(bundle, dim);
129: }
130: throw ALE::Exception("Invalid partition height");
131: };
134: // Partition a bundle on process 0 and scatter to all processes
135: static void scatterBundle(const Obj<bundle_type>& bundle, const int dim, const Obj<bundle_type>& bundleNew, const Obj<send_overlap_type>& sendOverlap, const Obj<recv_overlap_type>& recvOverlap, const std::string& partitioner, const int height = 0, const Obj<bundle_type>& subBundle = NULL, const Obj<bundle_type>& subBundleNew = NULL) {
136: if (height == 0) {
137: if (partitioner == "chaco") {
138: #ifdef PETSC_HAVE_CHACO
139: typedef typename ALE::New::Chaco::Partitioner<bundle_type> Partitioner;
140: typedef typename ALE::New::Partitioner<bundle_type> GenPartitioner;
141: typedef typename Partitioner::part_type part_type;
143: part_type *assignment = scatterBundle<Partitioner>(bundle, dim, bundleNew, sendOverlap, recvOverlap, height);
144: if (!subBundle.isNull() && !subBundleNew.isNull()) {
145: part_type *subAssignment = GenPartitioner::subordinatePartition(bundle, 1, subBundle, assignment);
146: const Obj<sieve_type>& sieve = subBundle->getSieve();
147: const Obj<sieve_type>& sieveNew = new Mesh::sieve_type(subBundle->comm(), subBundle->debug());
148: const int numCells = subBundle->heightStratum(height)->size();
150: subBundleNew->setSieve(sieveNew);
151: sieveCompletion::scatterSieve(subBundle, sieve, dim, sieveNew, sendOverlap, recvOverlap, height, numCells, subAssignment);
152: subBundleNew->stratify();
153: if (subAssignment != NULL) delete [] subAssignment;
154: }
155: if (assignment != NULL) delete [] assignment;
156: #else
157: throw ALE::Exception("Chaco is not installed. Reconfigure with the flag --download-chaco");
158: #endif
159: } else if (partitioner == "parmetis") {
160: #ifdef PETSC_HAVE_PARMETIS
161: typedef typename ALE::New::ParMetis::Partitioner<bundle_type> Partitioner;
162: typedef typename ALE::New::Partitioner<bundle_type> GenPartitioner;
163: typedef typename Partitioner::part_type part_type;
165: part_type *assignment = scatterBundle<Partitioner>(bundle, dim, bundleNew, sendOverlap, recvOverlap, height);
166: if (!subBundle.isNull() && !subBundleNew.isNull()) {
167: part_type *subAssignment = GenPartitioner::subordinatePartition(bundle, 1, subBundle, assignment);
168: const Obj<sieve_type>& sieve = subBundle->getSieve();
169: const Obj<sieve_type>& sieveNew = new Mesh::sieve_type(subBundle->comm(), subBundle->debug());
170: const int numCells = subBundle->heightStratum(height)->size();
172: subBundleNew->setSieve(sieveNew);
173: sieveCompletion::scatterSieve(subBundle, sieve, dim, sieveNew, sendOverlap, recvOverlap, height, numCells, subAssignment);
174: subBundleNew->stratify();
175: if (subAssignment != NULL) delete [] subAssignment;
176: }
177: if (assignment != NULL) delete [] assignment;
178: #else
179: throw ALE::Exception("ParMetis is not installed. Reconfigure with the flag --download-parmetis");
180: #endif
181: } else {
182: throw ALE::Exception("Unknown partitioner");
183: }
184: } else if (height == 1) {
185: if (partitioner == "zoltan") {
186: #ifdef PETSC_HAVE_ZOLTAN
187: typedef typename ALE::New::Zoltan::Partitioner<bundle_type> Partitioner;
188: typedef typename Partitioner::part_type part_type;
190: part_type *assignment = scatterBundle<Partitioner>(bundle, dim, bundleNew, sendOverlap, recvOverlap, height);
191: if (assignment != NULL) delete [] assignment;
192: #else
193: throw ALE::Exception("Zoltan is not installed. Reconfigure with the flag --download-zoltan");
194: #endif
195: } else if (partitioner == "parmetis") {
196: #ifdef PETSC_HAVE_PARMETIS
197: typedef typename ALE::New::ParMetis::Partitioner<bundle_type> Partitioner;
198: typedef typename Partitioner::part_type part_type;
200: part_type *assignment = scatterBundle<Partitioner>(bundle, dim, bundleNew, sendOverlap, recvOverlap, height);
201: if (assignment != NULL) delete [] assignment;
202: #else
203: throw ALE::Exception("ParMetis is not installed. Reconfigure with the flag --download-parmetis");
204: #endif
205: } else {
206: throw ALE::Exception("Unknown partitioner");
207: }
208: } else {
209: throw ALE::Exception("Invalid partition height");
210: }
211: };
212: template<typename Partitioner>
213: static typename Partitioner::part_type *scatterBundle(const Obj<bundle_type>& bundle, const int dim, const Obj<bundle_type>& bundleNew, const Obj<send_overlap_type>& sendOverlap, const Obj<recv_overlap_type>& recvOverlap, const int height = 0) {
214: typename Partitioner::part_type *assignment = createAssignment<Partitioner>(bundle, dim, sendOverlap, recvOverlap, height);
215: const Obj<sieve_type>& sieve = bundle->getSieve();
216: const Obj<sieve_type>& sieveNew = bundleNew->getSieve();
217: const int numPoints = bundle->heightStratum(height)->size();
219: sieveCompletion::scatterSieve(bundle, sieve, dim, sieveNew, sendOverlap, recvOverlap, height, numPoints, assignment);
220: bundleNew->stratify();
221: return assignment;
222: };
225: static Obj<ALE::Mesh> distributeMesh(const Obj<ALE::Mesh>& serialMesh, const int height = 0, const std::string& partitioner = "chaco") {
226: MPI_Comm comm = serialMesh->comm();
227: const int dim = serialMesh->getDimension();
228: Obj<ALE::Mesh> parallelMesh = new ALE::Mesh(comm, dim, serialMesh->debug());
229: const Obj<ALE::Mesh::sieve_type>& parallelSieve = new ALE::Mesh::sieve_type(comm, serialMesh->debug());
231: ALE_LOG_EVENT_BEGIN;
232: parallelMesh->setSieve(parallelSieve);
233: if (serialMesh->debug()) {serialMesh->view("Serial mesh");}
235: // Distribute cones
236: Obj<send_overlap_type> sendOverlap = new send_overlap_type(comm, serialMesh->debug());
237: Obj<recv_overlap_type> recvOverlap = new recv_overlap_type(comm, serialMesh->debug());
238: scatterBundle(serialMesh, dim, parallelMesh, sendOverlap, recvOverlap, partitioner, height);
239: parallelMesh->setDistSendOverlap(sendOverlap);
240: parallelMesh->setDistRecvOverlap(recvOverlap);
242: // Distribute labels
243: const typename bundle_type::labels_type& labels = serialMesh->getLabels();
245: for(typename bundle_type::labels_type::const_iterator l_iter = labels.begin(); l_iter != labels.end(); ++l_iter) {
246: if (parallelMesh->hasLabel(l_iter->first)) continue;
247: const Obj<typename bundle_type::label_type>& serialLabel = l_iter->second;
248: const Obj<typename bundle_type::label_type>& parallelLabel = parallelMesh->createLabel(l_iter->first);
249: // Create local label
250: const Obj<typename bundle_type::label_type::traits::baseSequence>& base = serialLabel->base();
252: for(typename bundle_type::label_type::traits::baseSequence::iterator b_iter = base->begin(); b_iter != base->end(); ++b_iter) {
253: if (parallelSieve->capContains(*b_iter) || parallelSieve->baseContains(*b_iter)) {
254: parallelLabel->addArrow(*serialLabel->cone(*b_iter)->begin(), *b_iter);
255: }
256: }
257: // Get remote labels
258: sieveCompletion::scatterCones(serialLabel, parallelLabel, sendOverlap, recvOverlap);
259: }
261: // Distribute sections
262: Obj<std::set<std::string> > sections = serialMesh->getRealSections();
264: for(std::set<std::string>::iterator name = sections->begin(); name != sections->end(); ++name) {
265: parallelMesh->setRealSection(*name, distributeSection(serialMesh->getRealSection(*name), parallelMesh, sendOverlap, recvOverlap));
266: }
267: sections = serialMesh->getIntSections();
268: for(std::set<std::string>::iterator name = sections->begin(); name != sections->end(); ++name) {
269: parallelMesh->setIntSection(*name, distributeSection(serialMesh->getIntSection(*name), parallelMesh, sendOverlap, recvOverlap));
270: }
271: sections = serialMesh->getArrowSections();
273: for(std::set<std::string>::iterator name = sections->begin(); name != sections->end(); ++name) {
274: parallelMesh->setArrowSection(*name, distributeArrowSection(serialMesh->getArrowSection(*name), serialMesh, parallelMesh, sendOverlap, recvOverlap));
275: }
276: if (parallelMesh->debug()) {parallelMesh->view("Parallel Mesh");}
277: ALE_LOG_EVENT_END;
278: return parallelMesh;
279: };
282: template<typename Section>
283: static void updateSectionLocal(const Obj<Section>& oldSection, const Obj<bundle_type>& newBundle, const Obj<Section>& newSection) {
284: const Obj<typename bundle_type::sieve_type>& newSieve = newBundle->getSieve();
285: const typename Section::atlas_type::chart_type& oldChart = oldSection->getChart();
287: for(typename Section::atlas_type::chart_type::const_iterator c_iter = oldChart.begin(); c_iter != oldChart.end(); ++c_iter) {
288: if (newSieve->capContains(*c_iter) || newSieve->baseContains(*c_iter)) {
289: newSection->setFiberDimension(*c_iter, oldSection->getFiberDimension(*c_iter));
290: }
291: }
292: newBundle->allocate(newSection);
293: const typename Section::atlas_type::chart_type& newChart = newSection->getChart();
295: for(typename Section::atlas_type::chart_type::const_iterator c_iter = newChart.begin(); c_iter != newChart.end(); ++c_iter) {
296: newSection->updatePointAll(*c_iter, oldSection->restrictPoint(*c_iter));
297: }
298: };
301: template<typename RecvSection, typename Section>
302: static void updateSectionRemote(const Obj<recv_overlap_type>& recvOverlap, const Obj<RecvSection>& recvSection, const Obj<bundle_type>& newBundle, const Obj<Section>& newSection) {
303: Obj<typename recv_overlap_type::traits::baseSequence> recvPoints = recvOverlap->base();
305: for(typename recv_overlap_type::traits::baseSequence::iterator r_iter = recvPoints->begin(); r_iter != recvPoints->end(); ++r_iter) {
306: const Obj<typename recv_overlap_type::traits::coneSequence>& recvPatches = recvOverlap->cone(*r_iter);
307: const typename recv_overlap_type::traits::coneSequence::iterator end = recvPatches->end();
309: for(typename recv_overlap_type::traits::coneSequence::iterator p_iter = recvPatches->begin(); p_iter != end; ++p_iter) {
310: newSection->addPoint(*r_iter, recvSection->getSection(*p_iter)->getFiberDimension(*r_iter));
311: }
312: }
313: newBundle->reallocate(newSection);
314: for(typename recv_overlap_type::traits::baseSequence::iterator r_iter = recvPoints->begin(); r_iter != recvPoints->end(); ++r_iter) {
315: const Obj<typename recv_overlap_type::traits::coneSequence>& recvPatches = recvOverlap->cone(*r_iter);
316: const typename recv_overlap_type::traits::coneSequence::iterator end = recvPatches->end();
318: for(typename recv_overlap_type::traits::coneSequence::iterator p_iter = recvPatches->begin(); p_iter != end; ++p_iter) {
319: if (recvSection->getSection(*p_iter)->getFiberDimension(*r_iter)) {
320: newSection->updatePointAll(*r_iter, recvSection->getSection(*p_iter)->restrictPoint(*r_iter));
321: }
322: }
323: }
324: };
327: template<typename Section>
328: static Obj<Section> distributeSection(const Obj<Section>& serialSection, const Obj<bundle_type>& parallelBundle, const Obj<send_overlap_type>& sendOverlap, const Obj<recv_overlap_type>& recvOverlap) {
329: if (serialSection->debug()) {
330: serialSection->view("Serial Section");
331: }
332: typedef ALE::Field<send_overlap_type, int, ALE::Section<point_type, typename Section::value_type> > send_section_type;
333: typedef ALE::Field<recv_overlap_type, int, ALE::Section<point_type, typename Section::value_type> > recv_section_type;
334: typedef ALE::New::SizeSection<Section> SectionSizer;
335: Obj<Section> parallelSection = new Section(serialSection->comm(), serialSection->debug());
336: const Obj<send_section_type> sendSection = new send_section_type(serialSection->comm(), serialSection->debug());
337: const Obj<recv_section_type> recvSection = new recv_section_type(serialSection->comm(), sendSection->getTag(), serialSection->debug());
338: const Obj<SectionSizer> sizer = new SectionSizer(serialSection);
340: updateSectionLocal(serialSection, parallelBundle, parallelSection);
341: sectionCompletion::completeSection(sendOverlap, recvOverlap, sizer, serialSection, sendSection, recvSection);
342: updateSectionRemote(recvOverlap, recvSection, parallelBundle, parallelSection);
343: if (parallelSection->debug()) {
344: parallelSection->view("Parallel Section");
345: }
346: return parallelSection;
347: };
350: template<typename Section>
351: static void updateArrowSectionLocal(const Obj<Section>& oldSection, const Obj<bundle_type>& newBundle, const Obj<Section>& newSection) {
352: const Obj<typename bundle_type::sieve_type>& newSieve = newBundle->getSieve();
353: const typename Section::atlas_type::chart_type& oldChart = oldSection->getChart();
355: for(typename Section::atlas_type::chart_type::const_iterator c_iter = oldChart.begin(); c_iter != oldChart.end(); ++c_iter) {
356: // Dmitry should provide a Sieve::contains(MinimalArrow) method
357: if (newSieve->capContains(c_iter->source) && newSieve->baseContains(c_iter->target)) {
358: newSection->setFiberDimension(*c_iter, oldSection->getFiberDimension(*c_iter));
359: }
360: }
361: //newBundle->allocate(newSection);
362: const typename Section::atlas_type::chart_type& newChart = newSection->getChart();
364: for(typename Section::atlas_type::chart_type::const_iterator c_iter = newChart.begin(); c_iter != newChart.end(); ++c_iter) {
365: newSection->updatePointAll(*c_iter, oldSection->restrictPoint(*c_iter));
366: }
367: };
370: template<typename RecvSection, typename Section>
371: static void updateArrowSectionRemote(const Obj<recv_overlap_type>& recvOverlap, const Obj<RecvSection>& recvSection, const Obj<bundle_type>& newBundle, const Obj<Section>& newSection) {
372: Obj<typename recv_overlap_type::traits::baseSequence> recvPoints = recvOverlap->base();
374: for(typename recv_overlap_type::traits::baseSequence::iterator r_iter = recvPoints->begin(); r_iter != recvPoints->end(); ++r_iter) {
375: const Obj<typename bundle_type::sieve_type::traits::coneSequence>& cone = newBundle->getSieve()->cone(*r_iter);
376: const typename bundle_type::sieve_type::traits::coneSequence::iterator end = cone->end();
378: for(typename bundle_type::sieve_type::traits::coneSequence::iterator c_iter = cone->begin(); c_iter != end; ++c_iter) {
379: newSection->setFiberDimension(typename Section::point_type(*c_iter, *r_iter), 1);
380: }
381: }
382: //newBundle->reallocate(newSection);
383: for(typename recv_overlap_type::traits::baseSequence::iterator r_iter = recvPoints->begin(); r_iter != recvPoints->end(); ++r_iter) {
384: const Obj<typename recv_overlap_type::traits::coneSequence>& recvPatches = recvOverlap->cone(*r_iter);
385: const typename recv_overlap_type::traits::coneSequence::iterator recvEnd = recvPatches->end();
387: for(typename recv_overlap_type::traits::coneSequence::iterator p_iter = recvPatches->begin(); p_iter != recvEnd; ++p_iter) {
388: const Obj<typename RecvSection::section_type>& section = recvSection->getSection(*p_iter);
390: if (section->getFiberDimension(*r_iter)) {
391: const Obj<typename bundle_type::sieve_type::traits::coneSequence>& cone = newBundle->getSieve()->cone(*r_iter);
392: const typename bundle_type::sieve_type::traits::coneSequence::iterator end = cone->end();
393: const typename RecvSection::value_type *values = section->restrictPoint(*r_iter);
394: int c = -1;
396: for(typename bundle_type::sieve_type::traits::coneSequence::iterator c_iter = cone->begin(); c_iter != end; ++c_iter) {
397: newSection->updatePoint(typename Section::point_type(*c_iter, *r_iter), &values[++c]);
398: }
399: }
400: }
401: }
402: };
405: template<typename Section>
406: static Obj<Section> distributeArrowSection(const Obj<Section>& serialSection, const Obj<bundle_type>& serialBundle, const Obj<bundle_type>& parallelBundle, const Obj<send_overlap_type>& sendOverlap, const Obj<recv_overlap_type>& recvOverlap) {
407: if (serialSection->debug()) {
408: serialSection->view("Serial ArrowSection");
409: }
410: typedef ALE::Field<send_overlap_type, int, ALE::Section<point_type, typename Section::value_type> > send_section_type;
411: typedef ALE::Field<recv_overlap_type, int, ALE::Section<point_type, typename Section::value_type> > recv_section_type;
412: typedef ALE::New::ConeSizeSection<bundle_type, sieve_type> SectionSizer;
413: typedef ALE::New::ArrowSection<sieve_type, Section> ArrowFiller;
414: Obj<Section> parallelSection = new Section(serialSection->comm(), serialSection->debug());
415: const Obj<send_section_type> sendSection = new send_section_type(serialSection->comm(), serialSection->debug());
416: const Obj<recv_section_type> recvSection = new recv_section_type(serialSection->comm(), sendSection->getTag(), serialSection->debug());
417: const Obj<SectionSizer> sizer = new SectionSizer(serialBundle, serialBundle->getSieve());
418: const Obj<ArrowFiller> filler = new ArrowFiller(serialBundle->getSieve(), serialSection);
420: updateArrowSectionLocal(serialSection, parallelBundle, parallelSection);
421: sectionCompletion::completeSection(sendOverlap, recvOverlap, sizer, filler, sendSection, recvSection);
422: updateArrowSectionRemote(recvOverlap, recvSection, parallelBundle, parallelSection);
423: if (parallelSection->debug()) {
424: parallelSection->view("Parallel ArrowSection");
425: }
426: return parallelSection;
427: };
428: static void unifyBundle(const Obj<bundle_type>& bundle, const int dim, const Obj<bundle_type>& bundleNew, const Obj<send_overlap_type>& sendOverlap, const Obj<recv_overlap_type>& recvOverlap) {
429: typedef int part_type;
430: const Obj<sieve_type>& sieve = bundle->getSieve();
431: const Obj<sieve_type>& sieveNew = bundleNew->getSieve();
432: const int rank = bundle->commRank();
433: const int debug = bundle->debug();
435: // 1) Form partition point overlap a priori
436: if (rank == 0) {
437: for(int p = 1; p < sieve->commSize(); p++) {
438: // The arrow is from remote partition point 0 on rank p to local partition point 0
439: recvOverlap->addCone(p, 0, 0);
440: }
441: } else {
442: // The arrow is from local partition point 0 to remote partition point 0 on rank 0
443: sendOverlap->addCone(0, 0, 0);
444: }
445: if (debug) {
446: sendOverlap->view("Send overlap for partition");
447: recvOverlap->view("Receive overlap for partition");
448: }
449: // 2) Partition the mesh
450: int numCells = bundle->heightStratum(0)->size();
451: part_type *assignment = new part_type[numCells];
453: for(int c = 0; c < numCells; ++c) {
454: assignment[c] = 0;
455: }
456: // 3) Scatter the sieve
457: sieveCompletion::scatterSieve(bundle, sieve, dim, sieveNew, sendOverlap, recvOverlap, 0, numCells, assignment);
458: bundleNew->stratify();
459: // 4) Cleanup
460: if (assignment != NULL) delete [] assignment;
461: };
464: static Obj<ALE::Mesh> unifyMesh(const Obj<ALE::Mesh>& parallelMesh) {
465: const int dim = parallelMesh->getDimension();
466: Obj<ALE::Mesh> serialMesh = new ALE::Mesh(parallelMesh->comm(), dim, parallelMesh->debug());
467: const Obj<ALE::Mesh::sieve_type>& serialSieve = new ALE::Mesh::sieve_type(parallelMesh->comm(), parallelMesh->debug());
469: ALE_LOG_EVENT_BEGIN;
470: serialMesh->setSieve(serialSieve);
471: if (parallelMesh->debug()) {
472: parallelMesh->view("Parallel topology");
473: }
475: // Unify cones
476: Obj<send_overlap_type> sendOverlap = new send_overlap_type(serialMesh->comm(), serialMesh->debug());
477: Obj<recv_overlap_type> recvOverlap = new recv_overlap_type(serialMesh->comm(), serialMesh->debug());
478: unifyBundle(parallelMesh, dim, serialMesh, sendOverlap, recvOverlap);
479: serialMesh->setDistSendOverlap(sendOverlap);
480: serialMesh->setDistRecvOverlap(recvOverlap);
482: // Unify labels
483: const typename bundle_type::labels_type& labels = parallelMesh->getLabels();
485: for(typename bundle_type::labels_type::const_iterator l_iter = labels.begin(); l_iter != labels.end(); ++l_iter) {
486: if (serialMesh->hasLabel(l_iter->first)) continue;
487: const Obj<typename bundle_type::label_type>& parallelLabel = l_iter->second;
488: const Obj<typename bundle_type::label_type>& serialLabel = serialMesh->createLabel(l_iter->first);
490: sieveCompletion::scatterCones(parallelLabel, serialLabel, sendOverlap, recvOverlap);
491: }
493: // Unify coordinates
494: Obj<std::set<std::string> > sections = parallelMesh->getRealSections();
496: for(std::set<std::string>::iterator name = sections->begin(); name != sections->end(); ++name) {
497: serialMesh->setRealSection(*name, distributeSection(parallelMesh->getRealSection(*name), serialMesh, sendOverlap, recvOverlap));
498: }
499: sections = parallelMesh->getIntSections();
500: for(std::set<std::string>::iterator name = sections->begin(); name != sections->end(); ++name) {
501: serialMesh->setIntSection(*name, distributeSection(parallelMesh->getIntSection(*name), serialMesh, sendOverlap, recvOverlap));
502: }
503: sections = parallelMesh->getArrowSections();
504: for(std::set<std::string>::iterator name = sections->begin(); name != sections->end(); ++name) {
505: serialMesh->setArrowSection(*name, distributeArrowSection(parallelMesh->getArrowSection(*name), parallelMesh, serialMesh, sendOverlap, recvOverlap));
506: }
507: if (serialMesh->debug()) {serialMesh->view("Serial Mesh");}
508: ALE_LOG_EVENT_END;
509: return serialMesh;
510: };
511: public: // Do not like these
514: // This is just crappy. We could introduce another phase to find out exactly what
515: // indices people do not have in the global order after communication
516: template<typename SendSection, typename RecvSection>
517: static void updateOverlap(const Obj<SendSection>& sendSection, const Obj<RecvSection>& recvSection, const Obj<send_overlap_type>& sendOverlap, const Obj<recv_overlap_type>& recvOverlap) {
518: const typename SendSection::sheaf_type& sendRanks = sendSection->getPatches();
519: const typename RecvSection::sheaf_type& recvRanks = recvSection->getPatches();
521: for(typename SendSection::sheaf_type::const_iterator p_iter = sendRanks.begin(); p_iter != sendRanks.end(); ++p_iter) {
522: const typename SendSection::patch_type& rank = p_iter->first;
523: const Obj<typename SendSection::section_type>& section = p_iter->second;
524: const typename SendSection::section_type::chart_type& chart = section->getChart();
526: for(typename SendSection::section_type::chart_type::iterator b_iter = chart.begin(); b_iter != chart.end(); ++b_iter) {
527: const typename SendSection::value_type *points = section->restrictPoint(*b_iter);
528: const int size = section->getFiberDimension(*b_iter);
530: for(int p = 0; p < size; p++) {
531: sendOverlap->addArrow(points[p], rank, points[p]);
532: }
533: }
534: }
535: for(typename RecvSection::sheaf_type::const_iterator p_iter = recvRanks.begin(); p_iter != recvRanks.end(); ++p_iter) {
536: const typename RecvSection::patch_type& rank = p_iter->first;
537: const Obj<typename RecvSection::section_type>& section = p_iter->second;
538: const typename RecvSection::section_type::chart_type& chart = section->getChart();
540: for(typename RecvSection::section_type::chart_type::iterator b_iter = chart.begin(); b_iter != chart.end(); ++b_iter) {
541: const typename RecvSection::value_type *points = section->restrictPoint(*b_iter);
542: const int size = section->getFiberDimension(*b_iter);
544: for(int p = 0; p < size; p++) {
545: recvOverlap->addArrow(rank, points[p], points[p]);
546: }
547: }
548: }
549: };
552: template<typename RecvSection>
553: static void updateSieve(const Obj<RecvSection>& recvSection, const Obj<sieve_type>& sieve) {
554: const typename RecvSection::sheaf_type& ranks = recvSection->getPatches();
556: for(typename RecvSection::sheaf_type::const_iterator p_iter = ranks.begin(); p_iter != ranks.end(); ++p_iter) {
557: const Obj<typename RecvSection::section_type>& section = p_iter->second;
558: const typename RecvSection::section_type::chart_type& chart = section->getChart();
560: for(typename RecvSection::section_type::chart_type::iterator b_iter = chart.begin(); b_iter != chart.end(); ++b_iter) {
561: const typename RecvSection::value_type *points = section->restrictPoint(*b_iter);
562: int size = section->getFiberDimension(*b_iter);
563: int c = 0;
565: for(int p = 0; p < size; p++) {
566: //sieve->addArrow(points[p], *b_iter, c++);
567: sieve->addArrow(points[p], *b_iter, c);
568: }
569: }
570: }
571: };
574: template<typename SendSection, typename RecvSection>
575: static void coneCompletion(const Obj<send_overlap_type>& sendOverlap, const Obj<recv_overlap_type>& recvOverlap, const Obj<bundle_type>& bundle, const Obj<SendSection>& sendSection, const Obj<RecvSection>& recvSection) {
576: if (sendOverlap->commSize() == 1) return;
577: // Distribute cones
578: const Obj<sieve_type>& sieve = bundle->getSieve();
579: const Obj<typename sieveCompletion::topology_type> secTopology = sieveCompletion::completion::createSendTopology(sendOverlap);
580: const Obj<typename sieveCompletion::cone_size_section> coneSizeSection = new typename sieveCompletion::cone_size_section(bundle, sieve);
581: const Obj<typename sieveCompletion::cone_section> coneSection = new typename sieveCompletion::cone_section(sieve);
582: sieveCompletion::completion::completeSection(sendOverlap, recvOverlap, coneSizeSection, coneSection, sendSection, recvSection);
583: // Update cones
584: updateSieve(recvSection, sieve);
585: };
588: template<typename Section>
589: static void completeSection(const Obj<bundle_type>& bundle, const Obj<Section>& section) {
590: typedef typename Distribution<bundle_type>::sieveCompletion sieveCompletion;
591: typedef typename bundle_type::send_overlap_type send_overlap_type;
592: typedef typename bundle_type::recv_overlap_type recv_overlap_type;
593: typedef typename Section::value_type value_type;
594: typedef typename ALE::Field<send_overlap_type, int, ALE::Section<point_type, value_type> > send_section_type;
595: typedef typename ALE::Field<recv_overlap_type, int, ALE::Section<point_type, value_type> > recv_section_type;
596: typedef ALE::New::SizeSection<Section> SectionSizer;
597: const int debug = section->debug();
599: bundle->constructOverlap();
600: const Obj<send_overlap_type> sendOverlap = bundle->getSendOverlap();
601: const Obj<recv_overlap_type> recvOverlap = bundle->getRecvOverlap();
602: const Obj<send_section_type> sendSection = new send_section_type(section->comm(), section->debug());
603: const Obj<recv_section_type> recvSection = new recv_section_type(section->comm(), sendSection->getTag(), section->debug());
604: const Obj<SectionSizer> sizer = new SectionSizer(section);
606: sectionCompletion::completeSection(sendOverlap, recvOverlap, sizer, section, sendSection, recvSection);
607: // Update section with remote data
608: const Obj<typename recv_overlap_type::traits::baseSequence> recvPoints = bundle->getRecvOverlap()->base();
610: for(typename recv_overlap_type::traits::baseSequence::iterator r_iter = recvPoints->begin(); r_iter != recvPoints->end(); ++r_iter) {
611: const Obj<typename recv_overlap_type::traits::coneSequence>& recvPatches = recvOverlap->cone(*r_iter);
612: const typename recv_overlap_type::traits::coneSequence::iterator end = recvPatches->end();
614: for(typename recv_overlap_type::traits::coneSequence::iterator p_iter = recvPatches->begin(); p_iter != end; ++p_iter) {
615: if (recvSection->getSection(*p_iter)->getFiberDimension(p_iter.color())) {
616: if (debug) {std::cout << "["<<section->commRank()<<"]Completed point " << *r_iter << std::endl;}
617: section->updateAddPoint(*r_iter, recvSection->getSection(*p_iter)->restrictPoint(p_iter.color()));
618: }
619: }
620: }
621: };
622: };
623: }
625: namespace ALECompat {
626: namespace New {
627: template<typename Topology_>
628: class Distribution {
629: public:
630: typedef Topology_ topology_type;
631: typedef typename topology_type::sieve_type sieve_type;
632: typedef ALECompat::New::Completion<Topology_, Mesh::sieve_type::point_type> sieveCompletion;
633: typedef ALECompat::New::SectionCompletion<Topology_, Mesh::real_section_type::value_type> sectionCompletion;
634: typedef typename sectionCompletion::send_overlap_type send_overlap_type;
635: typedef typename sectionCompletion::recv_overlap_type recv_overlap_type;
636: public:
639: // This is just crappy. WE could introduce another phase to find out exactly what
640: // indices people do not have in the global order after communication
641: template<typename SendSection, typename RecvSection>
642: static void updateOverlap(const Obj<SendSection>& sendSection, const Obj<RecvSection>& recvSection, const Obj<send_overlap_type>& sendOverlap, const Obj<recv_overlap_type>& recvOverlap) {
643: const typename SendSection::topology_type::sheaf_type& sendRanks = sendSection->getTopology()->getPatches();
644: const typename RecvSection::topology_type::sheaf_type& recvRanks = recvSection->getTopology()->getPatches();
646: for(typename SendSection::topology_type::sheaf_type::const_iterator p_iter = sendRanks.begin(); p_iter != sendRanks.end(); ++p_iter) {
647: int rank = p_iter->first;
648: const Obj<typename SendSection::topology_type::sieve_type::baseSequence>& base = p_iter->second->base();
650: for(typename SendSection::topology_type::sieve_type::baseSequence::iterator b_iter = base->begin(); b_iter != base->end(); ++b_iter) {
651: const typename SendSection::value_type *points = sendSection->restrict(rank, *b_iter);
652: int size = sendSection->getFiberDimension(rank, *b_iter);
654: for(int p = 0; p < size; p++) {
655: sendOverlap->addArrow(points[p], rank, points[p]);
656: }
657: }
658: }
659: for(typename RecvSection::topology_type::sheaf_type::const_iterator p_iter = recvRanks.begin(); p_iter != recvRanks.end(); ++p_iter) {
660: int rank = p_iter->first;
661: const Obj<typename RecvSection::topology_type::sieve_type::baseSequence>& base = p_iter->second->base();
663: for(typename RecvSection::topology_type::sieve_type::baseSequence::iterator b_iter = base->begin(); b_iter != base->end(); ++b_iter) {
664: const typename RecvSection::value_type *points = recvSection->restrict(rank, *b_iter);
665: int size = recvSection->getFiberDimension(rank, *b_iter);
667: for(int p = 0; p < size; p++) {
668: recvOverlap->addArrow(rank, points[p], points[p]);
669: }
670: }
671: }
672: };
673: static void createLabelOverlap(const Obj<send_overlap_type>& sendOverlap, const Obj<recv_overlap_type>& recvOverlap) {
674: };
677: template<typename Section>
678: static void updateSectionLocal(const Obj<Section>& oldSection, const Obj<Section>& newSection)
679: {
680: const typename Section::topology_type::sheaf_type& patches = newSection->getTopology()->getPatches();
682: for(typename Section::topology_type::sheaf_type::const_iterator p_iter = patches.begin(); p_iter != patches.end(); ++p_iter) {
683: const typename Section::patch_type& patch = p_iter->first;
684: const Obj<typename Section::topology_type::sieve_type>& newSieve = p_iter->second;
685: if (!oldSection->hasPatch(patch)) continue;
686: const typename Section::atlas_type::chart_type& oldChart = oldSection->getPatch(patch);
688: for(typename Section::atlas_type::chart_type::const_iterator c_iter = oldChart.begin(); c_iter != oldChart.end(); ++c_iter) {
689: if (newSieve->hasPoint(*c_iter)) {
690: newSection->setFiberDimension(patch, *c_iter, oldSection->getFiberDimension(patch, *c_iter));
691: }
692: }
693: }
694: newSection->allocate();
695: for(typename Section::topology_type::sheaf_type::const_iterator p_iter = patches.begin(); p_iter != patches.end(); ++p_iter) {
696: const typename Section::patch_type& patch = p_iter->first;
697: if (!oldSection->hasPatch(patch)) continue;
698: const typename Section::atlas_type::chart_type& newChart = newSection->getPatch(patch);
700: for(typename Section::atlas_type::chart_type::const_iterator c_iter = newChart.begin(); c_iter != newChart.end(); ++c_iter) {
701: newSection->updatePoint(patch, *c_iter, oldSection->restrictPoint(patch, *c_iter));
702: }
703: }
704: };
707: template<typename RecvSection, typename Section>
708: static void updateSectionRemote(const Obj<recv_overlap_type>& recvOverlap, const Obj<RecvSection>& recvSection, const Obj<Section>& newSection) {
709: const Mesh::real_section_type::patch_type patch = 0;
710: Obj<typename recv_overlap_type::traits::baseSequence> recvPoints = recvOverlap->base();
712: for(typename recv_overlap_type::traits::baseSequence::iterator r_iter = recvPoints->begin(); r_iter != recvPoints->end(); ++r_iter) {
713: const Obj<typename recv_overlap_type::traits::coneSequence>& recvPatches = recvOverlap->cone(*r_iter);
714: const typename recv_overlap_type::traits::coneSequence::iterator end = recvPatches->end();
716: for(typename recv_overlap_type::traits::coneSequence::iterator p_iter = recvPatches->begin(); p_iter != end; ++p_iter) {
717: newSection->addPoint(patch, *r_iter, recvSection->getFiberDimension(*p_iter, *r_iter));
718: }
719: }
720: newSection->reallocate();
721: for(typename recv_overlap_type::traits::baseSequence::iterator r_iter = recvPoints->begin(); r_iter != recvPoints->end(); ++r_iter) {
722: const Obj<typename recv_overlap_type::traits::coneSequence>& recvPatches = recvOverlap->cone(*r_iter);
723: const typename recv_overlap_type::traits::coneSequence::iterator end = recvPatches->end();
725: for(typename recv_overlap_type::traits::coneSequence::iterator p_iter = recvPatches->begin(); p_iter != end; ++p_iter) {
726: if (recvSection->getFiberDimension(*p_iter, *r_iter)) {
727: newSection->updatePoint(patch, *r_iter, recvSection->restrictPoint(*p_iter, *r_iter));
728: }
729: }
730: }
731: };
734: template<typename RecvSection>
735: static void updateSieve(const Obj<RecvSection>& recvSection, const Obj<topology_type>& topology) {
736: const typename RecvSection::patch_type patch = 0;
737: const typename RecvSection::topology_type::sheaf_type& ranks = recvSection->getTopology()->getPatches();
738: const Obj<typename topology_type::sieve_type>& sieve = topology->getPatch(patch);
740: for(typename RecvSection::topology_type::sheaf_type::const_iterator p_iter = ranks.begin(); p_iter != ranks.end(); ++p_iter) {
741: int rank = p_iter->first;
742: const Obj<typename RecvSection::topology_type::sieve_type::baseSequence>& base = p_iter->second->base();
744: for(typename RecvSection::topology_type::sieve_type::baseSequence::iterator b_iter = base->begin(); b_iter != base->end(); ++b_iter) {
745: const typename RecvSection::value_type *points = recvSection->restrict(rank, *b_iter);
746: int size = recvSection->getFiberDimension(rank, *b_iter);
747: int c = 0;
749: for(int p = 0; p < size; p++) {
750: //sieve->addArrow(points[p], *b_iter, c++);
751: sieve->addArrow(points[p], *b_iter, c);
752: }
753: }
754: }
755: };
758: template<typename SendSection, typename RecvSection>
759: static void coneCompletion(const Obj<send_overlap_type>& sendOverlap, const Obj<recv_overlap_type>& recvOverlap, const Obj<topology_type>& topology, const Obj<SendSection>& sendSection, const Obj<RecvSection>& recvSection) {
760: if (sendOverlap->commSize() == 1) return;
761: typedef typename ALECompat::New::SectionCompletion<topology_type, typename sieve_type::point_type> completion;
763: // Distribute cones
764: const typename topology_type::patch_type patch = 0;
765: const Obj<typename sieveCompletion::topology_type> secTopology = completion::createSendTopology(sendOverlap);
766: const Obj<typename sieveCompletion::cone_size_section> coneSizeSection = new typename sieveCompletion::cone_size_section(secTopology, topology, topology->getPatch(patch));
767: const Obj<typename sieveCompletion::cone_section> coneSection = new typename sieveCompletion::cone_section(secTopology, topology->getPatch(patch));
768: completion::completeSection(sendOverlap, recvOverlap, coneSizeSection, coneSection, sendSection, recvSection);
769: // Update cones
770: updateSieve(recvSection, topology);
771: };
774: static void createPartitionOverlap(const Obj<topology_type>& topology, const Obj<send_overlap_type>& sendOverlap, const Obj<recv_overlap_type>& recvOverlap) {
775: const Obj<send_overlap_type>& topSendOverlap = topology->getSendOverlap();
776: const Obj<recv_overlap_type>& topRecvOverlap = topology->getRecvOverlap();
777: const Obj<typename send_overlap_type::traits::baseSequence> base = topSendOverlap->base();
778: const Obj<typename recv_overlap_type::traits::capSequence> cap = topRecvOverlap->cap();
780: if (base->empty()) {
781: if (topology->commRank() == 0) {
782: for(int p = 1; p < topology->commSize(); p++) {
783: // The arrow is from local partition point p (source) to remote partition point p (color) on rank p (target)
784: sendOverlap->addCone(p, p, p);
785: }
786: }
787: } else {
788: for(typename send_overlap_type::traits::baseSequence::iterator b_iter = base->begin(); b_iter != base->end(); ++b_iter) {
789: const int& p = *b_iter;
790: // The arrow is from local partition point p (source) to remote partition point p (color) on rank p (target)
791: sendOverlap->addCone(p, p, p);
792: }
793: }
794: if (cap->empty()) {
795: if (topology->commRank() != 0) {
796: // The arrow is from local partition point rank (color) on rank 0 (source) to remote partition point rank (target)
797: recvOverlap->addCone(0, topology->commRank(), topology->commRank());
798: }
799: } else {
800: for(typename recv_overlap_type::traits::capSequence::iterator c_iter = cap->begin(); c_iter != cap->end(); ++c_iter) {
801: const int& p = *c_iter;
802: // The arrow is from local partition point rank (color) on rank p (source) to remote partition point rank (target)
803: recvOverlap->addCone(p, topology->commRank(), topology->commRank());
804: }
805: }
806: if (topology->debug()) {
807: sendOverlap->view("Initial send overlap");
808: recvOverlap->view("Initial receive overlap");
809: }
810: }
813: template<typename Partitioner>
814: static typename Partitioner::part_type *createAssignment(const Obj<topology_type>& topology, const int dim, const Obj<send_overlap_type>& sendOverlap, const Obj<recv_overlap_type>& recvOverlap) {
815: // 1) Form partition point overlap a priori
816: createPartitionOverlap(topology, sendOverlap, recvOverlap);
817: if (topology->debug()) {
818: sendOverlap->view("Send overlap for partition");
819: recvOverlap->view("Receive overlap for partition");
820: }
821: // 2) Partition the mesh
822: return Partitioner::partitionSieve(topology, dim);
823: };
826: // Partition a topology on process 0 and scatter to all processes
827: static void scatterTopology(const Obj<topology_type>& topology, const int dim, const Obj<topology_type>& topologyNew, const Obj<send_overlap_type>& sendOverlap, const Obj<recv_overlap_type>& recvOverlap, const std::string& partitioner, const Obj<topology_type>& subTopology = NULL, const Obj<topology_type>& subTopologyNew = NULL) {
828: if (partitioner == "chaco") {
829: #ifdef PETSC_HAVE_CHACO
830: typedef typename ALECompat::New::Chaco::Partitioner<topology_type> Partitioner;
831: typedef typename ALECompat::New::Partitioner<topology_type> GenPartitioner;
832: typedef typename Partitioner::part_type part_type;
834: part_type *assignment = scatterTopology<Partitioner>(topology, dim, topologyNew, sendOverlap, recvOverlap);
835: if (!subTopology.isNull() && !subTopologyNew.isNull()) {
836: part_type *subAssignment = GenPartitioner::subordinatePartition(topology, 1, subTopology, assignment);
837: Obj<send_overlap_type> subSendOverlap = new send_overlap_type(topology->comm(), topology->debug());
838: Obj<recv_overlap_type> subRecvOverlap = new recv_overlap_type(topology->comm(), topology->debug());
839: const typename topology_type::patch_type patch = 0;
840: const Obj<sieve_type>& sieve = subTopology->getPatch(patch);
841: const Obj<sieve_type>& sieveNew = new Mesh::sieve_type(subTopology->comm(), subTopology->debug());
842: const int numCells = subTopology->heightStratum(patch, 0)->size();
844: createPartitionOverlap(subTopology, subSendOverlap, subRecvOverlap);
845: subTopologyNew->setPatch(0, sieveNew);
846: sieveCompletion::scatterSieve(subTopology, sieve, dim, sieveNew, subSendOverlap, subRecvOverlap, numCells, subAssignment);
847: subTopologyNew->stratify();
848: subTopologyNew->setDistSendOverlap(subSendOverlap);
849: subTopologyNew->setDistRecvOverlap(subRecvOverlap);
850: if (subAssignment != NULL) delete [] subAssignment;
851: }
852: if (assignment != NULL) delete [] assignment;
853: #else
854: throw ALE::Exception("Chaco is not installed. Reconfigure with the flag --download-chaco");
855: #endif
856: } else if (partitioner == "parmetis") {
857: #ifdef PETSC_HAVE_PARMETIS
858: typedef typename ALECompat::New::ParMetis::Partitioner<topology_type> Partitioner;
859: typedef typename Partitioner::part_type part_type;
861: part_type *assignment = scatterTopology<Partitioner>(topology, dim, topologyNew, sendOverlap, recvOverlap);
862: if (assignment != NULL) delete [] assignment;
863: #else
864: throw ALE::Exception("ParMetis is not installed. Reconfigure with the flag --download-parmetis");
865: #endif
866: } else {
867: throw ALE::Exception("Unknown partitioner");
868: }
869: };
870: template<typename Partitioner>
871: static typename Partitioner::part_type *scatterTopology(const Obj<topology_type>& topology, const int dim, const Obj<topology_type>& topologyNew, const Obj<send_overlap_type>& sendOverlap, const Obj<recv_overlap_type>& recvOverlap) {
872: typename Partitioner::part_type *assignment = createAssignment<Partitioner>(topology, dim, sendOverlap, recvOverlap);
873: const typename topology_type::patch_type patch = 0;
874: const Obj<sieve_type>& sieve = topology->getPatch(patch);
875: const Obj<sieve_type>& sieveNew = topologyNew->getPatch(patch);
876: const int numCells = topology->heightStratum(patch, 0)->size();
878: sieveCompletion::scatterSieve(topology, sieve, dim, sieveNew, sendOverlap, recvOverlap, numCells, assignment);
879: topologyNew->stratify();
880: return assignment;
881: };
884: static Obj<Mesh> distributeMesh(const Obj<Mesh>& serialMesh, const std::string& partitioner = "chaco") {
885: Obj<Mesh> parallelMesh = new Mesh(serialMesh->comm(), serialMesh->getDimension(), serialMesh->debug());
886: const Obj<Mesh::topology_type>& serialTopology = serialMesh->getTopology();
887: const Obj<Mesh::topology_type>& parallelTopology = new Mesh::topology_type(serialMesh->comm(), serialMesh->debug());
888: const Obj<Mesh::topology_type>& tractionTopology = new Mesh::topology_type(serialMesh->comm(), serialMesh->debug());
889: const Obj<Mesh::sieve_type>& sieve = new Mesh::sieve_type(serialMesh->comm(), serialMesh->debug());
890: const int dim = serialMesh->getDimension();
892: //if (serialMesh->getDistributed()) return serialMesh;
893: ALE_LOG_EVENT_BEGIN;
894: parallelTopology->setPatch(0, sieve);
895: parallelMesh->setTopology(parallelTopology);
896: if (serialMesh->debug()) {
897: serialMesh->view("Serial topology");
898: }
900: // Distribute cones
901: Obj<send_overlap_type> sendOverlap = new send_overlap_type(serialTopology->comm(), serialTopology->debug());
902: Obj<recv_overlap_type> recvOverlap = new recv_overlap_type(serialTopology->comm(), serialTopology->debug());
903: if (serialMesh->hasRealSection("traction")) {
904: const Obj<Mesh::real_section_type>& traction = serialMesh->getRealSection("traction");
906: scatterTopology(serialTopology, dim, parallelTopology, sendOverlap, recvOverlap, partitioner, traction->getTopology(), tractionTopology);
907: } else {
908: scatterTopology(serialTopology, dim, parallelTopology, sendOverlap, recvOverlap, partitioner);
909: }
910: parallelTopology->setDistSendOverlap(sendOverlap);
911: parallelTopology->setDistRecvOverlap(recvOverlap);
913: // Distribute labels
914: const typename topology_type::labels_type& labels = serialTopology->getLabels();
916: for(typename topology_type::labels_type::const_iterator l_iter = labels.begin(); l_iter != labels.end(); ++l_iter) {
917: for(typename topology_type::label_type::const_iterator pl_iter = l_iter->second.begin(); pl_iter != l_iter->second.end(); ++pl_iter) {
918: if (parallelTopology->hasLabel(l_iter->first, pl_iter->first)) continue;
919: const Obj<typename topology_type::patch_label_type>& serialLabel = pl_iter->second;
920: const Obj<typename topology_type::patch_label_type>& parallelLabel = parallelTopology->createLabel(pl_iter->first, l_iter->first);
921: // Create local label
922: const Obj<typename topology_type::patch_label_type::traits::baseSequence>& base = serialLabel->base();
923: const Obj<typename topology_type::sieve_type>& parallelSieve = parallelTopology->getPatch(pl_iter->first);
925: for(typename topology_type::patch_label_type::traits::baseSequence::iterator b_iter = base->begin(); b_iter != base->end(); ++b_iter) {
926: if (parallelSieve->capContains(*b_iter) || parallelSieve->baseContains(*b_iter)) {
927: parallelLabel->addArrow(*serialLabel->cone(*b_iter)->begin(), *b_iter);
928: }
929: }
930: // Get remote labels
931: sieveCompletion::scatterCones(serialLabel, parallelLabel, sendOverlap, recvOverlap);
932: }
933: }
935: // Distribute sections
936: Obj<std::set<std::string> > sections = serialMesh->getRealSections();
938: for(std::set<std::string>::iterator name = sections->begin(); name != sections->end(); ++name) {
939: if (*name == "traction") {
940: parallelMesh->setRealSection(*name, distributeSection(serialMesh->getRealSection(*name), tractionTopology, sendOverlap, recvOverlap));
941: } else {
942: parallelMesh->setRealSection(*name, distributeSection(serialMesh->getRealSection(*name), parallelMesh->getTopology(), sendOverlap, recvOverlap));
943: }
944: }
945: sections = serialMesh->getIntSections();
946: for(std::set<std::string>::iterator name = sections->begin(); name != sections->end(); ++name) {
947: parallelMesh->setIntSection(*name, distributeSection(serialMesh->getIntSection(*name), parallelMesh->getTopology(), sendOverlap, recvOverlap));
948: }
949: sections = serialMesh->getPairSections();
950: for(std::set<std::string>::iterator name = sections->begin(); name != sections->end(); ++name) {
951: parallelMesh->setPairSection(*name, distributeSection(serialMesh->getPairSection(*name), parallelMesh->getTopology(), sendOverlap, recvOverlap));
952: }
954: // This is necessary since we create types (like PartitionSection) on a subset of processors
955: if (parallelMesh->debug()) {parallelMesh->view("Parallel Mesh");}
956: parallelMesh->setDistributed(true);
957: ALE_LOG_EVENT_END;
958: return parallelMesh;
959: };
962: template<typename Section>
963: static Obj<Section> distributeSection(const Obj<Section>& serialSection, const Obj<Mesh::topology_type>& parallelTopology, const Obj<send_overlap_type>& sendOverlap, const Obj<recv_overlap_type>& recvOverlap) {
964: if (serialSection->debug()) {
965: serialSection->view("Serial Section");
966: }
967: typedef OverlapValues<send_overlap_type, typename sieveCompletion::topology_type, typename Section::value_type> send_section_type;
968: typedef OverlapValues<recv_overlap_type, typename sieveCompletion::topology_type, typename Section::value_type> recv_section_type;
969: typedef SizeSection<Section> SectionSizer;
970: // TEST THIS! I think this is unnecessary
971: typedef PatchlessSection<Section> SectionFiller;
972: Obj<Section> parallelSection = new Section(parallelTopology);
973: const typename Section::patch_type patch = 0;
974: const Obj<send_section_type> sendSection = new send_section_type(serialSection->comm(), serialSection->debug());
975: const Obj<recv_section_type> recvSection = new recv_section_type(serialSection->comm(), sendSection->getTag(), serialSection->debug());
976: const Obj<SectionSizer> sizer = new SectionSizer(serialSection, patch);
977: const Obj<SectionFiller> filler = new SectionFiller(serialSection, patch);
979: updateSectionLocal(serialSection, parallelSection);
980: sectionCompletion::completeSection(sendOverlap, recvOverlap, sizer, filler, sendSection, recvSection);
981: updateSectionRemote(recvOverlap, recvSection, parallelSection);
982: if (parallelSection->debug()) {
983: parallelSection->view("Parallel Section");
984: }
985: return parallelSection;
986: };
989: template<typename Section>
990: static void completeSection(const Obj<Section>& section) {
991: typedef typename Section::topology_type topology_type;
992: typedef typename Distribution<topology_type>::sieveCompletion sieveCompletion;
993: typedef typename topology_type::send_overlap_type send_overlap_type;
994: typedef typename topology_type::recv_overlap_type recv_overlap_type;
995: typedef typename Section::value_type value_type;
996: typedef OverlapValues<send_overlap_type, typename sieveCompletion::topology_type, value_type> send_section_type;
997: typedef OverlapValues<recv_overlap_type, typename sieveCompletion::topology_type, value_type> recv_section_type;
998: typedef SizeSection<Section> SectionSizer;
999: typedef PatchlessSection<Section> SectionFiller;
1000: const Obj<topology_type>& topology = section->getTopology();
1001: const typename topology_type::patch_type patch = 0;
1002: const int debug = section->debug();
1003: topology->constructOverlap(patch);
1005: const Obj<send_overlap_type> sendOverlap = topology->getSendOverlap();
1006: const Obj<recv_overlap_type> recvOverlap = topology->getRecvOverlap();
1007: const Obj<send_section_type> sendSection = new send_section_type(section->comm(), section->debug());
1008: const Obj<recv_section_type> recvSection = new recv_section_type(section->comm(), sendSection->getTag(), section->debug());
1009: const Obj<SectionSizer> sizer = new SectionSizer(section, patch);
1010: const Obj<SectionFiller> filler = new SectionFiller(section, patch);
1012: sectionCompletion::completeSection(sendOverlap, recvOverlap, sizer, filler, sendSection, recvSection);
1013: // Update section with remote data
1014: const Obj<typename recv_overlap_type::traits::baseSequence> recvPoints = topology->getRecvOverlap()->base();
1016: for(typename recv_overlap_type::traits::baseSequence::iterator r_iter = recvPoints->begin(); r_iter != recvPoints->end(); ++r_iter) {
1017: const Obj<typename recv_overlap_type::traits::coneSequence>& recvPatches = recvOverlap->cone(*r_iter);
1018: const typename recv_overlap_type::traits::coneSequence::iterator end = recvPatches->end();
1020: for(typename recv_overlap_type::traits::coneSequence::iterator p_iter = recvPatches->begin(); p_iter != end; ++p_iter) {
1021: if (recvSection->getFiberDimension(*p_iter, p_iter.color())) {
1022: if (debug) {std::cout << "["<<section->commRank()<<"]Completed point " << *r_iter << std::endl;}
1023: section->updateAddPoint(patch, *r_iter, recvSection->restrictPoint(*p_iter, p_iter.color()));
1024: }
1025: }
1026: }
1027: };
1028: };
1029: }
1030: }
1031: #endif