*> \brief SSBEVD_2STAGE computes the eigenvalues and, optionally, the left and/or right eigenvectors for OTHER matrices
*
* @generated from dsbevd_2stage.f, fortran d -> s, Sat Nov 5 23:58:03 2016
*
* =========== DOCUMENTATION ===========
*
* Online html documentation available at
* http://www.netlib.org/lapack/explore-html/
*
*> \htmlonly
*> Download SSBEVD_2STAGE + dependencies
*>
*> [TGZ]
*>
*> [ZIP]
*>
*> [TXT]
*> \endhtmlonly
*
* Definition:
* ===========
*
* SUBROUTINE SSBEVD_2STAGE( JOBZ, UPLO, N, KD, AB, LDAB, W, Z, LDZ,
* WORK, LWORK, IWORK, LIWORK, INFO )
*
* IMPLICIT NONE
*
* .. Scalar Arguments ..
* CHARACTER JOBZ, UPLO
* INTEGER INFO, KD, LDAB, LDZ, LIWORK, LWORK, N
* ..
* .. Array Arguments ..
* INTEGER IWORK( * )
* REAL AB( LDAB, * ), W( * ), WORK( * ), Z( LDZ, * )
* ..
*
*
*> \par Purpose:
* =============
*>
*> \verbatim
*>
*> SSBEVD_2STAGE computes all the eigenvalues and, optionally, eigenvectors of
*> a real symmetric band matrix A using the 2stage technique for
*> the reduction to tridiagonal. If eigenvectors are desired, it uses
*> a divide and conquer algorithm.
*>
*> The divide and conquer algorithm makes very mild assumptions about
*> floating point arithmetic. It will work on machines with a guard
*> digit in add/subtract, or on those binary machines without guard
*> digits which subtract like the Cray X-MP, Cray Y-MP, Cray C-90, or
*> Cray-2. It could conceivably fail on hexadecimal or decimal machines
*> without guard digits, but we know of none.
*> \endverbatim
*
* Arguments:
* ==========
*
*> \param[in] JOBZ
*> \verbatim
*> JOBZ is CHARACTER*1
*> = 'N': Compute eigenvalues only;
*> = 'V': Compute eigenvalues and eigenvectors.
*> Not available in this release.
*> \endverbatim
*>
*> \param[in] UPLO
*> \verbatim
*> UPLO is CHARACTER*1
*> = 'U': Upper triangle of A is stored;
*> = 'L': Lower triangle of A is stored.
*> \endverbatim
*>
*> \param[in] N
*> \verbatim
*> N is INTEGER
*> The order of the matrix A. N >= 0.
*> \endverbatim
*>
*> \param[in] KD
*> \verbatim
*> KD is INTEGER
*> The number of superdiagonals of the matrix A if UPLO = 'U',
*> or the number of subdiagonals if UPLO = 'L'. KD >= 0.
*> \endverbatim
*>
*> \param[in,out] AB
*> \verbatim
*> AB is REAL array, dimension (LDAB, N)
*> On entry, the upper or lower triangle of the symmetric band
*> matrix A, stored in the first KD+1 rows of the array. The
*> j-th column of A is stored in the j-th column of the array AB
*> as follows:
*> if UPLO = 'U', AB(kd+1+i-j,j) = A(i,j) for max(1,j-kd)<=i<=j;
*> if UPLO = 'L', AB(1+i-j,j) = A(i,j) for j<=i<=min(n,j+kd).
*>
*> On exit, AB is overwritten by values generated during the
*> reduction to tridiagonal form. If UPLO = 'U', the first
*> superdiagonal and the diagonal of the tridiagonal matrix T
*> are returned in rows KD and KD+1 of AB, and if UPLO = 'L',
*> the diagonal and first subdiagonal of T are returned in the
*> first two rows of AB.
*> \endverbatim
*>
*> \param[in] LDAB
*> \verbatim
*> LDAB is INTEGER
*> The leading dimension of the array AB. LDAB >= KD + 1.
*> \endverbatim
*>
*> \param[out] W
*> \verbatim
*> W is REAL array, dimension (N)
*> If INFO = 0, the eigenvalues in ascending order.
*> \endverbatim
*>
*> \param[out] Z
*> \verbatim
*> Z is REAL array, dimension (LDZ, N)
*> If JOBZ = 'V', then if INFO = 0, Z contains the orthonormal
*> eigenvectors of the matrix A, with the i-th column of Z
*> holding the eigenvector associated with W(i).
*> If JOBZ = 'N', then Z is not referenced.
*> \endverbatim
*>
*> \param[in] LDZ
*> \verbatim
*> LDZ is INTEGER
*> The leading dimension of the array Z. LDZ >= 1, and if
*> JOBZ = 'V', LDZ >= max(1,N).
*> \endverbatim
*>
*> \param[out] WORK
*> \verbatim
*> WORK is REAL array, dimension LWORK
*> On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
*> \endverbatim
*>
*> \param[in] LWORK
*> \verbatim
*> LWORK is INTEGER
*> The length of the array WORK. LWORK >= 1, when N <= 1;
*> otherwise
*> If JOBZ = 'N' and N > 1, LWORK must be queried.
*> LWORK = MAX(1, dimension) where
*> dimension = (2KD+1)*N + KD*NTHREADS + N
*> where KD is the size of the band.
*> NTHREADS is the number of threads used when
*> openMP compilation is enabled, otherwise =1.
*> If JOBZ = 'V' and N > 1, LWORK must be queried. Not yet available.
*>
*> If LWORK = -1, then a workspace query is assumed; the routine
*> only calculates the optimal sizes of the WORK and IWORK
*> arrays, returns these values as the first entries of the WORK
*> and IWORK arrays, and no error message related to LWORK or
*> LIWORK is issued by XERBLA.
*> \endverbatim
*>
*> \param[out] IWORK
*> \verbatim
*> IWORK is INTEGER array, dimension (MAX(1,LIWORK))
*> On exit, if INFO = 0, IWORK(1) returns the optimal LIWORK.
*> \endverbatim
*>
*> \param[in] LIWORK
*> \verbatim
*> LIWORK is INTEGER
*> The dimension of the array IWORK.
*> If JOBZ = 'N' or N <= 1, LIWORK must be at least 1.
*> If JOBZ = 'V' and N > 2, LIWORK must be at least 3 + 5*N.
*>
*> If LIWORK = -1, then a workspace query is assumed; the
*> routine only calculates the optimal sizes of the WORK and
*> IWORK arrays, returns these values as the first entries of
*> the WORK and IWORK arrays, and no error message related to
*> LWORK or LIWORK is issued by XERBLA.
*> \endverbatim
*>
*> \param[out] INFO
*> \verbatim
*> INFO is INTEGER
*> = 0: successful exit
*> < 0: if INFO = -i, the i-th argument had an illegal value
*> > 0: if INFO = i, the algorithm failed to converge; i
*> off-diagonal elements of an intermediate tridiagonal
*> form did not converge to zero.
*> \endverbatim
*
* Authors:
* ========
*
*> \author Univ. of Tennessee
*> \author Univ. of California Berkeley
*> \author Univ. of Colorado Denver
*> \author NAG Ltd.
*
*> \ingroup realOTHEReigen
*
*> \par Further Details:
* =====================
*>
*> \verbatim
*>
*> All details about the 2stage techniques are available in:
*>
*> Azzam Haidar, Hatem Ltaief, and Jack Dongarra.
*> Parallel reduction to condensed forms for symmetric eigenvalue problems
*> using aggregated fine-grained and memory-aware kernels. In Proceedings
*> of 2011 International Conference for High Performance Computing,
*> Networking, Storage and Analysis (SC '11), New York, NY, USA,
*> Article 8 , 11 pages.
*> http://doi.acm.org/10.1145/2063384.2063394
*>
*> A. Haidar, J. Kurzak, P. Luszczek, 2013.
*> An improved parallel singular value algorithm and its implementation
*> for multicore hardware, In Proceedings of 2013 International Conference
*> for High Performance Computing, Networking, Storage and Analysis (SC '13).
*> Denver, Colorado, USA, 2013.
*> Article 90, 12 pages.
*> http://doi.acm.org/10.1145/2503210.2503292
*>
*> A. Haidar, R. Solca, S. Tomov, T. Schulthess and J. Dongarra.
*> A novel hybrid CPU-GPU generalized eigensolver for electronic structure
*> calculations based on fine-grained memory aware tasks.
*> International Journal of High Performance Computing Applications.
*> Volume 28 Issue 2, Pages 196-209, May 2014.
*> http://hpc.sagepub.com/content/28/2/196
*>
*> \endverbatim
*
* =====================================================================
SUBROUTINE SSBEVD_2STAGE( JOBZ, UPLO, N, KD, AB, LDAB, W, Z, LDZ,
$ WORK, LWORK, IWORK, LIWORK, INFO )
*
IMPLICIT NONE
*
* -- LAPACK driver routine --
* -- LAPACK is a software package provided by Univ. of Tennessee, --
* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
*
* .. Scalar Arguments ..
CHARACTER JOBZ, UPLO
INTEGER INFO, KD, LDAB, LDZ, LIWORK, LWORK, N
* ..
* .. Array Arguments ..
INTEGER IWORK( * )
REAL AB( LDAB, * ), W( * ), WORK( * ), Z( LDZ, * )
* ..
*
* =====================================================================
*
* .. Parameters ..
REAL ZERO, ONE
PARAMETER ( ZERO = 0.0E+0, ONE = 1.0E+0 )
* ..
* .. Local Scalars ..
LOGICAL LOWER, LQUERY, WANTZ
INTEGER IINFO, INDE, INDWK2, INDWRK, ISCALE, LIWMIN,
$ LLWORK, LWMIN, LHTRD, LWTRD, IB, INDHOUS,
$ LLWRK2
REAL ANRM, BIGNUM, EPS, RMAX, RMIN, SAFMIN, SIGMA,
$ SMLNUM
* ..
* .. External Functions ..
LOGICAL LSAME
INTEGER ILAENV2STAGE
REAL SLAMCH, SLANSB
EXTERNAL LSAME, SLAMCH, SLANSB, ILAENV2STAGE
* ..
* .. External Subroutines ..
EXTERNAL SGEMM, SLACPY, SLASCL, SSCAL, SSTEDC,
$ SSTERF, XERBLA, SSYTRD_SB2ST
* ..
* .. Intrinsic Functions ..
INTRINSIC SQRT
* ..
* .. Executable Statements ..
*
* Test the input parameters.
*
WANTZ = LSAME( JOBZ, 'V' )
LOWER = LSAME( UPLO, 'L' )
LQUERY = ( LWORK.EQ.-1 .OR. LIWORK.EQ.-1 )
*
INFO = 0
IF( N.LE.1 ) THEN
LIWMIN = 1
LWMIN = 1
ELSE
IB = ILAENV2STAGE( 2, 'SSYTRD_SB2ST', JOBZ, N, KD, -1, -1 )
LHTRD = ILAENV2STAGE( 3, 'SSYTRD_SB2ST', JOBZ, N, KD, IB, -1 )
LWTRD = ILAENV2STAGE( 4, 'SSYTRD_SB2ST', JOBZ, N, KD, IB, -1 )
IF( WANTZ ) THEN
LIWMIN = 3 + 5*N
LWMIN = 1 + 5*N + 2*N**2
ELSE
LIWMIN = 1
LWMIN = MAX( 2*N, N+LHTRD+LWTRD )
END IF
END IF
IF( .NOT.( LSAME( JOBZ, 'N' ) ) ) THEN
INFO = -1
ELSE IF( .NOT.( LOWER .OR. LSAME( UPLO, 'U' ) ) ) THEN
INFO = -2
ELSE IF( N.LT.0 ) THEN
INFO = -3
ELSE IF( KD.LT.0 ) THEN
INFO = -4
ELSE IF( LDAB.LT.KD+1 ) THEN
INFO = -6
ELSE IF( LDZ.LT.1 .OR. ( WANTZ .AND. LDZ.LT.N ) ) THEN
INFO = -9
END IF
*
IF( INFO.EQ.0 ) THEN
WORK( 1 ) = LWMIN
IWORK( 1 ) = LIWMIN
*
IF( LWORK.LT.LWMIN .AND. .NOT.LQUERY ) THEN
INFO = -11
ELSE IF( LIWORK.LT.LIWMIN .AND. .NOT.LQUERY ) THEN
INFO = -13
END IF
END IF
*
IF( INFO.NE.0 ) THEN
CALL XERBLA( 'SSBEVD_2STAGE', -INFO )
RETURN
ELSE IF( LQUERY ) THEN
RETURN
END IF
*
* Quick return if possible
*
IF( N.EQ.0 )
$ RETURN
*
IF( N.EQ.1 ) THEN
W( 1 ) = AB( 1, 1 )
IF( WANTZ )
$ Z( 1, 1 ) = ONE
RETURN
END IF
*
* Get machine constants.
*
SAFMIN = SLAMCH( 'Safe minimum' )
EPS = SLAMCH( 'Precision' )
SMLNUM = SAFMIN / EPS
BIGNUM = ONE / SMLNUM
RMIN = SQRT( SMLNUM )
RMAX = SQRT( BIGNUM )
*
* Scale matrix to allowable range, if necessary.
*
ANRM = SLANSB( 'M', UPLO, N, KD, AB, LDAB, WORK )
ISCALE = 0
IF( ANRM.GT.ZERO .AND. ANRM.LT.RMIN ) THEN
ISCALE = 1
SIGMA = RMIN / ANRM
ELSE IF( ANRM.GT.RMAX ) THEN
ISCALE = 1
SIGMA = RMAX / ANRM
END IF
IF( ISCALE.EQ.1 ) THEN
IF( LOWER ) THEN
CALL SLASCL( 'B', KD, KD, ONE, SIGMA, N, N, AB, LDAB, INFO )
ELSE
CALL SLASCL( 'Q', KD, KD, ONE, SIGMA, N, N, AB, LDAB, INFO )
END IF
END IF
*
* Call SSYTRD_SB2ST to reduce band symmetric matrix to tridiagonal form.
*
INDE = 1
INDHOUS = INDE + N
INDWRK = INDHOUS + LHTRD
LLWORK = LWORK - INDWRK + 1
INDWK2 = INDWRK + N*N
LLWRK2 = LWORK - INDWK2 + 1
*
CALL SSYTRD_SB2ST( "N", JOBZ, UPLO, N, KD, AB, LDAB, W,
$ WORK( INDE ), WORK( INDHOUS ), LHTRD,
$ WORK( INDWRK ), LLWORK, IINFO )
*
* For eigenvalues only, call SSTERF. For eigenvectors, call SSTEDC.
*
IF( .NOT.WANTZ ) THEN
CALL SSTERF( N, W, WORK( INDE ), INFO )
ELSE
CALL SSTEDC( 'I', N, W, WORK( INDE ), WORK( INDWRK ), N,
$ WORK( INDWK2 ), LLWRK2, IWORK, LIWORK, INFO )
CALL SGEMM( 'N', 'N', N, N, N, ONE, Z, LDZ, WORK( INDWRK ), N,
$ ZERO, WORK( INDWK2 ), N )
CALL SLACPY( 'A', N, N, WORK( INDWK2 ), N, Z, LDZ )
END IF
*
* If matrix was scaled, then rescale eigenvalues appropriately.
*
IF( ISCALE.EQ.1 )
$ CALL SSCAL( N, ONE / SIGMA, W, 1 )
*
WORK( 1 ) = LWMIN
IWORK( 1 ) = LIWMIN
RETURN
*
* End of SSBEVD_2STAGE
*
END