Previous topic

numpy.linalg.eig

Next topic

numpy.linalg.eigvals

numpy.linalg.eigh

numpy.linalg.eigh(a, UPLO='L')[source]

Return the eigenvalues and eigenvectors of a Hermitian or symmetric matrix.

Returns two objects, a 1-D array containing the eigenvalues of a, and a 2-D square array or matrix (depending on the input type) of the corresponding eigenvectors (in columns).

Parameters:
a : (…, M, M) array

Hermitian/Symmetric matrices whose eigenvalues and eigenvectors are to be computed.

UPLO : {‘L’, ‘U’}, optional

Specifies whether the calculation is done with the lower triangular part of a (‘L’, default) or the upper triangular part (‘U’). Irrespective of this value only the real parts of the diagonal will be considered in the computation to preserve the notion of a Hermitian matrix. It therefore follows that the imaginary part of the diagonal will always be treated as zero.

Returns:
w : (…, M) ndarray

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

v : {(…, M, M) ndarray, (…, M, M) matrix}

The column v[:, i] is the normalized eigenvector corresponding to the eigenvalue w[i]. Will return a matrix object if a is a matrix object.

Raises:
LinAlgError

If the eigenvalue computation does not converge.

See also

eigvalsh
eigenvalues of symmetric or Hermitian arrays.
eig
eigenvalues and right eigenvectors for non-symmetric arrays.
eigvals
eigenvalues of non-symmetric arrays.

Notes

New in version 1.8.0.

Broadcasting rules apply, see the numpy.linalg documentation for details.

The eigenvalues/eigenvectors are computed using LAPACK routines _syevd, _heevd

The eigenvalues of real symmetric or complex Hermitian matrices are always real. [1] The array v of (column) eigenvectors is unitary and a, w, and v satisfy the equations dot(a, v[:, i]) = w[i] * v[:, i].

References

[1](1, 2) G. Strang, Linear Algebra and Its Applications, 2nd Ed., Orlando, FL, Academic Press, Inc., 1980, pg. 222.

Examples

>>> from numpy import linalg as LA
>>> a = np.array([[1, -2j], [2j, 5]])
>>> a
array([[ 1.+0.j,  0.-2.j],
       [ 0.+2.j,  5.+0.j]])
>>> w, v = LA.eigh(a)
>>> w; v
array([ 0.17157288,  5.82842712])
array([[-0.92387953+0.j        , -0.38268343+0.j        ],
       [ 0.00000000+0.38268343j,  0.00000000-0.92387953j]])
>>> np.dot(a, v[:, 0]) - w[0] * v[:, 0] # verify 1st e-val/vec pair
array([2.77555756e-17 + 0.j, 0. + 1.38777878e-16j])
>>> np.dot(a, v[:, 1]) - w[1] * v[:, 1] # verify 2nd e-val/vec pair
array([ 0.+0.j,  0.+0.j])
>>> A = np.matrix(a) # what happens if input is a matrix object
>>> A
matrix([[ 1.+0.j,  0.-2.j],
        [ 0.+2.j,  5.+0.j]])
>>> w, v = LA.eigh(A)
>>> w; v
array([ 0.17157288,  5.82842712])
matrix([[-0.92387953+0.j        , -0.38268343+0.j        ],
        [ 0.00000000+0.38268343j,  0.00000000-0.92387953j]])
>>> # demonstrate the treatment of the imaginary part of the diagonal
>>> a = np.array([[5+2j, 9-2j], [0+2j, 2-1j]])
>>> a
array([[ 5.+2.j,  9.-2.j],
       [ 0.+2.j,  2.-1.j]])
>>> # with UPLO='L' this is numerically equivalent to using LA.eig() with:
>>> b = np.array([[5.+0.j, 0.-2.j], [0.+2.j, 2.-0.j]])
>>> b
array([[ 5.+0.j,  0.-2.j],
       [ 0.+2.j,  2.+0.j]])
>>> wa, va = LA.eigh(a)
>>> wb, vb = LA.eig(b)
>>> wa; wb
array([ 1.,  6.])
array([ 6.+0.j,  1.+0.j])
>>> va; vb
array([[-0.44721360-0.j        , -0.89442719+0.j        ],
       [ 0.00000000+0.89442719j,  0.00000000-0.4472136j ]])
array([[ 0.89442719+0.j       ,  0.00000000-0.4472136j],
       [ 0.00000000-0.4472136j,  0.89442719+0.j       ]])