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"
12 XBT_LOG_NEW_DEFAULT_SUBCATEGORY(ker_bmf, kernel, "Kernel BMF solver");
14 int sg_bmf_max_iterations = 1000; /* Change this with --cfg=bmf/max-iterations:VALUE */
20 AllocationGenerator::AllocationGenerator(Eigen::MatrixXd A) : A_(std::move(A)), alloc_(A_.cols(), 0)
22 // got a first valid allocation
23 for (size_t p = 0; p < alloc_.size(); p++) {
24 for (int r = 0; r < A_.rows(); r++) {
33 bool AllocationGenerator::next(std::vector<int>& next_alloc)
41 auto n_resources = A_.rows();
43 while (idx < alloc_.size()) {
44 alloc_[idx] = (alloc_[idx] + 1) % n_resources;
45 if (alloc_[idx] == 0) {
51 if (A_(alloc_[idx], idx) > 0) {
59 /*****************************************************************************/
61 BmfSolver::BmfSolver(Eigen::MatrixXd A, Eigen::MatrixXd maxA, Eigen::VectorXd C, std::vector<bool> shared,
64 , maxA_(std::move(maxA))
66 , C_shared_(std::move(shared))
67 , phi_(std::move(phi))
70 xbt_assert(max_iteration_ > 0,
71 "Invalid number of iterations for BMF solver. Please check your \"bmf/max-iterations\" configuration.");
72 xbt_assert(A_.cols() == maxA_.cols(), "Invalid number of cols in matrix A (%td) or maxA (%td)", A_.cols(),
74 xbt_assert(A_.cols() == phi_.size(), "Invalid size of phi vector (%td)", phi_.size());
75 xbt_assert(static_cast<long>(C_shared_.size()) == C_.size(), "Invalid size param shared (%zu)", C_shared_.size());
78 template <typename T> std::string BmfSolver::debug_eigen(const T& obj) const
80 std::stringstream debug;
85 template <typename C> std::string BmfSolver::debug_vector(const C& container) const
87 std::stringstream debug;
88 std::copy(container.begin(), container.end(),
89 std::ostream_iterator<typename std::remove_reference<decltype(container)>::type::value_type>(debug, " "));
93 std::string BmfSolver::debug_alloc(const allocation_map_t& alloc) const
95 std::stringstream debug;
96 for (const auto& e : alloc) {
97 debug << "{" + std::to_string(e.first) + ": [" + debug_vector(e.second) + "]}, ";
102 double BmfSolver::get_resource_capacity(int resource, const std::vector<int>& bounded_players) const
104 double capacity = C_[resource];
105 if (not C_shared_[resource])
108 for (int p : bounded_players) {
109 capacity -= A_(resource, p) * phi_[p];
111 return std::max(0.0, capacity);
114 std::vector<int> BmfSolver::alloc_map_to_vector(const allocation_map_t& alloc) const
116 std::vector<int> alloc_by_player(A_.cols(), -1);
117 for (const auto& it : alloc) {
118 for (auto p : it.second) {
119 alloc_by_player[p] = it.first;
122 return alloc_by_player;
125 std::vector<int> BmfSolver::get_bounded_players(const allocation_map_t& alloc) const
127 std::vector<int> bounded_players;
128 for (const auto& e : alloc) {
129 if (e.first == NO_RESOURCE) {
130 bounded_players.insert(bounded_players.end(), e.second.begin(), e.second.end());
133 return bounded_players;
136 Eigen::VectorXd BmfSolver::equilibrium(const allocation_map_t& alloc) const
138 auto n_players = A_.cols();
139 Eigen::MatrixXd A_p = Eigen::MatrixXd::Zero(n_players, n_players); // square matrix with number of players
140 Eigen::VectorXd C_p = Eigen::VectorXd::Zero(n_players);
143 auto bounded_players = get_bounded_players(alloc);
144 for (const auto& e : alloc) {
145 // add one row for the resource with A[r,]
146 int cur_resource = e.first;
147 /* bounded players, nothing to do */
148 if (cur_resource == NO_RESOURCE)
150 /* not shared resource, each player can receive the full capacity of the resource */
151 if (not C_shared_[cur_resource]) {
152 for (int i : e.second) {
153 C_p[row] = get_resource_capacity(cur_resource, bounded_players);
154 A_p(row, i) = A_(cur_resource, i);
160 /* shared resource: fairly share it between players */
161 A_p.row(row) = A_.row(cur_resource);
162 C_p[row] = get_resource_capacity(cur_resource, bounded_players);
164 if (e.second.size() > 1) {
165 // if 2 players have chosen the same resource
166 // they must have a fair sharing of this resource, adjust A_p and C_p accordingly
167 auto it = e.second.begin();
168 int i = *it; // first player
169 /* for each other player sharing this resource */
170 for (++it; it != e.second.end(); ++it) {
171 /* player i and k on this resource j: so maxA_ji*rho_i - maxA_jk*rho_k = 0 */
174 A_p(row, i) = maxA_(cur_resource, i);
175 A_p(row, k) = -maxA_(cur_resource, k);
180 /* clear players which are externally bounded */
181 for (int p : bounded_players) {
182 A_p.col(p).setZero();
185 XBT_DEBUG("A':\n%s", debug_eigen(A_p).c_str());
187 XBT_DEBUG("C':\n%s", debug_eigen(C_p).c_str());
188 /* PartialPivLU is much faster than FullPivLU but requires that the matrix is invertible
189 * FullPivLU however assures that it finds come solution even if the matrix is singular
190 * Ideally we would like to be optimist and try Partial and in case of error, go back
192 * However, this with isNaN doesn't work if compiler uses -Ofastmath. In our case,
193 * the icc compiler raises an error when compiling the code (comparison with NaN always evaluates to false in fast
194 * floating point modes).
195 * Eigen::VectorXd rho = Eigen::PartialPivLU<Eigen::MatrixXd>(A_p).solve(C_p);
196 * if (rho.array().isNaN().any()) {
197 * XBT_DEBUG("rho with nan values, falling back to FullPivLU, rho:\n%s", debug_eigen(rho).c_str());
198 * rho = Eigen::FullPivLU<Eigen::MatrixXd>(A_p).solve(C_p);
202 Eigen::VectorXd rho = Eigen::FullPivLU<Eigen::MatrixXd>(A_p).solve(C_p);
203 for (int p : bounded_players) {
209 bool BmfSolver::disturb_allocation(allocation_map_t& alloc, std::vector<int>& alloc_by_player)
211 while (gen_.next(alloc_by_player)) {
212 if (allocations_.find(alloc_by_player) == allocations_.end()) {
213 allocations_.clear();
214 allocations_.insert(alloc_by_player);
216 for (size_t p = 0; p < alloc_by_player.size(); p++) {
217 alloc[alloc_by_player[p]].insert(p);
225 bool BmfSolver::get_alloc(const Eigen::VectorXd& fair_sharing, const allocation_map_t& last_alloc,
226 allocation_map_t& alloc, bool initial)
229 for (int player_idx = 0; player_idx < A_.cols(); player_idx++) {
230 int selected_resource = NO_RESOURCE;
231 double bound = phi_[player_idx];
232 double min_share = (bound <= 0 || initial) ? -1 : bound;
233 for (int cnst_idx = 0; cnst_idx < A_.rows(); cnst_idx++) {
234 if (A_(cnst_idx, player_idx) <= 0.0)
237 double share = fair_sharing[cnst_idx] / A_(cnst_idx, player_idx);
238 if (min_share == -1 || share < min_share) {
239 selected_resource = cnst_idx;
243 alloc[selected_resource].insert(player_idx);
245 bool is_stable = (alloc == last_alloc);
249 std::vector<int> alloc_by_player = alloc_map_to_vector(alloc);
250 auto ret = allocations_.insert(alloc_by_player);
251 /* oops, allocation already tried, let's pertube it a bit */
252 if (not ret.second) {
253 XBT_DEBUG("Allocation already tried: %s", debug_alloc(alloc).c_str());
254 return disturb_allocation(alloc, alloc_by_player);
259 void BmfSolver::set_fair_sharing(const allocation_map_t& alloc, const Eigen::VectorXd& rho,
260 Eigen::VectorXd& fair_sharing) const
262 std::vector<int> bounded_players = get_bounded_players(alloc);
264 for (int r = 0; r < fair_sharing.size(); r++) {
265 auto it = alloc.find(r);
266 if (it != alloc.end()) { // resource selected by some player, fair share depends on rho
267 int player = *(it->second.begin()); // equilibrium assures that every player receives the same, use one of them to
268 // calculate the fair sharing for resource r
269 fair_sharing[r] = A_(r, player) * rho[player];
270 } else { // nobody selects this resource, fair_sharing depends on resource saturation
271 // resource r is saturated (A[r,*] * rho > C), divide it among players
272 double consumption_r = A_.row(r) * rho;
273 double_update(&consumption_r, C_[r], sg_maxmin_precision);
274 if (consumption_r > 0.0) {
275 auto n_players = (A_.row(r).array() > 0).count();
276 fair_sharing[r] = C_[r] / n_players;
278 fair_sharing[r] = get_resource_capacity(r, bounded_players);
284 bool BmfSolver::is_bmf(const Eigen::VectorXd& rho) const
288 // 1) the capacity of all resources is respected
289 Eigen::VectorXd shared(C_shared_.size());
290 for (int j = 0; j < shared.size(); j++)
291 shared[j] = C_shared_[j] ? 1.0 : 0.0;
293 Eigen::VectorXd remaining = (A_ * rho) - C_;
294 remaining = remaining.array() * shared.array(); // ignore non shared resources
295 bmf = bmf && (not std::any_of(remaining.data(), remaining.data() + remaining.size(),
296 [](double v) { return double_positive(v, sg_maxmin_precision); }));
298 // 3) every player receives maximum share in at least 1 saturated resource
299 // due to subflows, compare with the maximum consumption and not the A matrix
300 Eigen::MatrixXd usage =
301 maxA_.array().rowwise() * rho.transpose().array(); // usage_ji: indicates the usage of player i on resource j
303 XBT_DEBUG("Usage_ji considering max consumption:\n%s", debug_eigen(usage).c_str());
304 auto max_share = usage.rowwise().maxCoeff(); // max share for each resource j
306 // matrix_ji: boolean indicating player p has the maximum share at resource j
307 Eigen::MatrixXi player_max_share =
308 ((usage.array().colwise() - max_share.array()).abs() <= sg_maxmin_precision).cast<int>();
309 // but only saturated resources must be considered
310 Eigen::VectorXi saturated = (remaining.array().abs() <= sg_maxmin_precision).cast<int>();
311 XBT_DEBUG("Saturated_j resources:\n%s", debug_eigen(saturated).c_str());
312 player_max_share.array().colwise() *= saturated.array();
314 // just check if it has received at least it's bound
315 for (int p = 0; p < rho.size(); p++) {
316 if (double_equals(rho[p], phi_[p], sg_maxmin_precision)) {
317 player_max_share(0, p) = 1; // it doesn't really matter, just to say that it's a bmf
322 // 2) at least 1 resource is saturated
323 bmf = bmf && (saturated.array() == 1).any();
325 XBT_DEBUG("Player_ji usage of saturated resources:\n%s", debug_eigen(player_max_share).c_str());
326 // for all columns(players) it has to be the max at least in 1
327 bmf = bmf && (player_max_share.colwise().sum().array() >= 1).all();
331 Eigen::VectorXd BmfSolver::solve()
333 XBT_DEBUG("Starting BMF solver");
335 XBT_DEBUG("A:\n%s", debug_eigen(A_).c_str());
336 XBT_DEBUG("maxA:\n%s", debug_eigen(maxA_).c_str());
337 XBT_DEBUG("C:\n%s", debug_eigen(C_).c_str());
339 /* no flows to share, just returns */
344 auto fair_sharing = C_;
346 /* BMF allocation for each player (current and last one) stop when are equal */
347 allocation_map_t last_alloc;
348 allocation_map_t cur_alloc;
351 while (it < max_iteration_ && not get_alloc(fair_sharing, last_alloc, cur_alloc, it == 0)) {
352 last_alloc = cur_alloc;
353 XBT_DEBUG("BMF: iteration %d", it);
354 XBT_DEBUG("B (current allocation): %s", debug_alloc(cur_alloc).c_str());
356 // solve inv(A)*rho = C
357 rho = equilibrium(cur_alloc);
358 XBT_DEBUG("rho:\n%s", debug_eigen(rho).c_str());
360 // get fair sharing for each resource
361 set_fair_sharing(cur_alloc, rho, fair_sharing);
362 XBT_DEBUG("Fair sharing vector (per resource):\n%s", debug_eigen(fair_sharing).c_str());
364 // get new allocation for players
368 /* Not mandatory but a safe check to assure we have a proper solution */
369 if (not is_bmf(rho)) {
370 fprintf(stderr, "Unable to find a BMF allocation for your system.\n"
371 "You may try to increase the maximum number of iterations performed by BMF solver "
372 "(\"--cfg=bmf/max-iterations\").\n"
373 "Additionally, you could decrease numerical precision (\"--cfg=surf/precision\").\n");
374 fprintf(stderr, "Internal states (after %d iterations):\n", it);
375 fprintf(stderr, "A:\n%s\n", debug_eigen(A_).c_str());
376 fprintf(stderr, "maxA:\n%s\n", debug_eigen(maxA_).c_str());
377 fprintf(stderr, "C:\n%s\n", debug_eigen(C_).c_str());
378 fprintf(stderr, "C_shared:\n%s\n", debug_vector(C_shared_).c_str());
379 fprintf(stderr, "phi:\n%s\n", debug_eigen(phi_).c_str());
380 fprintf(stderr, "rho:\n%s\n", debug_eigen(rho).c_str());
384 XBT_DEBUG("BMF done after %d iterations", it);
388 /*****************************************************************************/
390 void BmfSystem::get_flows_data(Eigen::Index number_cnsts, Eigen::MatrixXd& A, Eigen::MatrixXd& maxA,
391 Eigen::VectorXd& phi)
393 A.resize(number_cnsts, variable_set.size());
395 maxA.resize(number_cnsts, variable_set.size());
397 phi.resize(variable_set.size());
400 for (Variable& var : variable_set) {
401 if (var.sharing_penalty_ <= 0)
404 bool linked = false; // variable is linked to some constraint (specially for selective_update)
405 for (const Element& elem : var.cnsts_) {
406 const boost::intrusive::list_member_hook<>& cnst_hook = selective_update_active
407 ? elem.constraint->modified_constraint_set_hook_
408 : elem.constraint->active_constraint_set_hook_;
409 if (not cnst_hook.is_linked())
411 /* active and linked variable, lets check its consumption */
413 double consumption = elem.consumption_weight;
414 if (consumption > 0) {
415 int cnst_idx = cnst2idx_[elem.constraint];
416 A(cnst_idx, var_idx) += consumption;
417 // a variable with double penalty must receive half share, so it max weight is greater
418 maxA(cnst_idx, var_idx) = std::max(maxA(cnst_idx, var_idx), elem.max_consumption_weight * var.sharing_penalty_);
422 /* skip variables not linked to any modified or active constraint */
426 phi[var_idx] = var.get_bound();
427 idx2Var_[var_idx] = &var;
430 var.value_ = 1; // assign something by default for tasks with 0 consumption
433 // resize matrix to active variables only
434 A.conservativeResize(Eigen::NoChange_t::NoChange, var_idx);
435 maxA.conservativeResize(Eigen::NoChange_t::NoChange, var_idx);
436 phi.conservativeResize(var_idx);
439 template <class CnstList>
440 void BmfSystem::get_constraint_data(const CnstList& cnst_list, Eigen::VectorXd& C, std::vector<bool>& shared)
442 C.resize(cnst_list.size());
443 shared.resize(cnst_list.size());
446 for (const Constraint& cnst : cnst_list) {
447 C(cnst_idx) = cnst.bound_;
448 if (cnst.get_sharing_policy() == Constraint::SharingPolicy::NONLINEAR && cnst.dyn_constraint_cb_) {
449 C(cnst_idx) = cnst.dyn_constraint_cb_(cnst.bound_, cnst.concurrency_current_);
450 if (not warned_nonlinear_) {
451 XBT_WARN("You are using dynamic constraint bound with parallel tasks and BMF model."
452 " The BMF solver assumes that all flows (and subflows) are always active and executing."
453 " This is quite pessimist, specially considering parallel tasks with small subflows."
454 " Analyze your results with caution.");
455 warned_nonlinear_ = true;
458 cnst2idx_[&cnst] = cnst_idx;
459 // FATPIPE links aren't really shared
460 shared[cnst_idx] = (cnst.sharing_policy_ != Constraint::SharingPolicy::FATPIPE);
465 void BmfSystem::solve()
468 if (selective_update_active)
469 bmf_solve(modified_constraint_set);
471 bmf_solve(active_constraint_set);
475 template <class CnstList> void BmfSystem::bmf_solve(const CnstList& cnst_list)
477 /* initialize players' weight and constraint matrices */
481 Eigen::MatrixXd maxA;
483 Eigen::VectorXd bounds;
484 std::vector<bool> shared;
485 get_constraint_data(cnst_list, C, shared);
486 get_flows_data(C.size(), A, maxA, bounds);
488 auto solver = BmfSolver(std::move(A), std::move(maxA), std::move(C), std::move(shared), std::move(bounds));
489 auto rho = solver.solve();
495 for (int i = 0; i < rho.size(); i++) {
496 idx2Var_[i]->value_ = rho[i];
503 } // namespace kernel
504 } // namespace simgrid