eig_banded — SciPy v1.15.2 Manual (original) (raw)

scipy.linalg.

scipy.linalg.eig_banded(a_band, lower=False, eigvals_only=False, overwrite_a_band=False, select='a', select_range=None, max_ev=0, check_finite=True)[source]#

Solve real symmetric or complex Hermitian band matrix eigenvalue problem.

Find eigenvalues w and optionally right eigenvectors v of a:

a v[:,i] = w[i] v[:,i] v.H v = identity

The matrix a is stored in a_band either in lower diagonal or upper diagonal ordered form:

a_band[u + i - j, j] == a[i,j] (if upper form; i <= j) a_band[ i - j, j] == a[i,j] (if lower form; i >= j)

where u is the number of bands above the diagonal.

Example of a_band (shape of a is (6,6), u=2):

upper form:

lower form: a00 a11 a22 a33 a44 a55 a10 a21 a32 a43 a54 * a20 a31 a42 a53 * *

Cells marked with * are not used.

Parameters:

a_band(u+1, M) array_like

The bands of the M by M matrix a.

lowerbool, optional

Is the matrix in the lower form. (Default is upper form)

eigvals_onlybool, optional

Compute only the eigenvalues and no eigenvectors. (Default: calculate also eigenvectors)

overwrite_a_bandbool, optional

Discard data in a_band (may enhance performance)

select{‘a’, ‘v’, ‘i’}, optional

Which eigenvalues to calculate

select calculated
‘a’ All eigenvalues
‘v’ Eigenvalues in the interval (min, max]
‘i’ Eigenvalues with indices min <= i <= max

select_range(min, max), optional

Range of selected eigenvalues

max_evint, optional

For select==’v’, maximum number of eigenvalues expected. For other values of select, has no meaning.

In doubt, leave this parameter untouched.

check_finitebool, optional

Whether to check that the input matrix contains 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.

Returns:

w(M,) ndarray

The eigenvalues, in ascending order, each repeated according to its multiplicity.

v(M, M) float or complex ndarray

The normalized eigenvector corresponding to the eigenvalue w[i] is the column v[:,i]. Only returned if eigvals_only=False.

Raises:

LinAlgError

If eigenvalue computation does not converge.

See also

eigvals_banded

eigenvalues for symmetric/Hermitian band matrices

eig

eigenvalues and right eigenvectors of general arrays.

eigh

eigenvalues and right eigenvectors for symmetric/Hermitian arrays

eigh_tridiagonal

eigenvalues and right eigenvectors for symmetric/Hermitian tridiagonal matrices

Examples

import numpy as np from scipy.linalg import eig_banded A = np.array([[1, 5, 2, 0], [5, 2, 5, 2], [2, 5, 3, 5], [0, 2, 5, 4]]) Ab = np.array([[1, 2, 3, 4], [5, 5, 5, 0], [2, 2, 0, 0]]) w, v = eig_banded(Ab, lower=True) np.allclose(A @ v - v @ np.diag(w), np.zeros((4, 4))) True w = eig_banded(Ab, lower=True, eigvals_only=True) w array([-4.26200532, -2.22987175, 3.95222349, 12.53965359])

Request only the eigenvalues between [-3, 4]

w, v = eig_banded(Ab, lower=True, select='v', select_range=[-3, 4]) w array([-2.22987175, 3.95222349])