def test_kopt_heuristic_single_identity(m): r"""Test k-opt heuristic single search algorithm with identity permutation.""" # create a random matrix A and random permutation of identity matrix a = np.random.uniform(-2.0, 2.0, (m, m)) p0 = np.eye(m) # find and check permutation for when B=A with guess p0=I perm, error = kopt_heuristic_single(lambda x: compute_error(a, a, x, x.T), p0, k=2) assert_equal(perm, np.eye(m)) assert_equal(error, 0.0) # find and check permutation for when B=A with guess p0 being swapped I p0[[m - 2, -1]] = p0[[-1, m - 2]] perm, error = kopt_heuristic_single(lambda x: compute_error(a, a, x, x.T), p0, k=2) assert_equal(perm, np.eye(m)) assert_equal(error, 0.0)
def test_kopt_heuristic_single(): r"""Test k-opt heuristic search algorithm.""" arr_a = np.array([[3, 6, 1, 0, 7], [4, 5, 2, 7, 6], [8, 6, 6, 1, 7], [4, 4, 7, 9, 4], [4, 8, 0, 3, 1]]) arr_b = np.array([[1, 8, 0, 4, 3], [6, 5, 2, 4, 7], [7, 6, 6, 8, 1], [7, 6, 1, 3, 0], [4, 4, 7, 4, 9]]) perm_guess = np.array([[0, 0, 1, 0, 0], [1, 0, 0, 0, 0], [0, 0, 0, 1, 0], [0, 0, 0, 0, 1], [0, 1, 0, 0, 0]]) perm_exact = np.array([[0, 0, 0, 1, 0], [0, 1, 0, 0, 0], [0, 0, 1, 0, 0], [0, 0, 0, 0, 1], [1, 0, 0, 0, 0]]) error_old = compute_error(arr_a, arr_b, perm_guess, perm_guess) perm, kopt_error = kopt_heuristic_single(arr_a, arr_b, error_old, perm_guess, 3, kopt_tol=1.e-8) assert_equal(perm, perm_exact) assert kopt_error == 0 # test the error exceptions assert_raises(ValueError, kopt_heuristic_single, arr_a, arr_b, error_old, perm_guess, 1, kopt_tol=1.e-8)
def test_kopt_heuristic_single_all_permutations(m): r"""Test k-opt heuristic single search algorithm going through all permutations.""" # create a random matrix A and random permutation of identity matrix a = np.random.uniform(-10.0, 10.0, (m, m)) p = np.random.permutation(np.eye(m)) # compute B = P^T A P b = np.linalg.multi_dot([p.T, a, p]) # find and check permutation perm, error = kopt_heuristic_single(lambda x: compute_error(a, b, x, x.T), np.eye(m), m) assert_equal(perm, p) assert_equal(error, 0.0)
def test_kopt_heuristic_single_k_permutations(m): r"""Test k-opt heuristic single search algorithm going upto k permutations.""" # create a random matrix A a = np.random.uniform(-10.0, 10.0, (m, m)) # create permutation matrix by swapping rows m-3 & -1 of identity matrix (this makes sures that # heuristic algorithm only finds the solution towards the end of its search) p = np.eye(m) p[[m - 3, -1]] = p[[-1, m - 3]] # compute B = P^T A P b = np.linalg.multi_dot([p.T, a, p]) # find and check permutation perm, error = kopt_heuristic_single(lambda x: compute_error(a, b, x, x.T), np.eye(m), k=2) assert_equal(perm, p) assert_equal(error, 0.0)
def permutation_2sided( a, b, single=True, method="kopt", guess_p1=None, guess_p2=None, pad=False, unpad_col=False, unpad_row=False, translate=False, scale=False, check_finite=True, options=None, weight=None, lapack_driver="gesvd", ): r"""Perform two-sided permutation Procrustes. Parameters ---------- a : ndarray The 2d-array :math:`\mathbf{A}` which is going to be transformed. b : ndarray The 2d-array :math:`\mathbf{B}` representing the reference matrix. single : bool, optional If `True`, the single-transformation Procrustes is performed to obtain :math:`\mathbf{P}`. If `False`, the two-transformations Procrustes is performed to obtain :math:`\mathbf{P}_1` and :math:`\mathbf{P}_2`. method : str, optional The method to solve for permutation matrices. For `single=False`, these include "flip-flop" and "k-opt" methods. For `single=True`, these include "approx-normal1", "approx-normal2", "approx-umeyama", "approx-umeyama-svd", "k-opt", "soft-assign", and "nmf". guess_p1 : np.ndarray, optional Guess for :math:`\mathbf{P}_1` matrix given as a 2D-array. This is only required for the two-transformations case specified by setting `single=False`. guess_p2 : np.ndarray, optional Guess for :math:`\mathbf{P}_2` matrix given as a 2D-array. pad : bool, optional Add zero rows (at the bottom) and/or columns (to the right-hand side) of matrices :math:`\mathbf{A}` and :math:`\mathbf{B}` so that they have the same shape. unpad_col : bool, optional If True, zero columns (with values less than 1.0e-8) on the right-hand side are removed. unpad_row : bool, optional If True, zero rows (with values less than 1.0e-8) at the bottom are removed. translate : bool, optional If True, both arrays are centered at origin (columns of the arrays will have mean zero). scale : bool, optional If True, both arrays are normalized with respect to the Frobenius norm, i.e., :math:`\text{Tr}\left[\mathbf{A}^\dagger\mathbf{A}\right] = 1` and :math:`\text{Tr}\left[\mathbf{B}^\dagger\mathbf{B}\right] = 1`. check_finite : bool, optional If True, convert the input to an array, checking for NaNs or Infs. options : dict, optional A dictionary of method options. weight : ndarray, optional The 1D-array representing the weights of each row of :math:`\mathbf{A}`. This defines the elements of the diagonal matrix :math:`\mathbf{W}` that is multiplied by :math:`\mathbf{A}` matrix, i.e., :math:`\mathbf{A} \rightarrow \mathbf{WA}`. lapack_driver : {'gesvd', 'gesdd'}, optional Whether to use the more efficient divide-and-conquer approach ('gesdd') or the more robust general rectangular approach ('gesvd') to compute the singular-value decomposition with `scipy.linalg.svd`. Returns ------- res : ProcrustesResult The Procrustes result represented as a class:`utils.ProcrustesResult` object. Notes ----- Given matrix :math:`\mathbf{A}_{n \times n}` and a reference :math:`\mathbf{B}_{n \times n}`, find a permutation of rows/columns of :math:`\mathbf{A}_{n \times n}` that makes it as close as possible to :math:`\mathbf{B}_{n \times n}`. I.e., .. math:: &\underbrace{\text{min}}_{\left\{\mathbf{P} \left| {p_{ij} \in \{0, 1\} \atop \sum_{i=1}^n p_{ij} = \sum_{j=1}^n p_{ij} = 1} \right. \right\}} \|\mathbf{P}^\dagger \mathbf{A} \mathbf{P} - \mathbf{B}\|_{F}^2\\ = &\underbrace{\text{min}}_{\left\{\mathbf{P} \left| {p_{ij} \in \{0, 1\} \atop \sum_{i=1}^n p_{ij} = \sum_{j=1}^n p_{ij} = 1} \right. \right\}} \text{Tr}\left[\left(\mathbf{P}^\dagger\mathbf{A}\mathbf{P} - \mathbf{B} \right)^\dagger \left(\mathbf{P}^\dagger\mathbf{A}\mathbf{P} - \mathbf{B} \right)\right] \\ = &\underbrace{\text{max}}_{\left\{\mathbf{P} \left| {p_{ij} \in \{0, 1\} \atop \sum_{i=1}^n p_{ij} = \sum_{j=1}^n p_{ij} = 1} \right. \right\}} \text{Tr}\left[\mathbf{P}^\dagger\mathbf{A}^\dagger\mathbf{P}\mathbf{B} \right]\\ Here, :math:`\mathbf{P}_{n \times n}` is the permutation matrix. Given an intial guess, the best local minimum can be obtained by the iterative procedure, .. math:: p_{ij}^{(n + 1)} = p_{ij}^{(n)} \sqrt{ \frac{2\left[\mathbf{T}^{(n)}\right]_{ij}}{\left[ \mathbf{P}^{(n)} \left( \left(\mathbf{P}^{(n)}\right)^T \mathbf{T} + \left( \left(\mathbf{P}^{(n)}\right)^T \mathbf{T} \right)^T \right) \right]_{ij}} } where, .. math:: \mathbf{T}^{(n)} = \mathbf{A} \mathbf{P}^{(n)} \mathbf{B} Using an initial guess, the iteration can stops when the change in :math:`d` is below the specified threshold, .. math:: d = \text{Tr} \left[\left(\mathbf{P}^{(n+1)} -\mathbf{P}^{(n)} \right)^T \left(\mathbf{P}^{(n+1)} -\mathbf{P}^{(n)} \right)\right] The outcome of the iterative procedure :math:`\mathbf{P}^{(\infty)}` is not a permutation matrix. So, the closest permutation can be found by setting ``refinement=True``. This uses :class:`procrustes.permutation.PermutationProcrustes` to find the closest permutation; that is, .. math:: \underbrace{\text{min}}_{\left\{\mathbf{P} \left| {p_{ij} \in \{0, 1\} \atop \sum_{i=1}^n p_{ij} = \sum_{j=1}^n p_{ij} = 1} \right. \right\}} \|\mathbf{P} - \mathbf{P}^{(\infty)}\|_{F}^2 = \underbrace{\text{max}}_{\left\{\mathbf{P} \left| {p_{ij} \in \{0, 1\} \atop \sum_{i=1}^n p_{ij} = \sum_{j=1}^n p_{ij} = 1} \right. \right\}} \text{Tr}\left[\mathbf{P}^\dagger\mathbf{P}^{(\infty)} \right] The answer to this problem is a heuristic solution for the matrix-matching problem that seems to be relatively accurate. **Initial Guess:** Two possible initial guesses are inferred from the Umeyama procedure. One can find either the closest permutation matrix to :math:`\mathbf{U}_\text{Umeyama}` or to :math:`\mathbf{U}_\text{Umeyama}^\text{approx.}`. Considering the :class:`procrustes.permutation.PermutationProcrustes`, the resulting permutation matrix can be specified as initial guess through ``guess=umeyama`` and ``guess=umeyama_approx``, which solves: .. math:: \underbrace{\text{max}}_{\left\{\mathbf{P} \left| {p_{ij} \in \{0, 1\} \atop \sum_{i=1}^n p_{ij} = \sum_{j=1}^n p_{ij} = 1} \right. \right\}} \text{Tr}\left[\mathbf{P}^\dagger\mathbf{U}_\text{Umeyama} \right] \\ \underbrace{\text{max}}_{\left\{\mathbf{P} \left| {p_{ij} \in \{0, 1\} \atop \sum_{i=1}^n p_{ij} = \sum_{j=1}^n p_{ij} = 1} \right. \right\}} \text{Tr}\left[\mathbf{P}^\dagger\mathbf{U}_\text{Umeyama}^\text{approx.} \right] Another choice is to start by solving a normal permutation Procrustes problem. In other words, write new matrices, :math:`\mathbf{A}^0` and :math:`\mathbf{B}^0`, with columns like, .. math:: \begin{bmatrix} a_{ii} \\ p \cdot \text{sgn}\left( a_{ij_\text{max}} \right) \underbrace{\text{max}}_{1 \le j \le n} \left(\left|a_{ij}\right|\right)\\ p^2 \cdot \text{sgn}\left( a_{ij_{\text{max}-1}} \right) \underbrace{\text{max}-1}_{1 \le j \le n} \left(\left|a_{ij}\right|\right)\\ \vdots \end{bmatrix} Here, :math:`\text{max}-1` denotes the second-largest absolute value of elements, :math:`\text{max}-2` is the third-largest abosule value of elements, etc. The matrices :math:`\mathbf{A}^0` and :math:`\mathbf{B}^0` have the diagonal elements of :math:`\mathbf{A}` and :math:`\mathbf{B}` in the first row, and below the first row has the largest off-diagonal element in row :math:`i`, the second-largest off-diagonal element, etc. The elements are weighted by a factor :math:`0 < p < 1`, so that the smaller elements are considered less important for matching. The matrices can be truncated after a few terms; for example, after the size of elements falls below some threshold. A reasonable choice would be to stop after :math:`\lfloor \frac{-2\ln 10}{\ln p} +1\rfloor` rows; this ensures that the size of the elements in the last row is less than 1% of those in the first off-diagonal row. There are obviously many different ways to construct the matrices :math:`\mathbf{A}^0` and :math:`\mathbf{B}^0`. Another, even better, method would be to try to encode not only what the off-diagonal elements are, but which element in the matrix they correspond to. One could do that by not only listing the diagonal elements, but also listing the associated off-diagonal element. I.e., the columns of :math:`\mathbf{A}^0` and :math:`\mathbf{B}^0` would be, .. math:: \begin{bmatrix} a_{ii} \\ p \cdot a_{j_\text{max} j_\text{max}} \\ p \cdot \text{sgn}\left( a_{ij_\text{max}} \right) \underbrace{\text{max}}_{1 \le j \le n} \left(\left|a_{ij}\right|\right)\\ p^2 \cdot a_{j_{\text{max}-1} j_{\text{max}-1}} \\ p^2 \cdot \text{sgn}\left( a_{ij_{\text{max}-1}} \right) \underbrace{\text{max}-1}_{1 \le j \le n} \left(\left|a_{ij}\right|\right)\\ \vdots \end{bmatrix} In this case, you would stop the procedure after :math:`m = \left\lfloor {\frac{{ - 4\ln 10}}{{\ln p}} + 1} \right \rfloor` rows. Then one uses the :class:`procrustes.permutation.PermutationProcrustes` to match the constructed matrices :math:`\mathbf{A}^0` and :math:`\mathbf{B}^0` instead of :math:`\mathbf{A}` and :math:`\mathbf{B}`. I.e., .. math:: \underbrace{\text{max}}_{\left\{\mathbf{P} \left| {p_{ij} \in \{0, 1\} \atop \sum_{i=1}^n p_{ij} = \sum_{j=1}^n p_{ij} = 1} \right. \right\}} \text{Tr}\left[\mathbf{P}^\dagger \left(\mathbf{A^0}^\dagger\mathbf{B^0}\right)\right] Please note that the "umeyama_approx" might give inaccurate permutation matrix. More specificity, this is a approximated Umeyama method. One example we can give is that when we compute the permutation matrix that transforms :math:`A` to :math:`B`, the "umeyama_approx" method can not give the exact permutation transformation matrix while "umeyama", "normal1" and "normal2" do. .. math:: A = \begin{bmatrix} 4 & 5 & -3 & 3 \\ 5 & 7 & 3 & -5 \\ -3 & 3 & 2 & 2 \\ 3 & -5 & 2 & 5 \\ \end{bmatrix} \\ B = \begin{bmatrix} 73 & 100 & 73 & -62 \\ 100 & 208 & -116 & 154 \\ 73 & -116 & 154 & 100 \\ -62 & 154 & 100 & 127 \\ \end{bmatrix} \\ References ---------- [1] C. Ding, T. Li and M. I. Jordan, "Nonnegative Matrix Factorization for Combinatorial Optimization: Spectral Clustering, Graph Matching, and Clique Finding," 2008 Eighth IEEE International Conference on Data Mining, Pisa, Italy, 2008, pp. 183-192, doi: 10.1109/ICDM.2008.130. [2] Papadimitriou, Pythagoras. "Parallel solution of SVD-related problems, with applications." PhD diss., University of Manchester, 1993. [3] S. Umeyama. An eigendecomposition approach toweighted graph matching problems. IEEE Trans. on Pattern Analysis and Machine Intelligence, 10:695 –703, 1988. """ # check single argument if not isinstance(single, bool): raise TypeError( f"Argument single is not a boolean! Given type={type(single)}") # check inputs new_a, new_b = setup_input_arrays(a, b, unpad_col, unpad_row, pad, translate, scale, check_finite, weight) # check that A & B are square in case of single transformation if single and new_a.shape[0] != new_a.shape[1]: raise ValueError( f"For single={single}, matrix A should be square but A.shape={new_a.shape}" "Check pad, unpad_col, and unpad_row arguments.") if single and new_b.shape[0] != new_b.shape[1]: raise ValueError( f"For single={single}, matrix B should be square but B.shape={new_b.shape}" "Check pad, unpad_col, and unpad_row arguments.") # print a statement if user-specified guess is not used if method.startswith("approx") and guess_p1 is not None: print( f"Method={method} does not use an initial guess, so guess_p1 is ignored!" ) if method.startswith("approx") and guess_p2 is not None: print( f"Method={method} does not use an initial guess, so guess_p2 is ignored!" ) # get the number of rows & columns of matrix A m, n = new_a.shape # assign & check initial guess for P1 if single and guess_p1 is not None: raise ValueError( f"For single={single}, P1 is transpose of P2, so guess_p1 should be None." ) if not single: if guess_p1 is None: guess_p1 = np.eye(m) if guess_p1.shape != (m, m): raise ValueError( f"Argument guess_p1 should be either None or a ({m}, {m}) array." ) # assign & check initial guess for P2 if guess_p2 is None: guess_p2 = np.eye(n) if guess_p2.shape != (n, n): raise ValueError( f"Argument guess_p2 should be either None or a ({n}, {n}) array.") # check options dictionary & assign default keys defaults = {"tol": 1.0e-8, "maxiter": 500, "k": 3} if options is not None: if not isinstance(options, dict): raise ValueError( f"Argument options should be a dictionary. Given type={type(options)}" ) if not all(k in defaults.keys() for k in options.keys()): raise ValueError( f"Argument options should only have {defaults.keys()} keys. " f"Given options contains {options.keys()} keys!") # update defaults dictionary to use the specified options defaults.update(options) # 2-sided permutation Procrustes with two transformations # ------------------------------------------------------- if not single: if method == "flip-flop": # compute permutations using flip-flop algorithm perm1, perm2, error = _permutation_2sided_2trans_flipflop( new_a, new_b, defaults["tol"], defaults["maxiter"], guess_p1, guess_p2) elif method == "k-opt": # compute permutations using k-opt heuristic search fun_error = lambda p1, p2: compute_error(new_a, new_b, p2, p1.T) perm1, perm2, error = kopt_heuristic_double(fun_error, p1=guess_p1, p2=guess_p2, k=defaults["k"]) else: raise ValueError( f"Method={method} not supported for single={single} transformation!" ) return ProcrustesResult(error=error, new_a=new_a, new_b=new_b, t=perm2, s=perm1) # 2-sided permutation Procrustes with one transformation # ------------------------------------------------------ # The (un)directed iterative procedure for finding the permutation matrix takes the square # root of the matrix entries, which can result in complex numbers if the entries are # negative. To avoid this, all matrix entries are shifted (by the smallest amount) to be # positive. This causes no change to the objective function, as it's a constant value # being added to all entries of a and b. shift = 1.0e-6 if np.min(new_a) < 0 or np.min(new_b) < 0: shift += abs(min(np.min(new_a), np.min(new_b))) # shift is a float, so even if new_a or new_b are ints, the positive matrices are floats # default shift is not zero to avoid division by zero later in the algorithm pos_a = new_a + shift pos_b = new_b + shift if method == "approx-normal1": tmp_a = _approx_permutation_2sided_1trans_normal1(a) tmp_b = _approx_permutation_2sided_1trans_normal1(b) perm = permutation(tmp_a, tmp_b).t elif method == "approx-normal2": tmp_a = _approx_permutation_2sided_1trans_normal2(a) tmp_b = _approx_permutation_2sided_1trans_normal2(b) perm = permutation(tmp_a, tmp_b).t elif method == "approx-umeyama": perm = _approx_permutation_2sided_1trans_umeyama(pos_a, pos_b) elif method == "approx-umeyama-svd": perm = _approx_permutation_2sided_1trans_umeyama_svd( a, b, lapack_driver) elif method == "k-opt": fun_error = lambda p: compute_error(pos_a, pos_b, p, p.T) perm, error = kopt_heuristic_single(fun_error, p0=guess_p2, k=defaults["k"]) elif method == "soft-assign": raise NotImplementedError elif method == "nmf": # check whether A & B are symmetric (within a relative & absolute tolerance) is_pos_a_symmetric = np.allclose(pos_a, pos_a.T, rtol=1.0e-05, atol=1.0e-08) is_pos_b_symmetric = np.allclose(pos_b, pos_b.T, rtol=1.0e-05, atol=1.0e-08) if is_pos_a_symmetric and is_pos_b_symmetric: # undirected graph matching problem (iterative procedure) perm = _permutation_2sided_1trans_undirected( pos_a, pos_b, guess_p2, defaults['tol'], defaults['maxiter']) else: # directed graph matching problem (iterative procedure) perm = _permutation_2sided_1trans_directed(pos_a, pos_b, guess_p2, defaults['tol'], defaults['maxiter']) else: raise ValueError( f"Method={method} not supported for single={single} transformation!" ) # some of the methods for 2-sided-1-transformation permutation procrustes does not produce a # permutation matrix. So, their output is treated like a guess, and the closest permutation # matrix is found using 1-sided permutation procrustes (where A=I & B=perm) # Even though this step is not needed for ALL methods (e.g. k-opt, normal1, & normal2), to # make the code simple, this step is performed for all methods as its cost is negligible. perm = permutation( np.eye(perm.shape[0]), perm, translate=False, scale=False, unpad_col=False, unpad_row=False, check_finite=True, ).t # compute error error = compute_error(new_a, new_b, t=perm, s=perm.T) return ProcrustesResult(error=error, new_a=new_a, new_b=new_b, t=perm, s=perm.T)
def softassign( a: np.ndarray, b: np.ndarray, pad: bool = True, translate: bool = False, scale: bool = False, unpad_col: bool = False, unpad_row: bool = False, check_finite: bool = True, weight: Optional[np.ndarray] = None, iteration_soft: int = 50, iteration_sink: int = 200, beta_r: float = 1.10, beta_f: float = 1.0e5, epsilon: float = 0.05, epsilon_soft: float = 1.0e-3, epsilon_sink: float = 1.0e-3, k: float = 0.15, gamma_scaler: float = 1.01, n_stop: int = 3, adapted: bool = True, beta_0: Optional[float] = None, m_guess: Optional[float] = None, iteration_anneal: Optional[int] = None, kopt: bool = False, kopt_k: int = 3, ) -> ProcrustesResult: r""" Find the transformation matrix for 2-sided permutation Procrustes with softassign algorithm. Parameters ---------- a : ndarray The 2D-array :math:`\mathbf{A}_{m \times n}` which is going to be transformed. b : ndarray The 2D-array :math:`\mathbf{B}_{m \times n}` representing the reference. pad : bool, optional Add zero rows (at the bottom) and/or columns (to the right-hand side) of matrices :math:`\mathbf{A}` and :math:`\mathbf{B}` so that they have the same shape. translate : bool, optional If True, both arrays are centered at origin (columns of the arrays will have mean zero). scale : bool, optional If True, both arrays are normalized with respect to the Frobenius norm, i.e., :math:`\text{Tr}\left[\mathbf{A}^\dagger\mathbf{A}\right] = 1` and :math:`\text{Tr}\left[\mathbf{B}^\dagger\mathbf{B}\right] = 1`. unpad_col : bool, optional If True, zero columns (with values less than 1.0e-8) on the right-hand side of the intial :math:`\mathbf{A}` and :math:`\mathbf{B}` matrices are removed. unpad_row : bool, optional If True, zero rows (with values less than 1.0e-8) at the bottom of the intial :math:`\mathbf{A}` and :math:`\mathbf{B}` matrices are removed. check_finite : bool, optional If true, convert the input to an array, checking for NaNs or Infs. Default=True. weight : ndarray, optional The 1D-array representing the weights of each row of :math:`\mathbf{A}`. This defines the elements of the diagonal matrix :math:`\mathbf{W}` that is multiplied by :math:`\mathbf{A}` matrix, i.e., :math:`\mathbf{A} \rightarrow \mathbf{WA}`. iteration_soft : int, optional Number of iterations for softassign loop. iteration_sink : int, optional Number of iterations for Sinkhorn loop. beta_r : float, optional Annealing rate which should greater than 1. beta_f : float, optional The final inverse temperature. epsilon : float, optional The tolerance value for annealing loop. epsilon_soft : float, optional The tolerance value used for softassign. epsilon_sink : float, optional The tolerance value used for Sinkhorn loop. If adapted version is used, it will use the adapted tolerance value for Sinkhorn instead. k : float, optional This parameter controls how much tighter the coverage threshold for the interior loop should be than the coverage threshold for the loops outside. It has be be within the integral :math:`(0,1)`. gamma_scaler : float, optional This parameter ensures the quadratic cost function including self-amplification positive define. n_stop : int, optional Number of running steps after the calculation converges in the relaxation procedure. adapted : bool, optional If adapted, this function will use the tighter covergence threshold for the interior loops. beta_0 : float, optional Initial inverse temperature. beta_f : float, optional Final inverse temperature. m_guess : ndarray, optional The initial guess of the doubly-stochastic matrix. iteration_anneal : int, optional Number of iterations for annealing loop. kopt : bool, optional If True, the k_opt heuristic search will be performed. kopt_k : int, optional Defines the oder of k-opt heuristic local search. For example, kopt_k=3 leads to a local search of 3 items and kopt_k=2 only searches for two items locally. weight : ndarray, optional The weighting matrix. Returns ------- res : ProcrustesResult The Procrustes result represented as a class:`utils.ProcrustesResult` object. Notes ----- Quadratic assignment problem (QAP) has played a very special but fundamental role in combinatorial optimization problems. The problem can be defined as a optimization problem to minimize the cost to assign a set of facilities to a set of locations. The cost is a function of the flow between the facilities and the geographical distances among various facilities. The objective function (also named loss function in machine learning) is defined as .. math:: E_{qap}(M, \mu, \nu) = - \frac{1}{2}\Sigma_{aibj}C_{ai;bj}M_{ai}M_{bj} + \Sigma_{a}{\mu}_a (\Sigma_i M_{ai} -1) \\ + \Sigma_i {\nu}_i (\Sigma_i M_{ai} -1) - \frac{\gamma}{2}\Sigma_{ai} {M_{ai}}^2 + \frac{1}{\beta} \Sigma_{ai} M_{ai}\log{M_{ai}} where :math:`C_{ai,bj}` is the benefit matrix, :math:`M` is the desired :math:`N \times N` permutation matrix. :math:`E` is the energy function which comes along with a self-amplification term with `\gamma`, two Lagrange parameters :math:`\mu` and :math:`\nu` for constrained optimization and :math:`M_{ai} \log{M_{ai}}` servers as a barrier function which ensures positivity of :math:`M_{ai}`. The inverse temperature :math:`\beta` is a deterministic annealing control parameter. Examples -------- >>> import numpy as np >>> array_a = np.array([[4, 5, 3, 3], [5, 7, 3, 5], ... [3, 3, 2, 2], [3, 5, 2, 5]]) # define a random matrix >>> perm = np.array([[0., 0., 1., 0.], [1., 0., 0., 0.], ... [0., 0., 0., 1.], [0., 1., 0., 0.]]) # define b by permuting array_a >>> b = np.dot(perm.T, np.dot(a, perm)) >>> new_a, new_b, M_ai, error = softassign(a,b,unpad_col=False,unpad_row=False) >>> M_ai # the permutation matrix array([[0., 0., 1., 0.], [1., 0., 0., 0.], [0., 0., 0., 1.], [0., 1., 0., 0.]]) >>> error # the error 0.0 """ # pylint: disable-msg=too-many-arguments # pylint: disable-msg=too-many-branches # todo: add linear_cost_func with default value 0 # Check beta_r if beta_r <= 1: raise ValueError("Argument beta_r cannot be less than 1.") new_a, new_b = setup_input_arrays(a, b, unpad_col, unpad_row, pad, translate, scale, check_finite, weight) # Check that A & B are square and that they match each other. if new_a.shape[0] != new_a.shape[1]: raise ValueError(f"Matrix A should be square but A.shape={new_a.shape}" "Check pad, unpad_col, and unpad_row arguments.") if new_b.shape[0] != new_b.shape[1]: raise ValueError(f"Matrix B should be square but B.shape={new_b.shape}" "Check pad, unpad_col, and unpad_row arguments.") if new_a.shape != new_b.shape: raise ValueError(f"New matrix A {new_a.shape} should match the new" f" matrix B shape {new_b.shape}.") # Initialization # Compute the benefit matrix array_c = np.kron(new_a, new_b) # Get the shape of A (B and the permutation matrix as well) row_num = new_a.shape[0] c_tensor = array_c.reshape(row_num, row_num, row_num, row_num) # Compute the beta_0 gamma = _compute_gamma(array_c, row_num, gamma_scaler) if beta_0 is None: c_gamma = array_c + gamma * (np.identity(row_num * row_num)) eival_gamma = np.amax(np.abs(np.linalg.eigvalsh(c_gamma))) beta_0 = gamma_scaler * max(1.0e-10, eival_gamma / row_num) beta_0 = 1 / beta_0 else: beta_0 *= row_num beta = beta_0 # We will use iteration_anneal if provided even if the final inverse temperature is specified # iteration_anneal is not None, beta_f can be None or not if iteration_anneal is not None: beta_f = beta_0 * np.power(beta_r, iteration_anneal) * row_num # iteration_anneal is None and beta_f is not None elif iteration_anneal is None and beta_f is not None: beta_f *= row_num # Both iteration_anneal and beta_f are None else: raise ValueError( "We must specify at least one of iteration_anneal and beta_f and " "specify only one is recommended.") # Initialization of m_ai # check shape of m_guess if m_guess is not None: if np.any(m_guess < 0): raise ValueError( "The initial guess of permutation matrix cannot contain any negative values." ) if m_guess.shape[0] == row_num and m_guess.shape[1] == row_num: array_m = m_guess else: warnings.warn( f"The shape of m_guess does not match ({row_num}, {row_num})." "Use random initial guess instead.") array_m = np.abs( np.random.normal(loc=1.0, scale=0.1, size=(row_num, row_num))) else: # m_relax_old = 1 / N + np.random.rand(N, N) array_m = np.abs( np.random.normal(loc=1.0, scale=0.1, size=(row_num, row_num))) array_m[array_m < 0] = 0 array_m = array_m / row_num nochange = 0 if adapted: epsilon_sink = epsilon_soft * k while beta < beta_f: # relaxation m_old_beta = deepcopy(array_m) # softassign loop for _ in np.arange(iteration_soft): m_old_soft = deepcopy(array_m) # Compute Z in relaxation step # C_gamma_tensor = C_gamma.reshape(N, N, N, N) # Z = -np.einsum('ijkl,jl->ik', C_gamma_tensor, M) # Z -= linear_cost_func array_z = np.einsum("aibj,bj->ai", c_tensor, array_m) array_z += gamma * array_m # soft_assign array_m = np.exp(beta * array_z) # Sinkhorn loop for _ in np.arange(iteration_sink): # Row normalization array_m = array_m / array_m.sum(axis=1, keepdims=1) # Column normalization array_m = array_m / array_m.sum(axis=0, keepdims=1) # Compute the delata_M_sink if np.amax(np.abs(array_m.sum(axis=1, keepdims=1) - 1)) < epsilon_sink: array_m = array_m / array_m.sum(axis=1, keepdims=1) break change_soft = np.amax(np.abs(array_m - m_old_soft)) # pylint: disable-msg=no-else-break if change_soft < epsilon_soft: break else: if adapted: epsilon_sink = change_soft * k else: continue change_annealing = np.amax(np.abs(array_m - m_old_beta)) if change_annealing < epsilon: nochange += 1 if nochange > n_stop: break else: nochange = 0 beta *= beta_r if adapted: epsilon_soft = change_soft * k epsilon_sink = epsilon_soft * k # Compute the error array_m = permutation(np.eye(array_m.shape[0]), array_m)["t"] # k-opt heuristic if kopt: fun_error = lambda p: compute_error(new_a, new_b, p, p.T) array_m, error = kopt_heuristic_single(fun_error, p0=array_m, k=kopt_k) else: error = compute_error(new_a, new_b, array_m, array_m.T) return ProcrustesResult(error=error, new_a=new_a, new_b=new_b, t=array_m, s=None)
def permutation_2sided(array_a, array_b, transform_mode="single", remove_zero_col=True, remove_zero_row=True, pad_mode="row-col", translate=False, scale=False, mode="normal1", check_finite=True, iteration=500, add_noise=False, tol=1.0e-8, kopt=False, kopt_k=3, kopt_tol=1.e-8, weight=None): r"""Double sided permutation Procrustes. Parameters ---------- array_a : ndarray The 2d-array :math:`\mathbf{A}_{m \times n}` which is going to be transformed. array_b : ndarray The 2d-array :math:`\mathbf{B}_{m \times n}` representing the reference. transform_mode : str, optional When transform_mode="single", it is the two-sided permutation Procrustes with one transformation. (1). If the input matrices (adjacency matrices) are symmetric within the threshold of 1.e-5, undirected graph matching algorithm will be applied. (2). If the input matrices (adjacency matrices) are asymmetric, the directed graph matching is applied. When transform_mode="double", the problem becomes two-sided permutation Procrustes with two transformations (denoted as regular two-sided permutation Procrustes here). An flip-flop algorithm taken from *Parallel solution of SVD-related problems, with applications, Pythagoras Papadimitriou, Ph.D. Thesis, University of Manchester, 1993* is used to solve the problem. Default="single". Otherwise, transform_mode="double", the two-sided permutation Procrustes with two transformations will be performed. Default="single_undirected". remove_zero_col : bool, optional If True, zero columns (values less than 1e-8) on the right side will be removed. Default= True. remove_zero_row : bool, optional If True, zero rows (values less than 1e-8) on the bottom will be removed. Default= True. pad_mode : str, optional Specifying how to pad the arrays, listed below. Default="row-col". - "row" The array with fewer rows is padded with zero rows so that both have the same number of rows. - "col" The array with fewer columns is padded with zero columns so that both have the same number of columns. - "row-col" The array with fewer rows is padded with zero rows, and the array with fewer columns is padded with zero columns, so that both have the same dimensions. This does not necessarily result in square arrays. - "square" The arrays are padded with zero rows and zero columns so that they are both squared arrays. The dimension of square array is specified based on the highest dimension, i.e. :math:`\text{max}(n_a, m_a, n_b, m_b)`. translate : bool, optional If True, both arrays are translated to be centered at origin, ie columns of the arrays will have mean zero. Default=False. scale : bool, optional If True, both arrays are normalized to one with respect to the Frobenius norm, ie :math:`Tr(A^T A) = 1`. Default=False. mode : string, optional Option for choosing the initial guess methods, including "normal1", "normal2", "umeyama" and "umeyama_approx". "umeyama_approx" is the approximated umeyama method. check_finite : bool, optional If true, convert the input to an array, checking for NaNs or Infs. Default=True. iteration : int, optional Maximum number for iterations. Default=500. add_noise : bool, optional Add small noise if the arrays are non-diagonalizable. Default=False. tol : float, optional The tolerance value used for updating the initial guess. Default=1.e-8. kopt : bool, optional If True, the k_opt heuristic search will be performed. Default=False. kopt_k : int, optional Defines the oder of k-opt heuristic local search. For example, kopt_k=3 leads to a local search of 3 items and kopt_k=2 only searches for two items locally. Default=3. kopt_tol : float, optional Tolerance value to check if k-opt heuristic converges. Default=1.e-8. Returns ------- res : ProcrustesResult Procrustes analysis result object. Attributes ---------- new_a : ndarray The transformed ndarray A. new_b : ndarray The transformed ndarray B. array_u : ndarray The optimum permutation transformation matrix. array_p : ndarray The optimum permutation transformation matrix when using double transform mode. array_q : ndarray The optimum permutation transformation matrix when using double transform mode. error : float Two-sided permutation Procrustes error. Notes ----- Given matrix :math:`\mathbf{A}_{n \times n}` and a reference :math:`\mathbf{B}_{n \times n}`, find a permutation of rows/columns of :math:`\mathbf{A}_{n \times n}` that makes it as close as possible to :math:`\mathbf{B}_{n \times n}`. I.e., .. math:: &\underbrace{\text{min}}_{\left\{\mathbf{P} \left| {p_{ij} \in \{0, 1\} \atop \sum_{i=1}^n p_{ij} = \sum_{j=1}^n p_{ij} = 1} \right. \right\}} \|\mathbf{P}^\dagger \mathbf{A} \mathbf{P} - \mathbf{B}\|_{F}^2\\ = &\underbrace{\text{min}}_{\left\{\mathbf{P} \left| {p_{ij} \in \{0, 1\} \atop \sum_{i=1}^n p_{ij} = \sum_{j=1}^n p_{ij} = 1} \right. \right\}} \text{Tr}\left[\left(\mathbf{P}^\dagger\mathbf{A}\mathbf{P} - \mathbf{B} \right)^\dagger \left(\mathbf{P}^\dagger\mathbf{A}\mathbf{P} - \mathbf{B} \right)\right] \\ = &\underbrace{\text{max}}_{\left\{\mathbf{P} \left| {p_{ij} \in \{0, 1\} \atop \sum_{i=1}^n p_{ij} = \sum_{j=1}^n p_{ij} = 1} \right. \right\}} \text{Tr}\left[\mathbf{P}^\dagger\mathbf{A}^\dagger\mathbf{P}\mathbf{B} \right]\\ Here, :math:`\mathbf{P}_{n \times n}` is the permutation matrix. Given an intial guess, the best local minimum can be obtained by the iterative procedure, .. math:: p_{ij}^{(n + 1)} = p_{ij}^{(n)} \sqrt{ \frac{2\left[\mathbf{T}^{(n)}\right]_{ij}}{\left[ \mathbf{P}^{(n)} \left( \left(\mathbf{P}^{(n)}\right)^T \mathbf{T} + \left( \left(\mathbf{P}^{(n)}\right)^T \mathbf{T} \right)^T \right) \right]_{ij}} } where, .. math:: \mathbf{T}^{(n)} = \mathbf{A} \mathbf{P}^{(n)} \mathbf{B} Using an initial guess, the iteration can stops when the change in :math:`d` is below the specified threshold, .. math:: d = \text{Tr} \left[\left(\mathbf{P}^{(n+1)} -\mathbf{P}^{(n)} \right)^T \left(\mathbf{P}^{(n+1)} -\mathbf{P}^{(n)} \right)\right] The outcome of the iterative procedure :math:`\mathbf{P}^{(\infty)}` is not a permutation matrix. So, the closest permutation can be found by setting ``refinement=True``. This uses :class:`procrustes.permutation.PermutationProcrustes` to find the closest permutation; that is, .. math:: \underbrace{\text{min}}_{\left\{\mathbf{P} \left| {p_{ij} \in \{0, 1\} \atop \sum_{i=1}^n p_{ij} = \sum_{j=1}^n p_{ij} = 1} \right. \right\}} \|\mathbf{P} - \mathbf{P}^{(\infty)}\|_{F}^2 = \underbrace{\text{max}}_{\left\{\mathbf{P} \left| {p_{ij} \in \{0, 1\} \atop \sum_{i=1}^n p_{ij} = \sum_{j=1}^n p_{ij} = 1} \right. \right\}} \text{Tr}\left[\mathbf{P}^\dagger\mathbf{P}^{(\infty)} \right] The answer to this problem is a heuristic solution for the matrix-matching problem that seems to be relatively accurate. **Initial Guess:** Two possible initial guesses are inferred from the Umeyama procedure. One can find either the closest permutation matrix to :math:`\mathbf{U}_\text{Umeyama}` or to :math:`\mathbf{U}_\text{Umeyama}^\text{approx.}`. Considering the :class:`procrustes.permutation.PermutationProcrustes`, the resulting permutation matrix can be specified as initial guess through ``guess=umeyama`` and ``guess=umeyama_approx``, which solves: .. math:: \underbrace{\text{max}}_{\left\{\mathbf{P} \left| {p_{ij} \in \{0, 1\} \atop \sum_{i=1}^n p_{ij} = \sum_{j=1}^n p_{ij} = 1} \right. \right\}} \text{Tr}\left[\mathbf{P}^\dagger\mathbf{U}_\text{Umeyama} \right] \\ \underbrace{\text{max}}_{\left\{\mathbf{P} \left| {p_{ij} \in \{0, 1\} \atop \sum_{i=1}^n p_{ij} = \sum_{j=1}^n p_{ij} = 1} \right. \right\}} \text{Tr}\left[\mathbf{P}^\dagger\mathbf{U}_\text{Umeyama}^\text{approx.} \right] Another choice is to start by solving a normal permutation Procrustes problem. In other words, write new matrices, :math:`\mathbf{A}^0` and :math:`\mathbf{B}^0`, with columns like, .. math:: \begin{bmatrix} a_{ii} \\ p \cdot \text{sgn}\left( a_{ij_\text{max}} \right) \underbrace{\text{max}}_{1 \le j \le n} \left(\left|a_{ij}\right|\right)\\ p^2 \cdot \text{sgn}\left( a_{ij_{\text{max}-1}} \right) \underbrace{\text{max}-1}_{1 \le j \le n} \left(\left|a_{ij}\right|\right)\\ \vdots \end{bmatrix} Here, :math:`\text{max}-1` denotes the second-largest absolute value of elements, :math:`\text{max}-2` is the third-largest abosule value of elements, etc. The matrices :math:`\mathbf{A}^0` and :math:`\mathbf{B}^0` have the diagonal elements of :math:`\mathbf{A}` and :math:`\mathbf{B}` in the first row, and below the first row has the largest off-diagonal element in row :math:`i`, the second-largest off-diagonal element, etc. The elements are weighted by a factor :math:`0 < p < 1`, so that the smaller elements are considered less important for matching. The matrices can be truncated after a few terms; for example, after the size of elements falls below some threshold. A reasonable choice would be to stop after :math:`\lfloor \frac{-2\ln 10}{\ln p} +1\rfloor` rows; this ensures that the size of the elements in the last row is less than 1% of those in the first off-diagonal row. There are obviously many different ways to construct the matrices :math:`\mathbf{A}^0` and :math:`\mathbf{B}^0`. Another, even better, method would be to try to encode not only what the off-diagonal elements are, but which element in the matrix they correspond to. One could do that by not only listing the diagonal elements, but also listing the associated off-diagonal element. I.e., the columns of :math:`\mathbf{A}^0` and :math:`\mathbf{B}^0` would be, .. math:: \begin{bmatrix} a_{ii} \\ p \cdot a_{j_\text{max} j_\text{max}} \\ p \cdot \text{sgn}\left( a_{ij_\text{max}} \right) \underbrace{\text{max}}_{1 \le j \le n} \left(\left|a_{ij}\right|\right)\\ p^2 \cdot a_{j_{\text{max}-1} j_{\text{max}-1}} \\ p^2 \cdot \text{sgn}\left( a_{ij_{\text{max}-1}} \right) \underbrace{\text{max}-1}_{1 \le j \le n} \left(\left|a_{ij}\right|\right)\\ \vdots \end{bmatrix} In this case, you would stop the procedure after :math:`m = \left\lfloor {\frac{{ - 4\ln 10}}{{\ln p}} + 1} \right \rfloor` rows. Then one uses the :class:`procrustes.permutation.PermutationProcrustes` to match the constructed matrices :math:`\mathbf{A}^0` and :math:`\mathbf{B}^0` instead of :math:`\mathbf{A}` and :math:`\mathbf{B}`. I.e., .. math:: \underbrace{\text{max}}_{\left\{\mathbf{P} \left| {p_{ij} \in \{0, 1\} \atop \sum_{i=1}^n p_{ij} = \sum_{j=1}^n p_{ij} = 1} \right. \right\}} \text{Tr}\left[\mathbf{P}^\dagger \left(\mathbf{A^0}^\dagger\mathbf{B^0}\right)\right] Please note that the "umeyama_approx" might give inaccurate permutation matrix. More specificity, this is a approximated Umeyama method. One example we can give is that when we compute the permutation matrix that transforms :math:`A` to :math:`B`, the "umeyama_approx" method can not give the exact permutation transformation matrix while "umeyama", "normal1" and "normal2" do. .. math:: A = \begin{bmatrix} 4 & 5 & -3 & 3 \\ 5 & 7 & 3 & -5 \\ -3 & 3 & 2 & 2 \\ 3 & -5 & 2 & 5 \\ \end{bmatrix} \\ B = \begin{bmatrix} 73 & 100 & 73 & -62 \\ 100 & 208 & -116 & 154 \\ 73 & -116 & 154 & 100 \\ -62 & 154 & 100 & 127 \\ \end{bmatrix} \\ """ # check inputs new_a, new_b = setup_input_arrays(array_a, array_b, remove_zero_col, remove_zero_row, pad_mode, translate, scale, check_finite, weight) # np.power() can not handle the negatives values # Try to convert the matrices to non-negative # shift the the matrices to avoid negative values # otherwise it will cause an error in the Eq. 28 in the research notes maximum = max(np.amax(np.abs(new_a)), np.amax(np.abs(new_b))) new_a_positive = new_a.astype(np.float) + maximum new_b_positive = new_b.astype(np.float) + maximum # Do single-transformation computation if requested transform_mode = transform_mode.lower() if transform_mode == "single": # algorithm for undirected graph matching problem # check if two matrices are symmetric within a relative tolerance and absolute tolerance. if np.allclose(new_a_positive, new_a_positive.T, rtol=1.e-05, atol=1.e-08) and \ np.allclose(new_b_positive, new_b_positive.T, rtol=1.e-05, atol=1.e-08): # the initial guess guess = _guess_initial_permutation_undirected( new_a_positive, new_b_positive, mode, add_noise) # Compute the permutation matrix by iterations array_u = _compute_transform(new_a_positive, new_b_positive, guess, tol, iteration) error = compute_error(new_a_positive, new_b_positive, array_u, array_u) # k-opt heuristic if kopt: array_u, error = kopt_heuristic_single(array_a=new_a_positive, array_b=new_b_positive, ref_error=error, perm=array_u, kopt_k=kopt_k, kopt_tol=kopt_tol) # algorithm for directed graph matching problem else: # the initial guess guess = _2sided_1trans_initial_guess_directed( new_a_positive, new_b_positive) # Compute the permutation matrix by iterations array_u = _compute_transform_directed(new_a_positive, new_b_positive, guess, tol, iteration) error = compute_error(new_a_positive, new_b_positive, array_u, array_u) # k-opt heuristic if kopt: array_u, error = kopt_heuristic_single(array_a=new_a_positive, array_b=new_b_positive, ref_error=error, perm=array_u, kopt_k=kopt_k, kopt_tol=kopt_tol) return ProcrustesResult(new_a=new_a, new_b=new_b, array_u=array_u, error=error) # Do regular computation elif transform_mode == "double": array_m = new_a array_n = new_b array_p, array_q, error = _2sided_regular(array_m, array_n, tol, iteration) # perform k-opt heuristic search twice if kopt: array_p, array_q, error = kopt_heuristic_double(array_m=array_m, array_n=array_n, ref_error=error, perm_p=array_p, perm_q=array_q, kopt_k=kopt_k, kopt_tol=kopt_tol) # return array_m, array_n, array_p, array_q, error return ProcrustesResult(new_a=new_a, new_b=new_b, array_p=array_p, array_q=array_q, error=error) else: raise ValueError( """Invalid transform_mode argument, use "single"or "double".""")