quimb.linalg.scipy_linalg ========================= .. py:module:: quimb.linalg.scipy_linalg .. autoapi-nested-parse:: Scipy based linear algebra. Attributes ---------- .. autoapisummary:: quimb.linalg.scipy_linalg.eigs_primme quimb.linalg.scipy_linalg.svds_primme Functions --------- .. autoapisummary:: quimb.linalg.scipy_linalg.maybe_sort_and_project quimb.linalg.scipy_linalg.eigs_scipy quimb.linalg.scipy_linalg.eigs_lobpcg quimb.linalg.scipy_linalg.svds_scipy Module Contents --------------- .. py:function:: maybe_sort_and_project(lk, vk, P, sort=True) .. py:function:: eigs_scipy(A, k, *, B=None, which=None, return_vecs=True, sigma=None, isherm=True, sort=True, P=None, tol=None, backend=None, **eigs_opts) Returns a few eigenpairs from a possibly sparse hermitian operator :param A: The operator to solve for. :type A: array_like, sparse_matrix, LinearOperator or quimb.Lazy :param k: Number of eigenpairs to return :type k: int :param B: If given, the RHS operator (which should be positive) defining a generalized eigen problem. :type B: array_like, sparse_matrix, LinearOperator or quimb.Lazy, optional :param which: where in spectrum to take eigenvalues from (see :func:`scipy.sparse.linalg.eigsh`). :type which: str, optional :param return_vecs: Whether to return the eigenvectors as well. :type return_vecs: bool, optional :param sigma: Shift, if targeting interior eigenpairs. :type sigma: float, optional :param isherm: Whether ``A`` is hermitian. :type isherm: bool, optional :param P: Perform the eigensolve in the subspace defined by this projector. :type P: array_like, sparse_matrix, LinearOperator or quimb.Lazy, optional :param sort: Whether to ensure the eigenvalues are sorted in ascending value. :type sort: bool, optional :param backend: Which backend to use. :type backend: None or 'primme', optional :param eigs_opts: Supplied to :func:`scipy.sparse.linalg.eigsh` or :func:`scipy.sparse.linalg.eigs`. :returns: * **lk** (*(k,) array*) -- The eigenvalues. * **vk** (*(m, k) array*) -- Corresponding eigenvectors (if ``return_vecs=True``). .. py:function:: eigs_lobpcg(A, k, *, B=None, v0=None, which=None, return_vecs=True, sigma=None, isherm=True, P=None, sort=True, **lobpcg_opts) Interface to scipy's lobpcg eigensolver, which can be good for generalized eigenproblems with matrix-free operators. Seems to a be a bit innacurate though (e.g. on the order of ~ 1e-6 for eigenvalues). Also only takes real, symmetric problems, targeting smallest eigenvalues (though scipy will soon have complex support, and its easy to add oneself). Note that the slepc eigensolver also has a lobpcg backend (``EPSType='lobpcg'``) which accepts complex input and is more accurate - though seems slower. :param A: The operator to solve for. :type A: array_like, sparse_matrix, LinearOperator or callable :param k: Number of eigenpairs to return :type k: int :param B: If given, the RHS operator (which should be positive) defining a generalized eigen problem. :type B: array_like, sparse_matrix, LinearOperator or callable, optional :param v0: The initial subspace to iterate with. :type v0: array_like (d, k), optional :param which: Find the smallest or largest eigenvalues. :type which: {'SA', 'LA'}, optional :param return_vecs: Whether to return the eigenvectors found. :type return_vecs: bool, optional :param P: Perform the eigensolve in the subspace defined by this projector. :type P: array_like, sparse_matrix, LinearOperator or callable, optional :param sort: Whether to ensure the eigenvalues are sorted in ascending value. :type sort: bool, optional :param lobpcg_opts: Supplied to :func:`scipy.sparse.linagl.lobpcg`. :returns: * **lk** (*array_like (k,)*) -- The eigenvalues. * **vk** (*array_like (d, k)*) -- The eigenvectors, if `return_vecs=True`. .. seealso:: :py:obj:`eigs_scipy`, :py:obj:`eigs_numpy`, :py:obj:`eigs_slepc` .. py:function:: svds_scipy(A, k=6, *, return_vecs=True, backend=None, **svds_opts) Compute a number of singular value pairs :param A: The operator to solve. :type A: (m, n) dense, sparse or linear operator. :param k: Number of requested singular values. :type k: int :param return_vecs: Whether to return the singular vectors. :type return_vecs: bool, optional :returns: * **U** (*(m, k) array*) -- Left singular vectors (if ``return_vecs=True``) as columns. * **s** (*(k,) array*) -- Singular values. * **VH** (*(k, n) array*) -- Right singular vectors (if ``return_vecs=True``) as rows. .. py:data:: eigs_primme .. py:data:: svds_primme