1 #ifndef LOAD_BALANCE_HPP
2 #define LOAD_BALANCE_HPP
31 MPI_Comm inter_period_comm;
42 void update_model_no_history(
const View<int*,CLayout,HostType>& current_partition,
const View<double*, HostType>& all_periods_timings){
44 for(
int i=0; i<all_periods_timings.size(); i++){
45 int node_offset = current_partition(i) - 1;
46 int next_node_offset = current_partition(i+1) - 1;
47 int nnodes = next_node_offset - node_offset;
49 double time_per_node = all_periods_timings(i)/nnodes;
50 for(
int j=0; j<nnodes; j++){
58 double new_accumulated_time = 0.0;
59 int ierr = GPTLget_wallclock(
region_name.c_str(), -1, &new_accumulated_time);
62 if(ierr != 0) new_accumulated_time = 0.0;
68 time_accumulated = new_accumulated_time;
80 inter_period_comm(inter_period_comm),
81 period_comm(period_comm),
84 update_method(update_method),
85 region_name(region_name)
93 MPI_Comm_size(inter_period_comm, &
n_periods);
100 period_comm(period_comm),
103 update_method(update_method),
104 region_name(region_name)
120 update_method(update_method),
121 region_name(region_name)
159 double largest_time = 0.0;
161 for(
int i=0; i<(proposed_partition.size()-1); i++){
162 double proc_time = 0.0;
163 for(
int i_node=proposed_partition(i)-1; i_node<proposed_partition(i+1)-1; i_node++){
166 if(proc_time>largest_time){
168 largest_time = proc_time;
181 void update_model(
const View<int*,CLayout,HostType>& current_partition){
191 MPI_Reduce(&time_spent_this_rank, &time_spent, 1, MPI_DOUBLE, MPI_MAX, 0, inter_period_comm);
194 time_spent = time_spent_this_rank;
201 MPI_Gather(&time_spent, 1, MPI_DOUBLE, all_periods_timings.data(), 1, MPI_DOUBLE, 0, period_comm);
203 all_periods_timings(0) = time_spent;
213 for(
int i=0; i<all_periods_timings.size(); i++){
227 exit_XGC(
"Error: Update method not available\n");
265 Kokkos::parallel_reduce(
"sum", Kokkos::RangePolicy<HostExSpace>(0,input.size()), [=](
const int i,
double& l_total){
273 bool greedily_fill_partition(
const View<double*,HostType>& weight,
const View<double*,HostType>& constraint1,
double target_weight_per_rank){
274 int nnode = weight.size();
278 double constraint1_in_this_rank = 0.0;
279 double weight_in_this_rank = 0.0;
282 bool assign_one_node_per_proc =
false;
283 for(
int i_node = 0; i_node<nnode; i_node++){
286 if(nnode-i_node==nproc-pid) assign_one_node_per_proc =
true;
289 bool rank_is_loaded =
false;
292 if(assign_one_node_per_proc) rank_is_loaded =
true;
298 if(weight_in_this_rank>=target_weight_per_rank) rank_is_loaded =
true;
304 constraint1_in_this_rank = 0.0;
305 weight_in_this_rank = 0.0;
307 if(pid==nproc-1)
break;
309 constraint1_in_this_rank += constraint1(i_node);
310 weight_in_this_rank += weight(i_node);
317 constraint1_in_this_rank = 0.0;
319 constraint1_in_this_rank += constraint1(i_node);
331 double largest_time = 0.0;
333 for(
int i=0; i<(partition.size()-1); i++){
334 double proc_time = 0.0;
335 for(
int i_node=partition(i)-1; i_node<partition(i+1)-1; i_node++){
336 proc_time += weight(i_node);
338 if(proc_time>largest_time){
340 largest_time = proc_time;
346 void one_weight_balance(
const View<double*,HostType>& weight,
const View<double*,CLayout,HostType> constraint1){
355 if(!meets_constraints){
359 if(!meets_constraints)
exit_XGC(
"\nUnexpected issue in load balance: constraint-based partition doesn't satisfy constraints\n");
363 if(
verbose) printf(
"The ideal distribution (%1.3e) does not satisfy constraints. This distribution (%1.3e) does.\n", ideal_weight_per_rank, upper_limit_weight_per_rank);
367 const double desired_precision = 0.01;
368 double desired_step_size = ideal_weight_per_rank*desired_precision;
369 double step_size = upper_limit_weight_per_rank - ideal_weight_per_rank;
370 double compromise_weight_per_rank = upper_limit_weight_per_rank;
371 if(
verbose) printf(
"\nEmploying a binary search to narrow down on the best option.");
372 while(step_size>desired_step_size){
376 if(meets_constraints){
378 compromise_weight_per_rank -= step_size;
381 compromise_weight_per_rank += step_size;
386 if(
verbose) printf(
"\n Stepped by %1.3e to %1.3e. The new partition does%s meet constraints.", step_size, compromise_weight_per_rank, meets_constraints?
"" :
"nt");
389 while(!meets_constraints){
390 compromise_weight_per_rank += step_size;
392 if(
verbose) printf(
"\n Stepped by %1.3e UP to %1.3e. The new partition does%s meet constraints.", step_size, compromise_weight_per_rank, meets_constraints?
"" :
"nt");
409 double observed_max_all_rgns_time = 0.0;
410 for(
int i=0; i<
regions.size(); i++){
411 observed_max_all_rgns_time +=
regions[i].get_observed_max_region_time();
415 double predicted_max_all_rgns_time = 0.0;
416 for(
int i=0; i<
regions.size(); i++){
420 double fractional_improvement = 1.0 - predicted_max_all_rgns_time/observed_max_all_rgns_time;
423 printf(
"\nDetermining whether to adopt the proposed partition:");
424 printf(
"\n Observed time with current partition was: %1.3e", observed_max_all_rgns_time);
425 printf(
"\n Predicted time with proposed partition, adjusting for historical undershoot (of Rgn 0: %1.3e): %1.3e",
regions[0].get_prediction_undershoot(), predicted_max_all_rgns_time);
426 printf(
"\n Adopted if Fractional improvement of adjusted prediction (%1.3e) exceeds specified threshold (%1.3e)\n", fractional_improvement,
threshold_to_rebalance);
445 void update_model(
const View<int*,CLayout,HostType>& current_partition){
446 for(
int i=0; i<
regions.size(); i++){
447 regions[i].update_model(current_partition);
452 for(
int i=0; i<
regions.size(); i++){
458 bool is_initialized=
true;
459 for(
int i=0; i<
regions.size(); i++){
460 is_initialized = (is_initialized &&
regions[i].get_model_is_initialized());
462 return is_initialized;
466 auto constraint1 = ptl_count;
468 if(
regions.size()!=1)
exit_XGC(
"Error: Load balancing is currently single-region only\n");
475 printf(
"\nNewly PROPOSED partition:\n");
499 printf(
"\nNEW PARTITION:\n");
510 const View<int*,CLayout,HostType>& old_partition){
513 TIMER(
"F0_REDISTRIBUTE",
514 f0_redistribute(plasma, pol_decomp, grid, magnetic_field, vgrid, old_partition) );
533 return (should_rebalance_int==1);
536 bool should_rebalance;
546 int should_rebalance_int = (should_rebalance ? 1 : 0);
550 return (should_rebalance_int==1);
561 double max_mem_redist_gb = 10.0;
562 std::string weighting_algorithm_str =
"Fortran";
565 max_mem_redist_gb = nlr.
get<
double>(
"max_mem_redist_gb", 10.0);
566 weighting_algorithm_str = nlr.
get<std::string>(
"weighting_algorithm",
"Fortran");
575 double max_n_ptl_on_rank = max_mem_redist_gb*1024*1024*1024/80.0;
577 if(weighting_algorithm_str==
"Fortran"){
579 }
else if(weighting_algorithm_str==
"ParticleBalance"){
581 }
else if(weighting_algorithm_str==
"SingleRegionBalance"){
584 exit_XGC(
"\nError: weighting_algorithm input not valid.");
591 std::string region_name =
"ipc2:PUSHE";
599 for(
int i=0; i<n_regions; i++){
601 regions.push_back(
LoadRegion(region_name, grid_nnode, pol_decomp.mpi.intpl_comm, pol_decomp.mpi.plane_comm, update_method,
verbose));
619 exit_XGC(
"\nError: ParticleConstraint is only available ConstraintOption.\n");
648 double avg_n_ptl = (double)(total_n_ptl)/pol_decomp.mpi.nranks;
650 double avg_n_ptl = total_n_ptl;
652 double max_ratio = max_n_ptl/avg_n_ptl;
653 if(
is_rank_zero()) printf(
" Species %d: (max/avg ptl per rank = %1.2f); total n_ptl = %lld\n", species.
idx, max_ratio, total_n_ptl);
664 TIMER(
"LOAD_BAL_UPDATE",
667 if(
is_rank_zero() &&
verbose) printf(
"Initializing timing model, no partition proposal yet.");
683 View<int*,CLayout,HostType> old_partition(
NoInit(
"old_partition"), pol_decomp.
gvid0_pid_h.layout());
684 Kokkos::deep_copy(old_partition, pol_decomp.
gvid0_pid_h);
687 TIMER(
"LOAD_BAL_SET_NEW",
688 set_new_partition(sml, grid, magnetic_field, vgrid, plasma, pol_decomp, weighting_algorithm) );
691 TIMER(
"LOAD_BAL_REDIST",
692 redistribute_load(sml, grid, magnetic_field, vgrid, plasma, pol_decomp, old_partition) );
View< double *, HostType > estimated_time_per_vertex
Definition: load_balance.hpp:22
static int GPTLstart(const char *name)
Definition: timer_macro.hpp:9
bool is_rank_zero()
Definition: globals.hpp:27
void one_weight_balance(const View< double *, HostType > &weight, const View< double *, CLayout, HostType > constraint1)
Definition: load_balance.hpp:346
MPI_Comm SML_COMM_WORLD
Definition: my_mpi.cpp:4
double threshold_to_rebalance
Definition: load_balance.hpp:258
T get(const string ¶m, const T default_val, int val_ind=0)
Definition: NamelistReader.hpp:373
double time_accumulated
Definition: load_balance.hpp:38
double predicted_max_region_time
Definition: load_balance.hpp:26
double observed_max_region_time
Definition: load_balance.hpp:27
double get_largest_predicted_time(const View< int *, CLayout, HostType > &partition, const View< double *, HostType > &weight) const
Definition: load_balance.hpp:330
Definition: velocity_grid.hpp:8
int n_unique_ranks
How many ranks are in a 'period' (in tokamaks, in a plane)
Definition: load_balance.hpp:35
subroutine plasma(grid, itr, p, dene_out, deni_out, Te_out, Ti_out, Vparai_out)
Calculate the plasma density, temperature, and parallel velocity for a point in triangle itr using pl...
Definition: neutral_totalf.F90:1235
Definition: plasma.hpp:94
Definition: NamelistReader.hpp:193
Definition: magnetic_field.hpp:12
void propose_new_partition(const Kokkos::View< double *, Kokkos::LayoutRight, HostType > &ptl_count, WeightingAlgorithm weighting_algorithm)
Definition: load_balance.hpp:465
int idx
Index in all_species.
Definition: species.hpp:77
bool verbose
Definition: load_balance.hpp:259
void for_all_nonadiabatic_species(F func, DevicePtlOpt device_ptl_opt=UseDevicePtl)
Definition: plasma.hpp:113
double prediction_undershoot
Definition: load_balance.hpp:28
double get_even_division(const View< double *, HostType > &input, int n) const
Definition: load_balance.hpp:263
void set_weights(double *ptl_count)
Kokkos::View< double **, Kokkos::LayoutRight, HostType > count_ptl_per_node_elec_main_ion(const Grid< DeviceType > &grid, const MagneticField< DeviceType > &magnetic_field, Plasma &plasma, bool kinetic_electrons)
Definition: count_ptl_per_node.cpp:29
void update_pol_decomp()
Definition: domain_decomposition.tpp:103
void update_model(const View< int *, CLayout, HostType > ¤t_partition)
Definition: load_balance.hpp:445
long long int get_total_n_ptl()
Definition: species.hpp:685
void f0_redistribute(Plasma &plasma, DomainDecomposition< DeviceType > &pol_decomp, const Grid< DeviceType > &grid, const MagneticField< DeviceType > &magnetic_field, const VelocityGrid &vgrid, const View< int *, CLayout, HostType > &old_partition)
Definition: f0_redistribute.cpp:8
void set_new_partition(const Simulation< DeviceType > &sml, const Grid< DeviceType > &grid, const MagneticField< DeviceType > &magnetic_field, const VelocityGrid &vgrid, Plasma &plasma, DomainDecomposition< DeviceType > &pol_decomp, WeightingAlgorithm weighting_algorithm)
Definition: load_balance.hpp:480
#define TIMER(N, F)
Definition: timer_macro.hpp:24
void calculate_load_imbalance()
bool will_rebalance(ReweightOption reweight_option, WeightingAlgorithm weighting_algorithm)
Definition: load_balance.hpp:522
void rebalance(const Simulation< DeviceType > &sml, const Grid< DeviceType > &grid, const MagneticField< DeviceType > &magnetic_field, const VelocityGrid &vgrid, Plasma &plasma, DomainDecomposition< DeviceType > &pol_decomp, ReweightOption reweight_option, WeightingAlgorithm weighting_algorithm=WeightingAlgorithm::Default)
Definition: load_balance.hpp:623
void use_namelist(const string &namelist)
Definition: NamelistReader.hpp:355
ConstraintOption
Definition: load_balance.hpp:241
double get_estimated_time_per_vertex(int i) const
Definition: load_balance.hpp:137
int my_period_rank
Definition: load_balance.hpp:36
void initialize_model()
Definition: load_balance.hpp:174
WeightingAlgorithm default_weighting_algorithm
Definition: load_balance.hpp:252
View< int *, CLayout, HostType > proposed_partition
Which processors get which vertices.
Definition: load_balance.hpp:256
double get_time_since_previous_call()
Definition: load_balance.hpp:56
bool get_model_is_initialized() const
Definition: load_balance.hpp:153
UpdateMethod
Definition: load_balance.hpp:16
std::string region_name
Definition: load_balance.hpp:23
bool verbose
Definition: load_balance.hpp:24
bool f0_grid
Whether to use f0 grid.
Definition: sml.hpp:60
std::vector< LoadRegion > regions
Definition: load_balance.hpp:254
bool namelist_present(const string &namelist)
Definition: NamelistReader.hpp:351
UpdateMethod update_method
Definition: load_balance.hpp:40
bool model_is_initialized
Definition: load_balance.hpp:25
bool electron_on
Use kinetic electrons.
Definition: sml.hpp:59
void redistribute_load(const Simulation< DeviceType > &sml, const Grid< DeviceType > &grid, const MagneticField< DeviceType > &magnetic_field, const VelocityGrid &vgrid, Plasma &plasma, DomainDecomposition< DeviceType > &pol_decomp, const View< int *, CLayout, HostType > &old_partition)
Definition: load_balance.hpp:508
double get_largest_predicted_time(const View< int *, CLayout, HostType > &proposed_partition) const
Definition: load_balance.hpp:158
void initialize_model()
Definition: load_balance.hpp:451
bool recommend_proposed_partition()
Definition: load_balance.hpp:407
void exit_XGC(std::string msg)
Definition: globals.hpp:37
bool greedily_fill_partition(const View< double *, HostType > &weight, const View< double *, HostType > &constraint1, double target_weight_per_rank)
Definition: load_balance.hpp:273
int n_periods
How many repeating periods there are; in tokamaks this is planes.
Definition: load_balance.hpp:34
Definition: magnetic_field.F90:1
Definition: load_balance.hpp:14
void print_new_partition()
Definition: load_balance.hpp:437
WeightingAlgorithm
Definition: load_balance.hpp:234
int assess_whether_to_rebalance_load()
LoadBalance(NLReader::NamelistReader &nlr, const DomainDecomposition< DeviceType > &pol_decomp)
Definition: load_balance.hpp:558
LoadRegion(std::string region_name, int n_vertices, UpdateMethod update_method, bool verbose)
Definition: load_balance.hpp:114
View< double *, HostType > get_estimated_time_per_vertex() const
Definition: load_balance.hpp:141
Definition: plasma.hpp:14
double constraint1_max
Definition: load_balance.hpp:261
View< double *, CLayout, HostType > count_all_ptl_per_node(const Grid< DeviceType > &grid, const MagneticField< DeviceType > &magnetic_field, Plasma &plasma)
Definition: count_ptl_per_node.cpp:71
ReweightOption
Definition: load_balance.hpp:245
void shift_all_species(Plasma &plasma, const Grid< DeviceType > &grid, const MagneticField< DeviceType > &magnetic_field, const DomainDecomposition< DeviceType > &pol_decomp, Shift::ShiftPh0 shift_ph0)
Definition: shift.cpp:301
void reset_timer()
Definition: load_balance.hpp:132
Definition: species.hpp:74
bool model_is_initialized()
Definition: load_balance.hpp:457
void update_model(const View< int *, CLayout, HostType > ¤t_partition)
Definition: load_balance.hpp:181
double get_observed_max_region_time() const
Definition: load_balance.hpp:149
Definition: load_balance.hpp:232
Kokkos::ViewAllocateWithoutInitializing NoInit
Definition: space_settings.hpp:68
void update_model_no_history(const View< int *, CLayout, HostType > ¤t_partition, const View< double *, HostType > &all_periods_timings)
Definition: load_balance.hpp:42
double get_prediction_undershoot() const
Definition: load_balance.hpp:145
int get_max_n_ptl()
Definition: species.hpp:696
static int GPTLstop(const char *name)
Definition: timer_macro.hpp:10
Kokkos::View< int *, Kokkos::LayoutRight, HostType > gvid0_pid_h
Which processors get which vertices (host)
Definition: domain_decomposition.hpp:45