OpenRadioss 2025.1.11
OpenRadioss project
Loading...
Searching...
No Matches
chbevx_2stage.f
Go to the documentation of this file.
1*> \brief <b> CHBEVX_2STAGE computes the eigenvalues and, optionally, the left and/or right eigenvectors for OTHER matrices</b>
2*
3* @generated from zhbevx_2stage.f, fortran z -> c, Sat Nov 5 23:18:22 2016
4*
5* =========== DOCUMENTATION ===========
6*
7* Online html documentation available at
8* http://www.netlib.org/lapack/explore-html/
9*
10*> \htmlonly
11*> Download CHBEVX_2STAGE + dependencies
12*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/chbevx_2stage.f">
13*> [TGZ]</a>
14*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/chbevx_2stage.f">
15*> [ZIP]</a>
16*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/chbevx_2stage.f">
17*> [TXT]</a>
18*> \endhtmlonly
19*
20* Definition:
21* ===========
22*
23* SUBROUTINE CHBEVX_2STAGE( JOBZ, RANGE, UPLO, N, KD, AB, LDAB,
24* Q, LDQ, VL, VU, IL, IU, ABSTOL, M, W,
25* Z, LDZ, WORK, LWORK, RWORK, IWORK,
26* IFAIL, INFO )
27*
28* IMPLICIT NONE
29*
30* .. Scalar Arguments ..
31* CHARACTER JOBZ, RANGE, UPLO
32* INTEGER IL, INFO, IU, KD, LDAB, LDQ, LDZ, M, N, LWORK
33* REAL ABSTOL, VL, VU
34* ..
35* .. Array Arguments ..
36* INTEGER IFAIL( * ), IWORK( * )
37* REAL RWORK( * ), W( * )
38* COMPLEX AB( LDAB, * ), Q( LDQ, * ), WORK( * ),
39* $ Z( LDZ, * )
40* ..
41*
42*
43*> \par Purpose:
44* =============
45*>
46*> \verbatim
47*>
48*> CHBEVX_2STAGE computes selected eigenvalues and, optionally, eigenvectors
49*> of a complex Hermitian band matrix A using the 2stage technique for
50*> the reduction to tridiagonal. Eigenvalues and eigenvectors
51*> can be selected by specifying either a range of values or a range of
52*> indices for the desired eigenvalues.
53*> \endverbatim
54*
55* Arguments:
56* ==========
57*
58*> \param[in] JOBZ
59*> \verbatim
60*> JOBZ is CHARACTER*1
61*> = 'N': Compute eigenvalues only;
62*> = 'V': Compute eigenvalues and eigenvectors.
63*> Not available in this release.
64*> \endverbatim
65*>
66*> \param[in] RANGE
67*> \verbatim
68*> RANGE is CHARACTER*1
69*> = 'A': all eigenvalues will be found;
70*> = 'V': all eigenvalues in the half-open interval (VL,VU]
71*> will be found;
72*> = 'I': the IL-th through IU-th eigenvalues will be found.
73*> \endverbatim
74*>
75*> \param[in] UPLO
76*> \verbatim
77*> UPLO is CHARACTER*1
78*> = 'U': Upper triangle of A is stored;
79*> = 'L': Lower triangle of A is stored.
80*> \endverbatim
81*>
82*> \param[in] N
83*> \verbatim
84*> N is INTEGER
85*> The order of the matrix A. N >= 0.
86*> \endverbatim
87*>
88*> \param[in] KD
89*> \verbatim
90*> KD is INTEGER
91*> The number of superdiagonals of the matrix A if UPLO = 'U',
92*> or the number of subdiagonals if UPLO = 'L'. KD >= 0.
93*> \endverbatim
94*>
95*> \param[in,out] AB
96*> \verbatim
97*> AB is COMPLEX array, dimension (LDAB, N)
98*> On entry, the upper or lower triangle of the Hermitian band
99*> matrix A, stored in the first KD+1 rows of the array. The
100*> j-th column of A is stored in the j-th column of the array AB
101*> as follows:
102*> if UPLO = 'U', AB(kd+1+i-j,j) = A(i,j) for max(1,j-kd)<=i<=j;
103*> if UPLO = 'L', AB(1+i-j,j) = A(i,j) for j<=i<=min(n,j+kd).
104*>
105*> On exit, AB is overwritten by values generated during the
106*> reduction to tridiagonal form.
107*> \endverbatim
108*>
109*> \param[in] LDAB
110*> \verbatim
111*> LDAB is INTEGER
112*> The leading dimension of the array AB. LDAB >= KD + 1.
113*> \endverbatim
114*>
115*> \param[out] Q
116*> \verbatim
117*> Q is COMPLEX array, dimension (LDQ, N)
118*> If JOBZ = 'V', the N-by-N unitary matrix used in the
119*> reduction to tridiagonal form.
120*> If JOBZ = 'N', the array Q is not referenced.
121*> \endverbatim
122*>
123*> \param[in] LDQ
124*> \verbatim
125*> LDQ is INTEGER
126*> The leading dimension of the array Q. If JOBZ = 'V', then
127*> LDQ >= max(1,N).
128*> \endverbatim
129*>
130*> \param[in] VL
131*> \verbatim
132*> VL is REAL
133*> If RANGE='V', the lower bound of the interval to
134*> be searched for eigenvalues. VL < VU.
135*> Not referenced if RANGE = 'A' or 'I'.
136*> \endverbatim
137*>
138*> \param[in] VU
139*> \verbatim
140*> VU is REAL
141*> If RANGE='V', the upper bound of the interval to
142*> be searched for eigenvalues. VL < VU.
143*> Not referenced if RANGE = 'A' or 'I'.
144*> \endverbatim
145*>
146*> \param[in] IL
147*> \verbatim
148*> IL is INTEGER
149*> If RANGE='I', the index of the
150*> smallest eigenvalue to be returned.
151*> 1 <= IL <= IU <= N, if N > 0; IL = 1 and IU = 0 if N = 0.
152*> Not referenced if RANGE = 'A' or 'V'.
153*> \endverbatim
154*>
155*> \param[in] IU
156*> \verbatim
157*> IU is INTEGER
158*> If RANGE='I', the index of the
159*> largest eigenvalue to be returned.
160*> 1 <= IL <= IU <= N, if N > 0; IL = 1 and IU = 0 if N = 0.
161*> Not referenced if RANGE = 'A' or 'V'.
162*> \endverbatim
163*>
164*> \param[in] ABSTOL
165*> \verbatim
166*> ABSTOL is REAL
167*> The absolute error tolerance for the eigenvalues.
168*> An approximate eigenvalue is accepted as converged
169*> when it is determined to lie in an interval [a,b]
170*> of width less than or equal to
171*>
172*> ABSTOL + EPS * max( |a|,|b| ) ,
173*>
174*> where EPS is the machine precision. If ABSTOL is less than
175*> or equal to zero, then EPS*|T| will be used in its place,
176*> where |T| is the 1-norm of the tridiagonal matrix obtained
177*> by reducing AB to tridiagonal form.
178*>
179*> Eigenvalues will be computed most accurately when ABSTOL is
180*> set to twice the underflow threshold 2*SLAMCH('S'), not zero.
181*> If this routine returns with INFO>0, indicating that some
182*> eigenvectors did not converge, try setting ABSTOL to
183*> 2*SLAMCH('S').
184*>
185*> See "Computing Small Singular Values of Bidiagonal Matrices
186*> with Guaranteed High Relative Accuracy," by Demmel and
187*> Kahan, LAPACK Working Note #3.
188*> \endverbatim
189*>
190*> \param[out] M
191*> \verbatim
192*> M is INTEGER
193*> The total number of eigenvalues found. 0 <= M <= N.
194*> If RANGE = 'A', M = N, and if RANGE = 'I', M = IU-IL+1.
195*> \endverbatim
196*>
197*> \param[out] W
198*> \verbatim
199*> W is REAL array, dimension (N)
200*> The first M elements contain the selected eigenvalues in
201*> ascending order.
202*> \endverbatim
203*>
204*> \param[out] Z
205*> \verbatim
206*> Z is COMPLEX array, dimension (LDZ, max(1,M))
207*> If JOBZ = 'V', then if INFO = 0, the first M columns of Z
208*> contain the orthonormal eigenvectors of the matrix A
209*> corresponding to the selected eigenvalues, with the i-th
210*> column of Z holding the eigenvector associated with W(i).
211*> If an eigenvector fails to converge, then that column of Z
212*> contains the latest approximation to the eigenvector, and the
213*> index of the eigenvector is returned in IFAIL.
214*> If JOBZ = 'N', then Z is not referenced.
215*> Note: the user must ensure that at least max(1,M) columns are
216*> supplied in the array Z; if RANGE = 'V', the exact value of M
217*> is not known in advance and an upper bound must be used.
218*> \endverbatim
219*>
220*> \param[in] LDZ
221*> \verbatim
222*> LDZ is INTEGER
223*> The leading dimension of the array Z. LDZ >= 1, and if
224*> JOBZ = 'V', LDZ >= max(1,N).
225*> \endverbatim
226*>
227*> \param[out] WORK
228*> \verbatim
229*> WORK is COMPLEX array, dimension (LWORK)
230*> \endverbatim
231*>
232*> \param[in] LWORK
233*> \verbatim
234*> LWORK is INTEGER
235*> The length of the array WORK. LWORK >= 1, when N <= 1;
236*> otherwise
237*> If JOBZ = 'N' and N > 1, LWORK must be queried.
238*> LWORK = MAX(1, dimension) where
239*> dimension = (2KD+1)*N + KD*NTHREADS
240*> where KD is the size of the band.
241*> NTHREADS is the number of threads used when
242*> openMP compilation is enabled, otherwise =1.
243*> If JOBZ = 'V' and N > 1, LWORK must be queried. Not yet available.
244*>
245*> If LWORK = -1, then a workspace query is assumed; the routine
246*> only calculates the optimal sizes of the WORK, RWORK and
247*> IWORK arrays, returns these values as the first entries of
248*> the WORK, RWORK and IWORK arrays, and no error message
249*> related to LWORK or LRWORK or LIWORK is issued by XERBLA.
250*> \endverbatim
251*>
252*> \param[out] RWORK
253*> \verbatim
254*> RWORK is REAL array, dimension (7*N)
255*> \endverbatim
256*>
257*> \param[out] IWORK
258*> \verbatim
259*> IWORK is INTEGER array, dimension (5*N)
260*> \endverbatim
261*>
262*> \param[out] IFAIL
263*> \verbatim
264*> IFAIL is INTEGER array, dimension (N)
265*> If JOBZ = 'V', then if INFO = 0, the first M elements of
266*> IFAIL are zero. If INFO > 0, then IFAIL contains the
267*> indices of the eigenvectors that failed to converge.
268*> If JOBZ = 'N', then IFAIL is not referenced.
269*> \endverbatim
270*>
271*> \param[out] INFO
272*> \verbatim
273*> INFO is INTEGER
274*> = 0: successful exit
275*> < 0: if INFO = -i, the i-th argument had an illegal value
276*> > 0: if INFO = i, then i eigenvectors failed to converge.
277*> Their indices are stored in array IFAIL.
278*> \endverbatim
279*
280* Authors:
281* ========
282*
283*> \author Univ. of Tennessee
284*> \author Univ. of California Berkeley
285*> \author Univ. of Colorado Denver
286*> \author NAG Ltd.
287*
288*> \ingroup complexOTHEReigen
289*
290*> \par Further Details:
291* =====================
292*>
293*> \verbatim
294*>
295*> All details about the 2stage techniques are available in:
296*>
297*> Azzam Haidar, Hatem Ltaief, and Jack Dongarra.
298*> Parallel reduction to condensed forms for symmetric eigenvalue problems
299*> using aggregated fine-grained and memory-aware kernels. In Proceedings
300*> of 2011 International Conference for High Performance Computing,
301*> Networking, Storage and Analysis (SC '11), New York, NY, USA,
302*> Article 8 , 11 pages.
303*> http://doi.acm.org/10.1145/2063384.2063394
304*>
305*> A. Haidar, J. Kurzak, P. Luszczek, 2013.
306*> An improved parallel singular value algorithm and its implementation
307*> for multicore hardware, In Proceedings of 2013 International Conference
308*> for High Performance Computing, Networking, Storage and Analysis (SC '13).
309*> Denver, Colorado, USA, 2013.
310*> Article 90, 12 pages.
311*> http://doi.acm.org/10.1145/2503210.2503292
312*>
313*> A. Haidar, R. Solca, S. Tomov, T. Schulthess and J. Dongarra.
314*> A novel hybrid CPU-GPU generalized eigensolver for electronic structure
315*> calculations based on fine-grained memory aware tasks.
316*> International Journal of High Performance Computing Applications.
317*> Volume 28 Issue 2, Pages 196-209, May 2014.
318*> http://hpc.sagepub.com/content/28/2/196
319*>
320*> \endverbatim
321*
322* =====================================================================
323 SUBROUTINE chbevx_2stage( JOBZ, RANGE, UPLO, N, KD, AB, LDAB,
324 $ Q, LDQ, VL, VU, IL, IU, ABSTOL, M, W,
325 $ Z, LDZ, WORK, LWORK, RWORK, IWORK,
326 $ IFAIL, INFO )
327*
328 IMPLICIT NONE
329*
330* -- LAPACK driver routine --
331* -- LAPACK is a software package provided by Univ. of Tennessee, --
332* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
333*
334* .. Scalar Arguments ..
335 CHARACTER JOBZ, RANGE, UPLO
336 INTEGER IL, INFO, IU, KD, LDAB, LDQ, LDZ, M, N, LWORK
337 REAL ABSTOL, VL, VU
338* ..
339* .. Array Arguments ..
340 INTEGER IFAIL( * ), IWORK( * )
341 REAL RWORK( * ), W( * )
342 COMPLEX AB( LDAB, * ), Q( LDQ, * ), WORK( * ),
343 $ z( ldz, * )
344* ..
345*
346* =====================================================================
347*
348* .. Parameters ..
349 REAL ZERO, ONE
350 PARAMETER ( ZERO = 0.0e0, one = 1.0e0 )
351 COMPLEX CZERO, CONE
352 parameter( czero = ( 0.0e0, 0.0e0 ),
353 $ cone = ( 1.0e0, 0.0e0 ) )
354* ..
355* .. Local Scalars ..
356 LOGICAL ALLEIG, INDEIG, LOWER, TEST, VALEIG, WANTZ,
357 $ LQUERY
358 CHARACTER ORDER
359 INTEGER I, IINFO, IMAX, INDD, INDE, INDEE, INDIBL,
360 $ INDISP, INDIWK, INDRWK, INDWRK, ISCALE, ITMP1,
361 $ llwork, lwmin, lhtrd, lwtrd, ib, indhous,
362 $ j, jj, nsplit
363 REAL ABSTLL, ANRM, BIGNUM, EPS, RMAX, RMIN, SAFMIN,
364 $ SIGMA, SMLNUM, TMP1, VLL, VUU
365 COMPLEX CTMP1
366* ..
367* .. External Functions ..
368 LOGICAL LSAME
369 INTEGER ILAENV2STAGE
370 REAL SLAMCH, CLANHB
371 EXTERNAL lsame, slamch, clanhb, ilaenv2stage
372* ..
373* .. External Subroutines ..
374 EXTERNAL scopy, sscal, sstebz, ssterf, xerbla, ccopy,
377* ..
378* .. Intrinsic Functions ..
379 INTRINSIC real, max, min, sqrt
380* ..
381* .. Executable Statements ..
382*
383* Test the input parameters.
384*
385 wantz = lsame( jobz, 'V' )
386 alleig = lsame( range, 'A' )
387 valeig = lsame( range, 'V' )
388 indeig = lsame( range, 'I' )
389 lower = lsame( uplo, 'L' )
390 lquery = ( lwork.EQ.-1 )
391*
392 info = 0
393 IF( .NOT.( lsame( jobz, 'N' ) ) ) THEN
394 info = -1
395 ELSE IF( .NOT.( alleig .OR. valeig .OR. indeig ) ) THEN
396 info = -2
397 ELSE IF( .NOT.( lower .OR. lsame( uplo, 'U' ) ) ) THEN
398 info = -3
399 ELSE IF( n.LT.0 ) THEN
400 info = -4
401 ELSE IF( kd.LT.0 ) THEN
402 info = -5
403 ELSE IF( ldab.LT.kd+1 ) THEN
404 info = -7
405 ELSE IF( wantz .AND. ldq.LT.max( 1, n ) ) THEN
406 info = -9
407 ELSE
408 IF( valeig ) THEN
409 IF( n.GT.0 .AND. vu.LE.vl )
410 $ info = -11
411 ELSE IF( indeig ) THEN
412 IF( il.LT.1 .OR. il.GT.max( 1, n ) ) THEN
413 info = -12
414 ELSE IF( iu.LT.min( n, il ) .OR. iu.GT.n ) THEN
415 info = -13
416 END IF
417 END IF
418 END IF
419 IF( info.EQ.0 ) THEN
420 IF( ldz.LT.1 .OR. ( wantz .AND. ldz.LT.n ) )
421 $ info = -18
422 END IF
423*
424 IF( info.EQ.0 ) THEN
425 IF( n.LE.1 ) THEN
426 lwmin = 1
427 work( 1 ) = lwmin
428 ELSE
429 ib = ilaenv2stage( 2, 'CHETRD_HB2ST', jobz,
430 $ n, kd, -1, -1 )
431 lhtrd = ilaenv2stage( 3, 'CHETRD_HB2ST', jobz,
432 $ n, kd, ib, -1 )
433 lwtrd = ilaenv2stage( 4, 'chetrd_hb2st', JOBZ,
434 $ N, KD, IB, -1 )
435 LWMIN = LHTRD + LWTRD
436 WORK( 1 ) = LWMIN
437 ENDIF
438*
439.LT..AND..NOT. IF( LWORKLWMIN LQUERY )
440 $ INFO = -20
441 END IF
442*
443.NE. IF( INFO0 ) THEN
444 CALL XERBLA( 'chbevx_2stage', -INFO )
445 RETURN
446 ELSE IF( LQUERY ) THEN
447 RETURN
448 END IF
449*
450* Quick return if possible
451*
452 M = 0
453.EQ. IF( N0 )
454 $ RETURN
455*
456.EQ. IF( N1 ) THEN
457 M = 1
458 IF( LOWER ) THEN
459 CTMP1 = AB( 1, 1 )
460 ELSE
461 CTMP1 = AB( KD+1, 1 )
462 END IF
463 TMP1 = REAL( CTMP1 )
464 IF( VALEIG ) THEN
465.NOT..LT..AND..GE. IF( ( VLTMP1 VUTMP1 ) )
466 $ M = 0
467 END IF
468.EQ. IF( M1 ) THEN
469 W( 1 ) = REAL( CTMP1 )
470 IF( WANTZ )
471 $ Z( 1, 1 ) = CONE
472 END IF
473 RETURN
474 END IF
475*
476* Get machine constants.
477*
478 SAFMIN = SLAMCH( 'safe minimum' )
479 EPS = SLAMCH( 'precision' )
480 SMLNUM = SAFMIN / EPS
481 BIGNUM = ONE / SMLNUM
482 RMIN = SQRT( SMLNUM )
483 RMAX = MIN( SQRT( BIGNUM ), ONE / SQRT( SQRT( SAFMIN ) ) )
484*
485* Scale matrix to allowable range, if necessary.
486*
487 ISCALE = 0
488 ABSTLL = ABSTOL
489 IF( VALEIG ) THEN
490 VLL = VL
491 VUU = VU
492 ELSE
493 VLL = ZERO
494 VUU = ZERO
495 END IF
496 ANRM = CLANHB( 'm', UPLO, N, KD, AB, LDAB, RWORK )
497.GT..AND..LT. IF( ANRMZERO ANRMRMIN ) THEN
498 ISCALE = 1
499 SIGMA = RMIN / ANRM
500.GT. ELSE IF( ANRMRMAX ) THEN
501 ISCALE = 1
502 SIGMA = RMAX / ANRM
503 END IF
504.EQ. IF( ISCALE1 ) THEN
505 IF( LOWER ) THEN
506 CALL CLASCL( 'b', KD, KD, ONE, SIGMA, N, N, AB, LDAB, INFO )
507 ELSE
508 CALL CLASCL( 'q', KD, KD, ONE, SIGMA, N, N, AB, LDAB, INFO )
509 END IF
510.GT. IF( ABSTOL0 )
511 $ ABSTLL = ABSTOL*SIGMA
512 IF( VALEIG ) THEN
513 VLL = VL*SIGMA
514 VUU = VU*SIGMA
515 END IF
516 END IF
517*
518* Call CHBTRD_HB2ST to reduce Hermitian band matrix to tridiagonal form.
519*
520 INDD = 1
521 INDE = INDD + N
522 INDRWK = INDE + N
523*
524 INDHOUS = 1
525 INDWRK = INDHOUS + LHTRD
526 LLWORK = LWORK - INDWRK + 1
527*
528 CALL CHETRD_HB2ST( 'n', JOBZ, UPLO, N, KD, AB, LDAB,
529 $ RWORK( INDD ), RWORK( INDE ), WORK( INDHOUS ),
530 $ LHTRD, WORK( INDWRK ), LLWORK, IINFO )
531*
532* If all eigenvalues are desired and ABSTOL is less than or equal
533* to zero, then call SSTERF or CSTEQR. If this fails for some
534* eigenvalue, then try SSTEBZ.
535*
536 TEST = .FALSE.
537 IF (INDEIG) THEN
538.EQ..AND..EQ. IF (IL1 IUN) THEN
539 TEST = .TRUE.
540 END IF
541 END IF
542.OR..AND..LE. IF ((ALLEIG TEST) (ABSTOLZERO)) THEN
543 CALL SCOPY( N, RWORK( INDD ), 1, W, 1 )
544 INDEE = INDRWK + 2*N
545.NOT. IF( WANTZ ) THEN
546 CALL SCOPY( N-1, RWORK( INDE ), 1, RWORK( INDEE ), 1 )
547 CALL SSTERF( N, W, RWORK( INDEE ), INFO )
548 ELSE
549 CALL CLACPY( 'a', N, N, Q, LDQ, Z, LDZ )
550 CALL SCOPY( N-1, RWORK( INDE ), 1, RWORK( INDEE ), 1 )
551 CALL CSTEQR( JOBZ, N, W, RWORK( INDEE ), Z, LDZ,
552 $ RWORK( INDRWK ), INFO )
553.EQ. IF( INFO0 ) THEN
554 DO 10 I = 1, N
555 IFAIL( I ) = 0
556 10 CONTINUE
557 END IF
558 END IF
559.EQ. IF( INFO0 ) THEN
560 M = N
561 GO TO 30
562 END IF
563 INFO = 0
564 END IF
565*
566* Otherwise, call SSTEBZ and, if eigenvectors are desired, CSTEIN.
567*
568 IF( WANTZ ) THEN
569 ORDER = 'b'
570 ELSE
571 ORDER = 'e'
572 END IF
573 INDIBL = 1
574 INDISP = INDIBL + N
575 INDIWK = INDISP + N
576 CALL SSTEBZ( RANGE, ORDER, N, VLL, VUU, IL, IU, ABSTLL,
577 $ RWORK( INDD ), RWORK( INDE ), M, NSPLIT, W,
578 $ IWORK( INDIBL ), IWORK( INDISP ), RWORK( INDRWK ),
579 $ IWORK( INDIWK ), INFO )
580*
581 IF( WANTZ ) THEN
582 CALL CSTEIN( N, RWORK( INDD ), RWORK( INDE ), M, W,
583 $ IWORK( INDIBL ), IWORK( INDISP ), Z, LDZ,
584 $ RWORK( INDRWK ), IWORK( INDIWK ), IFAIL, INFO )
585*
586* Apply unitary matrix used in reduction to tridiagonal
587* form to eigenvectors returned by CSTEIN.
588*
589 DO 20 J = 1, M
590 CALL CCOPY( N, Z( 1, J ), 1, WORK( 1 ), 1 )
591 CALL CGEMV( 'n', N, N, CONE, Q, LDQ, WORK, 1, CZERO,
592 $ Z( 1, J ), 1 )
593 20 CONTINUE
594 END IF
595*
596* If matrix was scaled, then rescale eigenvalues appropriately.
597*
598 30 CONTINUE
599.EQ. IF( ISCALE1 ) THEN
600.EQ. IF( INFO0 ) THEN
601 IMAX = M
602 ELSE
603 IMAX = INFO - 1
604 END IF
605 CALL SSCAL( IMAX, ONE / SIGMA, W, 1 )
606 END IF
607*
608* If eigenvalues are not in order, then sort them, along with
609* eigenvectors.
610*
611 IF( WANTZ ) THEN
612 DO 50 J = 1, M - 1
613 I = 0
614 TMP1 = W( J )
615 DO 40 JJ = J + 1, M
616.LT. IF( W( JJ )TMP1 ) THEN
617 I = JJ
618 TMP1 = W( JJ )
619 END IF
620 40 CONTINUE
621*
622.NE. IF( I0 ) THEN
623 ITMP1 = IWORK( INDIBL+I-1 )
624 W( I ) = W( J )
625 IWORK( INDIBL+I-1 ) = IWORK( INDIBL+J-1 )
626 W( J ) = TMP1
627 IWORK( INDIBL+J-1 ) = ITMP1
628 CALL CSWAP( N, Z( 1, I ), 1, Z( 1, J ), 1 )
629.NE. IF( INFO0 ) THEN
630 ITMP1 = IFAIL( I )
631 IFAIL( I ) = IFAIL( J )
632 IFAIL( J ) = ITMP1
633 END IF
634 END IF
635 50 CONTINUE
636 END IF
637*
638* Set WORK(1) to optimal workspace size.
639*
640 WORK( 1 ) = LWMIN
641*
642 RETURN
643*
644* End of CHBEVX_2STAGE
645*
646 END
subroutine sstebz(range, order, n, vl, vu, il, iu, abstol, d, e, m, nsplit, w, iblock, isplit, work, iwork, info)
SSTEBZ
Definition sstebz.f:273
subroutine ssterf(n, d, e, info)
SSTERF
Definition ssterf.f:86
subroutine xerbla(srname, info)
XERBLA
Definition xerbla.f:60
subroutine clascl(type, kl, ku, cfrom, cto, m, n, a, lda, info)
CLASCL multiplies a general rectangular matrix by a real scalar defined as cto/cfrom.
Definition clascl.f:143
subroutine clacpy(uplo, m, n, a, lda, b, ldb)
CLACPY copies all or part of one two-dimensional array to another.
Definition clacpy.f:103
subroutine cstein(n, d, e, m, w, iblock, isplit, z, ldz, work, iwork, ifail, info)
CSTEIN
Definition cstein.f:182
subroutine chetrd_hb2st(stage1, vect, uplo, n, kd, ab, ldab, d, e, hous, lhous, work, lwork, info)
CHETRD_HB2ST reduces a complex Hermitian band matrix A to real symmetric tridiagonal form T
subroutine csteqr(compz, n, d, e, z, ldz, work, info)
CSTEQR
Definition csteqr.f:132
subroutine chbevx_2stage(jobz, range, uplo, n, kd, ab, ldab, q, ldq, vl, vu, il, iu, abstol, m, w, z, ldz, work, lwork, rwork, iwork, ifail, info)
CHBEVX_2STAGE computes the eigenvalues and, optionally, the left and/or right eigenvectors for OTHER ...
subroutine cswap(n, cx, incx, cy, incy)
CSWAP
Definition cswap.f:81
subroutine ccopy(n, cx, incx, cy, incy)
CCOPY
Definition ccopy.f:81
subroutine cgemv(trans, m, n, alpha, a, lda, x, incx, beta, y, incy)
CGEMV
Definition cgemv.f:158
subroutine sscal(n, sa, sx, incx)
SSCAL
Definition sscal.f:79
subroutine scopy(n, sx, incx, sy, incy)
SCOPY
Definition scopy.f:82
#define min(a, b)
Definition macros.h:20
#define max(a, b)
Definition macros.h:21