Python API Reference#

Hyperbolic Householder transformations for Up- and Downdating Cholesky factorizations.

hyhound.update_cholesky(L: ndarray[dtype=float64, shape=(*, *), device='cpu', writable=False], A: ndarray[dtype=float64, shape=(*, *), device='cpu', writable=False]) tuple[numpy.ndarray[dtype=float64, shape=(*, *), order='F'], numpy.ndarray[dtype=float64, shape=(*, *), order='F'], numpy.ndarray[dtype=float64, shape=(*, *), order='F']]#
hyhound.update_cholesky(L: ndarray[dtype=float32, shape=(*, *), device='cpu', writable=False], A: ndarray[dtype=float32, shape=(*, *), device='cpu', writable=False]) tuple[numpy.ndarray[dtype=float32, shape=(*, *), order='F'], numpy.ndarray[dtype=float32, shape=(*, *), order='F'], numpy.ndarray[dtype=float32, shape=(*, *), order='F']]

Cholesky factorization update. Returns updated copies.

L̃L̃ᵀ + ÃÃᵀ = LLᵀ + AAᵀ

Parameters:
L(k × n), lower-trapezoidal

The original Cholesky factor.

A(k × m), rectangular

The update matrix.

Returns:
(k × n)

The updated Cholesky factor.

A_rem(k × m)

Contains the k-n bottom rows of the remaining update matrix Ã. The top n rows of à are zero (not stored explicitly). The top n rows of A_rem contain Householder reflectors.

W(r × n)

The upper triangular Householder representations generated during the Cholesky factorization update. Together with the top n rows of A_rem, this can be used to apply the block Householder transformation to other matrices. The number of rows depends on the block size and is architecture-dependent.

hyhound.downdate_cholesky(L: ndarray[dtype=float64, shape=(*, *), device='cpu', writable=False], A: ndarray[dtype=float64, shape=(*, *), device='cpu', writable=False]) tuple[numpy.ndarray[dtype=float64, shape=(*, *), order='F'], numpy.ndarray[dtype=float64, shape=(*, *), order='F'], numpy.ndarray[dtype=float64, shape=(*, *), order='F']]#
hyhound.downdate_cholesky(L: ndarray[dtype=float32, shape=(*, *), device='cpu', writable=False], A: ndarray[dtype=float32, shape=(*, *), device='cpu', writable=False]) tuple[numpy.ndarray[dtype=float32, shape=(*, *), order='F'], numpy.ndarray[dtype=float32, shape=(*, *), order='F'], numpy.ndarray[dtype=float32, shape=(*, *), order='F']]

Cholesky factorization downdate. Returns updated copies.

L̃L̃ᵀ - ÃÃᵀ = LLᵀ - AAᵀ

Parameters:
L(k × n), lower-trapezoidal

The original Cholesky factor.

A(k × m), rectangular

The downdate matrix.

Returns:
(k × n)

The updated Cholesky factor.

A_rem(k × m)

Contains the k-n bottom rows of the remaining downdate matrix Ã. The top n rows of à are zero (not stored explicitly). The top n rows of A_rem contain Householder reflectors.

W(r × n)

The upper triangular Householder representations generated during the Cholesky factorization update. Together with the top n rows of A_rem, this can be used to apply the block Householder transformation to other matrices. The number of rows depends on the block size and is architecture-dependent.

hyhound.update_cholesky_sign(L: ndarray[dtype=float64, shape=(*, *), device='cpu', writable=False], A: ndarray[dtype=float64, shape=(*, *), device='cpu', writable=False], signs: ndarray[dtype=float64, shape=(*), order='A', device='cpu', writable=False]) tuple[numpy.ndarray[dtype=float64, shape=(*, *), order='F'], numpy.ndarray[dtype=float64, shape=(*, *), order='F'], numpy.ndarray[dtype=float64, shape=(*, *), order='F']]#
hyhound.update_cholesky_sign(L: ndarray[dtype=float32, shape=(*, *), device='cpu', writable=False], A: ndarray[dtype=float32, shape=(*, *), device='cpu', writable=False], signs: ndarray[dtype=float32, shape=(*), order='A', device='cpu', writable=False]) tuple[numpy.ndarray[dtype=float32, shape=(*, *), order='F'], numpy.ndarray[dtype=float32, shape=(*, *), order='F'], numpy.ndarray[dtype=float32, shape=(*, *), order='F']]

Cholesky factorization update with signed columns. Returns updated copies.

L̃L̃ᵀ + ÃSÃᵀ = LLᵀ + ASAᵀ, where S = np.diag(np.copysign(np.ones(m), signs)) and signs contains ±0.

Parameters:
L(k × n), lower-trapezoidal

The original Cholesky factor.

A(k × m), rectangular

The update matrix.

signsm-vector

Signs that determine whether a column of A is added (+0) or removed (-0). Values other than ±0 are not allowed.

Returns:
(k × n)

The updated Cholesky factor.

A_rem(k × m)

Contains the k-n bottom rows of the remaining update matrix Ã. The top n rows of à are zero (not stored explicitly). The top n rows of A_rem contain Householder reflectors.

W(r × n)

The upper triangular Householder representations generated during the Cholesky factorization update. Together with the top n rows of A_rem, this can be used to apply the block Householder transformation to other matrices. The number of rows depends on the block size and is architecture-dependent.

hyhound.update_cholesky_diag(L: ndarray[dtype=float64, shape=(*, *), device='cpu', writable=False], A: ndarray[dtype=float64, shape=(*, *), device='cpu', writable=False], diag: ndarray[dtype=float64, shape=(*), order='A', device='cpu', writable=False]) tuple[numpy.ndarray[dtype=float64, shape=(*, *), order='F'], numpy.ndarray[dtype=float64, shape=(*, *), order='F'], numpy.ndarray[dtype=float64, shape=(*, *), order='F']]#
hyhound.update_cholesky_diag(L: ndarray[dtype=float32, shape=(*, *), device='cpu', writable=False], A: ndarray[dtype=float32, shape=(*, *), device='cpu', writable=False], diag: ndarray[dtype=float32, shape=(*), order='A', device='cpu', writable=False]) tuple[numpy.ndarray[dtype=float32, shape=(*, *), order='F'], numpy.ndarray[dtype=float32, shape=(*, *), order='F'], numpy.ndarray[dtype=float32, shape=(*, *), order='F']]

Cholesky factorization update with diagonal scaling. Returns updated copies.

L̃L̃ᵀ + ÃDÃᵀ = LLᵀ + ADAᵀ, where D = np.diag(diag).

Parameters:
L(k × n), lower-trapezoidal

The original Cholesky factor.

A(k × m), rectangular

The update matrix.

diagm-vector

Scale factors corresponding to the columns of A.

Returns:
(k × n)

The updated Cholesky factor.

A_rem(k × m)

Contains the k-n bottom rows of the remaining update matrix Ã. The top n rows of à are zero (not stored explicitly). The top n rows of A_rem contain Householder reflectors.

W(r × n)

The upper triangular Householder representations generated during the Cholesky factorization update. Together with the top n rows of A_rem, this can be used to apply the block Householder transformation to other matrices. The number of rows depends on the block size and is architecture-dependent.

hyhound.update_cholesky_inplace(L: ndarray[dtype=float64, shape=(*, *), order='F', device='cpu'], A: ndarray[dtype=float64, shape=(*, *), order='F', device='cpu']) None#
hyhound.update_cholesky_inplace(L: ndarray[dtype=float32, shape=(*, *), order='F', device='cpu'], A: ndarray[dtype=float32, shape=(*, *), order='F', device='cpu']) None

Cholesky factorization update. Overwrites its arguments.

L̃L̃ᵀ + ÃÃᵀ = LLᵀ + AAᵀ

Parameters:
L(k × n), lower-trapezoidal, Fortran order

On entry, the original Cholesky factor L. On exit, contains the updated Cholesky factor L̃.

A(k × m), rectangular, Fortran order

On entry, the update matrix A. On exit, contains the k-n bottom rows of the remaining update matrix à (the top n rows of à are zero).

hyhound.downdate_cholesky_inplace(L: ndarray[dtype=float64, shape=(*, *), order='F', device='cpu'], A: ndarray[dtype=float64, shape=(*, *), order='F', device='cpu']) None#
hyhound.downdate_cholesky_inplace(L: ndarray[dtype=float32, shape=(*, *), order='F', device='cpu'], A: ndarray[dtype=float32, shape=(*, *), order='F', device='cpu']) None

Cholesky factorization downdate. Overwrites its arguments.

L̃L̃ᵀ - ÃÃᵀ = LLᵀ - AAᵀ

Parameters:
L(k × n), lower-trapezoidal, Fortran order

On entry, the original Cholesky factor L. On exit, contains the updated Cholesky factor L̃.

A(k × m), rectangular, Fortran order

On entry, the downdate matrix A. On exit, contains the k-n bottom rows of the remaining downdate matrix à (the top n rows of à are zero).

hyhound.update_cholesky_sign_inplace(L: ndarray[dtype=float64, shape=(*, *), order='F', device='cpu'], A: ndarray[dtype=float64, shape=(*, *), order='F', device='cpu'], signs: ndarray[dtype=float64, shape=(*), order='A', device='cpu', writable=False]) None#
hyhound.update_cholesky_sign_inplace(L: ndarray[dtype=float32, shape=(*, *), order='F', device='cpu'], A: ndarray[dtype=float32, shape=(*, *), order='F', device='cpu'], signs: ndarray[dtype=float32, shape=(*), order='A', device='cpu', writable=False]) None

Cholesky factorization update with signed columns. Overwrites its arguments.

L̃L̃ᵀ + ÃSÃᵀ = LLᵀ + ASAᵀ, where S = np.diag(np.copysign(np.ones(m), signs)) and signs contains ±0.

Parameters:
L(k × n), lower-trapezoidal, Fortran order

On entry, the original Cholesky factor L. On exit, contains the updated Cholesky factor L̃.

A(k × m), rectangular, Fortran order

On entry, the update matrix A. On exit, contains the k-n bottom rows of the remaining update matrix à (the top n rows of à are zero).

signsm-vector

Signs that determine whether a column of A is added (+0) or removed (-0). Values other than ±0 are not allowed.

hyhound.update_cholesky_diag_inplace(L: ndarray[dtype=float64, shape=(*, *), order='F', device='cpu'], A: ndarray[dtype=float64, shape=(*, *), order='F', device='cpu'], diag: ndarray[dtype=float64, shape=(*), order='A', device='cpu', writable=False]) None#
hyhound.update_cholesky_diag_inplace(L: ndarray[dtype=float32, shape=(*, *), order='F', device='cpu'], A: ndarray[dtype=float32, shape=(*, *), order='F', device='cpu'], diag: ndarray[dtype=float32, shape=(*), order='A', device='cpu', writable=False]) None

Cholesky factorization update with diagonal scaling. Overwrites its arguments.

L̃L̃ᵀ + ÃDÃᵀ = LLᵀ + ADAᵀ, where D = np.diag(diag).

Parameters:
L(k × n), lower-trapezoidal, Fortran order

On entry, the original Cholesky factor L. On exit, contains the updated Cholesky factor L̃.

A(k × m), rectangular, Fortran order

On entry, the update matrix A. On exit, contains the k-n bottom rows of the remaining update matrix à (the top n rows of à are zero).

diagm-vector

Scale factors corresponding to the columns of A.

hyhound.update_apply_householder(L: ndarray[dtype=float64, shape=(*, *), device='cpu', writable=False], A: ndarray[dtype=float64, shape=(*, *), device='cpu', writable=False], B: ndarray[dtype=float64, shape=(*, *), order='F', device='cpu', writable=False], W: ndarray[dtype=float64, shape=(*, *), order='F', device='cpu', writable=False]) tuple[numpy.ndarray[dtype=float64, shape=(*, *), order='F'], numpy.ndarray[dtype=float64, shape=(*, *), order='F']]#
hyhound.update_apply_householder(L: ndarray[dtype=float32, shape=(*, *), device='cpu', writable=False], A: ndarray[dtype=float32, shape=(*, *), device='cpu', writable=False], B: ndarray[dtype=float32, shape=(*, *), order='F', device='cpu', writable=False], W: ndarray[dtype=float32, shape=(*, *), order='F', device='cpu', writable=False]) tuple[numpy.ndarray[dtype=float32, shape=(*, *), order='F'], numpy.ndarray[dtype=float32, shape=(*, *), order='F']]

Apply a block Householder transformation generated during a Cholesky factorization update. Returns updated copies.

(L̃ Ã) = (L A) Q

where Q is the block Householder transformation represented by W and B.

Parameters:
L(l × n), rectangular

Matrix to apply the transformation to.

A(l × m), rectangular

Matrix to apply the transformation to.

B(k × m), rectangular

The Householder reflector vectors generated during the Cholesky factorization update.

W(r × n)

The upper triangular Householder representations generated during the Cholesky factorization update.

Returns:
(l × n)

The updated matrix L.

Ã(l × m)

The updated matrix A.

hyhound.downdate_apply_householder(L: ndarray[dtype=float64, shape=(*, *), device='cpu', writable=False], A: ndarray[dtype=float64, shape=(*, *), device='cpu', writable=False], B: ndarray[dtype=float64, shape=(*, *), order='F', device='cpu', writable=False], W: ndarray[dtype=float64, shape=(*, *), order='F', device='cpu', writable=False]) tuple[numpy.ndarray[dtype=float64, shape=(*, *), order='F'], numpy.ndarray[dtype=float64, shape=(*, *), order='F']]#
hyhound.downdate_apply_householder(L: ndarray[dtype=float32, shape=(*, *), device='cpu', writable=False], A: ndarray[dtype=float32, shape=(*, *), device='cpu', writable=False], B: ndarray[dtype=float32, shape=(*, *), order='F', device='cpu', writable=False], W: ndarray[dtype=float32, shape=(*, *), order='F', device='cpu', writable=False]) tuple[numpy.ndarray[dtype=float32, shape=(*, *), order='F'], numpy.ndarray[dtype=float32, shape=(*, *), order='F']]

Apply a block Householder transformation generated during a Cholesky factorization downdate. Returns updated copies.

(L̃ Ã) = (L A) Q

where Q is the block Householder transformation represented by W and B.

Parameters:
L(l × n), rectangular

Matrix to apply the transformation to.

A(l × m), rectangular

Matrix to apply the transformation to.

B(k × m), rectangular

The Householder reflector vectors generated during the Cholesky factorization downdate.

W(r × n)

The upper triangular Householder representations generated during the Cholesky factorization downdate.

Returns:
(l × n)

The updated matrix L.

Ã(l × m)

The updated matrix A.

hyhound.update_apply_householder_sign(L: ndarray[dtype=float64, shape=(*, *), device='cpu', writable=False], A: ndarray[dtype=float64, shape=(*, *), device='cpu', writable=False], signs: ndarray[dtype=float64, shape=(*), order='A', device='cpu', writable=False], B: ndarray[dtype=float64, shape=(*, *), order='F', device='cpu', writable=False], W: ndarray[dtype=float64, shape=(*, *), order='F', device='cpu', writable=False]) tuple[numpy.ndarray[dtype=float64, shape=(*, *), order='F'], numpy.ndarray[dtype=float64, shape=(*, *), order='F']]#
hyhound.update_apply_householder_sign(L: ndarray[dtype=float32, shape=(*, *), device='cpu', writable=False], A: ndarray[dtype=float32, shape=(*, *), device='cpu', writable=False], signs: ndarray[dtype=float32, shape=(*), order='A', device='cpu', writable=False], B: ndarray[dtype=float32, shape=(*, *), order='F', device='cpu', writable=False], W: ndarray[dtype=float32, shape=(*, *), order='F', device='cpu', writable=False]) tuple[numpy.ndarray[dtype=float32, shape=(*, *), order='F'], numpy.ndarray[dtype=float32, shape=(*, *), order='F']]

Apply a block Householder transformation generated during a Cholesky factorization update with signed columns. Returns updated copies.

(L̃ Ã) = (L A) Q

where Q is the block Householder transformation represented by W, B and signs.

Parameters:
L(l × n), rectangular

Matrix to apply the transformation to.

A(l × m), rectangular

Matrix to apply the transformation to.

signsm-vector

Signs that determine whether a column of A was added (+0) or removed (-0). Values other than ±0 are not allowed.

B(k × m), rectangular

The Householder reflector vectors generated during the Cholesky factorization update.

W(r × n)

The upper triangular Householder representations generated during the Cholesky factorization update.

Returns:
(l × n)

The updated matrix L.

Ã(l × m)

The updated matrix A.

hyhound.update_apply_householder_diag(L: ndarray[dtype=float64, shape=(*, *), device='cpu', writable=False], A: ndarray[dtype=float64, shape=(*, *), device='cpu', writable=False], diag: ndarray[dtype=float64, shape=(*), order='A', device='cpu', writable=False], B: ndarray[dtype=float64, shape=(*, *), order='F', device='cpu', writable=False], W: ndarray[dtype=float64, shape=(*, *), order='F', device='cpu', writable=False]) tuple[numpy.ndarray[dtype=float64, shape=(*, *), order='F'], numpy.ndarray[dtype=float64, shape=(*, *), order='F']]#
hyhound.update_apply_householder_diag(L: ndarray[dtype=float32, shape=(*, *), device='cpu', writable=False], A: ndarray[dtype=float32, shape=(*, *), device='cpu', writable=False], diag: ndarray[dtype=float32, shape=(*), order='A', device='cpu', writable=False], B: ndarray[dtype=float32, shape=(*, *), order='F', device='cpu', writable=False], W: ndarray[dtype=float32, shape=(*, *), order='F', device='cpu', writable=False]) tuple[numpy.ndarray[dtype=float32, shape=(*, *), order='F'], numpy.ndarray[dtype=float32, shape=(*, *), order='F']]

Apply a block Householder transformation generated during a Cholesky factorization update with diagonal scaling. Returns updated copies.

(L̃ Ã) = (L A) Q

where Q is the block Householder transformation represented by W, B and diag.

Parameters:
L(l × n), rectangular

Matrix to apply the transformation to.

A(l × m), rectangular

Matrix to apply the transformation to.

diagm-vector

Scale factors corresponding to the columns of A used when generating the Householder transformation.

B(k × m), rectangular

The Householder reflector vectors generated during the Cholesky factorization update.

W(r × n)

The upper triangular Householder representations generated during the Cholesky factorization update.

Returns:
(l × n)

The updated matrix L.

Ã(l × m)

The updated matrix A.