1 /* Copyright (c) 2007-2022. 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 "xbt/config.hpp"
14 XBT_LOG_NEW_DEFAULT_SUBCATEGORY(ker_bmf, kernel, "Kernel BMF solver");
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);
20 simgrid::config::Flag<double> cfg_bmf_precision{"bmf/precision",
21 "Numerical precision used when computing resource sharing", 1E-12};
27 AllocationGenerator::AllocationGenerator(Eigen::MatrixXd A) : A_(std::move(A)), alloc_(A_.cols(), 0)
29 // got a first valid allocation
30 for (size_t p = 0; p < alloc_.size(); p++) {
31 for (int r = 0; r < A_.rows(); r++) {
40 bool AllocationGenerator::next(std::vector<int>& next_alloc)
48 auto n_resources = A_.rows();
50 while (idx < alloc_.size()) {
51 alloc_[idx] = (alloc_[idx] + 1) % n_resources;
52 if (alloc_[idx] == 0) {
58 if (A_(alloc_[idx], idx) > 0) {
66 /*****************************************************************************/
68 BmfSolver::BmfSolver(Eigen::MatrixXd A, Eigen::MatrixXd maxA, Eigen::VectorXd C, std::vector<bool> shared,
71 , maxA_(std::move(maxA))
73 , C_shared_(std::move(shared))
74 , phi_(std::move(phi))
76 , max_iteration_(cfg_bmf_max_iteration)
79 xbt_assert(max_iteration_ > 0,
80 "Invalid number of iterations for BMF solver. Please check your \"bmf/max-iterations\" configuration.");
81 xbt_assert(A_.cols() == maxA_.cols(), "Invalid number of cols in matrix A (%td) or maxA (%td)", A_.cols(),
83 xbt_assert(A_.cols() == phi_.size(), "Invalid size of phi vector (%td)", phi_.size());
84 xbt_assert(static_cast<long>(C_shared_.size()) == C_.size(), "Invalid size param shared (%zu)", C_shared_.size());
87 template <typename T> std::string BmfSolver::debug_eigen(const T& obj) const
89 std::stringstream debug;
94 template <typename C> std::string BmfSolver::debug_vector(const C& container) const
96 std::stringstream debug;
97 std::copy(container.begin(), container.end(),
98 std::ostream_iterator<typename std::remove_reference<decltype(container)>::type::value_type>(debug, " "));
102 std::string BmfSolver::debug_alloc(const allocation_map_t& alloc) const
104 std::stringstream debug;
105 for (const auto& e : alloc) {
106 debug << "{" + std::to_string(e.first) + ": [" + debug_vector(e.second) + "]}, ";
111 double BmfSolver::get_resource_capacity(int resource, const std::vector<int>& bounded_players) const
113 double capacity = C_[resource];
114 if (not C_shared_[resource])
117 for (int p : bounded_players) {
118 capacity -= A_(resource, p) * phi_[p];
120 return std::max(0.0, capacity);
123 double BmfSolver::get_maxmin_share(int resource, const std::vector<int>& bounded_players) const
125 auto n_players = (A_.row(resource).array() > 0).count() - bounded_players.size();
126 double capacity = get_resource_capacity(resource, bounded_players);
128 capacity /= n_players;
132 std::vector<int> BmfSolver::alloc_map_to_vector(const allocation_map_t& alloc) const
134 std::vector<int> alloc_by_player(A_.cols(), -1);
135 for (const auto& it : alloc) {
136 for (auto p : it.second) {
137 alloc_by_player[p] = it.first;
140 return alloc_by_player;
143 std::vector<int> BmfSolver::get_bounded_players(const allocation_map_t& alloc) const
145 std::vector<int> bounded_players;
146 for (const auto& e : alloc) {
147 if (e.first == NO_RESOURCE) {
148 bounded_players.insert(bounded_players.end(), e.second.begin(), e.second.end());
151 return bounded_players;
154 Eigen::VectorXd BmfSolver::equilibrium(const allocation_map_t& alloc) const
156 auto n_players = A_.cols();
157 Eigen::MatrixXd A_p = Eigen::MatrixXd::Zero(n_players, n_players); // square matrix with number of players
158 Eigen::VectorXd C_p = Eigen::VectorXd::Zero(n_players);
161 auto bounded_players = get_bounded_players(alloc);
162 for (const auto& e : alloc) {
163 // add one row for the resource with A[r,]
164 int cur_resource = e.first;
165 /* bounded players, nothing to do */
166 if (cur_resource == NO_RESOURCE)
168 /* not shared resource, each player can receive the full capacity of the resource */
169 if (not C_shared_[cur_resource]) {
170 for (int i : e.second) {
171 C_p[row] = get_resource_capacity(cur_resource, bounded_players);
172 A_p(row, i) = A_(cur_resource, i);
178 /* shared resource: fairly share it between players */
179 A_p.row(row) = A_.row(cur_resource);
180 C_p[row] = get_resource_capacity(cur_resource, bounded_players);
182 if (e.second.size() > 1) {
183 // if 2 players have chosen the same resource
184 // they must have a fair sharing of this resource, adjust A_p and C_p accordingly
185 auto it = e.second.begin();
186 int i = *it; // first player
187 /* for each other player sharing this resource */
188 for (++it; it != e.second.end(); ++it) {
189 /* player i and k on this resource j: so maxA_ji*rho_i - maxA_jk*rho_k = 0 */
192 A_p(row, i) = maxA_(cur_resource, i);
193 A_p(row, k) = -maxA_(cur_resource, k);
198 /* clear players which are externally bounded */
199 for (int p : bounded_players) {
200 A_p.col(p).setZero();
203 XBT_DEBUG("A':\n%s", debug_eigen(A_p).c_str());
205 XBT_DEBUG("C':\n%s", debug_eigen(C_p).c_str());
206 /* PartialPivLU is much faster than FullPivLU but requires that the matrix is invertible
207 * FullPivLU however assures that it finds come solution even if the matrix is singular
208 * Ideally we would like to be optimist and try Partial and in case of error, go back
210 * However, this with isNaN doesn't work if compiler uses -Ofastmath. In our case,
211 * the icc compiler raises an error when compiling the code (comparison with NaN always evaluates to false in fast
212 * floating point modes).
213 * Eigen::VectorXd rho = Eigen::PartialPivLU<Eigen::MatrixXd>(A_p).solve(C_p);
214 * if (rho.array().isNaN().any()) {
215 * XBT_DEBUG("rho with nan values, falling back to FullPivLU, rho:\n%s", debug_eigen(rho).c_str());
216 * rho = Eigen::FullPivLU<Eigen::MatrixXd>(A_p).solve(C_p);
220 Eigen::VectorXd rho = Eigen::FullPivLU<Eigen::MatrixXd>(A_p).solve(C_p);
221 for (int p : bounded_players) {
227 bool BmfSolver::disturb_allocation(allocation_map_t& alloc, std::vector<int>& alloc_by_player)
229 while (gen_.next(alloc_by_player)) {
230 if (allocations_.find(alloc_by_player) == allocations_.end()) {
231 allocations_.clear();
232 allocations_.insert(alloc_by_player);
234 for (size_t p = 0; p < alloc_by_player.size(); p++) {
235 alloc[alloc_by_player[p]].insert(p);
243 bool BmfSolver::get_alloc(const Eigen::VectorXd& fair_sharing, const allocation_map_t& last_alloc,
244 allocation_map_t& alloc, bool initial)
247 for (int player_idx = 0; player_idx < A_.cols(); player_idx++) {
248 int selected_resource = NO_RESOURCE;
250 /* the player's maximal rate is the minimum among all resources */
251 double min_rate = -1;
252 for (int cnst_idx = 0; cnst_idx < A_.rows(); cnst_idx++) {
253 if (A_(cnst_idx, player_idx) <= 0.0)
256 /* Note: the max_ may artificially increase the rate if priority < 0
257 * The equilibrium sets a rho which respects the C_ though */
258 double rate = fair_sharing[cnst_idx] / maxA_(cnst_idx, player_idx);
259 if (min_rate == -1 || double_positive(min_rate - rate, cfg_bmf_precision)) {
260 selected_resource = cnst_idx;
263 double bound = initial ? -1 : phi_[player_idx];
264 /* Given that the priority may artificially increase the rate,
265 * we need to check that the bound given by user respects the resource capacity C_ */
266 if (bound > 0 && bound * A_(cnst_idx, player_idx) < C_[cnst_idx] &&
267 double_positive(min_rate - bound, cfg_bmf_precision)) {
268 selected_resource = NO_RESOURCE;
272 alloc[selected_resource].insert(player_idx);
274 bool is_stable = (alloc == last_alloc);
278 std::vector<int> alloc_by_player = alloc_map_to_vector(alloc);
279 bool inserted = allocations_.insert(alloc_by_player).second;
280 /* oops, allocation already tried, let's pertube it a bit */
282 XBT_DEBUG("Allocation already tried: %s", debug_alloc(alloc).c_str());
283 return disturb_allocation(alloc, alloc_by_player);
288 void BmfSolver::set_fair_sharing(const allocation_map_t& alloc, const Eigen::VectorXd& rho,
289 Eigen::VectorXd& fair_sharing) const
291 std::vector<int> bounded_players = get_bounded_players(alloc);
293 for (int r = 0; r < fair_sharing.size(); r++) {
294 auto it = alloc.find(r);
295 if (it != alloc.end()) { // resource selected by some player, fair share depends on rho
296 double min_share = std::numeric_limits<double>::max();
297 for (int p : it->second) {
298 double share = A_(r, p) * rho[p];
299 min_share = std::min(min_share, share);
301 fair_sharing[r] = min_share;
302 } else { // nobody selects this resource, fair_sharing depends on resource saturation
303 // resource r is saturated (A[r,*] * rho > C), divide it among players
304 double consumption_r = A_.row(r) * rho;
305 double_update(&consumption_r, C_[r], cfg_bmf_precision);
306 if (consumption_r > 0.0) {
307 fair_sharing[r] = get_maxmin_share(r, bounded_players);
309 fair_sharing[r] = C_[r];
315 bool BmfSolver::is_bmf(const Eigen::VectorXd& rho) const
319 // 1) the capacity of all resources is respected
320 Eigen::VectorXd shared(C_shared_.size());
321 for (int j = 0; j < shared.size(); j++)
322 shared[j] = C_shared_[j] ? 1.0 : 0.0;
324 Eigen::VectorXd remaining = (A_ * rho) - C_;
325 remaining = remaining.array() * shared.array(); // ignore non shared resources
326 bmf = bmf && (not std::any_of(remaining.data(), remaining.data() + remaining.size(),
327 [](double v) { return double_positive(v, sg_maxmin_precision); }));
329 // 3) every player receives maximum share in at least 1 saturated resource
330 // due to subflows, compare with the maximum consumption and not the A matrix
331 Eigen::MatrixXd usage =
332 maxA_.array().rowwise() * rho.transpose().array(); // usage_ji: indicates the usage of player i on resource j
334 XBT_DEBUG("Usage_ji considering max consumption:\n%s", debug_eigen(usage).c_str());
335 auto max_share = usage.rowwise().maxCoeff(); // max share for each resource j
337 // matrix_ji: boolean indicating player p has the maximum share at resource j
338 Eigen::MatrixXi player_max_share =
339 ((usage.array().colwise() - max_share.array()).abs() <= sg_maxmin_precision).cast<int>();
340 // but only saturated resources must be considered
341 Eigen::VectorXi saturated = (remaining.array().abs() <= sg_maxmin_precision).cast<int>();
342 XBT_DEBUG("Saturated_j resources:\n%s", debug_eigen(saturated).c_str());
343 player_max_share.array().colwise() *= saturated.array();
345 // just check if it has received at least it's bound
346 for (int p = 0; p < rho.size(); p++) {
347 if (double_equals(rho[p], phi_[p], sg_maxmin_precision)) {
348 player_max_share(0, p) = 1; // it doesn't really matter, just to say that it's a bmf
353 // 2) at least 1 resource is saturated
354 bmf = bmf && (saturated.array() == 1).any();
356 XBT_DEBUG("Player_ji usage of saturated resources:\n%s", debug_eigen(player_max_share).c_str());
357 // for all columns(players) it has to be the max at least in 1
358 bmf = bmf && (player_max_share.colwise().sum().array() >= 1).all();
362 Eigen::VectorXd BmfSolver::solve()
364 XBT_DEBUG("Starting BMF solver");
366 XBT_DEBUG("A:\n%s", debug_eigen(A_).c_str());
367 XBT_DEBUG("maxA:\n%s", debug_eigen(maxA_).c_str());
368 XBT_DEBUG("C:\n%s", debug_eigen(C_).c_str());
369 XBT_DEBUG("phi:\n%s", debug_eigen(phi_).c_str());
371 /* no flows to share, just returns */
376 auto fair_sharing = C_;
378 /* BMF allocation for each player (current and last one) stop when are equal */
379 allocation_map_t last_alloc;
380 allocation_map_t cur_alloc;
383 while (it < max_iteration_ && not get_alloc(fair_sharing, last_alloc, cur_alloc, it == 0)) {
384 last_alloc = cur_alloc;
385 XBT_DEBUG("BMF: iteration %d", it);
386 XBT_DEBUG("B (current allocation): %s", debug_alloc(cur_alloc).c_str());
388 // solve inv(A)*rho = C
389 rho = equilibrium(cur_alloc);
390 XBT_DEBUG("rho:\n%s", debug_eigen(rho).c_str());
392 // get fair sharing for each resource
393 set_fair_sharing(cur_alloc, rho, fair_sharing);
394 XBT_DEBUG("Fair sharing vector (per resource):\n%s", debug_eigen(fair_sharing).c_str());
396 // get new allocation for players
400 /* Not mandatory but a safe check to assure we have a proper solution */
401 if (not is_bmf(rho)) {
402 fprintf(stderr, "Unable to find a BMF allocation for your system.\n"
403 "You may try to increase the maximum number of iterations performed by BMF solver "
404 "(\"--cfg=bmf/max-iterations\").\n"
405 "Additionally, you could adjust numerical precision (\"--cfg=bmf/precision\").\n");
406 fprintf(stderr, "Internal states (after %d iterations):\n", it);
407 fprintf(stderr, "A:\n%s\n", debug_eigen(A_).c_str());
408 fprintf(stderr, "maxA:\n%s\n", debug_eigen(maxA_).c_str());
409 fprintf(stderr, "C:\n%s\n", debug_eigen(C_).c_str());
410 fprintf(stderr, "C_shared:\n%s\n", debug_vector(C_shared_).c_str());
411 fprintf(stderr, "phi:\n%s\n", debug_eigen(phi_).c_str());
412 fprintf(stderr, "rho:\n%s\n", debug_eigen(rho).c_str());
416 XBT_DEBUG("BMF done after %d iterations", it);
420 /*****************************************************************************/
422 void BmfSystem::get_flows_data(Eigen::Index number_cnsts, Eigen::MatrixXd& A, Eigen::MatrixXd& maxA,
423 Eigen::VectorXd& phi)
425 A.resize(number_cnsts, variable_set.size());
427 maxA.resize(number_cnsts, variable_set.size());
429 phi.resize(variable_set.size());
432 for (Variable& var : variable_set) {
433 if (var.sharing_penalty_ <= 0)
436 bool linked = false; // variable is linked to some constraint (specially for selective_update)
437 for (const Element& elem : var.cnsts_) {
438 const boost::intrusive::list_member_hook<>& cnst_hook = selective_update_active
439 ? elem.constraint->modified_constraint_set_hook_
440 : elem.constraint->active_constraint_set_hook_;
441 if (not cnst_hook.is_linked())
443 /* active and linked variable, lets check its consumption */
445 double consumption = elem.consumption_weight;
446 if (consumption > 0) {
447 int cnst_idx = cnst2idx_[elem.constraint];
448 A(cnst_idx, var_idx) += consumption;
449 // a variable with double penalty must receive half share, so it max weight is greater
450 maxA(cnst_idx, var_idx) = std::max(maxA(cnst_idx, var_idx), elem.max_consumption_weight * var.sharing_penalty_);
454 /* skip variables not linked to any modified or active constraint */
458 phi[var_idx] = var.get_bound();
459 idx2Var_[var_idx] = &var;
462 var.value_ = 1; // assign something by default for tasks with 0 consumption
465 // resize matrix to active variables only
466 A.conservativeResize(Eigen::NoChange_t::NoChange, var_idx);
467 maxA.conservativeResize(Eigen::NoChange_t::NoChange, var_idx);
468 phi.conservativeResize(var_idx);
471 template <class CnstList>
472 void BmfSystem::get_constraint_data(const CnstList& cnst_list, Eigen::VectorXd& C, std::vector<bool>& shared)
474 C.resize(cnst_list.size());
475 shared.resize(cnst_list.size());
478 for (const Constraint& cnst : cnst_list) {
479 C(cnst_idx) = cnst.bound_;
480 if (cnst.get_sharing_policy() == Constraint::SharingPolicy::NONLINEAR && cnst.dyn_constraint_cb_) {
481 C(cnst_idx) = cnst.dyn_constraint_cb_(cnst.bound_, cnst.concurrency_current_);
483 cnst2idx_[&cnst] = cnst_idx;
484 // FATPIPE links aren't really shared
485 shared[cnst_idx] = (cnst.sharing_policy_ != Constraint::SharingPolicy::FATPIPE);
490 void BmfSystem::do_solve()
492 if (selective_update_active)
493 bmf_solve(modified_constraint_set);
495 bmf_solve(active_constraint_set);
498 template <class CnstList> void BmfSystem::bmf_solve(const CnstList& cnst_list)
503 Eigen::MatrixXd maxA;
505 Eigen::VectorXd bounds;
506 std::vector<bool> shared;
507 get_constraint_data(cnst_list, C, shared);
508 get_flows_data(C.size(), A, maxA, bounds);
510 auto solver = BmfSolver(std::move(A), std::move(maxA), std::move(C), std::move(shared), std::move(bounds));
511 auto rho = solver.solve();
517 for (int i = 0; i < rho.size(); i++) {
518 idx2Var_[i]->value_ = rho[i];
523 } // namespace kernel
524 } // namespace simgrid