Go to main content
Oracle Developer Studio 12.5 Man Pages

Exit Print View

Updated: June 2017
 
 

zstemr (3p)

Name

zstemr - compute selected eigenvalues and, optionally, eigenvectors of a real symmetric tridiagonal matrix T

Synopsis

SUBROUTINE ZSTEMR(JOBZ, RANGE, N, D, E, VL, VU, IL, IU, M, W,
Z, LDZ, NZC, ISUPPZ, TRYRAC, WORK, LWORK, IWORK, LIWORK, INFO)

CHARACTER*1 JOBZ, RANGE
LOGICAL TRYRAC
INTEGER IL, INFO, IU, LDZ, NZC, LIWORK, LWORK, M, N
DOUBLE PRECISION VL, VU
INTEGER ISUPPZ(*), IWORK(*)
DOUBLE PRECISION D(*), E(*), W(*), WORK(*)
COMPLEX*16 Z(LDZ, *)

SUBROUTINE ZSTEMR_64(JOBZ, RANGE, N, D, E, VL, VU, IL, IU, M, W,
Z, LDZ, NZC, ISUPPZ, TRYRAC, WORK, LWORK, IWORK, LIWORK, INFO)

CHARACTER*1 JOBZ, RANGE
LOGICAL TRYRAC
INTEGER*8 IL, INFO, IU, LDZ, NZC, LIWORK, LWORK, M, N
INTEGER*8 ISUPPZ(*), IWORK(*)
DOUBLE PRECISION VL, VU
DOUBLE PRECISION D(*), E(*), W(*), WORK(*)
COMPLEX*16 Z(LDZ, *)




F95 INTERFACE
SUBROUTINE STEMR(JOBZ, RANGE, N, D, E, VL, VU, IL, IU, M, W,
Z, LDZ, NZC, ISUPPZ, TRYRAC, WORK, LWORK, IWORK, LIWORK, INFO)

CHARACTER(LEN=1) :: JOBZ, RANGE
LOGICAL TRYRAC
INTEGER :: IL, INFO, IU, LDZ, NZC, LIWORK, LWORK, M, N
INTEGER, DIMENSION(:) :: ISUPPZ, IWORK
DOUBLE PRECISION :: VL, VU
DOUBLE PRECISION, DIMENSION(:) :: D, E, W, WORK
COMPLEX*16, DIMENSION(:,:) :: Z

SUBROUTINE STEMR_64(JOBZ, RANGE, N, D, E, VL, VU, IL, IU, M, W,
Z, LDZ, NZC, ISUPPZ, TRYRAC, WORK, LWORK, IWORK, LIWORK, INFO)

CHARACTER(LEN=1) :: JOBZ, RANGE
LOGICAL TRYRAC
INTEGER(8) :: IL, INFO, IU, LDZ, NZC, LIWORK, LWORK, M, N
INTEGER(8), DIMENSION(:) :: ISUPPZ, IWORK
DOUBLE PRECISION :: VL, VU
DOUBLE PRECISION, DIMENSION(:) :: D, E, W, WORK
COMPLEX*16, DIMENSION(:,:) :: Z




C INTERFACE
#include <sunperf.h>

void zstemr(char jobz, char range, int n, double *d, double *e,  double
vl,  double vu, int il, int iu, int *m, double *w, doublecom-
plex *z, int ldz, int nzc,  int  *isuppz,  int  *tryrac,  int
*info);

void  zstemr_64(char  jobz,  char  range, long n, double *d, double *e,
double vl, double vu, long il, long iu, long *m,  double  *w,
doublecomplex  *z,  long  ldz,  long  nzc, long *isuppz, long
*tryrac, long *info);

Description

Oracle Solaris Studio Performance Library                           zstemr(3P)



NAME
       zstemr -  compute selected eigenvalues and, optionally, eigenvectors of
       a real symmetric tridiagonal matrix T


SYNOPSIS
       SUBROUTINE ZSTEMR(JOBZ, RANGE, N, D, E, VL, VU, IL, IU, M, W,
             Z, LDZ, NZC, ISUPPZ, TRYRAC, WORK, LWORK, IWORK, LIWORK, INFO)

       CHARACTER*1 JOBZ, RANGE
       LOGICAL TRYRAC
       INTEGER IL, INFO, IU, LDZ, NZC, LIWORK, LWORK, M, N
       DOUBLE PRECISION VL, VU
       INTEGER ISUPPZ(*), IWORK(*)
       DOUBLE PRECISION D(*), E(*), W(*), WORK(*)
       COMPLEX*16 Z(LDZ, *)

       SUBROUTINE ZSTEMR_64(JOBZ, RANGE, N, D, E, VL, VU, IL, IU, M, W,
             Z, LDZ, NZC, ISUPPZ, TRYRAC, WORK, LWORK, IWORK, LIWORK, INFO)

       CHARACTER*1 JOBZ, RANGE
       LOGICAL TRYRAC
       INTEGER*8 IL, INFO, IU, LDZ, NZC, LIWORK, LWORK, M, N
       INTEGER*8 ISUPPZ(*), IWORK(*)
       DOUBLE PRECISION VL, VU
       DOUBLE PRECISION D(*), E(*), W(*), WORK(*)
       COMPLEX*16 Z(LDZ, *)




   F95 INTERFACE
       SUBROUTINE STEMR(JOBZ, RANGE, N, D, E, VL, VU, IL, IU, M, W,
              Z, LDZ, NZC, ISUPPZ, TRYRAC, WORK, LWORK, IWORK, LIWORK, INFO)

       CHARACTER(LEN=1) :: JOBZ, RANGE
       LOGICAL TRYRAC
       INTEGER :: IL, INFO, IU, LDZ, NZC, LIWORK, LWORK, M, N
       INTEGER, DIMENSION(:) :: ISUPPZ, IWORK
       DOUBLE PRECISION :: VL, VU
       DOUBLE PRECISION, DIMENSION(:) :: D, E, W, WORK
       COMPLEX*16, DIMENSION(:,:) :: Z

       SUBROUTINE STEMR_64(JOBZ, RANGE, N, D, E, VL, VU, IL, IU, M, W,
              Z, LDZ, NZC, ISUPPZ, TRYRAC, WORK, LWORK, IWORK, LIWORK, INFO)

       CHARACTER(LEN=1) :: JOBZ, RANGE
       LOGICAL TRYRAC
       INTEGER(8) :: IL, INFO, IU, LDZ, NZC, LIWORK, LWORK, M, N
       INTEGER(8), DIMENSION(:) :: ISUPPZ, IWORK
       DOUBLE PRECISION :: VL, VU
       DOUBLE PRECISION, DIMENSION(:) :: D, E, W, WORK
       COMPLEX*16, DIMENSION(:,:) :: Z




   C INTERFACE
       #include <sunperf.h>

       void zstemr(char jobz, char range, int n, double *d, double *e,  double
                 vl,  double vu, int il, int iu, int *m, double *w, doublecom-
                 plex *z, int ldz, int nzc,  int  *isuppz,  int  *tryrac,  int
                 *info);

       void  zstemr_64(char  jobz,  char  range, long n, double *d, double *e,
                 double vl, double vu, long il, long iu, long *m,  double  *w,
                 doublecomplex  *z,  long  ldz,  long  nzc, long *isuppz, long
                 *tryrac, long *info);



PURPOSE
       ZSTEMR computes selected eigenvalues and, optionally, eigenvectors of a
       real  symmetric  tridiagonal  matrix T. Any such unreduced matrix has a
       well defined set of pairwise different  real  eigenvalues,  the  corre-
       sponding real eigenvectors are pairwise orthogonal.

       The spectrum may be computed either completely or partially by specify-
       ing either an interval (VL,VU] or a range  of  indices  IL:IU  for  the
       desired eigenvalues.

       Depending  on  the  number  of  desired eigenvalues, these are computed
       either by bisection  or  the  dqds  algorithm.  Numerically  orthogonal
       eigenvectors  are  computed by the use of various suitable L D L^T fac-
       torizations near clusters of close eigenvalues (referred  to  as  RRRs,
       Relatively Robust Representations). An informal sketch of the algorithm
       follows.

       For each unreduced block (submatrix) of T,
          (a) Compute T - sigma I  = L D L^T, so that L and D
              define all the wanted eigenvalues to high relative accuracy.
              This means that small relative changes in the entries of D and L
              cause only small relative changes in the eigenvalues and
              eigenvectors. The standard (unfactored) representation of the
              tridiagonal matrix T does not have this property in general.
          (b) Compute the eigenvalues to suitable accuracy.
              If the eigenvectors are desired, the algorithm attains full
              accuracy of the computed eigenvalues only right before
              the  corresponding vectors have to be computed, see steps c) and
       d).
          (c) For each cluster of close eigenvalues, select a new
              shift close to the cluster, find a new factorization, and refine
              the shifted eigenvalues to suitable accuracy.
          (d) For each eigenvalue with a large enough relative separation com-
       pute
              the  corresponding  eigenvector  by  forming  a  rank  revealing
       twisted
              factorization. Go back to (c) for any clusters that remain.

       For more details, see:
       -  Inderjit  S. Dhillon and Beresford N. Parlett: "Multiple representa-
       tions
         to compute orthogonal eigenvectors of  symmetric  tridiagonal  matri-
       ces,"
         Linear  Algebra  and its Applications, 387(1), pp. 1-28, August 2004.
       - Inderjit Dhillon and Beresford Parlett: "Orthogonal Eigenvectors and
         Relative Gaps," SIAM Journal on  Matrix  Analysis  and  Applications,
       Vol. 25,
         2004.  Also LAPACK Working Note 154.
       - Inderjit Dhillon: "A new O(n^2) algorithm for the symmetric
         tridiagonal eigenvalue/eigenvector problem",
         Computer Science Division Technical Report No. UCB/CSD-97-971,
         UC Berkeley, May 1997.

       Notes:
       1.   ZSTEMR works only on machines which follow IEEE-754 floating-point
       standard in their handling of infinities and NaNs.   This  permits  the
       use of efficient inner loops avoiding a check for zero divisors.
       2.  LAPACK routines can be used to reduce a complex Hermitian matrix to
       real symmetric tridiagonal form.  (Any  complex  Hermitian  tridiagonal
       matrix  has real values on its diagonal and potentially complex numbers
       on its off-diagonals. By applying a similarity transform with an appro-
       priate diagonal matrix, the complex Hermitian matrix can be transformed
       into a real symmetric matrix and complex  arithmetic  can  be  entirely
       avoided.)   While  the  eigenvectors  of the real symmetric tridiagonal
       matrix are real, the eigenvectors of original complex Hermitean  matrix
       have  complex  entries  in general.  Since LAPACK drivers overwrite the
       matrix data with the eigenvectors, ZSTEMR accepts complex workspace  to
       facilitate interoperability with CUNMTR or CUPMTR.


ARGUMENTS
       JOBZ (input) CHARACTER*1
                 = 'N':  Compute eigenvalues only;
                 = 'V':  Compute eigenvalues and eigenvectors.


       RANGE (input) CHARACTER*1
                 = 'A': all eigenvalues will be found.
                 = 'V': all eigenvalues in the half-open interval (VL,VU] will
                 be found.  = 'I': the IL-th through IU-th eigenvalues will be
                 found.


       N (input) INTEGER
                 The order of the matrix.  N >= 0.


       D (input/output) DOUBLE PRECISION array, dimension (N)
                 On  entry,  the n diagonal elements of the tridiagonal matrix
                 T. On exit, D is overwritten.


       E (input/output) DOUBLE PRECISION array, dimension (N)
                 On entry, the (n-1) subdiagonal elements of  the  tridiagonal
                 matrix T in elements 1 to N-1 of E; E(N) need not be set.  On
                 exit, E is overwritten.


       VL (input) INTEGER
                 If RANGE='V', the lower and upper bounds of the  interval  to
                 be  searched  for  eigenvalues.  VL  < VU.  Not referenced if
                 RANGE = 'A' or 'I'.


       VU (input) INTEGER
                 See the description of VL.


       IL (input) INTEGER
                 If RANGE='I', the indices (in ascending order) of the  small-
                 est and largest eigenvalues to be returned.  1 <= IL <= IU <=
                 N, if N > 0; IL = 1 and IU = 0 if N = 0.  Not  referenced  if
                 RANGE = 'A' or 'V'.


       IU (input) INTEGER
                 See the description of IL.


       M (output) INTEGER
                 The  total  number  of  eigenvalues  found.  0 <= M <= N.  If
                 RANGE = 'A', M = N, and if RANGE = 'I', M = IU-IL+1.


       W (output) DOUBLE PRECISION array, dimension (N)
                 The first M elements  contain  the  selected  eigenvalues  in
                 ascending order.


       Z (output) COMPLEX*16 array, dimension (LDZ, max(1,M) )
                 If  JOBZ  =  'V',  then if INFO = 0, the first M columns of Z
                 contain the orthonormal eigenvectors of the matrix  T  corre-
                 sponding to the selected eigenvalues, with the i-th column of
                 Z holding the eigenvector associated with W(i).   If  JOBZ  =
                 'N',  then  Z  is not referenced.  Note: the user must ensure
                 that at least max(1,M) columns are supplied in the  array  Z;
                 if  RANGE = 'V', the exact value of M is not known in advance
                 and can be computed with a workspace query by setting  NZC  =
                 -1, see below.


       LDZ (input) INTEGER
                 The  leading dimension of the array Z.  LDZ >= 1, and if JOBZ
                 = 'V', LDZ >= max(1,N).

                 NZC     (input) INTEGER The number of eigenvectors to be held
                 in the array Z.
                 If RANGE = 'A', then NZC >= max(1,N).
                 If  RANGE  =  'V',  then  NZC >= the number of eigenvalues in
                 (VL,VU].
                 If RANGE = 'I', then NZC >= IU-IL+1.
                 If NZC = -1, then a workspace query is assumed;  the  routine
                 calculates  the  number  of  columns  of the array Z that are
                 needed to hold the eigenvectors.  This value is  returned  as
                 the  first entry of the Z array, and no error message related
                 to NZC is issued by XERBLA.


       ISUPPZ (output) INTEGER ARRAY, dimension ( 2*max(1,M) )
                 The support of the eigenvectors in Z, i.e., the indices indi-
                 cating  the  nonzero  elements  in Z. The i-th eigenvector is
                 nonzero only in elements ISUPPZ( 2*i-1 ) through ISUPPZ(  2*i
                 ).   This  is  relevant in the case when the matrix is split.
                 ISUPPZ is only accessed when JOBZ is 'V' and N > 0.


       TRYRAC  (input/output) LOGICAL
                 If TRYRAC.EQ..TRUE., indicates that  the  code  should  check
                 whether  the  tridiagonal  matrix  defines its eigenvalues to
                 high relative accuracy.  If so, the code uses  relative-accu-
                 racy  preserving  algorithms  that  might  be  (a bit) slower
                 depending on the matrix.  If the matrix does not  define  its
                 eigenvalues to high relative accuracy, the code can uses pos-
                 sibly faster algorithms.  If TRYRAC.EQ..FALSE., the  code  is
                 not required to guarantee relatively accurate eigenvalues and
                 can use the fastest possible techniques.  On exit,  a  .TRUE.
                 TRYRAC  will  be set to .FALSE. if the matrix does not define
                 its eigenvalues to high relative accuracy.


       WORK (workspace/output) DOUBLE PRECISION array, dimension (LWORK)
                 On exit, if INFO = 0, WORK(1) returns the optimal (and  mini-
                 mal) LWORK.


       LWORK (input)
                 The  dimension  of  the  array WORK.  LWORK >= max(1,18*N) if
                 JOBZ = 'V', and LWORK >= max(1,12*N) if JOBZ = 'N'.  If LWORK
                 =  -1,  then  a  workspace query is assumed; the routine only
                 calculates the optimal size of the WORK array,  returns  this
                 value as the first entry of the WORK array, and no error mes-
                 sage related to LWORK is issued by XERBLA.


       IWORK (workspace/output) INTEGER array, dimension (LIWORK)
                 On exit, if INFO = 0, IWORK(1) returns the optimal LIWORK.


       LIWORK (input) INTEGER
                 The dimension of the array IWORK.  LIWORK >=  max(1,10*N)  if
                 the  eigenvectors  are  desired,  and LIWORK >= max(1,8*N) if
                 only the eigenvalues are to be computed.   If  LIWORK  =  -1,
                 then  a  workspace  query is assumed; the routine only calcu-
                 lates the optimal size of the IWORK array, returns this value
                 as  the  first entry of the IWORK array, and no error message
                 related to LIWORK is issued by XERBLA.


       INFO (output)
                 = 0:  successful exit
                 < 0:  if INFO = -i, the i-th argument had an illegal value
                 > 0:  if INFO = 1, internal error in SLARRE, if  INFO  =  2X,
                 internal error in CLARRV.  Here, the digit X = ABS( IINFO ) <
                 10, where IINFO is the nonzero error code returned by  SLARRE
                 or CLARRV, respectively.


FURTHER DETAILS
       Based on contributions by
          Beresford Parlett, University of California, Berkeley, USA
          Jim Demmel, University of California, Berkeley, USA
          Inderjit Dhillon, University of Texas, Austin, USA
          Osni Marques, LBNL/NERSC, USA
          Christof Voemel, University of California, Berkeley, USA




                                  7 Nov 2015                        zstemr(3P)