Preconditioners

pymatting.ichol

🔗

Signature

ichol( A, discard_threshold=1e-4, shifts=[0.0, 1e-4, 1e-3, 1e-2, 0.1, 0.5, 1.0, 10.0, 100, 1e3, 1e4, 1e5], max_nnz=int(4e9 / 16), relative_discard_threshold=0.0, diag_keep_discarded=True)

Function Description

Implements the thresholded incomplete Cholesky decomposition. For reference, a dense ichol implementation with relative threshold would look like this:
L = np.tril(A) for j in range(n): col = L[j:, j] col -= np.sum(L[j, :j] * L[j:, :j], axis=1) discard_mask = abs(col[1:]) < relative_discard_threshold * np.sum(np.abs(A[j:, j])) col[1:][discard_mask] = 0 col[0] **= 0.5 col[1:] /= col[0]

Parameters

  • A (scipy.sparse.csc_matrix)
      Matrix for which the preconditioner should be calculated
  • discard_threshold (float)
      Values having an absolute value smaller than this threshold will be discarded while calculating the Cholesky decompositions
  • shifts (array of floats)
      Values to try for regularizing the matrix of interest in case it is not positive definite after discarding the small values
  • max_nnz (int)
      Maximum number of non-zero entries in the Cholesky decomposition. Defaults to 250 million, which should usually be around 4 GB.
  • relative_discard_threshold (float)
      Values with an absolute value of less than relative_discard_threshold * sum(abs(A[j:, j])) will be discarded.
  • diag_keep_discarded (bool)
      Whether to update the diagonal with the discarded values. Usually better if True.

Returns

  • chol (CholeskyDecomposition)
      Preconditioner or solver object.

Raises

  • ValueError
      If inappropriate parameter values were passed

Example

>>> from pymatting import * >>> import numpy as np >>> from scipy.sparse import csc_matrix >>> A = np.array([[2.0, 3.0], [3.0, 5.0]]) >>> cholesky_decomposition = ichol(csc_matrix(A)) >>> cholesky_decomposition(np.array([1.0, 2.0])) array([-1., 1.])

pymatting.CholeskyDecomposition.L

🔗

Signature

L(self)

Function Description

Returns the Cholesky factor

Returns

  • L (scipy.sparse.csc_matrix)
      Cholesky factor

pymatting.jacobi

🔗

Signature

jacobi(A)

Function Description

Compute the Jacobi preconditioner function for the matrix A.

Parameters

  • A (np.array)
      Input matrix to compute the Jacobi preconditioner for.

Returns

  • precondition_matvec (function)
      Function which applies the Jacobi preconditioner to a vector

Example

>>> from pymatting import * >>> import numpy as np >>> A = np.array([[2, 3], [3, 5]]) >>> preconditioner = jacobi(A) >>> preconditioner(np.array([1, 2])) array([0.5, 0.4])

pymatting.vcycle

🔗

Signature

vcycle( A, shape, num_pre_iter=1, num_post_iter=1, omega=0.8, direct_solve_size=64, cache=None)

Function Description

Implements the V-Cycle preconditioner. The V-Cycle solver was recommended by [LW14] to solve the alpha matting problem.

Parameters

  • A (numpy.ndarray)
      Input matrix
  • shape (tuple of ints)
      Describing the height and width of the image
  • num_pre_iter (int)
      Number of Jacobi iterations before each V-Cycle, defaults to 1
  • num_post_iter (int)
      Number of Jacobi iterations after each V-Cycle, defaults to 1
  • omega (float)
      Weight parameter for the Jacobi method. If method fails to converge, try different values.

Returns

  • precondition (function)
      Function which applies the V-Cycle preconditioner to a vector

Example

>>> from pymatting import * >>> import numpy as np >>> from scipy.sparse import csc_matrix >>> A = np.array([[2, 3], [3, 5]]) >>> preconditioner = vcycle(A, (2, 2)) >>> preconditioner(np.array([1, 2])) array([-1., 1.])