quadratic_assignment(method=’faq’)#

scipy.optimize.quadratic_assignment(A, B, method='faq', options=None)

Solve the quadratic assignment problem (approximately).

This function solves the Quadratic Assignment Problem (QAP) and the Graph Matching Problem (GMP) using the Fast Approximate QAP Algorithm (FAQ) [1].

Quadratic assignment solves problems of the following form:

\[\begin{split}\min_P & \ {\ \text{trace}(A^T P B P^T)}\\ \mbox{s.t. } & {P \ \epsilon \ \mathcal{P}}\\\end{split}\]

where \(\mathcal{P}\) is the set of all permutation matrices, and \(A\) and \(B\) are square matrices.

Graph matching tries to maximize the same objective function. This algorithm can be thought of as finding the alignment of the nodes of two graphs that minimizes the number of induced edge disagreements, or, in the case of weighted graphs, the sum of squared edge weight differences.

Note that the quadratic assignment problem is NP-hard. The results given here are approximations and are not guaranteed to be optimal.

Parameters:
A2-D array, square

The square matrix \(A\) in the objective function above.

B2-D array, square

The square matrix \(B\) in the objective function above.

methodstr in {‘faq’, ‘2opt’} (default: ‘faq’)

The algorithm used to solve the problem. This is the method-specific documentation for ‘faq’. ‘2opt’ is also available.

Returns:
resOptimizeResult

OptimizeResult containing the following fields.

col_ind1-D array

Column indices corresponding to the best permutation found of the nodes of B.

funfloat

The objective value of the solution.

nitint

The number of Frank-Wolfe iterations performed.

See also

For documentation for the rest of the parameters, see scipy.optimize.quadratic_assignment

Options:
——-
maximizebool (default: False)

Maximizes the objective function if True.

partial_match2-D array of integers, optional (default: None)

Fixes part of the matching. Also known as a “seed” [2].

Each row of partial_match specifies a pair of matched nodes: node partial_match[i, 0] of A is matched to node partial_match[i, 1] of B. The array has shape (m, 2), where m is not greater than the number of nodes, \(n\).

rng{None, int, numpy.random.Generator,

If seed is None (or np.random), the numpy.random.RandomState singleton is used. If seed is an int, a new RandomState instance is used, seeded with seed. If seed is already a Generator or RandomState instance then that instance is used.

P02-D array, “barycenter”, or “randomized” (default: “barycenter”)

Initial position. Must be a doubly-stochastic matrix [3].

If the initial position is an array, it must be a doubly stochastic matrix of size \(m' \times m'\) where \(m' = n - m\).

If "barycenter" (default), the initial position is the barycenter of the Birkhoff polytope (the space of doubly stochastic matrices). This is a \(m' \times m'\) matrix with all entries equal to \(1 / m'\).

If "randomized" the initial search position is \(P_0 = (J + K) / 2\), where \(J\) is the barycenter and \(K\) is a random doubly stochastic matrix.

shuffle_inputbool (default: False)

Set to True to resolve degenerate gradients randomly. For non-degenerate gradients this option has no effect.

maxiterint, positive (default: 30)

Integer specifying the max number of Frank-Wolfe iterations performed.

tolfloat (default: 0.03)

Tolerance for termination. Frank-Wolfe iteration terminates when \(\frac{||P_{i}-P_{i+1}||_F}{\sqrt{m')}} \leq tol\), where \(i\) is the iteration number.

Notes

The algorithm may be sensitive to the initial permutation matrix (or search “position”) due to the possibility of several local minima within the feasible region. A barycenter initialization is more likely to result in a better solution than a single random initialization. However, calling quadratic_assignment several times with different random initializations may result in a better optimum at the cost of longer total execution time.

References

[1]

J.T. Vogelstein, J.M. Conroy, V. Lyzinski, L.J. Podrazik, S.G. Kratzer, E.T. Harley, D.E. Fishkind, R.J. Vogelstein, and C.E. Priebe, “Fast approximate quadratic programming for graph matching,” PLOS one, vol. 10, no. 4, p. e0121002, 2015, DOI:10.1371/journal.pone.0121002

[2]

D. Fishkind, S. Adali, H. Patsolic, L. Meng, D. Singh, V. Lyzinski, C. Priebe, “Seeded graph matching”, Pattern Recognit. 87 (2019): 203-215, DOI:10.1016/j.patcog.2018.09.014

[3]

“Doubly stochastic Matrix,” Wikipedia. https://en.wikipedia.org/wiki/Doubly_stochastic_matrix

Examples

As mentioned above, a barycenter initialization often results in a better solution than a single random initialization.

>>> from numpy.random import default_rng
>>> rng = default_rng()
>>> n = 15
>>> A = rng.random((n, n))
>>> B = rng.random((n, n))
>>> res = quadratic_assignment(A, B)  # FAQ is default method
>>> print(res.fun)
46.871483385480545  # may vary
>>> options = {"P0": "randomized"}  # use randomized initialization
>>> res = quadratic_assignment(A, B, options=options)
>>> print(res.fun)
47.224831071310625 # may vary

However, consider running from several randomized initializations and keeping the best result.

>>> res = min([quadratic_assignment(A, B, options=options)
...            for i in range(30)], key=lambda x: x.fun)
>>> print(res.fun)
46.671852533681516 # may vary

The ‘2-opt’ method can be used to further refine the results.

>>> options = {"partial_guess": np.array([np.arange(n), res.col_ind]).T}
>>> res = quadratic_assignment(A, B, method="2opt", options=options)
>>> print(res.fun)
46.47160735721583 # may vary