1 #ifndef LOAD_BALANCE_HPP
2 #define LOAD_BALANCE_HPP
14 extern "C" void set_weights(
int* gvid0_pid,
double* ptl_count,
double* f0_node_cost);
34 MPI_Comm inter_period_comm;
45 void update_model_no_history(
const View<int*,CLayout,HostType>& current_partition,
const View<double*, HostType>& all_periods_timings){
47 for(
int i=0; i<all_periods_timings.size(); i++){
48 int node_offset = current_partition(i) - 1;
49 int next_node_offset = current_partition(i+1) - 1;
50 int nnodes = next_node_offset - node_offset;
52 double time_per_node = all_periods_timings(i)/nnodes;
53 for(
int j=0; j<nnodes; j++){
61 double new_accumulated_time = 0.0;
65 double accumulated_single_timer = 0;
66 int ierr = GPTLget_wallclock(
timer_names[i].c_str(), -1, &accumulated_single_timer);
69 if(ierr != 0) accumulated_single_timer = 0.0;
71 new_accumulated_time += accumulated_single_timer;
78 time_accumulated = new_accumulated_time;
97 inter_period_comm(inter_period_comm),
98 period_comm(period_comm),
101 update_method(update_method),
102 region_name(region_name),
103 timer_names(timer_names)
114 MPI_Comm_size(inter_period_comm, &
n_periods);
145 double largest_time = 0.0;
147 for(
int i=0; i<(proposed_partition.size()-1); i++){
148 double proc_time = 0.0;
149 for(
int i_node=proposed_partition(i)-1; i_node<proposed_partition(i+1)-1; i_node++){
152 if(proc_time>largest_time){
154 largest_time = proc_time;
167 void update_model(
const View<int*,CLayout,HostType>& current_partition){
177 MPI_Reduce(&time_spent_this_rank, &time_spent, 1, MPI_DOUBLE, MPI_MAX, 0, inter_period_comm);
180 time_spent = time_spent_this_rank;
187 MPI_Gather(&time_spent, 1, MPI_DOUBLE, all_periods_timings.data(), 1, MPI_DOUBLE, 0, period_comm);
197 for(
int i=0; i<all_periods_timings.size(); i++){
211 exit_XGC(
"Error: Update method not available\n");
249 Kokkos::parallel_reduce(
"sum", Kokkos::RangePolicy<HostExSpace>(0,input.size()), [=](
const int i,
double& l_total){
257 bool greedily_fill_partition(
const View<double*,HostType>& weight,
const View<double*,HostType>& constraint1,
double target_weight_per_rank){
258 int nnode = weight.size();
262 double constraint1_in_this_rank = 0.0;
263 double weight_in_this_rank = 0.0;
266 bool assign_one_node_per_proc =
false;
267 for(
int i_node = 0; i_node<nnode; i_node++){
270 if(nnode-i_node==nproc-pid) assign_one_node_per_proc =
true;
273 bool rank_is_loaded =
false;
276 if(assign_one_node_per_proc) rank_is_loaded =
true;
282 if(weight_in_this_rank>=target_weight_per_rank) rank_is_loaded =
true;
288 constraint1_in_this_rank = 0.0;
289 weight_in_this_rank = 0.0;
291 if(pid==nproc-1)
break;
293 constraint1_in_this_rank += constraint1(i_node);
294 weight_in_this_rank += weight(i_node);
301 constraint1_in_this_rank = 0.0;
303 constraint1_in_this_rank += constraint1(i_node);
315 double largest_time = 0.0;
317 for(
int i=0; i<(partition.size()-1); i++){
318 double proc_time = 0.0;
319 for(
int i_node=partition(i)-1; i_node<partition(i+1)-1; i_node++){
320 proc_time += weight(i_node);
322 if(proc_time>largest_time){
324 largest_time = proc_time;
330 void one_weight_balance(
const View<double*,HostType>& weight,
const View<double*,CLayout,HostType> constraint1){
339 if(!meets_constraints){
343 if(!meets_constraints)
exit_XGC(
"\nUnexpected issue in load balance: constraint-based partition doesn't satisfy constraints\n");
347 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);
351 const double desired_precision = 0.01;
352 double desired_step_size = ideal_weight_per_rank*desired_precision;
353 double step_size = upper_limit_weight_per_rank - ideal_weight_per_rank;
354 double compromise_weight_per_rank = upper_limit_weight_per_rank;
355 if(
verbose) printf(
"\nEmploying a binary search to narrow down on the best option.");
356 while(step_size>desired_step_size){
360 if(meets_constraints){
362 compromise_weight_per_rank -= step_size;
365 compromise_weight_per_rank += step_size;
370 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");
373 while(!meets_constraints){
374 compromise_weight_per_rank += step_size;
376 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");
393 double observed_max_all_rgns_time = 0.0;
394 for(
int i=0; i<
regions.size(); i++){
395 observed_max_all_rgns_time +=
regions[i].get_observed_max_region_time();
399 double predicted_max_all_rgns_time = 0.0;
400 for(
int i=0; i<
regions.size(); i++){
404 double fractional_improvement = 1.0 - predicted_max_all_rgns_time/observed_max_all_rgns_time;
407 printf(
"\nDetermining whether to adopt the proposed partition:");
408 printf(
"\n Observed time with current partition was: %1.3e", observed_max_all_rgns_time);
409 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);
410 printf(
"\n Adopted if Fractional improvement of adjusted prediction (%1.3e) exceeds specified threshold (%1.3e)\n", fractional_improvement,
threshold_to_rebalance);
429 void update_model(
const View<int*,CLayout,HostType>& current_partition){
430 for(
int i=0; i<
regions.size(); i++){
431 regions[i].update_model(current_partition);
436 for(
int i=0; i<
regions.size(); i++){
442 bool is_initialized=
true;
443 for(
int i=0; i<
regions.size(); i++){
444 is_initialized = (is_initialized &&
regions[i].get_model_is_initialized());
446 return is_initialized;
450 auto constraint1 = ptl_count;
452 if(
regions.size()!=1)
exit_XGC(
"Error: Load balancing is currently single-region only\n");
459 printf(
"\nNewly PROPOSED partition:\n");
489 printf(
"\nNEW PARTITION:\n");
500 const View<int*,CLayout,HostType>& old_partition){
503 TIMER(
"F0_REDISTRIBUTE",
504 f0_redistribute(plasma, pol_decomp, grid, magnetic_field, vgrid, old_partition) );
527 return (should_rebalance_int==1);
530 bool should_rebalance;
540 int should_rebalance_int = (should_rebalance ? 1 : 0);
544 return (should_rebalance_int==1);
555 double max_mem_redist_gb = 10.0;
556 std::string weighting_algorithm_str =
"Fortran";
559 max_mem_redist_gb = nlr.
get<
double>(
"max_mem_redist_gb", 10.0);
560 weighting_algorithm_str = nlr.
get<std::string>(
"weighting_algorithm",
"Fortran");
569 double max_n_ptl_on_rank = max_mem_redist_gb*1024*1024*1024/80.0;
571 if(weighting_algorithm_str==
"Fortran"){
573 }
else if(weighting_algorithm_str==
"ParticleBalance"){
575 }
else if(weighting_algorithm_str==
"SingleRegionBalance"){
578 exit_XGC(
"\nError: weighting_algorithm input not valid.");
601 exit_XGC(
"\nError: ParticleConstraint is only available ConstraintOption.\n");
618 if(weighting_algorithm == WeightingAlgorithm::Default){
619 weighting_algorithm = default_weighting_algorithm;
622 if(default_weighting_algorithm==WeightingAlgorithm::Fortran){
624 weighting_algorithm = default_weighting_algorithm;
628 if(weighting_algorithm==WeightingAlgorithm::Fortran) printf(
"\nLoad balance called with Fortran algorithm\n");
629 else if(weighting_algorithm==WeightingAlgorithm::ParticleBalance) printf(
"\nLoad balance called with Ptl algorithm\n");
630 else if(weighting_algorithm==WeightingAlgorithm::SingleRegionBalance) printf(
"\nLoad balance called with SingleRegion algorithm\n");
634 if(weighting_algorithm!=WeightingAlgorithm::Fortran){
639 double avg_n_ptl = (double)(total_n_ptl)/pol_decomp.mpi.nranks;
640 double max_ratio = max_n_ptl/avg_n_ptl;
641 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);
649 if(weighting_algorithm!=WeightingAlgorithm::ParticleBalance){
650 if(model_is_initialized()){
652 TIMER(
"LOAD_BAL_UPDATE",
655 if(
is_rank_zero() && verbose) printf(
"Initializing timing model, no partition proposal yet.");
664 if (
is_rank_zero()) propose_new_partition(ptl_count, weighting_algorithm);
670 if(will_rebalance(reweight_option, weighting_algorithm, f0_cost)){
672 View<int*,CLayout,HostType> old_partition(
NoInit(
"old_partition"), pol_decomp.
gvid0_pid_h.layout());
673 Kokkos::deep_copy(old_partition, pol_decomp.
gvid0_pid_h);
676 TIMER(
"LOAD_BAL_SET_NEW",
677 set_new_partition(sml, grid, magnetic_field, vgrid, plasma, pol_decomp, weighting_algorithm) );
680 TIMER(
"LOAD_BAL_REDIST",
681 redistribute_load(sml, grid, magnetic_field, vgrid, plasma, pol_decomp, old_partition) );
View< double *, HostType > estimated_time_per_vertex
Definition: load_balance.hpp:24
static int GPTLstart(const char *name)
Definition: timer_macro.hpp:9
bool will_rebalance(ReweightOption reweight_option, WeightingAlgorithm weighting_algorithm, double f0_cost)
Definition: load_balance.hpp:512
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:330
MPI_Comm SML_COMM_WORLD
Definition: my_mpi.cpp:4
double threshold_to_rebalance
Definition: load_balance.hpp:242
T get(const string ¶m, const T default_val, int val_ind=0)
Definition: NamelistReader.hpp:373
double time_accumulated
Definition: load_balance.hpp:41
double predicted_max_region_time
Definition: load_balance.hpp:29
double observed_max_region_time
Definition: load_balance.hpp:30
double get_largest_predicted_time(const View< int *, CLayout, HostType > &partition, const View< double *, HostType > &weight) const
Definition: load_balance.hpp:314
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:38
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:1238
Definition: plasma.hpp:109
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:449
int idx
Index in all_species.
Definition: species.hpp:78
bool verbose
Definition: load_balance.hpp:243
void for_all_nonadiabatic_species(F func, DevicePtlOpt device_ptl_opt=UseDevicePtl)
Definition: plasma.hpp:128
double prediction_undershoot
Definition: load_balance.hpp:31
double get_even_division(const View< double *, HostType > &input, int n) const
Definition: load_balance.hpp:247
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:109
void update_model(const View< int *, CLayout, HostType > ¤t_partition)
Definition: load_balance.hpp:429
long long int get_total_n_ptl()
Definition: species.hpp:671
void f0_redistribute(Plasma &plasma, const 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:184
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:464
#define TIMER(N, F)
Definition: timer_macro.hpp:24
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:606
void use_namelist(const string &namelist)
Definition: NamelistReader.hpp:355
ConstraintOption
Definition: load_balance.hpp:225
double get_estimated_time_per_vertex(int i) const
Definition: load_balance.hpp:123
int my_period_rank
Definition: load_balance.hpp:39
int nnodes_on_plane
Number of nodes on local plane.
Definition: domain_decomposition.hpp:43
void calculate_load_imbalance(double f0_cost)
void initialize_model()
Definition: load_balance.hpp:160
WeightingAlgorithm default_weighting_algorithm
Definition: load_balance.hpp:236
View< int *, CLayout, HostType > proposed_partition
Which processors get which vertices.
Definition: load_balance.hpp:240
double get_time_since_previous_call()
Definition: load_balance.hpp:59
bool get_model_is_initialized() const
Definition: load_balance.hpp:139
UpdateMethod
Definition: load_balance.hpp:18
std::string region_name
Definition: load_balance.hpp:25
bool verbose
Definition: load_balance.hpp:27
bool f0_grid
Whether to use f0 grid.
Definition: sml.hpp:65
std::vector< LoadRegion > regions
Definition: load_balance.hpp:238
bool namelist_present(const string &namelist)
Definition: NamelistReader.hpp:351
UpdateMethod update_method
Definition: load_balance.hpp:43
bool model_is_initialized
Definition: load_balance.hpp:28
bool electron_on
Use kinetic electrons.
Definition: sml.hpp:64
void reset_cost_trackers()
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:498
double get_largest_predicted_time(const View< int *, CLayout, HostType > &proposed_partition) const
Definition: load_balance.hpp:144
void initialize_model()
Definition: load_balance.hpp:435
T::value_type sum_view(const T &view)
Definition: view_arithmetic.hpp:56
bool recommend_proposed_partition()
Definition: load_balance.hpp:391
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:257
int n_periods
How many repeating periods there are; in tokamaks this is planes.
Definition: load_balance.hpp:37
Definition: magnetic_field.F90:1
Definition: load_balance.hpp:16
void print_new_partition()
Definition: load_balance.hpp:421
WeightingAlgorithm
Definition: load_balance.hpp:218
int assess_whether_to_rebalance_load()
LoadBalance(NLReader::NamelistReader &nlr, const DomainDecomposition< DeviceType > &pol_decomp)
Definition: load_balance.hpp:552
View< double *, HostType > get_estimated_time_per_vertex() const
Definition: load_balance.hpp:127
Definition: plasma.hpp:14
double constraint1_max
Definition: load_balance.hpp:245
void set_weights(int *gvid0_pid, double *ptl_count, double *f0_node_cost)
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:229
View< double *, CLayout, HostType > f0_node_cost
Definition: plasma.hpp:38
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:118
Definition: species.hpp:75
bool model_is_initialized()
Definition: load_balance.hpp:441
void update_model(const View< int *, CLayout, HostType > ¤t_partition)
Definition: load_balance.hpp:167
std::vector< std::string > timer_names
Definition: load_balance.hpp:26
bool pol_decomp
Use poloidal decomposition.
Definition: domain_decomposition.hpp:36
double get_observed_max_region_time() const
Definition: load_balance.hpp:135
Definition: load_balance.hpp:216
void touch_timers()
Definition: load_balance.hpp:83
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:45
double get_prediction_undershoot() const
Definition: load_balance.hpp:131
int get_max_n_ptl()
Definition: species.hpp:682
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:53