quimb.linalg.autoblock ====================== .. py:module:: quimb.linalg.autoblock .. autoapi-nested-parse:: Numba accelerated functions for finding charge sectors and subselecting submatrices. Functions --------- .. autoapisummary:: quimb.linalg.autoblock.get_nz quimb.linalg.autoblock.compute_blocks quimb.linalg.autoblock.subselect quimb.linalg.autoblock.subselect_set quimb.linalg.autoblock._eigh_autoblocked quimb.linalg.autoblock._eigvalsh_autoblocked quimb.linalg.autoblock.eigensystem_autoblocked Module Contents --------------- .. py:function:: get_nz(A) .. py:function:: compute_blocks(ix, jx, d) Find the charge sectors (blocks in matrix terms) given element coordinates ``ix`` and ``jx`` and total size ``d``. :param ix: The row coordinates of non-zero elements. :type ix: array of int :param jx: The column coordinates of non-zero elements. :type jx: array of int :param d: The total size of the operator. :type d: int :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)``. :rtype: list[list[int]] .. rubric:: 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]] .. py:function:: subselect(A, p) Select only the intersection of rows and columns of ``A`` matching the basis indices ``p``. Faster than double numpy slicing. :param A: Dense matrix to select from. :type A: 2D-array :param p: The basis indices to select. :type p: sequence of int :returns: **B** -- The matrix, of size ``(len(p), len(p))``. :rtype: 2D-array .. rubric:: 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]]) .. py:function:: subselect_set(A, B, p) Set only the intersection of rows and colums of ``A`` matching the basis indices ``p`` to ``B``. :param A: The matrix to set elements in. :type A: array with shape (d, d) :param B: The matrix to set elements from. :type B: array with shape (dp, dp) :param p: The basis indices. :type p: sequence of size dp .. rubric:: 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]]) .. py:function:: _eigh_autoblocked(A, sort=True) .. py:function:: _eigvalsh_autoblocked(A, sort=True) .. py:function:: eigensystem_autoblocked(A, sort=True, return_vecs=True, isherm=True) 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``. :param A: The operator to eigen-decompose. :type A: array_like :param sort: Whether to sort into ascending order, default True. :type sort: bool, optional :param isherm: Whether ``A`` is hermitian, default True. :type isherm: bool, optional :param return_vecs: Whether to return the eigenvectors, default True. :type return_vecs: bool, optional :returns: * **evals** (*1D-array*) -- The eigenvalues. * **evecs** (*qarray*) -- If ``return_vecs=True``, the eigenvectors.