symmray.linalg¶
Interface to linear algebra submodule functions.
Functions¶
|
Hermitian eigen-decomposition of an assumed hermitian symmray array. |
|
Truncated hermitian eigen-decomposition of an assumed hermitian symmray |
|
|
|
QR decomposition of a symmray array. |
|
Stabilized QR decomposition of a symmray array, returning results in an |
|
LQ decomposition of a symmray array. |
|
Stabilized LQ decomposition of a symmray array, returning results in an |
|
|
|
Singular value decomposition of a symmray array, returning all singular |
|
Truncated singular value decomposition of a symmray array. |
|
Truncated singular value decomposition of a symmray array, using |
|
Truncated singular value decomposition of a symmray array, using |
|
Cholesky decomposition of an assumed positive-definite symmray array. |
|
Cholesky decomposition with optional diagonal regularization, |
|
LQ decomposition via Cholesky factorization of |
|
QR decomposition via Cholesky factorization of |
Module Contents¶
- symmray.linalg.eigh(x: symmray.array_common.ArrayCommon, *args, **kwargs)[source]¶
Hermitian eigen-decomposition of an assumed hermitian symmray array.
- Returns:
w (VectorCommon) – The eigenvalues as a vector.
u (ArrayCommon) – The array of eigenvectors.
- symmray.linalg.eigh_truncated(x: symmray.array_common.ArrayCommon, *args, **kwargs)[source]¶
Truncated hermitian eigen-decomposition of an assumed hermitian symmray array.
- Parameters:
cutoff (float, optional) – Absolute eigenvalue cutoff threshold.
cutoff_mode (int or str, optional) –
How to perform the truncation:
1 or ‘abs’: trim values below
cutoff2 or ‘rel’: trim values below
s[0] * cutoff3 or ‘sum2’: trim s.t.
sum(s_trim**2) < cutoff.4 or ‘rsum2’: trim s.t.
sum(s_trim**2) < sum(s**2) * cutoff.5 or ‘sum1’: trim s.t.
sum(s_trim**1) < cutoff.6 or ‘rsum1’: trim s.t.
sum(s_trim**1) < sum(s**1) * cutoff.
max_bond (int) – An explicit maximum bond dimension, use -1 for none.
absorb ({-1, 0, 1, None}) –
How to absorb the eigenvalues.
-1 or ‘left’: absorb into the left factor (U).
0 or ‘both’: absorb the square root into both factors.
1 or ‘right’: absorb into the right factor (VH).
None: do not absorb, return eigenvalues as a BlockVector.
renorm ({0, 1}) – Whether to renormalize the eigenvalues (depends on cutoff_mode).
- Returns:
u (ArrayCommon) – The abelian array of left eigenvectors.
w (VectorCommon or None) – The vector of eigenvalues, or None if absorbed.
uh (ArrayCommon) – The abelian array of right eigenvectors.
- symmray.linalg.norm(x: symmray.array_common.ArrayCommon, *args, **kwargs)[source]¶
- symmray.linalg.qr(x: symmray.array_common.ArrayCommon, *args, **kwargs)[source]¶
QR decomposition of a symmray array.
- symmray.linalg.qr_stabilized(x: symmray.array_common.ArrayCommon, *args, **kwargs)[source]¶
Stabilized QR decomposition of a symmray array, returning results in an SVD-like
(Q, None, R)format for compatibility with tensor network split drivers.
- symmray.linalg.lq(x: symmray.array_common.ArrayCommon, *args, **kwargs)[source]¶
LQ decomposition of a symmray array.
- symmray.linalg.lq_stabilized(x: symmray.array_common.ArrayCommon, *args, **kwargs)[source]¶
Stabilized LQ decomposition of a symmray array, returning results in an SVD-like
(L, None, Q)format for compatibility with tensor network split drivers.
- symmray.linalg.solve(x: symmray.array_common.ArrayCommon, *args, **kwargs)[source]¶
- symmray.linalg.svd(x: symmray.array_common.ArrayCommon, *args, **kwargs) tuple[symmray.array_common.ArrayCommon, symmray.vector_common.VectorCommon, symmray.array_common.ArrayCommon][source]¶
Singular value decomposition of a symmray array, returning all singular values and vectors without truncation.
- symmray.linalg.svd_truncated(x: symmray.array_common.ArrayCommon, *args, **kwargs)[source]¶
Truncated singular value decomposition of a symmray array.
- Parameters:
cutoff (float, optional) – Singular value cutoff threshold.
cutoff_mode (int or str, optional) –
How to perform the truncation:
1 or ‘abs’: trim values below
cutoff2 or ‘rel’: trim values below
s[0] * cutoff3 or ‘sum2’: trim s.t.
sum(s_trim**2) < cutoff.4 or ‘rsum2’: trim s.t.
sum(s_trim**2) < sum(s**2) * cutoff.5 or ‘sum1’: trim s.t.
sum(s_trim**1) < cutoff.6 or ‘rsum1’: trim s.t.
sum(s_trim**1) < sum(s**1) * cutoff.
max_bond (int) – An explicit maximum bond dimension, use -1 for none.
absorb ({-1, 0, 1, None}) –
How to absorb the singular values.
-1 or ‘left’: absorb into the left factor (U).
0 or ‘both’: absorb the square root into both factors.
1 or ‘right’: absorb into the right factor (VH).
None: do not absorb, return singular values as a BlockVector.
renorm ({0, 1}) – Whether to renormalize the singular values (depends on cutoff_mode).
- Returns:
u (ArrayCommon) – The abelian array of left singular vectors.
s (VectorCommon or None) – The vector of singular values, or None if absorbed.
vh (ArrayCommon) – The abelian array of right singular vectors.
- symmray.linalg.svd_rand_truncated(x: symmray.array_common.ArrayCommon, *args, **kwargs)[source]¶
Truncated singular value decomposition of a symmray array, using randomized projection. This is efficient for low-rank approximations when a target
max_bondis known.- Parameters:
max_bond (int) – Target rank / maximum bond dimension.
absorb ({-1, 0, 1, None}) –
How to absorb the singular values.
-1 or ‘left’: absorb into the left factor (U).
0 or ‘both’: absorb the square root into both factors.
1 or ‘right’: absorb into the right factor (VH).
None: do not absorb, return singular values.
oversample (int, optional) – Extra sketch dimensions for accuracy. Default is 10.
num_iterations (int, optional) – Number of power iterations for accuracy. Default is 2.
seed (int, Generator or None, optional) – Random seed or generator for reproducibility.
- Returns:
u (ArrayCommon or None) – The array of left singular vectors.
s (VectorCommon or None) – The singular values, or None if absorbed.
vh (ArrayCommon or None) – The array of right singular vectors.
- symmray.linalg.svd_via_eig_truncated(x: symmray.array_common.ArrayCommon, *args, **kwargs)[source]¶
Truncated singular value decomposition of a symmray array, using eigen-decomposition of the gram (xdag @ x or x @ xdag) matrix. This can be faster, but also incur a loss of precision due to the squaring.
- Parameters:
cutoff (float, optional) – Singular value cutoff threshold.
cutoff_mode (int or str, optional) –
How to perform the truncation:
1 or ‘abs’: trim values below
cutoff2 or ‘rel’: trim values below
s[0] * cutoff3 or ‘sum2’: trim s.t.
sum(s_trim**2) < cutoff.4 or ‘rsum2’: trim s.t.
sum(s_trim**2) < sum(s**2) * cutoff.5 or ‘sum1’: trim s.t.
sum(s_trim**1) < cutoff.6 or ‘rsum1’: trim s.t.
sum(s_trim**1) < sum(s**1) * cutoff.
max_bond (int) – An explicit maximum bond dimension, use -1 for none.
absorb ({-1, 0, 1, None}) –
How to absorb the singular values.
-1 or ‘left’: absorb into the left factor (U).
0 or ‘both’: absorb the square root into both factors.
1 or ‘right’: absorb into the right factor (VH).
None: do not absorb, return singular values as a BlockVector.
renorm ({0, 1}) – Whether to renormalize the singular values (depends on cutoff_mode).
- Returns:
u (ArrayCommon) – The abelian array of left singular vectors.
s (VectorCommon or None) – The vector of singular values, or None if absorbed.
vh (ArrayCommon) – The abelian array of right singular vectors.
- symmray.linalg.cholesky(x: symmray.array_common.ArrayCommon, *args, **kwargs)[source]¶
Cholesky decomposition of an assumed positive-definite symmray array.
- Parameters:
x (ArrayCommon) – The 2D block-symmetric array to decompose.
upper (bool, optional) – Whether to return the upper triangular Cholesky factor. Default is False, returning the lower triangular factor.
- Returns:
l_or_r – The Cholesky factor. Lower triangular if
upper=False, upper triangular ifupper=True.- Return type:
- symmray.linalg.cholesky_regularized(x: symmray.array_common.ArrayCommon, *args, **kwargs)[source]¶
Cholesky decomposition with optional diagonal regularization, returning results in an SVD-like
(left, None, right)format for compatibility with tensor network split drivers.- Parameters:
x (ArrayCommon) – The 2D block-symmetric array to decompose. Must be positive (semi-)definite.
absorb ({-12, 0, 12}, optional) –
How to return the factors:
0('both'): return(L, None, L^H).-12('lsqrt'): return(L, None, None).12('rsqrt'): return(None, None, L^H).
shift (float, optional) – Diagonal regularization shift. If True or negative, auto-compute from dtype machine epsilon. The shift is always applied as a relative shift scaled by the trace of each block. Default is True.
- Returns:
left (ArrayCommon or None) – The lower Cholesky factor, or None.
s (None) – Always None (no singular values).
right (ArrayCommon or None) – The conjugate transpose of the Cholesky factor, or None.
- symmray.linalg.lq_via_cholesky(x: symmray.array_common.ArrayCommon, *args, **kwargs)[source]¶
LQ decomposition via Cholesky factorization of
x @ x^H.Computes
x = L @ QwhereLis lower triangular andQis isometric.- Parameters:
x (ArrayCommon) – The 2D block-symmetric array to decompose.
absorb ({-1, -10, -11} or str, optional) –
How to return the factors:
-1('left'): return(L, None, Q).-10('lfactor'): return(L, None, None).-11('rorthog'): return(None, None, Q).
shift (float, optional) – Diagonal regularization shift. If True or negative, auto-compute from dtype machine epsilon. The shift is always applied as a relative shift scaled by the trace of each block. Default is True.
solve_triangular (bool, optional) – Whether to use triangular solve (faster) or general solve to compute Q. Default is True.
- Returns:
L (ArrayCommon or None) – The lower triangular factor.
s (None) – Always None.
Q (ArrayCommon or None) – The isometric factor.
- symmray.linalg.qr_via_cholesky(x: symmray.array_common.ArrayCommon, *args, **kwargs)[source]¶
QR decomposition via Cholesky factorization of
x^H @ x.Computes
x = Q @ RwhereQis isometric andRis upper triangular.- Parameters:
x (ArrayCommon) – The 2D block-symmetric array to decompose.
absorb ({1, 11, 10} or str, optional) –
How to return the factors:
1('right'): return(Q, None, R).11('rfactor'): return(None, None, R).10('lorthog'): return(Q, None, None).
shift (float, optional) – Diagonal regularization shift. If True or negative, auto-compute from dtype machine epsilon. The shift is always applied as a relative shift scaled by the trace of each block. Default is True.
solve_triangular (bool, optional) – Whether to use triangular solve (faster) or general solve. Default is True.
- Returns:
Q (ArrayCommon or None) – The isometric factor.
s (None) – Always None.
R (ArrayCommon or None) – The upper triangular factor.