lemon/capacity_scaling.h
author convert-repo
Thu, 04 Mar 2010 19:34:55 +0000
changeset 2659 611ced85018b
parent 2623 90defb96ee61
permissions -rw-r--r--
update tags
     1 /* -*- C++ -*-
     2  *
     3  * This file is a part of LEMON, a generic C++ optimization library
     4  *
     5  * Copyright (C) 2003-2008
     6  * Egervary Jeno Kombinatorikus Optimalizalasi Kutatocsoport
     7  * (Egervary Research Group on Combinatorial Optimization, EGRES).
     8  *
     9  * Permission to use, modify and distribute this software is granted
    10  * provided that this copyright notice appears in all copies. For
    11  * precise terms see the accompanying LICENSE file.
    12  *
    13  * This software is provided "AS IS" with no warranty of any kind,
    14  * express or implied, and with no claim as to its suitability for any
    15  * purpose.
    16  *
    17  */
    18 
    19 #ifndef LEMON_CAPACITY_SCALING_H
    20 #define LEMON_CAPACITY_SCALING_H
    21 
    22 /// \ingroup min_cost_flow
    23 ///
    24 /// \file
    25 /// \brief Capacity scaling algorithm for finding a minimum cost flow.
    26 
    27 #include <vector>
    28 #include <lemon/bin_heap.h>
    29 
    30 namespace lemon {
    31 
    32   /// \addtogroup min_cost_flow
    33   /// @{
    34 
    35   /// \brief Implementation of the capacity scaling algorithm for
    36   /// finding a minimum cost flow.
    37   ///
    38   /// \ref CapacityScaling implements the capacity scaling version
    39   /// of the successive shortest path algorithm for finding a minimum
    40   /// cost flow.
    41   ///
    42   /// \tparam Graph The directed graph type the algorithm runs on.
    43   /// \tparam LowerMap The type of the lower bound map.
    44   /// \tparam CapacityMap The type of the capacity (upper bound) map.
    45   /// \tparam CostMap The type of the cost (length) map.
    46   /// \tparam SupplyMap The type of the supply map.
    47   ///
    48   /// \warning
    49   /// - Edge capacities and costs should be \e non-negative \e integers.
    50   /// - Supply values should be \e signed \e integers.
    51   /// - The value types of the maps should be convertible to each other.
    52   /// - \c CostMap::Value must be signed type.
    53   ///
    54   /// \author Peter Kovacs
    55   template < typename Graph,
    56              typename LowerMap = typename Graph::template EdgeMap<int>,
    57              typename CapacityMap = typename Graph::template EdgeMap<int>,
    58              typename CostMap = typename Graph::template EdgeMap<int>,
    59              typename SupplyMap = typename Graph::template NodeMap<int> >
    60   class CapacityScaling
    61   {
    62     GRAPH_TYPEDEFS(typename Graph);
    63 
    64     typedef typename CapacityMap::Value Capacity;
    65     typedef typename CostMap::Value Cost;
    66     typedef typename SupplyMap::Value Supply;
    67     typedef typename Graph::template EdgeMap<Capacity> CapacityEdgeMap;
    68     typedef typename Graph::template NodeMap<Supply> SupplyNodeMap;
    69     typedef typename Graph::template NodeMap<Edge> PredMap;
    70 
    71   public:
    72 
    73     /// The type of the flow map.
    74     typedef typename Graph::template EdgeMap<Capacity> FlowMap;
    75     /// The type of the potential map.
    76     typedef typename Graph::template NodeMap<Cost> PotentialMap;
    77 
    78   private:
    79 
    80     /// \brief Special implementation of the \ref Dijkstra algorithm
    81     /// for finding shortest paths in the residual network.
    82     ///
    83     /// \ref ResidualDijkstra is a special implementation of the
    84     /// \ref Dijkstra algorithm for finding shortest paths in the
    85     /// residual network of the graph with respect to the reduced edge
    86     /// costs and modifying the node potentials according to the
    87     /// distance of the nodes.
    88     class ResidualDijkstra
    89     {
    90       typedef typename Graph::template NodeMap<int> HeapCrossRef;
    91       typedef BinHeap<Cost, HeapCrossRef> Heap;
    92 
    93     private:
    94 
    95       // The directed graph the algorithm runs on
    96       const Graph &_graph;
    97 
    98       // The main maps
    99       const FlowMap &_flow;
   100       const CapacityEdgeMap &_res_cap;
   101       const CostMap &_cost;
   102       const SupplyNodeMap &_excess;
   103       PotentialMap &_potential;
   104 
   105       // The distance map
   106       PotentialMap _dist;
   107       // The pred edge map
   108       PredMap &_pred;
   109       // The processed (i.e. permanently labeled) nodes
   110       std::vector<Node> _proc_nodes;
   111 
   112     public:
   113 
   114       /// Constructor.
   115       ResidualDijkstra( const Graph &graph,
   116                         const FlowMap &flow,
   117                         const CapacityEdgeMap &res_cap,
   118                         const CostMap &cost,
   119                         const SupplyMap &excess,
   120                         PotentialMap &potential,
   121                         PredMap &pred ) :
   122         _graph(graph), _flow(flow), _res_cap(res_cap), _cost(cost),
   123         _excess(excess), _potential(potential), _dist(graph),
   124         _pred(pred)
   125       {}
   126 
   127       /// Run the algorithm from the given source node.
   128       Node run(Node s, Capacity delta = 1) {
   129         HeapCrossRef heap_cross_ref(_graph, Heap::PRE_HEAP);
   130         Heap heap(heap_cross_ref);
   131         heap.push(s, 0);
   132         _pred[s] = INVALID;
   133         _proc_nodes.clear();
   134 
   135         // Processing nodes
   136         while (!heap.empty() && _excess[heap.top()] > -delta) {
   137           Node u = heap.top(), v;
   138           Cost d = heap.prio() + _potential[u], nd;
   139           _dist[u] = heap.prio();
   140           heap.pop();
   141           _proc_nodes.push_back(u);
   142 
   143           // Traversing outgoing edges
   144           for (OutEdgeIt e(_graph, u); e != INVALID; ++e) {
   145             if (_res_cap[e] >= delta) {
   146               v = _graph.target(e);
   147               switch(heap.state(v)) {
   148               case Heap::PRE_HEAP:
   149                 heap.push(v, d + _cost[e] - _potential[v]);
   150                 _pred[v] = e;
   151                 break;
   152               case Heap::IN_HEAP:
   153                 nd = d + _cost[e] - _potential[v];
   154                 if (nd < heap[v]) {
   155                   heap.decrease(v, nd);
   156                   _pred[v] = e;
   157                 }
   158                 break;
   159               case Heap::POST_HEAP:
   160                 break;
   161               }
   162             }
   163           }
   164 
   165           // Traversing incoming edges
   166           for (InEdgeIt e(_graph, u); e != INVALID; ++e) {
   167             if (_flow[e] >= delta) {
   168               v = _graph.source(e);
   169               switch(heap.state(v)) {
   170               case Heap::PRE_HEAP:
   171                 heap.push(v, d - _cost[e] - _potential[v]);
   172                 _pred[v] = e;
   173                 break;
   174               case Heap::IN_HEAP:
   175                 nd = d - _cost[e] - _potential[v];
   176                 if (nd < heap[v]) {
   177                   heap.decrease(v, nd);
   178                   _pred[v] = e;
   179                 }
   180                 break;
   181               case Heap::POST_HEAP:
   182                 break;
   183               }
   184             }
   185           }
   186         }
   187         if (heap.empty()) return INVALID;
   188 
   189         // Updating potentials of processed nodes
   190         Node t = heap.top();
   191         Cost t_dist = heap.prio();
   192         for (int i = 0; i < int(_proc_nodes.size()); ++i)
   193           _potential[_proc_nodes[i]] += _dist[_proc_nodes[i]] - t_dist;
   194 
   195         return t;
   196       }
   197 
   198     }; //class ResidualDijkstra
   199 
   200   private:
   201 
   202     // The directed graph the algorithm runs on
   203     const Graph &_graph;
   204     // The original lower bound map
   205     const LowerMap *_lower;
   206     // The modified capacity map
   207     CapacityEdgeMap _capacity;
   208     // The original cost map
   209     const CostMap &_cost;
   210     // The modified supply map
   211     SupplyNodeMap _supply;
   212     bool _valid_supply;
   213 
   214     // Edge map of the current flow
   215     FlowMap *_flow;
   216     bool _local_flow;
   217     // Node map of the current potentials
   218     PotentialMap *_potential;
   219     bool _local_potential;
   220 
   221     // The residual capacity map
   222     CapacityEdgeMap _res_cap;
   223     // The excess map
   224     SupplyNodeMap _excess;
   225     // The excess nodes (i.e. nodes with positive excess)
   226     std::vector<Node> _excess_nodes;
   227     // The deficit nodes (i.e. nodes with negative excess)
   228     std::vector<Node> _deficit_nodes;
   229 
   230     // The delta parameter used for capacity scaling
   231     Capacity _delta;
   232     // The maximum number of phases
   233     int _phase_num;
   234 
   235     // The pred edge map
   236     PredMap _pred;
   237     // Implementation of the Dijkstra algorithm for finding augmenting
   238     // shortest paths in the residual network
   239     ResidualDijkstra *_dijkstra;
   240 
   241   public:
   242 
   243     /// \brief General constructor (with lower bounds).
   244     ///
   245     /// General constructor (with lower bounds).
   246     ///
   247     /// \param graph The directed graph the algorithm runs on.
   248     /// \param lower The lower bounds of the edges.
   249     /// \param capacity The capacities (upper bounds) of the edges.
   250     /// \param cost The cost (length) values of the edges.
   251     /// \param supply The supply values of the nodes (signed).
   252     CapacityScaling( const Graph &graph,
   253                      const LowerMap &lower,
   254                      const CapacityMap &capacity,
   255                      const CostMap &cost,
   256                      const SupplyMap &supply ) :
   257       _graph(graph), _lower(&lower), _capacity(capacity), _cost(cost),
   258       _supply(supply), _flow(NULL), _local_flow(false),
   259       _potential(NULL), _local_potential(false),
   260       _res_cap(capacity), _excess(supply), _pred(graph), _dijkstra(NULL)
   261     {
   262       // Check the sum of supply values
   263       Supply sum = 0;
   264       for (NodeIt n(_graph); n != INVALID; ++n) sum += _supply[n];
   265       _valid_supply = sum == 0;
   266 
   267       // Remove non-zero lower bounds
   268       typename LowerMap::Value lcap;
   269       for (EdgeIt e(_graph); e != INVALID; ++e) {
   270         if ((lcap = lower[e]) != 0) {
   271           _capacity[e] -= lcap;
   272           _res_cap[e] -= lcap;
   273           _supply[_graph.source(e)] -= lcap;
   274           _supply[_graph.target(e)] += lcap;
   275           _excess[_graph.source(e)] -= lcap;
   276           _excess[_graph.target(e)] += lcap;
   277         }
   278       }
   279     }
   280 
   281     /// \brief General constructor (without lower bounds).
   282     ///
   283     /// General constructor (without lower bounds).
   284     ///
   285     /// \param graph The directed graph the algorithm runs on.
   286     /// \param capacity The capacities (upper bounds) of the edges.
   287     /// \param cost The cost (length) values of the edges.
   288     /// \param supply The supply values of the nodes (signed).
   289     CapacityScaling( const Graph &graph,
   290                      const CapacityMap &capacity,
   291                      const CostMap &cost,
   292                      const SupplyMap &supply ) :
   293       _graph(graph), _lower(NULL), _capacity(capacity), _cost(cost),
   294       _supply(supply), _flow(NULL), _local_flow(false),
   295       _potential(NULL), _local_potential(false),
   296       _res_cap(capacity), _excess(supply), _pred(graph), _dijkstra(NULL)
   297     {
   298       // Check the sum of supply values
   299       Supply sum = 0;
   300       for (NodeIt n(_graph); n != INVALID; ++n) sum += _supply[n];
   301       _valid_supply = sum == 0;
   302     }
   303 
   304     /// \brief Simple constructor (with lower bounds).
   305     ///
   306     /// Simple constructor (with lower bounds).
   307     ///
   308     /// \param graph The directed graph the algorithm runs on.
   309     /// \param lower The lower bounds of the edges.
   310     /// \param capacity The capacities (upper bounds) of the edges.
   311     /// \param cost The cost (length) values of the edges.
   312     /// \param s The source node.
   313     /// \param t The target node.
   314     /// \param flow_value The required amount of flow from node \c s
   315     /// to node \c t (i.e. the supply of \c s and the demand of \c t).
   316     CapacityScaling( const Graph &graph,
   317                      const LowerMap &lower,
   318                      const CapacityMap &capacity,
   319                      const CostMap &cost,
   320                      Node s, Node t,
   321                      Supply flow_value ) :
   322       _graph(graph), _lower(&lower), _capacity(capacity), _cost(cost),
   323       _supply(graph, 0), _flow(NULL), _local_flow(false),
   324       _potential(NULL), _local_potential(false),
   325       _res_cap(capacity), _excess(graph, 0), _pred(graph), _dijkstra(NULL)
   326     {
   327       // Remove non-zero lower bounds
   328       _supply[s] = _excess[s] =  flow_value;
   329       _supply[t] = _excess[t] = -flow_value;
   330       typename LowerMap::Value lcap;
   331       for (EdgeIt e(_graph); e != INVALID; ++e) {
   332         if ((lcap = lower[e]) != 0) {
   333           _capacity[e] -= lcap;
   334           _res_cap[e] -= lcap;
   335           _supply[_graph.source(e)] -= lcap;
   336           _supply[_graph.target(e)] += lcap;
   337           _excess[_graph.source(e)] -= lcap;
   338           _excess[_graph.target(e)] += lcap;
   339         }
   340       }
   341       _valid_supply = true;
   342     }
   343 
   344     /// \brief Simple constructor (without lower bounds).
   345     ///
   346     /// Simple constructor (without lower bounds).
   347     ///
   348     /// \param graph The directed graph the algorithm runs on.
   349     /// \param capacity The capacities (upper bounds) of the edges.
   350     /// \param cost The cost (length) values of the edges.
   351     /// \param s The source node.
   352     /// \param t The target node.
   353     /// \param flow_value The required amount of flow from node \c s
   354     /// to node \c t (i.e. the supply of \c s and the demand of \c t).
   355     CapacityScaling( const Graph &graph,
   356                      const CapacityMap &capacity,
   357                      const CostMap &cost,
   358                      Node s, Node t,
   359                      Supply flow_value ) :
   360       _graph(graph), _lower(NULL), _capacity(capacity), _cost(cost),
   361       _supply(graph, 0), _flow(NULL), _local_flow(false),
   362       _potential(NULL), _local_potential(false),
   363       _res_cap(capacity), _excess(graph, 0), _pred(graph), _dijkstra(NULL)
   364     {
   365       _supply[s] = _excess[s] =  flow_value;
   366       _supply[t] = _excess[t] = -flow_value;
   367       _valid_supply = true;
   368     }
   369 
   370     /// Destructor.
   371     ~CapacityScaling() {
   372       if (_local_flow) delete _flow;
   373       if (_local_potential) delete _potential;
   374       delete _dijkstra;
   375     }
   376 
   377     /// \brief Set the flow map.
   378     ///
   379     /// Set the flow map.
   380     ///
   381     /// \return \c (*this)
   382     CapacityScaling& flowMap(FlowMap &map) {
   383       if (_local_flow) {
   384         delete _flow;
   385         _local_flow = false;
   386       }
   387       _flow = &map;
   388       return *this;
   389     }
   390 
   391     /// \brief Set the potential map.
   392     ///
   393     /// Set the potential map.
   394     ///
   395     /// \return \c (*this)
   396     CapacityScaling& potentialMap(PotentialMap &map) {
   397       if (_local_potential) {
   398         delete _potential;
   399         _local_potential = false;
   400       }
   401       _potential = &map;
   402       return *this;
   403     }
   404 
   405     /// \name Execution control
   406 
   407     /// @{
   408 
   409     /// \brief Run the algorithm.
   410     ///
   411     /// This function runs the algorithm.
   412     ///
   413     /// \param scaling Enable or disable capacity scaling.
   414     /// If the maximum edge capacity and/or the amount of total supply
   415     /// is rather small, the algorithm could be slightly faster without
   416     /// scaling.
   417     ///
   418     /// \return \c true if a feasible flow can be found.
   419     bool run(bool scaling = true) {
   420       return init(scaling) && start();
   421     }
   422 
   423     /// @}
   424 
   425     /// \name Query Functions
   426     /// The results of the algorithm can be obtained using these
   427     /// functions.\n
   428     /// \ref lemon::CapacityScaling::run() "run()" must be called before
   429     /// using them.
   430 
   431     /// @{
   432 
   433     /// \brief Return a const reference to the edge map storing the
   434     /// found flow.
   435     ///
   436     /// Return a const reference to the edge map storing the found flow.
   437     ///
   438     /// \pre \ref run() must be called before using this function.
   439     const FlowMap& flowMap() const {
   440       return *_flow;
   441     }
   442 
   443     /// \brief Return a const reference to the node map storing the
   444     /// found potentials (the dual solution).
   445     ///
   446     /// Return a const reference to the node map storing the found
   447     /// potentials (the dual solution).
   448     ///
   449     /// \pre \ref run() must be called before using this function.
   450     const PotentialMap& potentialMap() const {
   451       return *_potential;
   452     }
   453 
   454     /// \brief Return the flow on the given edge.
   455     ///
   456     /// Return the flow on the given edge.
   457     ///
   458     /// \pre \ref run() must be called before using this function.
   459     Capacity flow(const Edge& edge) const {
   460       return (*_flow)[edge];
   461     }
   462 
   463     /// \brief Return the potential of the given node.
   464     ///
   465     /// Return the potential of the given node.
   466     ///
   467     /// \pre \ref run() must be called before using this function.
   468     Cost potential(const Node& node) const {
   469       return (*_potential)[node];
   470     }
   471 
   472     /// \brief Return the total cost of the found flow.
   473     ///
   474     /// Return the total cost of the found flow. The complexity of the
   475     /// function is \f$ O(e) \f$.
   476     ///
   477     /// \pre \ref run() must be called before using this function.
   478     Cost totalCost() const {
   479       Cost c = 0;
   480       for (EdgeIt e(_graph); e != INVALID; ++e)
   481         c += (*_flow)[e] * _cost[e];
   482       return c;
   483     }
   484 
   485     /// @}
   486 
   487   private:
   488 
   489     /// Initialize the algorithm.
   490     bool init(bool scaling) {
   491       if (!_valid_supply) return false;
   492 
   493       // Initializing maps
   494       if (!_flow) {
   495         _flow = new FlowMap(_graph);
   496         _local_flow = true;
   497       }
   498       if (!_potential) {
   499         _potential = new PotentialMap(_graph);
   500         _local_potential = true;
   501       }
   502       for (EdgeIt e(_graph); e != INVALID; ++e) (*_flow)[e] = 0;
   503       for (NodeIt n(_graph); n != INVALID; ++n) (*_potential)[n] = 0;
   504 
   505       _dijkstra = new ResidualDijkstra( _graph, *_flow, _res_cap, _cost,
   506                                         _excess, *_potential, _pred );
   507 
   508       // Initializing delta value
   509       if (scaling) {
   510         // With scaling
   511         Supply max_sup = 0, max_dem = 0;
   512         for (NodeIt n(_graph); n != INVALID; ++n) {
   513           if ( _supply[n] > max_sup) max_sup =  _supply[n];
   514           if (-_supply[n] > max_dem) max_dem = -_supply[n];
   515         }
   516         Capacity max_cap = 0;
   517         for (EdgeIt e(_graph); e != INVALID; ++e) {
   518           if (_capacity[e] > max_cap) max_cap = _capacity[e];
   519         }
   520         max_sup = std::min(std::min(max_sup, max_dem), max_cap);
   521         _phase_num = 0;
   522         for (_delta = 1; 2 * _delta <= max_sup; _delta *= 2)
   523           ++_phase_num;
   524       } else {
   525         // Without scaling
   526         _delta = 1;
   527       }
   528 
   529       return true;
   530     }
   531 
   532     bool start() {
   533       if (_delta > 1)
   534         return startWithScaling();
   535       else
   536         return startWithoutScaling();
   537     }
   538 
   539     /// Execute the capacity scaling algorithm.
   540     bool startWithScaling() {
   541       // Processing capacity scaling phases
   542       Node s, t;
   543       int phase_cnt = 0;
   544       int factor = 4;
   545       while (true) {
   546         // Saturating all edges not satisfying the optimality condition
   547         for (EdgeIt e(_graph); e != INVALID; ++e) {
   548           Node u = _graph.source(e), v = _graph.target(e);
   549           Cost c = _cost[e] + (*_potential)[u] - (*_potential)[v];
   550           if (c < 0 && _res_cap[e] >= _delta) {
   551             _excess[u] -= _res_cap[e];
   552             _excess[v] += _res_cap[e];
   553             (*_flow)[e] = _capacity[e];
   554             _res_cap[e] = 0;
   555           }
   556           else if (c > 0 && (*_flow)[e] >= _delta) {
   557             _excess[u] += (*_flow)[e];
   558             _excess[v] -= (*_flow)[e];
   559             (*_flow)[e] = 0;
   560             _res_cap[e] = _capacity[e];
   561           }
   562         }
   563 
   564         // Finding excess nodes and deficit nodes
   565         _excess_nodes.clear();
   566         _deficit_nodes.clear();
   567         for (NodeIt n(_graph); n != INVALID; ++n) {
   568           if (_excess[n] >=  _delta) _excess_nodes.push_back(n);
   569           if (_excess[n] <= -_delta) _deficit_nodes.push_back(n);
   570         }
   571         int next_node = 0, next_def_node = 0;
   572 
   573         // Finding augmenting shortest paths
   574         while (next_node < int(_excess_nodes.size())) {
   575           // Checking deficit nodes
   576           if (_delta > 1) {
   577             bool delta_deficit = false;
   578             for ( ; next_def_node < int(_deficit_nodes.size());
   579                     ++next_def_node ) {
   580               if (_excess[_deficit_nodes[next_def_node]] <= -_delta) {
   581                 delta_deficit = true;
   582                 break;
   583               }
   584             }
   585             if (!delta_deficit) break;
   586           }
   587 
   588           // Running Dijkstra
   589           s = _excess_nodes[next_node];
   590           if ((t = _dijkstra->run(s, _delta)) == INVALID) {
   591             if (_delta > 1) {
   592               ++next_node;
   593               continue;
   594             }
   595             return false;
   596           }
   597 
   598           // Augmenting along a shortest path from s to t.
   599           Capacity d = std::min(_excess[s], -_excess[t]);
   600           Node u = t;
   601           Edge e;
   602           if (d > _delta) {
   603             while ((e = _pred[u]) != INVALID) {
   604               Capacity rc;
   605               if (u == _graph.target(e)) {
   606                 rc = _res_cap[e];
   607                 u = _graph.source(e);
   608               } else {
   609                 rc = (*_flow)[e];
   610                 u = _graph.target(e);
   611               }
   612               if (rc < d) d = rc;
   613             }
   614           }
   615           u = t;
   616           while ((e = _pred[u]) != INVALID) {
   617             if (u == _graph.target(e)) {
   618               (*_flow)[e] += d;
   619               _res_cap[e] -= d;
   620               u = _graph.source(e);
   621             } else {
   622               (*_flow)[e] -= d;
   623               _res_cap[e] += d;
   624               u = _graph.target(e);
   625             }
   626           }
   627           _excess[s] -= d;
   628           _excess[t] += d;
   629 
   630           if (_excess[s] < _delta) ++next_node;
   631         }
   632 
   633         if (_delta == 1) break;
   634         if (++phase_cnt > _phase_num / 4) factor = 2;
   635         _delta = _delta <= factor ? 1 : _delta / factor;
   636       }
   637 
   638       // Handling non-zero lower bounds
   639       if (_lower) {
   640         for (EdgeIt e(_graph); e != INVALID; ++e)
   641           (*_flow)[e] += (*_lower)[e];
   642       }
   643       return true;
   644     }
   645 
   646     /// Execute the successive shortest path algorithm.
   647     bool startWithoutScaling() {
   648       // Finding excess nodes
   649       for (NodeIt n(_graph); n != INVALID; ++n)
   650         if (_excess[n] > 0) _excess_nodes.push_back(n);
   651       if (_excess_nodes.size() == 0) return true;
   652       int next_node = 0;
   653 
   654       // Finding shortest paths
   655       Node s, t;
   656       while ( _excess[_excess_nodes[next_node]] > 0 ||
   657               ++next_node < int(_excess_nodes.size()) )
   658       {
   659         // Running Dijkstra
   660         s = _excess_nodes[next_node];
   661         if ((t = _dijkstra->run(s)) == INVALID) return false;
   662 
   663         // Augmenting along a shortest path from s to t
   664         Capacity d = std::min(_excess[s], -_excess[t]);
   665         Node u = t;
   666         Edge e;
   667         if (d > 1) {
   668           while ((e = _pred[u]) != INVALID) {
   669             Capacity rc;
   670             if (u == _graph.target(e)) {
   671               rc = _res_cap[e];
   672               u = _graph.source(e);
   673             } else {
   674               rc = (*_flow)[e];
   675               u = _graph.target(e);
   676             }
   677             if (rc < d) d = rc;
   678           }
   679         }
   680         u = t;
   681         while ((e = _pred[u]) != INVALID) {
   682           if (u == _graph.target(e)) {
   683             (*_flow)[e] += d;
   684             _res_cap[e] -= d;
   685             u = _graph.source(e);
   686           } else {
   687             (*_flow)[e] -= d;
   688             _res_cap[e] += d;
   689             u = _graph.target(e);
   690           }
   691         }
   692         _excess[s] -= d;
   693         _excess[t] += d;
   694       }
   695 
   696       // Handling non-zero lower bounds
   697       if (_lower) {
   698         for (EdgeIt e(_graph); e != INVALID; ++e)
   699           (*_flow)[e] += (*_lower)[e];
   700       }
   701       return true;
   702     }
   703 
   704   }; //class CapacityScaling
   705 
   706   ///@}
   707 
   708 } //namespace lemon
   709 
   710 #endif //LEMON_CAPACITY_SCALING_H