SciPy

scipy.sparse.csgraph.maximum_flow

scipy.sparse.csgraph.maximum_flow(csgraph, source, sink)

Maximize the flow between two vertices in a graph.

New in version 1.4.0.

Parameters
csgraphcsr_matrix

The square matrix representing a directed graph whose (i, j)’th entry is an integer representing the capacity of the edge between vertices i and j.

sourceint

The source vertex from which the flow flows.

sinkint

The sink vertex to which the flow flows.

Returns
resMaximumFlowResult

A maximum flow represented by a MaximumFlowResult which includes the value of the flow in flow_value, and the residual graph in residual.

Raises
TypeError:

if the input graph is not in CSR format.

ValueError:

if the capacity values are not integers, or the source or sink are out of bounds.

Notes

This solves the maximum flow problem on a given directed weighted graph: A flow associates to every edge a value, also called a flow, less than the capacity of the edge, so that for every vertex (apart from the source and the sink vertices), the total incoming flow is equal to the total outgoing flow. The value of a flow is the sum of the flow of all edges leaving the source vertex, and the maximum flow problem consists of finding a flow whose value is maximal.

By the max-flow min-cut theorem, the maximal value of the flow is also the total weight of the edges in a minimum cut.

To solve the problem, we use the Edmonds–Karp algorithm. [1] This particular implementation strives to exploit sparsity. Its time complexity is \(O(VE^2)\) and its space complexity is \(O(E)\).

The maximum flow problem is usually defined with real valued capacities, but we require that all capacities are integral to ensure convergence. When dealing with rational capacities, or capacities belonging to \(x\mathbb{Q}\) for some fixed \(x \in \mathbb{R}\), it is possible to reduce the problem to the integral case by scaling all capacities accordingly.

Solving a maximum-flow problem can be used for example for graph cuts optimization in computer vision [3].

References

1

Edmonds, J. and Karp, R. M. Theoretical improvements in algorithmic efficiency for network flow problems. 1972. Journal of the ACM. 19 (2): pp. 248-264

2

Cormen, T. H. and Leiserson, C. E. and Rivest, R. L. and Stein C. Introduction to Algorithms. Second Edition. 2001. MIT Press.

3

https://en.wikipedia.org/wiki/Graph_cuts_in_computer_vision

Examples

Perhaps the simplest flow problem is that of a graph of only two vertices with an edge from source (0) to sink (1):

(0) --5--> (1)

Here, the maximum flow is simply the capacity of the edge:

>>> from scipy.sparse import csr_matrix
>>> from scipy.sparse.csgraph import maximum_flow
>>> graph = csr_matrix([[0, 5], [0, 0]])
>>> maximum_flow(graph, 0, 1).flow_value
5

If, on the other hand, there is a bottleneck between source and sink, that can reduce the maximum flow:

(0) --5--> (1) --3--> (2)
>>> graph = csr_matrix([[0, 5, 0], [0, 0, 3], [0, 0, 0]])
>>> maximum_flow(graph, 0, 2).flow_value
3

A less trivial example is given in [2], Chapter 26.1:

>>> graph = csr_matrix([[0, 16, 13,  0,  0,  0],
...                     [0, 10,  0, 12,  0,  0],
...                     [0,  4,  0,  0, 14,  0],
...                     [0,  0,  9,  0,  0, 20],
...                     [0,  0,  0,  7,  0,  4],
...                     [0,  0,  0,  0,  0,  0]])
>>> maximum_flow(graph, 0, 5).flow_value
23

It is possible to reduce the problem of finding a maximum matching in a bipartite graph to a maximum flow problem: Let \(G = ((U, V), E)\) be a bipartite graph. Then, add to the graph a source vertex with edges to every vertex in \(U\) and a sink vertex with edges from every vertex in \(V\). Finally, give every edge in the resulting graph a capacity of 1. Then, a maximum flow in the new graph gives a maximum matching in the original graph consisting of the edges in \(E\) whose flow is positive.

Assume that the edges are represented by a \(\lvert U \rvert \times \lvert V \rvert\) matrix in CSR format whose \((i, j)\)’th entry is 1 if there is an edge from \(i \in U\) to \(j \in V\) and 0 otherwise; that is, the input is of the form required by maximum_bipartite_matching. Then the CSR representation of the graph constructed above contains this matrix as a block. Here’s an example:

>>> graph = csr_matrix([[0, 1, 0, 1], [1, 0, 1, 0], [0, 1, 1, 0]])
>>> print(graph.toarray())
[[0 1 0 1]
 [1 0 1 0]
 [0 1 1 0]]
>>> i, j = graph.shape
>>> n = graph.nnz
>>> indptr = np.concatenate([[0],
...                          graph.indptr + i,
...                          np.arange(n + i + 1, n + i + j + 1),
...                          [n + i + j]])
>>> indices = np.concatenate([np.arange(1, i + 1),
...                           graph.indices + i + 1,
...                           np.repeat(i + j + 1, j)])
>>> data = np.ones(n + i + j, dtype=int)
>>>
>>> graph_flow = csr_matrix((data, indices, indptr))
>>> print(graph_flow.toarray())
[[0 1 1 1 0 0 0 0 0]
 [0 0 0 0 0 1 0 1 0]
 [0 0 0 0 1 0 1 0 0]
 [0 0 0 0 0 1 1 0 0]
 [0 0 0 0 0 0 0 0 1]
 [0 0 0 0 0 0 0 0 1]
 [0 0 0 0 0 0 0 0 1]
 [0 0 0 0 0 0 0 0 1]
 [0 0 0 0 0 0 0 0 0]]

At this point, we can find the maximum flow between the added sink and the added source and the desired matching can be obtained by restricting the residual graph to the block corresponding to the original graph:

>>> flow = maximum_flow(graph_flow, 0, i+j+1)
>>> matching = flow.residual[1:i+1, i+1:i+j+1]
>>> print(matching.toarray())
[[0 1 0 0]
 [1 0 0 0]
 [0 0 1 0]]

This tells us that the first, second, and third vertex in \(U\) are matched with the second, first, and third vertex in \(V\) respectively.

While this solves the maximum bipartite matching problem in general, note that algorithms specialized to that problem, such as maximum_bipartite_matching, will generally perform better.

This approach can also be used to solve various common generalizations of the maximum bipartite matching problem. If, for instance, some vertices can be matched with more than one other vertex, this may be handled by modifying the capacities of the new graph appropriately.

Previous topic

scipy.sparse.csgraph.reverse_cuthill_mckee

Next topic

scipy.sparse.csgraph.maximum_bipartite_matching