Logo AND Algorithmique Numérique Distribuée

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