Google OR-Tools v9.12
a fast and portable software suite for combinatorial optimization
Loading...
Searching...
No Matches
generic_max_flow.h
Go to the documentation of this file.
1// Copyright 2010-2025 Google LLC
2// Licensed under the Apache License, Version 2.0 (the "License");
3// you may not use this file except in compliance with the License.
4// You may obtain a copy of the License at
5//
6// http://www.apache.org/licenses/LICENSE-2.0
7//
8// Unless required by applicable law or agreed to in writing, software
9// distributed under the License is distributed on an "AS IS" BASIS,
10// WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
11// See the License for the specific language governing permissions and
12// limitations under the License.
13
14// An implementation of a push-relabel algorithm for the max flow problem.
15//
16// In the following, we consider a graph G = (V,E,s,t) where V denotes the set
17// of nodes (vertices) in the graph, E denotes the set of arcs (edges). s and t
18// denote distinguished nodes in G called source and target. n = |V| denotes the
19// number of nodes in the graph, and m = |E| denotes the number of arcs in the
20// graph.
21//
22// Each arc (v,w) is associated a capacity c(v,w).
23//
24// A flow is a function from E to R such that:
25//
26// a) f(v,w) <= c(v,w) for all (v,w) in E (capacity constraint.)
27//
28// b) f(v,w) = -f(w,v) for all (v,w) in E (flow antisymmetry constraint.)
29//
30// c) sum on v f(v,w) = 0 (flow conservation.)
31//
32// The goal of this algorithm is to find the maximum flow from s to t, i.e.
33// for example to maximize sum v f(s,v).
34//
35// The starting reference for this class of algorithms is:
36// A.V. Goldberg and R.E. Tarjan. A new approach to the maximum flow problem.
37// ACM Symposium on Theory of Computing, pp. 136-146.
38// http://portal.acm.org/citation.cfm?id=12144.
39//
40// The basic idea of the algorithm is to handle preflows instead of flows,
41// and to refine preflows until a maximum flow is obtained.
42// A preflow is like a flow, except that the inflow can be larger than the
43// outflow. If it is the case at a given node v, it is said that there is an
44// excess at node v, and inflow = outflow + excess.
45//
46// More formally, a preflow is a function f such that:
47//
48// 1) f(v,w) <= c(v,w) for all (v,w) in E (capacity constraint). c(v,w) is a
49// value representing the maximum capacity for arc (v,w).
50//
51// 2) f(v,w) = -f(w,v) for all (v,w) in E (flow antisymmetry constraint)
52//
53// 3) excess(v) = sum on u f(u,v) >= 0 is the excess at node v, the
54// algebraic sum of all the incoming preflows at this node.
55//
56// Each node has an associated "height", in addition to its excess. The
57// height of the source is defined to be equal to n, and cannot change. The
58// height of the target is defined to be zero, and cannot change either. The
59// height of all the other nodes is initialized at zero and is updated during
60// the algorithm (see below). For those who want to know the details, the height
61// of a node, corresponds to a reduced cost, and this enables one to prove that
62// the algorithm actually computes the max flow. Note that the height of a node
63// can be initialized to the distance to the target node in terms of number of
64// nodes. This has not been tried in this implementation.
65//
66// A node v is said to be *active* if excess(v) > 0.
67//
68// In this case the following operations can be applied to it:
69//
70// - if there are *admissible* incident arcs, i.e. arcs which are not saturated,
71// and whose head's height is lower than the height of the active node
72// considered, a PushFlow operation can be applied. It consists in sending as
73// much flow as both the excess at the node and the capacity of the arc
74// permit.
75// - if there are no admissible arcs, the active node considered is relabeled,
76// i.e. its height is increased to 1 + the minimum height of its neighboring
77// nodes on admissible arcs.
78// This is implemented in Discharge, which itself calls PushFlow and Relabel.
79//
80// Before running Discharge, it is necessary to initialize the algorithm with a
81// preflow. This is done in InitializePreflow, which saturates all the arcs
82// leaving the source node, and sets the excess at the heads of those arcs
83// accordingly.
84//
85// The algorithm terminates when there are no remaining active nodes, i.e. all
86// the excesses at all nodes are equal to zero. In this case, a maximum flow is
87// obtained.
88//
89// The complexity of this algorithm depends amongst other things on the choice
90// of the next active node. It has been shown, for example in:
91// L. Tuncel, "On the Complexity of Preflow-Push Algorithms for Maximum-Flow
92// Problems", Algorithmica 11(4): 353-359 (1994).
93// and
94// J. Cheriyan and K. Mehlhorn, "An analysis of the highest-level selection rule
95// in the preflow-push max-flow algorithm", Information processing letters,
96// 69(5):239-242 (1999).
97// http://www.math.uwaterloo.ca/~jcheriya/PS_files/me3.0.ps
98//
99// ...that choosing the active node with the highest level yields a
100// complexity of O(n^2 * sqrt(m)).
101//
102// TODO(user): implement the above active node choice rule.
103//
104// This has been validated experimentally in:
105// R.K. Ahuja, M. Kodialam, A.K. Mishra, and J.B. Orlin, "Computational
106// Investigations of Maximum Flow Algorithms", EJOR 97:509-542(1997).
107// http://jorlin.scripts.mit.edu/docs/publications/58-comput%20investigations%20of.pdf.
108//
109//
110// TODO(user): an alternative would be to evaluate:
111// A.V. Goldberg, "The Partial Augment-Relabel Algorithm for the Maximum Flow
112// Problem.” In Proceedings of Algorithms ESA, LNCS 5193:466-477, Springer 2008.
113// http://www.springerlink.com/index/5535k2j1mt646338.pdf
114//
115// An interesting general reference on network flows is:
116// R. K. Ahuja, T. L. Magnanti, J. B. Orlin, "Network Flows: Theory, Algorithms,
117// and Applications," Prentice Hall, 1993, ISBN: 978-0136175490,
118// http://www.amazon.com/dp/013617549X
119//
120// Keywords: Push-relabel, max-flow, network, graph, Goldberg, Tarjan, Dinic,
121// Dinitz.
122
123#ifndef OR_TOOLS_GRAPH_GENERIC_MAX_FLOW_H_
124#define OR_TOOLS_GRAPH_GENERIC_MAX_FLOW_H_
125
126#include <algorithm>
127#include <cstdint>
128#include <limits>
129#include <memory>
130#include <string>
131#include <utility>
132#include <vector>
133
134#include "absl/strings/string_view.h"
135#include "ortools/base/logging.h"
137#include "ortools/graph/flow_problem.pb.h"
138#include "ortools/util/stats.h"
139#include "ortools/util/zvector.h"
140
141namespace operations_research {
142
143// Specific but efficient priority queue implementation. The priority type must
144// be an integer. The queue allows to retrieve the element with highest priority
145// but only allows pushes with a priority greater or equal to the highest
146// priority in the queue minus one. All operations are in O(1) and the memory is
147// in O(num elements in the queue). Elements with the same priority are
148// retrieved with LIFO order.
149//
150// Note(user): As far as I know, this is an original idea and is the only code
151// that use this in the Maximum Flow context. Papers usually refer to an
152// height-indexed array of simple linked lists of active node with the same
153// height. Even worse, sometimes they use double-linked list to allow arbitrary
154// height update in order to detect missing height (used for the Gap heuristic).
155// But this can actually be implemented a lot more efficiently by just
156// maintaining the height distribution of all the node in the graph.
157template <typename Element, typename IntegerPriority>
159 public:
160 PriorityQueueWithRestrictedPush() : even_queue_(), odd_queue_() {}
161
162#ifndef SWIG
163 // This type is neither copyable nor movable.
165 delete;
167 const PriorityQueueWithRestrictedPush&) = delete;
168#endif
169
170 // Is the queue empty?
171 bool IsEmpty() const;
172
173 // Clears the queue.
174 void Clear();
175
176 // Push a new element in the queue. Its priority must be greater or equal to
177 // the highest priority present in the queue, minus one. This condition is
178 // DCHECKed, and violating it yields erroneous queue behavior in NDEBUG mode.
179 void Push(Element element, IntegerPriority priority);
180
181 // Returns the element with highest priority and remove it from the queue.
182 // IsEmpty() must be false, this condition is DCHECKed.
183 Element Pop();
184
185 private:
186 // Helper function to get the last element of a vector and pop it.
187 Element PopBack(std::vector<std::pair<Element, IntegerPriority> >* queue);
188
189 // This is the heart of the algorithm. basically we split the elements by
190 // parity of their priority and the precondition on the Push() ensures that
191 // both vectors are always sorted by increasing priority.
192 std::vector<std::pair<Element, IntegerPriority> > even_queue_;
193 std::vector<std::pair<Element, IntegerPriority> > odd_queue_;
194};
195
196// We want an enum for the Status of a max flow run, and we want this
197// enum to be scoped under GenericMaxFlow<>. Unfortunately, swig
198// doesn't handle templated enums very well, so we need a base,
199// untemplated class to hold it.
201 public:
202 enum Status {
203 NOT_SOLVED, // The problem was not solved, or its data were edited.
204 OPTIMAL, // Solve() was called and found an optimal solution.
205 INT_OVERFLOW, // There is a feasible flow > max possible flow.
206
207 // TODO(user): These are no longer used. Remove.
210 };
211};
212
213// Generic MaxFlow that works on graphs with the notion of "reverse" arcs.
214// That is, the graph is directed, and each arc 'tail -> head' must be
215// associated to an unique reverse arc going in the opposite direction
216// 'head -> tail'. We must also have reverse[reverse[arc]] = arc.
217//
218// This works with all the reverse arc graphs from 'util/graph/graph.h' and
219// uses the API defined there.
220//
221// We actually support two kind of graphs with "reverse" arcs depending on the
222// value of Graph::kHasNegativeReverseArcs:
223// - If it is true, the arcs from the "user graph" are called "direct" arcs and
224// are assumed to be indexed in [0, num_arcs). These are the only ones that
225// can have a positive capacity. All the "reverse/opposite" of these arcs
226// will always have negative indices in [-num_arcs, 0) and a capacity of
227// zero.
228// - If it is false, then we only have "direct" arcs in [0, num_arcs), the
229// reverse of each arc must lives in the same space, and both an arc and its
230// reverse can have a positive initial capacity. This can lead to twice fewer
231// arcs and a faster algo if the "user graph" as a lot of reverse arcs
232// already.
233//
234// Flow types and overflow: To optimize memory usage and speed we have two
235// integer types representing flow quantities.
236// - ArcFlowType (can be unsigned) is used to represent the flow per arc. The
237// type should be big enough to hold the initial capacity of an arc + its
238// reverse. So if all your initial capacity <=
239// std::numeric_limit<ArcFlowType>::max() / 2, you are safe.
240// - FlowSumType (must be signed) is used to hold the sum of a flow going
241// through a node or through the full network. Ideally, the sum of all the
242// capacity of the arcs leaving the source should fit in a FlowSumType. If
243// this is the case, you will get the better performance, and always an
244// OPTIMAL result. Otherwise, the code might be slower, and you will get
245// OPTIMAL only if the max_flow that can be sent through the network fit a
246// FlowSumType, you will get the INT_OVERFLOW status if more flow than that
247// can be sent.
248template <typename Graph, typename ArcFlowT = int64_t,
249 typename FlowSumT = int64_t>
251 public:
252 typedef typename Graph::NodeIndex NodeIndex;
253 typedef typename Graph::ArcIndex ArcIndex;
254 typedef ArcFlowT ArcFlowType;
255 typedef FlowSumT FlowSumType;
256 static_assert(sizeof(FlowSumType) >= sizeof(ArcFlowType));
257 static_assert(std::is_signed_v<FlowSumType>);
258
259 // The height of a node never excess 2 times the number of node, so we
260 // use the same type as a Node index.
262
263 // Initialize a MaxFlow instance on the given graph. The graph does not need
264 // to be fully built yet, but its capacity reservation are used to initialize
265 // the memory of this class. source and sink must also be valid node of
266 // graph.
268
269#ifndef SWIG
270 // This type is neither copyable nor movable.
273#endif
274
275 // Returns the graph associated to the current object.
276 const Graph* graph() const { return graph_; }
277
278 // Returns the status of last call to Solve(). NOT_SOLVED is returned if
279 // Solve() has never been called or if the problem has been modified in such a
280 // way that the previous solution becomes invalid.
281 Status status() const { return status_; }
282
283 // Returns the index of the node corresponding to the source of the network.
285
286 // Returns the index of the node corresponding to the sink of the network.
287 NodeIndex GetSinkNodeIndex() const { return sink_; }
288
289 // Sets the capacity for arc to new_capacity.
290 //
291 // Note that this will be ignored for self-arc, so do not be surprised if you
292 // get zero when reading the capacity of a self-arc back.
293 void SetArcCapacity(ArcIndex arc, ArcFlowType new_capacity);
294
295 // Returns true if a maximum flow was solved.
296 bool Solve();
297
298 // Returns the total flow found by the algorithm.
300
301 // Returns the current flow on the given arc.
302 //
303 // Note that on a couple (arc, opposite_arc) the flow only goes in one
304 // direction (where it is positive) and the other direction will have the
305 // negation of that flow.
306 //
307 // Tricky: This returns a FlowSumType in order to be able to distinguish
308 // negative flow from positive flow when ArcFlowType is unsigned.
309 // TODO(user): We could change the API to avoid that if needed.
311 if constexpr (Graph::kHasNegativeReverseArcs) {
312 if (IsArcDirect(arc)) {
313 return static_cast<FlowSumType>(residual_arc_capacity_[Opposite(arc)]);
314 }
315 return -static_cast<FlowSumType>(residual_arc_capacity_[arc]);
316 }
317 return static_cast<FlowSumType>(initial_capacity_[arc]) -
318 static_cast<FlowSumType>(residual_arc_capacity_[arc]);
319 }
320
321 // Returns the initial capacity of an arc.
323 if constexpr (Graph::kHasNegativeReverseArcs) {
324 if (!IsArcDirect(arc)) return 0;
325 return residual_arc_capacity_[arc] +
327 }
328 return initial_capacity_[arc];
329 }
330
331 // Returns the nodes reachable from the source in the residual graph, the
332 // outgoing arcs of this set form a minimum cut.
333 void GetSourceSideMinCut(std::vector<NodeIndex>* result);
334
335 // Returns the nodes that can reach the sink in the residual graph, the
336 // outgoing arcs of this set form a minimum cut. Note that if this is the
337 // complement of GetNodeReachableFromSource(), then the min-cut is unique.
338 //
339 // TODO(user): In the two-phases algorithm, we can get this minimum cut
340 // without doing the second phase. Add an option for this if there is a need
341 // to, note that the second phase is pretty fast so the gain will be small.
342 void GetSinkSideMinCut(std::vector<NodeIndex>* result);
343
344 // Returns true if there exists a path from the source to the sink with
345 // remaining capacity. This allows us to easily check at the end that the flow
346 // we computed is indeed optimal (provided that all the conditions tested by
347 // CheckResult() also hold).
349
350 // Returns the protocol buffer representation of the current problem.
351 FlowModelProto CreateFlowModel();
352
353 protected:
354 // Checks whether the result is valid, i.e. that node excesses are all equal
355 // to zero (we have a flow) and that residual capacities are all non-negative
356 // or zero.
357 bool CheckResult() const;
358
359 // Returns true if arc is admissible.
361 const NodeHeight* potentials) const {
362 DCHECK_EQ(tail, Tail(arc));
363 return residual_arc_capacity_[arc] > 0 &&
364 potentials[tail] == potentials[Head(arc)] + 1;
365 }
366
367 // Returns true if node is active, i.e. if its excess is positive and it
368 // is neither the source or the sink of the graph.
369 bool IsActive(NodeIndex node) const {
370 return (node != source_) && (node != sink_) && (node_excess_[node] > 0);
371 }
372
373 // Returns true if a precondition for Relabel is met, i.e. the outgoing arcs
374 // of node are all either saturated or the heights of their heads are greater
375 // or equal to the height of node.
377
378 // Returns context concatenated with information about arc
379 // in a human-friendly way.
380 std::string DebugString(absl::string_view context, ArcIndex arc) const;
381
382 // Initializes the container active_nodes_.
384
385 // Get the first element from the active node container.
389
390 // Push element to the active node container.
391 void PushActiveNode(const NodeIndex& node) {
392 active_node_by_height_.Push(node, node_potential_[node]);
393 }
394
395 // Check the emptiness of the container.
397
398 // Performs optimization step.
399 void Refine();
401
402 // Discharges an active node node by saturating its admissible adjacent arcs,
403 // if any, and by relabelling it when it becomes inactive.
405
406 // Initializes the preflow to a state that enables to run Refine.
408
409 // Clears the flow excess at each node by pushing the flow back to the source:
410 // - Do a depth-first search from the source in the direct graph to cancel
411 // flow cycles.
412 // - Then, return flow excess along the depth-first search tree (by pushing
413 // the flow in the reverse dfs topological order).
414 // The theoretical complexity is O(mn), but it is a lot faster in practice.
416
417 // Computes the best possible node potential given the current flow using a
418 // reverse breadth-first search from the sink in the reverse residual graph.
419 // This is an implementation of the global update heuristic mentioned in many
420 // max-flow papers. See for instance: B.V. Cherkassky, A.V. Goldberg, "On
421 // implementing push-relabel methods for the maximum flow problem",
422 // Algorithmica, 19:390-410, 1997.
423 // ftp://reports.stanford.edu/pub/cstr/reports/cs/tr/94/1523/CS-TR-94-1523.pdf
425
426 // Tries to saturate all the outgoing arcs from the source that can reach the
427 // sink. Most of the time, we can do that in one go, except when more flow
428 // than kMaxFlowSum can be pushed out of the source in which case we have
429 // to be careful. Returns true if some flow was pushed.
431
432 // Pushes flow on arc, i.e. consumes flow on residual_arc_capacity_[arc],
433 // and consumes -flow on residual_arc_capacity_[Opposite(arc)]. Updates
434 // node_excess_ at the tail and head of arc accordingly.
435 void PushFlow(ArcFlowType flow, NodeIndex tail, ArcIndex arc);
436
437 // Similar to PushFlow(), but this will always push
438 // std::min(node_excess_[tail], residual_arc_capacity_[arc]).
439 //
440 // We also use a cached node_excess_.data() as this is used in hot loops and
441 // allow to remove bound checking and fetching the pointer which is constant.
443 FlowSumType* node_excess);
444
445 // Relabels a node, i.e. increases its height by the minimum necessary amount.
446 // This version of Relabel is relaxed in a way such that if an admissible arc
447 // exists at the current node height, then the node is not relabeled. This
448 // enables us to deal with wrong values of first_admissible_arc_[node] when
449 // updating it is too costly.
450 void Relabel(NodeIndex node);
451
452 // Handy member functions to make the code more compact.
453 NodeIndex Head(ArcIndex arc) const { return graph_->Head(arc); }
454 NodeIndex Tail(ArcIndex arc) const { return graph_->Tail(arc); }
456 bool IsArcDirect(ArcIndex arc) const;
457 bool IsArcValid(ArcIndex arc) const;
458
459 // Returns the set of nodes reachable from start in the residual graph or in
460 // the reverse residual graph (if reverse is true).
461 template <bool reverse>
462 void ComputeReachableNodes(NodeIndex start, std::vector<NodeIndex>* result);
463
464 // Maximum manageable flow.
465 static constexpr FlowSumType kMaxFlowSum =
466 std::numeric_limits<FlowSumType>::max();
467
468 // A pointer to the graph passed as argument.
469 const Graph* graph_;
470
471 // An array representing the excess for each node in graph_.
472 // With the current algorithm this is always positive except at the source.
473 // TODO(user): If needed we could tweak that to allow unsigned FlowSumType.
474 std::unique_ptr<FlowSumType[]> node_excess_;
475
476 // An array representing the height function for each node in graph_. For a
477 // given node, this is a lower bound on the shortest path length from this
478 // node to the sink in the residual network. The height of a node always goes
479 // up during the course of a Solve().
480 //
481 // Since initially we saturate all the outgoing arcs of the source, we can
482 // never reach the sink from the source in the residual graph. Initially we
483 // set the height of the source to n (the number of node of the graph) and it
484 // never changes. If a node as an height >= n, then this node can't reach the
485 // sink and its height minus n is a lower bound on the shortest path length
486 // from this node to the source in the residual graph.
487 std::unique_ptr<NodeHeight[]> node_potential_;
488
489 // An array representing the residual_capacity for each arc in graph_.
490 // Residual capacities enable one to represent the capacity and flow for all
491 // arcs in the graph in the following manner.
492 // For all arc, residual_arc_capacity_[arc] = capacity[arc] - flow[arc]
493 // Moreover, for reverse arcs, capacity[arc] = 0 by definition,
494 // Also flow[Opposite(arc)] = -flow[arc] by definition.
495 // Therefore:
496 // - for a direct arc:
497 // flow[arc] = 0 - flow[Opposite(arc)]
498 // = capacity[Opposite(arc)] - flow[Opposite(arc)]
499 // = residual_arc_capacity_[Opposite(arc)]
500 // - for a reverse arc:
501 // flow[arc] = -residual_arc_capacity_[arc]
502 // Using these facts enables one to only maintain residual_arc_capacity_,
503 // instead of both capacity and flow, for each direct and indirect arc. This
504 // reduces the amount of memory for this information by a factor 2.
506
507 // The initial capacity as set by SetArcCapacity(), unused if
508 // `Graph::kNegativeReverseArcs`, as we can always recover the initial
509 // capacity from the residual capacities of an arc and its reverse.
510 std::unique_ptr<ArcFlowType[]> initial_capacity_;
511
512 // An array representing the first admissible arc for each node in graph_.
513 std::unique_ptr<ArcIndex[]> first_admissible_arc_;
514
515 // A priority queue used for managing active nodes in the algorithm. It allows
516 // to select the active node with highest height before each Discharge().
517 // Moreover, since all pushes from this node will be to nodes with height
518 // greater or equal to the initial discharged node height minus one, the
519 // PriorityQueueWithRestrictedPush is a perfect fit.
521
522 // The index of the source node in graph_.
524
525 // The index of the sink node in graph_.
527
528 // The status of the problem.
530
531 // BFS queue used by the GlobalUpdate() function. We do not use a C++ queue
532 // because we need access to the vector for different optimizations.
533 std::vector<bool> node_in_bfs_queue_;
534 std::vector<NodeIndex> bfs_queue_;
535
536 // Statistics about this class.
538};
539
540template <typename Element, typename IntegerPriority>
542 const {
543 return even_queue_.empty() && odd_queue_.empty();
544}
545
546template <typename Element, typename IntegerPriority>
548 even_queue_.clear();
549 odd_queue_.clear();
550}
551
552template <typename Element, typename IntegerPriority>
554 Element element, IntegerPriority priority) {
555 // Since users may rely on it, we DCHECK the exact condition.
556 DCHECK(even_queue_.empty() || priority >= even_queue_.back().second - 1);
557 DCHECK(odd_queue_.empty() || priority >= odd_queue_.back().second - 1);
558
559 // Note that the DCHECK() below are less restrictive than the ones above but
560 // check a necessary and sufficient condition for the priority queue to behave
561 // as expected.
562 if (priority & 1) {
563 DCHECK(odd_queue_.empty() || priority >= odd_queue_.back().second);
564 odd_queue_.push_back(std::make_pair(element, priority));
565 } else {
566 DCHECK(even_queue_.empty() || priority >= even_queue_.back().second);
567 even_queue_.push_back(std::make_pair(element, priority));
568 }
569}
570
571template <typename Element, typename IntegerPriority>
573 DCHECK(!IsEmpty());
574 if (even_queue_.empty()) return PopBack(&odd_queue_);
575 if (odd_queue_.empty()) return PopBack(&even_queue_);
576 if (odd_queue_.back().second > even_queue_.back().second) {
577 return PopBack(&odd_queue_);
578 } else {
579 return PopBack(&even_queue_);
580 }
581}
582
583template <typename Element, typename IntegerPriority>
584Element PriorityQueueWithRestrictedPush<Element, IntegerPriority>::PopBack(
585 std::vector<std::pair<Element, IntegerPriority> >* queue) {
586 DCHECK(!queue->empty());
587 Element element = queue->back().first;
588 queue->pop_back();
589 return element;
590}
591
592template <typename Graph, typename ArcFlowT, typename FlowSumT>
594 NodeIndex source,
595 NodeIndex sink)
596 : graph_(graph),
598 source_(source),
599 sink_(sink),
600 stats_("MaxFlow") {
602 DCHECK(graph->IsNodeValid(source));
603 DCHECK(graph->IsNodeValid(sink));
604 const NodeIndex max_num_nodes = graph_->node_capacity();
605 if (max_num_nodes > 0) {
606 // We will initialize them in InitializePreflow(), so no need for memset.
607 //
608 // TODO(user): Unfortunately std/absl::make_unique_for_overwrite is not
609 // yet available in open-source.
610 node_excess_ = std::make_unique<FlowSumType[]>(max_num_nodes);
611 node_potential_ = std::make_unique<NodeHeight[]>(max_num_nodes);
612 first_admissible_arc_ = std::make_unique<ArcIndex[]>(max_num_nodes);
613 bfs_queue_.reserve(max_num_nodes);
614 }
615 const ArcIndex max_num_arcs = graph_->arc_capacity();
616 if (max_num_arcs > 0) {
617 if constexpr (Graph::kHasNegativeReverseArcs) {
618 residual_arc_capacity_.Reserve(-max_num_arcs, max_num_arcs - 1);
619 } else {
620 // We will need to store the initial capacity in this case.
621 initial_capacity_ = std::make_unique<ArcFlowType[]>(max_num_arcs);
622 residual_arc_capacity_.Reserve(0, max_num_arcs - 1);
623 }
624 residual_arc_capacity_.SetAll(0);
625 }
626}
627
628template <typename Graph, typename ArcFlowT, typename FlowSumT>
630 ArcIndex arc, ArcFlowType new_capacity) {
632 DCHECK_LE(0, new_capacity);
633 DCHECK(IsArcDirect(arc));
634
635 // This serves no purpose from a max-flow point of view, so it is safer to
636 // just leave the capacity of all self-arc at zero.
637 if (Head(arc) == Tail(arc)) return;
638
640 residual_arc_capacity_[arc] = new_capacity;
641
642 // Since the class might have already be used, we need to make sure we clear
643 // any previous state on this arc and its reverse.
644 if constexpr (Graph::kHasNegativeReverseArcs) {
646 } else {
647 initial_capacity_[arc] = new_capacity;
648 DCHECK_LE(initial_capacity_[arc], std::numeric_limits<ArcFlowType>::max() -
650 }
651}
652
653template <typename Graph, typename ArcFlowT, typename FlowSumT>
658
659template <typename Graph, typename ArcFlowT, typename FlowSumT>
661 std::vector<NodeIndex>* result) {
663}
664
665template <typename Graph, typename ArcFlowT, typename FlowSumT>
669 LOG(DFATAL) << "-node_excess_[source_] = " << -node_excess_[source_]
670 << " != node_excess_[sink_] = " << node_excess_[sink_];
671 return false;
672 }
673 for (NodeIndex node = 0; node < graph_->num_nodes(); ++node) {
674 if (node != source_ && node != sink_) {
675 if (node_excess_[node] != 0) {
676 LOG(DFATAL) << "node_excess_[" << node << "] = " << node_excess_[node]
677 << " != 0";
678 return false;
679 }
680 }
681 }
682 for (ArcIndex arc = 0; arc < graph_->num_arcs(); ++arc) {
683 const ArcIndex opposite = Opposite(arc);
684 const ArcFlowType direct_capacity = residual_arc_capacity_[arc];
685 const ArcFlowType opposite_capacity = residual_arc_capacity_[opposite];
686 if (direct_capacity < 0) {
687 LOG(DFATAL) << "residual_arc_capacity_[" << arc
688 << "] = " << direct_capacity << " < 0";
689 return false;
690 }
691 if (opposite_capacity < 0) {
692 LOG(DFATAL) << "residual_arc_capacity_[" << opposite
693 << "] = " << opposite_capacity << " < 0";
694 return false;
695 }
696 // The initial capacity of the direct arcs is non-negative.
697 if (direct_capacity + opposite_capacity < 0) {
698 LOG(DFATAL) << "initial capacity [" << arc
699 << "] = " << direct_capacity + opposite_capacity << " < 0";
700 return false;
701 }
702 }
703
705 LOG(DFATAL) << "The algorithm terminated, but the flow is not maximal!";
706 return false;
707 }
708
709 return true;
710}
711
712template <typename Graph, typename ArcFlowT, typename FlowSumT>
715
716 // We simply compute the reachability from the source in the residual graph.
717 const NodeIndex num_nodes = graph_->num_nodes();
718 std::vector<bool> is_reached(num_nodes, false);
719 std::vector<NodeIndex> to_process;
720
721 to_process.push_back(source_);
722 is_reached[source_] = true;
723 while (!to_process.empty()) {
724 const NodeIndex node = to_process.back();
725 to_process.pop_back();
726 for (const ArcIndex arc : graph_->OutgoingOrOppositeIncomingArcs(node)) {
727 if (residual_arc_capacity_[arc] > 0) {
728 const NodeIndex head = graph_->Head(arc);
729 if (!is_reached[head]) {
730 is_reached[head] = true;
731 to_process.push_back(head);
732 }
733 }
734 }
735 }
736 return is_reached[sink_];
737}
738
739template <typename Graph, typename ArcFlowT, typename FlowSumT>
741 NodeIndex node) const {
742 DCHECK(IsActive(node));
743 for (const ArcIndex arc : graph_->OutgoingOrOppositeIncomingArcs(node)) {
744 DCHECK(!IsAdmissible(node, arc, node_potential_.data()))
745 << DebugString("CheckRelabelPrecondition:", arc);
746 }
747 return true;
748}
749
750template <typename Graph, typename ArcFlowT, typename FlowSumT>
752 absl::string_view context, ArcIndex arc) const {
753 const NodeIndex tail = Tail(arc);
754 const NodeIndex head = Head(arc);
755 return absl::StrFormat(
756 "%s Arc %d, from %d to %d, "
757 "Residual capacity = %d, "
758 "Residual capacity for reverse arc = %d, "
759 "Height(tail) = %d, Height(head) = %d, "
760 "Excess(tail) = %d, Excess(head) = %d",
761 context, arc, tail, head, residual_arc_capacity_[arc],
763 node_potential_[head], node_excess_[tail], node_excess_[head]);
764}
765
766template <typename Graph, typename ArcFlowT, typename FlowSumT>
770
771 // Deal with the case when source_ or sink_ is not inside graph_.
772 // Since they are both specified independently of the graph, we do need to
773 // take care of this corner case.
774 const NodeIndex num_nodes = graph_->num_nodes();
775 if (sink_ >= num_nodes || source_ >= num_nodes) {
776 // Behave like a normal graph where source_ and sink_ are disconnected.
777 // Note that the arc flow is set to 0 by InitializePreflow().
779 return true;
780 }
781
783
785 DCHECK(CheckResult());
786
788 // In this case, we are sure that the flow is > kMaxFlowSum.
790 }
791 IF_STATS_ENABLED(VLOG(1) << stats_.StatString());
792 return true;
793}
794
795template <typename Graph, typename ArcFlowT, typename FlowSumT>
798
799 // TODO(user): Ebert graph has an issue with nodes with no arcs, so we
800 // use max_num_nodes here to resize vectors.
801 const NodeIndex num_nodes = graph_->num_nodes();
802 const NodeIndex max_num_nodes = graph_->node_capacity();
803
804 // InitializePreflow() clears the whole flow that could have been computed
805 // by a previous Solve(). This is not optimal in terms of complexity.
806 //
807 // TODO(user): find a way to make the re-solving incremental (not an obvious
808 // task, and there has not been a lot of literature on the subject.)
809 std::fill(node_excess_.get(), node_excess_.get() + max_num_nodes, 0);
810
811 // Restart from a clear state with no flow, and an initial arc capacity.
812 const ArcIndex num_arcs = graph_->num_arcs();
813 if constexpr (Graph::kHasNegativeReverseArcs) {
814 for (ArcIndex arc = 0; arc < num_arcs; ++arc) {
815 const ArcIndex opposite_arc = Opposite(arc);
817 residual_arc_capacity_[opposite_arc] = 0;
818 }
819 } else {
820 for (ArcIndex arc = 0; arc < num_arcs; ++arc) {
822 }
823 }
824
825 // All the initial heights are zero except for the source whose height is
826 // equal to the number of nodes and will never change during the algorithm.
827 std::fill(node_potential_.get(), node_potential_.get() + max_num_nodes, 0);
828 node_potential_[source_] = num_nodes;
829
830 // Initially we set first_admissible_arc_ to the first arc in the iteration.
831 for (NodeIndex node = 0; node < num_nodes; ++node) {
833 for (const ArcIndex arc : graph_->OutgoingOrOppositeIncomingArcs(node)) {
834 first_admissible_arc_[node] = arc;
835 break;
836 }
837 }
838}
839
840// Note(user): Calling this function will break the property on the node
841// potentials because of the way we cancel flow on cycle. However, we only call
842// that at the end of the algorithm, or just before a GlobalUpdate() that will
843// restore the precondition on the node potentials.
844template <typename Graph, typename ArcFlowT, typename FlowSumT>
847 const NodeIndex num_nodes = graph_->num_nodes();
848
849 // We implement a variation of Tarjan's strongly connected component algorithm
850 // to detect cycles published in: Tarjan, R. E. (1972), "Depth-first search
851 // and linear graph algorithms", SIAM Journal on Computing. A description can
852 // also be found in wikipedia.
853
854 // Stored nodes are settled nodes already stored in the
855 // reverse_topological_order (except the sink_ that we do not actually store).
856 std::vector<bool> stored(num_nodes, false);
857 stored[sink_] = true;
858
859 // The visited nodes that are not yet stored are all the nodes from the
860 // source_ to the current node in the current dfs branch.
861 std::vector<bool> visited(num_nodes, false);
862 visited[sink_] = true;
863
864 // Stack of arcs to explore in the dfs search.
865 // The current node is Head(arc_stack.back()).
866 std::vector<ArcIndex> arc_stack;
867
868 // Increasing list of indices into the arc_stack that correspond to the list
869 // of arcs in the current dfs branch from the source_ to the current node.
870 std::vector<int> index_branch;
871
872 // Node in reverse_topological_order in the final dfs tree.
873 std::vector<NodeIndex> reverse_topological_order;
874
875 // We start by pushing all the outgoing arcs from the source on the stack to
876 // avoid special conditions in the code. As a result, source_ will not be
877 // stored in reverse_topological_order, and this is what we want.
878 for (const ArcIndex arc : graph_->OutgoingArcs(source_)) {
879 if (Flow(arc) > 0) {
880 arc_stack.push_back(arc);
881 }
882 }
883 visited[source_] = true;
884
885 // Start the dfs on the subgraph formed by the direct arcs with positive flow.
886 while (!arc_stack.empty()) {
887 DCHECK_GT(Flow(arc_stack.back()), 0)
888 << arc_stack.size() - 1 << " arc " << arc_stack.back() << " "
889 << Tail(arc_stack.back()) << "->" << Head(arc_stack.back());
890 const NodeIndex node = Head(arc_stack.back());
891
892 // If the node is visited, it means we have explored all its arcs and we
893 // have just backtracked in the dfs. Store it if it is not already stored
894 // and process the next arc on the stack.
895 if (visited[node]) {
896 if (!stored[node]) {
897 stored[node] = true;
898 reverse_topological_order.push_back(node);
899 DCHECK(!index_branch.empty());
900 index_branch.pop_back();
901 }
902 arc_stack.pop_back();
903 continue;
904 }
905
906 // The node is a new unexplored node, add all its outgoing arcs with
907 // positive flow to the stack and go deeper in the dfs.
908 DCHECK(!stored[node]);
909 DCHECK(index_branch.empty() ||
910 (arc_stack.size() - 1 > index_branch.back()));
911 visited[node] = true;
912 index_branch.push_back(arc_stack.size() - 1);
913 for (const ArcIndex arc : graph_->OutgoingArcs(node)) {
914 const FlowSumType flow = Flow(arc);
915 const NodeIndex head = Head(arc);
916 if (flow > 0 && !stored[head]) {
917 if (!visited[head]) {
918 arc_stack.push_back(arc);
919 } else {
920 // There is a cycle.
921 // Find the first index to consider,
922 // arc_stack[index_branch[cycle_begin]] will be the first arc on the
923 // cycle.
924 int cycle_begin = index_branch.size();
925 while (cycle_begin > 0 &&
926 Head(arc_stack[index_branch[cycle_begin - 1]]) != head) {
927 --cycle_begin;
928 }
929
930 // Compute the maximum flow that can be canceled on the cycle and the
931 // min index such that arc_stack[index_branch[i]] will be saturated.
932 ArcFlowType flow_on_cycle = static_cast<ArcFlowType>(flow);
933 int first_saturated_index = index_branch.size();
934 for (int i = index_branch.size() - 1; i >= cycle_begin; --i) {
935 const ArcIndex arc_on_cycle = arc_stack[index_branch[i]];
936 if (const ArcFlowType arc_flow = Flow(arc_on_cycle);
937 arc_flow <= flow_on_cycle) {
938 flow_on_cycle = arc_flow;
939 first_saturated_index = i;
940 }
941 }
942
943 // This is just here for a DCHECK() below.
944 const FlowSumType excess = node_excess_[head];
945
946 // Cancel the flow on the cycle, and set visited[node] = false for
947 // the node that will be backtracked over.
948 PushFlow(-flow_on_cycle, node, arc);
949 for (int i = index_branch.size() - 1; i >= cycle_begin; --i) {
950 const ArcIndex arc_on_cycle = arc_stack[index_branch[i]];
951 PushFlow(-flow_on_cycle, Tail(arc_on_cycle), arc_on_cycle);
952 if (i >= first_saturated_index) {
953 DCHECK(visited[Head(arc_on_cycle)]);
954 visited[Head(arc_on_cycle)] = false;
955 } else {
956 DCHECK_GT(Flow(arc_on_cycle), 0);
957 }
958 }
959
960 // This is a simple check that the flow was pushed properly.
961 DCHECK_EQ(excess, node_excess_[head]);
962
963 // Backtrack the dfs just before index_branch[first_saturated_index].
964 // If the current node is still active, there is nothing to do.
965 if (first_saturated_index < index_branch.size()) {
966 arc_stack.resize(index_branch[first_saturated_index]);
967 index_branch.resize(first_saturated_index);
968
969 // We backtracked over the current node, so there is no need to
970 // continue looping over its arcs.
971 break;
972 }
973 }
974 }
975 }
976 }
977 DCHECK(arc_stack.empty());
978 DCHECK(index_branch.empty());
979
980 // Return the flow to the source. Note that the sink_ and the source_ are not
981 // stored in reverse_topological_order.
982 for (int i = 0; i < reverse_topological_order.size(); i++) {
983 const NodeIndex node = reverse_topological_order[i];
984 if (node_excess_[node] == 0) continue;
985 for (const ArcIndex arc : graph_->OutgoingOrOppositeIncomingArcs(node)) {
986 const FlowSumType flow = Flow(arc);
987 if (flow < 0) {
988 DCHECK_GT(residual_arc_capacity_[arc], 0);
989 const ArcFlowType to_push =
990 static_cast<ArcFlowType>(std::min(node_excess_[node], -flow));
991 PushFlow(to_push, node, arc);
992 if (node_excess_[node] == 0) break;
993 }
994 }
995 DCHECK_EQ(0, node_excess_[node]);
996 }
997 DCHECK_EQ(-node_excess_[source_], node_excess_[sink_]);
998}
999
1000template <typename Graph, typename ArcFlowT, typename FlowSumT>
1003 bfs_queue_.clear();
1004 int queue_index = 0;
1005 const NodeIndex num_nodes = graph_->num_nodes();
1006 node_in_bfs_queue_.assign(num_nodes, false);
1007 node_in_bfs_queue_[sink_] = true;
1008
1009 // We cache this as this is a hot-loop and we don't want bound checking.
1010 FlowSumType* const node_excess = node_excess_.get();
1011
1012 // We do a BFS in the reverse residual graph, starting from the sink.
1013 // Because all the arcs from the source are saturated (except in
1014 // presence of integer overflow), the source cannot reach the sink in the
1015 // residual graph. If there is possible overflow and the source is reachable,
1016 // we still do not want to relabel it, so we start with the source marked.
1018 bfs_queue_.push_back(sink_);
1019
1020 while (queue_index != bfs_queue_.size()) {
1021 const NodeIndex node = bfs_queue_[queue_index];
1022 ++queue_index;
1023 const NodeIndex candidate_distance = node_potential_[node] + 1;
1024 for (const ArcIndex arc : graph_->OutgoingOrOppositeIncomingArcs(node)) {
1025 const NodeIndex head = Head(arc);
1026
1027 // Skip the arc if the height of head was already set to the correct
1028 // value (Remember we are doing reverse BFS).
1029 if (node_in_bfs_queue_[head]) continue;
1030
1031 // TODO(user): By using more memory we can speed this up quite a bit by
1032 // avoiding to take the opposite arc here, too options:
1033 // - if (residual_arc_capacity_[arc] != arc_capacity_[arc])
1034 // - if (opposite_arc_is_admissible_[arc]) // need updates.
1035 // Experiment with the first option shows more than 10% gain on this
1036 // function running time, which is the bottleneck on many instances.
1037 const ArcIndex opposite_arc = Opposite(arc);
1038 if (residual_arc_capacity_[opposite_arc] > 0) {
1039 // Note(user): We used to have a DCHECK_GE(candidate_distance,
1040 // node_potential_[head]); which is always true except in the case
1041 // where we can push more than kMaxFlowSum out of the source. The
1042 // problem comes from the fact that in this case, we call
1043 // PushFlowExcessBackToSource() in the middle of the algorithm. The
1044 // later call will break the properties of the node potential. Note
1045 // however, that this function will recompute a good node potential
1046 // for all the nodes and thus fix the issue.
1047
1048 // If head is active, we can steal some or all of its excess.
1049 // This brings a huge gain on some problems.
1050 // Note(user): I haven't seen this anywhere in the literature.
1051 // TODO(user): Investigate more and maybe write a publication :)
1052 if (node_excess[head] > 0) {
1053 PushAsMuchFlowAsPossible(head, opposite_arc, node_excess);
1054
1055 // If the arc became saturated, it is no longer in the residual
1056 // graph, so we do not need to consider head at this time.
1057 if (residual_arc_capacity_[opposite_arc] == 0) continue;
1058 }
1059
1060 // Note that there is no need to touch first_admissible_arc_[node]
1061 // because of the relaxed Relabel() we use.
1062 node_potential_[head] = candidate_distance;
1063 node_in_bfs_queue_[head] = true;
1064 bfs_queue_.push_back(head);
1065 }
1066 }
1067 }
1068
1069 // At the end of the search, some nodes may not be in the bfs_queue_. Such
1070 // nodes cannot reach the sink_ or source_ in the residual graph, so there is
1071 // no point trying to push flow toward them. We obtain this effect by setting
1072 // their height to something unreachable.
1073 //
1074 // Note that this also prevents cycling due to our anti-overflow procedure.
1075 // For instance, suppose there is an edge s -> n outgoing from the source. If
1076 // node n has no other connection and some excess, we will push the flow back
1077 // to the source, but if we don't update the height of n
1078 // SaturateOutgoingArcsFromSource() will push the flow to n again.
1079 // TODO(user): This is another argument for another anti-overflow algorithm.
1080 for (NodeIndex node = 0; node < num_nodes; ++node) {
1081 if (!node_in_bfs_queue_[node]) {
1082 node_potential_[node] = 2 * num_nodes - 1;
1083 }
1084 }
1085
1086 // Reset the active nodes. Doing it like this pushes the nodes in increasing
1087 // order of height. Note that bfs_queue_[0] is the sink_ so we skip it.
1089 for (int i = 1; i < bfs_queue_.size(); ++i) {
1090 const NodeIndex node = bfs_queue_[i];
1091 if (node_excess_[node] > 0) {
1092 DCHECK(IsActive(node));
1093 PushActiveNode(node);
1094 }
1095 }
1096}
1097
1098template <typename Graph, typename ArcFlowT, typename FlowSumT>
1099bool GenericMaxFlow<Graph, ArcFlowT,
1102 const NodeIndex num_nodes = graph_->num_nodes();
1103
1104 // If sink_ or source_ already have kMaxFlowSum, then there is no
1105 // point pushing more flow since it will cause an integer overflow.
1106 if (node_excess_[sink_] == kMaxFlowSum) return false;
1107 if (node_excess_[source_] == -kMaxFlowSum) return false;
1108
1109 bool flow_pushed = false;
1110 for (const ArcIndex arc : graph_->OutgoingArcs(source_)) {
1111 const ArcFlowType flow = residual_arc_capacity_[arc];
1112
1113 // This is a special IsAdmissible() condition for the source.
1114 if (flow == 0 || node_potential_[Head(arc)] >= num_nodes) continue;
1115
1116 // We are careful in case the sum of the flow out of the source is greater
1117 // than kMaxFlowSum to avoid overflow.
1118 const FlowSumType current_flow_out_of_source = -node_excess_[source_];
1119 DCHECK_GE(flow, 0) << flow;
1120 DCHECK_GE(current_flow_out_of_source, 0) << current_flow_out_of_source;
1121 const FlowSumType capped_flow = kMaxFlowSum - current_flow_out_of_source;
1122 if (capped_flow < static_cast<FlowSumType>(flow)) {
1123 // We push as much flow as we can so the current flow on the network will
1124 // be kMaxFlowSum.
1125
1126 // Since at the beginning of the function, current_flow_out_of_source
1127 // was different from kMaxFlowSum, we are sure to have pushed some
1128 // flow before if capped_flow is 0.
1129 if (capped_flow == 0) return true;
1130 PushFlow(static_cast<ArcFlowType>(capped_flow), source_, arc);
1131 return true;
1132 }
1133 PushFlow(flow, source_, arc);
1134 flow_pushed = true;
1135 }
1136 DCHECK_LE(node_excess_[source_], 0);
1137 return flow_pushed;
1138}
1139
1140template <typename Graph, typename ArcFlowT, typename FlowSumT>
1142 NodeIndex tail,
1143 ArcIndex arc) {
1145
1146 // Update the residual capacity of the arc and its opposite arc.
1147 DCHECK_NE(flow, 0);
1148 residual_arc_capacity_[arc] -= flow;
1149 residual_arc_capacity_[Opposite(arc)] += flow;
1150 DCHECK_GE(residual_arc_capacity_[arc], 0);
1151 DCHECK_GE(residual_arc_capacity_[Opposite(arc)], 0);
1152
1153 // Update the excesses at the tail and head of the arc.
1154 //
1155 // Note(user): node_excess_ should be always greater than or equal to 0 except
1156 // for the source where it should always be smaller than or equal to 0. Note
1157 // however that we cannot check this because when we cancel the flow on a
1158 // cycle in PushFlowExcessBackToSource(), we may break this invariant during
1159 // the operation even if it is still valid at the end.
1160 node_excess_[tail] -= static_cast<FlowSumType>(flow);
1161 node_excess_[Head(arc)] += static_cast<FlowSumType>(flow);
1162}
1163
1164template <typename Graph, typename ArcFlowT, typename FlowSumT>
1166 NodeIndex tail, ArcIndex arc, FlowSumType* node_excess) {
1168
1169 // We either push all node_excess or saturate the arc by consuming all its
1170 // residual capacity.
1171 const ArcFlowType flow = static_cast<ArcFlowType>(
1172 std::min(static_cast<FlowSumType>(residual_arc_capacity_[arc]),
1173 node_excess[tail]));
1174
1175 DCHECK_NE(flow, 0);
1176 residual_arc_capacity_[arc] -= flow;
1177 residual_arc_capacity_[Opposite(arc)] += flow;
1178 DCHECK_GE(residual_arc_capacity_[arc], 0);
1179 DCHECK_GE(residual_arc_capacity_[Opposite(arc)], 0);
1180
1181 node_excess[tail] -= static_cast<FlowSumType>(flow);
1182 node_excess[Head(arc)] += static_cast<FlowSumType>(flow);
1183}
1184
1185template <typename Graph, typename ArcFlowT, typename FlowSumT>
1186void GenericMaxFlow<Graph, ArcFlowT,
1190 const NodeIndex num_nodes = graph_->num_nodes();
1191 for (NodeIndex node = 0; node < num_nodes; ++node) {
1192 if (IsActive(node)) {
1193 // This node cannot reach the sink in the residual graph, we don't need to
1194 // consider it anymore, and we will just push its potential excess back to
1195 // the source in PushFlowExcessBackToSource() at the end.
1196 if (node_potential_[node] >= num_nodes) continue;
1197 PushActiveNode(node);
1198 }
1199 }
1200}
1201
1202template <typename Graph, typename ArcFlowT, typename FlowSumT>
1205
1206 const NodeIndex num_nodes = graph_->num_nodes();
1207 std::vector<int> skip_active_node;
1208
1209 // Usually SaturateOutgoingArcsFromSource() will saturate all the arcs from
1210 // the source in one go, and we will loop just once. But in case we can push
1211 // more than kMaxFlowSum out of the source the loop is as follow:
1212 // - Push up to kMaxFlowSum out of the source on the admissible outgoing
1213 // arcs. Stop if no flow was pushed.
1214 // - Compute the current max-flow. This will push some flow back to the
1215 // source and render more outgoing arcs from the source not admissible.
1216 //
1217 // TODO(user): This may not be the most efficient algorithm if we need to loop
1218 // many times. An alternative may be to handle the source like the other nodes
1219 // in the algorithm, initially putting an excess of kMaxFlowSum on it,
1220 // and making the source active like any other node with positive excess. To
1221 // investigate.
1223 int num_skipped;
1224 do {
1225 num_skipped = 0;
1226 skip_active_node.assign(num_nodes, 0);
1227 skip_active_node[sink_] = 2;
1228 skip_active_node[source_] = 2;
1229 GlobalUpdate();
1230 while (!IsEmptyActiveNodeContainer()) {
1232 if (skip_active_node[node] > 1) {
1233 if (node != sink_ && node != source_) ++num_skipped;
1234 continue;
1235 }
1236 const NodeIndex old_height = node_potential_[node];
1237 Discharge(node);
1238
1239 // The idea behind this is that if a node height augments by more than
1240 // one, then it is likely to push flow back the way it came. This can
1241 // lead to very costly loops. A bad case is: source -> n1 -> n2 and n2
1242 // just recently isolated from the sink. Then n2 will push flow back to
1243 // n1, and n1 to n2 and so on. The height of each node will increase by
1244 // steps of two until the height of the source is reached, which can
1245 // take a long time. If the chain is longer, the situation is even
1246 // worse. The behavior of this heuristic is related to the Gap
1247 // heuristic.
1248 //
1249 // Note that the global update will fix all such cases efficiently. So
1250 // the idea is to discharge the active node as much as possible, and
1251 // then do a global update.
1252 //
1253 // We skip a node when this condition was true 2 times to avoid doing a
1254 // global update too frequently.
1255 if (node_potential_[node] > old_height + 1) {
1256 ++skip_active_node[node];
1257 }
1258 }
1259 } while (num_skipped > 0);
1260
1261 // We use a two-phase algorithm:
1262 // 1/ Only deal with nodes that can reach the sink. At the end we know the
1263 // value of the maximum flow and we have a min-cut.
1264 // 2/ Call PushFlowExcessBackToSource() to obtain a max-flow. This is
1265 // usually a lot faster than the first phase.
1267 }
1268}
1269
1270template <typename Graph, typename ArcFlowT, typename FlowSumT>
1272 const NodeIndex node) {
1274 const NodeIndex num_nodes = graph_->num_nodes();
1275
1276 // We cache this as this is a hot-loop and we don't want bound checking.
1277 FlowSumType* const node_excess = node_excess_.get();
1278 NodeHeight* const node_potentials = node_potential_.get();
1279 ArcIndex* const first_admissible_arc = first_admissible_arc_.get();
1280
1281 while (true) {
1282 DCHECK(IsActive(node));
1283 for (const ArcIndex arc :
1284 graph_->OutgoingOrOppositeIncomingArcsStartingFrom(
1285 node, first_admissible_arc_[node])) {
1286 if (IsAdmissible(node, arc, node_potentials)) {
1287 DCHECK(IsActive(node));
1288 const NodeIndex head = Head(arc);
1289 if (node_excess[head] == 0) {
1290 // The push below will make the node active for sure. Note that we may
1291 // push the sink_, but that is handled properly in Refine().
1292 PushActiveNode(head);
1293 }
1294 PushAsMuchFlowAsPossible(node, arc, node_excess);
1295 if (node_excess[node] == 0) {
1296 first_admissible_arc[node] = arc; // arc may still be admissible.
1297 return;
1298 }
1299 }
1300 }
1301 Relabel(node);
1302
1303 // This node can no longer reach the sink, skip until
1304 // PushFlowExcessBackToSource().
1305 if (node_potentials[node] >= num_nodes) break;
1306 }
1307}
1308
1309template <typename Graph, typename ArcFlowT, typename FlowSumT>
1312 // Because we use a relaxed version, this is no longer true if the
1313 // first_admissible_arc_[node] was not actually the first arc!
1314 // DCHECK(CheckRelabelPrecondition(node));
1315 NodeHeight min_height = std::numeric_limits<NodeHeight>::max();
1316 ArcIndex first_admissible_arc = Graph::kNilArc;
1317 for (const ArcIndex arc : graph_->OutgoingOrOppositeIncomingArcs(node)) {
1318 if (residual_arc_capacity_[arc] > 0) {
1319 // Update min_height only for arcs with available capacity.
1320 NodeHeight head_height = node_potential_[Head(arc)];
1321 if (head_height < min_height) {
1322 min_height = head_height;
1323 first_admissible_arc = arc;
1324
1325 // We found an admissible arc at the current height, just stop there.
1326 // This is the true first_admissible_arc_[node].
1327 if (min_height + 1 == node_potential_[node]) break;
1328 }
1329 }
1330 }
1331 DCHECK_NE(first_admissible_arc, Graph::kNilArc);
1332 node_potential_[node] = min_height + 1;
1333
1334 // Note that after a Relabel(), the loop will continue in Discharge(), and
1335 // we are sure that all the arcs before first_admissible_arc are not
1336 // admissible since their height is > min_height.
1337 first_admissible_arc_[node] = first_admissible_arc;
1338}
1339
1340template <typename Graph, typename ArcFlowT, typename FlowSumT>
1342 ArcIndex arc) const {
1343 return graph_->OppositeArc(arc);
1344}
1345
1346template <typename Graph, typename ArcFlowT, typename FlowSumT>
1348 ArcIndex arc) const {
1349 return IsArcValid(arc) && arc >= 0;
1350}
1351
1352template <typename Graph, typename ArcFlowT, typename FlowSumT>
1354 return graph_->IsArcValid(arc);
1355}
1356
1357template <typename Graph, typename ArcFlowT, typename FlowSumT>
1358template <bool reverse>
1360 NodeIndex start, std::vector<NodeIndex>* result) {
1361 // If start is not a valid node index, it can reach only itself.
1362 // Note(user): This is needed because source and sink are given independently
1363 // of the graph and sometimes before it is even constructed.
1364 const NodeIndex num_nodes = graph_->num_nodes();
1365 if (start >= num_nodes) {
1366 result->clear();
1367 result->push_back(start);
1368 return;
1369 }
1370 bfs_queue_.clear();
1371 node_in_bfs_queue_.assign(num_nodes, false);
1372
1373 int queue_index = 0;
1374 bfs_queue_.push_back(start);
1375 node_in_bfs_queue_[start] = true;
1376 while (queue_index != bfs_queue_.size()) {
1377 const NodeIndex node = bfs_queue_[queue_index];
1378 ++queue_index;
1379 for (const ArcIndex arc : graph_->OutgoingOrOppositeIncomingArcs(node)) {
1380 const NodeIndex head = Head(arc);
1381 if (node_in_bfs_queue_[head]) continue;
1382 if (residual_arc_capacity_[reverse ? Opposite(arc) : arc] == 0) continue;
1383 node_in_bfs_queue_[head] = true;
1384 bfs_queue_.push_back(head);
1385 }
1386 }
1387 *result = bfs_queue_;
1388}
1389
1390template <typename Graph, typename ArcFlowT, typename FlowSumT>
1392 FlowModelProto model;
1393 model.set_problem_type(FlowModelProto::MAX_FLOW);
1394 for (int n = 0; n < graph_->num_nodes(); ++n) {
1395 FlowNodeProto* node = model.add_nodes();
1396 node->set_id(n);
1397 if (n == source_) node->set_supply(1);
1398 if (n == sink_) node->set_supply(-1);
1399 }
1400 for (int a = 0; a < graph_->num_arcs(); ++a) {
1401 FlowArcProto* arc = model.add_arcs();
1402 arc->set_tail(graph_->Tail(a));
1403 arc->set_head(graph_->Head(a));
1404 arc->set_capacity(Capacity(a));
1405 }
1406 return model;
1407}
1408
1409} // namespace operations_research
1410
1411#endif // OR_TOOLS_GRAPH_GENERIC_MAX_FLOW_H_
void InitializePreflow()
Initializes the preflow to a state that enables to run Refine.
const Graph * graph_
A pointer to the graph passed as argument.
std::unique_ptr< NodeHeight[]> node_potential_
PriorityQueueWithRestrictedPush< NodeIndex, NodeHeight > active_node_by_height_
bool IsAdmissible(NodeIndex tail, ArcIndex arc, const NodeHeight *potentials) const
Returns true if arc is admissible.
std::unique_ptr< ArcIndex[]> first_admissible_arc_
An array representing the first admissible arc for each node in graph_.
void GetSourceSideMinCut(std::vector< NodeIndex > *result)
void PushFlow(ArcFlowType flow, NodeIndex tail, ArcIndex arc)
NodeIndex GetSourceNodeIndex() const
Returns the index of the node corresponding to the source of the network.
ArcIndex Opposite(ArcIndex arc) const
bool CheckRelabelPrecondition(NodeIndex node) const
StatsGroup stats_
Statistics about this class.
NodeIndex sink_
The index of the sink node in graph_.
void InitializeActiveNodeContainer()
Initializes the container active_nodes_.
std::unique_ptr< ArcFlowType[]> initial_capacity_
GenericMaxFlow & operator=(const GenericMaxFlow &)=delete
const util::ReverseArcStaticGraph< NodeIndex, ArcIndex > * graph() const
FlowSumType Flow(ArcIndex arc) const
NodeIndex Tail(ArcIndex arc) const
std::vector< NodeIndex > bfs_queue_
FlowModelProto CreateFlowModel()
Returns the protocol buffer representation of the current problem.
std::unique_ptr< FlowSumType[]> node_excess_
void PushActiveNode(const NodeIndex &node)
Push element to the active node container.
bool Solve()
Returns true if a maximum flow was solved.
void GetSinkSideMinCut(std::vector< NodeIndex > *result)
NodeIndex source_
The index of the source node in graph_.
void ComputeReachableNodes(NodeIndex start, std::vector< NodeIndex > *result)
void PushAsMuchFlowAsPossible(NodeIndex tail, ArcIndex arc, FlowSumType *node_excess)
NodeIndex Head(ArcIndex arc) const
Handy member functions to make the code more compact.
ArcFlowType Capacity(ArcIndex arc) const
Returns the initial capacity of an arc.
GenericMaxFlow(const Graph *graph, NodeIndex source, NodeIndex sink)
void Refine()
Performs optimization step.
bool IsActive(NodeIndex node) const
static constexpr FlowSumType kMaxFlowSum
Maximum manageable flow.
std::string DebugString(absl::string_view context, ArcIndex arc) const
NodeIndex GetSinkNodeIndex() const
Returns the index of the node corresponding to the sink of the network.
ZVector< ArcFlowType > residual_arc_capacity_
NodeIndex GetAndRemoveFirstActiveNode()
Get the first element from the active node container.
FlowSumType GetOptimalFlow() const
Returns the total flow found by the algorithm.
bool IsEmptyActiveNodeContainer()
Check the emptiness of the container.
Status status_
The status of the problem.
void SetArcCapacity(ArcIndex arc, ArcFlowType new_capacity)
GenericMaxFlow(const GenericMaxFlow &)=delete
This type is neither copyable nor movable.
void Push(Element element, IntegerPriority priority)
PriorityQueueWithRestrictedPush & operator=(const PriorityQueueWithRestrictedPush &)=delete
PriorityQueueWithRestrictedPush(const PriorityQueueWithRestrictedPush &)=delete
This type is neither copyable nor movable.
Base class to print a nice summary of a group of statistics.
Definition stats.h:128
static constexpr bool kHasNegativeReverseArcs
Definition graph.h:201
In SWIG mode, we don't want anything besides these top-level includes.
util::ReverseArcStaticGraph Graph
Type of graph to use.
#define IF_STATS_ENABLED(instructions)
Definition stats.h:417
#define SCOPED_TIME_STAT(stats)
Definition stats.h:418