1 /* Copyright (c) 2007-2023. The SimGrid Team. All rights reserved. */
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. */
6 #include "src/kernel/lmm/bmf.hpp"
7 #include "simgrid/math_utils.h"
14 XBT_LOG_NEW_DEFAULT_SUBCATEGORY(ker_bmf, kernel, "Kernel BMF solver");
16 namespace simgrid::kernel::lmm {
18 AllocationGenerator::AllocationGenerator(Eigen::MatrixXd A) : A_(std::move(A)), alloc_(A_.cols(), 0)
20 // got a first valid allocation
21 for (size_t p = 0; p < alloc_.size(); p++) {
22 for (int r = 0; r < A_.rows(); r++) {
31 bool AllocationGenerator::next(std::vector<int>& next_alloc)
39 auto n_resources = A_.rows();
41 while (idx < alloc_.size()) {
42 alloc_[idx] = (alloc_[idx] + 1) % n_resources;
43 if (alloc_[idx] == 0) {
49 if (A_(alloc_[idx], idx) > 0) {
57 /*****************************************************************************/
59 BmfSolver::BmfSolver(Eigen::MatrixXd A, Eigen::MatrixXd maxA, Eigen::VectorXd C, std::vector<bool> shared,
62 , maxA_(std::move(maxA))
64 , C_shared_(std::move(shared))
65 , phi_(std::move(phi))
69 xbt_assert(max_iteration_ > 0,
70 "Invalid number of iterations for BMF solver. Please check your \"bmf/max-iterations\" configuration.");
71 xbt_assert(A_.cols() == maxA_.cols(), "Invalid number of cols in matrix A (%td) or maxA (%td)", A_.cols(),
73 xbt_assert(A_.cols() == phi_.size(), "Invalid size of phi vector (%td)", phi_.size());
74 xbt_assert(static_cast<long>(C_shared_.size()) == C_.size(), "Invalid size param shared (%zu)", C_shared_.size());
77 template <typename T> std::string BmfSolver::debug_eigen(const T& obj) const
79 std::stringstream debug;
84 template <typename C> std::string BmfSolver::debug_vector(const C& container) const
86 std::stringstream debug;
87 std::copy(container.begin(), container.end(),
88 std::ostream_iterator<typename std::remove_reference<decltype(container)>::type::value_type>(debug, " "));
92 std::string BmfSolver::debug_alloc(const allocation_map_t& alloc) const
94 std::stringstream debug;
95 for (const auto& [resource, players] : alloc) {
96 debug << "{" + std::to_string(resource) + ": [" + debug_vector(players) + "]}, ";
101 double BmfSolver::get_resource_capacity(int resource, const std::vector<int>& bounded_players) const
103 double capacity = C_[resource];
104 if (not C_shared_[resource])
107 for (int p : bounded_players) {
108 capacity -= A_(resource, p) * phi_[p];
110 return std::max(0.0, capacity);
113 double BmfSolver::get_maxmin_share(int resource, const std::vector<int>& bounded_players) const
115 auto n_players = (A_.row(resource).array() > 0).count() - bounded_players.size();
116 double capacity = get_resource_capacity(resource, bounded_players);
118 capacity /= n_players;
122 std::vector<int> BmfSolver::alloc_map_to_vector(const allocation_map_t& alloc) const
124 std::vector<int> alloc_by_player(A_.cols(), -1);
125 for (const auto& [resource, players] : alloc) {
126 for (auto p : players) {
127 alloc_by_player[p] = resource;
130 return alloc_by_player;
133 std::vector<int> BmfSolver::get_bounded_players(const allocation_map_t& alloc) const
135 std::vector<int> bounded_players;
136 for (const auto& [resource, players] : alloc) {
137 if (resource == NO_RESOURCE) {
138 bounded_players.insert(bounded_players.end(), players.begin(), players.end());
141 return bounded_players;
144 Eigen::VectorXd BmfSolver::equilibrium(const allocation_map_t& alloc) const
146 auto n_players = A_.cols();
147 Eigen::MatrixXd A_p = Eigen::MatrixXd::Zero(n_players, n_players); // square matrix with number of players
148 Eigen::VectorXd C_p = Eigen::VectorXd::Zero(n_players);
151 auto bounded_players = get_bounded_players(alloc);
152 for (const auto& [resource, players] : alloc) {
153 // add one row for the resource with A[r,]
154 /* bounded players, nothing to do */
155 if (resource == NO_RESOURCE)
157 /* not shared resource, each player can receive the full capacity of the resource */
158 if (not C_shared_[resource]) {
159 for (int i : players) {
160 C_p[row] = get_resource_capacity(resource, bounded_players);
161 A_p(row, i) = A_(resource, i);
167 /* shared resource: fairly share it between players */
168 A_p.row(row) = A_.row(resource);
169 C_p[row] = get_resource_capacity(resource, bounded_players);
171 if (players.size() > 1) {
172 // if 2 players have chosen the same resource
173 // they must have a fair sharing of this resource, adjust A_p and C_p accordingly
174 auto it = players.begin();
175 int i = *it; // first player
176 /* for each other player sharing this resource */
177 for (++it; it != players.end(); ++it) {
178 /* player i and k on this resource j: so maxA_ji*rho_i - maxA_jk*rho_k = 0 */
181 A_p(row, i) = maxA_(resource, i);
182 A_p(row, k) = -maxA_(resource, k);
187 /* clear players which are externally bounded */
188 for (int p : bounded_players) {
189 A_p.col(p).setZero();
192 XBT_DEBUG("A':\n%s", debug_eigen(A_p).c_str());
194 XBT_DEBUG("C':\n%s", debug_eigen(C_p).c_str());
195 /* PartialPivLU is much faster than FullPivLU but requires that the matrix is invertible
196 * FullPivLU however assures that it finds come solution even if the matrix is singular
197 * Ideally we would like to be optimist and try Partial and in case of error, go back
199 * However, this with isNaN doesn't work if compiler uses -Ofastmath. In our case,
200 * the icc compiler raises an error when compiling the code (comparison with NaN always evaluates to false in fast
201 * floating point modes).
202 * Eigen::VectorXd rho = Eigen::PartialPivLU<Eigen::MatrixXd>(A_p).solve(C_p);
203 * if (rho.array().isNaN().any()) {
204 * XBT_DEBUG("rho with nan values, falling back to FullPivLU, rho:\n%s", debug_eigen(rho).c_str());
205 * rho = Eigen::FullPivLU<Eigen::MatrixXd>(A_p).solve(C_p);
209 Eigen::VectorXd rho = Eigen::FullPivLU<Eigen::MatrixXd>(A_p).solve(C_p);
210 for (int p : bounded_players) {
216 bool BmfSolver::disturb_allocation(allocation_map_t& alloc, std::vector<int>& alloc_by_player)
218 while (gen_.next(alloc_by_player)) {
219 if (allocations_.find(alloc_by_player) == allocations_.end()) {
220 allocations_.clear();
221 allocations_.insert(alloc_by_player);
223 for (size_t p = 0; p < alloc_by_player.size(); p++) {
224 alloc[alloc_by_player[p]].insert(p);
232 bool BmfSolver::get_alloc(const Eigen::VectorXd& fair_sharing, const allocation_map_t& last_alloc,
233 allocation_map_t& alloc, bool initial)
236 for (int player_idx = 0; player_idx < A_.cols(); player_idx++) {
237 int selected_resource = NO_RESOURCE;
239 /* the player's maximal rate is the minimum among all resources */
240 double min_rate = -1;
241 for (int cnst_idx = 0; cnst_idx < A_.rows(); cnst_idx++) {
242 if (A_(cnst_idx, player_idx) <= 0.0)
245 /* Note: the max_ may artificially increase the rate if priority < 0
246 * The equilibrium sets a rho which respects the C_ though */
247 if (double rate = fair_sharing[cnst_idx] / maxA_(cnst_idx, player_idx);
248 min_rate == -1 || double_positive(min_rate - rate, cfg_bmf_precision)) {
249 selected_resource = cnst_idx;
252 /* Given that the priority may artificially increase the rate,
253 * we need to check that the bound given by user respects the resource capacity C_ */
254 if (double bound = initial ? -1 : phi_[player_idx]; bound > 0 &&
255 bound * A_(cnst_idx, player_idx) < C_[cnst_idx] &&
256 double_positive(min_rate - bound, cfg_bmf_precision)) {
257 selected_resource = NO_RESOURCE;
261 alloc[selected_resource].insert(player_idx);
263 if (alloc == last_alloc) // considered stable
266 if (auto alloc_by_player = alloc_map_to_vector(alloc); not allocations_.insert(alloc_by_player).second) {
267 /* oops, allocation already tried, let's pertube it a bit */
268 XBT_DEBUG("Allocation already tried: %s", debug_alloc(alloc).c_str());
269 return disturb_allocation(alloc, alloc_by_player);
274 void BmfSolver::set_fair_sharing(const allocation_map_t& alloc, const Eigen::VectorXd& rho,
275 Eigen::VectorXd& fair_sharing) const
277 std::vector<int> bounded_players = get_bounded_players(alloc);
279 for (int r = 0; r < fair_sharing.size(); r++) {
280 auto it = alloc.find(r);
281 if (it != alloc.end()) { // resource selected by some player, fair share depends on rho
282 double min_share = std::numeric_limits<double>::max();
283 for (int p : it->second) {
284 double share = A_(r, p) * rho[p];
285 min_share = std::min(min_share, share);
287 fair_sharing[r] = min_share;
288 } else { // nobody selects this resource, fair_sharing depends on resource saturation
289 // resource r is saturated (A[r,*] * rho > C), divide it among players
290 double consumption_r = A_.row(r) * rho;
291 double_update(&consumption_r, C_[r], cfg_bmf_precision);
292 if (consumption_r > 0.0) {
293 fair_sharing[r] = get_maxmin_share(r, bounded_players);
295 fair_sharing[r] = C_[r];
301 bool BmfSolver::is_bmf(const Eigen::VectorXd& rho) const
305 // 1) the capacity of all resources is respected
306 Eigen::VectorXd shared(C_shared_.size());
307 for (int j = 0; j < shared.size(); j++)
308 shared[j] = C_shared_[j] ? 1.0 : 0.0;
310 Eigen::VectorXd remaining = (A_ * rho) - C_;
311 remaining = remaining.array() * shared.array(); // ignore non shared resources
312 bmf = bmf && (not std::any_of(remaining.data(), remaining.data() + remaining.size(),
313 [](double v) { return double_positive(v, sg_precision_workamount); }));
315 // 3) every player receives maximum share in at least 1 saturated resource
316 // due to subflows, compare with the maximum consumption and not the A matrix
317 Eigen::MatrixXd usage =
318 maxA_.array().rowwise() * rho.transpose().array(); // usage_ji: indicates the usage of player i on resource j
320 XBT_DEBUG("Usage_ji considering max consumption:\n%s", debug_eigen(usage).c_str());
321 auto max_share = usage.rowwise().maxCoeff(); // max share for each resource j
323 // matrix_ji: boolean indicating player p has the maximum share at resource j
324 Eigen::MatrixXi player_max_share =
325 ((usage.array().colwise() - max_share.array()).abs() <= sg_precision_workamount).cast<int>();
326 // but only saturated resources must be considered
327 Eigen::VectorXi saturated = (remaining.array().abs() <= sg_precision_workamount).cast<int>();
328 XBT_DEBUG("Saturated_j resources:\n%s", debug_eigen(saturated).c_str());
329 player_max_share.array().colwise() *= saturated.array();
331 // just check if it has received at least it's bound
332 for (int p = 0; p < rho.size(); p++) {
333 if (double_equals(rho[p], phi_[p], sg_precision_workamount)) {
334 player_max_share(0, p) = 1; // it doesn't really matter, just to say that it's a bmf
339 // 2) at least 1 resource is saturated
340 bmf = bmf && (saturated.array() == 1).any();
342 XBT_DEBUG("Player_ji usage of saturated resources:\n%s", debug_eigen(player_max_share).c_str());
343 // for all columns(players) it has to be the max at least in 1
344 bmf = bmf && (player_max_share.colwise().sum().array() >= 1).all();
348 Eigen::VectorXd BmfSolver::solve()
350 XBT_DEBUG("Starting BMF solver");
352 XBT_DEBUG("A:\n%s", debug_eigen(A_).c_str());
353 XBT_DEBUG("maxA:\n%s", debug_eigen(maxA_).c_str());
354 XBT_DEBUG("C:\n%s", debug_eigen(C_).c_str());
355 XBT_DEBUG("phi:\n%s", debug_eigen(phi_).c_str());
357 /* no flows to share, just returns */
362 auto fair_sharing = C_;
364 /* BMF allocation for each player (current and last one) stop when are equal */
365 allocation_map_t last_alloc;
366 allocation_map_t cur_alloc;
369 while (it < max_iteration_ && not get_alloc(fair_sharing, last_alloc, cur_alloc, it == 0)) {
370 last_alloc = cur_alloc;
371 XBT_DEBUG("BMF: iteration %d", it);
372 XBT_DEBUG("B (current allocation): %s", debug_alloc(cur_alloc).c_str());
374 // solve inv(A)*rho = C
375 rho = equilibrium(cur_alloc);
376 XBT_DEBUG("rho:\n%s", debug_eigen(rho).c_str());
378 // get fair sharing for each resource
379 set_fair_sharing(cur_alloc, rho, fair_sharing);
380 XBT_DEBUG("Fair sharing vector (per resource):\n%s", debug_eigen(fair_sharing).c_str());
382 // get new allocation for players
386 /* Not mandatory but a safe check to assure we have a proper solution */
387 if (not is_bmf(rho)) {
388 fprintf(stderr, "Unable to find a BMF allocation for your system.\n"
389 "You may try to increase the maximum number of iterations performed by BMF solver "
390 "(\"--cfg=bmf/max-iterations\").\n"
391 "Additionally, you could adjust numerical precision (\"--cfg=bmf/precision\").\n");
392 fprintf(stderr, "Internal states (after %d iterations):\n", it);
393 fprintf(stderr, "A:\n%s\n", debug_eigen(A_).c_str());
394 fprintf(stderr, "maxA:\n%s\n", debug_eigen(maxA_).c_str());
395 fprintf(stderr, "C:\n%s\n", debug_eigen(C_).c_str());
396 fprintf(stderr, "C_shared:\n%s\n", debug_vector(C_shared_).c_str());
397 fprintf(stderr, "phi:\n%s\n", debug_eigen(phi_).c_str());
398 fprintf(stderr, "rho:\n%s\n", debug_eigen(rho).c_str());
402 XBT_DEBUG("BMF done after %d iterations", it);
406 /*****************************************************************************/
408 void BmfSystem::get_flows_data(Eigen::Index number_cnsts, Eigen::MatrixXd& A, Eigen::MatrixXd& maxA,
409 Eigen::VectorXd& phi)
411 A.resize(number_cnsts, variable_set.size());
413 maxA.resize(number_cnsts, variable_set.size());
415 phi.resize(variable_set.size());
418 for (Variable& var : variable_set) {
419 if (var.sharing_penalty_ <= 0)
422 bool linked = false; // variable is linked to some constraint (specially for selective_update)
423 for (const Element& elem : var.cnsts_) {
424 if (const auto& cnst_hook = selective_update_active ? elem.constraint->modified_constraint_set_hook_
425 : elem.constraint->active_constraint_set_hook_;
426 not cnst_hook.is_linked())
428 /* active and linked variable, lets check its consumption */
430 double consumption = elem.consumption_weight;
431 if (consumption > 0) {
432 int cnst_idx = cnst2idx_[elem.constraint];
433 A(cnst_idx, var_idx) += consumption;
434 // a variable with double penalty must receive half share, so it max weight is greater
435 maxA(cnst_idx, var_idx) = std::max(maxA(cnst_idx, var_idx), elem.max_consumption_weight * var.sharing_penalty_);
439 /* skip variables not linked to any modified or active constraint */
443 phi[var_idx] = var.get_bound();
444 idx2Var_[var_idx] = &var;
447 var.value_ = 1; // assign something by default for tasks with 0 consumption
450 // resize matrix to active variables only
451 A.conservativeResize(Eigen::NoChange_t::NoChange, var_idx);
452 maxA.conservativeResize(Eigen::NoChange_t::NoChange, var_idx);
453 phi.conservativeResize(var_idx);
456 template <class CnstList>
457 void BmfSystem::get_constraint_data(const CnstList& cnst_list, Eigen::VectorXd& C, std::vector<bool>& shared)
459 C.resize(cnst_list.size());
460 shared.resize(cnst_list.size());
463 for (const Constraint& cnst : cnst_list) {
464 C(cnst_idx) = cnst.bound_;
465 if (cnst.get_sharing_policy() == Constraint::SharingPolicy::NONLINEAR && cnst.dyn_constraint_cb_) {
466 C(cnst_idx) = cnst.dyn_constraint_cb_(cnst.bound_, cnst.concurrency_current_);
468 cnst2idx_[&cnst] = cnst_idx;
469 // FATPIPE links aren't really shared
470 shared[cnst_idx] = (cnst.sharing_policy_ != Constraint::SharingPolicy::FATPIPE);
475 void BmfSystem::do_solve()
477 if (selective_update_active)
478 bmf_solve(modified_constraint_set);
480 bmf_solve(active_constraint_set);
483 template <class CnstList> void BmfSystem::bmf_solve(const CnstList& cnst_list)
488 Eigen::MatrixXd maxA;
490 Eigen::VectorXd bounds;
491 std::vector<bool> shared;
492 get_constraint_data(cnst_list, C, shared);
493 get_flows_data(C.size(), A, maxA, bounds);
495 auto solver = BmfSolver(std::move(A), std::move(maxA), std::move(C), std::move(shared), std::move(bounds));
496 auto rho = solver.solve();
502 for (int i = 0; i < rho.size(); i++) {
503 idx2Var_[i]->value_ = rho[i];
507 } // namespace simgrid::kernel::lmm