symmray.array_common¶
Methods that apply to all Abelian symmetric arrays, regardless of backend.
Classes¶
Common base class for all abelian arrays with symmetry, fermionic or |
Functions¶
|
|
|
|
|
Match |
|
Match |
|
Attach squeezed singleton axes to a neighbouring kept dimension. |
|
Convert target entries to unfuse/fuse/expand operation arguments. |
|
Given a current block sparse shape |
|
Return a tuple of the elements in |
|
Parse the axes argument for single integer and also negative indices. |
|
Decide whether to keep labels from operands when combining two arrays. |
Module Contents¶
- symmray.array_common._match_reshape_forward(shape, subshapes, newshape, i_start, i_stop, j_start, j_stop)[source]¶
Match
shape[i_start:i_stop]tonewshape[j_start:j_stop]from left to right, returning target entries and the first unconsumed input axis.
- symmray.array_common._match_reshape_reverse(shape, subshapes, newshape, i_start, i_stop, j_start, j_stop)[source]¶
Match
shape[i_start:i_stop]tonewshape[j_start:j_stop]from right to left, returning target entries and the first consumed input axis.
- symmray.array_common._absorb_squeeze_entries(entries)[source]¶
Attach squeezed singleton axes to a neighbouring kept dimension.
- symmray.array_common._entries_to_reshape_args(entries, ndim_old)[source]¶
Convert target entries to unfuse/fuse/expand operation arguments.
- symmray.array_common.calc_reshape_args(shape, newshape, subshapes)[source]¶
Given a current block sparse shape
shapea target shapenewshapeand current sub index sizessubshapes(i.e. previously fused dimensions) compute the sequence of axes to unfuse, fuse and expand to reshape the array.A single
-1entry innewshapeis interpreted structurally: it consumes whichever input axes are not matched by the rest of the target shape (fusing them if there is more than one). This avoids the division-by-size resolution which can be ambiguous for sparse arrays whose fused-axis sizes do not equal the product of their subshape.- Parameters:
- Returns:
axs_unfuse (tuple[int]) – The axes to unfuse.
axs_fuse (tuple[tuple[tuple[int]]]) – The axes (after unfusing) to fuse, grouped by contiguous groups.
axs_expand (tuple[int]) – The axes (after unfusing and fusing) to expand.
- class symmray.array_common.ArrayCommon[source]¶
Common base class for all abelian arrays with symmetry, fermionic or not, sparse or flat etc.
- property symmetry: symmray.symmetries.Symmetry¶
The symmetry object of the array.
- property indices: tuple[symmray.index_common.Index, Ellipsis]¶
The indices of the array.
- classmethod get_class_symmetry(symmetry=None) symmray.symmetries.Symmetry[source]¶
- property T: ArrayCommon¶
The transpose of the block array.
- _dagger_abelian(inplace=False) ArrayCommon[source]¶
Conjugate transpose or adjoint, implements the .H property.
- property H: ArrayCommon¶
- fuse(*axes_groups, expand_empty=True, inplace=False, **kwargs) ArrayCommon[source]¶
Fuse the given group or groups of axes. The new fused axes will be inserted at the minimum index of any fused axis (even if it is not in the first group). For example,
x.fuse([5, 3], [7, 2, 6])will produce an array with axes like:groups inserted at axis 2, removed beyond that. ......<-- (0, 1, g0, g1, 4, 8, ...) | | | g1=(7, 2, 6) g0=(5, 3)
The fused axes will carry subindex information, which can be used to automatically unfuse them back into their original components. Depending on expand_empty, any empty groups can be expanded to new singlet dimensions, or simply ignored.
- Parameters:
axes_groups (Sequence[Sequence[int]]) – The axes to fuse. Each group of axes will be fused into a single axis.
expand_empty (bool, optional) – Whether to expand empty groups into new axes.
mode ("auto", "insert", "concat", optional) – The method to use for fusing. “insert” creates the new fused blocks and insert the subblocks inplace. “concat” recursively concatenates the subblocks, which can be slightly slower but is more compatible with e.g. autodiff. “auto” will use “insert” if the backend is numpy, otherwise “concat”.
inplace (bool, optional) – Whether to perform the operation inplace or return a new array.
kwargs (dict, optional) – Additional keyword arguments to pass to the core fusing method.
- Return type:
- unfuse_all(inplace=False) ArrayCommon[source]¶
Unfuse all indices that carry subindex information, likely from a fusing operation.
- Parameters:
inplace (bool, optional) – Whether to perform the operation inplace or return a new array.
- Return type:
- reshape(newshape, inplace=False) ArrayCommon[source]¶
Reshape this abelian array to
newshape, assuming it can be done by any mix of fusing, unfusing, and expanding new axes.Restrictions and complications vs normal array reshaping arise from the fact that
only previously fused axes can be unfused, and their total size may not be the product of the individual sizes due to sparsity.
the indices carry information beyond size and how they are grouped potentially matters, relevant for singleton dimensions.
Accordingly the approach here is as follows:
Unfuse any axes that match the new shape.
If there are singleton dimensions that don’t appear in the new shape, (i.e. are being ‘squeezed’) these are grouped with the axis to the their left to then be fused. If they are already left-most, they are grouped with the right.
Fuse any groups of axes required to match the new shape. Adjacent groups are fused simultaneously for efficiency.
Expand new axes required to match singlet dimensions in the new shape. By default these will have zero charge and dual-ness iherited from whichever axis is to their left, or right if they are the left-most axis already.
To avoid the effective grouping of ‘squeezed’ axes you can explicitly squeeze them before reshaping. Similarly use
expand_dimsto add new axes with specific charges and dual-ness.newshapemay contain at most one-1entry, which is resolved structurally: it absorbs whichever input axes are not matched by the rest of the target shape (fusing them if there is more than one). This differs from numpy’s division-based resolution and avoids ambiguity when fused-index sizes do not equal the product of their subshape due to sparsity.- Parameters:
- Return type:
See also
fuse,unfuse,squeeze,expand_dims
- eigh_truncated(cutoff=-1.0, cutoff_mode='rsum2', max_bond=-1, absorb=0, renorm=False, positive=False, **kwargs) tuple[ArrayCommon, symmray.vector_common.VectorCommon, ArrayCommon][source]¶
Truncated hermitian eigen-decomposition of this assumed hermitian 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 VectorCommon.
renorm ({0, 1}) – Whether to renormalize the eigenvalues (depends on
cutoff_mode).positive (bool, optional) – If True, assume all eigenvalues are positive for a faster sort. By default False.
- Returns:
u (ArrayCommon) – The array of left eigenvectors.
w (VectorCommon or None) – The vector of eigenvalues, or None if absorbed.
uh (ArrayCommon) – The array of right eigenvectors.
- svd(**kwargs) tuple[ArrayCommon, symmray.vector_common.VectorCommon, ArrayCommon][source]¶
Singular value decomposition of this array.
- Returns:
u (ArrayCommon) – The left singular vectors.
s (VectorCommon) – The singular values.
vh (ArrayCommon) – The right singular vectors (hermitian transposed).
- svd_truncated(cutoff=0.0, cutoff_mode='rsum2', max_bond=None, absorb='both', renorm=False, **kwargs) tuple[ArrayCommon, symmray.vector_common.VectorCommon, ArrayCommon][source]¶
Truncated singular value decomposition of this 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 VectorCommon.
renorm ({0, 1}) – Whether to renormalize the singular values (depends on cutoff_mode).
- Returns:
u (ArrayCommon or None) – The abelian array of left singular vectors.
s (VectorCommon or None) – The vector of singular values, or None if absorbed.
vh (ArrayCommon or None) – The abelian array of right singular vectors.
- svd_via_eig(**kwargs)[source]¶
Singular value decomposition of this array, using eigen-decomposition of the Gram matrix (xdag @ x or x @ xdag). This can be faster, but also incur a loss of precision due to the squaring.
- svd_via_eig_truncated(cutoff=0.0, cutoff_mode='rsum2', max_bond=None, absorb='both', renorm=False, **kwargs)[source]¶
Truncated singular value decomposition of this array, using eigen-decomposition of the gram matrix (xdag @ x or x @ xdag). 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 VectorCommon.
renorm ({0, 1}) – Whether to renormalize singular values (depends on cutoff_mode).
- Returns:
u (ArrayCommon or None) – The abelian array of left singular vectors.
s (VectorCommon or None) – The vector of singular values, or None if absorbed.
vh (ArrayCommon or None) – The abelian array of right singular vectors.
- svd_rand_truncated(max_bond, absorb=0, oversample=10, num_iterations=2, seed=None, **kwargs) tuple[ArrayCommon, symmray.vector_common.VectorCommon, ArrayCommon][source]¶
Truncated singular value decomposition of this 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.
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 abelian array of left singular vectors.
s (VectorCommon or None) – The singular values, or None if absorbed.
vh (ArrayCommon or None) – The abelian array of right singular vectors.
- qr(absorb='right', stabilized=False, **kwargs) tuple[ArrayCommon, ArrayCommon][source]¶
QR decomposition.
- Parameters:
x (ArrayCommon) – The array to decompose.
absorb (int or str, optional) –
Absorption mode:
Absorb.U_sVH(1, ‘right’): return(Q, None, R).Absorb.sVH(11, ‘rfactor’): return(None, None, R).Absorb.U(10, ‘lorthog’): return(Q, None, None).
stabilized (bool, optional) – Whether to use a stabilized QR decomposition, that is, ensure positive diagonal elements in the R factor. Default is False.
- Returns:
q (ArrayCommon) – The orthogonal matrix.
r (ArrayCommon) – The right factor matrix.
- qr_via_cholesky(absorb='right', shift=True, solve_triangular=True, **kwargs) tuple[ArrayCommon, None, ArrayCommon][source]¶
QR decomposition of a this array via Cholesky factorization.
By default (absorb=”right”) computes
x = Q @ RwhereQis isometric andRis upper triangular. LQ-like decompositions can be obtained with absorb=”left” for example.- Parameters:
absorb (int or str, optional) –
Absorption mode:
Absorb.U_sVH(1, ‘right’): return(Q, None, R).Absorb.sVH(11, ‘rfactor’): return(None, None, R).Absorb.U(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.
- lq(absorb='left', stabilized=False, **kwargs) tuple[ArrayCommon, ArrayCommon][source]¶
LQ decomposition.
- Parameters:
stabilized (bool, optional) – Whether to use a stabilized LQ decomposition, that is, with positive diagonal elements in the L factor. Default is False.
absorb (str or int, optional) –
Which factors to return in the output, by default ‘Us_VH’ (both). Options are:
”left” or “Us_VH”: return (L, Q)
”lfactor” or “Us”: return (L, None)
”rorthog” or “VH”: return (None, Q)
- Returns:
l (ArrayCommon) – The lower triangular matrix.
q (ArrayCommon) – The orthogonal matrix.
- lq_via_cholesky(absorb='left', shift=True, solve_triangular=True, **kwargs) tuple[ArrayCommon, None, ArrayCommon][source]¶
LQ decomposition of a this array via Cholesky factorization of the Gram matrix
x @ x^H.By default (absorb=”left”) computes
x = L @ QwhereLis lower triangular andQis isometric.- Parameters:
absorb (int or str, optional) –
Absorption mode:
Absorb.Us_VH(-1, ‘left’): return(L, None, Q).Absorb.Us(-10, ‘lfactor’): return(L, None, None).Absorb.VH(-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.array_common.without(it, remove)[source]¶
Return a tuple of the elements in
itwith those at positionsremoveremoved.