Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
[project-description] Fix extraction of the ns-3 version.
[simgrid.git] / src / kernel / lmm / bmf.cpp
1 /* Copyright (c) 2007-2022. The SimGrid Team. All rights reserved.          */
2
3 /* This program is free software; you can redistribute it and/or modify it
4  * under the terms of the license (GNU LGPL) which comes with this package. */
5
6 #include "src/kernel/lmm/bmf.hpp"
7 #include "xbt/config.hpp"
8
9 #include <Eigen/LU>
10 #include <iostream>
11 #include <numeric>
12 #include <sstream>
13
14 XBT_LOG_NEW_DEFAULT_SUBCATEGORY(ker_bmf, kernel, "Kernel BMF solver");
15
16 simgrid::config::Flag<int>
17     cfg_bmf_max_iteration("bmf/max-iterations",
18                           "Maximum number of steps to be performed while searching for a BMF allocation", 1000);
19
20 simgrid::config::Flag<double> cfg_bmf_precision{"bmf/precision",
21                                                 "Numerical precision used when computing resource sharing", 1E-12};
22
23 namespace simgrid::kernel::lmm {
24
25 AllocationGenerator::AllocationGenerator(Eigen::MatrixXd A) : A_(std::move(A)), alloc_(A_.cols(), 0)
26 {
27   // got a first valid allocation
28   for (size_t p = 0; p < alloc_.size(); p++) {
29     for (int r = 0; r < A_.rows(); r++) {
30       if (A_(r, p) > 0) {
31         alloc_[p] = r;
32         break;
33       }
34     }
35   }
36 }
37
38 bool AllocationGenerator::next(std::vector<int>& next_alloc)
39 {
40   if (first_) {
41     next_alloc = alloc_;
42     first_     = false;
43     return true;
44   }
45
46   auto n_resources = A_.rows();
47   size_t idx      = 0;
48   while (idx < alloc_.size()) {
49     alloc_[idx] = (alloc_[idx] + 1) % n_resources;
50     if (alloc_[idx] == 0) {
51       idx++;
52       continue;
53     } else {
54       idx = 0;
55     }
56     if (A_(alloc_[idx], idx) > 0) {
57       next_alloc = alloc_;
58       return true;
59     }
60   }
61   return false;
62 }
63
64 /*****************************************************************************/
65
66 BmfSolver::BmfSolver(Eigen::MatrixXd A, Eigen::MatrixXd maxA, Eigen::VectorXd C, std::vector<bool> shared,
67                      Eigen::VectorXd phi)
68     : A_(std::move(A))
69     , maxA_(std::move(maxA))
70     , C_(std::move(C))
71     , C_shared_(std::move(shared))
72     , phi_(std::move(phi))
73     , gen_(A_)
74     , max_iteration_(cfg_bmf_max_iteration)
75
76 {
77   xbt_assert(max_iteration_ > 0,
78              "Invalid number of iterations for BMF solver. Please check your \"bmf/max-iterations\" configuration.");
79   xbt_assert(A_.cols() == maxA_.cols(), "Invalid number of cols in matrix A (%td) or maxA (%td)", A_.cols(),
80              maxA_.cols());
81   xbt_assert(A_.cols() == phi_.size(), "Invalid size of phi vector (%td)", phi_.size());
82   xbt_assert(static_cast<long>(C_shared_.size()) == C_.size(), "Invalid size param shared (%zu)", C_shared_.size());
83 }
84
85 template <typename T> std::string BmfSolver::debug_eigen(const T& obj) const
86 {
87   std::stringstream debug;
88   debug << obj;
89   return debug.str();
90 }
91
92 template <typename C> std::string BmfSolver::debug_vector(const C& container) const
93 {
94   std::stringstream debug;
95   std::copy(container.begin(), container.end(),
96             std::ostream_iterator<typename std::remove_reference<decltype(container)>::type::value_type>(debug, " "));
97   return debug.str();
98 }
99
100 std::string BmfSolver::debug_alloc(const allocation_map_t& alloc) const
101 {
102   std::stringstream debug;
103   for (const auto& [resource, players] : alloc) {
104     debug << "{" + std::to_string(resource) + ": [" + debug_vector(players) + "]}, ";
105   }
106   return debug.str();
107 }
108
109 double BmfSolver::get_resource_capacity(int resource, const std::vector<int>& bounded_players) const
110 {
111   double capacity = C_[resource];
112   if (not C_shared_[resource])
113     return capacity;
114
115   for (int p : bounded_players) {
116     capacity -= A_(resource, p) * phi_[p];
117   }
118   return std::max(0.0, capacity);
119 }
120
121 double BmfSolver::get_maxmin_share(int resource, const std::vector<int>& bounded_players) const
122 {
123   auto n_players  = (A_.row(resource).array() > 0).count() - bounded_players.size();
124   double capacity = get_resource_capacity(resource, bounded_players);
125   if (n_players > 0)
126     capacity /= n_players;
127   return capacity;
128 }
129
130 std::vector<int> BmfSolver::alloc_map_to_vector(const allocation_map_t& alloc) const
131 {
132   std::vector<int> alloc_by_player(A_.cols(), -1);
133   for (const auto& [resource, players] : alloc) {
134     for (auto p : players) {
135       alloc_by_player[p] = resource;
136     }
137   }
138   return alloc_by_player;
139 }
140
141 std::vector<int> BmfSolver::get_bounded_players(const allocation_map_t& alloc) const
142 {
143   std::vector<int> bounded_players;
144   for (const auto& [resource, players] : alloc) {
145     if (resource == NO_RESOURCE) {
146       bounded_players.insert(bounded_players.end(), players.begin(), players.end());
147     }
148   }
149   return bounded_players;
150 }
151
152 Eigen::VectorXd BmfSolver::equilibrium(const allocation_map_t& alloc) const
153 {
154   auto n_players      = A_.cols();
155   Eigen::MatrixXd A_p = Eigen::MatrixXd::Zero(n_players, n_players); // square matrix with number of players
156   Eigen::VectorXd C_p = Eigen::VectorXd::Zero(n_players);
157
158   int row = 0;
159   auto bounded_players = get_bounded_players(alloc);
160   for (const auto& [resource, players] : alloc) {
161     // add one row for the resource with A[r,]
162     /* bounded players, nothing to do */
163     if (resource == NO_RESOURCE)
164       continue;
165     /* not shared resource, each player can receive the full capacity of the resource */
166     if (not C_shared_[resource]) {
167       for (int i : players) {
168         C_p[row]    = get_resource_capacity(resource, bounded_players);
169         A_p(row, i) = A_(resource, i);
170         row++;
171       }
172       continue;
173     }
174
175     /* shared resource: fairly share it between players */
176     A_p.row(row) = A_.row(resource);
177     C_p[row]     = get_resource_capacity(resource, bounded_players);
178     row++;
179     if (players.size() > 1) {
180       // if 2 players have chosen the same resource
181       // they must have a fair sharing of this resource, adjust A_p and C_p accordingly
182       auto it = players.begin();
183       int i   = *it; // first player
184       /* for each other player sharing this resource */
185       for (++it; it != players.end(); ++it) {
186         /* player i and k on this resource j: so maxA_ji*rho_i - maxA_jk*rho_k = 0 */
187         int k       = *it;
188         C_p[row]    = 0;
189         A_p(row, i) = maxA_(resource, i);
190         A_p(row, k) = -maxA_(resource, k);
191         row++;
192       }
193     }
194   }
195   /* clear players which are externally bounded */
196   for (int p : bounded_players) {
197     A_p.col(p).setZero();
198   }
199
200   XBT_DEBUG("A':\n%s", debug_eigen(A_p).c_str());
201
202   XBT_DEBUG("C':\n%s", debug_eigen(C_p).c_str());
203   /* PartialPivLU is much faster than FullPivLU but requires that the matrix is invertible
204    * FullPivLU however assures that it finds come solution even if the matrix is singular
205    * Ideally we would like to be optimist and try Partial and in case of error, go back
206    * to FullPivLU.
207    * However, this with isNaN doesn't work if compiler uses -Ofastmath. In our case,
208    * the icc compiler raises an error when compiling the code (comparison with NaN always evaluates to false in fast
209    * floating point modes).
210    * Eigen::VectorXd rho = Eigen::PartialPivLU<Eigen::MatrixXd>(A_p).solve(C_p);
211    * if (rho.array().isNaN().any()) {
212    *   XBT_DEBUG("rho with nan values, falling back to FullPivLU, rho:\n%s", debug_eigen(rho).c_str());
213    *   rho = Eigen::FullPivLU<Eigen::MatrixXd>(A_p).solve(C_p);
214    * }
215    */
216
217   Eigen::VectorXd rho = Eigen::FullPivLU<Eigen::MatrixXd>(A_p).solve(C_p);
218   for (int p : bounded_players) {
219     rho[p] = phi_[p];
220   }
221   return rho;
222 }
223
224 bool BmfSolver::disturb_allocation(allocation_map_t& alloc, std::vector<int>& alloc_by_player)
225 {
226   while (gen_.next(alloc_by_player)) {
227     if (allocations_.find(alloc_by_player) == allocations_.end()) {
228       allocations_.clear();
229       allocations_.insert(alloc_by_player);
230       alloc.clear();
231       for (size_t p = 0; p < alloc_by_player.size(); p++) {
232         alloc[alloc_by_player[p]].insert(p);
233       }
234       return false;
235     }
236   }
237   return true;
238 }
239
240 bool BmfSolver::get_alloc(const Eigen::VectorXd& fair_sharing, const allocation_map_t& last_alloc,
241                           allocation_map_t& alloc, bool initial)
242 {
243   alloc.clear();
244   for (int player_idx = 0; player_idx < A_.cols(); player_idx++) {
245     int selected_resource = NO_RESOURCE;
246
247     /* the player's maximal rate is the minimum among all resources */
248     double min_rate = -1;
249     for (int cnst_idx = 0; cnst_idx < A_.rows(); cnst_idx++) {
250       if (A_(cnst_idx, player_idx) <= 0.0)
251         continue;
252
253       /* Note: the max_ may artificially increase the rate if priority < 0
254        * The equilibrium sets a rho which respects the C_ though */
255       if (double rate = fair_sharing[cnst_idx] / maxA_(cnst_idx, player_idx);
256           min_rate == -1 || double_positive(min_rate - rate, cfg_bmf_precision)) {
257         selected_resource = cnst_idx;
258         min_rate          = rate;
259       }
260       /* Given that the priority may artificially increase the rate,
261        * we need to check that the bound given by user respects the resource capacity C_ */
262       if (double bound = initial ? -1 : phi_[player_idx]; bound > 0 &&
263                                                           bound * A_(cnst_idx, player_idx) < C_[cnst_idx] &&
264                                                           double_positive(min_rate - bound, cfg_bmf_precision)) {
265         selected_resource = NO_RESOURCE;
266         min_rate          = bound;
267       }
268     }
269     alloc[selected_resource].insert(player_idx);
270   }
271   if (alloc == last_alloc) // considered stable
272     return true;
273
274   if (auto alloc_by_player = alloc_map_to_vector(alloc); not allocations_.insert(alloc_by_player).second) {
275     /* oops, allocation already tried, let's pertube it a bit */
276     XBT_DEBUG("Allocation already tried: %s", debug_alloc(alloc).c_str());
277     return disturb_allocation(alloc, alloc_by_player);
278   }
279   return false;
280 }
281
282 void BmfSolver::set_fair_sharing(const allocation_map_t& alloc, const Eigen::VectorXd& rho,
283                                  Eigen::VectorXd& fair_sharing) const
284 {
285   std::vector<int> bounded_players = get_bounded_players(alloc);
286
287   for (int r = 0; r < fair_sharing.size(); r++) {
288     auto it = alloc.find(r);
289     if (it != alloc.end()) { // resource selected by some player, fair share depends on rho
290       double min_share = std::numeric_limits<double>::max();
291       for (int p : it->second) {
292         double share = A_(r, p) * rho[p];
293         min_share    = std::min(min_share, share);
294       }
295       fair_sharing[r] = min_share;
296     } else { // nobody selects this resource, fair_sharing depends on resource saturation
297       // resource r is saturated (A[r,*] * rho > C), divide it among players
298       double consumption_r = A_.row(r) * rho;
299       double_update(&consumption_r, C_[r], cfg_bmf_precision);
300       if (consumption_r > 0.0) {
301         fair_sharing[r] = get_maxmin_share(r, bounded_players);
302       } else {
303         fair_sharing[r] = C_[r];
304       }
305     }
306   }
307 }
308
309 bool BmfSolver::is_bmf(const Eigen::VectorXd& rho) const
310 {
311   bool bmf = true;
312
313   // 1) the capacity of all resources is respected
314   Eigen::VectorXd shared(C_shared_.size());
315   for (int j = 0; j < shared.size(); j++)
316     shared[j] = C_shared_[j] ? 1.0 : 0.0;
317
318   Eigen::VectorXd remaining = (A_ * rho) - C_;
319   remaining                 = remaining.array() * shared.array(); // ignore non shared resources
320   bmf                       = bmf && (not std::any_of(remaining.data(), remaining.data() + remaining.size(),
321                                 [](double v) { return double_positive(v, sg_maxmin_precision); }));
322
323   // 3) every player receives maximum share in at least 1 saturated resource
324   // due to subflows, compare with the maximum consumption and not the A matrix
325   Eigen::MatrixXd usage =
326       maxA_.array().rowwise() * rho.transpose().array(); // usage_ji: indicates the usage of player i on resource j
327
328   XBT_DEBUG("Usage_ji considering max consumption:\n%s", debug_eigen(usage).c_str());
329   auto max_share = usage.rowwise().maxCoeff(); // max share for each resource j
330
331   // matrix_ji: boolean indicating player p has the maximum share at resource j
332   Eigen::MatrixXi player_max_share =
333       ((usage.array().colwise() - max_share.array()).abs() <= sg_maxmin_precision).cast<int>();
334   // but only saturated resources must be considered
335   Eigen::VectorXi saturated = (remaining.array().abs() <= sg_maxmin_precision).cast<int>();
336   XBT_DEBUG("Saturated_j resources:\n%s", debug_eigen(saturated).c_str());
337   player_max_share.array().colwise() *= saturated.array();
338
339   // just check if it has received at least it's bound
340   for (int p = 0; p < rho.size(); p++) {
341     if (double_equals(rho[p], phi_[p], sg_maxmin_precision)) {
342       player_max_share(0, p) = 1; // it doesn't really matter, just to say that it's a bmf
343       saturated[0]           = 1;
344     }
345   }
346
347   // 2) at least 1 resource is saturated
348   bmf = bmf && (saturated.array() == 1).any();
349
350   XBT_DEBUG("Player_ji usage of saturated resources:\n%s", debug_eigen(player_max_share).c_str());
351   // for all columns(players) it has to be the max at least in 1
352   bmf = bmf && (player_max_share.colwise().sum().array() >= 1).all();
353   return bmf;
354 }
355
356 Eigen::VectorXd BmfSolver::solve()
357 {
358   XBT_DEBUG("Starting BMF solver");
359
360   XBT_DEBUG("A:\n%s", debug_eigen(A_).c_str());
361   XBT_DEBUG("maxA:\n%s", debug_eigen(maxA_).c_str());
362   XBT_DEBUG("C:\n%s", debug_eigen(C_).c_str());
363   XBT_DEBUG("phi:\n%s", debug_eigen(phi_).c_str());
364
365   /* no flows to share, just returns */
366   if (A_.cols() == 0)
367     return {};
368
369   int it            = 0;
370   auto fair_sharing = C_;
371
372   /* BMF allocation for each player (current and last one) stop when are equal */
373   allocation_map_t last_alloc;
374   allocation_map_t cur_alloc;
375   Eigen::VectorXd rho;
376
377   while (it < max_iteration_ && not get_alloc(fair_sharing, last_alloc, cur_alloc, it == 0)) {
378     last_alloc = cur_alloc;
379     XBT_DEBUG("BMF: iteration %d", it);
380     XBT_DEBUG("B (current allocation): %s", debug_alloc(cur_alloc).c_str());
381
382     // solve inv(A)*rho = C
383     rho = equilibrium(cur_alloc);
384     XBT_DEBUG("rho:\n%s", debug_eigen(rho).c_str());
385
386     // get fair sharing for each resource
387     set_fair_sharing(cur_alloc, rho, fair_sharing);
388     XBT_DEBUG("Fair sharing vector (per resource):\n%s", debug_eigen(fair_sharing).c_str());
389
390     // get new allocation for players
391     it++;
392   }
393
394   /* Not mandatory but a safe check to assure we have a proper solution */
395   if (not is_bmf(rho)) {
396     fprintf(stderr, "Unable to find a BMF allocation for your system.\n"
397                     "You may try to increase the maximum number of iterations performed by BMF solver "
398                     "(\"--cfg=bmf/max-iterations\").\n"
399                     "Additionally, you could adjust numerical precision (\"--cfg=bmf/precision\").\n");
400     fprintf(stderr, "Internal states (after %d iterations):\n", it);
401     fprintf(stderr, "A:\n%s\n", debug_eigen(A_).c_str());
402     fprintf(stderr, "maxA:\n%s\n", debug_eigen(maxA_).c_str());
403     fprintf(stderr, "C:\n%s\n", debug_eigen(C_).c_str());
404     fprintf(stderr, "C_shared:\n%s\n", debug_vector(C_shared_).c_str());
405     fprintf(stderr, "phi:\n%s\n", debug_eigen(phi_).c_str());
406     fprintf(stderr, "rho:\n%s\n", debug_eigen(rho).c_str());
407     xbt_abort();
408   }
409
410   XBT_DEBUG("BMF done after %d iterations", it);
411   return rho;
412 }
413
414 /*****************************************************************************/
415
416 void BmfSystem::get_flows_data(Eigen::Index number_cnsts, Eigen::MatrixXd& A, Eigen::MatrixXd& maxA,
417                                Eigen::VectorXd& phi)
418 {
419   A.resize(number_cnsts, variable_set.size());
420   A.setZero();
421   maxA.resize(number_cnsts, variable_set.size());
422   maxA.setZero();
423   phi.resize(variable_set.size());
424
425   int var_idx = 0;
426   for (Variable& var : variable_set) {
427     if (var.sharing_penalty_ <= 0)
428       continue;
429     bool active = false;
430     bool linked = false; // variable is linked to some constraint (specially for selective_update)
431     for (const Element& elem : var.cnsts_) {
432       if (const auto& cnst_hook = selective_update_active ? elem.constraint->modified_constraint_set_hook_
433                                                           : elem.constraint->active_constraint_set_hook_;
434           not cnst_hook.is_linked())
435         continue;
436       /* active and linked variable, lets check its consumption */
437       linked             = true;
438       double consumption = elem.consumption_weight;
439       if (consumption > 0) {
440         int cnst_idx = cnst2idx_[elem.constraint];
441         A(cnst_idx, var_idx) += consumption;
442         // a variable with double penalty must receive half share, so it max weight is greater
443         maxA(cnst_idx, var_idx) = std::max(maxA(cnst_idx, var_idx), elem.max_consumption_weight * var.sharing_penalty_);
444         active                  = true;
445       }
446     }
447     /* skip variables not linked to any modified or active constraint */
448     if (not linked)
449       continue;
450     if (active) {
451       phi[var_idx]      = var.get_bound();
452       idx2Var_[var_idx] = &var;
453       var_idx++;
454     } else {
455       var.value_ = 1; // assign something by default for tasks with 0 consumption
456     }
457   }
458   // resize matrix to active variables only
459   A.conservativeResize(Eigen::NoChange_t::NoChange, var_idx);
460   maxA.conservativeResize(Eigen::NoChange_t::NoChange, var_idx);
461   phi.conservativeResize(var_idx);
462 }
463
464 template <class CnstList>
465 void BmfSystem::get_constraint_data(const CnstList& cnst_list, Eigen::VectorXd& C, std::vector<bool>& shared)
466 {
467   C.resize(cnst_list.size());
468   shared.resize(cnst_list.size());
469   cnst2idx_.clear();
470   int cnst_idx = 0;
471   for (const Constraint& cnst : cnst_list) {
472     C(cnst_idx)      = cnst.bound_;
473     if (cnst.get_sharing_policy() == Constraint::SharingPolicy::NONLINEAR && cnst.dyn_constraint_cb_) {
474       C(cnst_idx) = cnst.dyn_constraint_cb_(cnst.bound_, cnst.concurrency_current_);
475     }
476     cnst2idx_[&cnst] = cnst_idx;
477     // FATPIPE links aren't really shared
478     shared[cnst_idx] = (cnst.sharing_policy_ != Constraint::SharingPolicy::FATPIPE);
479     cnst_idx++;
480   }
481 }
482
483 void BmfSystem::do_solve()
484 {
485   if (selective_update_active)
486     bmf_solve(modified_constraint_set);
487   else
488     bmf_solve(active_constraint_set);
489 }
490
491 template <class CnstList> void BmfSystem::bmf_solve(const CnstList& cnst_list)
492 {
493   idx2Var_.clear();
494   cnst2idx_.clear();
495   Eigen::MatrixXd A;
496   Eigen::MatrixXd maxA;
497   Eigen::VectorXd C;
498   Eigen::VectorXd bounds;
499   std::vector<bool> shared;
500   get_constraint_data(cnst_list, C, shared);
501   get_flows_data(C.size(), A, maxA, bounds);
502
503   auto solver = BmfSolver(std::move(A), std::move(maxA), std::move(C), std::move(shared), std::move(bounds));
504   auto rho    = solver.solve();
505
506   if (rho.size() == 0)
507     return;
508
509   /* setting rhos */
510   for (int i = 0; i < rho.size(); i++) {
511     idx2Var_[i]->value_ = rho[i];
512   }
513 }
514
515 } // namespace simgrid::kernel::lmm