Annotation of rpl/lapack/lapack/dgsvj0.f, revision 1.7

1.7     ! bertrand    1: *> \brief \b DGSVJ0
        !             2: *
        !             3: *  =========== DOCUMENTATION ===========
        !             4: *
        !             5: * Online html documentation available at 
        !             6: *            http://www.netlib.org/lapack/explore-html/ 
        !             7: *
        !             8: *> \htmlonly
        !             9: *> Download DGSVJ0 + dependencies 
        !            10: *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dgsvj0.f"> 
        !            11: *> [TGZ]</a> 
        !            12: *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dgsvj0.f"> 
        !            13: *> [ZIP]</a> 
        !            14: *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dgsvj0.f"> 
        !            15: *> [TXT]</a>
        !            16: *> \endhtmlonly 
        !            17: *
        !            18: *  Definition:
        !            19: *  ===========
        !            20: *
        !            21: *       SUBROUTINE DGSVJ0( JOBV, M, N, A, LDA, D, SVA, MV, V, LDV, EPS,
        !            22: *                          SFMIN, TOL, NSWEEP, WORK, LWORK, INFO )
        !            23: * 
        !            24: *       .. Scalar Arguments ..
        !            25: *       INTEGER            INFO, LDA, LDV, LWORK, M, MV, N, NSWEEP
        !            26: *       DOUBLE PRECISION   EPS, SFMIN, TOL
        !            27: *       CHARACTER*1        JOBV
        !            28: *       ..
        !            29: *       .. Array Arguments ..
        !            30: *       DOUBLE PRECISION   A( LDA, * ), SVA( N ), D( N ), V( LDV, * ),
        !            31: *      $                   WORK( LWORK )
        !            32: *       ..
        !            33: *  
        !            34: *
        !            35: *> \par Purpose:
        !            36: *  =============
        !            37: *>
        !            38: *> \verbatim
        !            39: *>
        !            40: *> DGSVJ0 is called from DGESVJ as a pre-processor and that is its main
        !            41: *> purpose. It applies Jacobi rotations in the same way as DGESVJ does, but
        !            42: *> it does not check convergence (stopping criterion). Few tuning
        !            43: *> parameters (marked by [TP]) are available for the implementer.
        !            44: *> \endverbatim
        !            45: *
        !            46: *  Arguments:
        !            47: *  ==========
        !            48: *
        !            49: *> \param[in] JOBV
        !            50: *> \verbatim
        !            51: *>          JOBV is CHARACTER*1
        !            52: *>          Specifies whether the output from this procedure is used
        !            53: *>          to compute the matrix V:
        !            54: *>          = 'V': the product of the Jacobi rotations is accumulated
        !            55: *>                 by postmulyiplying the N-by-N array V.
        !            56: *>                (See the description of V.)
        !            57: *>          = 'A': the product of the Jacobi rotations is accumulated
        !            58: *>                 by postmulyiplying the MV-by-N array V.
        !            59: *>                (See the descriptions of MV and V.)
        !            60: *>          = 'N': the Jacobi rotations are not accumulated.
        !            61: *> \endverbatim
        !            62: *>
        !            63: *> \param[in] M
        !            64: *> \verbatim
        !            65: *>          M is INTEGER
        !            66: *>          The number of rows of the input matrix A.  M >= 0.
        !            67: *> \endverbatim
        !            68: *>
        !            69: *> \param[in] N
        !            70: *> \verbatim
        !            71: *>          N is INTEGER
        !            72: *>          The number of columns of the input matrix A.
        !            73: *>          M >= N >= 0.
        !            74: *> \endverbatim
        !            75: *>
        !            76: *> \param[in,out] A
        !            77: *> \verbatim
        !            78: *>          A is DOUBLE PRECISION array, dimension (LDA,N)
        !            79: *>          On entry, M-by-N matrix A, such that A*diag(D) represents
        !            80: *>          the input matrix.
        !            81: *>          On exit,
        !            82: *>          A_onexit * D_onexit represents the input matrix A*diag(D)
        !            83: *>          post-multiplied by a sequence of Jacobi rotations, where the
        !            84: *>          rotation threshold and the total number of sweeps are given in
        !            85: *>          TOL and NSWEEP, respectively.
        !            86: *>          (See the descriptions of D, TOL and NSWEEP.)
        !            87: *> \endverbatim
        !            88: *>
        !            89: *> \param[in] LDA
        !            90: *> \verbatim
        !            91: *>          LDA is INTEGER
        !            92: *>          The leading dimension of the array A.  LDA >= max(1,M).
        !            93: *> \endverbatim
        !            94: *>
        !            95: *> \param[in,out] D
        !            96: *> \verbatim
        !            97: *>          D is DOUBLE PRECISION array, dimension (N)
        !            98: *>          The array D accumulates the scaling factors from the fast scaled
        !            99: *>          Jacobi rotations.
        !           100: *>          On entry, A*diag(D) represents the input matrix.
        !           101: *>          On exit, A_onexit*diag(D_onexit) represents the input matrix
        !           102: *>          post-multiplied by a sequence of Jacobi rotations, where the
        !           103: *>          rotation threshold and the total number of sweeps are given in
        !           104: *>          TOL and NSWEEP, respectively.
        !           105: *>          (See the descriptions of A, TOL and NSWEEP.)
        !           106: *> \endverbatim
        !           107: *>
        !           108: *> \param[in,out] SVA
        !           109: *> \verbatim
        !           110: *>          SVA is DOUBLE PRECISION array, dimension (N)
        !           111: *>          On entry, SVA contains the Euclidean norms of the columns of
        !           112: *>          the matrix A*diag(D).
        !           113: *>          On exit, SVA contains the Euclidean norms of the columns of
        !           114: *>          the matrix onexit*diag(D_onexit).
        !           115: *> \endverbatim
        !           116: *>
        !           117: *> \param[in] MV
        !           118: *> \verbatim
        !           119: *>          MV is INTEGER
        !           120: *>          If JOBV .EQ. 'A', then MV rows of V are post-multipled by a
        !           121: *>                           sequence of Jacobi rotations.
        !           122: *>          If JOBV = 'N',   then MV is not referenced.
        !           123: *> \endverbatim
        !           124: *>
        !           125: *> \param[in,out] V
        !           126: *> \verbatim
        !           127: *>          V is DOUBLE PRECISION array, dimension (LDV,N)
        !           128: *>          If JOBV .EQ. 'V' then N rows of V are post-multipled by a
        !           129: *>                           sequence of Jacobi rotations.
        !           130: *>          If JOBV .EQ. 'A' then MV rows of V are post-multipled by a
        !           131: *>                           sequence of Jacobi rotations.
        !           132: *>          If JOBV = 'N',   then V is not referenced.
        !           133: *> \endverbatim
        !           134: *>
        !           135: *> \param[in] LDV
        !           136: *> \verbatim
        !           137: *>          LDV is INTEGER
        !           138: *>          The leading dimension of the array V,  LDV >= 1.
        !           139: *>          If JOBV = 'V', LDV .GE. N.
        !           140: *>          If JOBV = 'A', LDV .GE. MV.
        !           141: *> \endverbatim
        !           142: *>
        !           143: *> \param[in] EPS
        !           144: *> \verbatim
        !           145: *>          EPS is DOUBLE PRECISION
        !           146: *>          EPS = DLAMCH('Epsilon')
        !           147: *> \endverbatim
        !           148: *>
        !           149: *> \param[in] SFMIN
        !           150: *> \verbatim
        !           151: *>          SFMIN is DOUBLE PRECISION
        !           152: *>          SFMIN = DLAMCH('Safe Minimum')
        !           153: *> \endverbatim
        !           154: *>
        !           155: *> \param[in] TOL
        !           156: *> \verbatim
        !           157: *>          TOL is DOUBLE PRECISION
        !           158: *>          TOL is the threshold for Jacobi rotations. For a pair
        !           159: *>          A(:,p), A(:,q) of pivot columns, the Jacobi rotation is
        !           160: *>          applied only if DABS(COS(angle(A(:,p),A(:,q)))) .GT. TOL.
        !           161: *> \endverbatim
        !           162: *>
        !           163: *> \param[in] NSWEEP
        !           164: *> \verbatim
        !           165: *>          NSWEEP is INTEGER
        !           166: *>          NSWEEP is the number of sweeps of Jacobi rotations to be
        !           167: *>          performed.
        !           168: *> \endverbatim
        !           169: *>
        !           170: *> \param[out] WORK
        !           171: *> \verbatim
        !           172: *>          WORK is DOUBLE PRECISION array, dimension (LWORK)
        !           173: *> \endverbatim
        !           174: *>
        !           175: *> \param[in] LWORK
        !           176: *> \verbatim
        !           177: *>          LWORK is INTEGER
        !           178: *>          LWORK is the dimension of WORK. LWORK .GE. M.
        !           179: *> \endverbatim
        !           180: *>
        !           181: *> \param[out] INFO
        !           182: *> \verbatim
        !           183: *>          INFO is INTEGER
        !           184: *>          = 0 : successful exit.
        !           185: *>          < 0 : if INFO = -i, then the i-th argument had an illegal value
        !           186: *> \endverbatim
        !           187: *
        !           188: *  Authors:
        !           189: *  ========
        !           190: *
        !           191: *> \author Univ. of Tennessee 
        !           192: *> \author Univ. of California Berkeley 
        !           193: *> \author Univ. of Colorado Denver 
        !           194: *> \author NAG Ltd. 
        !           195: *
        !           196: *> \date November 2011
        !           197: *
        !           198: *> \ingroup doubleOTHERcomputational
        !           199: *
        !           200: *> \par Further Details:
        !           201: *  =====================
        !           202: *>
        !           203: *> DGSVJ0 is used just to enable DGESVJ to call a simplified version of
        !           204: *> itself to work on a submatrix of the original matrix.
        !           205: *>
        !           206: *> \par Contributors:
        !           207: *  ==================
        !           208: *>
        !           209: *> Zlatko Drmac (Zagreb, Croatia) and Kresimir Veselic (Hagen, Germany)
        !           210: *>
        !           211: *> \par Bugs, Examples and Comments:
        !           212: *  =================================
        !           213: *>
        !           214: *> Please report all bugs and send interesting test examples and comments to
        !           215: *> drmac@math.hr. Thank you.
        !           216: *
        !           217: *  =====================================================================
1.1       bertrand  218:       SUBROUTINE DGSVJ0( JOBV, M, N, A, LDA, D, SVA, MV, V, LDV, EPS,
1.6       bertrand  219:      $                   SFMIN, TOL, NSWEEP, WORK, LWORK, INFO )
1.1       bertrand  220: *
1.7     ! bertrand  221: *  -- LAPACK computational routine (version 3.4.0) --
1.1       bertrand  222: *  -- LAPACK is a software package provided by Univ. of Tennessee,    --
                    223: *  -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
1.7     ! bertrand  224: *     November 2011
1.1       bertrand  225: *
                    226: *     .. Scalar Arguments ..
                    227:       INTEGER            INFO, LDA, LDV, LWORK, M, MV, N, NSWEEP
                    228:       DOUBLE PRECISION   EPS, SFMIN, TOL
                    229:       CHARACTER*1        JOBV
                    230: *     ..
                    231: *     .. Array Arguments ..
                    232:       DOUBLE PRECISION   A( LDA, * ), SVA( N ), D( N ), V( LDV, * ),
1.6       bertrand  233:      $                   WORK( LWORK )
1.1       bertrand  234: *     ..
                    235: *
                    236: *  =====================================================================
                    237: *
                    238: *     .. Local Parameters ..
                    239:       DOUBLE PRECISION   ZERO, HALF, ONE, TWO
                    240:       PARAMETER          ( ZERO = 0.0D0, HALF = 0.5D0, ONE = 1.0D0,
1.6       bertrand  241:      $                   TWO = 2.0D0 )
1.1       bertrand  242: *     ..
                    243: *     .. Local Scalars ..
                    244:       DOUBLE PRECISION   AAPP, AAPP0, AAPQ, AAQQ, APOAQ, AQOAP, BIG,
1.6       bertrand  245:      $                   BIGTHETA, CS, MXAAPQ, MXSINJ, ROOTBIG, ROOTEPS,
                    246:      $                   ROOTSFMIN, ROOTTOL, SMALL, SN, T, TEMP1, THETA,
                    247:      $                   THSIGN
1.1       bertrand  248:       INTEGER            BLSKIP, EMPTSW, i, ibr, IERR, igl, IJBLSK, ir1,
1.6       bertrand  249:      $                   ISWROT, jbc, jgl, KBL, LKAHEAD, MVL, NBL,
                    250:      $                   NOTROT, p, PSKIPPED, q, ROWSKIP, SWBAND
1.1       bertrand  251:       LOGICAL            APPLV, ROTOK, RSVEC
                    252: *     ..
                    253: *     .. Local Arrays ..
                    254:       DOUBLE PRECISION   FASTR( 5 )
                    255: *     ..
                    256: *     .. Intrinsic Functions ..
                    257:       INTRINSIC          DABS, DMAX1, DBLE, MIN0, DSIGN, DSQRT
                    258: *     ..
                    259: *     .. External Functions ..
                    260:       DOUBLE PRECISION   DDOT, DNRM2
                    261:       INTEGER            IDAMAX
                    262:       LOGICAL            LSAME
                    263:       EXTERNAL           IDAMAX, LSAME, DDOT, DNRM2
                    264: *     ..
                    265: *     .. External Subroutines ..
                    266:       EXTERNAL           DAXPY, DCOPY, DLASCL, DLASSQ, DROTM, DSWAP
                    267: *     ..
                    268: *     .. Executable Statements ..
                    269: *
1.6       bertrand  270: *     Test the input parameters.
                    271: *
1.1       bertrand  272:       APPLV = LSAME( JOBV, 'A' )
                    273:       RSVEC = LSAME( JOBV, 'V' )
                    274:       IF( .NOT.( RSVEC .OR. APPLV .OR. LSAME( JOBV, 'N' ) ) ) THEN
                    275:          INFO = -1
                    276:       ELSE IF( M.LT.0 ) THEN
                    277:          INFO = -2
                    278:       ELSE IF( ( N.LT.0 ) .OR. ( N.GT.M ) ) THEN
                    279:          INFO = -3
                    280:       ELSE IF( LDA.LT.M ) THEN
                    281:          INFO = -5
1.4       bertrand  282:       ELSE IF( ( RSVEC.OR.APPLV ) .AND. ( MV.LT.0 ) ) THEN
1.1       bertrand  283:          INFO = -8
1.4       bertrand  284:       ELSE IF( ( RSVEC.AND.( LDV.LT.N ) ).OR. 
1.6       bertrand  285:      $         ( APPLV.AND.( LDV.LT.MV ) ) ) THEN
1.1       bertrand  286:          INFO = -10
                    287:       ELSE IF( TOL.LE.EPS ) THEN
                    288:          INFO = -13
                    289:       ELSE IF( NSWEEP.LT.0 ) THEN
                    290:          INFO = -14
                    291:       ELSE IF( LWORK.LT.M ) THEN
                    292:          INFO = -16
                    293:       ELSE
                    294:          INFO = 0
                    295:       END IF
                    296: *
                    297: *     #:(
                    298:       IF( INFO.NE.0 ) THEN
                    299:          CALL XERBLA( 'DGSVJ0', -INFO )
                    300:          RETURN
                    301:       END IF
                    302: *
                    303:       IF( RSVEC ) THEN
                    304:          MVL = N
                    305:       ELSE IF( APPLV ) THEN
                    306:          MVL = MV
                    307:       END IF
                    308:       RSVEC = RSVEC .OR. APPLV
                    309: 
                    310:       ROOTEPS = DSQRT( EPS )
                    311:       ROOTSFMIN = DSQRT( SFMIN )
                    312:       SMALL = SFMIN / EPS
                    313:       BIG = ONE / SFMIN
                    314:       ROOTBIG = ONE / ROOTSFMIN
                    315:       BIGTHETA = ONE / ROOTEPS
                    316:       ROOTTOL = DSQRT( TOL )
                    317: *
                    318: *     -#- Row-cyclic Jacobi SVD algorithm with column pivoting -#-
                    319: *
                    320:       EMPTSW = ( N*( N-1 ) ) / 2
                    321:       NOTROT = 0
                    322:       FASTR( 1 ) = ZERO
                    323: *
                    324: *     -#- Row-cyclic pivot strategy with de Rijk's pivoting -#-
                    325: *
                    326: 
                    327:       SWBAND = 0
                    328: *[TP] SWBAND is a tuning parameter. It is meaningful and effective
                    329: *     if SGESVJ is used as a computational routine in the preconditioned
                    330: *     Jacobi SVD algorithm SGESVJ. For sweeps i=1:SWBAND the procedure
                    331: *     ......
                    332: 
                    333:       KBL = MIN0( 8, N )
                    334: *[TP] KBL is a tuning parameter that defines the tile size in the
                    335: *     tiling of the p-q loops of pivot pairs. In general, an optimal
                    336: *     value of KBL depends on the matrix dimensions and on the
                    337: *     parameters of the computer's memory.
                    338: *
                    339:       NBL = N / KBL
                    340:       IF( ( NBL*KBL ).NE.N )NBL = NBL + 1
                    341: 
                    342:       BLSKIP = ( KBL**2 ) + 1
                    343: *[TP] BLKSKIP is a tuning parameter that depends on SWBAND and KBL.
                    344: 
                    345:       ROWSKIP = MIN0( 5, KBL )
                    346: *[TP] ROWSKIP is a tuning parameter.
                    347: 
                    348:       LKAHEAD = 1
                    349: *[TP] LKAHEAD is a tuning parameter.
                    350:       SWBAND = 0
                    351:       PSKIPPED = 0
                    352: *
                    353:       DO 1993 i = 1, NSWEEP
                    354: *     .. go go go ...
                    355: *
                    356:          MXAAPQ = ZERO
                    357:          MXSINJ = ZERO
                    358:          ISWROT = 0
                    359: *
                    360:          NOTROT = 0
                    361:          PSKIPPED = 0
                    362: *
                    363:          DO 2000 ibr = 1, NBL
                    364: 
                    365:             igl = ( ibr-1 )*KBL + 1
                    366: *
                    367:             DO 1002 ir1 = 0, MIN0( LKAHEAD, NBL-ibr )
                    368: *
                    369:                igl = igl + ir1*KBL
                    370: *
                    371:                DO 2001 p = igl, MIN0( igl+KBL-1, N-1 )
                    372: 
                    373: *     .. de Rijk's pivoting
                    374:                   q = IDAMAX( N-p+1, SVA( p ), 1 ) + p - 1
                    375:                   IF( p.NE.q ) THEN
                    376:                      CALL DSWAP( M, A( 1, p ), 1, A( 1, q ), 1 )
                    377:                      IF( RSVEC )CALL DSWAP( MVL, V( 1, p ), 1,
1.6       bertrand  378:      $                                      V( 1, q ), 1 )
1.1       bertrand  379:                      TEMP1 = SVA( p )
                    380:                      SVA( p ) = SVA( q )
                    381:                      SVA( q ) = TEMP1
                    382:                      TEMP1 = D( p )
                    383:                      D( p ) = D( q )
                    384:                      D( q ) = TEMP1
                    385:                   END IF
                    386: *
                    387:                   IF( ir1.EQ.0 ) THEN
                    388: *
                    389: *        Column norms are periodically updated by explicit
                    390: *        norm computation.
                    391: *        Caveat:
                    392: *        Some BLAS implementations compute DNRM2(M,A(1,p),1)
                    393: *        as DSQRT(DDOT(M,A(1,p),1,A(1,p),1)), which may result in
                    394: *        overflow for ||A(:,p)||_2 > DSQRT(overflow_threshold), and
                    395: *        undeflow for ||A(:,p)||_2 < DSQRT(underflow_threshold).
                    396: *        Hence, DNRM2 cannot be trusted, not even in the case when
                    397: *        the true norm is far from the under(over)flow boundaries.
                    398: *        If properly implemented DNRM2 is available, the IF-THEN-ELSE
                    399: *        below should read "AAPP = DNRM2( M, A(1,p), 1 ) * D(p)".
                    400: *
                    401:                      IF( ( SVA( p ).LT.ROOTBIG ) .AND.
1.6       bertrand  402:      $                   ( SVA( p ).GT.ROOTSFMIN ) ) THEN
1.1       bertrand  403:                         SVA( p ) = DNRM2( M, A( 1, p ), 1 )*D( p )
                    404:                      ELSE
                    405:                         TEMP1 = ZERO
1.4       bertrand  406:                         AAPP = ONE
1.1       bertrand  407:                         CALL DLASSQ( M, A( 1, p ), 1, TEMP1, AAPP )
                    408:                         SVA( p ) = TEMP1*DSQRT( AAPP )*D( p )
                    409:                      END IF
                    410:                      AAPP = SVA( p )
                    411:                   ELSE
                    412:                      AAPP = SVA( p )
                    413:                   END IF
                    414: 
                    415: *
                    416:                   IF( AAPP.GT.ZERO ) THEN
                    417: *
                    418:                      PSKIPPED = 0
                    419: *
                    420:                      DO 2002 q = p + 1, MIN0( igl+KBL-1, N )
                    421: *
                    422:                         AAQQ = SVA( q )
                    423: 
                    424:                         IF( AAQQ.GT.ZERO ) THEN
                    425: *
                    426:                            AAPP0 = AAPP
                    427:                            IF( AAQQ.GE.ONE ) THEN
                    428:                               ROTOK = ( SMALL*AAPP ).LE.AAQQ
                    429:                               IF( AAPP.LT.( BIG / AAQQ ) ) THEN
                    430:                                  AAPQ = ( DDOT( M, A( 1, p ), 1, A( 1,
1.6       bertrand  431:      $                                  q ), 1 )*D( p )*D( q ) / AAQQ )
                    432:      $                                  / AAPP
1.1       bertrand  433:                               ELSE
                    434:                                  CALL DCOPY( M, A( 1, p ), 1, WORK, 1 )
                    435:                                  CALL DLASCL( 'G', 0, 0, AAPP, D( p ),
1.6       bertrand  436:      $                                        M, 1, WORK, LDA, IERR )
1.1       bertrand  437:                                  AAPQ = DDOT( M, WORK, 1, A( 1, q ),
1.6       bertrand  438:      $                                  1 )*D( q ) / AAQQ
1.1       bertrand  439:                               END IF
                    440:                            ELSE
                    441:                               ROTOK = AAPP.LE.( AAQQ / SMALL )
                    442:                               IF( AAPP.GT.( SMALL / AAQQ ) ) THEN
                    443:                                  AAPQ = ( DDOT( M, A( 1, p ), 1, A( 1,
1.6       bertrand  444:      $                                  q ), 1 )*D( p )*D( q ) / AAQQ )
                    445:      $                                  / AAPP
1.1       bertrand  446:                               ELSE
                    447:                                  CALL DCOPY( M, A( 1, q ), 1, WORK, 1 )
                    448:                                  CALL DLASCL( 'G', 0, 0, AAQQ, D( q ),
1.6       bertrand  449:      $                                        M, 1, WORK, LDA, IERR )
1.1       bertrand  450:                                  AAPQ = DDOT( M, WORK, 1, A( 1, p ),
1.6       bertrand  451:      $                                  1 )*D( p ) / AAPP
1.1       bertrand  452:                               END IF
                    453:                            END IF
                    454: *
                    455:                            MXAAPQ = DMAX1( MXAAPQ, DABS( AAPQ ) )
                    456: *
                    457: *        TO rotate or NOT to rotate, THAT is the question ...
                    458: *
                    459:                            IF( DABS( AAPQ ).GT.TOL ) THEN
                    460: *
                    461: *           .. rotate
                    462: *           ROTATED = ROTATED + ONE
                    463: *
                    464:                               IF( ir1.EQ.0 ) THEN
                    465:                                  NOTROT = 0
                    466:                                  PSKIPPED = 0
                    467:                                  ISWROT = ISWROT + 1
                    468:                               END IF
                    469: *
                    470:                               IF( ROTOK ) THEN
                    471: *
                    472:                                  AQOAP = AAQQ / AAPP
                    473:                                  APOAQ = AAPP / AAQQ
1.6       bertrand  474:                                  THETA = -HALF*DABS( AQOAP-APOAQ )/AAPQ
1.1       bertrand  475: *
                    476:                                  IF( DABS( THETA ).GT.BIGTHETA ) THEN
                    477: *
                    478:                                     T = HALF / THETA
                    479:                                     FASTR( 3 ) = T*D( p ) / D( q )
                    480:                                     FASTR( 4 ) = -T*D( q ) / D( p )
                    481:                                     CALL DROTM( M, A( 1, p ), 1,
1.6       bertrand  482:      $                                          A( 1, q ), 1, FASTR )
1.1       bertrand  483:                                     IF( RSVEC )CALL DROTM( MVL,
1.6       bertrand  484:      $                                              V( 1, p ), 1,
                    485:      $                                              V( 1, q ), 1,
                    486:      $                                              FASTR )
1.1       bertrand  487:                                     SVA( q ) = AAQQ*DSQRT( DMAX1( ZERO,
1.6       bertrand  488:      $                                         ONE+T*APOAQ*AAPQ ) )
1.4       bertrand  489:                                     AAPP = AAPP*DSQRT( DMAX1( ZERO, 
1.6       bertrand  490:      $                                     ONE-T*AQOAP*AAPQ ) )
1.1       bertrand  491:                                     MXSINJ = DMAX1( MXSINJ, DABS( T ) )
                    492: *
                    493:                                  ELSE
                    494: *
                    495: *                 .. choose correct signum for THETA and rotate
                    496: *
                    497:                                     THSIGN = -DSIGN( ONE, AAPQ )
                    498:                                     T = ONE / ( THETA+THSIGN*
1.6       bertrand  499:      $                                  DSQRT( ONE+THETA*THETA ) )
1.1       bertrand  500:                                     CS = DSQRT( ONE / ( ONE+T*T ) )
                    501:                                     SN = T*CS
                    502: *
                    503:                                     MXSINJ = DMAX1( MXSINJ, DABS( SN ) )
                    504:                                     SVA( q ) = AAQQ*DSQRT( DMAX1( ZERO,
1.6       bertrand  505:      $                                         ONE+T*APOAQ*AAPQ ) )
1.1       bertrand  506:                                     AAPP = AAPP*DSQRT( DMAX1( ZERO,
1.6       bertrand  507:      $                                     ONE-T*AQOAP*AAPQ ) )
1.1       bertrand  508: *
                    509:                                     APOAQ = D( p ) / D( q )
                    510:                                     AQOAP = D( q ) / D( p )
                    511:                                     IF( D( p ).GE.ONE ) THEN
                    512:                                        IF( D( q ).GE.ONE ) THEN
                    513:                                           FASTR( 3 ) = T*APOAQ
                    514:                                           FASTR( 4 ) = -T*AQOAP
                    515:                                           D( p ) = D( p )*CS
                    516:                                           D( q ) = D( q )*CS
                    517:                                           CALL DROTM( M, A( 1, p ), 1,
1.6       bertrand  518:      $                                                A( 1, q ), 1,
                    519:      $                                                FASTR )
1.1       bertrand  520:                                           IF( RSVEC )CALL DROTM( MVL,
1.6       bertrand  521:      $                                        V( 1, p ), 1, V( 1, q ),
                    522:      $                                        1, FASTR )
1.1       bertrand  523:                                        ELSE
                    524:                                           CALL DAXPY( M, -T*AQOAP,
1.6       bertrand  525:      $                                                A( 1, q ), 1,
                    526:      $                                                A( 1, p ), 1 )
1.1       bertrand  527:                                           CALL DAXPY( M, CS*SN*APOAQ,
1.6       bertrand  528:      $                                                A( 1, p ), 1,
                    529:      $                                                A( 1, q ), 1 )
1.1       bertrand  530:                                           D( p ) = D( p )*CS
                    531:                                           D( q ) = D( q ) / CS
                    532:                                           IF( RSVEC ) THEN
                    533:                                              CALL DAXPY( MVL, -T*AQOAP,
1.6       bertrand  534:      $                                                   V( 1, q ), 1,
                    535:      $                                                   V( 1, p ), 1 )
1.1       bertrand  536:                                              CALL DAXPY( MVL,
1.6       bertrand  537:      $                                                   CS*SN*APOAQ,
                    538:      $                                                   V( 1, p ), 1,
                    539:      $                                                   V( 1, q ), 1 )
1.1       bertrand  540:                                           END IF
                    541:                                        END IF
                    542:                                     ELSE
                    543:                                        IF( D( q ).GE.ONE ) THEN
                    544:                                           CALL DAXPY( M, T*APOAQ,
1.6       bertrand  545:      $                                                A( 1, p ), 1,
                    546:      $                                                A( 1, q ), 1 )
1.1       bertrand  547:                                           CALL DAXPY( M, -CS*SN*AQOAP,
1.6       bertrand  548:      $                                                A( 1, q ), 1,
                    549:      $                                                A( 1, p ), 1 )
1.1       bertrand  550:                                           D( p ) = D( p ) / CS
                    551:                                           D( q ) = D( q )*CS
                    552:                                           IF( RSVEC ) THEN
                    553:                                              CALL DAXPY( MVL, T*APOAQ,
1.6       bertrand  554:      $                                                   V( 1, p ), 1,
                    555:      $                                                   V( 1, q ), 1 )
1.1       bertrand  556:                                              CALL DAXPY( MVL,
1.6       bertrand  557:      $                                                   -CS*SN*AQOAP,
                    558:      $                                                   V( 1, q ), 1,
                    559:      $                                                   V( 1, p ), 1 )
1.1       bertrand  560:                                           END IF
                    561:                                        ELSE
                    562:                                           IF( D( p ).GE.D( q ) ) THEN
                    563:                                              CALL DAXPY( M, -T*AQOAP,
1.6       bertrand  564:      $                                                   A( 1, q ), 1,
                    565:      $                                                   A( 1, p ), 1 )
1.1       bertrand  566:                                              CALL DAXPY( M, CS*SN*APOAQ,
1.6       bertrand  567:      $                                                   A( 1, p ), 1,
                    568:      $                                                   A( 1, q ), 1 )
1.1       bertrand  569:                                              D( p ) = D( p )*CS
                    570:                                              D( q ) = D( q ) / CS
                    571:                                              IF( RSVEC ) THEN
                    572:                                                 CALL DAXPY( MVL,
1.6       bertrand  573:      $                                               -T*AQOAP,
                    574:      $                                               V( 1, q ), 1,
                    575:      $                                               V( 1, p ), 1 )
1.1       bertrand  576:                                                 CALL DAXPY( MVL,
1.6       bertrand  577:      $                                               CS*SN*APOAQ,
                    578:      $                                               V( 1, p ), 1,
                    579:      $                                               V( 1, q ), 1 )
1.1       bertrand  580:                                              END IF
                    581:                                           ELSE
                    582:                                              CALL DAXPY( M, T*APOAQ,
1.6       bertrand  583:      $                                                   A( 1, p ), 1,
                    584:      $                                                   A( 1, q ), 1 )
1.1       bertrand  585:                                              CALL DAXPY( M,
1.6       bertrand  586:      $                                                   -CS*SN*AQOAP,
                    587:      $                                                   A( 1, q ), 1,
                    588:      $                                                   A( 1, p ), 1 )
1.1       bertrand  589:                                              D( p ) = D( p ) / CS
                    590:                                              D( q ) = D( q )*CS
                    591:                                              IF( RSVEC ) THEN
                    592:                                                 CALL DAXPY( MVL,
1.6       bertrand  593:      $                                               T*APOAQ, V( 1, p ),
                    594:      $                                               1, V( 1, q ), 1 )
1.1       bertrand  595:                                                 CALL DAXPY( MVL,
1.6       bertrand  596:      $                                               -CS*SN*AQOAP,
                    597:      $                                               V( 1, q ), 1,
                    598:      $                                               V( 1, p ), 1 )
1.1       bertrand  599:                                              END IF
                    600:                                           END IF
                    601:                                        END IF
                    602:                                     END IF
                    603:                                  END IF
                    604: *
                    605:                               ELSE
                    606: *              .. have to use modified Gram-Schmidt like transformation
                    607:                                  CALL DCOPY( M, A( 1, p ), 1, WORK, 1 )
                    608:                                  CALL DLASCL( 'G', 0, 0, AAPP, ONE, M,
1.6       bertrand  609:      $                                        1, WORK, LDA, IERR )
1.1       bertrand  610:                                  CALL DLASCL( 'G', 0, 0, AAQQ, ONE, M,
1.6       bertrand  611:      $                                        1, A( 1, q ), LDA, IERR )
1.1       bertrand  612:                                  TEMP1 = -AAPQ*D( p ) / D( q )
                    613:                                  CALL DAXPY( M, TEMP1, WORK, 1,
1.6       bertrand  614:      $                                       A( 1, q ), 1 )
1.1       bertrand  615:                                  CALL DLASCL( 'G', 0, 0, ONE, AAQQ, M,
1.6       bertrand  616:      $                                        1, A( 1, q ), LDA, IERR )
1.1       bertrand  617:                                  SVA( q ) = AAQQ*DSQRT( DMAX1( ZERO,
1.6       bertrand  618:      $                                      ONE-AAPQ*AAPQ ) )
1.1       bertrand  619:                                  MXSINJ = DMAX1( MXSINJ, SFMIN )
                    620:                               END IF
                    621: *           END IF ROTOK THEN ... ELSE
                    622: *
                    623: *           In the case of cancellation in updating SVA(q), SVA(p)
                    624: *           recompute SVA(q), SVA(p).
                    625:                               IF( ( SVA( q ) / AAQQ )**2.LE.ROOTEPS )
1.6       bertrand  626:      $                            THEN
1.1       bertrand  627:                                  IF( ( AAQQ.LT.ROOTBIG ) .AND.
1.6       bertrand  628:      $                               ( AAQQ.GT.ROOTSFMIN ) ) THEN
1.1       bertrand  629:                                     SVA( q ) = DNRM2( M, A( 1, q ), 1 )*
1.6       bertrand  630:      $                                         D( q )
1.1       bertrand  631:                                  ELSE
                    632:                                     T = ZERO
1.4       bertrand  633:                                     AAQQ = ONE
1.1       bertrand  634:                                     CALL DLASSQ( M, A( 1, q ), 1, T,
1.6       bertrand  635:      $                                           AAQQ )
1.1       bertrand  636:                                     SVA( q ) = T*DSQRT( AAQQ )*D( q )
                    637:                                  END IF
                    638:                               END IF
                    639:                               IF( ( AAPP / AAPP0 ).LE.ROOTEPS ) THEN
                    640:                                  IF( ( AAPP.LT.ROOTBIG ) .AND.
1.6       bertrand  641:      $                               ( AAPP.GT.ROOTSFMIN ) ) THEN
1.1       bertrand  642:                                     AAPP = DNRM2( M, A( 1, p ), 1 )*
1.6       bertrand  643:      $                                     D( p )
1.1       bertrand  644:                                  ELSE
                    645:                                     T = ZERO
1.4       bertrand  646:                                     AAPP = ONE
1.1       bertrand  647:                                     CALL DLASSQ( M, A( 1, p ), 1, T,
1.6       bertrand  648:      $                                           AAPP )
1.1       bertrand  649:                                     AAPP = T*DSQRT( AAPP )*D( p )
                    650:                                  END IF
                    651:                                  SVA( p ) = AAPP
                    652:                               END IF
                    653: *
                    654:                            ELSE
                    655: *        A(:,p) and A(:,q) already numerically orthogonal
                    656:                               IF( ir1.EQ.0 )NOTROT = NOTROT + 1
                    657:                               PSKIPPED = PSKIPPED + 1
                    658:                            END IF
                    659:                         ELSE
                    660: *        A(:,q) is zero column
                    661:                            IF( ir1.EQ.0 )NOTROT = NOTROT + 1
                    662:                            PSKIPPED = PSKIPPED + 1
                    663:                         END IF
                    664: *
                    665:                         IF( ( i.LE.SWBAND ) .AND.
1.6       bertrand  666:      $                      ( PSKIPPED.GT.ROWSKIP ) ) THEN
1.1       bertrand  667:                            IF( ir1.EQ.0 )AAPP = -AAPP
                    668:                            NOTROT = 0
                    669:                            GO TO 2103
                    670:                         END IF
                    671: *
                    672:  2002                CONTINUE
                    673: *     END q-LOOP
                    674: *
                    675:  2103                CONTINUE
                    676: *     bailed out of q-loop
                    677: 
                    678:                      SVA( p ) = AAPP
                    679: 
                    680:                   ELSE
                    681:                      SVA( p ) = AAPP
                    682:                      IF( ( ir1.EQ.0 ) .AND. ( AAPP.EQ.ZERO ) )
1.6       bertrand  683:      $                   NOTROT = NOTROT + MIN0( igl+KBL-1, N ) - p
1.1       bertrand  684:                   END IF
                    685: *
                    686:  2001          CONTINUE
                    687: *     end of the p-loop
                    688: *     end of doing the block ( ibr, ibr )
                    689:  1002       CONTINUE
                    690: *     end of ir1-loop
                    691: *
                    692: *........................................................
                    693: * ... go to the off diagonal blocks
                    694: *
                    695:             igl = ( ibr-1 )*KBL + 1
                    696: *
                    697:             DO 2010 jbc = ibr + 1, NBL
                    698: *
                    699:                jgl = ( jbc-1 )*KBL + 1
                    700: *
                    701: *        doing the block at ( ibr, jbc )
                    702: *
                    703:                IJBLSK = 0
                    704:                DO 2100 p = igl, MIN0( igl+KBL-1, N )
                    705: *
                    706:                   AAPP = SVA( p )
                    707: *
                    708:                   IF( AAPP.GT.ZERO ) THEN
                    709: *
                    710:                      PSKIPPED = 0
                    711: *
                    712:                      DO 2200 q = jgl, MIN0( jgl+KBL-1, N )
                    713: *
                    714:                         AAQQ = SVA( q )
                    715: *
                    716:                         IF( AAQQ.GT.ZERO ) THEN
                    717:                            AAPP0 = AAPP
                    718: *
                    719: *     -#- M x 2 Jacobi SVD -#-
                    720: *
                    721: *        -#- Safe Gram matrix computation -#-
                    722: *
                    723:                            IF( AAQQ.GE.ONE ) THEN
                    724:                               IF( AAPP.GE.AAQQ ) THEN
                    725:                                  ROTOK = ( SMALL*AAPP ).LE.AAQQ
                    726:                               ELSE
                    727:                                  ROTOK = ( SMALL*AAQQ ).LE.AAPP
                    728:                               END IF
                    729:                               IF( AAPP.LT.( BIG / AAQQ ) ) THEN
                    730:                                  AAPQ = ( DDOT( M, A( 1, p ), 1, A( 1,
1.6       bertrand  731:      $                                  q ), 1 )*D( p )*D( q ) / AAQQ )
                    732:      $                                  / AAPP
1.1       bertrand  733:                               ELSE
                    734:                                  CALL DCOPY( M, A( 1, p ), 1, WORK, 1 )
                    735:                                  CALL DLASCL( 'G', 0, 0, AAPP, D( p ),
1.6       bertrand  736:      $                                        M, 1, WORK, LDA, IERR )
1.1       bertrand  737:                                  AAPQ = DDOT( M, WORK, 1, A( 1, q ),
1.6       bertrand  738:      $                                  1 )*D( q ) / AAQQ
1.1       bertrand  739:                               END IF
                    740:                            ELSE
                    741:                               IF( AAPP.GE.AAQQ ) THEN
                    742:                                  ROTOK = AAPP.LE.( AAQQ / SMALL )
                    743:                               ELSE
                    744:                                  ROTOK = AAQQ.LE.( AAPP / SMALL )
                    745:                               END IF
                    746:                               IF( AAPP.GT.( SMALL / AAQQ ) ) THEN
                    747:                                  AAPQ = ( DDOT( M, A( 1, p ), 1, A( 1,
1.6       bertrand  748:      $                                  q ), 1 )*D( p )*D( q ) / AAQQ )
                    749:      $                                  / AAPP
1.1       bertrand  750:                               ELSE
                    751:                                  CALL DCOPY( M, A( 1, q ), 1, WORK, 1 )
                    752:                                  CALL DLASCL( 'G', 0, 0, AAQQ, D( q ),
1.6       bertrand  753:      $                                        M, 1, WORK, LDA, IERR )
1.1       bertrand  754:                                  AAPQ = DDOT( M, WORK, 1, A( 1, p ),
1.6       bertrand  755:      $                                  1 )*D( p ) / AAPP
1.1       bertrand  756:                               END IF
                    757:                            END IF
                    758: *
                    759:                            MXAAPQ = DMAX1( MXAAPQ, DABS( AAPQ ) )
                    760: *
                    761: *        TO rotate or NOT to rotate, THAT is the question ...
                    762: *
                    763:                            IF( DABS( AAPQ ).GT.TOL ) THEN
                    764:                               NOTROT = 0
                    765: *           ROTATED  = ROTATED + 1
                    766:                               PSKIPPED = 0
                    767:                               ISWROT = ISWROT + 1
                    768: *
                    769:                               IF( ROTOK ) THEN
                    770: *
                    771:                                  AQOAP = AAQQ / AAPP
                    772:                                  APOAQ = AAPP / AAQQ
1.6       bertrand  773:                                  THETA = -HALF*DABS( AQOAP-APOAQ )/AAPQ
1.1       bertrand  774:                                  IF( AAQQ.GT.AAPP0 )THETA = -THETA
                    775: *
                    776:                                  IF( DABS( THETA ).GT.BIGTHETA ) THEN
                    777:                                     T = HALF / THETA
                    778:                                     FASTR( 3 ) = T*D( p ) / D( q )
                    779:                                     FASTR( 4 ) = -T*D( q ) / D( p )
                    780:                                     CALL DROTM( M, A( 1, p ), 1,
1.6       bertrand  781:      $                                          A( 1, q ), 1, FASTR )
1.1       bertrand  782:                                     IF( RSVEC )CALL DROTM( MVL,
1.6       bertrand  783:      $                                              V( 1, p ), 1,
                    784:      $                                              V( 1, q ), 1,
                    785:      $                                              FASTR )
1.1       bertrand  786:                                     SVA( q ) = AAQQ*DSQRT( DMAX1( ZERO,
1.6       bertrand  787:      $                                         ONE+T*APOAQ*AAPQ ) )
1.1       bertrand  788:                                     AAPP = AAPP*DSQRT( DMAX1( ZERO,
1.6       bertrand  789:      $                                     ONE-T*AQOAP*AAPQ ) )
1.1       bertrand  790:                                     MXSINJ = DMAX1( MXSINJ, DABS( T ) )
                    791:                                  ELSE
                    792: *
                    793: *                 .. choose correct signum for THETA and rotate
                    794: *
                    795:                                     THSIGN = -DSIGN( ONE, AAPQ )
                    796:                                     IF( AAQQ.GT.AAPP0 )THSIGN = -THSIGN
                    797:                                     T = ONE / ( THETA+THSIGN*
1.6       bertrand  798:      $                                  DSQRT( ONE+THETA*THETA ) )
1.1       bertrand  799:                                     CS = DSQRT( ONE / ( ONE+T*T ) )
                    800:                                     SN = T*CS
                    801:                                     MXSINJ = DMAX1( MXSINJ, DABS( SN ) )
                    802:                                     SVA( q ) = AAQQ*DSQRT( DMAX1( ZERO,
1.6       bertrand  803:      $                                         ONE+T*APOAQ*AAPQ ) )
1.4       bertrand  804:                                     AAPP = AAPP*DSQRT( DMAX1( ZERO, 
1.6       bertrand  805:      $                                     ONE-T*AQOAP*AAPQ ) )
1.1       bertrand  806: *
                    807:                                     APOAQ = D( p ) / D( q )
                    808:                                     AQOAP = D( q ) / D( p )
                    809:                                     IF( D( p ).GE.ONE ) THEN
                    810: *
                    811:                                        IF( D( q ).GE.ONE ) THEN
                    812:                                           FASTR( 3 ) = T*APOAQ
                    813:                                           FASTR( 4 ) = -T*AQOAP
                    814:                                           D( p ) = D( p )*CS
                    815:                                           D( q ) = D( q )*CS
                    816:                                           CALL DROTM( M, A( 1, p ), 1,
1.6       bertrand  817:      $                                                A( 1, q ), 1,
                    818:      $                                                FASTR )
1.1       bertrand  819:                                           IF( RSVEC )CALL DROTM( MVL,
1.6       bertrand  820:      $                                        V( 1, p ), 1, V( 1, q ),
                    821:      $                                        1, FASTR )
1.1       bertrand  822:                                        ELSE
                    823:                                           CALL DAXPY( M, -T*AQOAP,
1.6       bertrand  824:      $                                                A( 1, q ), 1,
                    825:      $                                                A( 1, p ), 1 )
1.1       bertrand  826:                                           CALL DAXPY( M, CS*SN*APOAQ,
1.6       bertrand  827:      $                                                A( 1, p ), 1,
                    828:      $                                                A( 1, q ), 1 )
1.1       bertrand  829:                                           IF( RSVEC ) THEN
                    830:                                              CALL DAXPY( MVL, -T*AQOAP,
1.6       bertrand  831:      $                                                   V( 1, q ), 1,
                    832:      $                                                   V( 1, p ), 1 )
1.1       bertrand  833:                                              CALL DAXPY( MVL,
1.6       bertrand  834:      $                                                   CS*SN*APOAQ,
                    835:      $                                                   V( 1, p ), 1,
                    836:      $                                                   V( 1, q ), 1 )
1.1       bertrand  837:                                           END IF
                    838:                                           D( p ) = D( p )*CS
                    839:                                           D( q ) = D( q ) / CS
                    840:                                        END IF
                    841:                                     ELSE
                    842:                                        IF( D( q ).GE.ONE ) THEN
                    843:                                           CALL DAXPY( M, T*APOAQ,
1.6       bertrand  844:      $                                                A( 1, p ), 1,
                    845:      $                                                A( 1, q ), 1 )
1.1       bertrand  846:                                           CALL DAXPY( M, -CS*SN*AQOAP,
1.6       bertrand  847:      $                                                A( 1, q ), 1,
                    848:      $                                                A( 1, p ), 1 )
1.1       bertrand  849:                                           IF( RSVEC ) THEN
                    850:                                              CALL DAXPY( MVL, T*APOAQ,
1.6       bertrand  851:      $                                                   V( 1, p ), 1,
                    852:      $                                                   V( 1, q ), 1 )
1.1       bertrand  853:                                              CALL DAXPY( MVL,
1.6       bertrand  854:      $                                                   -CS*SN*AQOAP,
                    855:      $                                                   V( 1, q ), 1,
                    856:      $                                                   V( 1, p ), 1 )
1.1       bertrand  857:                                           END IF
                    858:                                           D( p ) = D( p ) / CS
                    859:                                           D( q ) = D( q )*CS
                    860:                                        ELSE
                    861:                                           IF( D( p ).GE.D( q ) ) THEN
                    862:                                              CALL DAXPY( M, -T*AQOAP,
1.6       bertrand  863:      $                                                   A( 1, q ), 1,
                    864:      $                                                   A( 1, p ), 1 )
1.1       bertrand  865:                                              CALL DAXPY( M, CS*SN*APOAQ,
1.6       bertrand  866:      $                                                   A( 1, p ), 1,
                    867:      $                                                   A( 1, q ), 1 )
1.1       bertrand  868:                                              D( p ) = D( p )*CS
                    869:                                              D( q ) = D( q ) / CS
                    870:                                              IF( RSVEC ) THEN
                    871:                                                 CALL DAXPY( MVL,
1.6       bertrand  872:      $                                               -T*AQOAP,
                    873:      $                                               V( 1, q ), 1,
                    874:      $                                               V( 1, p ), 1 )
1.1       bertrand  875:                                                 CALL DAXPY( MVL,
1.6       bertrand  876:      $                                               CS*SN*APOAQ,
                    877:      $                                               V( 1, p ), 1,
                    878:      $                                               V( 1, q ), 1 )
1.1       bertrand  879:                                              END IF
                    880:                                           ELSE
                    881:                                              CALL DAXPY( M, T*APOAQ,
1.6       bertrand  882:      $                                                   A( 1, p ), 1,
                    883:      $                                                   A( 1, q ), 1 )
1.1       bertrand  884:                                              CALL DAXPY( M,
1.6       bertrand  885:      $                                                   -CS*SN*AQOAP,
                    886:      $                                                   A( 1, q ), 1,
                    887:      $                                                   A( 1, p ), 1 )
1.1       bertrand  888:                                              D( p ) = D( p ) / CS
                    889:                                              D( q ) = D( q )*CS
                    890:                                              IF( RSVEC ) THEN
                    891:                                                 CALL DAXPY( MVL,
1.6       bertrand  892:      $                                               T*APOAQ, V( 1, p ),
                    893:      $                                               1, V( 1, q ), 1 )
1.1       bertrand  894:                                                 CALL DAXPY( MVL,
1.6       bertrand  895:      $                                               -CS*SN*AQOAP,
                    896:      $                                               V( 1, q ), 1,
                    897:      $                                               V( 1, p ), 1 )
1.1       bertrand  898:                                              END IF
                    899:                                           END IF
                    900:                                        END IF
                    901:                                     END IF
                    902:                                  END IF
                    903: *
                    904:                               ELSE
                    905:                                  IF( AAPP.GT.AAQQ ) THEN
                    906:                                     CALL DCOPY( M, A( 1, p ), 1, WORK,
1.6       bertrand  907:      $                                          1 )
1.1       bertrand  908:                                     CALL DLASCL( 'G', 0, 0, AAPP, ONE,
1.6       bertrand  909:      $                                           M, 1, WORK, LDA, IERR )
1.1       bertrand  910:                                     CALL DLASCL( 'G', 0, 0, AAQQ, ONE,
1.6       bertrand  911:      $                                           M, 1, A( 1, q ), LDA,
                    912:      $                                           IERR )
1.1       bertrand  913:                                     TEMP1 = -AAPQ*D( p ) / D( q )
                    914:                                     CALL DAXPY( M, TEMP1, WORK, 1,
1.6       bertrand  915:      $                                          A( 1, q ), 1 )
1.1       bertrand  916:                                     CALL DLASCL( 'G', 0, 0, ONE, AAQQ,
1.6       bertrand  917:      $                                           M, 1, A( 1, q ), LDA,
                    918:      $                                           IERR )
1.1       bertrand  919:                                     SVA( q ) = AAQQ*DSQRT( DMAX1( ZERO,
1.6       bertrand  920:      $                                         ONE-AAPQ*AAPQ ) )
1.1       bertrand  921:                                     MXSINJ = DMAX1( MXSINJ, SFMIN )
                    922:                                  ELSE
                    923:                                     CALL DCOPY( M, A( 1, q ), 1, WORK,
1.6       bertrand  924:      $                                          1 )
1.1       bertrand  925:                                     CALL DLASCL( 'G', 0, 0, AAQQ, ONE,
1.6       bertrand  926:      $                                           M, 1, WORK, LDA, IERR )
1.1       bertrand  927:                                     CALL DLASCL( 'G', 0, 0, AAPP, ONE,
1.6       bertrand  928:      $                                           M, 1, A( 1, p ), LDA,
                    929:      $                                           IERR )
1.1       bertrand  930:                                     TEMP1 = -AAPQ*D( q ) / D( p )
                    931:                                     CALL DAXPY( M, TEMP1, WORK, 1,
1.6       bertrand  932:      $                                          A( 1, p ), 1 )
1.1       bertrand  933:                                     CALL DLASCL( 'G', 0, 0, ONE, AAPP,
1.6       bertrand  934:      $                                           M, 1, A( 1, p ), LDA,
                    935:      $                                           IERR )
1.1       bertrand  936:                                     SVA( p ) = AAPP*DSQRT( DMAX1( ZERO,
1.6       bertrand  937:      $                                         ONE-AAPQ*AAPQ ) )
1.1       bertrand  938:                                     MXSINJ = DMAX1( MXSINJ, SFMIN )
                    939:                                  END IF
                    940:                               END IF
                    941: *           END IF ROTOK THEN ... ELSE
                    942: *
                    943: *           In the case of cancellation in updating SVA(q)
                    944: *           .. recompute SVA(q)
                    945:                               IF( ( SVA( q ) / AAQQ )**2.LE.ROOTEPS )
1.6       bertrand  946:      $                            THEN
1.1       bertrand  947:                                  IF( ( AAQQ.LT.ROOTBIG ) .AND.
1.6       bertrand  948:      $                               ( AAQQ.GT.ROOTSFMIN ) ) THEN
1.1       bertrand  949:                                     SVA( q ) = DNRM2( M, A( 1, q ), 1 )*
1.6       bertrand  950:      $                                         D( q )
1.1       bertrand  951:                                  ELSE
                    952:                                     T = ZERO
1.4       bertrand  953:                                     AAQQ = ONE
1.1       bertrand  954:                                     CALL DLASSQ( M, A( 1, q ), 1, T,
1.6       bertrand  955:      $                                           AAQQ )
1.1       bertrand  956:                                     SVA( q ) = T*DSQRT( AAQQ )*D( q )
                    957:                                  END IF
                    958:                               END IF
                    959:                               IF( ( AAPP / AAPP0 )**2.LE.ROOTEPS ) THEN
                    960:                                  IF( ( AAPP.LT.ROOTBIG ) .AND.
1.6       bertrand  961:      $                               ( AAPP.GT.ROOTSFMIN ) ) THEN
1.1       bertrand  962:                                     AAPP = DNRM2( M, A( 1, p ), 1 )*
1.6       bertrand  963:      $                                     D( p )
1.1       bertrand  964:                                  ELSE
                    965:                                     T = ZERO
1.4       bertrand  966:                                     AAPP = ONE
1.1       bertrand  967:                                     CALL DLASSQ( M, A( 1, p ), 1, T,
1.6       bertrand  968:      $                                           AAPP )
1.1       bertrand  969:                                     AAPP = T*DSQRT( AAPP )*D( p )
                    970:                                  END IF
                    971:                                  SVA( p ) = AAPP
                    972:                               END IF
                    973: *              end of OK rotation
                    974:                            ELSE
                    975:                               NOTROT = NOTROT + 1
                    976:                               PSKIPPED = PSKIPPED + 1
                    977:                               IJBLSK = IJBLSK + 1
                    978:                            END IF
                    979:                         ELSE
                    980:                            NOTROT = NOTROT + 1
                    981:                            PSKIPPED = PSKIPPED + 1
                    982:                            IJBLSK = IJBLSK + 1
                    983:                         END IF
                    984: *
                    985:                         IF( ( i.LE.SWBAND ) .AND. ( IJBLSK.GE.BLSKIP ) )
1.6       bertrand  986:      $                      THEN
1.1       bertrand  987:                            SVA( p ) = AAPP
                    988:                            NOTROT = 0
                    989:                            GO TO 2011
                    990:                         END IF
                    991:                         IF( ( i.LE.SWBAND ) .AND.
1.6       bertrand  992:      $                      ( PSKIPPED.GT.ROWSKIP ) ) THEN
1.1       bertrand  993:                            AAPP = -AAPP
                    994:                            NOTROT = 0
                    995:                            GO TO 2203
                    996:                         END IF
                    997: *
                    998:  2200                CONTINUE
                    999: *        end of the q-loop
                   1000:  2203                CONTINUE
                   1001: *
                   1002:                      SVA( p ) = AAPP
                   1003: *
                   1004:                   ELSE
                   1005:                      IF( AAPP.EQ.ZERO )NOTROT = NOTROT +
1.6       bertrand 1006:      $                   MIN0( jgl+KBL-1, N ) - jgl + 1
1.1       bertrand 1007:                      IF( AAPP.LT.ZERO )NOTROT = 0
                   1008:                   END IF
                   1009: 
                   1010:  2100          CONTINUE
                   1011: *     end of the p-loop
                   1012:  2010       CONTINUE
                   1013: *     end of the jbc-loop
                   1014:  2011       CONTINUE
                   1015: *2011 bailed out of the jbc-loop
                   1016:             DO 2012 p = igl, MIN0( igl+KBL-1, N )
                   1017:                SVA( p ) = DABS( SVA( p ) )
                   1018:  2012       CONTINUE
                   1019: *
                   1020:  2000    CONTINUE
                   1021: *2000 :: end of the ibr-loop
                   1022: *
                   1023: *     .. update SVA(N)
                   1024:          IF( ( SVA( N ).LT.ROOTBIG ) .AND. ( SVA( N ).GT.ROOTSFMIN ) )
1.6       bertrand 1025:      $       THEN
1.1       bertrand 1026:             SVA( N ) = DNRM2( M, A( 1, N ), 1 )*D( N )
                   1027:          ELSE
                   1028:             T = ZERO
1.4       bertrand 1029:             AAPP = ONE
1.1       bertrand 1030:             CALL DLASSQ( M, A( 1, N ), 1, T, AAPP )
                   1031:             SVA( N ) = T*DSQRT( AAPP )*D( N )
                   1032:          END IF
                   1033: *
                   1034: *     Additional steering devices
                   1035: *
                   1036:          IF( ( i.LT.SWBAND ) .AND. ( ( MXAAPQ.LE.ROOTTOL ) .OR.
1.6       bertrand 1037:      $       ( ISWROT.LE.N ) ) )SWBAND = i
1.1       bertrand 1038: *
                   1039:          IF( ( i.GT.SWBAND+1 ) .AND. ( MXAAPQ.LT.DBLE( N )*TOL ) .AND.
1.6       bertrand 1040:      $       ( DBLE( N )*MXAAPQ*MXSINJ.LT.TOL ) ) THEN
1.1       bertrand 1041:             GO TO 1994
                   1042:          END IF
                   1043: *
                   1044:          IF( NOTROT.GE.EMPTSW )GO TO 1994
                   1045: 
                   1046:  1993 CONTINUE
                   1047: *     end i=1:NSWEEP loop
                   1048: * #:) Reaching this point means that the procedure has comleted the given
                   1049: *     number of iterations.
                   1050:       INFO = NSWEEP - 1
                   1051:       GO TO 1995
                   1052:  1994 CONTINUE
                   1053: * #:) Reaching this point means that during the i-th sweep all pivots were
                   1054: *     below the given tolerance, causing early exit.
                   1055: *
                   1056:       INFO = 0
                   1057: * #:) INFO = 0 confirms successful iterations.
                   1058:  1995 CONTINUE
                   1059: *
                   1060: *     Sort the vector D.
                   1061:       DO 5991 p = 1, N - 1
                   1062:          q = IDAMAX( N-p+1, SVA( p ), 1 ) + p - 1
                   1063:          IF( p.NE.q ) THEN
                   1064:             TEMP1 = SVA( p )
                   1065:             SVA( p ) = SVA( q )
                   1066:             SVA( q ) = TEMP1
                   1067:             TEMP1 = D( p )
                   1068:             D( p ) = D( q )
                   1069:             D( q ) = TEMP1
                   1070:             CALL DSWAP( M, A( 1, p ), 1, A( 1, q ), 1 )
                   1071:             IF( RSVEC )CALL DSWAP( MVL, V( 1, p ), 1, V( 1, q ), 1 )
                   1072:          END IF
                   1073:  5991 CONTINUE
                   1074: *
                   1075:       RETURN
                   1076: *     ..
                   1077: *     .. END OF DGSVJ0
                   1078: *     ..
                   1079:       END

CVSweb interface <joel.bertrand@systella.fr>