scipy.linalg.eigvalsh¶
-
scipy.linalg.
eigvalsh
(a, b=None, lower=True, overwrite_a=False, overwrite_b=False, turbo=True, eigvals=None, type=1, check_finite=True, subset_by_index=None, subset_by_value=None, driver=None)[source]¶ Solves a standard or generalized eigenvalue problem for a complex Hermitian or real symmetric matrix.
Find eigenvalues array
w
of arraya
, whereb
is positive definite such that for every eigenvalue λ (i-th entry of w) and its eigenvector vi (i-th column of v) satisfies:a @ vi = λ * b @ vi vi.conj().T @ a @ vi = λ vi.conj().T @ b @ vi = 1
In the standard problem, b is assumed to be the identity matrix.
- Parameters
- a(M, M) array_like
A complex Hermitian or real symmetric matrix whose eigenvalues will be computed.
- b(M, M) array_like, optional
A complex Hermitian or real symmetric definite positive matrix in. If omitted, identity matrix is assumed.
- lowerbool, optional
Whether the pertinent array data is taken from the lower or upper triangle of a. (Default: lower)
- eigvals_onlybool, optional
Whether to calculate only eigenvalues and no eigenvectors. (Default: both are calculated)
- subset_by_indexiterable, optional
If provided, this two-element iterable defines the start and the end indices of the desired eigenvalues (ascending order and 0-indexed). To return only the second smallest to fifth smallest eigenvalues,
[1, 4]
is used.[n-3, n-1]
returns the largest three. Only available with “evr”, “evx”, and “gvx” drivers. The entries are directly converted to integers viaint()
.- subset_by_valueiterable, optional
If provided, this two-element iterable defines the half-open interval
(a, b]
that, if any, only the eigenvalues between these values are returned. Only available with “evr”, “evx”, and “gvx” drivers. Usenp.inf
for the unconstrained ends.- driver: str, optional
Defines which LAPACK driver should be used. Valid options are “ev”, “evd”, “evr”, “evx” for standard problems and “gv”, “gvd”, “gvx” for generalized (where b is not None) problems. See the Notes section of
scipy.linalg.eigh
.- typeint, optional
For the generalized problems, this keyword specifies the problem type to be solved for
w
andv
(only takes 1, 2, 3 as possible inputs):1 => a @ v = w @ b @ v 2 => a @ b @ v = w @ v 3 => b @ a @ v = w @ v
This keyword is ignored for standard problems.
- overwrite_abool, optional
Whether to overwrite data in
a
(may improve performance). Default is False.- overwrite_bbool, optional
Whether to overwrite data in
b
(may improve performance). Default is False.- check_finitebool, optional
Whether to check that the input matrices contain only finite numbers. Disabling may give a performance gain, but may result in problems (crashes, non-termination) if the inputs do contain infinities or NaNs.
- turbobool, optional
Deprecated by ``driver=gvd`` option. Has no significant effect for eigenvalue computations since no eigenvectors are requested.
..Deprecated in v1.5.0
- eigvalstuple (lo, hi), optional
Deprecated by ``subset_by_index`` keyword. Indexes of the smallest and largest (in ascending order) eigenvalues and corresponding eigenvectors to be returned: 0 <= lo <= hi <= M-1. If omitted, all eigenvalues and eigenvectors are returned.
- Returns
- w(N,) ndarray
The
N
(1<=N<=M
) selected eigenvalues, in ascending order, each repeated according to its multiplicity.
- Raises
- LinAlgError
If eigenvalue computation does not converge, an error occurred, or b matrix is not definite positive. Note that if input matrices are not symmetric or Hermitian, no error will be reported but results will be wrong.
See also
eigh
eigenvalues and right eigenvectors for symmetric/Hermitian arrays
eigvals
eigenvalues of general arrays
eigvals_banded
eigenvalues for symmetric/Hermitian band matrices
eigvalsh_tridiagonal
eigenvalues of symmetric/Hermitian tridiagonal matrices
Notes
This function does not check the input array for being Hermitian/symmetric in order to allow for representing arrays with only their upper/lower triangular parts.
This function serves as a one-liner shorthand for
scipy.linalg.eigh
with the optioneigvals_only=True
to get the eigenvalues and not the eigenvectors. Here it is kept as a legacy convenience. It might be beneficial to use the main function to have full control and to be a bit more pythonic.Examples
For more examples see
scipy.linalg.eigh
.>>> from scipy.linalg import eigvalsh >>> A = np.array([[6, 3, 1, 5], [3, 0, 5, 1], [1, 5, 6, 2], [5, 1, 2, 2]]) >>> w = eigvalsh(A) >>> w array([-3.74637491, -0.76263923, 6.08502336, 12.42399079])