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

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

CVSweb interface <joel.bertrand@systella.fr>