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<bool> cfg_bmf_selective_update{
21 "bmf/selective-update", "Update the constraint set propagating recursively to others constraints (off by default)",
28 AllocationGenerator::AllocationGenerator(Eigen::MatrixXd A) : A_(std::move(A)), alloc_(A_.cols(), 0)
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++) {
41 bool AllocationGenerator::next(std::vector<int>& next_alloc)
49 auto n_resources = A_.rows();
51 while (idx < alloc_.size()) {
52 alloc_[idx] = (alloc_[idx] + 1) % n_resources;
53 if (alloc_[idx] == 0) {
59 if (A_(alloc_[idx], idx) > 0) {
67 /*****************************************************************************/
69 BmfSolver::BmfSolver(Eigen::MatrixXd A, Eigen::MatrixXd maxA, Eigen::VectorXd C, std::vector<bool> shared,
72 , maxA_(std::move(maxA))
74 , C_shared_(std::move(shared))
75 , phi_(std::move(phi))
77 , max_iteration_(cfg_bmf_max_iteration)
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(),
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());
88 template <typename T> std::string BmfSolver::debug_eigen(const T& obj) const
90 std::stringstream debug;
95 template <typename C> std::string BmfSolver::debug_vector(const C& container) const
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, " "));
103 std::string BmfSolver::debug_alloc(const allocation_map_t& alloc) const
105 std::stringstream debug;
106 for (const auto& e : alloc) {
107 debug << "{" + std::to_string(e.first) + ": [" + debug_vector(e.second) + "]}, ";
112 double BmfSolver::get_resource_capacity(int resource, const std::vector<int>& bounded_players) const
114 double capacity = C_[resource];
115 if (not C_shared_[resource])
118 for (int p : bounded_players) {
119 capacity -= A_(resource, p) * phi_[p];
121 return std::max(0.0, capacity);
124 double BmfSolver::get_maxmin_share(int resource) const
126 auto n_players = (A_.row(resource).array() > 0).count();
127 return C_[resource] / n_players;
130 std::vector<int> BmfSolver::alloc_map_to_vector(const allocation_map_t& alloc) const
132 std::vector<int> alloc_by_player(A_.cols(), -1);
133 for (const auto& it : alloc) {
134 for (auto p : it.second) {
135 alloc_by_player[p] = it.first;
138 return alloc_by_player;
141 std::vector<int> BmfSolver::get_bounded_players(const allocation_map_t& alloc) const
143 std::vector<int> bounded_players;
144 for (const auto& e : alloc) {
145 if (e.first == NO_RESOURCE) {
146 bounded_players.insert(bounded_players.end(), e.second.begin(), e.second.end());
149 return bounded_players;
152 Eigen::VectorXd BmfSolver::equilibrium(const allocation_map_t& alloc) const
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);
159 auto bounded_players = get_bounded_players(alloc);
160 for (const auto& e : alloc) {
161 // add one row for the resource with A[r,]
162 int cur_resource = e.first;
163 /* bounded players, nothing to do */
164 if (cur_resource == NO_RESOURCE)
166 /* not shared resource, each player can receive the full capacity of the resource */
167 if (not C_shared_[cur_resource]) {
168 for (int i : e.second) {
169 C_p[row] = get_resource_capacity(cur_resource, bounded_players);
170 A_p(row, i) = A_(cur_resource, i);
176 /* shared resource: fairly share it between players */
177 A_p.row(row) = A_.row(cur_resource);
178 C_p[row] = get_resource_capacity(cur_resource, bounded_players);
180 if (e.second.size() > 1) {
181 // if 2 players have chosen the same resource
182 // they must have a fair sharing of this resource, adjust A_p and C_p accordingly
183 auto it = e.second.begin();
184 int i = *it; // first player
185 /* for each other player sharing this resource */
186 for (++it; it != e.second.end(); ++it) {
187 /* player i and k on this resource j: so maxA_ji*rho_i - maxA_jk*rho_k = 0 */
190 A_p(row, i) = maxA_(cur_resource, i);
191 A_p(row, k) = -maxA_(cur_resource, k);
196 /* clear players which are externally bounded */
197 for (int p : bounded_players) {
198 A_p.col(p).setZero();
201 XBT_DEBUG("A':\n%s", debug_eigen(A_p).c_str());
203 XBT_DEBUG("C':\n%s", debug_eigen(C_p).c_str());
204 /* PartialPivLU is much faster than FullPivLU but requires that the matrix is invertible
205 * FullPivLU however assures that it finds come solution even if the matrix is singular
206 * Ideally we would like to be optimist and try Partial and in case of error, go back
208 * However, this with isNaN doesn't work if compiler uses -Ofastmath. In our case,
209 * the icc compiler raises an error when compiling the code (comparison with NaN always evaluates to false in fast
210 * floating point modes).
211 * Eigen::VectorXd rho = Eigen::PartialPivLU<Eigen::MatrixXd>(A_p).solve(C_p);
212 * if (rho.array().isNaN().any()) {
213 * XBT_DEBUG("rho with nan values, falling back to FullPivLU, rho:\n%s", debug_eigen(rho).c_str());
214 * rho = Eigen::FullPivLU<Eigen::MatrixXd>(A_p).solve(C_p);
218 Eigen::VectorXd rho = Eigen::FullPivLU<Eigen::MatrixXd>(A_p).solve(C_p);
219 for (int p : bounded_players) {
225 bool BmfSolver::disturb_allocation(allocation_map_t& alloc, std::vector<int>& alloc_by_player)
227 while (gen_.next(alloc_by_player)) {
228 if (allocations_.find(alloc_by_player) == allocations_.end()) {
229 allocations_.clear();
230 allocations_.insert(alloc_by_player);
232 for (size_t p = 0; p < alloc_by_player.size(); p++) {
233 alloc[alloc_by_player[p]].insert(p);
241 bool BmfSolver::get_alloc(const Eigen::VectorXd& fair_sharing, const allocation_map_t& last_alloc,
242 allocation_map_t& alloc, bool initial)
245 for (int player_idx = 0; player_idx < A_.cols(); player_idx++) {
246 int selected_resource = NO_RESOURCE;
247 double bound = phi_[player_idx];
248 /* the player's maximal rate is the minimum among all resources */
249 double min_rate = (bound <= 0 || initial) ? -1 : bound;
250 for (int cnst_idx = 0; cnst_idx < A_.rows(); cnst_idx++) {
251 if (A_(cnst_idx, player_idx) <= 0.0)
254 double rate = fair_sharing[cnst_idx] / maxA_(cnst_idx, player_idx);
255 if (min_rate == -1 || rate < min_rate) {
256 selected_resource = cnst_idx;
260 alloc[selected_resource].insert(player_idx);
262 bool is_stable = (alloc == last_alloc);
266 std::vector<int> alloc_by_player = alloc_map_to_vector(alloc);
267 auto ret = allocations_.insert(alloc_by_player);
268 /* oops, allocation already tried, let's pertube it a bit */
269 if (not ret.second) {
270 XBT_DEBUG("Allocation already tried: %s", debug_alloc(alloc).c_str());
271 return disturb_allocation(alloc, alloc_by_player);
276 void BmfSolver::set_fair_sharing(const allocation_map_t& alloc, const Eigen::VectorXd& rho,
277 Eigen::VectorXd& fair_sharing) const
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 int player = *(it->second.begin()); // equilibrium assures that every player receives the same, use one of them to
283 // calculate the fair sharing for resource r
284 fair_sharing[r] = maxA_(r, player) * rho[player];
285 } else { // nobody selects this resource, fair_sharing depends on resource saturation
286 // resource r is saturated (A[r,*] * rho > C), divide it among players
287 double consumption_r = A_.row(r) * rho;
288 double_update(&consumption_r, C_[r], sg_maxmin_precision);
289 if (consumption_r > 0.0) {
290 fair_sharing[r] = get_maxmin_share(r);
292 fair_sharing[r] = C_[r];
298 bool BmfSolver::is_bmf(const Eigen::VectorXd& rho) const
302 // 1) the capacity of all resources is respected
303 Eigen::VectorXd shared(C_shared_.size());
304 for (int j = 0; j < shared.size(); j++)
305 shared[j] = C_shared_[j] ? 1.0 : 0.0;
307 Eigen::VectorXd remaining = (A_ * rho) - C_;
308 remaining = remaining.array() * shared.array(); // ignore non shared resources
309 bmf = bmf && (not std::any_of(remaining.data(), remaining.data() + remaining.size(),
310 [](double v) { return double_positive(v, sg_maxmin_precision); }));
312 // 3) every player receives maximum share in at least 1 saturated resource
313 // due to subflows, compare with the maximum consumption and not the A matrix
314 Eigen::MatrixXd usage =
315 maxA_.array().rowwise() * rho.transpose().array(); // usage_ji: indicates the usage of player i on resource j
317 XBT_DEBUG("Usage_ji considering max consumption:\n%s", debug_eigen(usage).c_str());
318 auto max_share = usage.rowwise().maxCoeff(); // max share for each resource j
320 // matrix_ji: boolean indicating player p has the maximum share at resource j
321 Eigen::MatrixXi player_max_share =
322 ((usage.array().colwise() - max_share.array()).abs() <= sg_maxmin_precision).cast<int>();
323 // but only saturated resources must be considered
324 Eigen::VectorXi saturated = (remaining.array().abs() <= sg_maxmin_precision).cast<int>();
325 XBT_DEBUG("Saturated_j resources:\n%s", debug_eigen(saturated).c_str());
326 player_max_share.array().colwise() *= saturated.array();
328 // just check if it has received at least it's bound
329 for (int p = 0; p < rho.size(); p++) {
330 if (double_equals(rho[p], phi_[p], sg_maxmin_precision)) {
331 player_max_share(0, p) = 1; // it doesn't really matter, just to say that it's a bmf
336 // 2) at least 1 resource is saturated
337 bmf = bmf && (saturated.array() == 1).any();
339 XBT_DEBUG("Player_ji usage of saturated resources:\n%s", debug_eigen(player_max_share).c_str());
340 // for all columns(players) it has to be the max at least in 1
341 bmf = bmf && (player_max_share.colwise().sum().array() >= 1).all();
345 Eigen::VectorXd BmfSolver::solve()
347 XBT_DEBUG("Starting BMF solver");
349 XBT_DEBUG("A:\n%s", debug_eigen(A_).c_str());
350 XBT_DEBUG("maxA:\n%s", debug_eigen(maxA_).c_str());
351 XBT_DEBUG("C:\n%s", debug_eigen(C_).c_str());
352 XBT_DEBUG("phi:\n%s", debug_eigen(phi_).c_str());
354 /* no flows to share, just returns */
359 auto fair_sharing = C_;
361 /* BMF allocation for each player (current and last one) stop when are equal */
362 allocation_map_t last_alloc;
363 allocation_map_t cur_alloc;
366 while (it < max_iteration_ && not get_alloc(fair_sharing, last_alloc, cur_alloc, it == 0)) {
367 last_alloc = cur_alloc;
368 XBT_DEBUG("BMF: iteration %d", it);
369 XBT_DEBUG("B (current allocation): %s", debug_alloc(cur_alloc).c_str());
371 // solve inv(A)*rho = C
372 rho = equilibrium(cur_alloc);
373 XBT_DEBUG("rho:\n%s", debug_eigen(rho).c_str());
375 // get fair sharing for each resource
376 set_fair_sharing(cur_alloc, rho, fair_sharing);
377 XBT_DEBUG("Fair sharing vector (per resource):\n%s", debug_eigen(fair_sharing).c_str());
379 // get new allocation for players
383 /* Not mandatory but a safe check to assure we have a proper solution */
384 if (not is_bmf(rho)) {
385 fprintf(stderr, "Unable to find a BMF allocation for your system.\n"
386 "You may try to increase the maximum number of iterations performed by BMF solver "
387 "(\"--cfg=bmf/max-iterations\").\n"
388 "Additionally, you could decrease numerical precision (\"--cfg=surf/precision\").\n");
389 fprintf(stderr, "Internal states (after %d iterations):\n", it);
390 fprintf(stderr, "A:\n%s\n", debug_eigen(A_).c_str());
391 fprintf(stderr, "maxA:\n%s\n", debug_eigen(maxA_).c_str());
392 fprintf(stderr, "C:\n%s\n", debug_eigen(C_).c_str());
393 fprintf(stderr, "C_shared:\n%s\n", debug_vector(C_shared_).c_str());
394 fprintf(stderr, "phi:\n%s\n", debug_eigen(phi_).c_str());
395 fprintf(stderr, "rho:\n%s\n", debug_eigen(rho).c_str());
399 XBT_DEBUG("BMF done after %d iterations", it);
403 /*****************************************************************************/
405 void BmfSystem::get_flows_data(Eigen::Index number_cnsts, Eigen::MatrixXd& A, Eigen::MatrixXd& maxA,
406 Eigen::VectorXd& phi)
408 A.resize(number_cnsts, variable_set.size());
410 maxA.resize(number_cnsts, variable_set.size());
412 phi.resize(variable_set.size());
415 for (Variable& var : variable_set) {
416 if (var.sharing_penalty_ <= 0)
419 bool linked = false; // variable is linked to some constraint (specially for selective_update)
420 for (const Element& elem : var.cnsts_) {
421 const boost::intrusive::list_member_hook<>& cnst_hook = selective_update_active
422 ? elem.constraint->modified_constraint_set_hook_
423 : elem.constraint->active_constraint_set_hook_;
424 if (not cnst_hook.is_linked())
426 /* active and linked variable, lets check its consumption */
428 double consumption = elem.consumption_weight;
429 if (consumption > 0) {
430 int cnst_idx = cnst2idx_[elem.constraint];
431 A(cnst_idx, var_idx) += consumption;
432 // a variable with double penalty must receive half share, so it max weight is greater
433 maxA(cnst_idx, var_idx) = std::max(maxA(cnst_idx, var_idx), elem.max_consumption_weight * var.sharing_penalty_);
437 /* skip variables not linked to any modified or active constraint */
441 phi[var_idx] = var.get_bound();
442 idx2Var_[var_idx] = &var;
445 var.value_ = 1; // assign something by default for tasks with 0 consumption
448 // resize matrix to active variables only
449 A.conservativeResize(Eigen::NoChange_t::NoChange, var_idx);
450 maxA.conservativeResize(Eigen::NoChange_t::NoChange, var_idx);
451 phi.conservativeResize(var_idx);
454 template <class CnstList>
455 void BmfSystem::get_constraint_data(const CnstList& cnst_list, Eigen::VectorXd& C, std::vector<bool>& shared)
457 C.resize(cnst_list.size());
458 shared.resize(cnst_list.size());
461 for (const Constraint& cnst : cnst_list) {
462 C(cnst_idx) = cnst.bound_;
463 if (cnst.get_sharing_policy() == Constraint::SharingPolicy::NONLINEAR && cnst.dyn_constraint_cb_) {
464 C(cnst_idx) = cnst.dyn_constraint_cb_(cnst.bound_, cnst.concurrency_current_);
466 cnst2idx_[&cnst] = cnst_idx;
467 // FATPIPE links aren't really shared
468 shared[cnst_idx] = (cnst.sharing_policy_ != Constraint::SharingPolicy::FATPIPE);
473 void BmfSystem::do_solve()
475 if (selective_update_active)
476 bmf_solve(modified_constraint_set);
478 bmf_solve(active_constraint_set);
481 template <class CnstList> void BmfSystem::bmf_solve(const CnstList& cnst_list)
483 /* initialize players' weight and constraint matrices */
487 Eigen::MatrixXd maxA;
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);
494 auto solver = BmfSolver(std::move(A), std::move(maxA), std::move(C), std::move(shared), std::move(bounds));
495 auto rho = solver.solve();
501 for (int i = 0; i < rho.size(); i++) {
502 idx2Var_[i]->value_ = rho[i];
507 } // namespace kernel
508 } // namespace simgrid