lemon/capacity_scaling.h
author kpeter
Sun, 05 Oct 2008 13:37:17 +0000
changeset 2620 8f41a3129746
parent 2589 1bbb28acb8c9
child 2623 90defb96ee61
permissions -rw-r--r--
Doc improvements
     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(graph), _cost(cost),
   258       _supply(graph), _flow(0), _local_flow(false),
   259       _potential(0), _local_potential(false),
   260       _res_cap(graph), _excess(graph), _pred(graph)
   261     {
   262       // Removing non-zero lower bounds
   263       _capacity = subMap(capacity, lower);
   264       _res_cap = _capacity;
   265       Supply sum = 0;
   266       for (NodeIt n(_graph); n != INVALID; ++n) {
   267         Supply s = supply[n];
   268         for (InEdgeIt e(_graph, n); e != INVALID; ++e)
   269           s += lower[e];
   270         for (OutEdgeIt e(_graph, n); e != INVALID; ++e)
   271           s -= lower[e];
   272         _supply[n] = s;
   273         sum += s;
   274       }
   275       _valid_supply = sum == 0;
   276     }
   277 
   278     /// \brief General constructor (without lower bounds).
   279     ///
   280     /// General constructor (without lower bounds).
   281     ///
   282     /// \param graph The directed graph the algorithm runs on.
   283     /// \param capacity The capacities (upper bounds) of the edges.
   284     /// \param cost The cost (length) values of the edges.
   285     /// \param supply The supply values of the nodes (signed).
   286     CapacityScaling( const Graph &graph,
   287                      const CapacityMap &capacity,
   288                      const CostMap &cost,
   289                      const SupplyMap &supply ) :
   290       _graph(graph), _lower(NULL), _capacity(capacity), _cost(cost),
   291       _supply(supply), _flow(0), _local_flow(false),
   292       _potential(0), _local_potential(false),
   293       _res_cap(capacity), _excess(graph), _pred(graph)
   294     {
   295       // Checking the sum of supply values
   296       Supply sum = 0;
   297       for (NodeIt n(_graph); n != INVALID; ++n) sum += _supply[n];
   298       _valid_supply = sum == 0;
   299     }
   300 
   301     /// \brief Simple constructor (with lower bounds).
   302     ///
   303     /// Simple constructor (with lower bounds).
   304     ///
   305     /// \param graph The directed graph the algorithm runs on.
   306     /// \param lower The lower bounds of the edges.
   307     /// \param capacity The capacities (upper bounds) of the edges.
   308     /// \param cost The cost (length) values of the edges.
   309     /// \param s The source node.
   310     /// \param t The target node.
   311     /// \param flow_value The required amount of flow from node \c s
   312     /// to node \c t (i.e. the supply of \c s and the demand of \c t).
   313     CapacityScaling( const Graph &graph,
   314                      const LowerMap &lower,
   315                      const CapacityMap &capacity,
   316                      const CostMap &cost,
   317                      Node s, Node t,
   318                      Supply flow_value ) :
   319       _graph(graph), _lower(&lower), _capacity(graph), _cost(cost),
   320       _supply(graph), _flow(0), _local_flow(false),
   321       _potential(0), _local_potential(false),
   322       _res_cap(graph), _excess(graph), _pred(graph)
   323     {
   324       // Removing non-zero lower bounds
   325       _capacity = subMap(capacity, lower);
   326       _res_cap = _capacity;
   327       for (NodeIt n(_graph); n != INVALID; ++n) {
   328         Supply sum = 0;
   329         if (n == s) sum =  flow_value;
   330         if (n == t) sum = -flow_value;
   331         for (InEdgeIt e(_graph, n); e != INVALID; ++e)
   332           sum += lower[e];
   333         for (OutEdgeIt e(_graph, n); e != INVALID; ++e)
   334           sum -= lower[e];
   335         _supply[n] = sum;
   336       }
   337       _valid_supply = true;
   338     }
   339 
   340     /// \brief Simple constructor (without lower bounds).
   341     ///
   342     /// Simple constructor (without lower bounds).
   343     ///
   344     /// \param graph The directed graph the algorithm runs on.
   345     /// \param capacity The capacities (upper bounds) of the edges.
   346     /// \param cost The cost (length) values of the edges.
   347     /// \param s The source node.
   348     /// \param t The target node.
   349     /// \param flow_value The required amount of flow from node \c s
   350     /// to node \c t (i.e. the supply of \c s and the demand of \c t).
   351     CapacityScaling( const Graph &graph,
   352                      const CapacityMap &capacity,
   353                      const CostMap &cost,
   354                      Node s, Node t,
   355                      Supply flow_value ) :
   356       _graph(graph), _lower(NULL), _capacity(capacity), _cost(cost),
   357       _supply(graph, 0), _flow(0), _local_flow(false),
   358       _potential(0), _local_potential(false),
   359       _res_cap(capacity), _excess(graph), _pred(graph)
   360     {
   361       _supply[s] =  flow_value;
   362       _supply[t] = -flow_value;
   363       _valid_supply = true;
   364     }
   365 
   366     /// Destructor.
   367     ~CapacityScaling() {
   368       if (_local_flow) delete _flow;
   369       if (_local_potential) delete _potential;
   370       delete _dijkstra;
   371     }
   372 
   373     /// \brief Set the flow map.
   374     ///
   375     /// Set the flow map.
   376     ///
   377     /// \return \c (*this)
   378     CapacityScaling& flowMap(FlowMap &map) {
   379       if (_local_flow) {
   380         delete _flow;
   381         _local_flow = false;
   382       }
   383       _flow = &map;
   384       return *this;
   385     }
   386 
   387     /// \brief Set the potential map.
   388     ///
   389     /// Set the potential map.
   390     ///
   391     /// \return \c (*this)
   392     CapacityScaling& potentialMap(PotentialMap &map) {
   393       if (_local_potential) {
   394         delete _potential;
   395         _local_potential = false;
   396       }
   397       _potential = &map;
   398       return *this;
   399     }
   400 
   401     /// \name Execution control
   402 
   403     /// @{
   404 
   405     /// \brief Run the algorithm.
   406     ///
   407     /// This function runs the algorithm.
   408     ///
   409     /// \param scaling Enable or disable capacity scaling.
   410     /// If the maximum edge capacity and/or the amount of total supply
   411     /// is rather small, the algorithm could be slightly faster without
   412     /// scaling.
   413     ///
   414     /// \return \c true if a feasible flow can be found.
   415     bool run(bool scaling = true) {
   416       return init(scaling) && start();
   417     }
   418 
   419     /// @}
   420 
   421     /// \name Query Functions
   422     /// The results of the algorithm can be obtained using these
   423     /// functions.\n
   424     /// \ref lemon::CapacityScaling::run() "run()" must be called before
   425     /// using them.
   426 
   427     /// @{
   428 
   429     /// \brief Return a const reference to the edge map storing the
   430     /// found flow.
   431     ///
   432     /// Return a const reference to the edge map storing the found flow.
   433     ///
   434     /// \pre \ref run() must be called before using this function.
   435     const FlowMap& flowMap() const {
   436       return *_flow;
   437     }
   438 
   439     /// \brief Return a const reference to the node map storing the
   440     /// found potentials (the dual solution).
   441     ///
   442     /// Return a const reference to the node map storing the found
   443     /// potentials (the dual solution).
   444     ///
   445     /// \pre \ref run() must be called before using this function.
   446     const PotentialMap& potentialMap() const {
   447       return *_potential;
   448     }
   449 
   450     /// \brief Return the flow on the given edge.
   451     ///
   452     /// Return the flow on the given edge.
   453     ///
   454     /// \pre \ref run() must be called before using this function.
   455     Capacity flow(const Edge& edge) const {
   456       return (*_flow)[edge];
   457     }
   458 
   459     /// \brief Return the potential of the given node.
   460     ///
   461     /// Return the potential of the given node.
   462     ///
   463     /// \pre \ref run() must be called before using this function.
   464     Cost potential(const Node& node) const {
   465       return (*_potential)[node];
   466     }
   467 
   468     /// \brief Return the total cost of the found flow.
   469     ///
   470     /// Return the total cost of the found flow. The complexity of the
   471     /// function is \f$ O(e) \f$.
   472     ///
   473     /// \pre \ref run() must be called before using this function.
   474     Cost totalCost() const {
   475       Cost c = 0;
   476       for (EdgeIt e(_graph); e != INVALID; ++e)
   477         c += (*_flow)[e] * _cost[e];
   478       return c;
   479     }
   480 
   481     /// @}
   482 
   483   private:
   484 
   485     /// Initialize the algorithm.
   486     bool init(bool scaling) {
   487       if (!_valid_supply) return false;
   488 
   489       // Initializing maps
   490       if (!_flow) {
   491         _flow = new FlowMap(_graph);
   492         _local_flow = true;
   493       }
   494       if (!_potential) {
   495         _potential = new PotentialMap(_graph);
   496         _local_potential = true;
   497       }
   498       for (EdgeIt e(_graph); e != INVALID; ++e) (*_flow)[e] = 0;
   499       for (NodeIt n(_graph); n != INVALID; ++n) (*_potential)[n] = 0;
   500       _excess = _supply;
   501 
   502       _dijkstra = new ResidualDijkstra( _graph, *_flow, _res_cap, _cost,
   503                                         _excess, *_potential, _pred );
   504 
   505       // Initializing delta value
   506       if (scaling) {
   507         // With scaling
   508         Supply max_sup = 0, max_dem = 0;
   509         for (NodeIt n(_graph); n != INVALID; ++n) {
   510           if ( _supply[n] > max_sup) max_sup =  _supply[n];
   511           if (-_supply[n] > max_dem) max_dem = -_supply[n];
   512         }
   513         Capacity max_cap = 0;
   514         for (EdgeIt e(_graph); e != INVALID; ++e) {
   515           if (_capacity[e] > max_cap) max_cap = _capacity[e];
   516         }
   517         max_sup = std::min(std::min(max_sup, max_dem), max_cap);
   518         _phase_num = 0;
   519         for (_delta = 1; 2 * _delta <= max_sup; _delta *= 2)
   520           ++_phase_num;
   521       } else {
   522         // Without scaling
   523         _delta = 1;
   524       }
   525 
   526       return true;
   527     }
   528 
   529     bool start() {
   530       if (_delta > 1)
   531         return startWithScaling();
   532       else
   533         return startWithoutScaling();
   534     }
   535 
   536     /// Execute the capacity scaling algorithm.
   537     bool startWithScaling() {
   538       // Processing capacity scaling phases
   539       Node s, t;
   540       int phase_cnt = 0;
   541       int factor = 4;
   542       while (true) {
   543         // Saturating all edges not satisfying the optimality condition
   544         for (EdgeIt e(_graph); e != INVALID; ++e) {
   545           Node u = _graph.source(e), v = _graph.target(e);
   546           Cost c = _cost[e] + (*_potential)[u] - (*_potential)[v];
   547           if (c < 0 && _res_cap[e] >= _delta) {
   548             _excess[u] -= _res_cap[e];
   549             _excess[v] += _res_cap[e];
   550             (*_flow)[e] = _capacity[e];
   551             _res_cap[e] = 0;
   552           }
   553           else if (c > 0 && (*_flow)[e] >= _delta) {
   554             _excess[u] += (*_flow)[e];
   555             _excess[v] -= (*_flow)[e];
   556             (*_flow)[e] = 0;
   557             _res_cap[e] = _capacity[e];
   558           }
   559         }
   560 
   561         // Finding excess nodes and deficit nodes
   562         _excess_nodes.clear();
   563         _deficit_nodes.clear();
   564         for (NodeIt n(_graph); n != INVALID; ++n) {
   565           if (_excess[n] >=  _delta) _excess_nodes.push_back(n);
   566           if (_excess[n] <= -_delta) _deficit_nodes.push_back(n);
   567         }
   568         int next_node = 0, next_def_node = 0;
   569 
   570         // Finding augmenting shortest paths
   571         while (next_node < int(_excess_nodes.size())) {
   572           // Checking deficit nodes
   573           if (_delta > 1) {
   574             bool delta_deficit = false;
   575             for ( ; next_def_node < int(_deficit_nodes.size());
   576                     ++next_def_node ) {
   577               if (_excess[_deficit_nodes[next_def_node]] <= -_delta) {
   578                 delta_deficit = true;
   579                 break;
   580               }
   581             }
   582             if (!delta_deficit) break;
   583           }
   584 
   585           // Running Dijkstra
   586           s = _excess_nodes[next_node];
   587           if ((t = _dijkstra->run(s, _delta)) == INVALID) {
   588             if (_delta > 1) {
   589               ++next_node;
   590               continue;
   591             }
   592             return false;
   593           }
   594 
   595           // Augmenting along a shortest path from s to t.
   596           Capacity d = std::min(_excess[s], -_excess[t]);
   597           Node u = t;
   598           Edge e;
   599           if (d > _delta) {
   600             while ((e = _pred[u]) != INVALID) {
   601               Capacity rc;
   602               if (u == _graph.target(e)) {
   603                 rc = _res_cap[e];
   604                 u = _graph.source(e);
   605               } else {
   606                 rc = (*_flow)[e];
   607                 u = _graph.target(e);
   608               }
   609               if (rc < d) d = rc;
   610             }
   611           }
   612           u = t;
   613           while ((e = _pred[u]) != INVALID) {
   614             if (u == _graph.target(e)) {
   615               (*_flow)[e] += d;
   616               _res_cap[e] -= d;
   617               u = _graph.source(e);
   618             } else {
   619               (*_flow)[e] -= d;
   620               _res_cap[e] += d;
   621               u = _graph.target(e);
   622             }
   623           }
   624           _excess[s] -= d;
   625           _excess[t] += d;
   626 
   627           if (_excess[s] < _delta) ++next_node;
   628         }
   629 
   630         if (_delta == 1) break;
   631         if (++phase_cnt > _phase_num / 4) factor = 2;
   632         _delta = _delta <= factor ? 1 : _delta / factor;
   633       }
   634 
   635       // Handling non-zero lower bounds
   636       if (_lower) {
   637         for (EdgeIt e(_graph); e != INVALID; ++e)
   638           (*_flow)[e] += (*_lower)[e];
   639       }
   640       return true;
   641     }
   642 
   643     /// Execute the successive shortest path algorithm.
   644     bool startWithoutScaling() {
   645       // Finding excess nodes
   646       for (NodeIt n(_graph); n != INVALID; ++n)
   647         if (_excess[n] > 0) _excess_nodes.push_back(n);
   648       if (_excess_nodes.size() == 0) return true;
   649       int next_node = 0;
   650 
   651       // Finding shortest paths
   652       Node s, t;
   653       while ( _excess[_excess_nodes[next_node]] > 0 ||
   654               ++next_node < int(_excess_nodes.size()) )
   655       {
   656         // Running Dijkstra
   657         s = _excess_nodes[next_node];
   658         if ((t = _dijkstra->run(s)) == INVALID) return false;
   659 
   660         // Augmenting along a shortest path from s to t
   661         Capacity d = std::min(_excess[s], -_excess[t]);
   662         Node u = t;
   663         Edge e;
   664         if (d > 1) {
   665           while ((e = _pred[u]) != INVALID) {
   666             Capacity rc;
   667             if (u == _graph.target(e)) {
   668               rc = _res_cap[e];
   669               u = _graph.source(e);
   670             } else {
   671               rc = (*_flow)[e];
   672               u = _graph.target(e);
   673             }
   674             if (rc < d) d = rc;
   675           }
   676         }
   677         u = t;
   678         while ((e = _pred[u]) != INVALID) {
   679           if (u == _graph.target(e)) {
   680             (*_flow)[e] += d;
   681             _res_cap[e] -= d;
   682             u = _graph.source(e);
   683           } else {
   684             (*_flow)[e] -= d;
   685             _res_cap[e] += d;
   686             u = _graph.target(e);
   687           }
   688         }
   689         _excess[s] -= d;
   690         _excess[t] += d;
   691       }
   692 
   693       // Handling non-zero lower bounds
   694       if (_lower) {
   695         for (EdgeIt e(_graph); e != INVALID; ++e)
   696           (*_flow)[e] += (*_lower)[e];
   697       }
   698       return true;
   699     }
   700 
   701   }; //class CapacityScaling
   702 
   703   ///@}
   704 
   705 } //namespace lemon
   706 
   707 #endif //LEMON_CAPACITY_SCALING_H