hsein#

Functions

void shsein(
    const char*          side,
    const char*          eigsrc,
    const char*          initv,
          INT*  restrict select,
    const INT            n,
    const f32*  restrict H,
    const INT            ldh,
          f32*  restrict wr,
    const f32*  restrict wi,
          f32*  restrict VL,
    const INT            ldvl,
          f32*  restrict VR,
    const INT            ldvr,
    const INT            mm,
          INT*           m,
          f32*  restrict work,
          INT*  restrict ifaill,
          INT*  restrict ifailr,
          INT*           info
);
void shsein(const char *side, const char *eigsrc, const char *initv, INT *restrict select, const INT n, const f32 *restrict H, const INT ldh, f32 *restrict wr, const f32 *restrict wi, f32 *restrict VL, const INT ldvl, f32 *restrict VR, const INT ldvr, const INT mm, INT *m, f32 *restrict work, INT *restrict ifaill, INT *restrict ifailr, INT *info)#

SHSEIN uses inverse iteration to find specified right and/or left eigenvectors of a real upper Hessenberg matrix H.

The right eigenvector x and the left eigenvector y of the matrix H corresponding to an eigenvalue w are defined by:

         H * x = w * x,     y**h * H = w * y**h
where y**h denotes the conjugate transpose of the vector y.

Each eigenvector is normalized so that the element of largest magnitude has magnitude 1; here the magnitude of a complex number (x,y) is taken to be |x|+|y|.

Parameters

in
side

‘R’: compute right eigenvectors only; ‘L’: compute left eigenvectors only; ‘B’: compute both right and left eigenvectors.

in
eigsrc

‘Q’: eigenvalues were found using SHSEQR, so matrix splitting can be exploited; ‘N’: no assumptions on correspondence between eigenvalues and diagonal blocks.

in
initv

‘N’: no initial vectors supplied; ‘U’: user-supplied initial vectors in VL and/or VR.

inout
select

Array of dimension n. Specifies which eigenvectors to compute. select[j] = 1 (nonzero) selects the eigenvalue wr[j] + i*wi[j].

in
n

The order of the matrix H (n >= 0).

in
H

Upper Hessenberg matrix H. Array of dimension (ldh, n).

in
ldh

The leading dimension of H (ldh >= max(1,n)).

inout
wr

Array of dimension n. Real parts of eigenvalues. May be perturbed on exit.

in
wi

Array of dimension n. Imaginary parts of eigenvalues.

inout
VL

Left eigenvectors. Array of dimension (ldvl, mm).

in
ldvl

Leading dimension of VL (ldvl >= 1; ldvl >= n if computing left eigenvectors).

inout
VR

Right eigenvectors. Array of dimension (ldvr, mm).

in
ldvr

Leading dimension of VR (ldvr >= 1; ldvr >= n if computing right eigenvectors).

in
mm

Number of columns in VL and/or VR (mm >= m).

out
m

Number of columns required to store the eigenvectors.

out
work

Workspace array of dimension (n+2)*n.

out
ifaill

Array of dimension mm. Convergence status for left eigenvectors (-1 if converged, k (0-based) if eigenvector corresponding to eigenvalue k failed).

out
ifailr

Array of dimension mm. Convergence status for right eigenvectors (-1 if converged, k (0-based) if failed).

out
info

  • = 0: successful exit

  • < 0: if info = -i, the i-th argument had an illegal value

  • > 0: info is the number of eigenvectors which failed to converge.

Functions

void dhsein(
    const char*          side,
    const char*          eigsrc,
    const char*          initv,
          INT*  restrict select,
    const INT            n,
    const f64*  restrict H,
    const INT            ldh,
          f64*  restrict wr,
    const f64*  restrict wi,
          f64*  restrict VL,
    const INT            ldvl,
          f64*  restrict VR,
    const INT            ldvr,
    const INT            mm,
          INT*           m,
          f64*  restrict work,
          INT*  restrict ifaill,
          INT*  restrict ifailr,
          INT*           info
);
void dhsein(const char *side, const char *eigsrc, const char *initv, INT *restrict select, const INT n, const f64 *restrict H, const INT ldh, f64 *restrict wr, const f64 *restrict wi, f64 *restrict VL, const INT ldvl, f64 *restrict VR, const INT ldvr, const INT mm, INT *m, f64 *restrict work, INT *restrict ifaill, INT *restrict ifailr, INT *info)#

DHSEIN uses inverse iteration to find specified right and/or left eigenvectors of a real upper Hessenberg matrix H.

The right eigenvector x and the left eigenvector y of the matrix H corresponding to an eigenvalue w are defined by:

         H * x = w * x,     y**h * H = w * y**h
where y**h denotes the conjugate transpose of the vector y.

Each eigenvector is normalized so that the element of largest magnitude has magnitude 1; here the magnitude of a complex number (x,y) is taken to be |x|+|y|.

Parameters

in
side

‘R’: compute right eigenvectors only; ‘L’: compute left eigenvectors only; ‘B’: compute both right and left eigenvectors.

in
eigsrc

‘Q’: eigenvalues were found using DHSEQR, so matrix splitting can be exploited; ‘N’: no assumptions on correspondence between eigenvalues and diagonal blocks.

in
initv

‘N’: no initial vectors supplied; ‘U’: user-supplied initial vectors in VL and/or VR.

inout
select

Array of dimension n. Specifies which eigenvectors to compute. select[j] = 1 (nonzero) selects the eigenvalue wr[j] + i*wi[j].

in
n

The order of the matrix H (n >= 0).

in
H

Upper Hessenberg matrix H. Array of dimension (ldh, n).

in
ldh

The leading dimension of H (ldh >= max(1,n)).

inout
wr

Array of dimension n. Real parts of eigenvalues. May be perturbed on exit.

in
wi

Array of dimension n. Imaginary parts of eigenvalues.

inout
VL

Left eigenvectors. Array of dimension (ldvl, mm).

in
ldvl

Leading dimension of VL (ldvl >= 1; ldvl >= n if computing left eigenvectors).

inout
VR

Right eigenvectors. Array of dimension (ldvr, mm).

in
ldvr

Leading dimension of VR (ldvr >= 1; ldvr >= n if computing right eigenvectors).

in
mm

Number of columns in VL and/or VR (mm >= m).

out
m

Number of columns required to store the eigenvectors.

out
work

Workspace array of dimension (n+2)*n.

out
ifaill

Array of dimension mm. Convergence status for left eigenvectors (-1 if converged, k (0-based) if eigenvector corresponding to eigenvalue k failed).

out
ifailr

Array of dimension mm. Convergence status for right eigenvectors (-1 if converged, k (0-based) if failed).

out
info

  • = 0: successful exit

  • < 0: if info = -i, the i-th argument had an illegal value

  • > 0: info is the number of eigenvectors which failed to converge.

Functions

void chsein(
    const char*          side,
    const char*          eigsrc,
    const char*          initv,
    const INT*  restrict select,
    const INT            n,
    const c64*  restrict H,
    const INT            ldh,
          c64*  restrict W,
          c64*  restrict VL,
    const INT            ldvl,
          c64*  restrict VR,
    const INT            ldvr,
    const INT            mm,
          INT*           m,
          c64*  restrict work,
          f32*  restrict rwork,
          INT*  restrict ifaill,
          INT*  restrict ifailr,
          INT*           info
);
void chsein(const char *side, const char *eigsrc, const char *initv, const INT *restrict select, const INT n, const c64 *restrict H, const INT ldh, c64 *restrict W, c64 *restrict VL, const INT ldvl, c64 *restrict VR, const INT ldvr, const INT mm, INT *m, c64 *restrict work, f32 *restrict rwork, INT *restrict ifaill, INT *restrict ifailr, INT *info)#

CHSEIN uses inverse iteration to find specified right and/or left eigenvectors of a complex upper Hessenberg matrix H.

The right eigenvector x and the left eigenvector y of the matrix H corresponding to an eigenvalue w are defined by:

         H * x = w * x,     y**h * H = w * y**h
where y**h denotes the conjugate transpose of the vector y.

Each eigenvector is normalized so that the element of largest magnitude has magnitude 1; here the magnitude of a complex number (x,y) is taken to be |x|+|y|.

Parameters

in
side

‘R’: compute right eigenvectors only; ‘L’: compute left eigenvectors only; ‘B’: compute both right and left eigenvectors.

in
eigsrc

‘Q’: eigenvalues were found using CHSEQR, so matrix splitting can be exploited; ‘N’: no assumptions on correspondence between eigenvalues and diagonal blocks.

in
initv

‘N’: no initial vectors supplied; ‘U’: user-supplied initial vectors in VL and/or VR.

in
select

Array of dimension n. Specifies which eigenvectors to compute. select[j] nonzero selects the eigenvector corresponding to W(j).

in
n

The order of the matrix H (n >= 0).

in
H

Upper Hessenberg matrix. Single complex array, dimension (ldh, n). If a NaN is detected in H, the routine will return with info=-6.

in
ldh

The leading dimension of H (ldh >= max(1,n)).

inout
W

Single complex array, dimension (n). On entry, the eigenvalues of H. On exit, the real parts of W may have been altered since close eigenvalues are perturbed slightly in searching for independent eigenvectors.

inout
VL

Left eigenvectors. Single complex array, dimension (ldvl, mm).

in
ldvl

Leading dimension of VL (ldvl >= 1; ldvl >= n if computing left eigenvectors).

inout
VR

Right eigenvectors. Single complex array, dimension (ldvr, mm).

in
ldvr

Leading dimension of VR (ldvr >= 1; ldvr >= n if computing right eigenvectors).

in
mm

Number of columns in VL and/or VR (mm >= m).

out
m

Number of columns required to store the eigenvectors.

out
work

Single complex workspace array, dimension (n*n).

out
rwork

Single precision array, dimension (n).

out
ifaill

Integer array, dimension (mm). Convergence status for left eigenvectors (-1 if converged, k (0-based) if eigenvector corresponding to eigenvalue k failed). Not referenced if side = ‘R’.

out
ifailr

Integer array, dimension (mm). Convergence status for right eigenvectors (-1 if converged, k (0-based) if eigenvector corresponding to eigenvalue k failed). Not referenced if side = ‘L’.

out
info

  • = 0: successful exit

  • < 0: if info = -i, the i-th argument had an illegal value

  • > 0: info is the number of eigenvectors which failed to converge.

Functions

void zhsein(
    const char*          side,
    const char*          eigsrc,
    const char*          initv,
    const INT*  restrict select,
    const INT            n,
    const c128* restrict H,
    const INT            ldh,
          c128* restrict W,
          c128* restrict VL,
    const INT            ldvl,
          c128* restrict VR,
    const INT            ldvr,
    const INT            mm,
          INT*           m,
          c128* restrict work,
          f64*  restrict rwork,
          INT*  restrict ifaill,
          INT*  restrict ifailr,
          INT*           info
);
void zhsein(const char *side, const char *eigsrc, const char *initv, const INT *restrict select, const INT n, const c128 *restrict H, const INT ldh, c128 *restrict W, c128 *restrict VL, const INT ldvl, c128 *restrict VR, const INT ldvr, const INT mm, INT *m, c128 *restrict work, f64 *restrict rwork, INT *restrict ifaill, INT *restrict ifailr, INT *info)#

ZHSEIN uses inverse iteration to find specified right and/or left eigenvectors of a complex upper Hessenberg matrix H.

The right eigenvector x and the left eigenvector y of the matrix H corresponding to an eigenvalue w are defined by:

         H * x = w * x,     y**h * H = w * y**h
where y**h denotes the conjugate transpose of the vector y.

Each eigenvector is normalized so that the element of largest magnitude has magnitude 1; here the magnitude of a complex number (x,y) is taken to be |x|+|y|.

Parameters

in
side

‘R’: compute right eigenvectors only; ‘L’: compute left eigenvectors only; ‘B’: compute both right and left eigenvectors.

in
eigsrc

‘Q’: eigenvalues were found using ZHSEQR, so matrix splitting can be exploited; ‘N’: no assumptions on correspondence between eigenvalues and diagonal blocks.

in
initv

‘N’: no initial vectors supplied; ‘U’: user-supplied initial vectors in VL and/or VR.

in
select

Array of dimension n. Specifies which eigenvectors to compute. select[j] nonzero selects the eigenvector corresponding to W(j).

in
n

The order of the matrix H (n >= 0).

in
H

Upper Hessenberg matrix. Double complex array, dimension (ldh, n). If a NaN is detected in H, the routine will return with info=-6.

in
ldh

The leading dimension of H (ldh >= max(1,n)).

inout
W

Double complex array, dimension (n). On entry, the eigenvalues of H. On exit, the real parts of W may have been altered since close eigenvalues are perturbed slightly in searching for independent eigenvectors.

inout
VL

Left eigenvectors. Double complex array, dimension (ldvl, mm).

in
ldvl

Leading dimension of VL (ldvl >= 1; ldvl >= n if computing left eigenvectors).

inout
VR

Right eigenvectors. Double complex array, dimension (ldvr, mm).

in
ldvr

Leading dimension of VR (ldvr >= 1; ldvr >= n if computing right eigenvectors).

in
mm

Number of columns in VL and/or VR (mm >= m).

out
m

Number of columns required to store the eigenvectors.

out
work

Double complex workspace array, dimension (n*n).

out
rwork

Double precision array, dimension (n).

out
ifaill

Integer array, dimension (mm). Convergence status for left eigenvectors (-1 if converged, k (0-based) if eigenvector corresponding to eigenvalue k failed). Not referenced if side = ‘R’.

out
ifailr

Integer array, dimension (mm). Convergence status for right eigenvectors (-1 if converged, k (0-based) if eigenvector corresponding to eigenvalue k failed). Not referenced if side = ‘L’.

out
info

  • = 0: successful exit

  • < 0: if info = -i, the i-th argument had an illegal value

  • > 0: info is the number of eigenvectors which failed to converge.