LAPACK 3.12.1
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches

◆ dstebz()

subroutine dstebz ( character range,
character order,
integer n,
double precision vl,
double precision vu,
integer il,
integer iu,
double precision abstol,
double precision, dimension( * ) d,
double precision, dimension( * ) e,
integer m,
integer nsplit,
double precision, dimension( * ) w,
integer, dimension( * ) iblock,
integer, dimension( * ) isplit,
double precision, dimension( * ) work,
integer, dimension( * ) iwork,
integer info )

DSTEBZ

Download DSTEBZ + dependencies [TGZ] [ZIP] [TXT]

Purpose:
!>
!> DSTEBZ computes the eigenvalues of a symmetric tridiagonal
!> matrix T.  The user may ask for all eigenvalues, all eigenvalues
!> in the half-open interval (VL, VU], or the IL-th through IU-th
!> eigenvalues.
!>
!> To avoid overflow, the matrix must be scaled so that its
!> largest element is no greater than overflow**(1/2) * underflow**(1/4) in absolute value, and for greatest
!> accuracy, it should not be much smaller than that.
!>
!> See W. Kahan , Report CS41, Computer Science Dept., Stanford
!> University, July 21, 1966.
!> 
Parameters
[in]RANGE
!>          RANGE is 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 (of the
!>                           entire matrix) will be found.
!> 
[in]ORDER
!>          ORDER is CHARACTER*1
!>          = 'B': () the eigenvalues will be grouped by
!>                              split-off block (see IBLOCK, ISPLIT) and
!>                              ordered from smallest to largest within
!>                              the block.
!>          = 'E': ()
!>                              the eigenvalues for the entire matrix
!>                              will be ordered from smallest to
!>                              largest.
!> 
[in]N
!>          N is INTEGER
!>          The order of the tridiagonal matrix T.  N >= 0.
!> 
[in]VL
!>          VL is DOUBLE PRECISION
!>
!>          If RANGE='V', the lower bound of the interval to
!>          be searched for eigenvalues.  Eigenvalues less than or equal
!>          to VL, or greater than VU, will not be returned.  VL < VU.
!>          Not referenced if RANGE = 'A' or 'I'.
!> 
[in]VU
!>          VU is DOUBLE PRECISION
!>
!>          If RANGE='V', the upper bound of the interval to
!>          be searched for eigenvalues.  Eigenvalues less than or equal
!>          to VL, or greater than VU, will not be returned.  VL < VU.
!>          Not referenced if RANGE = 'A' or 'I'.
!> 
[in]IL
!>          IL is INTEGER
!>
!>          If RANGE='I', the index of the
!>          smallest eigenvalue 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'.
!> 
[in]IU
!>          IU is INTEGER
!>
!>          If RANGE='I', the index of the
!>          largest eigenvalue 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'.
!> 
[in]ABSTOL
!>          ABSTOL is DOUBLE PRECISION
!>          The absolute tolerance for the eigenvalues.  An eigenvalue
!>          (or cluster) is considered to be located if it has been
!>          determined to lie in an interval whose width is ABSTOL or
!>          less.  If ABSTOL is less than or equal to zero, then ULP*|T|
!>          will be used, where |T| means the 1-norm of T.
!>
!>          Eigenvalues will be computed most accurately when ABSTOL is
!>          set to twice the underflow threshold 2*DLAMCH('S'), not zero.
!> 
[in]D
!>          D is DOUBLE PRECISION array, dimension (N)
!>          The n diagonal elements of the tridiagonal matrix T.
!> 
[in]E
!>          E is DOUBLE PRECISION array, dimension (N-1)
!>          The (n-1) off-diagonal elements of the tridiagonal matrix T.
!> 
[out]M
!>          M is INTEGER
!>          The actual number of eigenvalues found. 0 <= M <= N.
!>          (See also the description of INFO=2,3.)
!> 
[out]NSPLIT
!>          NSPLIT is INTEGER
!>          The number of diagonal blocks in the matrix T.
!>          1 <= NSPLIT <= N.
!> 
[out]W
!>          W is DOUBLE PRECISION array, dimension (N)
!>          On exit, the first M elements of W will contain the
!>          eigenvalues.  (DSTEBZ may use the remaining N-M elements as
!>          workspace.)
!> 
[out]IBLOCK
!>          IBLOCK is INTEGER array, dimension (N)
!>          At each row/column j where E(j) is zero or small, the
!>          matrix T is considered to split into a block diagonal
!>          matrix.  On exit, if INFO = 0, IBLOCK(i) specifies to which
!>          block (from 1 to the number of blocks) the eigenvalue W(i)
!>          belongs.  (DSTEBZ may use the remaining N-M elements as
!>          workspace.)
!> 
[out]ISPLIT
!>          ISPLIT is INTEGER array, dimension (N)
!>          The splitting points, at which T breaks up into submatrices.
!>          The first submatrix consists of rows/columns 1 to ISPLIT(1),
!>          the second of rows/columns ISPLIT(1)+1 through ISPLIT(2),
!>          etc., and the NSPLIT-th consists of rows/columns
!>          ISPLIT(NSPLIT-1)+1 through ISPLIT(NSPLIT)=N.
!>          (Only the first NSPLIT elements will actually be used, but
!>          since the user cannot know a priori what value NSPLIT will
!>          have, N words must be reserved for ISPLIT.)
!> 
[out]WORK
!>          WORK is DOUBLE PRECISION array, dimension (4*N)
!> 
[out]IWORK
!>          IWORK is INTEGER array, dimension (3*N)
!> 
[out]INFO
!>          INFO is INTEGER
!>          = 0:  successful exit
!>          < 0:  if INFO = -i, the i-th argument had an illegal value
!>          > 0:  some or all of the eigenvalues failed to converge or
!>                were not computed:
!>                =1 or 3: Bisection failed to converge for some
!>                        eigenvalues; these eigenvalues are flagged by a
!>                        negative block number.  The effect is that the
!>                        eigenvalues may not be as accurate as the
!>                        absolute and relative tolerances.  This is
!>                        generally caused by unexpectedly inaccurate
!>                        arithmetic.
!>                =2 or 3: RANGE='I' only: Not all of the eigenvalues
!>                        IL:IU were found.
!>                        Effect: M < IU+1-IL
!>                        Cause:  non-monotonic arithmetic, causing the
!>                                Sturm sequence to be non-monotonic.
!>                        Cure:   recalculate, using RANGE='A', and pick
!>                                out eigenvalues IL:IU.  In some cases,
!>                                increasing the PARAMETER  may
!>                                make things work.
!>                = 4:    RANGE='I', and the Gershgorin interval
!>                        initially used was too small.  No eigenvalues
!>                        were computed.
!>                        Probable cause: your machine has sloppy
!>                                        floating-point arithmetic.
!>                        Cure: Increase the PARAMETER ,
!>                              recompile, and try again.
!> 
Internal Parameters:
!>  RELFAC  DOUBLE PRECISION, default = 2.0e0
!>          The relative tolerance.  An interval (a,b] lies within
!>           if  b-a < RELFAC*ulp*max(|a|,|b|),
!>          where  is the machine precision (distance from 1 to
!>          the next larger floating point number.)
!>
!>  FUDGE   DOUBLE PRECISION, default = 2
!>          A  to widen the Gershgorin intervals.  Ideally,
!>          a value of 1 should work, but on machines with sloppy
!>          arithmetic, this needs to be larger.  The default for
!>          publicly released versions should be large enough to handle
!>          the worst machine around.  Note that this has no effect
!>          on accuracy of the solution.
!> 
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.

Definition at line 268 of file dstebz.f.

272*
273* -- LAPACK computational routine --
274* -- LAPACK is a software package provided by Univ. of Tennessee, --
275* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
276*
277* .. Scalar Arguments ..
278 CHARACTER ORDER, RANGE
279 INTEGER IL, INFO, IU, M, N, NSPLIT
280 DOUBLE PRECISION ABSTOL, VL, VU
281* ..
282* .. Array Arguments ..
283 INTEGER IBLOCK( * ), ISPLIT( * ), IWORK( * )
284 DOUBLE PRECISION D( * ), E( * ), W( * ), WORK( * )
285* ..
286*
287* =====================================================================
288*
289* .. Parameters ..
290 DOUBLE PRECISION ZERO, ONE, TWO, HALF
291 parameter( zero = 0.0d0, one = 1.0d0, two = 2.0d0,
292 $ half = 1.0d0 / two )
293 DOUBLE PRECISION FUDGE, RELFAC
294 parameter( fudge = 2.1d0, relfac = 2.0d0 )
295* ..
296* .. Local Scalars ..
297 LOGICAL NCNVRG, TOOFEW
298 INTEGER IB, IBEGIN, IDISCL, IDISCU, IE, IEND, IINFO,
299 $ IM, IN, IOFF, IORDER, IOUT, IRANGE, ITMAX,
300 $ ITMP1, IW, IWOFF, J, JB, JDISC, JE, NB, NWL,
301 $ NWU
302 DOUBLE PRECISION ATOLI, BNORM, GL, GU, PIVMIN, RTOLI, SAFEMN,
303 $ TMP1, TMP2, TNORM, ULP, WKILL, WL, WLU, WU, WUL
304* ..
305* .. Local Arrays ..
306 INTEGER IDUMMA( 1 )
307* ..
308* .. External Functions ..
309 LOGICAL LSAME
310 INTEGER ILAENV
311 DOUBLE PRECISION DLAMCH
312 EXTERNAL lsame, ilaenv, dlamch
313* ..
314* .. External Subroutines ..
315 EXTERNAL dlaebz, xerbla
316* ..
317* .. Intrinsic Functions ..
318 INTRINSIC abs, int, log, max, min, sqrt
319* ..
320* .. Executable Statements ..
321*
322 info = 0
323*
324* Decode RANGE
325*
326 IF( lsame( range, 'A' ) ) THEN
327 irange = 1
328 ELSE IF( lsame( range, 'V' ) ) THEN
329 irange = 2
330 ELSE IF( lsame( range, 'I' ) ) THEN
331 irange = 3
332 ELSE
333 irange = 0
334 END IF
335*
336* Decode ORDER
337*
338 IF( lsame( order, 'B' ) ) THEN
339 iorder = 2
340 ELSE IF( lsame( order, 'E' ) ) THEN
341 iorder = 1
342 ELSE
343 iorder = 0
344 END IF
345*
346* Check for Errors
347*
348 IF( irange.LE.0 ) THEN
349 info = -1
350 ELSE IF( iorder.LE.0 ) THEN
351 info = -2
352 ELSE IF( n.LT.0 ) THEN
353 info = -3
354 ELSE IF( irange.EQ.2 ) THEN
355 IF( vl.GE.vu )
356 $ info = -5
357 ELSE IF( irange.EQ.3 .AND. ( il.LT.1 .OR. il.GT.max( 1, n ) ) )
358 $ THEN
359 info = -6
360 ELSE IF( irange.EQ.3 .AND. ( iu.LT.min( n, il ) .OR. iu.GT.n ) )
361 $ THEN
362 info = -7
363 END IF
364*
365 IF( info.NE.0 ) THEN
366 CALL xerbla( 'DSTEBZ', -info )
367 RETURN
368 END IF
369*
370* Initialize error flags
371*
372 info = 0
373 ncnvrg = .false.
374 toofew = .false.
375*
376* Quick return if possible
377*
378 m = 0
379 IF( n.EQ.0 )
380 $ RETURN
381*
382* Simplifications:
383*
384 IF( irange.EQ.3 .AND. il.EQ.1 .AND. iu.EQ.n )
385 $ irange = 1
386*
387* Get machine constants
388* NB is the minimum vector length for vector bisection, or 0
389* if only scalar is to be done.
390*
391 safemn = dlamch( 'S' )
392 ulp = dlamch( 'P' )
393 rtoli = ulp*relfac
394 nb = ilaenv( 1, 'DSTEBZ', ' ', n, -1, -1, -1 )
395 IF( nb.LE.1 )
396 $ nb = 0
397*
398* Special Case when N=1
399*
400 IF( n.EQ.1 ) THEN
401 nsplit = 1
402 isplit( 1 ) = 1
403 IF( irange.EQ.2 .AND. ( vl.GE.d( 1 ) .OR. vu.LT.d( 1 ) ) ) THEN
404 m = 0
405 ELSE
406 w( 1 ) = d( 1 )
407 iblock( 1 ) = 1
408 m = 1
409 END IF
410 RETURN
411 END IF
412*
413* Compute Splitting Points
414*
415 nsplit = 1
416 work( n ) = zero
417 pivmin = one
418*
419 DO 10 j = 2, n
420 tmp1 = e( j-1 )**2
421 IF( abs( d( j )*d( j-1 ) )*ulp**2+safemn.GT.tmp1 ) THEN
422 isplit( nsplit ) = j - 1
423 nsplit = nsplit + 1
424 work( j-1 ) = zero
425 ELSE
426 work( j-1 ) = tmp1
427 pivmin = max( pivmin, tmp1 )
428 END IF
429 10 CONTINUE
430 isplit( nsplit ) = n
431 pivmin = pivmin*safemn
432*
433* Compute Interval and ATOLI
434*
435 IF( irange.EQ.3 ) THEN
436*
437* RANGE='I': Compute the interval containing eigenvalues
438* IL through IU.
439*
440* Compute Gershgorin interval for entire (split) matrix
441* and use it as the initial interval
442*
443 gu = d( 1 )
444 gl = d( 1 )
445 tmp1 = zero
446*
447 DO 20 j = 1, n - 1
448 tmp2 = sqrt( work( j ) )
449 gu = max( gu, d( j )+tmp1+tmp2 )
450 gl = min( gl, d( j )-tmp1-tmp2 )
451 tmp1 = tmp2
452 20 CONTINUE
453*
454 gu = max( gu, d( n )+tmp1 )
455 gl = min( gl, d( n )-tmp1 )
456 tnorm = max( abs( gl ), abs( gu ) )
457 gl = gl - fudge*tnorm*ulp*n - fudge*two*pivmin
458 gu = gu + fudge*tnorm*ulp*n + fudge*pivmin
459*
460* Compute Iteration parameters
461*
462 itmax = int( ( log( tnorm+pivmin )-log( pivmin ) ) /
463 $ log( two ) ) + 2
464 IF( abstol.LE.zero ) THEN
465 atoli = ulp*tnorm
466 ELSE
467 atoli = abstol
468 END IF
469*
470 work( n+1 ) = gl
471 work( n+2 ) = gl
472 work( n+3 ) = gu
473 work( n+4 ) = gu
474 work( n+5 ) = gl
475 work( n+6 ) = gu
476 iwork( 1 ) = -1
477 iwork( 2 ) = -1
478 iwork( 3 ) = n + 1
479 iwork( 4 ) = n + 1
480 iwork( 5 ) = il - 1
481 iwork( 6 ) = iu
482*
483 CALL dlaebz( 3, itmax, n, 2, 2, nb, atoli, rtoli, pivmin, d,
484 $ e,
485 $ work, iwork( 5 ), work( n+1 ), work( n+5 ), iout,
486 $ iwork, w, iblock, iinfo )
487*
488 IF( iwork( 6 ).EQ.iu ) THEN
489 wl = work( n+1 )
490 wlu = work( n+3 )
491 nwl = iwork( 1 )
492 wu = work( n+4 )
493 wul = work( n+2 )
494 nwu = iwork( 4 )
495 ELSE
496 wl = work( n+2 )
497 wlu = work( n+4 )
498 nwl = iwork( 2 )
499 wu = work( n+3 )
500 wul = work( n+1 )
501 nwu = iwork( 3 )
502 END IF
503*
504 IF( nwl.LT.0 .OR. nwl.GE.n .OR. nwu.LT.1 .OR. nwu.GT.n ) THEN
505 info = 4
506 RETURN
507 END IF
508 ELSE
509*
510* RANGE='A' or 'V' -- Set ATOLI
511*
512 tnorm = max( abs( d( 1 ) )+abs( e( 1 ) ),
513 $ abs( d( n ) )+abs( e( n-1 ) ) )
514*
515 DO 30 j = 2, n - 1
516 tnorm = max( tnorm, abs( d( j ) )+abs( e( j-1 ) )+
517 $ abs( e( j ) ) )
518 30 CONTINUE
519*
520 IF( abstol.LE.zero ) THEN
521 atoli = ulp*tnorm
522 ELSE
523 atoli = abstol
524 END IF
525*
526 IF( irange.EQ.2 ) THEN
527 wl = vl
528 wu = vu
529 ELSE
530 wl = zero
531 wu = zero
532 END IF
533 END IF
534*
535* Find Eigenvalues -- Loop Over Blocks and recompute NWL and NWU.
536* NWL accumulates the number of eigenvalues .le. WL,
537* NWU accumulates the number of eigenvalues .le. WU
538*
539 m = 0
540 iend = 0
541 info = 0
542 nwl = 0
543 nwu = 0
544*
545 DO 70 jb = 1, nsplit
546 ioff = iend
547 ibegin = ioff + 1
548 iend = isplit( jb )
549 in = iend - ioff
550*
551 IF( in.EQ.1 ) THEN
552*
553* Special Case -- IN=1
554*
555 IF( irange.EQ.1 .OR. wl.GE.d( ibegin )-pivmin )
556 $ nwl = nwl + 1
557 IF( irange.EQ.1 .OR. wu.GE.d( ibegin )-pivmin )
558 $ nwu = nwu + 1
559 IF( irange.EQ.1 .OR. ( wl.LT.d( ibegin )-pivmin .AND. wu.GE.
560 $ d( ibegin )-pivmin ) ) THEN
561 m = m + 1
562 w( m ) = d( ibegin )
563 iblock( m ) = jb
564 END IF
565 ELSE
566*
567* General Case -- IN > 1
568*
569* Compute Gershgorin Interval
570* and use it as the initial interval
571*
572 gu = d( ibegin )
573 gl = d( ibegin )
574 tmp1 = zero
575*
576 DO 40 j = ibegin, iend - 1
577 tmp2 = abs( e( j ) )
578 gu = max( gu, d( j )+tmp1+tmp2 )
579 gl = min( gl, d( j )-tmp1-tmp2 )
580 tmp1 = tmp2
581 40 CONTINUE
582*
583 gu = max( gu, d( iend )+tmp1 )
584 gl = min( gl, d( iend )-tmp1 )
585 bnorm = max( abs( gl ), abs( gu ) )
586 gl = gl - fudge*bnorm*ulp*in - fudge*pivmin
587 gu = gu + fudge*bnorm*ulp*in + fudge*pivmin
588*
589* Compute ATOLI for the current submatrix
590*
591 IF( abstol.LE.zero ) THEN
592 atoli = ulp*max( abs( gl ), abs( gu ) )
593 ELSE
594 atoli = abstol
595 END IF
596*
597 IF( irange.GT.1 ) THEN
598 IF( gu.LT.wl ) THEN
599 nwl = nwl + in
600 nwu = nwu + in
601 GO TO 70
602 END IF
603 gl = max( gl, wl )
604 gu = min( gu, wu )
605 IF( gl.GE.gu )
606 $ GO TO 70
607 END IF
608*
609* Set Up Initial Interval
610*
611 work( n+1 ) = gl
612 work( n+in+1 ) = gu
613 CALL dlaebz( 1, 0, in, in, 1, nb, atoli, rtoli, pivmin,
614 $ d( ibegin ), e( ibegin ), work( ibegin ),
615 $ idumma, work( n+1 ), work( n+2*in+1 ), im,
616 $ iwork, w( m+1 ), iblock( m+1 ), iinfo )
617*
618 nwl = nwl + iwork( 1 )
619 nwu = nwu + iwork( in+1 )
620 iwoff = m - iwork( 1 )
621*
622* Compute Eigenvalues
623*
624 itmax = int( ( log( gu-gl+pivmin )-log( pivmin ) ) /
625 $ log( two ) ) + 2
626 CALL dlaebz( 2, itmax, in, in, 1, nb, atoli, rtoli,
627 $ pivmin,
628 $ d( ibegin ), e( ibegin ), work( ibegin ),
629 $ idumma, work( n+1 ), work( n+2*in+1 ), iout,
630 $ iwork, w( m+1 ), iblock( m+1 ), iinfo )
631*
632* Copy Eigenvalues Into W and IBLOCK
633* Use -JB for block number for unconverged eigenvalues.
634*
635 DO 60 j = 1, iout
636 tmp1 = half*( work( j+n )+work( j+in+n ) )
637*
638* Flag non-convergence.
639*
640 IF( j.GT.iout-iinfo ) THEN
641 ncnvrg = .true.
642 ib = -jb
643 ELSE
644 ib = jb
645 END IF
646 DO 50 je = iwork( j ) + 1 + iwoff,
647 $ iwork( j+in ) + iwoff
648 w( je ) = tmp1
649 iblock( je ) = ib
650 50 CONTINUE
651 60 CONTINUE
652*
653 m = m + im
654 END IF
655 70 CONTINUE
656*
657* If RANGE='I', then (WL,WU) contains eigenvalues NWL+1,...,NWU
658* If NWL+1 < IL or NWU > IU, discard extra eigenvalues.
659*
660 IF( irange.EQ.3 ) THEN
661 im = 0
662 idiscl = il - 1 - nwl
663 idiscu = nwu - iu
664*
665 IF( idiscl.GT.0 .OR. idiscu.GT.0 ) THEN
666 DO 80 je = 1, m
667 IF( w( je ).LE.wlu .AND. idiscl.GT.0 ) THEN
668 idiscl = idiscl - 1
669 ELSE IF( w( je ).GE.wul .AND. idiscu.GT.0 ) THEN
670 idiscu = idiscu - 1
671 ELSE
672 im = im + 1
673 w( im ) = w( je )
674 iblock( im ) = iblock( je )
675 END IF
676 80 CONTINUE
677 m = im
678 END IF
679 IF( idiscl.GT.0 .OR. idiscu.GT.0 ) THEN
680*
681* Code to deal with effects of bad arithmetic:
682* Some low eigenvalues to be discarded are not in (WL,WLU],
683* or high eigenvalues to be discarded are not in (WUL,WU]
684* so just kill off the smallest IDISCL/largest IDISCU
685* eigenvalues, by simply finding the smallest/largest
686* eigenvalue(s).
687*
688* (If N(w) is monotone non-decreasing, this should never
689* happen.)
690*
691 IF( idiscl.GT.0 ) THEN
692 wkill = wu
693 DO 100 jdisc = 1, idiscl
694 iw = 0
695 DO 90 je = 1, m
696 IF( iblock( je ).NE.0 .AND.
697 $ ( w( je ).LT.wkill .OR. iw.EQ.0 ) ) THEN
698 iw = je
699 wkill = w( je )
700 END IF
701 90 CONTINUE
702 iblock( iw ) = 0
703 100 CONTINUE
704 END IF
705 IF( idiscu.GT.0 ) THEN
706*
707 wkill = wl
708 DO 120 jdisc = 1, idiscu
709 iw = 0
710 DO 110 je = 1, m
711 IF( iblock( je ).NE.0 .AND.
712 $ ( w( je ).GT.wkill .OR. iw.EQ.0 ) ) THEN
713 iw = je
714 wkill = w( je )
715 END IF
716 110 CONTINUE
717 iblock( iw ) = 0
718 120 CONTINUE
719 END IF
720 im = 0
721 DO 130 je = 1, m
722 IF( iblock( je ).NE.0 ) THEN
723 im = im + 1
724 w( im ) = w( je )
725 iblock( im ) = iblock( je )
726 END IF
727 130 CONTINUE
728 m = im
729 END IF
730 IF( idiscl.LT.0 .OR. idiscu.LT.0 ) THEN
731 toofew = .true.
732 END IF
733 END IF
734*
735* If ORDER='B', do nothing -- the eigenvalues are already sorted
736* by block.
737* If ORDER='E', sort the eigenvalues from smallest to largest
738*
739 IF( iorder.EQ.1 .AND. nsplit.GT.1 ) THEN
740 DO 150 je = 1, m - 1
741 ie = 0
742 tmp1 = w( je )
743 DO 140 j = je + 1, m
744 IF( w( j ).LT.tmp1 ) THEN
745 ie = j
746 tmp1 = w( j )
747 END IF
748 140 CONTINUE
749*
750 IF( ie.NE.0 ) THEN
751 itmp1 = iblock( ie )
752 w( ie ) = w( je )
753 iblock( ie ) = iblock( je )
754 w( je ) = tmp1
755 iblock( je ) = itmp1
756 END IF
757 150 CONTINUE
758 END IF
759*
760 info = 0
761 IF( ncnvrg )
762 $ info = info + 1
763 IF( toofew )
764 $ info = info + 2
765 RETURN
766*
767* End of DSTEBZ
768*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
integer function ilaenv(ispec, name, opts, n1, n2, n3, n4)
ILAENV
Definition ilaenv.f:160
subroutine dlaebz(ijob, nitmax, n, mmax, minp, nbmin, abstol, reltol, pivmin, d, e, e2, nval, ab, c, mout, nab, work, iwork, info)
DLAEBZ computes the number of eigenvalues of a real symmetric tridiagonal matrix which are less than ...
Definition dlaebz.f:317
double precision function dlamch(cmach)
DLAMCH
Definition dlamch.f:69
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
Here is the call graph for this function:
Here is the caller graph for this function: