File:  [local] / rpl / lapack / lapack / dlaqz0.f
Revision 1.1: download - view: text, annotated - select for diffs - revision graph
Mon Aug 7 08:55:29 2023 UTC (9 months, 1 week ago) by bertrand
Branches: MAIN
CVS tags: rpl-4_1_35, rpl-4_1_34, HEAD
Ajout de fichiers de lapack 3.11

    1: *> \brief \b DLAQZ0
    2: *
    3: *  =========== DOCUMENTATION ===========
    4: *
    5: * Online html documentation available at
    6: *            http://www.netlib.org/lapack/explore-html/
    7: *
    8: *> \htmlonly
    9: *> Download DLAQZ0 + dependencies
   10: *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dlaqz0.f">
   11: *> [TGZ]</a>
   12: *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dlaqz0.f">
   13: *> [ZIP]</a>
   14: *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dlaqz0.f">
   15: *> [TXT]</a>
   16: *> \endhtmlonly
   17: *
   18: *  Definition:
   19: *  ===========
   20: *
   21: *     RECURSIVE SUBROUTINE DLAQZ0( WANTS, WANTQ, WANTZ, N, ILO, IHI, A,
   22: *    $                             LDA, B, LDB, ALPHAR, ALPHAI, BETA,
   23: *    $                             Q, LDQ, Z, LDZ, WORK, LWORK, REC,
   24: *    $                             INFO )
   25: *     IMPLICIT NONE
   26: *
   27: *     Arguments
   28: *     CHARACTER, INTENT( IN ) :: WANTS, WANTQ, WANTZ
   29: *     INTEGER, INTENT( IN ) :: N, ILO, IHI, LDA, LDB, LDQ, LDZ, LWORK,
   30: *    $         REC
   31: *
   32: *     INTEGER, INTENT( OUT ) :: INFO
   33: *
   34: *     DOUBLE PRECISION, INTENT( INOUT ) :: A( LDA, * ), B( LDB, * ),
   35: *    $                  Q( LDQ, * ), Z( LDZ, * ), ALPHAR( * ),
   36: *    $                  ALPHAI( * ), BETA( * ), WORK( * )
   37: *      ..
   38: *
   39: *
   40: *> \par Purpose:
   41: *  =============
   42: *>
   43: *> \verbatim
   44: *>
   45: *> DLAQZ0 computes the eigenvalues of a real matrix pair (H,T),
   46: *> where H is an upper Hessenberg matrix and T is upper triangular,
   47: *> using the double-shift QZ method.
   48: *> Matrix pairs of this type are produced by the reduction to
   49: *> generalized upper Hessenberg form of a real matrix pair (A,B):
   50: *>
   51: *>    A = Q1*H*Z1**T,  B = Q1*T*Z1**T,
   52: *>
   53: *> as computed by DGGHRD.
   54: *>
   55: *> If JOB='S', then the Hessenberg-triangular pair (H,T) is
   56: *> also reduced to generalized Schur form,
   57: *>
   58: *>    H = Q*S*Z**T,  T = Q*P*Z**T,
   59: *>
   60: *> where Q and Z are orthogonal matrices, P is an upper triangular
   61: *> matrix, and S is a quasi-triangular matrix with 1-by-1 and 2-by-2
   62: *> diagonal blocks.
   63: *>
   64: *> The 1-by-1 blocks correspond to real eigenvalues of the matrix pair
   65: *> (H,T) and the 2-by-2 blocks correspond to complex conjugate pairs of
   66: *> eigenvalues.
   67: *>
   68: *> Additionally, the 2-by-2 upper triangular diagonal blocks of P
   69: *> corresponding to 2-by-2 blocks of S are reduced to positive diagonal
   70: *> form, i.e., if S(j+1,j) is non-zero, then P(j+1,j) = P(j,j+1) = 0,
   71: *> P(j,j) > 0, and P(j+1,j+1) > 0.
   72: *>
   73: *> Optionally, the orthogonal matrix Q from the generalized Schur
   74: *> factorization may be postmultiplied into an input matrix Q1, and the
   75: *> orthogonal matrix Z may be postmultiplied into an input matrix Z1.
   76: *> If Q1 and Z1 are the orthogonal matrices from DGGHRD that reduced
   77: *> the matrix pair (A,B) to generalized upper Hessenberg form, then the
   78: *> output matrices Q1*Q and Z1*Z are the orthogonal factors from the
   79: *> generalized Schur factorization of (A,B):
   80: *>
   81: *>    A = (Q1*Q)*S*(Z1*Z)**T,  B = (Q1*Q)*P*(Z1*Z)**T.
   82: *>
   83: *> To avoid overflow, eigenvalues of the matrix pair (H,T) (equivalently,
   84: *> of (A,B)) are computed as a pair of values (alpha,beta), where alpha is
   85: *> complex and beta real.
   86: *> If beta is nonzero, lambda = alpha / beta is an eigenvalue of the
   87: *> generalized nonsymmetric eigenvalue problem (GNEP)
   88: *>    A*x = lambda*B*x
   89: *> and if alpha is nonzero, mu = beta / alpha is an eigenvalue of the
   90: *> alternate form of the GNEP
   91: *>    mu*A*y = B*y.
   92: *> Real eigenvalues can be read directly from the generalized Schur
   93: *> form:
   94: *>   alpha = S(i,i), beta = P(i,i).
   95: *>
   96: *> Ref: C.B. Moler & G.W. Stewart, "An Algorithm for Generalized Matrix
   97: *>      Eigenvalue Problems", SIAM J. Numer. Anal., 10(1973),
   98: *>      pp. 241--256.
   99: *>
  100: *> Ref: B. Kagstrom, D. Kressner, "Multishift Variants of the QZ
  101: *>      Algorithm with Aggressive Early Deflation", SIAM J. Numer.
  102: *>      Anal., 29(2006), pp. 199--227.
  103: *>
  104: *> Ref: T. Steel, D. Camps, K. Meerbergen, R. Vandebril "A multishift,
  105: *>      multipole rational QZ method with agressive early deflation"
  106: *> \endverbatim
  107: *
  108: *  Arguments:
  109: *  ==========
  110: *
  111: *> \param[in] WANTS
  112: *> \verbatim
  113: *>          WANTS is CHARACTER*1
  114: *>          = 'E': Compute eigenvalues only;
  115: *>          = 'S': Compute eigenvalues and the Schur form.
  116: *> \endverbatim
  117: *>
  118: *> \param[in] WANTQ
  119: *> \verbatim
  120: *>          WANTQ is CHARACTER*1
  121: *>          = 'N': Left Schur vectors (Q) are not computed;
  122: *>          = 'I': Q is initialized to the unit matrix and the matrix Q
  123: *>                 of left Schur vectors of (A,B) is returned;
  124: *>          = 'V': Q must contain an orthogonal matrix Q1 on entry and
  125: *>                 the product Q1*Q is returned.
  126: *> \endverbatim
  127: *>
  128: *> \param[in] WANTZ
  129: *> \verbatim
  130: *>          WANTZ is CHARACTER*1
  131: *>          = 'N': Right Schur vectors (Z) are not computed;
  132: *>          = 'I': Z is initialized to the unit matrix and the matrix Z
  133: *>                 of right Schur vectors of (A,B) is returned;
  134: *>          = 'V': Z must contain an orthogonal matrix Z1 on entry and
  135: *>                 the product Z1*Z is returned.
  136: *> \endverbatim
  137: *>
  138: *> \param[in] N
  139: *> \verbatim
  140: *>          N is INTEGER
  141: *>          The order of the matrices A, B, Q, and Z.  N >= 0.
  142: *> \endverbatim
  143: *>
  144: *> \param[in] ILO
  145: *> \verbatim
  146: *>          ILO is INTEGER
  147: *> \endverbatim
  148: *>
  149: *> \param[in] IHI
  150: *> \verbatim
  151: *>          IHI is INTEGER
  152: *>          ILO and IHI mark the rows and columns of A which are in
  153: *>          Hessenberg form.  It is assumed that A is already upper
  154: *>          triangular in rows and columns 1:ILO-1 and IHI+1:N.
  155: *>          If N > 0, 1 <= ILO <= IHI <= N; if N = 0, ILO=1 and IHI=0.
  156: *> \endverbatim
  157: *>
  158: *> \param[in,out] A
  159: *> \verbatim
  160: *>          A is DOUBLE PRECISION array, dimension (LDA, N)
  161: *>          On entry, the N-by-N upper Hessenberg matrix A.
  162: *>          On exit, if JOB = 'S', A contains the upper quasi-triangular
  163: *>          matrix S from the generalized Schur factorization.
  164: *>          If JOB = 'E', the diagonal blocks of A match those of S, but
  165: *>          the rest of A is unspecified.
  166: *> \endverbatim
  167: *>
  168: *> \param[in] LDA
  169: *> \verbatim
  170: *>          LDA is INTEGER
  171: *>          The leading dimension of the array A.  LDA >= max( 1, N ).
  172: *> \endverbatim
  173: *>
  174: *> \param[in,out] B
  175: *> \verbatim
  176: *>          B is DOUBLE PRECISION array, dimension (LDB, N)
  177: *>          On entry, the N-by-N upper triangular matrix B.
  178: *>          On exit, if JOB = 'S', B contains the upper triangular
  179: *>          matrix P from the generalized Schur factorization;
  180: *>          2-by-2 diagonal blocks of P corresponding to 2-by-2 blocks of S
  181: *>          are reduced to positive diagonal form, i.e., if A(j+1,j) is
  182: *>          non-zero, then B(j+1,j) = B(j,j+1) = 0, B(j,j) > 0, and
  183: *>          B(j+1,j+1) > 0.
  184: *>          If JOB = 'E', the diagonal blocks of B match those of P, but
  185: *>          the rest of B is unspecified.
  186: *> \endverbatim
  187: *>
  188: *> \param[in] LDB
  189: *> \verbatim
  190: *>          LDB is INTEGER
  191: *>          The leading dimension of the array B.  LDB >= max( 1, N ).
  192: *> \endverbatim
  193: *>
  194: *> \param[out] ALPHAR
  195: *> \verbatim
  196: *>          ALPHAR is DOUBLE PRECISION array, dimension (N)
  197: *>          The real parts of each scalar alpha defining an eigenvalue
  198: *>          of GNEP.
  199: *> \endverbatim
  200: *>
  201: *> \param[out] ALPHAI
  202: *> \verbatim
  203: *>          ALPHAI is DOUBLE PRECISION array, dimension (N)
  204: *>          The imaginary parts of each scalar alpha defining an
  205: *>          eigenvalue of GNEP.
  206: *>          If ALPHAI(j) is zero, then the j-th eigenvalue is real; if
  207: *>          positive, then the j-th and (j+1)-st eigenvalues are a
  208: *>          complex conjugate pair, with ALPHAI(j+1) = -ALPHAI(j).
  209: *> \endverbatim
  210: *>
  211: *> \param[out] BETA
  212: *> \verbatim
  213: *>          BETA is DOUBLE PRECISION array, dimension (N)
  214: *>          The scalars beta that define the eigenvalues of GNEP.
  215: *>          Together, the quantities alpha = (ALPHAR(j),ALPHAI(j)) and
  216: *>          beta = BETA(j) represent the j-th eigenvalue of the matrix
  217: *>          pair (A,B), in one of the forms lambda = alpha/beta or
  218: *>          mu = beta/alpha.  Since either lambda or mu may overflow,
  219: *>          they should not, in general, be computed.
  220: *> \endverbatim
  221: *>
  222: *> \param[in,out] Q
  223: *> \verbatim
  224: *>          Q is DOUBLE PRECISION array, dimension (LDQ, N)
  225: *>          On entry, if COMPQ = 'V', the orthogonal matrix Q1 used in
  226: *>          the reduction of (A,B) to generalized Hessenberg form.
  227: *>          On exit, if COMPQ = 'I', the orthogonal matrix of left Schur
  228: *>          vectors of (A,B), and if COMPQ = 'V', the orthogonal matrix
  229: *>          of left Schur vectors of (A,B).
  230: *>          Not referenced if COMPQ = 'N'.
  231: *> \endverbatim
  232: *>
  233: *> \param[in] LDQ
  234: *> \verbatim
  235: *>          LDQ is INTEGER
  236: *>          The leading dimension of the array Q.  LDQ >= 1.
  237: *>          If COMPQ='V' or 'I', then LDQ >= N.
  238: *> \endverbatim
  239: *>
  240: *> \param[in,out] Z
  241: *> \verbatim
  242: *>          Z is DOUBLE PRECISION array, dimension (LDZ, N)
  243: *>          On entry, if COMPZ = 'V', the orthogonal matrix Z1 used in
  244: *>          the reduction of (A,B) to generalized Hessenberg form.
  245: *>          On exit, if COMPZ = 'I', the orthogonal matrix of
  246: *>          right Schur vectors of (H,T), and if COMPZ = 'V', the
  247: *>          orthogonal matrix of right Schur vectors of (A,B).
  248: *>          Not referenced if COMPZ = 'N'.
  249: *> \endverbatim
  250: *>
  251: *> \param[in] LDZ
  252: *> \verbatim
  253: *>          LDZ is INTEGER
  254: *>          The leading dimension of the array Z.  LDZ >= 1.
  255: *>          If COMPZ='V' or 'I', then LDZ >= N.
  256: *> \endverbatim
  257: *>
  258: *> \param[out] WORK
  259: *> \verbatim
  260: *>          WORK is DOUBLE PRECISION array, dimension (MAX(1,LWORK))
  261: *>          On exit, if INFO >= 0, WORK(1) returns the optimal LWORK.
  262: *> \endverbatim
  263: *>
  264: *> \param[in] LWORK
  265: *> \verbatim
  266: *>          LWORK is INTEGER
  267: *>          The dimension of the array WORK.  LWORK >= max(1,N).
  268: *>
  269: *>          If LWORK = -1, then a workspace query is assumed; the routine
  270: *>          only calculates the optimal size of the WORK array, returns
  271: *>          this value as the first entry of the WORK array, and no error
  272: *>          message related to LWORK is issued by XERBLA.
  273: *> \endverbatim
  274: *>
  275: *> \param[in] REC
  276: *> \verbatim
  277: *>          REC is INTEGER
  278: *>             REC indicates the current recursion level. Should be set
  279: *>             to 0 on first call.
  280: *> \endverbatim
  281: *>
  282: *> \param[out] INFO
  283: *> \verbatim
  284: *>          INFO is INTEGER
  285: *>          = 0: successful exit
  286: *>          < 0: if INFO = -i, the i-th argument had an illegal value
  287: *>          = 1,...,N: the QZ iteration did not converge.  (A,B) is not
  288: *>                     in Schur form, but ALPHAR(i), ALPHAI(i), and
  289: *>                     BETA(i), i=INFO+1,...,N should be correct.
  290: *> \endverbatim
  291: *
  292: *  Authors:
  293: *  ========
  294: *
  295: *> \author Thijs Steel, KU Leuven
  296: *
  297: *> \date May 2020
  298: *
  299: *> \ingroup doubleGEcomputational
  300: *>
  301: *  =====================================================================
  302:       RECURSIVE SUBROUTINE DLAQZ0( WANTS, WANTQ, WANTZ, N, ILO, IHI, A,
  303:      $                             LDA, B, LDB, ALPHAR, ALPHAI, BETA,
  304:      $                             Q, LDQ, Z, LDZ, WORK, LWORK, REC,
  305:      $                             INFO )
  306:       IMPLICIT NONE
  307: 
  308: *     Arguments
  309:       CHARACTER, INTENT( IN ) :: WANTS, WANTQ, WANTZ
  310:       INTEGER, INTENT( IN ) :: N, ILO, IHI, LDA, LDB, LDQ, LDZ, LWORK,
  311:      $         REC
  312: 
  313:       INTEGER, INTENT( OUT ) :: INFO
  314: 
  315:       DOUBLE PRECISION, INTENT( INOUT ) :: A( LDA, * ), B( LDB, * ),
  316:      $                  Q( LDQ, * ), Z( LDZ, * ), ALPHAR( * ),
  317:      $                  ALPHAI( * ), BETA( * ), WORK( * )
  318: 
  319: *     Parameters
  320:       DOUBLE PRECISION :: ZERO, ONE, HALF
  321:       PARAMETER( ZERO = 0.0D0, ONE = 1.0D0, HALF = 0.5D0 )
  322: 
  323: *     Local scalars
  324:       DOUBLE PRECISION :: SMLNUM, ULP, ESHIFT, SAFMIN, SAFMAX, C1, S1,
  325:      $                    TEMP, SWAP, BNORM, BTOL
  326:       INTEGER :: ISTART, ISTOP, IITER, MAXIT, ISTART2, K, LD, NSHIFTS,
  327:      $           NBLOCK, NW, NMIN, NIBBLE, N_UNDEFLATED, N_DEFLATED,
  328:      $           NS, SWEEP_INFO, SHIFTPOS, LWORKREQ, K2, ISTARTM,
  329:      $           ISTOPM, IWANTS, IWANTQ, IWANTZ, NORM_INFO, AED_INFO,
  330:      $           NWR, NBR, NSR, ITEMP1, ITEMP2, RCOST, I
  331:       LOGICAL :: ILSCHUR, ILQ, ILZ
  332:       CHARACTER :: JBCMPZ*3
  333: 
  334: *     External Functions
  335:       EXTERNAL :: XERBLA, DHGEQZ, DLASET, DLAQZ3, DLAQZ4, DLABAD,
  336:      $            DLARTG, DROT
  337:       DOUBLE PRECISION, EXTERNAL :: DLAMCH, DLANHS
  338:       LOGICAL, EXTERNAL :: LSAME
  339:       INTEGER, EXTERNAL :: ILAENV
  340: 
  341: *
  342: *     Decode wantS,wantQ,wantZ
  343: *      
  344:       IF( LSAME( WANTS, 'E' ) ) THEN
  345:          ILSCHUR = .FALSE.
  346:          IWANTS = 1
  347:       ELSE IF( LSAME( WANTS, 'S' ) ) THEN
  348:          ILSCHUR = .TRUE.
  349:          IWANTS = 2
  350:       ELSE
  351:          IWANTS = 0
  352:       END IF
  353: 
  354:       IF( LSAME( WANTQ, 'N' ) ) THEN
  355:          ILQ = .FALSE.
  356:          IWANTQ = 1
  357:       ELSE IF( LSAME( WANTQ, 'V' ) ) THEN
  358:          ILQ = .TRUE.
  359:          IWANTQ = 2
  360:       ELSE IF( LSAME( WANTQ, 'I' ) ) THEN
  361:          ILQ = .TRUE.
  362:          IWANTQ = 3
  363:       ELSE
  364:          IWANTQ = 0
  365:       END IF
  366: 
  367:       IF( LSAME( WANTZ, 'N' ) ) THEN
  368:          ILZ = .FALSE.
  369:          IWANTZ = 1
  370:       ELSE IF( LSAME( WANTZ, 'V' ) ) THEN
  371:          ILZ = .TRUE.
  372:          IWANTZ = 2
  373:       ELSE IF( LSAME( WANTZ, 'I' ) ) THEN
  374:          ILZ = .TRUE.
  375:          IWANTZ = 3
  376:       ELSE
  377:          IWANTZ = 0
  378:       END IF
  379: *
  380: *     Check Argument Values
  381: *
  382:       INFO = 0
  383:       IF( IWANTS.EQ.0 ) THEN
  384:          INFO = -1
  385:       ELSE IF( IWANTQ.EQ.0 ) THEN
  386:          INFO = -2
  387:       ELSE IF( IWANTZ.EQ.0 ) THEN
  388:          INFO = -3
  389:       ELSE IF( N.LT.0 ) THEN
  390:          INFO = -4
  391:       ELSE IF( ILO.LT.1 ) THEN
  392:          INFO = -5
  393:       ELSE IF( IHI.GT.N .OR. IHI.LT.ILO-1 ) THEN
  394:          INFO = -6
  395:       ELSE IF( LDA.LT.N ) THEN
  396:          INFO = -8
  397:       ELSE IF( LDB.LT.N ) THEN
  398:          INFO = -10
  399:       ELSE IF( LDQ.LT.1 .OR. ( ILQ .AND. LDQ.LT.N ) ) THEN
  400:          INFO = -15
  401:       ELSE IF( LDZ.LT.1 .OR. ( ILZ .AND. LDZ.LT.N ) ) THEN
  402:          INFO = -17
  403:       END IF
  404:       IF( INFO.NE.0 ) THEN
  405:          CALL XERBLA( 'DLAQZ0', -INFO )
  406:          RETURN
  407:       END IF
  408:    
  409: *
  410: *     Quick return if possible
  411: *
  412:       IF( N.LE.0 ) THEN
  413:          WORK( 1 ) = DBLE( 1 )
  414:          RETURN
  415:       END IF
  416: 
  417: *
  418: *     Get the parameters
  419: *
  420:       JBCMPZ( 1:1 ) = WANTS
  421:       JBCMPZ( 2:2 ) = WANTQ
  422:       JBCMPZ( 3:3 ) = WANTZ
  423: 
  424:       NMIN = ILAENV( 12, 'DLAQZ0', JBCMPZ, N, ILO, IHI, LWORK )
  425: 
  426:       NWR = ILAENV( 13, 'DLAQZ0', JBCMPZ, N, ILO, IHI, LWORK )
  427:       NWR = MAX( 2, NWR )
  428:       NWR = MIN( IHI-ILO+1, ( N-1 ) / 3, NWR )
  429: 
  430:       NIBBLE = ILAENV( 14, 'DLAQZ0', JBCMPZ, N, ILO, IHI, LWORK )
  431:       
  432:       NSR = ILAENV( 15, 'DLAQZ0', JBCMPZ, N, ILO, IHI, LWORK )
  433:       NSR = MIN( NSR, ( N+6 ) / 9, IHI-ILO )
  434:       NSR = MAX( 2, NSR-MOD( NSR, 2 ) )
  435: 
  436:       RCOST = ILAENV( 17, 'DLAQZ0', JBCMPZ, N, ILO, IHI, LWORK )
  437:       ITEMP1 = INT( NSR/SQRT( 1+2*NSR/( DBLE( RCOST )/100*N ) ) )
  438:       ITEMP1 = ( ( ITEMP1-1 )/4 )*4+4
  439:       NBR = NSR+ITEMP1
  440: 
  441:       IF( N .LT. NMIN .OR. REC .GE. 2 ) THEN
  442:          CALL DHGEQZ( WANTS, WANTQ, WANTZ, N, ILO, IHI, A, LDA, B, LDB,
  443:      $                ALPHAR, ALPHAI, BETA, Q, LDQ, Z, LDZ, WORK,
  444:      $                LWORK, INFO )
  445:          RETURN
  446:       END IF
  447: 
  448: *
  449: *     Find out required workspace
  450: *
  451: 
  452: *     Workspace query to dlaqz3
  453:       NW = MAX( NWR, NMIN )
  454:       CALL DLAQZ3( ILSCHUR, ILQ, ILZ, N, ILO, IHI, NW, A, LDA, B, LDB,
  455:      $             Q, LDQ, Z, LDZ, N_UNDEFLATED, N_DEFLATED, ALPHAR,
  456:      $             ALPHAI, BETA, WORK, NW, WORK, NW, WORK, -1, REC,
  457:      $             AED_INFO )
  458:       ITEMP1 = INT( WORK( 1 ) )
  459: *     Workspace query to dlaqz4
  460:       CALL DLAQZ4( ILSCHUR, ILQ, ILZ, N, ILO, IHI, NSR, NBR, ALPHAR,
  461:      $             ALPHAI, BETA, A, LDA, B, LDB, Q, LDQ, Z, LDZ, WORK,
  462:      $             NBR, WORK, NBR, WORK, -1, SWEEP_INFO )
  463:       ITEMP2 = INT( WORK( 1 ) )
  464: 
  465:       LWORKREQ = MAX( ITEMP1+2*NW**2, ITEMP2+2*NBR**2 )
  466:       IF ( LWORK .EQ.-1 ) THEN
  467:          WORK( 1 ) = DBLE( LWORKREQ )
  468:          RETURN
  469:       ELSE IF ( LWORK .LT. LWORKREQ ) THEN
  470:          INFO = -19
  471:       END IF
  472:       IF( INFO.NE.0 ) THEN
  473:          CALL XERBLA( 'DLAQZ0', INFO )
  474:          RETURN
  475:       END IF
  476: *
  477: *     Initialize Q and Z
  478: *
  479:       IF( IWANTQ.EQ.3 ) CALL DLASET( 'FULL', N, N, ZERO, ONE, Q, LDQ )
  480:       IF( IWANTZ.EQ.3 ) CALL DLASET( 'FULL', N, N, ZERO, ONE, Z, LDZ )
  481: 
  482: *     Get machine constants
  483:       SAFMIN = DLAMCH( 'SAFE MINIMUM' )
  484:       SAFMAX = ONE/SAFMIN
  485:       CALL DLABAD( SAFMIN, SAFMAX )
  486:       ULP = DLAMCH( 'PRECISION' )
  487:       SMLNUM = SAFMIN*( DBLE( N )/ULP )
  488: 
  489:       BNORM = DLANHS( 'F', IHI-ILO+1, B( ILO, ILO ), LDB, WORK )
  490:       BTOL = MAX( SAFMIN, ULP*BNORM )
  491: 
  492:       ISTART = ILO
  493:       ISTOP = IHI
  494:       MAXIT = 3*( IHI-ILO+1 )
  495:       LD = 0
  496:  
  497:       DO IITER = 1, MAXIT
  498:          IF( IITER .GE. MAXIT ) THEN
  499:             INFO = ISTOP+1
  500:             GOTO 80
  501:          END IF
  502:          IF ( ISTART+1 .GE. ISTOP ) THEN
  503:             ISTOP = ISTART
  504:             EXIT
  505:          END IF
  506: 
  507: *        Check deflations at the end
  508:          IF ( ABS( A( ISTOP-1, ISTOP-2 ) ) .LE. MAX( SMLNUM,
  509:      $      ULP*( ABS( A( ISTOP-1, ISTOP-1 ) )+ABS( A( ISTOP-2,
  510:      $      ISTOP-2 ) ) ) ) ) THEN
  511:             A( ISTOP-1, ISTOP-2 ) = ZERO
  512:             ISTOP = ISTOP-2
  513:             LD = 0
  514:             ESHIFT = ZERO
  515:          ELSE IF ( ABS( A( ISTOP, ISTOP-1 ) ) .LE. MAX( SMLNUM,
  516:      $      ULP*( ABS( A( ISTOP, ISTOP ) )+ABS( A( ISTOP-1,
  517:      $      ISTOP-1 ) ) ) ) ) THEN
  518:             A( ISTOP, ISTOP-1 ) = ZERO
  519:             ISTOP = ISTOP-1
  520:             LD = 0
  521:             ESHIFT = ZERO
  522:          END IF
  523: *        Check deflations at the start
  524:          IF ( ABS( A( ISTART+2, ISTART+1 ) ) .LE. MAX( SMLNUM,
  525:      $      ULP*( ABS( A( ISTART+1, ISTART+1 ) )+ABS( A( ISTART+2,
  526:      $      ISTART+2 ) ) ) ) ) THEN
  527:             A( ISTART+2, ISTART+1 ) = ZERO
  528:             ISTART = ISTART+2
  529:             LD = 0
  530:             ESHIFT = ZERO
  531:          ELSE IF ( ABS( A( ISTART+1, ISTART ) ) .LE. MAX( SMLNUM,
  532:      $      ULP*( ABS( A( ISTART, ISTART ) )+ABS( A( ISTART+1,
  533:      $      ISTART+1 ) ) ) ) ) THEN
  534:             A( ISTART+1, ISTART ) = ZERO
  535:             ISTART = ISTART+1
  536:             LD = 0
  537:             ESHIFT = ZERO
  538:          END IF
  539: 
  540:          IF ( ISTART+1 .GE. ISTOP ) THEN
  541:             EXIT
  542:          END IF
  543: 
  544: *        Check interior deflations
  545:          ISTART2 = ISTART
  546:          DO K = ISTOP, ISTART+1, -1
  547:             IF ( ABS( A( K, K-1 ) ) .LE. MAX( SMLNUM, ULP*( ABS( A( K,
  548:      $         K ) )+ABS( A( K-1, K-1 ) ) ) ) ) THEN
  549:                A( K, K-1 ) = ZERO
  550:                ISTART2 = K
  551:                EXIT
  552:             END IF
  553:          END DO
  554: 
  555: *        Get range to apply rotations to
  556:          IF ( ILSCHUR ) THEN
  557:             ISTARTM = 1
  558:             ISTOPM = N
  559:          ELSE
  560:             ISTARTM = ISTART2
  561:             ISTOPM = ISTOP
  562:          END IF
  563: 
  564: *        Check infinite eigenvalues, this is done without blocking so might
  565: *        slow down the method when many infinite eigenvalues are present
  566:          K = ISTOP
  567:          DO WHILE ( K.GE.ISTART2 )
  568: 
  569:             IF( ABS( B( K, K ) ) .LT. BTOL ) THEN
  570: *              A diagonal element of B is negligable, move it
  571: *              to the top and deflate it
  572:                
  573:                DO K2 = K, ISTART2+1, -1
  574:                   CALL DLARTG( B( K2-1, K2 ), B( K2-1, K2-1 ), C1, S1,
  575:      $                         TEMP )
  576:                   B( K2-1, K2 ) = TEMP
  577:                   B( K2-1, K2-1 ) = ZERO
  578: 
  579:                   CALL DROT( K2-2-ISTARTM+1, B( ISTARTM, K2 ), 1,
  580:      $                       B( ISTARTM, K2-1 ), 1, C1, S1 )
  581:                   CALL DROT( MIN( K2+1, ISTOP )-ISTARTM+1, A( ISTARTM,
  582:      $                       K2 ), 1, A( ISTARTM, K2-1 ), 1, C1, S1 )
  583:                   IF ( ILZ ) THEN
  584:                      CALL DROT( N, Z( 1, K2 ), 1, Z( 1, K2-1 ), 1, C1,
  585:      $                          S1 )
  586:                   END IF
  587: 
  588:                   IF( K2.LT.ISTOP ) THEN
  589:                      CALL DLARTG( A( K2, K2-1 ), A( K2+1, K2-1 ), C1,
  590:      $                            S1, TEMP )
  591:                      A( K2, K2-1 ) = TEMP
  592:                      A( K2+1, K2-1 ) = ZERO
  593: 
  594:                      CALL DROT( ISTOPM-K2+1, A( K2, K2 ), LDA, A( K2+1,
  595:      $                          K2 ), LDA, C1, S1 )
  596:                      CALL DROT( ISTOPM-K2+1, B( K2, K2 ), LDB, B( K2+1,
  597:      $                          K2 ), LDB, C1, S1 )
  598:                      IF( ILQ ) THEN
  599:                         CALL DROT( N, Q( 1, K2 ), 1, Q( 1, K2+1 ), 1,
  600:      $                             C1, S1 )
  601:                      END IF
  602:                   END IF
  603: 
  604:                END DO
  605: 
  606:                IF( ISTART2.LT.ISTOP )THEN
  607:                   CALL DLARTG( A( ISTART2, ISTART2 ), A( ISTART2+1,
  608:      $                         ISTART2 ), C1, S1, TEMP )
  609:                   A( ISTART2, ISTART2 ) = TEMP
  610:                   A( ISTART2+1, ISTART2 ) = ZERO
  611: 
  612:                   CALL DROT( ISTOPM-( ISTART2+1 )+1, A( ISTART2,
  613:      $                       ISTART2+1 ), LDA, A( ISTART2+1,
  614:      $                       ISTART2+1 ), LDA, C1, S1 )
  615:                   CALL DROT( ISTOPM-( ISTART2+1 )+1, B( ISTART2,
  616:      $                       ISTART2+1 ), LDB, B( ISTART2+1,
  617:      $                       ISTART2+1 ), LDB, C1, S1 )
  618:                   IF( ILQ ) THEN
  619:                      CALL DROT( N, Q( 1, ISTART2 ), 1, Q( 1,
  620:      $                          ISTART2+1 ), 1, C1, S1 )
  621:                   END IF
  622:                END IF
  623: 
  624:                ISTART2 = ISTART2+1
  625:    
  626:             END IF
  627:             K = K-1
  628:          END DO
  629: 
  630: *        istart2 now points to the top of the bottom right
  631: *        unreduced Hessenberg block
  632:          IF ( ISTART2 .GE. ISTOP ) THEN
  633:             ISTOP = ISTART2-1
  634:             LD = 0
  635:             ESHIFT = ZERO
  636:             CYCLE
  637:          END IF
  638: 
  639:          NW = NWR
  640:          NSHIFTS = NSR
  641:          NBLOCK = NBR
  642: 
  643:          IF ( ISTOP-ISTART2+1 .LT. NMIN ) THEN
  644: *           Setting nw to the size of the subblock will make AED deflate
  645: *           all the eigenvalues. This is slightly more efficient than just
  646: *           using DHGEQZ because the off diagonal part gets updated via BLAS.
  647:             IF ( ISTOP-ISTART+1 .LT. NMIN ) THEN
  648:                NW = ISTOP-ISTART+1
  649:                ISTART2 = ISTART
  650:             ELSE
  651:                NW = ISTOP-ISTART2+1
  652:             END IF
  653:          END IF
  654: 
  655: *
  656: *        Time for AED
  657: *
  658:          CALL DLAQZ3( ILSCHUR, ILQ, ILZ, N, ISTART2, ISTOP, NW, A, LDA,
  659:      $                B, LDB, Q, LDQ, Z, LDZ, N_UNDEFLATED, N_DEFLATED,
  660:      $                ALPHAR, ALPHAI, BETA, WORK, NW, WORK( NW**2+1 ),
  661:      $                NW, WORK( 2*NW**2+1 ), LWORK-2*NW**2, REC,
  662:      $                AED_INFO )
  663: 
  664:          IF ( N_DEFLATED > 0 ) THEN
  665:             ISTOP = ISTOP-N_DEFLATED
  666:             LD = 0
  667:             ESHIFT = ZERO
  668:          END IF
  669: 
  670:          IF ( 100*N_DEFLATED > NIBBLE*( N_DEFLATED+N_UNDEFLATED ) .OR.
  671:      $      ISTOP-ISTART2+1 .LT. NMIN ) THEN
  672: *           AED has uncovered many eigenvalues. Skip a QZ sweep and run
  673: *           AED again.
  674:             CYCLE
  675:          END IF
  676: 
  677:          LD = LD+1
  678: 
  679:          NS = MIN( NSHIFTS, ISTOP-ISTART2 )
  680:          NS = MIN( NS, N_UNDEFLATED )
  681:          SHIFTPOS = ISTOP-N_DEFLATED-N_UNDEFLATED+1
  682: *
  683: *        Shuffle shifts to put double shifts in front
  684: *        This ensures that we don't split up a double shift
  685: *
  686:          DO I = SHIFTPOS, SHIFTPOS+N_UNDEFLATED-1, 2
  687:             IF( ALPHAI( I ).NE.-ALPHAI( I+1 ) ) THEN
  688: *
  689:                SWAP = ALPHAR( I )
  690:                ALPHAR( I ) = ALPHAR( I+1 )
  691:                ALPHAR( I+1 ) = ALPHAR( I+2 )
  692:                ALPHAR( I+2 ) = SWAP
  693: 
  694:                SWAP = ALPHAI( I )
  695:                ALPHAI( I ) = ALPHAI( I+1 )
  696:                ALPHAI( I+1 ) = ALPHAI( I+2 )
  697:                ALPHAI( I+2 ) = SWAP
  698:                
  699:                SWAP = BETA( I )
  700:                BETA( I ) = BETA( I+1 )
  701:                BETA( I+1 ) = BETA( I+2 )
  702:                BETA( I+2 ) = SWAP
  703:             END IF
  704:          END DO
  705: 
  706:          IF ( MOD( LD, 6 ) .EQ. 0 ) THEN
  707:   708: *           Exceptional shift.  Chosen for no particularly good reason.
  709: *
  710:             IF( ( DBLE( MAXIT )*SAFMIN )*ABS( A( ISTOP,
  711:      $         ISTOP-1 ) ).LT.ABS( A( ISTOP-1, ISTOP-1 ) ) ) THEN
  712:                ESHIFT = A( ISTOP, ISTOP-1 )/B( ISTOP-1, ISTOP-1 )
  713:             ELSE
  714:                ESHIFT = ESHIFT+ONE/( SAFMIN*DBLE( MAXIT ) )
  715:             END IF
  716:             ALPHAR( SHIFTPOS ) = ONE
  717:             ALPHAR( SHIFTPOS+1 ) = ZERO
  718:             ALPHAI( SHIFTPOS ) = ZERO
  719:             ALPHAI( SHIFTPOS+1 ) = ZERO
  720:             BETA( SHIFTPOS ) = ESHIFT
  721:             BETA( SHIFTPOS+1 ) = ESHIFT
  722:             NS = 2
  723:          END IF
  724: 
  725: *
  726: *        Time for a QZ sweep
  727: *
  728:          CALL DLAQZ4( ILSCHUR, ILQ, ILZ, N, ISTART2, ISTOP, NS, NBLOCK,
  729:      $                ALPHAR( SHIFTPOS ), ALPHAI( SHIFTPOS ),
  730:      $                BETA( SHIFTPOS ), A, LDA, B, LDB, Q, LDQ, Z, LDZ,
  731:      $                WORK, NBLOCK, WORK( NBLOCK**2+1 ), NBLOCK,
  732:      $                WORK( 2*NBLOCK**2+1 ), LWORK-2*NBLOCK**2,
  733:      $                SWEEP_INFO )
  734: 
  735:       END DO
  736: 
  737: *
  738: *     Call DHGEQZ to normalize the eigenvalue blocks and set the eigenvalues
  739: *     If all the eigenvalues have been found, DHGEQZ will not do any iterations
  740: *     and only normalize the blocks. In case of a rare convergence failure,
  741: *     the single shift might perform better.
  742: *
  743:    80 CALL DHGEQZ( WANTS, WANTQ, WANTZ, N, ILO, IHI, A, LDA, B, LDB,
  744:      $             ALPHAR, ALPHAI, BETA, Q, LDQ, Z, LDZ, WORK, LWORK,
  745:      $             NORM_INFO )
  746:       
  747:       INFO = NORM_INFO
  748: 
  749:       END SUBROUTINE

CVSweb interface <joel.bertrand@systella.fr>