/* -*- mode: C++; indent-tabs-mode: nil; -*-
* This file is a part of LEMON, a generic C++ optimization library.
* Copyright (C) 2003-2009
* Egervary Jeno Kombinatorikus Optimalizalasi Kutatocsoport
* (Egervary Research Group on Combinatorial Optimization, EGRES).
* Permission to use, modify and distribute this software is granted
* provided that this copyright notice appears in all copies. For
* precise terms see the accompanying LICENSE file.
* This software is provided "AS IS" with no warranty of any kind,
* express or implied, and with no claim as to its suitability for any
#ifndef LEMON_NETWORK_SIMPLEX_H
#define LEMON_NETWORK_SIMPLEX_H
/// \ingroup min_cost_flow
/// \brief Network simplex algorithm for finding a minimum cost flow.
/// \addtogroup min_cost_flow
/// \brief Implementation of the primal network simplex algorithm
/// for finding a \ref min_cost_flow "minimum cost flow".
/// \ref NetworkSimplex implements the primal network simplex algorithm
/// for finding a \ref min_cost_flow "minimum cost flow".
/// \tparam Digraph The digraph type the algorithm runs on.
/// \tparam LowerMap The type of the lower bound map.
/// \tparam CapacityMap The type of the capacity (upper bound) map.
/// \tparam CostMap The type of the cost (length) map.
/// \tparam SupplyMap The type of the supply map.
/// - Arc capacities and costs should be \e non-negative \e integers.
/// - Supply values should be \e signed \e integers.
/// - The value types of the maps should be convertible to each other.
/// - \c CostMap::Value must be signed type.
/// \note \ref NetworkSimplex provides five different pivot rule
/// implementations that significantly affect the efficiency of the
/// By default "Block Search" pivot rule is used, which proved to be
/// by far the most efficient according to our benchmark tests.
/// However another pivot rule can be selected using \ref run()
/// function with the proper parameter.
template < typename Digraph,
template < typename Digraph,
typename LowerMap = typename Digraph::template ArcMap<int>,
typename CapacityMap = typename Digraph::template ArcMap<int>,
typename CostMap = typename Digraph::template ArcMap<int>,
typename SupplyMap = typename Digraph::template NodeMap<int> >
TEMPLATE_DIGRAPH_TYPEDEFS(Digraph);
typedef typename CapacityMap::Value Capacity;
typedef typename CostMap::Value Cost;
typedef typename SupplyMap::Value Supply;
typedef std::vector<Arc> ArcVector;
typedef std::vector<Node> NodeVector;
typedef std::vector<int> IntVector;
typedef std::vector<bool> BoolVector;
typedef std::vector<Capacity> CapacityVector;
typedef std::vector<Cost> CostVector;
typedef std::vector<Supply> SupplyVector;
/// The type of the flow map
typedef typename Digraph::template ArcMap<Capacity> FlowMap;
/// The type of the potential map
typedef typename Digraph::template NodeMap<Cost> PotentialMap;
/// Enum type for selecting the pivot rule used by \ref run()
// State constants for arcs
// References for the original data
const Digraph &_orig_graph;
const LowerMap *_orig_lower;
const CapacityMap &_orig_cap;
const CostMap &_orig_cost;
const SupplyMap *_orig_supply;
Capacity _orig_flow_value;
PotentialMap *_potential_result;
// Data structures for storing the graph
// The number of nodes and arcs in the original graph
// Node and arc maps for the spanning tree structure
// The entering arc in the current pivot iteration
// Temporary data used in the current pivot iteration
int join, u_in, v_in, u_out, v_out;
int right, first, second, last;
int stem, par_stem, new_stem;
/// \brief Implementation of the "First Eligible" pivot rule for the
/// \ref NetworkSimplex "network simplex" algorithm.
/// This class implements the "First Eligible" pivot rule
/// for the \ref NetworkSimplex "network simplex" algorithm.
/// For more information see \ref NetworkSimplex::run().
class FirstEligiblePivotRule
// References to the NetworkSimplex class
const IntVector &_source;
const IntVector &_target;
FirstEligiblePivotRule(NetworkSimplex &ns) :
_arc(ns._arc), _source(ns._source), _target(ns._target),
_cost(ns._cost), _state(ns._state), _pi(ns._pi),
_in_arc(ns._in_arc), _arc_num(ns._arc_num), _next_arc(0)
/// Find next entering arc
for (int e = _next_arc; e < _arc_num; ++e) {
c = _state[e] * (_cost[e] + _pi[_source[e]] - _pi[_target[e]]);
for (int e = 0; e < _next_arc; ++e) {
c = _state[e] * (_cost[e] + _pi[_source[e]] - _pi[_target[e]]);
}; //class FirstEligiblePivotRule
/// \brief Implementation of the "Best Eligible" pivot rule for the
/// \ref NetworkSimplex "network simplex" algorithm.
/// This class implements the "Best Eligible" pivot rule
/// for the \ref NetworkSimplex "network simplex" algorithm.
/// For more information see \ref NetworkSimplex::run().
class BestEligiblePivotRule
// References to the NetworkSimplex class
const IntVector &_source;
const IntVector &_target;
BestEligiblePivotRule(NetworkSimplex &ns) :
_arc(ns._arc), _source(ns._source), _target(ns._target),
_cost(ns._cost), _state(ns._state), _pi(ns._pi),
_in_arc(ns._in_arc), _arc_num(ns._arc_num)
/// Find next entering arc
for (int e = 0; e < _arc_num; ++e) {
c = _state[e] * (_cost[e] + _pi[_source[e]] - _pi[_target[e]]);
}; //class BestEligiblePivotRule
/// \brief Implementation of the "Block Search" pivot rule for the
/// \ref NetworkSimplex "network simplex" algorithm.
/// This class implements the "Block Search" pivot rule
/// for the \ref NetworkSimplex "network simplex" algorithm.
/// For more information see \ref NetworkSimplex::run().
class BlockSearchPivotRule
// References to the NetworkSimplex class
const IntVector &_source;
const IntVector &_target;
BlockSearchPivotRule(NetworkSimplex &ns) :
_arc(ns._arc), _source(ns._source), _target(ns._target),
_cost(ns._cost), _state(ns._state), _pi(ns._pi),
_in_arc(ns._in_arc), _arc_num(ns._arc_num), _next_arc(0)
// The main parameters of the pivot rule
const double BLOCK_SIZE_FACTOR = 2.0;
const int MIN_BLOCK_SIZE = 10;
_block_size = std::max( int(BLOCK_SIZE_FACTOR * sqrt(_arc_num)),
/// Find next entering arc
int e, min_arc = _next_arc;
for (e = _next_arc; e < _arc_num; ++e) {
c = _state[e] * (_cost[e] + _pi[_source[e]] - _pi[_target[e]]);
if (min == 0 || cnt > 0) {
for (e = 0; e < _next_arc; ++e) {
c = _state[e] * (_cost[e] + _pi[_source[e]] - _pi[_target[e]]);
if (min >= 0) return false;
}; //class BlockSearchPivotRule
/// \brief Implementation of the "Candidate List" pivot rule for the
/// \ref NetworkSimplex "network simplex" algorithm.
/// This class implements the "Candidate List" pivot rule
/// for the \ref NetworkSimplex "network simplex" algorithm.
/// For more information see \ref NetworkSimplex::run().
class CandidateListPivotRule
// References to the NetworkSimplex class
const IntVector &_source;
const IntVector &_target;
int _list_length, _minor_limit;
int _curr_length, _minor_count;
CandidateListPivotRule(NetworkSimplex &ns) :
_arc(ns._arc), _source(ns._source), _target(ns._target),
_cost(ns._cost), _state(ns._state), _pi(ns._pi),
_in_arc(ns._in_arc), _arc_num(ns._arc_num), _next_arc(0)
// The main parameters of the pivot rule
const double LIST_LENGTH_FACTOR = 1.0;
const int MIN_LIST_LENGTH = 10;
const double MINOR_LIMIT_FACTOR = 0.1;
const int MIN_MINOR_LIMIT = 3;
_list_length = std::max( int(LIST_LENGTH_FACTOR * sqrt(_arc_num)),
_minor_limit = std::max( int(MINOR_LIMIT_FACTOR * _list_length),
_curr_length = _minor_count = 0;
_candidates.resize(_list_length);
/// Find next entering arc
int e, min_arc = _next_arc;
if (_curr_length > 0 && _minor_count < _minor_limit) {
// Minor iteration: select the best eligible arc from the
// current candidate list
for (int i = 0; i < _curr_length; ++i) {
c = _state[e] * (_cost[e] + _pi[_source[e]] - _pi[_target[e]]);
_candidates[i--] = _candidates[--_curr_length];
// Major iteration: build a new candidate list
for (e = _next_arc; e < _arc_num; ++e) {
c = _state[e] * (_cost[e] + _pi[_source[e]] - _pi[_target[e]]);
_candidates[_curr_length++] = e;
if (_curr_length == _list_length) break;
if (_curr_length < _list_length) {
for (e = 0; e < _next_arc; ++e) {
c = _state[e] * (_cost[e] + _pi[_source[e]] - _pi[_target[e]]);
_candidates[_curr_length++] = e;
if (_curr_length == _list_length) break;
if (_curr_length == 0) return false;
}; //class CandidateListPivotRule
/// \brief Implementation of the "Altering Candidate List" pivot rule
/// for the \ref NetworkSimplex "network simplex" algorithm.
/// This class implements the "Altering Candidate List" pivot rule
/// for the \ref NetworkSimplex "network simplex" algorithm.
/// For more information see \ref NetworkSimplex::run().
class AlteringListPivotRule
// References to the NetworkSimplex class
const IntVector &_source;
const IntVector &_target;
int _block_size, _head_length, _curr_length;
// Functor class to compare arcs during sort of the candidate list
SortFunc(const CostVector &map) : _map(map) {}
bool operator()(int left, int right) {
return _map[left] > _map[right];
AlteringListPivotRule(NetworkSimplex &ns) :
_arc(ns._arc), _source(ns._source), _target(ns._target),
_cost(ns._cost), _state(ns._state), _pi(ns._pi),
_in_arc(ns._in_arc), _arc_num(ns._arc_num),
_next_arc(0), _cand_cost(ns._arc_num), _sort_func(_cand_cost)
// The main parameters of the pivot rule
const double BLOCK_SIZE_FACTOR = 1.5;
const int MIN_BLOCK_SIZE = 10;
const double HEAD_LENGTH_FACTOR = 0.1;
const int MIN_HEAD_LENGTH = 3;
_block_size = std::max( int(BLOCK_SIZE_FACTOR * sqrt(_arc_num)),
_head_length = std::max( int(HEAD_LENGTH_FACTOR * _block_size),
_candidates.resize(_head_length + _block_size);
/// Find next entering arc
// Check the current candidate list
for (int i = 0; i < _curr_length; ++i) {
_cand_cost[e] = _state[e] *
(_cost[e] + _pi[_source[e]] - _pi[_target[e]]);
if (_cand_cost[e] >= 0) {
_candidates[i--] = _candidates[--_curr_length];
int limit = _head_length;
for (int e = _next_arc; e < _arc_num; ++e) {
_cand_cost[e] = _state[e] *
(_cost[e] + _pi[_source[e]] - _pi[_target[e]]);
_candidates[_curr_length++] = e;
if (_curr_length > limit) break;
if (_curr_length <= limit) {
for (int e = 0; e < _next_arc; ++e) {
_cand_cost[e] = _state[e] *
(_cost[e] + _pi[_source[e]] - _pi[_target[e]]);
_candidates[_curr_length++] = e;
if (_curr_length > limit) break;
if (_curr_length == 0) return false;
_next_arc = last_edge + 1;
// Make heap of the candidate list (approximating a partial sort)
make_heap( _candidates.begin(), _candidates.begin() + _curr_length,
// Pop the first element of the heap
_in_arc = _candidates[0];
pop_heap( _candidates.begin(), _candidates.begin() + _curr_length,
_curr_length = std::min(_head_length, _curr_length - 1);
}; //class AlteringListPivotRule
/// \brief General constructor (with lower bounds).
/// General constructor (with lower bounds).
/// \param digraph The digraph the algorithm runs on.
/// \param lower The lower bounds of the arcs.
/// \param capacity The capacities (upper bounds) of the arcs.
/// \param cost The cost (length) values of the arcs.
/// \param supply The supply values of the nodes (signed).
NetworkSimplex( const Digraph &digraph,
const CapacityMap &capacity,
const SupplyMap &supply ) :
_orig_graph(digraph), _orig_lower(&lower), _orig_cap(capacity),
_orig_cost(cost), _orig_supply(&supply),
_flow_result(NULL), _potential_result(NULL),
_local_flow(false), _local_potential(false),
/// \brief General constructor (without lower bounds).
/// General constructor (without lower bounds).
/// \param digraph The digraph the algorithm runs on.
/// \param capacity The capacities (upper bounds) of the arcs.
/// \param cost The cost (length) values of the arcs.
/// \param supply The supply values of the nodes (signed).
NetworkSimplex( const Digraph &digraph,
const CapacityMap &capacity,
const SupplyMap &supply ) :
_orig_graph(digraph), _orig_lower(NULL), _orig_cap(capacity),
_orig_cost(cost), _orig_supply(&supply),
_flow_result(NULL), _potential_result(NULL),
_local_flow(false), _local_potential(false),
/// \brief Simple constructor (with lower bounds).
/// Simple constructor (with lower bounds).
/// \param digraph The digraph the algorithm runs on.
/// \param lower The lower bounds of the arcs.
/// \param capacity The capacities (upper bounds) of the arcs.
/// \param cost The cost (length) values of the arcs.
/// \param s The source node.
/// \param t The target node.
/// \param flow_value The required amount of flow from node \c s
/// to node \c t (i.e. the supply of \c s and the demand of \c t).
NetworkSimplex( const Digraph &digraph,
const CapacityMap &capacity,
_orig_graph(digraph), _orig_lower(&lower), _orig_cap(capacity),
_orig_cost(cost), _orig_supply(NULL),
_orig_source(s), _orig_target(t), _orig_flow_value(flow_value),
_flow_result(NULL), _potential_result(NULL),
_local_flow(false), _local_potential(false),
/// \brief Simple constructor (without lower bounds).
/// Simple constructor (without lower bounds).
/// \param digraph The digraph the algorithm runs on.
/// \param capacity The capacities (upper bounds) of the arcs.
/// \param cost The cost (length) values of the arcs.
/// \param s The source node.
/// \param t The target node.
/// \param flow_value The required amount of flow from node \c s
/// to node \c t (i.e. the supply of \c s and the demand of \c t).
NetworkSimplex( const Digraph &digraph,
const CapacityMap &capacity,
_orig_graph(digraph), _orig_lower(NULL), _orig_cap(capacity),
_orig_cost(cost), _orig_supply(NULL),
_orig_source(s), _orig_target(t), _orig_flow_value(flow_value),
_flow_result(NULL), _potential_result(NULL),
_local_flow(false), _local_potential(false),
if (_local_flow) delete _flow_result;
if (_local_potential) delete _potential_result;
/// \brief Set the flow map.
/// This function sets the flow map.
/// \return <tt>(*this)</tt>
NetworkSimplex& flowMap(FlowMap &map) {
/// \brief Set the potential map.
/// This function sets the potential map.
/// \return <tt>(*this)</tt>
NetworkSimplex& potentialMap(PotentialMap &map) {
delete _potential_result;
_local_potential = false;
_potential_result = ↦
/// \name Execution control
/// The algorithm can be executed using the
/// \ref lemon::NetworkSimplex::run() "run()" function.
/// \brief Run the algorithm.
/// This function runs the algorithm.
/// \param pivot_rule The pivot rule that is used during the
/// The available pivot rules:
/// - FIRST_ELIGIBLE_PIVOT The next eligible arc is selected in
/// a wraparound fashion in every iteration
/// (\ref FirstEligiblePivotRule).
/// - BEST_ELIGIBLE_PIVOT The best eligible arc is selected in
/// every iteration (\ref BestEligiblePivotRule).
/// - BLOCK_SEARCH_PIVOT A specified number of arcs are examined in
/// every iteration in a wraparound fashion and the best eligible
/// arc is selected from this block (\ref BlockSearchPivotRule).
/// - CANDIDATE_LIST_PIVOT In a major iteration a candidate list is
/// built from eligible arcs in a wraparound fashion and in the
/// following minor iterations the best eligible arc is selected
/// from this list (\ref CandidateListPivotRule).
/// - ALTERING_LIST_PIVOT It is a modified version of the
/// "Candidate List" pivot rule. It keeps only the several best
/// eligible arcs from the former candidate list and extends this
/// list in every iteration (\ref AlteringListPivotRule).
/// According to our comprehensive benchmark tests the "Block Search"
/// pivot rule proved to be the fastest and the most robust on
/// various test inputs. Thus it is the default option.
/// \return \c true if a feasible flow can be found.
bool run(PivotRuleEnum pivot_rule = BLOCK_SEARCH_PIVOT) {
return init() && start(pivot_rule);
/// \name Query Functions
/// The results of the algorithm can be obtained using these
/// \ref lemon::NetworkSimplex::run() "run()" must be called before
/// \brief Return a const reference to the flow map.
/// This function returns a const reference to an arc map storing
/// \pre \ref run() must be called before using this function.
const FlowMap& flowMap() const {
/// \brief Return a const reference to the potential map
/// This function returns a const reference to a node map storing
/// the found potentials (the dual solution).
/// \pre \ref run() must be called before using this function.
const PotentialMap& potentialMap() const {
return *_potential_result;
/// \brief Return the flow on the given arc.
/// This function returns the flow on the given arc.
/// \pre \ref run() must be called before using this function.
Capacity flow(const Arc& arc) const {
return (*_flow_result)[arc];
/// \brief Return the potential of the given node.
/// This function returns the potential of the given node.
/// \pre \ref run() must be called before using this function.
Cost potential(const Node& node) const {
return (*_potential_result)[node];
/// \brief Return the total cost of the found flow.
/// This function returns the total cost of the found flow.
/// The complexity of the function is \f$ O(e) \f$.
/// \pre \ref run() must be called before using this function.
for (ArcIt e(_orig_graph); e != INVALID; ++e)
c += (*_flow_result)[e] * _orig_cost[e];
// Initialize internal data structures
// Initialize result maps
_flow_result = new FlowMap(_orig_graph);
if (!_potential_result) {
_potential_result = new PotentialMap(_orig_graph);
_node_num = countNodes(_orig_graph);
_arc_num = countArcs(_orig_graph);
int all_node_num = _node_num + 1;
int all_edge_num = _arc_num + _node_num;
_node.reserve(_node_num);
_source.resize(all_edge_num);
_target.resize(all_edge_num);
_cap.resize(all_edge_num);
_cost.resize(all_edge_num);
_supply.resize(all_node_num);
_flow.resize(all_edge_num, 0);
_pi.resize(all_node_num, 0);
_depth.resize(all_node_num);
_parent.resize(all_node_num);
_pred.resize(all_node_num);
_thread.resize(all_node_num);
_forward.resize(all_node_num);
_state.resize(all_edge_num, STATE_LOWER);
// Initialize node related data
bool valid_supply = true;
for (NodeIt n(_orig_graph); n != INVALID; ++n, ++i) {
_supply[i] = (*_orig_supply)[n];
valid_supply = (sum == 0);
for (NodeIt n(_orig_graph); n != INVALID; ++n, ++i) {
_supply[_node_id[_orig_source]] = _orig_flow_value;
_supply[_node_id[_orig_target]] = -_orig_flow_value;
if (!valid_supply) return false;
// Set data for the artificial root node
// Store the arcs in a mixed order
int k = std::max(int(sqrt(_arc_num)), 10);
for (ArcIt e(_orig_graph); e != INVALID; ++e) {
if ((i += k) >= _arc_num) i = (i % k) + 1;
for (int i = 0; i != _arc_num; ++i) {
_source[i] = _node_id[_orig_graph.source(e)];
_target[i] = _node_id[_orig_graph.target(e)];
_cost[i] = _orig_cost[e];
// Remove non-zero lower bounds
for (int i = 0; i != _arc_num; ++i) {
Capacity c = (*_orig_lower)[_arc[i]];
_supply[_source[i]] -= c;
_supply[_target[i]] += c;
// Add artificial arcs and initialize the spanning tree data structure
Cost max_cost = std::numeric_limits<Cost>::max() / 4;
Capacity max_cap = std::numeric_limits<Capacity>::max();
for (int u = 0, e = _arc_num; u != _node_num; ++u, ++e) {
int u = _source[_in_arc];
int v = _target[_in_arc];
while (_depth[u] > _depth[v]) u = _parent[u];
while (_depth[v] > _depth[u]) v = _parent[v];
// Find the leaving arc of the cycle and returns true if the
// leaving arc is not the same as the entering arc
// Initialize first and second nodes according to the direction
if (_state[_in_arc] == STATE_LOWER) {
first = _source[_in_arc];
second = _target[_in_arc];
first = _target[_in_arc];
second = _source[_in_arc];
// Search the cycle along the path form the first node to the root
for (int u = first; u != join; u = _parent[u]) {
d = _forward[u] ? _flow[e] : _cap[e] - _flow[e];
// Search the cycle along the path form the second node to the root
for (int u = second; u != join; u = _parent[u]) {
d = _forward[u] ? _cap[e] - _flow[e] : _flow[e];
// Change _flow and _state vectors
void changeFlow(bool change) {
// Augment along the cycle
Capacity val = _state[_in_arc] * delta;
for (int u = _source[_in_arc]; u != join; u = _parent[u]) {
_flow[_pred[u]] += _forward[u] ? -val : val;
for (int u = _target[_in_arc]; u != join; u = _parent[u]) {
_flow[_pred[u]] += _forward[u] ? val : -val;
// Update the state of the entering and leaving arcs
_state[_in_arc] = STATE_TREE;
(_flow[_pred[u_out]] == 0) ? STATE_LOWER : STATE_UPPER;
_state[_in_arc] = -_state[_in_arc];
// Update _thread and _parent vectors
void updateThreadParent() {
// Handle the case when join and v_out coincide
for (u = join; u != u_in && u != v_in; u = _thread[u]) ;
while (_thread[u] != u_out) u = _thread[u];
// Find the last successor of u_in (u) and the node after it (right)
// according to the thread index
for (u = u_in; _depth[_thread[u]] > _depth[u_in]; u = _thread[u]) ;
if (_thread[v_in] == u_out) {
for (last = u; _depth[last] > _depth[u_out]; last = _thread[last]) ;
if (last == u_out) last = _thread[last];
else last = _thread[v_in];
_thread[v_in] = stem = u_in;
_thread[u] = new_stem = _parent[stem];
// Find the node just before the stem node (u) according to
// the original thread index
for (u = new_stem; _thread[u] != stem; u = _thread[u]) ;
// Change the parent node of stem and shift stem and par_stem nodes
_parent[stem] = par_stem;
// Find the last successor of stem (u) and the node after it (right)
// according to the thread index
for (u = stem; _depth[_thread[u]] > _depth[stem]; u = _thread[u]) ;
_parent[u_out] = par_stem;
if (join == v_out && par_first) {
if (first != v_in) _thread[first] = right;
for (u = v_out; _thread[u] != u_out; u = _thread[u]) ;
// Update _pred and _forward vectors
_forward[u] = !_forward[v];
_forward[u_in] = (u_in == _source[_in_arc]);
// Update _depth and _potential vectors
void updateDepthPotential() {
_depth[u_in] = _depth[v_in] + 1;
Cost sigma = _forward[u_in] ?
_pi[v_in] - _pi[u_in] - _cost[_pred[u_in]] :
_pi[v_in] - _pi[u_in] + _cost[_pred[u_in]];
for(int u = _thread[u_in]; _parent[u] != -1; u = _thread[u]) {
_depth[u] = _depth[_parent[u]] + 1;
if (_depth[u] <= _depth[u_in]) break;
bool start(PivotRuleEnum pivot_rule) {
// Select the pivot rule implementation
case FIRST_ELIGIBLE_PIVOT:
return start<FirstEligiblePivotRule>();
case BEST_ELIGIBLE_PIVOT:
return start<BestEligiblePivotRule>();
return start<BlockSearchPivotRule>();
case CANDIDATE_LIST_PIVOT:
return start<CandidateListPivotRule>();
case ALTERING_LIST_PIVOT:
return start<AlteringListPivotRule>();
template<class PivotRuleImplementation>
PivotRuleImplementation pivot(*this);
// Execute the network simplex algorithm
while (pivot.findEnteringArc()) {
bool change = findLeavingArc();
// Check if the flow amount equals zero on all the artificial arcs
for (int e = _arc_num; e != _arc_num + _node_num; ++e) {
if (_flow[e] > 0) return false;
// Copy flow values to _flow_result
for (int i = 0; i != _arc_num; ++i) {
(*_flow_result)[e] = (*_orig_lower)[e] + _flow[i];
for (int i = 0; i != _arc_num; ++i) {
(*_flow_result)[_arc[i]] = _flow[i];
// Copy potential values to _potential_result
for (int i = 0; i != _node_num; ++i) {
(*_potential_result)[_node[i]] = _pi[i];
}; //class NetworkSimplex
#endif //LEMON_NETWORK_SIMPLEX_H