lemon/cancel_and_tighten.h
author convert-repo
Thu, 04 Mar 2010 19:34:55 +0000
changeset 2659 611ced85018b
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_CANCEL_AND_TIGHTEN_H
    20 #define LEMON_CANCEL_AND_TIGHTEN_H
    21 
    22 /// \ingroup min_cost_flow
    23 ///
    24 /// \file
    25 /// \brief Cancel and Tighten algorithm for finding a minimum cost flow.
    26 
    27 #include <vector>
    28 
    29 #include <lemon/circulation.h>
    30 #include <lemon/bellman_ford.h>
    31 #include <lemon/min_mean_cycle.h>
    32 #include <lemon/graph_adaptor.h>
    33 #include <lemon/tolerance.h>
    34 #include <lemon/math.h>
    35 
    36 #include <lemon/static_graph.h>
    37 
    38 namespace lemon {
    39 
    40   /// \addtogroup min_cost_flow
    41   /// @{
    42 
    43   /// \brief Implementation of the Cancel and Tighten algorithm for
    44   /// finding a minimum cost flow.
    45   ///
    46   /// \ref CancelAndTighten implements the Cancel and Tighten algorithm for
    47   /// finding a minimum cost flow.
    48   ///
    49   /// \tparam Graph The directed graph type the algorithm runs on.
    50   /// \tparam LowerMap The type of the lower bound map.
    51   /// \tparam CapacityMap The type of the capacity (upper bound) map.
    52   /// \tparam CostMap The type of the cost (length) map.
    53   /// \tparam SupplyMap The type of the supply map.
    54   ///
    55   /// \warning
    56   /// - Edge capacities and costs should be \e non-negative \e integers.
    57   /// - Supply values should be \e signed \e integers.
    58   /// - The value types of the maps should be convertible to each other.
    59   /// - \c CostMap::Value must be signed type.
    60   ///
    61   /// \author Peter Kovacs
    62   template < typename Graph,
    63              typename LowerMap = typename Graph::template EdgeMap<int>,
    64              typename CapacityMap = typename Graph::template EdgeMap<int>,
    65              typename CostMap = typename Graph::template EdgeMap<int>,
    66              typename SupplyMap = typename Graph::template NodeMap<int> >
    67   class CancelAndTighten
    68   {
    69     GRAPH_TYPEDEFS(typename Graph);
    70 
    71     typedef typename CapacityMap::Value Capacity;
    72     typedef typename CostMap::Value Cost;
    73     typedef typename SupplyMap::Value Supply;
    74     typedef typename Graph::template EdgeMap<Capacity> CapacityEdgeMap;
    75     typedef typename Graph::template NodeMap<Supply> SupplyNodeMap;
    76 
    77     typedef ResGraphAdaptor< const Graph, Capacity,
    78                              CapacityEdgeMap, CapacityEdgeMap > ResGraph;
    79 
    80   public:
    81 
    82     /// The type of the flow map.
    83     typedef typename Graph::template EdgeMap<Capacity> FlowMap;
    84     /// The type of the potential map.
    85     typedef typename Graph::template NodeMap<Cost> PotentialMap;
    86 
    87   private:
    88 
    89     /// \brief Map adaptor class for handling residual edge costs.
    90     ///
    91     /// Map adaptor class for handling residual edge costs.
    92     class ResidualCostMap : public MapBase<typename ResGraph::Edge, Cost>
    93     {
    94       typedef typename ResGraph::Edge Edge;
    95       
    96     private:
    97 
    98       const CostMap &_cost_map;
    99 
   100     public:
   101 
   102       ///\e
   103       ResidualCostMap(const CostMap &cost_map) : _cost_map(cost_map) {}
   104 
   105       ///\e
   106       Cost operator[](const Edge &e) const {
   107         return ResGraph::forward(e) ? _cost_map[e] : -_cost_map[e];
   108       }
   109 
   110     }; //class ResidualCostMap
   111 
   112     /// \brief Map adaptor class for handling reduced edge costs.
   113     ///
   114     /// Map adaptor class for handling reduced edge costs.
   115     class ReducedCostMap : public MapBase<Edge, Cost>
   116     {
   117     private:
   118 
   119       const Graph &_gr;
   120       const CostMap &_cost_map;
   121       const PotentialMap &_pot_map;
   122 
   123     public:
   124 
   125       ///\e
   126       ReducedCostMap( const Graph &gr,
   127                       const CostMap &cost_map,
   128                       const PotentialMap &pot_map ) :
   129         _gr(gr), _cost_map(cost_map), _pot_map(pot_map) {}
   130 
   131       ///\e
   132       inline Cost operator[](const Edge &e) const {
   133         return _cost_map[e] + _pot_map[_gr.source(e)]
   134                             - _pot_map[_gr.target(e)];
   135       }
   136 
   137     }; //class ReducedCostMap
   138 
   139     struct BFOperationTraits {
   140       static double zero() { return 0; }
   141 
   142       static double infinity() {
   143         return std::numeric_limits<double>::infinity();
   144       }
   145 
   146       static double plus(const double& left, const double& right) {
   147         return left + right;
   148       }
   149 
   150       static bool less(const double& left, const double& right) {
   151         return left + 1e-6 < right;
   152       }
   153     }; // class BFOperationTraits
   154 
   155   private:
   156 
   157     // The directed graph the algorithm runs on
   158     const Graph &_graph;
   159     // The original lower bound map
   160     const LowerMap *_lower;
   161     // The modified capacity map
   162     CapacityEdgeMap _capacity;
   163     // The original cost map
   164     const CostMap &_cost;
   165     // The modified supply map
   166     SupplyNodeMap _supply;
   167     bool _valid_supply;
   168 
   169     // Edge map of the current flow
   170     FlowMap *_flow;
   171     bool _local_flow;
   172     // Node map of the current potentials
   173     PotentialMap *_potential;
   174     bool _local_potential;
   175 
   176     // The residual graph
   177     ResGraph *_res_graph;
   178     // The residual cost map
   179     ResidualCostMap _res_cost;
   180 
   181   public:
   182 
   183     /// \brief General constructor (with lower bounds).
   184     ///
   185     /// General constructor (with lower bounds).
   186     ///
   187     /// \param graph The directed graph the algorithm runs on.
   188     /// \param lower The lower bounds of the edges.
   189     /// \param capacity The capacities (upper bounds) of the edges.
   190     /// \param cost The cost (length) values of the edges.
   191     /// \param supply The supply values of the nodes (signed).
   192     CancelAndTighten( const Graph &graph,
   193                       const LowerMap &lower,
   194                       const CapacityMap &capacity,
   195                       const CostMap &cost,
   196                       const SupplyMap &supply ) :
   197       _graph(graph), _lower(&lower), _capacity(capacity), _cost(cost),
   198       _supply(supply), _flow(NULL), _local_flow(false),
   199       _potential(NULL), _local_potential(false),
   200       _res_graph(NULL), _res_cost(_cost)
   201     {
   202       // Check the sum of supply values
   203       Supply sum = 0;
   204       for (NodeIt n(_graph); n != INVALID; ++n) sum += _supply[n];
   205       _valid_supply = sum == 0;
   206 
   207       // Remove non-zero lower bounds
   208       for (EdgeIt e(_graph); e != INVALID; ++e) {
   209         if (lower[e] != 0) {
   210           _capacity[e] -= lower[e];
   211           _supply[_graph.source(e)] -= lower[e];
   212           _supply[_graph.target(e)] += lower[e];
   213         }
   214       }
   215     }
   216 
   217     /// \brief General constructor (without lower bounds).
   218     ///
   219     /// General constructor (without lower bounds).
   220     ///
   221     /// \param graph The directed graph the algorithm runs on.
   222     /// \param capacity The capacities (upper bounds) of the edges.
   223     /// \param cost The cost (length) values of the edges.
   224     /// \param supply The supply values of the nodes (signed).
   225     CancelAndTighten( const Graph &graph,
   226                       const CapacityMap &capacity,
   227                       const CostMap &cost,
   228                       const SupplyMap &supply ) :
   229       _graph(graph), _lower(NULL), _capacity(capacity), _cost(cost),
   230       _supply(supply), _flow(NULL), _local_flow(false),
   231       _potential(NULL), _local_potential(false),
   232       _res_graph(NULL), _res_cost(_cost)
   233     {
   234       // Check the sum of supply values
   235       Supply sum = 0;
   236       for (NodeIt n(_graph); n != INVALID; ++n) sum += _supply[n];
   237       _valid_supply = sum == 0;
   238     }
   239 
   240     /// \brief Simple constructor (with lower bounds).
   241     ///
   242     /// Simple constructor (with lower bounds).
   243     ///
   244     /// \param graph The directed graph the algorithm runs on.
   245     /// \param lower The lower bounds of the edges.
   246     /// \param capacity The capacities (upper bounds) of the edges.
   247     /// \param cost The cost (length) values of the edges.
   248     /// \param s The source node.
   249     /// \param t The target node.
   250     /// \param flow_value The required amount of flow from node \c s
   251     /// to node \c t (i.e. the supply of \c s and the demand of \c t).
   252     CancelAndTighten( const Graph &graph,
   253                       const LowerMap &lower,
   254                       const CapacityMap &capacity,
   255                       const CostMap &cost,
   256                       Node s, Node t,
   257                       Supply flow_value ) :
   258       _graph(graph), _lower(&lower), _capacity(capacity), _cost(cost),
   259       _supply(graph, 0), _flow(NULL), _local_flow(false),
   260       _potential(NULL), _local_potential(false),
   261       _res_graph(NULL), _res_cost(_cost)
   262     {
   263       // Remove non-zero lower bounds
   264       _supply[s] =  flow_value;
   265       _supply[t] = -flow_value;
   266       for (EdgeIt e(_graph); e != INVALID; ++e) {
   267         if (lower[e] != 0) {
   268           _capacity[e] -= lower[e];
   269           _supply[_graph.source(e)] -= lower[e];
   270           _supply[_graph.target(e)] += lower[e];
   271         }
   272       }
   273       _valid_supply = true;
   274     }
   275 
   276     /// \brief Simple constructor (without lower bounds).
   277     ///
   278     /// Simple constructor (without lower bounds).
   279     ///
   280     /// \param graph The directed graph the algorithm runs on.
   281     /// \param capacity The capacities (upper bounds) of the edges.
   282     /// \param cost The cost (length) values of the edges.
   283     /// \param s The source node.
   284     /// \param t The target node.
   285     /// \param flow_value The required amount of flow from node \c s
   286     /// to node \c t (i.e. the supply of \c s and the demand of \c t).
   287     CancelAndTighten( const Graph &graph,
   288                       const CapacityMap &capacity,
   289                       const CostMap &cost,
   290                       Node s, Node t,
   291                       Supply flow_value ) :
   292       _graph(graph), _lower(NULL), _capacity(capacity), _cost(cost),
   293       _supply(graph, 0), _flow(NULL), _local_flow(false),
   294       _potential(NULL), _local_potential(false),
   295       _res_graph(NULL), _res_cost(_cost)
   296     {
   297       _supply[s] =  flow_value;
   298       _supply[t] = -flow_value;
   299       _valid_supply = true;
   300     }
   301 
   302     /// Destructor.
   303     ~CancelAndTighten() {
   304       if (_local_flow) delete _flow;
   305       if (_local_potential) delete _potential;
   306       delete _res_graph;
   307     }
   308 
   309     /// \brief Set the flow map.
   310     ///
   311     /// Set the flow map.
   312     ///
   313     /// \return \c (*this)
   314     CancelAndTighten& flowMap(FlowMap &map) {
   315       if (_local_flow) {
   316         delete _flow;
   317         _local_flow = false;
   318       }
   319       _flow = &map;
   320       return *this;
   321     }
   322 
   323     /// \brief Set the potential map.
   324     ///
   325     /// Set the potential map.
   326     ///
   327     /// \return \c (*this)
   328     CancelAndTighten& potentialMap(PotentialMap &map) {
   329       if (_local_potential) {
   330         delete _potential;
   331         _local_potential = false;
   332       }
   333       _potential = &map;
   334       return *this;
   335     }
   336 
   337     /// \name Execution control
   338 
   339     /// @{
   340 
   341     /// \brief Run the algorithm.
   342     ///
   343     /// Run the algorithm.
   344     ///
   345     /// \return \c true if a feasible flow can be found.
   346     bool run() {
   347       return init() && start();
   348     }
   349 
   350     /// @}
   351 
   352     /// \name Query Functions
   353     /// The result of the algorithm can be obtained using these
   354     /// functions.\n
   355     /// \ref lemon::CancelAndTighten::run() "run()" must be called before
   356     /// using them.
   357 
   358     /// @{
   359 
   360     /// \brief Return a const reference to the edge map storing the
   361     /// found flow.
   362     ///
   363     /// Return a const reference to the edge map storing the found flow.
   364     ///
   365     /// \pre \ref run() must be called before using this function.
   366     const FlowMap& flowMap() const {
   367       return *_flow;
   368     }
   369 
   370     /// \brief Return a const reference to the node map storing the
   371     /// found potentials (the dual solution).
   372     ///
   373     /// Return a const reference to the node map storing the found
   374     /// potentials (the dual solution).
   375     ///
   376     /// \pre \ref run() must be called before using this function.
   377     const PotentialMap& potentialMap() const {
   378       return *_potential;
   379     }
   380 
   381     /// \brief Return the flow on the given edge.
   382     ///
   383     /// Return the flow on the given edge.
   384     ///
   385     /// \pre \ref run() must be called before using this function.
   386     Capacity flow(const Edge& edge) const {
   387       return (*_flow)[edge];
   388     }
   389 
   390     /// \brief Return the potential of the given node.
   391     ///
   392     /// Return the potential of the given node.
   393     ///
   394     /// \pre \ref run() must be called before using this function.
   395     Cost potential(const Node& node) const {
   396       return (*_potential)[node];
   397     }
   398 
   399     /// \brief Return the total cost of the found flow.
   400     ///
   401     /// Return the total cost of the found flow. The complexity of the
   402     /// function is \f$ O(e) \f$.
   403     ///
   404     /// \pre \ref run() must be called before using this function.
   405     Cost totalCost() const {
   406       Cost c = 0;
   407       for (EdgeIt e(_graph); e != INVALID; ++e)
   408         c += (*_flow)[e] * _cost[e];
   409       return c;
   410     }
   411 
   412     /// @}
   413 
   414   private:
   415 
   416     /// Initialize the algorithm.
   417     bool init() {
   418       if (!_valid_supply) return false;
   419 
   420       // Initialize flow and potential maps
   421       if (!_flow) {
   422         _flow = new FlowMap(_graph);
   423         _local_flow = true;
   424       }
   425       if (!_potential) {
   426         _potential = new PotentialMap(_graph);
   427         _local_potential = true;
   428       }
   429 
   430       _res_graph = new ResGraph(_graph, _capacity, *_flow);
   431 
   432       // Find a feasible flow using Circulation
   433       Circulation< Graph, ConstMap<Edge, Capacity>,
   434                    CapacityEdgeMap, SupplyMap >
   435         circulation( _graph, constMap<Edge>(Capacity(0)),
   436                      _capacity, _supply );
   437       return circulation.flowMap(*_flow).run();
   438     }
   439 
   440     bool start() {
   441       const double LIMIT_FACTOR = 0.01;
   442       const int MIN_LIMIT = 3;
   443 
   444       typedef typename Graph::template NodeMap<double> FloatPotentialMap;
   445       typedef typename Graph::template NodeMap<int> LevelMap;
   446       typedef typename Graph::template NodeMap<bool> BoolNodeMap;
   447       typedef typename Graph::template NodeMap<Node> PredNodeMap;
   448       typedef typename Graph::template NodeMap<Edge> PredEdgeMap;
   449       typedef typename ResGraph::template EdgeMap<double> ResShiftCostMap;
   450       FloatPotentialMap pi(_graph);
   451       LevelMap level(_graph);
   452       BoolNodeMap reached(_graph);
   453       BoolNodeMap processed(_graph);
   454       PredNodeMap pred_node(_graph);
   455       PredEdgeMap pred_edge(_graph);
   456       int node_num = countNodes(_graph);
   457       typedef std::pair<Edge, bool> pair;
   458       std::vector<pair> stack(node_num);
   459       std::vector<Node> proc_vector(node_num);
   460       ResShiftCostMap shift_cost(*_res_graph);
   461 
   462       Tolerance<double> tol;
   463       tol.epsilon(1e-6);
   464 
   465       Timer t1, t2, t3;
   466       t1.reset();
   467       t2.reset();
   468       t3.reset();
   469 
   470       // Initialize epsilon and the node potentials
   471       double epsilon = 0;
   472       for (EdgeIt e(_graph); e != INVALID; ++e) {
   473         if (_capacity[e] - (*_flow)[e] > 0 && _cost[e] < -epsilon)
   474           epsilon = -_cost[e];
   475         else if ((*_flow)[e] > 0 && _cost[e] > epsilon)
   476           epsilon = _cost[e];
   477       }
   478       for (NodeIt v(_graph); v != INVALID; ++v) {
   479         pi[v] = 0;
   480       }
   481 
   482       // Start phases
   483       int limit = int(LIMIT_FACTOR * node_num);
   484       if (limit < MIN_LIMIT) limit = MIN_LIMIT;
   485       int iter = limit;
   486       while (epsilon * node_num >= 1) {
   487         t1.start();
   488         // Find and cancel cycles in the admissible graph using DFS
   489         for (NodeIt n(_graph); n != INVALID; ++n) {
   490           reached[n] = false;
   491           processed[n] = false;
   492         }
   493         int stack_head = -1;
   494         int proc_head = -1;
   495 
   496         for (NodeIt start(_graph); start != INVALID; ++start) {
   497           if (reached[start]) continue;
   498 
   499           // New start node
   500           reached[start] = true;
   501           pred_edge[start] = INVALID;
   502           pred_node[start] = INVALID;
   503 
   504           // Find the first admissible residual outgoing edge
   505           double p = pi[start];
   506           Edge e;
   507           _graph.firstOut(e, start);
   508           while ( e != INVALID && (_capacity[e] - (*_flow)[e] == 0 ||
   509                   !tol.negative(_cost[e] + p - pi[_graph.target(e)])) )
   510             _graph.nextOut(e);
   511           if (e != INVALID) {
   512             stack[++stack_head] = pair(e, true);
   513             goto next_step_1;
   514           }
   515           _graph.firstIn(e, start);
   516           while ( e != INVALID && ((*_flow)[e] == 0 ||
   517                   !tol.negative(-_cost[e] + p - pi[_graph.source(e)])) )
   518             _graph.nextIn(e);
   519           if (e != INVALID) {
   520             stack[++stack_head] = pair(e, false);
   521             goto next_step_1;
   522           }
   523           processed[start] = true;
   524           proc_vector[++proc_head] = start;
   525           continue;
   526         next_step_1:
   527 
   528           while (stack_head >= 0) {
   529             Edge se = stack[stack_head].first;
   530             bool sf = stack[stack_head].second;
   531             Node u, v;
   532             if (sf) {
   533               u = _graph.source(se);
   534               v = _graph.target(se);
   535             } else {
   536               u = _graph.target(se);
   537               v = _graph.source(se);
   538             }
   539 
   540             if (!reached[v]) {
   541               // A new node is reached
   542               reached[v] = true;
   543               pred_node[v] = u;
   544               pred_edge[v] = se;
   545               // Find the first admissible residual outgoing edge
   546               double p = pi[v];
   547               Edge e;
   548               _graph.firstOut(e, v);
   549               while ( e != INVALID && (_capacity[e] - (*_flow)[e] == 0 ||
   550                       !tol.negative(_cost[e] + p - pi[_graph.target(e)])) )
   551                 _graph.nextOut(e);
   552               if (e != INVALID) {
   553                 stack[++stack_head] = pair(e, true);
   554                 goto next_step_2;
   555               }
   556               _graph.firstIn(e, v);
   557               while ( e != INVALID && ((*_flow)[e] == 0 ||
   558                       !tol.negative(-_cost[e] + p - pi[_graph.source(e)])) )
   559                 _graph.nextIn(e);
   560               stack[++stack_head] = pair(e, false);
   561             next_step_2: ;
   562             } else {
   563               if (!processed[v]) {
   564                 // A cycle is found
   565                 Node n, w = u;
   566                 Capacity d, delta = sf ? _capacity[se] - (*_flow)[se] :
   567                                          (*_flow)[se];
   568                 for (n = u; n != v; n = pred_node[n]) {
   569                   d = _graph.target(pred_edge[n]) == n ?
   570                       _capacity[pred_edge[n]] - (*_flow)[pred_edge[n]] :
   571                       (*_flow)[pred_edge[n]];
   572                   if (d <= delta) {
   573                     delta = d;
   574                     w = pred_node[n];
   575                   }
   576                 }
   577 
   578 /*
   579                 std::cout << "CYCLE FOUND: ";
   580                 if (sf)
   581                   std::cout << _cost[se] + pi[_graph.source(se)] - pi[_graph.target(se)];
   582                 else
   583                   std::cout << _graph.id(se) << ":" << -(_cost[se] + pi[_graph.source(se)] - pi[_graph.target(se)]);
   584                 for (n = u; n != v; n = pred_node[n]) {
   585                   if (_graph.target(pred_edge[n]) == n)
   586                     std::cout << " " << _cost[pred_edge[n]] + pi[_graph.source(pred_edge[n])] - pi[_graph.target(pred_edge[n])];
   587                   else
   588                     std::cout << " " << -(_cost[pred_edge[n]] + pi[_graph.source(pred_edge[n])] - pi[_graph.target(pred_edge[n])]);
   589                 }
   590                 std::cout << "\n";
   591 */
   592                 // Augment along the cycle
   593                 (*_flow)[se] = sf ? (*_flow)[se] + delta :
   594                                     (*_flow)[se] - delta;
   595                 for (n = u; n != v; n = pred_node[n]) {
   596                   if (_graph.target(pred_edge[n]) == n)
   597                     (*_flow)[pred_edge[n]] += delta;
   598                   else
   599                     (*_flow)[pred_edge[n]] -= delta;
   600                 }
   601                 for (n = u; stack_head > 0 && n != w; n = pred_node[n]) {
   602                   --stack_head;
   603                   reached[n] = false;
   604                 }
   605                 u = w;
   606               }
   607               v = u;
   608 
   609               // Find the next admissible residual outgoing edge
   610               double p = pi[v];
   611               Edge e = stack[stack_head].first;
   612               if (!stack[stack_head].second) {
   613                 _graph.nextIn(e);
   614                 goto in_edge_3;
   615               }
   616               _graph.nextOut(e);
   617               while ( e != INVALID && (_capacity[e] - (*_flow)[e] == 0 ||
   618                       !tol.negative(_cost[e] + p - pi[_graph.target(e)])) )
   619                 _graph.nextOut(e);
   620               if (e != INVALID) {
   621                 stack[stack_head] = pair(e, true);
   622                 goto next_step_3;
   623               }
   624               _graph.firstIn(e, v);
   625             in_edge_3:
   626               while ( e != INVALID && ((*_flow)[e] == 0 ||
   627                       !tol.negative(-_cost[e] + p - pi[_graph.source(e)])) )
   628                 _graph.nextIn(e);
   629               stack[stack_head] = pair(e, false);
   630             next_step_3: ;
   631             }
   632 
   633             while (stack_head >= 0 && stack[stack_head].first == INVALID) {
   634               processed[v] = true;
   635               proc_vector[++proc_head] = v;
   636               if (--stack_head >= 0) {
   637                 v = stack[stack_head].second ?
   638                     _graph.source(stack[stack_head].first) :
   639                     _graph.target(stack[stack_head].first);
   640                 // Find the next admissible residual outgoing edge
   641                 double p = pi[v];
   642                 Edge e = stack[stack_head].first;
   643                 if (!stack[stack_head].second) {
   644                   _graph.nextIn(e);
   645                   goto in_edge_4;
   646                 }
   647                 _graph.nextOut(e);
   648                 while ( e != INVALID && (_capacity[e] - (*_flow)[e] == 0 ||
   649                         !tol.negative(_cost[e] + p - pi[_graph.target(e)])) )
   650                   _graph.nextOut(e);
   651                 if (e != INVALID) {
   652                   stack[stack_head] = pair(e, true);
   653                   goto next_step_4;
   654                 }
   655                 _graph.firstIn(e, v);
   656               in_edge_4:
   657                 while ( e != INVALID && ((*_flow)[e] == 0 ||
   658                         !tol.negative(-_cost[e] + p - pi[_graph.source(e)])) )
   659                   _graph.nextIn(e);
   660                 stack[stack_head] = pair(e, false);
   661               next_step_4: ;
   662               }
   663             }
   664           }
   665         }
   666         t1.stop();
   667 
   668         // Tighten potentials and epsilon
   669         if (--iter > 0) {
   670           // Compute levels
   671           t2.start();
   672           for (int i = proc_head; i >= 0; --i) {
   673             Node v = proc_vector[i];
   674             double p = pi[v];
   675             int l = 0;
   676             for (InEdgeIt e(_graph, v); e != INVALID; ++e) {
   677               Node u = _graph.source(e);
   678               if ( _capacity[e] - (*_flow)[e] > 0 &&
   679                    tol.negative(_cost[e] + pi[u] - p) &&
   680                    level[u] + 1 > l ) l = level[u] + 1;
   681             }
   682             for (OutEdgeIt e(_graph, v); e != INVALID; ++e) {
   683               Node u = _graph.target(e);
   684               if ( (*_flow)[e] > 0 &&
   685                    tol.negative(-_cost[e] + pi[u] - p) &&
   686                    level[u] + 1 > l ) l = level[u] + 1;
   687             }
   688             level[v] = l;
   689           }
   690 
   691           // Modify potentials
   692           double p, q = -1;
   693           for (EdgeIt e(_graph); e != INVALID; ++e) {
   694             Node u = _graph.source(e);
   695             Node v = _graph.target(e);
   696             if (_capacity[e] - (*_flow)[e] > 0 && level[u] - level[v] > 0) {
   697               p = (_cost[e] + pi[u] - pi[v] + epsilon) /
   698                   (level[u] - level[v] + 1);
   699               if (q < 0 || p < q) q = p;
   700             }
   701             else if ((*_flow)[e] > 0 && level[v] - level[u] > 0) {
   702               p = (-_cost[e] - pi[u] + pi[v] + epsilon) /
   703                   (level[v] - level[u] + 1);
   704               if (q < 0 || p < q) q = p;
   705             }
   706           }
   707           for (NodeIt v(_graph); v != INVALID; ++v) {
   708             pi[v] -= q * level[v];
   709           }
   710 
   711           // Modify epsilon
   712           epsilon = 0;
   713           for (EdgeIt e(_graph); e != INVALID; ++e) {
   714             double curr = _cost[e] + pi[_graph.source(e)]
   715                                    - pi[_graph.target(e)];
   716             if (_capacity[e] - (*_flow)[e] > 0 && curr < -epsilon)
   717               epsilon = -curr;
   718             else if ((*_flow)[e] > 0 && curr > epsilon)
   719               epsilon = curr;
   720           }
   721           t2.stop();
   722         } else {
   723           // Set epsilon to the minimum cycle mean
   724           t3.start();
   725 
   726 /**/
   727           StaticGraph static_graph;
   728           typename ResGraph::template NodeMap<typename StaticGraph::Node> node_ref(*_res_graph);
   729           typename ResGraph::template EdgeMap<typename StaticGraph::Edge> edge_ref(*_res_graph);
   730           static_graph.build(*_res_graph, node_ref, edge_ref);
   731           typename StaticGraph::template NodeMap<double> static_pi(static_graph);
   732           typename StaticGraph::template EdgeMap<double> static_cost(static_graph);
   733 
   734           for (typename ResGraph::EdgeIt e(*_res_graph); e != INVALID; ++e)
   735             static_cost[edge_ref[e]] = _res_cost[e];
   736 
   737           MinMeanCycle<StaticGraph, typename StaticGraph::template EdgeMap<double> >
   738             mmc(static_graph, static_cost);
   739           mmc.init();
   740           mmc.findMinMean();
   741           epsilon = -mmc.cycleMean();
   742 /**/
   743 
   744 /*
   745           MinMeanCycle<ResGraph, ResidualCostMap> mmc(*_res_graph, _res_cost);
   746           mmc.init();
   747           mmc.findMinMean();
   748           epsilon = -mmc.cycleMean();
   749 */
   750 
   751           // Compute feasible potentials for the current epsilon
   752           for (typename StaticGraph::EdgeIt e(static_graph); e != INVALID; ++e)
   753             static_cost[e] += epsilon;
   754           typename BellmanFord<StaticGraph, typename StaticGraph::template EdgeMap<double> >::
   755             template DefDistMap<typename StaticGraph::template NodeMap<double> >::
   756             template DefOperationTraits<BFOperationTraits>::Create
   757               bf(static_graph, static_cost);
   758           bf.distMap(static_pi).init(0);
   759           bf.start();
   760           for (NodeIt n(_graph); n != INVALID; ++n)
   761             pi[n] = static_pi[node_ref[n]];
   762           
   763 /*
   764           for (typename ResGraph::EdgeIt e(*_res_graph); e != INVALID; ++e)
   765             shift_cost[e] = _res_cost[e] + epsilon;
   766           typename BellmanFord<ResGraph, ResShiftCostMap>::
   767             template DefDistMap<FloatPotentialMap>::
   768             template DefOperationTraits<BFOperationTraits>::Create
   769               bf(*_res_graph, shift_cost);
   770           bf.distMap(pi).init(0);
   771           bf.start();
   772 */
   773 
   774           iter = limit;
   775           t3.stop();
   776         }
   777       }
   778 
   779 //      std::cout << t1.realTime() << " " << t2.realTime() << " " << t3.realTime() << "\n";
   780 
   781       // Handle non-zero lower bounds
   782       if (_lower) {
   783         for (EdgeIt e(_graph); e != INVALID; ++e)
   784           (*_flow)[e] += (*_lower)[e];
   785       }
   786       return true;
   787     }
   788 
   789   }; //class CancelAndTighten
   790 
   791   ///@}
   792 
   793 } //namespace lemon
   794 
   795 #endif //LEMON_CANCEL_AND_TIGHTEN_H