Logo AND Algorithmique Numérique Distribuée

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