cholupdates.rank_1.downdate

cholupdates.rank_1.downdate(L, v, check_diag=True, overwrite_L=False, overwrite_v=False, method='seeger', **method_kwargs)

Update a Cholesky factorization after subtraction of a symmetric positive semidefinite rank-1 matrix.

In other words, given \(A = L L^T \in \mathbb{R}^{N \times N}\) and \(v \in \mathbb{R}^N\), compute \(L'\) such that

\[A' := A - v v^T = L' L'^T.\]
Parameters
  • L ((N, N) numpy.ndarray) – Lower-triangular Cholesky factor of \(A\). Must have a non-zero diagonal. The entries in the strict upper-triangular part of L can contain arbitrary values, since the algorithm neither reads from nor writes to this part of the matrix. This behavior is useful when using the Cholesky factors returned by scipy.linalg.cho_factor() which contain arbitrary values on the irrelevant triangular part of the matrix.

  • v ((N,) numpy.ndarray) – The vector \(v\) with shape (N,) defining the symmetric rank-1 update matrix \(v v^T\).

  • check_diag (bool) – If set to True, the function will check whether the diagonal of the given Cholesky factor L is non-zero and raise a ValueError if this is not the case. Setting check_diag to False can be used to speed up computations if it is clear that the Cholesky factor can not have zeros on its diagonal. Caution: If this argument is set to False and the Cholesky factor does contain zeros on its diagonal, the behavior of the function will be undefined.

  • overwrite_L (bool) – If set to True, the function may overwrite the array L with the upper Cholesky factor \(L'\) of \(A'\), i.e. the result is computed in-place. Passing False here ensures that the array L is not modified.

  • overwrite_v (bool) – If set to True, the function may reuse the array v as an internal computation buffer, which will modify v. Passing False here ensures that the array v is not modified.

  • method (str) –

    Algorithm to be used to compute the updated Cholesky factor. Must be one of

    Defaults to “seeger”.

  • method_kwargs (Dict[str, Any]) – Additional keyword arguments which will be passed to the function selected by method.

Returns

Lower triangular Cholesky factor \(L'\) of \(A - v v^T\). The diagonal entries of this matrix are guaranteed to be positive. The strict upper-triangular part of this matrix will contain the values from the upper-triangular part of L.

Return type

(N, N) numpy.ndarray, dtype=L.dtype

Raises
  • ValueError – If L does not have shape (N, N) for some N.

  • ValueError – If v does not have shape (N,), while L has shape (N, N).

  • numpy.linalg.LinAlgError – If the diagonal of L contains zeros and check_diag is set to True.

  • numpy.linalg.LinAlgError – If the downdate results in a matrix \(L'\), which is not positive definite.

  • Exception – Any exception raised by the function specified by method.

See also

cholupdates.rank_1.update

A similar function which performs a symmetric rank 1 update instead of a downdate.

Examples

Consider the following matrix-vector pair

>>> A = np.array([[ 7.77338976,  1.27256923,  1.58075291],
...               [ 1.27256923,  8.29126934,  0.80466256],
...               [ 1.58075291,  0.80466256, 13.65749896]])
>>> v = np.array([1.60994441, 0.21482681, 0.78780241])

We want to compute the lower-triangular Cholesky factor L_dd of

>>> A_dd = A - np.outer(v, v)
>>> A_dd
array([[ 5.18146876,  0.92671001,  0.31243482],
       [ 0.92671001,  8.24511878,  0.63542148],
       [ 0.31243482,  0.63542148, 13.03686632]])

We assume that the lower-triangular Cholesky factor L of A is given

>>> import scipy.linalg
>>> L = scipy.linalg.cholesky(A, lower=True)
>>> L
array([[2.78807994, 0.        , 0.        ],
       [0.45643212, 2.84305101, 0.        ],
       [0.56696829, 0.19200501, 3.64680408]])

The function cholupdates.rank_1.update() computes L_dd efficiently from L and v

>>> import cholupdates
>>> L_dd = cholupdates.rank_1.downdate(L, v, method="seeger")
>>> L_dd
array([[2.27628398, 0.        , 0.        ],
       [0.40711529, 2.8424243 , 0.        ],
       [0.13725652, 0.20389013, 3.6022848 ]])
>>> np.allclose(L_dd @ L_dd.T, A_dd)
True

We could also compute L_dd by applying a Cholesky factorization algorithm directly to A_dd (which is however less efficient)

>>> L_dd_cho = cholupdates.rank_1.downdate(L, v, method="cho_factor")
>>> L_dd_cho
array([[2.27628398, 0.        , 0.        ],
       [0.40711529, 2.8424243 , 0.        ],
       [0.13725652, 0.20389013, 3.6022848 ]])