quimb.linalg.autoblock¶
Numba accelerated functions for finding charge sectors and subselecting submatrices.
Functions¶
|
|
|
Find the charge sectors (blocks in matrix terms) given element |
|
Select only the intersection of rows and columns of |
|
Set only the intersection of rows and colums of |
|
|
|
|
|
Perform Hermitian eigen-decomposition, automatically identifying and |
Module Contents¶
- quimb.linalg.autoblock.compute_blocks(ix, jx, d)[source]¶
Find the charge sectors (blocks in matrix terms) given element coordinates
ix
andjx
and total sized
.- Parameters:
- Returns:
sectors – The list of charge sectors. Each element is itself a sorted list of the basis numbers that make up that sector. The permutation that would block diagonalize the operator is then
np.concatenate(sectors)
.- Return type:
Examples
>>> H = ham_hubbard_hardcore(4, sparse=True) >>> ix, jx = H.nonzero() >>> d = H.shape[0] >>> sectors = compute_blocks(ix, jx, d) >>> sectors [[0], [1, 2, 4, 8], [3, 5, 6, 9, 10, 12], [7, 11, 13, 14], [15]]
- quimb.linalg.autoblock.subselect(A, p)[source]¶
Select only the intersection of rows and columns of
A
matching the basis indicesp
. Faster than double numpy slicing.- Parameters:
A (2D-array) – Dense matrix to select from.
p (sequence of int) – The basis indices to select.
- Returns:
B – The matrix, of size
(len(p), len(p))
.- Return type:
2D-array
Examples
>>> A = np.arange(25).reshape(5, 5) >>> A array([[ 0, 1, 2, 3, 4], [ 5, 6, 7, 8, 9], [10, 11, 12, 13, 14], [15, 16, 17, 18, 19], [20, 21, 22, 23, 24]])
>>> subselect(A, [1, 3]) array([[ 6, 8], [16, 18]])
- quimb.linalg.autoblock.subselect_set(A, B, p)[source]¶
Set only the intersection of rows and colums of
A
matching the basis indicesp
toB
.- Parameters:
A (array with shape (d, d)) – The matrix to set elements in.
B (array with shape (dp, dp)) – The matrix to set elements from.
p (sequence of size dp) – The basis indices.
Examples
>>> A = np.zeros((5, 5)) >>> B = np.random.randn(3, 3) >>> p = [0, 2, 4] >>> subselect_set(A, B, p) array([[-0.31888218, 0. , 0.39293245, 0. , 0.21822712], [ 0. , 0. , 0. , 0. , 0. ], [ 0.66674486, 0. , 1.03388035, 0. , 1.7319345 ], [ 0. , 0. , 0. , 0. , 0. ], [-0.94542733, 0. , -0.37211882, 0. , 0.51951555]])
- quimb.linalg.autoblock.eigensystem_autoblocked(A, sort=True, return_vecs=True, isherm=True)[source]¶
Perform Hermitian eigen-decomposition, automatically identifying and exploiting symmetries appearing in the current basis as block diagonals formed via permutation of rows and columns. The whole process is accelerated using
numba
.- Parameters:
- Returns:
evals (1D-array) – The eigenvalues.
evecs (qarray) – If
return_vecs=True
, the eigenvectors.