FEAT 3
Finite Element Analysis Toolbox
Loading...
Searching...
No Matches
parti_domain_control_base.hpp
1// FEAT3: Finite Element Analysis Toolbox, Version 3
2// Copyright (C) 2010 by Stefan Turek & the FEAT group
3// FEAT3 is released under the GNU General Public License version 3,
4// see the file 'copyright.txt' in the top level directory for details.
5
6#pragma once
7
9
10#include <kernel/util/dist.hpp>
11#include <kernel/util/dist_file_io.hpp>
12#include <kernel/runtime.hpp>
13#include <kernel/util/simple_arg_parser.hpp>
14#include <kernel/util/property_map.hpp>
15#include <kernel/util/statistics.hpp>
16#include <kernel/geometry/mesh_node.hpp>
17#include <kernel/geometry/partition_set.hpp>
18#include <kernel/geometry/parti_2lvl.hpp>
19#include <kernel/geometry/parti_iterative.hpp>
20#include <kernel/geometry/parti_parmetis.hpp>
21#include <kernel/geometry/parti_zoltan.hpp>
22
23#include <control/domain/domain_control.hpp>
24
25namespace FEAT
26{
27 namespace Control
28 {
29 namespace Domain
30 {
40 template<typename DomainLevel_>
42 public Control::Domain::DomainControl<DomainLevel_>
43 {
44 public:
48 using typename BaseClass::LevelType;
50 using typename BaseClass::LayerType;
52 using typename BaseClass::MeshType;
54 using typename BaseClass::AtlasType;
58 typedef typename LevelType::PartType MeshPartType;
59
62
67 {
68 public:
70 int layer;
78 int desired_level_max, desired_level_min;
79
80 int progeny_group, progeny_child;
81 int progeny_first, progeny_count;
82 Dist::Comm progeny_comm;
83
94
95 Ancestor() :
96 layer(0), layer_p(0), num_procs(0), num_parts(0),
97 desired_level_max(0), desired_level_min(0),
98 progeny_group(0), progeny_child(0),
99 progeny_first(0), progeny_count(0),
100 progeny_comm(),
101 parti_found(false),
102 parti_info(),
103 parti_level(0),
104 parti_apriori(false),
106 {
107 }
108 }; // class Ancestor
109
110 protected:
115
124
127
129 std::deque<std::pair<int,int>> _desired_levels;
130
132 std::deque<std::pair<int,int>> _chosen_levels;
133
140
142 std::deque<Ancestor> _ancestry;
143
144 public:
154 explicit PartiDomainControlBase(const Dist::Comm& comm_, bool support_multi_layered) :
155 BaseClass(comm_),
156 _was_created(false),
157 _adapt_mode(Geometry::AdaptMode::chart),
158 _allow_parti_metis(false),
159 _allow_parti_genetic(false), // this one is exotic
160 _allow_parti_naive(true),
161 _allow_parti_zoltan(false),
162 _support_multi_layered(support_multi_layered),
167 _ancestry()
168 {
169 }
170
173 {
174 }
175
188 virtual bool parse_args(SimpleArgParser& args)
189 {
190 XASSERTM(!_was_created, "This function has to be called before domain control creation!");
191 // Try to parse --parti-type <types...>
192 {
193 auto it = args.query("parti-type");
194 if(it != nullptr)
195 {
196 // disallow all strategies by default
197 this->_reset_parti_types();
198
199 // loop over all allowed strategies
200 for(const auto& t : it->second)
201 {
202 // give derived classes a chance first
203 if(!this->_parse_parti_type(t))
204 {
205 this->_comm.print("ERROR: unknown partitioner type '" + t + "'");
206 return false;
207 }
208 }
209 }
210 }
211
212 // parse --parti-rank-elems <count>
213 args.parse("parti-rank-elems", _required_elems_per_rank);
214
215 // parse --parti-genetic-time <time-init> <time-mutate>
216 args.parse("parti-genetic-time", _genetic_time_init, _genetic_time_mutate);
217
218 // okay
219 return true;
220 }
221
232 virtual bool parse_property_map(const PropertyMap& pmap)
233 {
234 XASSERTM(!_was_created, "This function has to be called before domain control creation!");
235 auto parti_type_p = pmap.query("parti-type");
236 if(parti_type_p.second)
237 {
238 // disallow all strategies by default
239 this->_reset_parti_types();
240
241 std::deque<String> allowed_partitioners = parti_type_p.first.split_by_whitespaces();
242
243 for(const auto& t : allowed_partitioners)
244 {
245 if(!this->_parse_parti_type(t))
246 {
247 this->_comm.print("ERROR: unknown partitioner type '" + t + "'");
248 return false;
249 }
250 }
251 }
252
253 auto parti_rank_elems_p = pmap.query("parti-rank-elems");
254 if(parti_rank_elems_p.second)
255 {
256 if(!parti_rank_elems_p.first.parse(_required_elems_per_rank))
257 {
258 this->_comm.print("ERROR: Failed to parse 'parti-rank-elems'");
259 return false;
260 }
261 }
262
263 auto genetic_time_init_p = pmap.query("parti-genetic-time-init");
264 if(genetic_time_init_p.second)
265 {
266 if(!genetic_time_init_p.first.parse(_genetic_time_init))
267 {
268 this->_comm.print("ERROR: Failed to parse 'parti-genetic-time-init'");
269 return false;
270 }
271 }
272
273 auto genetic_time_mutate_p = pmap.query("parti-genetic-time-mutate");
274 if(genetic_time_mutate_p.second)
275 {
276 if(!genetic_time_mutate_p.first.parse(_genetic_time_mutate))
277 {
278 this->_comm.print("ERROR: Failed to parse 'parti-genetic-time-mutate'");
279 return false;
280 }
281 }
282
283 return true;
284 }
285
293 {
294 XASSERTM(!_was_created, "This function has to be called before domain control creation!");
295 _adapt_mode = adapt_mode;
296 }
297
305 {
306 return _adapt_mode;
307 }
308
315 void set_desired_levels(const String& slvls)
316 {
317 std::deque<String> sl = slvls.split_by_whitespaces();
319 }
320
329 virtual void set_desired_levels(const std::deque<String>& slvls)
330 {
331 XASSERTM(!_was_created, "This function has to be called before domain control creation!");
332
333 const int nranks = this->_comm.size();
334 std::deque<String> sv;
335
336 for(std::size_t i(0); i < slvls.size(); ++i)
337 {
338 int ilvl = -1, nprocs = -1;
339 sv = slvls.at(i).split_by_string(":");
340
341 if((sv.size() < std::size_t(1)) || (sv.size() > std::size_t(2)))
342 throw ParseError("Invalid input format", slvls.at(i), "an int-pair 'level:patches'");
343
344 if(!sv.front().parse(ilvl))
345 throw ParseError("Failed to parse level index", slvls.at(i), "an integer");
346 if((sv.size() > std::size_t(1)) && !sv.back().parse(nprocs))
347 throw ParseError("Failed to parse process count" , slvls.at(i), "an integer");
348
349 // level must be non-negative
350 if(ilvl < 0)
351 throw ParseError("Invalid negative level index", slvls.at(i), "a non-negative level index");
352
353 // first level index?
354 if(i == std::size_t(0))
355 {
356 // if the process count is given, it must match the communicator size
357 if((nprocs >= 0) && (nprocs != this->_comm.size()))
358 throw ParseError("Invalid number of processes for global level: '" + slvls.at(i) +
359 "', expected " + stringify(this->_comm.size()) + " but got " + stringify(nprocs));
360 _desired_levels.push_back(std::make_pair(ilvl, nranks));
361 continue;
362 }
363
364 // make sure the level parameter is non-ascending
365 if(_desired_levels.back().first < ilvl)
366 throw ParseError("Invalid non-descending level index: '" + slvls.at(i) +
367 "', expected <= " + stringify(_desired_levels.back().first) + " but got " + stringify(ilvl));
368
369 // ignore process count for last level
370 if((i + 1) == slvls.size())
371 nprocs = 0;
372 // if the number of processes is identical to the previous one, simply ignore this parameter
373 else if(_desired_levels.back().second == nprocs)
374 continue;
375 // the process count must not be ascending
376 else if(_desired_levels.back().second < nprocs)
377 throw ParseError("Invalid non-descending process count: '" + slvls.at(i) +
378 "', expected < " + stringify(_desired_levels.back().second) + " but got " + stringify(nprocs));
379 // the previous process count must be a multiple of the new one
380 else if(_desired_levels.back().second % nprocs != 0)
381 throw ParseError("Invalid indivisible process count: '" + slvls.at(i) +
382 "', expected a divisor of " + stringify(_desired_levels.back().second) + " but got " + stringify(nprocs));
383
384 // push the level-proc pair
385 _desired_levels.push_back(std::make_pair(ilvl, nprocs));
386 }
387 }
388
399 virtual void set_desired_levels(int lvl_max, int lvl_min = -1)
400 {
401 XASSERTM(!_was_created, "This function has to be called before domain control creation!");
402
403 // single-layered hierarchy
404 int _desired_level_max = lvl_max;
405 int _desired_level_min = (lvl_min >= 0 ? lvl_min : lvl_max + lvl_min + 1);
406
407 XASSERTM(_desired_level_max >= 0, "Invalid level-max");
408 XASSERTM(_desired_level_min >= 0, "Invalid level-min");
409 XASSERTM(_desired_level_max >= _desired_level_min, "Invalid level-min/max combination");
410
411 _desired_levels.emplace_back(std::make_pair(_desired_level_max, this->_comm.size()));
412 _desired_levels.emplace_back(std::make_pair(_desired_level_min, 0));
413 }
414
429 virtual void set_desired_levels(int lvl_max, int lvl_med, int lvl_min)
430 {
431 XASSERTM(!_was_created, "This function has to be called before domain control creation!");
432
433 // double-layered hierarchy
434 int _desired_level_max = lvl_max;
435 int _desired_level_med = (lvl_med >= 0 ? lvl_med : lvl_max + lvl_med + 1);
436 int _desired_level_min = (lvl_min >= 0 ? lvl_min : lvl_med + lvl_min + 1);
437
438 XASSERTM(_desired_level_max >= 0, "Invalid level-max");
439 XASSERTM(_desired_level_med >= 0, "Invalid level-med");
440 XASSERTM(_desired_level_min >= 0, "Invalid level-min");
441
442 // level_max must be strictly greater than level_med
443 XASSERTM(_desired_level_max > _desired_level_med, "Invalid level-max/med combination");
444
445 // level_med must be greater or equal level_min
446 XASSERTM(_desired_level_med >= _desired_level_min, "Invalid level-med/min combination");
447
448 _desired_levels.emplace_back(std::make_pair(_desired_level_max, this->_comm.size()));
449 _desired_levels.emplace_back(std::make_pair(_desired_level_med, 1));
450 _desired_levels.emplace_back(std::make_pair(_desired_level_min, 0));
451 }
452
459 {
460 return _desired_levels.front().first;
461 }
462
469 {
470 return _desired_levels.back().first;
471 }
472
481 {
482 String s;
483 for(std::size_t i(0); (i + 1) < _desired_levels.size(); ++i)
484 {
485 s += stringify(_desired_levels.at(i).first);
486 s += ":";
487 s += stringify(_desired_levels.at(i).second);
488 s += " ";
489 }
490 s += stringify(_desired_levels.back().first);
491 return s;
492 }
493
502 {
503 String s;
504 for(std::size_t i(0); (i + 1) < _chosen_levels.size(); ++i)
505 {
506 s += stringify(_chosen_levels.at(i).first);
507 s += ":";
508 s += stringify(_chosen_levels.at(i).second);
509 s += " ";
510 }
511 s += stringify(_chosen_levels.back().first);
512 return s;
513 }
514
528 {
529 String s;
530
531 for(auto it = this->_ancestry.rbegin(); it != this->_ancestry.rend(); ++it)
532 {
533 if(it != this->_ancestry.rbegin())
534 s += "\n";
535 s += it->parti_info;
536 s += " on level ";
537 s += stringify(it->parti_level);
538 s += " for ";
539 s += stringify(it->num_parts);
540 s += (it->num_parts > 1 ? " patches on " : " patch on ");
541 s += stringify(it->num_procs);
542 s += (it->num_procs > 1 ? " processes" : " process");
543 }
544
545 return s;
546 }
547
553 const std::deque<std::pair<int, int>> get_desired_levels() const
554 {
555 return this->_desired_levels;
556 }
557
563 const std::deque<std::pair<int, int>> get_chosen_levels() const
564 {
565 return this->_chosen_levels;
566 }
567
572 {
573 String s;
574 for(auto it = this->_ancestry.begin(); it != this->_ancestry.end(); ++it)
575 {
576 if(it != this->_ancestry.begin())
577 s += " | ";
578 s += stringify((*it).num_procs);
579 s += ":";
580 s += stringify((*it).num_parts);
581 s += "[";
582 s += stringify((*it).progeny_group).pad_front(2);
583 s += "+";
584 s += stringify((*it).progeny_child);
585 s += ":";
586 s += stringify((*it).progeny_first).pad_front(2);
587 s += "+";
588 s += stringify((*it).progeny_count);
589 s += ".";
590 s += stringify((*it).progeny_comm.rank());
591 s += "]";
592 s += "(";
593 s += stringify((*it).desired_level_max);
594 s += ":";
595 s += stringify((*it).desired_level_min);
596 s += ")";
597 s += ((*it).layer >= 0 ? "*" : " ");
598 s += ((*it).layer_p >= 0 ? ">" : " ");
599 }
600 return s;
601 }
602
606 const std::deque<Ancestor>& get_ancestry() const
607 {
608 return this->_ancestry;
609 }
610
614 std::vector<char> serialize_partitioning() const
615 {
616 BinaryStream bs;
617 typedef std::uint64_t u64;
618
619 // serialization magic number: "F3PaDoCo"
620 const std::uint64_t magic = 0x6F436F4461503346;
621
622 XASSERT(this->_num_global_layers >= this->_ancestry.size());
623
624 if(this->_comm.rank() == 0)
625 {
626 // dump magic
627 bs.write((const char*)&magic, 8u);
628
629 // write finest level
630 const u64 lev = u64(this->_virt_levels.front()->get_level_index());
631 bs.write((const char*)&lev, 8u);
632
633 // serialize ancestry/layers
634 const u64 na = u64(this->_ancestry.size());
635 bs.write((const char*)&na, 8u);
636
637 // write number of ranks per layer in reverse order
638 for(auto it = this->_ancestry.rbegin(); it != this->_ancestry.rend(); ++it)
639 {
640 // write ranks
641 const u64 np = u64(it->num_procs);
642 const u64 pl = u64(it->parti_level);
643 bs.write((const char*)&np, 8u);
644 bs.write((const char*)&pl, 8u);
645 }
646 }
647
648 // loop over all ancestry graphs
649 for(auto it = this->_ancestry.rbegin(); it != this->_ancestry.rend(); ++it)
650 {
651 if(it == this->_ancestry.rbegin())
652 {
653 std::vector<char> buf = it->parti_graph.serialize();
654 bs.write(buf.data(), std::streamsize(buf.size()));
655 continue;
656 }
657
658 // get layer index
659 if(it->layer_p < 0)
660 continue;
661
662 // get parent layer communicator
663 const std::size_t ilp = std::size_t(it->layer_p);
664 XASSERT(ilp < this->_layers.size());
665 const Dist::Comm& comm_p = this->_layers.at(ilp)->comm();
666
667 // serialize graph
668 std::vector<char> buf = it->parti_graph.serialize();
669
670 // choose maximum size
671 u64 buf_size = buf.size();
672
673 // gather individual buffer sizes on rank 0
674 std::vector<u64> all_sizes;
675 if(comm_p.rank() == 0)
676 all_sizes.resize(std::size_t(comm_p.size()));
677 comm_p.gather(&buf_size, 1u, all_sizes.data(), 1u, 0);
678
679 // allreduce maximum size
680 comm_p.allreduce(&buf_size, &buf_size, std::size_t(1), Dist::op_max);
681
682 // adjust buffer size
683 if(buf_size > u64(buf.size()))
684 buf.resize(buf_size);
685
686 // on rank 0, allocate common buffer
687 std::vector<char> com_buf;
688 if(comm_p.rank() == 0)
689 com_buf.resize(std::size_t(buf_size * u64(comm_p.size())));
690
691 // gather all buffers on rank 0
692 comm_p.gather(buf.data(), buf_size, com_buf.data(), buf_size, 0);
693
694 // write each individual buffer
695 if(comm_p.rank() == 0)
696 {
697 char* x = com_buf.data();
698 for(u64 k(0); k < u64(comm_p.size()); ++k)
699 bs.write(&x[k*buf_size], std::streamsize(all_sizes[k]));
700 }
701 }
702
703 return bs.container();
704 }
705
706 protected:
713 virtual void _reset_parti_types()
714 {
717 }
718
727 virtual bool _parse_parti_type(const String& type)
728 {
729 if((type.compare_no_case("metis") == 0) || (type.compare_no_case("parmetis") == 0))
730 return _allow_parti_metis = true;
731 else if(type.compare_no_case("genetic") == 0)
732 return _allow_parti_genetic = true;
733 else if(type.compare_no_case("naive") == 0)
734 return _allow_parti_naive = true;
735 else if(type.compare_no_case("zoltan") == 0)
736 return _allow_parti_zoltan = true;
737 else
738 return false;
739 }
740
751 virtual std::vector<WeightType> _compute_weights(Ancestor& DOXY(ancestor), const MeshNodeType& DOXY(base_mesh_node))
752 {
753 return std::vector<WeightType>();
754 }
755
760 {
761 // for more than 2 desired levels, call _create_ancestry_scattered
762 XASSERTM(this->_desired_levels.size() <= std::size_t(2), "multi-layered control is desired here");
763
764 // allocate and create ancestry
765 this->_ancestry.resize(std::size_t(1));
766 Ancestor& ancestor = this->_ancestry.front();
767
768 // set the layer index to 0
769 ancestor.layer = 0;
770
771 // set the parent layer index to -1
772 ancestor.layer_p = -1;
773
774 // set the total number of processes for this layer
775 ancestor.num_procs = this->_comm.size();
776
777 // set the total number of partitions per progeny group
778 ancestor.num_parts = this->_comm.size();
779
780 // set desired maximum level
781 ancestor.desired_level_max = this->_desired_levels.front().first;
782
783 // set desired minimum level (may be = level_max)
784 ancestor.desired_level_min = this->_desired_levels.back().first;
785
786 // set the progeny group
787 ancestor.progeny_group = 0;
788 ancestor.progeny_child = this->_comm.rank();
789
790 // create the progeny communicator
791 ancestor.progeny_count = this->_comm.size();
792 ancestor.progeny_first = 0;
793 ancestor.progeny_comm = this->_comm.comm_dup();
794 }
795
796 // Note: all following member functions are only required for parallel builds,
797 // so we enclose them in the following #if-block to reduce compile times.
798
799#if defined(FEAT_HAVE_MPI) || defined(DOXYGEN)
800
805 {
806 // create and push global layer
807 this->push_layer(std::make_shared<LayerType>(this->_comm.comm_dup(), 0));
808
809 // loop over all desired layers
810 for(std::size_t ilay(1u); (ilay+1u) < this->_desired_levels.size(); ++ilay)
811 {
812 // get child layer
813 std::shared_ptr<LayerType> layer_c = this->_layers.back();
814
815 // get child layer communicator
816 const Dist::Comm& comm_c = layer_c->comm();
817
818 // get number of processes in child comm
819 const int nprocs_c = comm_c.size();
820
821 // get desired number of processes for parent comm
822 const int nprocs_p = this->_desired_levels.at(ilay).second;
823 XASSERT(nprocs_p < nprocs_c);
824 XASSERT(nprocs_p > 0);
825 XASSERT(nprocs_c % nprocs_p == 0); // same number of children for each parent
826
827 // compute number of siblings = children per parent
828 const int num_sibs = nprocs_c / nprocs_p;
829
830 // get my rank in child comm
831 const int my_rank_c = comm_c.rank();
832
833 // compute my parent's rank in child comm
834 const int parent_rank_c = (my_rank_c / num_sibs) * num_sibs;
835
836 // create our sibling communicator and set the parent rank
837 layer_c->set_parent(comm_c.comm_create_range_incl(num_sibs, parent_rank_c), 0);
838
839 // next, create the actual parent communicator
840 Dist::Comm comm_p = comm_c.comm_create_range_incl(nprocs_p, 0, num_sibs);
841
842 // Are we the parent?
843 if(parent_rank_c == my_rank_c)
844 {
845 // make sure we have a valid communicator
846 XASSERT(!comm_p.is_null());
847
848 // push the parent layer
849 this->push_layer(std::make_shared<LayerType>(std::move(comm_p), int(ilay)));
850 }
851 else
852 {
853 // We are not a parent, so we must have received a null communicator
854 XASSERT(comm_p.is_null());
855
856 // Exit the loop, as we are not part of the party anymore...
857 break;
858 }
859 }
860 }
861
866 {
867 // the layers must have been created already
868 XASSERT(!this->_layers.empty());
869
870 // create a deque with the number of processes for each layer
871 std::deque<int> num_procs;
872 for(auto it = this->_desired_levels.begin(); it != this->_desired_levels.end(); ++it)
873 {
874 if(it->second > 1)
875 num_procs.push_back(it->second);
876 }
877 // manually add the base-layer
878 num_procs.push_back(1);
879
880 // allocate and create ancestry
881 this->_ancestry.resize(num_procs.size()-std::size_t(1));
882
883 // set up the ancestry
884 const int main_rank = this->_comm.rank();
885 const int main_size = this->_comm.size();
886 for(std::size_t i(0); i < this->_ancestry.size(); ++i)
887 {
888 // get our ancestor info
889 Ancestor& ancestor = this->_ancestry.at(i);
890
891 // set the layer index (or -1, if this process is not part of that layer)
892 ancestor.layer = (i < this->_layers.size() ? int(i) : -1);
893
894 // set the parent layer index (or -1, if this process is not in the parent layer)
895 ancestor.layer_p = ((i+1u) < this->_layers.size() ? int(i)+1 : -1);
896
897 // set the total number of processes for this layer
898 ancestor.num_procs = num_procs.at(i);
899
900 // set the total number of partitions per progeny group
901 ancestor.num_parts = num_procs.at(i) / num_procs.at(i+1u);
902
903 // set desired maximum level
904 XASSERT(i < this->_desired_levels.size());
905 ancestor.desired_level_max = this->_desired_levels.at(i).first;
906
907 // set desired minimum level
908 if((i+1u) < this->_desired_levels.size())
909 ancestor.desired_level_min = this->_desired_levels.at(i+1).first;
910 else
911 ancestor.desired_level_min = ancestor.desired_level_max;
912
913 // set the progeny group
914 ancestor.progeny_group = ((main_rank * num_procs.at(i+1u)) / main_size) * ancestor.num_parts;
915 ancestor.progeny_child = ((main_rank * num_procs.at(i)) / main_size) % ancestor.num_parts;
916
917 // create the progeny communicator
918 ancestor.progeny_count = main_size / num_procs.at(i+1u);
919 ancestor.progeny_first = (main_rank / ancestor.progeny_count) * ancestor.progeny_count;
920 ancestor.progeny_comm = this->_comm.comm_create_range_incl(ancestor.progeny_count, ancestor.progeny_first);
921 }
922 }
923
942 const Ancestor& ancestor,
943 const MeshNodeType& base_mesh_node,
944 MeshNodeType& patch_mesh_node,
945 std::vector<int>& neighbor_ranks)
946 {
947 // get the map of the base-mesh halos
948 const std::map<int, std::unique_ptr<MeshPartType>>& base_halo_map = base_mesh_node.get_halo_map();
949
950 // if the base mesh has no halos, then we can jump out of here
951 if(base_halo_map.empty())
952 return;
953
954 // get number of halos
955 const std::size_t num_halos = base_halo_map.size();
956
957 // create a halo splitter
958 Geometry::PatchHaloSplitter<MeshType> halo_splitter(*base_mesh_node.get_mesh(),
959 *base_mesh_node.get_patch(ancestor.progeny_child));
960
961 // add each base-mesh halo to our halo splitter and store the resulting split data size
962 std::vector<int> halo_ranks;
963 std::vector<std::size_t> halo_sizes;
964 for(auto it = base_halo_map.begin(); it != base_halo_map.end(); ++it)
965 {
966 // store halo rank
967 halo_ranks.push_back(it->first);
968
969 // try to split the halo and store the resulting data size
970 halo_sizes.push_back(halo_splitter.add_halo(it->first, *it->second));
971 }
972
973 // This vector will receive the split halo data from all our potential neighbor processes
974 std::vector<Index> halo_recv_data;
975 std::vector<std::vector<Index>> halo_send_data;
976
977 /* ******************************************************************************************************* */
978 /* ******************************************************************************************************* */
979 // PHASE I: collect halo data from our siblings
980 /* ******************************************************************************************************* */
981 /* ******************************************************************************************************* */
982 if(ancestor.layer >= 0)
983 {
984 XASSERT(this->_layers.size() > std::size_t(ancestor.layer));
985 const LayerType& layer = *this->_layers.at(std::size_t(ancestor.layer));
986 const Dist::Comm& sibling_comm = *layer.sibling_comm_ptr();
987 XASSERT(!sibling_comm.is_null());
988 const std::size_t num_sibls = std::size_t(sibling_comm.size());
989
990 // get the rank of this process in various communicators
991 const int sibl_rank = sibling_comm.rank();
992 const int layer_rank = layer.comm().rank();
993
994 // gather halo infos of all siblings at sibling rank 0
995 std::vector<std::size_t> sibl_halo_sizes(sibl_rank == 0 ? num_halos*num_sibls : 0u);
996 sibling_comm.gather(halo_sizes.data(), halo_sizes.size(), sibl_halo_sizes.data(), halo_sizes.size(), 0);
997
998 if(sibl_rank > 0)
999 {
1000 std::vector<std::vector<Index>> halo_split_data(num_halos);
1001 Dist::RequestVector send_reqs(num_halos);
1002
1003 // serialize all split halos
1004 for(std::size_t i(0); i < num_halos; ++i)
1005 {
1006 // skip empty halos
1007 if(halo_sizes.at(i) == Index(0))
1008 continue;
1009
1010 // serialize split halo data
1011 halo_split_data.at(i) = halo_splitter.serialize_split_halo(halo_ranks[i], layer_rank);
1012 XASSERT(halo_split_data.at(i).size() == halo_sizes.at(i));
1013
1014 // send split data over to our parent process
1015 send_reqs[i] = sibling_comm.isend(halo_split_data.at(i).data(), halo_sizes.at(i), 0);
1016 }
1017
1018 // wait for all pending sends to finish
1019 send_reqs.wait_all();
1020 }
1021 else // if(sibl_rank == 0)
1022 {
1023 halo_send_data.resize(num_halos);
1024
1025 // compute halo send data sizes
1026 for(std::size_t i(0); i < num_halos; ++i)
1027 {
1028 // determine the number of child processes for this halo
1029 // as well as the required send buffer size
1030 Index num_halo_childs = 0u;
1031 std::size_t buffer_size = 0u, offset = 0u;
1032 for(std::size_t j(0); j < num_sibls; ++j)
1033 {
1034 // update required buffer size
1035 buffer_size += sibl_halo_sizes.at(j*num_halos + i);
1036
1037 // this sibling is a child if its split halo is not empty
1038 if(sibl_halo_sizes.at(j*num_halos + i) > std::size_t(0))
1039 ++num_halo_childs;
1040 }
1041
1042 // increase buffer size to store child data offsets
1043 buffer_size += num_halo_childs + Index(1);
1044
1045 // allocate send buffer and get a pointer to its data array
1046 halo_send_data.at(i).resize(buffer_size);
1047 Index* halo_buffer = halo_send_data.at(i).data();
1048
1049 // store child count as first entry of buffer
1050 halo_buffer[0u] = num_halo_childs;
1051
1052 // initialize offset for first child
1053 offset = num_halo_childs + Index(1);
1054 Index coi = 0u; // child offset index
1055
1056 // collect my own split halo
1057 if(sibl_halo_sizes.at(i) > Index(0))
1058 {
1059 std::vector<Index> my_data(halo_splitter.serialize_split_halo(halo_ranks[i], layer_rank));
1060 std::size_t data_size = sibl_halo_sizes.at(i);
1061 halo_buffer[++coi] = Index(offset);
1062 for(std::size_t k(0); k < data_size; ++k)
1063 halo_buffer[offset+k] = my_data[k];
1064 offset += my_data.size();
1065 }
1066
1067 Dist::RequestVector sibl_recv_reqs(num_sibls);
1068
1069 // collect the other siblings
1070 for(std::size_t j(1); j < num_sibls; ++j)
1071 {
1072 // skips all sibling with empty halos
1073 std::size_t data_size = sibl_halo_sizes.at(j*num_halos + i);
1074 if(data_size == std::size_t(0))
1075 continue;
1076 XASSERT(offset+data_size <= buffer_size);
1077
1078 // store offset for this child
1079 halo_buffer[++coi] = Index(offset);
1080
1081 // receive serialized data from this sibling
1082 sibl_recv_reqs[j] = sibling_comm.irecv(&halo_buffer[offset], data_size, int(j));
1083 offset += data_size;
1084 }
1085 XASSERT(offset == buffer_size);
1086
1087 // wait for all pending receives to finish
1088 sibl_recv_reqs.wait_all();
1089 }
1090 }
1091 } // if(ancestor.layer >= 0)
1092
1093 /* ******************************************************************************************************* */
1094 /* ******************************************************************************************************* */
1095 // PHASE II: exchange split halos over parent layer communicator
1096 /* ******************************************************************************************************* */
1097 /* ******************************************************************************************************* */
1098
1099 if(ancestor.layer_p >= 0)
1100 {
1101 XASSERT(!halo_send_data.empty());
1102
1103 // get the parent layer communicator
1104 const Dist::Comm& parent_comm = this->_layers.at(std::size_t(ancestor.layer_p))->comm();
1105
1106 std::vector<std::size_t> halo_recv_sizes(num_halos), halo_send_sizes(num_halos);
1107 Dist::RequestVector halo_recv_reqs(num_halos), halo_send_reqs(num_halos);
1108
1109 // exchange halo send data sizes
1110 for(std::size_t i(0); i < num_halos; ++i)
1111 {
1112 // get halo send size
1113 halo_send_sizes.at(i) = halo_send_data.at(i).size();
1114
1115 // post receive requests
1116 halo_recv_reqs[i] = parent_comm.irecv(&halo_recv_sizes[i], std::size_t(1), halo_ranks[i]);
1117
1118 // post send requests
1119 halo_send_reqs[i] = parent_comm.isend(&halo_send_sizes[i], std::size_t(1), halo_ranks[i]);
1120 }
1121
1122 // wait for sends and receives to finish
1123 halo_recv_reqs.wait_all();
1124 halo_send_reqs.wait_all();
1125
1126 // compute total receive buffer size
1127 std::size_t recv_buf_size = num_halos + std::size_t(1);
1128 for(std::size_t i(0); i < num_halos; ++i)
1129 recv_buf_size += halo_recv_sizes[i];
1130
1131 // allocate receive buffer
1132 halo_recv_data.resize(recv_buf_size);
1133
1134 // set up receive data pointers
1135 halo_recv_data[0] = Index(num_halos) + Index(1);
1136 for(std::size_t i(0); i < num_halos; ++i)
1137 halo_recv_data[i+1u] = Index(halo_recv_data[i] + halo_recv_sizes[i]);
1138
1139 // allocate receive buffers and post receives
1140 for(std::size_t i(0); i < num_halos; ++i)
1141 {
1142 // resize buffer and post receive
1143 halo_recv_reqs[i] = parent_comm.irecv(&halo_recv_data[halo_recv_data[i]], halo_recv_sizes.at(i), halo_ranks[i]);
1144
1145 // post send of actual halo buffer
1146 halo_send_reqs[i] = parent_comm.isend(halo_send_data.at(i).data(), halo_send_sizes.at(i), halo_ranks[i]);
1147 }
1148
1149 // wait for sends and receives to finish
1150 halo_recv_reqs.wait_all();
1151 halo_send_reqs.wait_all();
1152 } // if(ancestor.layer_p >= 0)
1153
1154 /* ******************************************************************************************************* */
1155 /* ******************************************************************************************************* */
1156 // PHASE III: broadcast halo receive data over progeny comm
1157 /* ******************************************************************************************************* */
1158 /* ******************************************************************************************************* */
1159
1160 {
1161 // broadcast receive buffer size
1162 std::size_t recv_data_size = halo_recv_data.size();
1163 ancestor.progeny_comm.bcast(&recv_data_size, std::size_t(1), 0);
1164
1165 // allocate buffer
1166 if(ancestor.progeny_comm.rank() != 0)
1167 {
1168 XASSERT(halo_recv_data.empty()); // should be empty until now
1169 halo_recv_data.resize(recv_data_size);
1170 }
1171 else
1172 {
1173 XASSERT(!halo_recv_data.empty()); // must not be empty
1174 }
1175
1176 // broadcast buffer
1177 ancestor.progeny_comm.bcast(halo_recv_data.data(), recv_data_size, 0);
1178 }
1179
1180 /*{
1181 String s;
1182 for(std::size_t i(0); i < num_halos; ++i)
1183 {
1184 s += stringify(halo_ranks[i]);
1185 s += " >";
1186 for(Index j(halo_recv_data[i]); j < halo_recv_data[i+1]; ++j)
1187 (s += " ") += stringify(halo_recv_data[j]);
1188 s += "\n";
1189 }
1190 this->_comm.allprint(s);
1191 }*/
1192
1193 /* ******************************************************************************************************* */
1194 /* ******************************************************************************************************* */
1195 // PHASE IV: intersect received halo splits
1196 /* ******************************************************************************************************* */
1197 /* ******************************************************************************************************* */
1198
1199 for(std::size_t i(0); i < num_halos; ++i)
1200 {
1201 // get the offset of the first date for this halo
1202 const Index offset = halo_recv_data.at(i);
1203
1204 // get the number of child processes for this halo
1205 const Index num_childs = halo_recv_data.at(offset);
1206
1207 // loop over all child processes
1208 for(Index j(0); j < num_childs; ++j)
1209 {
1210 // compute the offset of this child's data within the large buffer
1211 const Index buffer_offset = offset + halo_recv_data.at(offset+j+1u);
1212
1213 // intersect with our other halo
1214 if(!halo_splitter.intersect_split_halo(halo_ranks[i], halo_recv_data, buffer_offset))
1215 continue; // no intersection
1216
1217 // get the new neighbor's rank
1218 const int neighbor_rank = int(halo_recv_data.at(buffer_offset));
1219
1220 // add the new neighbor to our list
1221 neighbor_ranks.push_back(neighbor_rank);
1222
1223 // create new mesh-part
1224 patch_mesh_node.add_halo(neighbor_rank, halo_splitter.make_unique());
1225 }
1226 }
1227 }
1228
1254 virtual bool _check_parti(Ancestor& ancestor, const MeshNodeType& mesh_node, bool DOXY(is_base_layer))
1255 {
1256 // We need to determine the required partitioning level for a-posteriori partitioning.
1257 // For this, first determine the factor by which the number of elements increases for each refinement:
1258 const Index factor = Index(Geometry::Intern::StandardRefinementTraits<typename BaseClass::ShapeType, BaseClass::ShapeType::dimension>::count);
1259
1260 // compute minimum number of elements for a-posteriori partitioning
1261 const Index min_elems = Index(ancestor.num_parts) * Index(_required_elems_per_rank);
1262
1263 // Okay, get the number of elements on base-mesh level 0
1264 Index num_elems = mesh_node.get_mesh()->get_num_elements();
1265
1266 // compute refinement level on which we have enough elements
1267 int level = 0;
1268 for(; num_elems < min_elems; ++level)
1269 {
1270 num_elems *= factor;
1271 }
1272
1273 // assume that we did not find a partitioning yet
1274 ancestor.parti_found = false;
1275
1276 // finally, the user may have ancestor a greater level for partitioning:
1277 ancestor.parti_level = Math::max(ancestor.parti_level, level);
1278
1279 // No a-priori partitioning; try a-posteriori partitioner later
1280 return false;
1281 }
1282
1295 virtual bool _apply_parti(Ancestor& ancestor, /*const*/ MeshNodeType& base_mesh_node)
1296 {
1297 // compute weights
1298 std::vector<WeightType> weights = this->_compute_weights(ancestor, base_mesh_node);
1299
1300 // First of all, check whether an a-priori partitioning was selected;
1301 // if so, then we do not have to apply any a-posteriori partitioner
1302 if(weights.empty() && ancestor.parti_apriori)
1303 return true;
1304
1305 // ensure that the mesh has enough elements
1306 //XASSERT(ancestor.num_parts <= int(base_mesh_node.get_mesh()->get_num_elements()));
1307 if(base_mesh_node.get_mesh()->get_num_elements() < Index(ancestor.num_parts))
1308 return false;
1309
1310 // try the various a-posteriori partitioners
1311 if(this->_apply_parti_metis(ancestor, base_mesh_node, weights))
1312 return true;
1313 if(this->_apply_parti_zoltan(ancestor, base_mesh_node, weights))
1314 return true;
1315 if(this->_apply_parti_genetic(ancestor, base_mesh_node, weights))
1316 return true;
1317 if(this->_apply_parti_naive(ancestor, base_mesh_node, weights))
1318 return true;
1319
1320 // we should not arrive here...
1321 return false;
1322 }
1323
1339 bool _apply_parti_zoltan(Ancestor& ancestor, const MeshNodeType& base_mesh_node, const std::vector<WeightType>& weights)
1340 {
1341#ifdef FEAT_HAVE_ZOLTAN
1342 // is this even allowed?
1343 if(!this->_allow_parti_zoltan)
1344 return false;
1345
1346 // create partitioner on the corresponding progeny communicator
1347 Geometry::PartiZoltan partitioner(ancestor.progeny_comm);
1348
1349 // get the corresponding adjacency graph for partitioning
1350 constexpr int shape_dim = MeshNodeType::MeshType::ShapeType::dimension;
1351 const auto& faces_at_elem = base_mesh_node.get_mesh()->template get_index_set<shape_dim, 0>();
1352
1353 // call the partitioner
1354 if(!partitioner.execute(faces_at_elem, Index(ancestor.num_parts), weights))
1355 return false;
1356
1357 // create partition graph
1358 ancestor.parti_graph = partitioner.build_elems_at_rank();
1359
1360 // set info string
1361 ancestor.parti_info = String("Applied Zoltan partitioner");
1362
1363 // partitioning found
1364 ancestor.parti_found = true;
1365
1366 // okay
1367 return true;
1368#else // FEAT_HAVE_ZOLTAN
1369 (void)ancestor;
1370 (void)base_mesh_node;
1371 (void)weights;
1372 return false;
1373#endif // FEAT_HAVE_ZOLTAN
1374 }
1375
1391 bool _apply_parti_metis(Ancestor& ancestor, const MeshNodeType& base_mesh_node, const std::vector<WeightType>& weights)
1392 {
1393#ifdef FEAT_HAVE_PARMETIS
1394 // is this even allowed?
1395 if(!this->_allow_parti_metis)
1396 return false;
1397
1398 // create partitioner on the corresponding progeny communicator
1399 Geometry::PartiParMETIS partitioner(ancestor.progeny_comm);
1400
1401 // get the corresponding adjacency graph for partitioning
1402 constexpr int shape_dim = MeshNodeType::MeshType::ShapeType::dimension;
1403 const auto& faces_at_elem = base_mesh_node.get_mesh()->template get_index_set<shape_dim, 0>();
1404
1405 // get vertices-at-element and vertex set
1406 const auto& verts_at_elem = base_mesh_node.get_mesh()->template get_index_set<shape_dim, 0>();
1407 const auto& vertex_set = base_mesh_node.get_mesh()->get_vertex_set();
1408
1409 // call the partitioner
1410 if(!partitioner.execute(faces_at_elem, verts_at_elem, vertex_set, ancestor.num_parts, weights))
1411 return false;
1412
1413 // create partition graph
1414 ancestor.parti_graph = partitioner.build_elems_at_rank();
1415
1416 // set info string
1417 ancestor.parti_info = String("Applied ParMETIS partitioner");
1418
1419 // partitioning found
1420 ancestor.parti_found = true;
1421
1422 // okay
1423 return true;
1424#else
1425 (void)ancestor;
1426 (void)base_mesh_node;
1427 (void)weights;
1428 return false;
1429#endif // FEAT_HAVE_PARMETIS
1430 }
1431
1444 bool _apply_parti_genetic(Ancestor& ancestor, /*const*/ MeshNodeType& base_mesh_node, const std::vector<WeightType>& weights)
1445 {
1446 // is this even allowed?
1447 if(!this->_allow_parti_genetic || !weights.empty())
1448 return false;
1449
1450 // create a genetic partitioner
1452 *base_mesh_node.get_mesh(),
1453 ancestor.progeny_comm,
1454 (Index)ancestor.num_parts,
1455 this->_genetic_time_init,
1456 this->_genetic_time_mutate);
1457
1458 // create elems-at-rank graph
1459 ancestor.parti_graph = partitioner.build_elems_at_rank();
1460
1461 // set info string
1462 ancestor.parti_info = String("Applied genetic partitioner");
1463
1464 // partitioning found
1465 ancestor.parti_found = true;
1466
1467 // okay
1468 return true;
1469 }
1470
1486 bool _apply_parti_naive(Ancestor& ancestor, const MeshNodeType& mesh_node, const std::vector<WeightType>& weights)
1487 {
1488 // is this even allowed?
1489 if(!this->_allow_parti_naive)
1490 return false;
1491
1492 const Index num_parts = Index(ancestor.num_parts);
1493 const Index num_elems = mesh_node.get_mesh()->get_num_elements();
1494 XASSERTM(num_parts <= num_elems, "Base-Mesh does not have enough elements");
1495
1496 // create elems-at-rank graph
1497 ancestor.parti_graph = Adjacency::Graph(num_parts, num_elems, num_elems);
1498 Index* ptr = ancestor.parti_graph.get_domain_ptr();
1499 Index* idx = ancestor.parti_graph.get_image_idx();
1500
1501 // weighted partitioning?
1502 if(!weights.empty())
1503 {
1504 // compute total number of non-zero weight elements
1505 WeightType weight_sum(0.0);
1506 for(auto w : weights)
1507 {
1508 weight_sum += w;
1509 }
1510
1511 ptr[0] = Index(0);
1512 WeightType weight_count(0.);
1513 for(Index i(1), last(0); i < num_parts; ++i)
1514 {
1515 // add as many elements as we need to achieve the target weight for this partition
1516 WeightType target_weight = weight_sum * WeightType(i) / WeightType(num_parts);
1517 for(; (weight_count < target_weight) && (last < num_elems); ++last)
1518 weight_count += weights[last];
1519
1520 ptr[i] = last;
1521 }
1522 ptr[num_parts] = num_elems;
1523 }
1524 else
1525 {
1526 // unweighted partitioning
1527 ptr[0] = Index(0);
1528 for(Index i(1); i < num_parts; ++i)
1529 ptr[i] = (i*num_elems) / num_parts;
1530 ptr[num_parts] = num_elems;
1531 }
1532
1533 for(Index j(0); j < num_elems; ++j)
1534 idx[j] = j;
1535
1536 // set info string
1537 ancestor.parti_info = String("Applied naive partitioner");
1538
1539 // partitioning found
1540 ancestor.parti_found = true;
1541
1542 // okay
1543 return true;
1544 }
1545#endif // defined(FEAT_HAVE_MPI) || defined(DOXYGEN)
1546 }; // class HierarchDomainControl<...>
1547 } // namespace Domain
1548 } // namespace Control
1549} // namespace FEAT
#define XASSERT(expr)
Assertion macro definition.
Definition: assertion.hpp:262
#define XASSERTM(expr, msg)
Assertion macro definition with custom message.
Definition: assertion.hpp:263
FEAT Kernel base header.
Adjacency Graph implementation.
Definition: graph.hpp:34
Index * get_domain_ptr()
Returns the domain pointer array.
Definition: graph.hpp:359
Index * get_image_idx()
Returns the image node index array.
Definition: graph.hpp:374
Binary Stream class.
std::vector< char > & container()
Returns a reference to the internal vector container.
Domain control base-class template.
bool parti_apriori
specifies whether the chosen partitioning is an a-priori partitioning strategy
int desired_level_max
the desired minimum and maximum refinement levels for this layer
int parti_level
the refinement level on which the patch is to be partitioned
String parti_info
a string containing some information about the chosen partitioning
bool parti_found
specifies whether a partitioning was found
int num_procs
the number of processes that participate in this layer
Adjacency::Graph parti_graph
this is the actual elements-at-rank partitioning graph
int layer_p
the index of the parent layer or -1, if this process is not part of the parent layer
int layer
the index of the layer that this ancestor object belongs to
int num_parts
the number of partitions for each patch of the parent layer
Base-Class for Hierarchical Domain Control implementations.
const std::deque< std::pair< int, int > > get_chosen_levels() const
Returns the deque of chosen refinement levels.
virtual String format_desired_levels() const
Returns the desired levels formatted as a parsable string.
int get_desired_level_min() const
Returns the desired minimum refinement level.
PartiDomainControlBase(const Dist::Comm &comm_, bool support_multi_layered)
Constructor.
virtual std::vector< WeightType > _compute_weights(Ancestor &ancestor, const MeshNodeType &base_mesh_node)
Computes the element partitioning weights.
virtual void _reset_parti_types()
Resets/disables all partitioner types.
void set_adapt_mode(Geometry::AdaptMode adapt_mode)
Sets the adapt-mode for refinement.
String dump_ancestry() const
Debugging function: Returns a string containing encoded ancestry information.
void set_desired_levels(const String &slvls)
Sets the desired levels for the partitioned hierarchy.
const std::deque< std::pair< int, int > > get_desired_levels() const
Returns the deque of desired refinement levels.
std::deque< std::pair< int, int > > _desired_levels
desired level deque
int _required_elems_per_rank
required number of elements per rank for a-posteriori partitioning
std::deque< std::pair< int, int > > _chosen_levels
chosen level deque
virtual void _create_ancestry_scattered()
Creates the layers for a multi-layered domain control in a scattered fashion.
Real WeightType
weight type for partitioner element weights; always Real
double _genetic_time_mutate
time for genetic partitioner mutation
std::deque< Ancestor > _ancestry
the partition ancestry deque
bool _apply_parti_zoltan(Ancestor &ancestor, const MeshNodeType &base_mesh_node, const std::vector< WeightType > &weights)
Applies the Zoltan partitioner onto the base-mesh.
bool _apply_parti_naive(Ancestor &ancestor, const MeshNodeType &mesh_node, const std::vector< WeightType > &weights)
Applies the naive partitioner onto the base-mesh.
bool _apply_parti_metis(Ancestor &ancestor, const MeshNodeType &base_mesh_node, const std::vector< WeightType > &weights)
Applies the METIS partitioner onto the base-mesh.
virtual bool parse_args(SimpleArgParser &args)
Parses the partitioner options from an argument parser.
bool _apply_parti_genetic(Ancestor &ancestor, MeshNodeType &base_mesh_node, const std::vector< WeightType > &weights)
Applies the genetic partitioner onto the base-mesh.
Geometry::RootMeshNode< MeshType > MeshNodeType
our root mesh node type
const std::deque< Ancestor > & get_ancestry() const
Returns a const reference to the ancestry.
std::vector< char > serialize_partitioning() const
Serializes the partitioning information into a binary buffer.
virtual bool _parse_parti_type(const String &type)
Parses a partitioner type.
virtual String format_chosen_levels() const
Returns the chosen levels formatted as a parsable string.
virtual void set_desired_levels(int lvl_max, int lvl_med, int lvl_min)
Sets the desired refinement levels for a double-layered hierarchy.
LevelType::PartType MeshPartType
our mesh-part type
Geometry::AdaptMode get_adapt_mode() const
Gets the adapt-mode for refinement.
virtual void _create_ancestry_single()
Creates the ancestry for a single layer (or a single process)
Geometry::AdaptMode _adapt_mode
the adapt mode for refinement
int get_desired_level_max() const
Returns the desired maximum refinement level.
virtual void set_desired_levels(const std::deque< String > &slvls)
Sets the desired levels for the partitioned hierarchy.
bool _was_created
specifies whether the domain control was already created
virtual String get_chosen_parti_info() const
Returns an informative string about the chosen partitioning.
virtual bool parse_property_map(const PropertyMap &pmap)
Parses the partitioner options from a PropertyMap.
double _genetic_time_init
time for genetic partitioner initialization
virtual bool _check_parti(Ancestor &ancestor, const MeshNodeType &mesh_node, bool is_base_layer)
Checks for an appropriate partitioning strategy.
bool _support_multi_layered
support multi-layered hierarchy?
virtual bool _apply_parti(Ancestor &ancestor, MeshNodeType &base_mesh_node)
Applies an a-posteriori partitioner.
Control::Domain::DomainControl< DomainLevel_ > BaseClass
Our base class.
virtual void set_desired_levels(int lvl_max, int lvl_min=-1)
Sets the desired refinement levels for a single-layered hierarchy.
virtual void _create_multi_layers_scattered()
Creates the layers for a multi-layered domain control in a scattered fashion.
virtual void _split_basemesh_halos(const Ancestor &ancestor, const MeshNodeType &base_mesh_node, MeshNodeType &patch_mesh_node, std::vector< int > &neighbor_ranks)
Splits the base-mesh halos and computes the inter-patch-mesh halos.
Communicator class.
Definition: dist.hpp:1349
void bcast(void *buffer, std::size_t count, const Datatype &datatype, int root) const
Blocking broadcast.
Definition: dist.cpp:541
Comm comm_create_range_incl(int count, int first=0, int stride=1) const
Creates a new sub-communicator from a strided range of ranks.
Definition: dist.cpp:472
void allreduce(const void *sendbuf, void *recvbuf, std::size_t count, const Datatype &datatype, const Operation &op) const
Blocking All-Reduce.
Definition: dist.cpp:655
int size() const
Returns the size of this communicator.
Definition: dist.hpp:1506
Request irecv(void *buffer, std::size_t count, const Datatype &datatype, int source, int tag=0) const
Nonblocking Receive.
Definition: dist.cpp:716
Request isend(const void *buffer, std::size_t count, const Datatype &datatype, int dest, int tag=0) const
Nonblocking Send.
Definition: dist.cpp:704
MPI_Comm comm
our MPI communicator handle
Definition: dist.hpp:1353
void gather(const void *sendbuf, std::size_t sendcount, const Datatype &sendtype, void *recvbuf, std::size_t recvcount, const Datatype &recvtype, int root) const
Blocking gather.
Definition: dist.cpp:553
Comm comm_dup() const
Creates a copy of this communicator.
Definition: dist.cpp:459
int rank() const
Returns the rank of this process in this communicator.
Definition: dist.hpp:1494
bool is_null() const
Checks whether this communicator is a null communicator.
Definition: dist.cpp:454
void print(std::ostream &os, const String &msg, int root=0) const
Prints a message line to an output stream.
Definition: dist.cpp:782
Communication Request vector class.
Definition: dist.hpp:640
void wait_all()
Blocks until all active requests are fulfilled.
Definition: dist.cpp:324
MeshType * get_mesh()
Returns the mesh of this node.
Definition: mesh_node.hpp:225
Iterative-Partitioner class template declaration.
ParMETIS mesh/graph partitioner backend.
bool execute(const IndexSet< nf_ > &faces_at_elem, const IndexSet< nv_ > &verts_at_elem, const VertexSet_ &vertices, const Index num_parts, const std::vector< Real > &weights)
Executes the ParMETIS graph partitioner.
Adjacency::Graph build_elems_at_rank() const
Builds and returns the elements-at-rank graph representing the partitioning.
Zoltan hypergraph partitioner backend.
bool execute(const IndexSet< n_ > &faces_at_elem, const Index num_parts, const std::vector< Real > &weights)
Executes the Zoltan partitioner.
Adjacency::Graph build_elems_at_rank() const
Builds and returns the elements-at-rank graph representing the partitioning.
Base-Mesh Patch Halo splitter.
Root mesh node class template.
Definition: mesh_node.hpp:748
const MeshPartType * get_patch(int rank) const
Returns a patch meshpart for a given child.
Definition: mesh_node.hpp:962
const std::map< int, std::unique_ptr< MeshPartType > > & get_halo_map() const
Definition: mesh_node.hpp:896
void add_halo(int rank, std::unique_ptr< MeshPartType > halo_part)
Adds a halo mesh part to this mesh node.
Definition: mesh_node.hpp:874
Class for parser related errors.
Definition: exception.hpp:132
A class organizing a tree of key-value pairs.
std::pair< String, bool > query(String key_path) const
Queries a value by its key path.
Simple argument parser implementation.
const std::pair< int, std::deque< String > > * query(const String &option) const
Query the parameters of an option.
int parse(const String &option, Prms_ &&... prms) const
Parses the parameters of an option.
String class implementation.
Definition: string.hpp:46
std::deque< String > split_by_whitespaces() const
Splits the string by white-spaces.
Definition: string.hpp:518
int compare_no_case(const String &other) const
Compares two strings without regard to case.
Definition: string.hpp:679
String pad_front(size_type len, char c=' ') const
Pads the front of the string up to a desired length.
Definition: string.hpp:392
const Operation op_max(MPI_MAX)
Operation wrapper for MPI_MAX.
Definition: dist.hpp:273
AdaptMode
Adapt mode enumeration.
Definition: mesh_node.hpp:53
T_ max(T_ a, T_ b)
Returns the maximum of two values.
Definition: math.hpp:137
FEAT namespace.
Definition: adjactor.hpp:12
double Real
Real data type.
String stringify(const T_ &item)
Converts an item into a String.
Definition: string.hpp:944
std::uint64_t Index
Index data type.