Annotation of rpl/lapack/lapack/dgesvj.f, revision 1.12

1.7       bertrand    1: *> \brief \b DGESVJ
                      2: *
                      3: *  =========== DOCUMENTATION ===========
                      4: *
                      5: * Online html documentation available at 
                      6: *            http://www.netlib.org/lapack/explore-html/ 
                      7: *
                      8: *> \htmlonly
                      9: *> Download DGESVJ + dependencies 
                     10: *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dgesvj.f"> 
                     11: *> [TGZ]</a> 
                     12: *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dgesvj.f"> 
                     13: *> [ZIP]</a> 
                     14: *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dgesvj.f"> 
                     15: *> [TXT]</a>
                     16: *> \endhtmlonly 
                     17: *
                     18: *  Definition:
                     19: *  ===========
                     20: *
                     21: *       SUBROUTINE DGESVJ( JOBA, JOBU, JOBV, M, N, A, LDA, SVA, MV, V,
                     22: *                          LDV, WORK, LWORK, INFO )
                     23: * 
                     24: *       .. Scalar Arguments ..
                     25: *       INTEGER            INFO, LDA, LDV, LWORK, M, MV, N
                     26: *       CHARACTER*1        JOBA, JOBU, JOBV
                     27: *       ..
                     28: *       .. Array Arguments ..
                     29: *       DOUBLE PRECISION   A( LDA, * ), SVA( N ), V( LDV, * ),
                     30: *      $                   WORK( LWORK )
                     31: *       ..
                     32: *  
                     33: *
                     34: *> \par Purpose:
                     35: *  =============
                     36: *>
                     37: *> \verbatim
                     38: *>
                     39: *> DGESVJ computes the singular value decomposition (SVD) of a real
                     40: *> M-by-N matrix A, where M >= N. The SVD of A is written as
                     41: *>                                    [++]   [xx]   [x0]   [xx]
                     42: *>              A = U * SIGMA * V^t,  [++] = [xx] * [ox] * [xx]
                     43: *>                                    [++]   [xx]
                     44: *> where SIGMA is an N-by-N diagonal matrix, U is an M-by-N orthonormal
                     45: *> matrix, and V is an N-by-N orthogonal matrix. The diagonal elements
                     46: *> of SIGMA are the singular values of A. The columns of U and V are the
                     47: *> left and the right singular vectors of A, respectively.
                     48: *> \endverbatim
                     49: *
                     50: *  Arguments:
                     51: *  ==========
                     52: *
                     53: *> \param[in] JOBA
                     54: *> \verbatim
                     55: *>          JOBA is CHARACTER* 1
                     56: *>          Specifies the structure of A.
                     57: *>          = 'L': The input matrix A is lower triangular;
                     58: *>          = 'U': The input matrix A is upper triangular;
                     59: *>          = 'G': The input matrix A is general M-by-N matrix, M >= N.
                     60: *> \endverbatim
                     61: *>
                     62: *> \param[in] JOBU
                     63: *> \verbatim
                     64: *>          JOBU is CHARACTER*1
                     65: *>          Specifies whether to compute the left singular vectors
                     66: *>          (columns of U):
                     67: *>          = 'U': The left singular vectors corresponding to the nonzero
                     68: *>                 singular values are computed and returned in the leading
                     69: *>                 columns of A. See more details in the description of A.
                     70: *>                 The default numerical orthogonality threshold is set to
                     71: *>                 approximately TOL=CTOL*EPS, CTOL=DSQRT(M), EPS=DLAMCH('E').
                     72: *>          = 'C': Analogous to JOBU='U', except that user can control the
                     73: *>                 level of numerical orthogonality of the computed left
                     74: *>                 singular vectors. TOL can be set to TOL = CTOL*EPS, where
                     75: *>                 CTOL is given on input in the array WORK.
                     76: *>                 No CTOL smaller than ONE is allowed. CTOL greater
                     77: *>                 than 1 / EPS is meaningless. The option 'C'
                     78: *>                 can be used if M*EPS is satisfactory orthogonality
                     79: *>                 of the computed left singular vectors, so CTOL=M could
                     80: *>                 save few sweeps of Jacobi rotations.
                     81: *>                 See the descriptions of A and WORK(1).
                     82: *>          = 'N': The matrix U is not computed. However, see the
                     83: *>                 description of A.
                     84: *> \endverbatim
                     85: *>
                     86: *> \param[in] JOBV
                     87: *> \verbatim
                     88: *>          JOBV is CHARACTER*1
                     89: *>          Specifies whether to compute the right singular vectors, that
                     90: *>          is, the matrix V:
                     91: *>          = 'V' : the matrix V is computed and returned in the array V
                     92: *>          = 'A' : the Jacobi rotations are applied to the MV-by-N
                     93: *>                  array V. In other words, the right singular vector
                     94: *>                  matrix V is not computed explicitly, instead it is
                     95: *>                  applied to an MV-by-N matrix initially stored in the
                     96: *>                  first MV rows of V.
                     97: *>          = 'N' : the matrix V is not computed and the array V is not
                     98: *>                  referenced
                     99: *> \endverbatim
                    100: *>
                    101: *> \param[in] M
                    102: *> \verbatim
                    103: *>          M is INTEGER
                    104: *>          The number of rows of the input matrix A. 1/DLAMCH('E') > M >= 0.  
                    105: *> \endverbatim
                    106: *>
                    107: *> \param[in] N
                    108: *> \verbatim
                    109: *>          N is INTEGER
                    110: *>          The number of columns of the input matrix A.
                    111: *>          M >= N >= 0.
                    112: *> \endverbatim
                    113: *>
                    114: *> \param[in,out] A
                    115: *> \verbatim
                    116: *>          A is DOUBLE PRECISION array, dimension (LDA,N)
                    117: *>          On entry, the M-by-N matrix A.
                    118: *>          On exit :
                    119: *>          If JOBU .EQ. 'U' .OR. JOBU .EQ. 'C' :
                    120: *>                 If INFO .EQ. 0 :
                    121: *>                 RANKA orthonormal columns of U are returned in the
                    122: *>                 leading RANKA columns of the array A. Here RANKA <= N
                    123: *>                 is the number of computed singular values of A that are
                    124: *>                 above the underflow threshold DLAMCH('S'). The singular
                    125: *>                 vectors corresponding to underflowed or zero singular
                    126: *>                 values are not computed. The value of RANKA is returned
                    127: *>                 in the array WORK as RANKA=NINT(WORK(2)). Also see the
                    128: *>                 descriptions of SVA and WORK. The computed columns of U
                    129: *>                 are mutually numerically orthogonal up to approximately
                    130: *>                 TOL=DSQRT(M)*EPS (default); or TOL=CTOL*EPS (JOBU.EQ.'C'),
                    131: *>                 see the description of JOBU.
                    132: *>                 If INFO .GT. 0 :
                    133: *>                 the procedure DGESVJ did not converge in the given number
                    134: *>                 of iterations (sweeps). In that case, the computed
                    135: *>                 columns of U may not be orthogonal up to TOL. The output
                    136: *>                 U (stored in A), SIGMA (given by the computed singular
                    137: *>                 values in SVA(1:N)) and V is still a decomposition of the
                    138: *>                 input matrix A in the sense that the residual
                    139: *>                 ||A-SCALE*U*SIGMA*V^T||_2 / ||A||_2 is small.
                    140: *>
                    141: *>          If JOBU .EQ. 'N' :
                    142: *>                 If INFO .EQ. 0 :
                    143: *>                 Note that the left singular vectors are 'for free' in the
                    144: *>                 one-sided Jacobi SVD algorithm. However, if only the
                    145: *>                 singular values are needed, the level of numerical
                    146: *>                 orthogonality of U is not an issue and iterations are
                    147: *>                 stopped when the columns of the iterated matrix are
                    148: *>                 numerically orthogonal up to approximately M*EPS. Thus,
                    149: *>                 on exit, A contains the columns of U scaled with the
                    150: *>                 corresponding singular values.
                    151: *>                 If INFO .GT. 0 :
                    152: *>                 the procedure DGESVJ did not converge in the given number
                    153: *>                 of iterations (sweeps).
                    154: *> \endverbatim
                    155: *>
                    156: *> \param[in] LDA
                    157: *> \verbatim
                    158: *>          LDA is INTEGER
                    159: *>          The leading dimension of the array A.  LDA >= max(1,M).
                    160: *> \endverbatim
                    161: *>
                    162: *> \param[out] SVA
                    163: *> \verbatim
                    164: *>          SVA is DOUBLE PRECISION array, dimension (N)
                    165: *>          On exit :
                    166: *>          If INFO .EQ. 0 :
                    167: *>          depending on the value SCALE = WORK(1), we have:
                    168: *>                 If SCALE .EQ. ONE :
                    169: *>                 SVA(1:N) contains the computed singular values of A.
                    170: *>                 During the computation SVA contains the Euclidean column
                    171: *>                 norms of the iterated matrices in the array A.
                    172: *>                 If SCALE .NE. ONE :
                    173: *>                 The singular values of A are SCALE*SVA(1:N), and this
                    174: *>                 factored representation is due to the fact that some of the
                    175: *>                 singular values of A might underflow or overflow.
                    176: *>          If INFO .GT. 0 :
                    177: *>          the procedure DGESVJ did not converge in the given number of
                    178: *>          iterations (sweeps) and SCALE*SVA(1:N) may not be accurate.
                    179: *> \endverbatim
                    180: *>
                    181: *> \param[in] MV
                    182: *> \verbatim
                    183: *>          MV is INTEGER
                    184: *>          If JOBV .EQ. 'A', then the product of Jacobi rotations in DGESVJ
                    185: *>          is applied to the first MV rows of V. See the description of JOBV.
                    186: *> \endverbatim
                    187: *>
                    188: *> \param[in,out] V
                    189: *> \verbatim
                    190: *>          V is DOUBLE PRECISION array, dimension (LDV,N)
                    191: *>          If JOBV = 'V', then V contains on exit the N-by-N matrix of
                    192: *>                         the right singular vectors;
                    193: *>          If JOBV = 'A', then V contains the product of the computed right
                    194: *>                         singular vector matrix and the initial matrix in
                    195: *>                         the array V.
                    196: *>          If JOBV = 'N', then V is not referenced.
                    197: *> \endverbatim
                    198: *>
                    199: *> \param[in] LDV
                    200: *> \verbatim
                    201: *>          LDV is INTEGER
                    202: *>          The leading dimension of the array V, LDV .GE. 1.
                    203: *>          If JOBV .EQ. 'V', then LDV .GE. max(1,N).
                    204: *>          If JOBV .EQ. 'A', then LDV .GE. max(1,MV) .
                    205: *> \endverbatim
                    206: *>
                    207: *> \param[in,out] WORK
                    208: *> \verbatim
                    209: *>          WORK is DOUBLE PRECISION array, dimension max(4,M+N).
                    210: *>          On entry :
                    211: *>          If JOBU .EQ. 'C' :
                    212: *>          WORK(1) = CTOL, where CTOL defines the threshold for convergence.
                    213: *>                    The process stops if all columns of A are mutually
                    214: *>                    orthogonal up to CTOL*EPS, EPS=DLAMCH('E').
                    215: *>                    It is required that CTOL >= ONE, i.e. it is not
                    216: *>                    allowed to force the routine to obtain orthogonality
                    217: *>                    below EPS.
                    218: *>          On exit :
                    219: *>          WORK(1) = SCALE is the scaling factor such that SCALE*SVA(1:N)
                    220: *>                    are the computed singular values of A.
                    221: *>                    (See description of SVA().)
                    222: *>          WORK(2) = NINT(WORK(2)) is the number of the computed nonzero
                    223: *>                    singular values.
                    224: *>          WORK(3) = NINT(WORK(3)) is the number of the computed singular
                    225: *>                    values that are larger than the underflow threshold.
                    226: *>          WORK(4) = NINT(WORK(4)) is the number of sweeps of Jacobi
                    227: *>                    rotations needed for numerical convergence.
                    228: *>          WORK(5) = max_{i.NE.j} |COS(A(:,i),A(:,j))| in the last sweep.
                    229: *>                    This is useful information in cases when DGESVJ did
                    230: *>                    not converge, as it can be used to estimate whether
                    231: *>                    the output is stil useful and for post festum analysis.
                    232: *>          WORK(6) = the largest absolute value over all sines of the
                    233: *>                    Jacobi rotation angles in the last sweep. It can be
                    234: *>                    useful for a post festum analysis.
                    235: *> \endverbatim
                    236: *>
                    237: *> \param[in] LWORK
                    238: *> \verbatim
                    239: *>          LWORK is INTEGER
                    240: *>          length of WORK, WORK >= MAX(6,M+N)
                    241: *> \endverbatim
                    242: *>
                    243: *> \param[out] INFO
                    244: *> \verbatim
                    245: *>          INFO is INTEGER
                    246: *>          = 0 : successful exit.
                    247: *>          < 0 : if INFO = -i, then the i-th argument had an illegal value
                    248: *>          > 0 : DGESVJ did not converge in the maximal allowed number (30)
                    249: *>                of sweeps. The output may still be useful. See the
                    250: *>                description of WORK.
                    251: *> \endverbatim
                    252: *
                    253: *  Authors:
                    254: *  ========
                    255: *
                    256: *> \author Univ. of Tennessee 
                    257: *> \author Univ. of California Berkeley 
                    258: *> \author Univ. of Colorado Denver 
                    259: *> \author NAG Ltd. 
                    260: *
1.11      bertrand  261: *> \date September 2012
1.7       bertrand  262: *
                    263: *> \ingroup doubleGEcomputational
                    264: *
                    265: *> \par Further Details:
                    266: *  =====================
                    267: *>
                    268: *> \verbatim
                    269: *>
                    270: *>  The orthogonal N-by-N matrix V is obtained as a product of Jacobi plane
                    271: *>  rotations. The rotations are implemented as fast scaled rotations of
                    272: *>  Anda and Park [1]. In the case of underflow of the Jacobi angle, a
                    273: *>  modified Jacobi transformation of Drmac [4] is used. Pivot strategy uses
                    274: *>  column interchanges of de Rijk [2]. The relative accuracy of the computed
                    275: *>  singular values and the accuracy of the computed singular vectors (in
                    276: *>  angle metric) is as guaranteed by the theory of Demmel and Veselic [3].
                    277: *>  The condition number that determines the accuracy in the full rank case
                    278: *>  is essentially min_{D=diag} kappa(A*D), where kappa(.) is the
                    279: *>  spectral condition number. The best performance of this Jacobi SVD
                    280: *>  procedure is achieved if used in an  accelerated version of Drmac and
                    281: *>  Veselic [5,6], and it is the kernel routine in the SIGMA library [7].
                    282: *>  Some tunning parameters (marked with [TP]) are available for the
                    283: *>  implementer.
                    284: *>  The computational range for the nonzero singular values is the  machine
                    285: *>  number interval ( UNDERFLOW , OVERFLOW ). In extreme cases, even
                    286: *>  denormalized singular values can be computed with the corresponding
                    287: *>  gradual loss of accurate digits.
                    288: *> \endverbatim
                    289: *
                    290: *> \par Contributors:
                    291: *  ==================
                    292: *>
                    293: *> \verbatim
                    294: *>
                    295: *>  ============
                    296: *>
                    297: *>  Zlatko Drmac (Zagreb, Croatia) and Kresimir Veselic (Hagen, Germany)
                    298: *> \endverbatim
                    299: *
                    300: *> \par References:
                    301: *  ================
                    302: *>
                    303: *> \verbatim
                    304: *>
                    305: *> [1] A. A. Anda and H. Park: Fast plane rotations with dynamic scaling.
                    306: *>     SIAM J. matrix Anal. Appl., Vol. 15 (1994), pp. 162-174.
                    307: *> [2] P. P. M. De Rijk: A one-sided Jacobi algorithm for computing the
                    308: *>     singular value decomposition on a vector computer.
                    309: *>     SIAM J. Sci. Stat. Comp., Vol. 10 (1998), pp. 359-371.
                    310: *> [3] J. Demmel and K. Veselic: Jacobi method is more accurate than QR.
                    311: *> [4] Z. Drmac: Implementation of Jacobi rotations for accurate singular
                    312: *>     value computation in floating point arithmetic.
                    313: *>     SIAM J. Sci. Comp., Vol. 18 (1997), pp. 1200-1222.
                    314: *> [5] Z. Drmac and K. Veselic: New fast and accurate Jacobi SVD algorithm I.
                    315: *>     SIAM J. Matrix Anal. Appl. Vol. 35, No. 2 (2008), pp. 1322-1342.
                    316: *>     LAPACK Working note 169.
                    317: *> [6] Z. Drmac and K. Veselic: New fast and accurate Jacobi SVD algorithm II.
                    318: *>     SIAM J. Matrix Anal. Appl. Vol. 35, No. 2 (2008), pp. 1343-1362.
                    319: *>     LAPACK Working note 170.
                    320: *> [7] Z. Drmac: SIGMA - mathematical software library for accurate SVD, PSV,
                    321: *>     QSVD, (H,K)-SVD computations.
                    322: *>     Department of Mathematics, University of Zagreb, 2008.
                    323: *> \endverbatim
                    324: *
                    325: *>  \par Bugs, examples and comments:
                    326: *   =================================
                    327: *>
                    328: *> \verbatim
                    329: *>  ===========================
                    330: *>  Please report all bugs and send interesting test examples and comments to
                    331: *>  drmac@math.hr. Thank you.
                    332: *> \endverbatim
                    333: *>
                    334: *  =====================================================================
1.1       bertrand  335:       SUBROUTINE DGESVJ( JOBA, JOBU, JOBV, M, N, A, LDA, SVA, MV, V,
1.6       bertrand  336:      $                   LDV, WORK, LWORK, INFO )
1.1       bertrand  337: *
1.11      bertrand  338: *  -- LAPACK computational routine (version 3.4.2) --
1.1       bertrand  339: *  -- LAPACK is a software package provided by Univ. of Tennessee,    --
                    340: *  -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
1.11      bertrand  341: *     September 2012
1.1       bertrand  342: *
                    343: *     .. Scalar Arguments ..
                    344:       INTEGER            INFO, LDA, LDV, LWORK, M, MV, N
                    345:       CHARACTER*1        JOBA, JOBU, JOBV
                    346: *     ..
                    347: *     .. Array Arguments ..
                    348:       DOUBLE PRECISION   A( LDA, * ), SVA( N ), V( LDV, * ),
1.6       bertrand  349:      $                   WORK( LWORK )
1.1       bertrand  350: *     ..
                    351: *
                    352: *  =====================================================================
                    353: *
                    354: *     .. Local Parameters ..
1.9       bertrand  355:       DOUBLE PRECISION   ZERO, HALF, ONE
                    356:       PARAMETER          ( ZERO = 0.0D0, HALF = 0.5D0, ONE = 1.0D0)
1.1       bertrand  357:       INTEGER            NSWEEP
                    358:       PARAMETER          ( NSWEEP = 30 )
                    359: *     ..
                    360: *     .. Local Scalars ..
                    361:       DOUBLE PRECISION   AAPP, AAPP0, AAPQ, AAQQ, APOAQ, AQOAP, BIG,
1.6       bertrand  362:      $                   BIGTHETA, CS, CTOL, EPSLN, LARGE, MXAAPQ,
                    363:      $                   MXSINJ, ROOTBIG, ROOTEPS, ROOTSFMIN, ROOTTOL,
                    364:      $                   SKL, SFMIN, SMALL, SN, T, TEMP1, THETA,
                    365:      $                   THSIGN, TOL
1.1       bertrand  366:       INTEGER            BLSKIP, EMPTSW, i, ibr, IERR, igl, IJBLSK, ir1,
1.6       bertrand  367:      $                   ISWROT, jbc, jgl, KBL, LKAHEAD, MVL, N2, N34,
                    368:      $                   N4, NBL, NOTROT, p, PSKIPPED, q, ROWSKIP,
                    369:      $                   SWBAND
1.1       bertrand  370:       LOGICAL            APPLV, GOSCALE, LOWER, LSVEC, NOSCALE, ROTOK,
1.6       bertrand  371:      $                   RSVEC, UCTOL, UPPER
1.1       bertrand  372: *     ..
                    373: *     .. Local Arrays ..
                    374:       DOUBLE PRECISION   FASTR( 5 )
                    375: *     ..
                    376: *     .. Intrinsic Functions ..
                    377:       INTRINSIC          DABS, DMAX1, DMIN1, DBLE, MIN0, DSIGN, DSQRT
                    378: *     ..
                    379: *     .. External Functions ..
                    380: *     ..
                    381: *     from BLAS
                    382:       DOUBLE PRECISION   DDOT, DNRM2
                    383:       EXTERNAL           DDOT, DNRM2
                    384:       INTEGER            IDAMAX
                    385:       EXTERNAL           IDAMAX
                    386: *     from LAPACK
                    387:       DOUBLE PRECISION   DLAMCH
                    388:       EXTERNAL           DLAMCH
                    389:       LOGICAL            LSAME
                    390:       EXTERNAL           LSAME
                    391: *     ..
                    392: *     .. External Subroutines ..
                    393: *     ..
                    394: *     from BLAS
                    395:       EXTERNAL           DAXPY, DCOPY, DROTM, DSCAL, DSWAP
                    396: *     from LAPACK
                    397:       EXTERNAL           DLASCL, DLASET, DLASSQ, XERBLA
                    398: *
                    399:       EXTERNAL           DGSVJ0, DGSVJ1
                    400: *     ..
                    401: *     .. Executable Statements ..
                    402: *
                    403: *     Test the input arguments
                    404: *
                    405:       LSVEC = LSAME( JOBU, 'U' )
                    406:       UCTOL = LSAME( JOBU, 'C' )
                    407:       RSVEC = LSAME( JOBV, 'V' )
                    408:       APPLV = LSAME( JOBV, 'A' )
                    409:       UPPER = LSAME( JOBA, 'U' )
                    410:       LOWER = LSAME( JOBA, 'L' )
                    411: *
                    412:       IF( .NOT.( UPPER .OR. LOWER .OR. LSAME( JOBA, 'G' ) ) ) THEN
                    413:          INFO = -1
                    414:       ELSE IF( .NOT.( LSVEC .OR. UCTOL .OR. LSAME( JOBU, 'N' ) ) ) THEN
                    415:          INFO = -2
                    416:       ELSE IF( .NOT.( RSVEC .OR. APPLV .OR. LSAME( JOBV, 'N' ) ) ) THEN
                    417:          INFO = -3
                    418:       ELSE IF( M.LT.0 ) THEN
                    419:          INFO = -4
                    420:       ELSE IF( ( N.LT.0 ) .OR. ( N.GT.M ) ) THEN
                    421:          INFO = -5
                    422:       ELSE IF( LDA.LT.M ) THEN
                    423:          INFO = -7
                    424:       ELSE IF( MV.LT.0 ) THEN
                    425:          INFO = -9
                    426:       ELSE IF( ( RSVEC .AND. ( LDV.LT.N ) ) .OR.
1.6       bertrand  427:      $         ( APPLV .AND. ( LDV.LT.MV ) ) ) THEN
1.1       bertrand  428:          INFO = -11
                    429:       ELSE IF( UCTOL .AND. ( WORK( 1 ).LE.ONE ) ) THEN
                    430:          INFO = -12
                    431:       ELSE IF( LWORK.LT.MAX0( M+N, 6 ) ) THEN
                    432:          INFO = -13
                    433:       ELSE
                    434:          INFO = 0
                    435:       END IF
                    436: *
                    437: *     #:(
                    438:       IF( INFO.NE.0 ) THEN
                    439:          CALL XERBLA( 'DGESVJ', -INFO )
                    440:          RETURN
                    441:       END IF
                    442: *
                    443: * #:) Quick return for void matrix
                    444: *
                    445:       IF( ( M.EQ.0 ) .OR. ( N.EQ.0 ) )RETURN
                    446: *
                    447: *     Set numerical parameters
                    448: *     The stopping criterion for Jacobi rotations is
                    449: *
                    450: *     max_{i<>j}|A(:,i)^T * A(:,j)|/(||A(:,i)||*||A(:,j)||) < CTOL*EPS
                    451: *
                    452: *     where EPS is the round-off and CTOL is defined as follows:
                    453: *
                    454:       IF( UCTOL ) THEN
                    455: *        ... user controlled
                    456:          CTOL = WORK( 1 )
                    457:       ELSE
                    458: *        ... default
                    459:          IF( LSVEC .OR. RSVEC .OR. APPLV ) THEN
                    460:             CTOL = DSQRT( DBLE( M ) )
                    461:          ELSE
                    462:             CTOL = DBLE( M )
                    463:          END IF
                    464:       END IF
                    465: *     ... and the machine dependent parameters are
                    466: *[!]  (Make sure that DLAMCH() works properly on the target machine.)
                    467: *
1.4       bertrand  468:       EPSLN = DLAMCH( 'Epsilon' )
                    469:       ROOTEPS = DSQRT( EPSLN )
1.1       bertrand  470:       SFMIN = DLAMCH( 'SafeMinimum' )
                    471:       ROOTSFMIN = DSQRT( SFMIN )
1.4       bertrand  472:       SMALL = SFMIN / EPSLN
1.1       bertrand  473:       BIG = DLAMCH( 'Overflow' )
                    474: *     BIG         = ONE    / SFMIN
                    475:       ROOTBIG = ONE / ROOTSFMIN
                    476:       LARGE = BIG / DSQRT( DBLE( M*N ) )
                    477:       BIGTHETA = ONE / ROOTEPS
                    478: *
1.4       bertrand  479:       TOL = CTOL*EPSLN
1.1       bertrand  480:       ROOTTOL = DSQRT( TOL )
                    481: *
1.4       bertrand  482:       IF( DBLE( M )*EPSLN.GE.ONE ) THEN
1.6       bertrand  483:          INFO = -4
1.1       bertrand  484:          CALL XERBLA( 'DGESVJ', -INFO )
                    485:          RETURN
                    486:       END IF
                    487: *
                    488: *     Initialize the right singular vector matrix.
                    489: *
                    490:       IF( RSVEC ) THEN
                    491:          MVL = N
                    492:          CALL DLASET( 'A', MVL, N, ZERO, ONE, V, LDV )
                    493:       ELSE IF( APPLV ) THEN
                    494:          MVL = MV
                    495:       END IF
                    496:       RSVEC = RSVEC .OR. APPLV
                    497: *
                    498: *     Initialize SVA( 1:N ) = ( ||A e_i||_2, i = 1:N )
                    499: *(!)  If necessary, scale A to protect the largest singular value
                    500: *     from overflow. It is possible that saving the largest singular
                    501: *     value destroys the information about the small ones.
                    502: *     This initial scaling is almost minimal in the sense that the
                    503: *     goal is to make sure that no column norm overflows, and that
                    504: *     DSQRT(N)*max_i SVA(i) does not overflow. If INFinite entries
                    505: *     in A are detected, the procedure returns with INFO=-6.
                    506: *
1.4       bertrand  507:       SKL= ONE / DSQRT( DBLE( M )*DBLE( N ) )
1.1       bertrand  508:       NOSCALE = .TRUE.
                    509:       GOSCALE = .TRUE.
                    510: *
                    511:       IF( LOWER ) THEN
                    512: *        the input matrix is M-by-N lower triangular (trapezoidal)
                    513:          DO 1874 p = 1, N
                    514:             AAPP = ZERO
1.4       bertrand  515:             AAQQ = ONE
1.1       bertrand  516:             CALL DLASSQ( M-p+1, A( p, p ), 1, AAPP, AAQQ )
                    517:             IF( AAPP.GT.BIG ) THEN
                    518:                INFO = -6
                    519:                CALL XERBLA( 'DGESVJ', -INFO )
                    520:                RETURN
                    521:             END IF
                    522:             AAQQ = DSQRT( AAQQ )
                    523:             IF( ( AAPP.LT.( BIG / AAQQ ) ) .AND. NOSCALE ) THEN
                    524:                SVA( p ) = AAPP*AAQQ
                    525:             ELSE
                    526:                NOSCALE = .FALSE.
1.4       bertrand  527:                SVA( p ) = AAPP*( AAQQ*SKL)
1.1       bertrand  528:                IF( GOSCALE ) THEN
                    529:                   GOSCALE = .FALSE.
                    530:                   DO 1873 q = 1, p - 1
1.4       bertrand  531:                      SVA( q ) = SVA( q )*SKL
1.1       bertrand  532:  1873             CONTINUE
                    533:                END IF
                    534:             END IF
                    535:  1874    CONTINUE
                    536:       ELSE IF( UPPER ) THEN
                    537: *        the input matrix is M-by-N upper triangular (trapezoidal)
                    538:          DO 2874 p = 1, N
                    539:             AAPP = ZERO
1.4       bertrand  540:             AAQQ = ONE
1.1       bertrand  541:             CALL DLASSQ( p, A( 1, p ), 1, AAPP, AAQQ )
                    542:             IF( AAPP.GT.BIG ) THEN
                    543:                INFO = -6
                    544:                CALL XERBLA( 'DGESVJ', -INFO )
                    545:                RETURN
                    546:             END IF
                    547:             AAQQ = DSQRT( AAQQ )
                    548:             IF( ( AAPP.LT.( BIG / AAQQ ) ) .AND. NOSCALE ) THEN
                    549:                SVA( p ) = AAPP*AAQQ
                    550:             ELSE
                    551:                NOSCALE = .FALSE.
1.4       bertrand  552:                SVA( p ) = AAPP*( AAQQ*SKL)
1.1       bertrand  553:                IF( GOSCALE ) THEN
                    554:                   GOSCALE = .FALSE.
                    555:                   DO 2873 q = 1, p - 1
1.4       bertrand  556:                      SVA( q ) = SVA( q )*SKL
1.1       bertrand  557:  2873             CONTINUE
                    558:                END IF
                    559:             END IF
                    560:  2874    CONTINUE
                    561:       ELSE
                    562: *        the input matrix is M-by-N general dense
                    563:          DO 3874 p = 1, N
                    564:             AAPP = ZERO
1.4       bertrand  565:             AAQQ = ONE
1.1       bertrand  566:             CALL DLASSQ( M, A( 1, p ), 1, AAPP, AAQQ )
                    567:             IF( AAPP.GT.BIG ) THEN
                    568:                INFO = -6
                    569:                CALL XERBLA( 'DGESVJ', -INFO )
                    570:                RETURN
                    571:             END IF
                    572:             AAQQ = DSQRT( AAQQ )
                    573:             IF( ( AAPP.LT.( BIG / AAQQ ) ) .AND. NOSCALE ) THEN
                    574:                SVA( p ) = AAPP*AAQQ
                    575:             ELSE
                    576:                NOSCALE = .FALSE.
1.4       bertrand  577:                SVA( p ) = AAPP*( AAQQ*SKL)
1.1       bertrand  578:                IF( GOSCALE ) THEN
                    579:                   GOSCALE = .FALSE.
                    580:                   DO 3873 q = 1, p - 1
1.4       bertrand  581:                      SVA( q ) = SVA( q )*SKL
1.1       bertrand  582:  3873             CONTINUE
                    583:                END IF
                    584:             END IF
                    585:  3874    CONTINUE
                    586:       END IF
                    587: *
1.4       bertrand  588:       IF( NOSCALE )SKL= ONE
1.1       bertrand  589: *
                    590: *     Move the smaller part of the spectrum from the underflow threshold
                    591: *(!)  Start by determining the position of the nonzero entries of the
                    592: *     array SVA() relative to ( SFMIN, BIG ).
                    593: *
                    594:       AAPP = ZERO
                    595:       AAQQ = BIG
                    596:       DO 4781 p = 1, N
                    597:          IF( SVA( p ).NE.ZERO )AAQQ = DMIN1( AAQQ, SVA( p ) )
                    598:          AAPP = DMAX1( AAPP, SVA( p ) )
                    599:  4781 CONTINUE
                    600: *
                    601: * #:) Quick return for zero matrix
                    602: *
                    603:       IF( AAPP.EQ.ZERO ) THEN
                    604:          IF( LSVEC )CALL DLASET( 'G', M, N, ZERO, ONE, A, LDA )
                    605:          WORK( 1 ) = ONE
                    606:          WORK( 2 ) = ZERO
                    607:          WORK( 3 ) = ZERO
                    608:          WORK( 4 ) = ZERO
                    609:          WORK( 5 ) = ZERO
                    610:          WORK( 6 ) = ZERO
                    611:          RETURN
                    612:       END IF
                    613: *
                    614: * #:) Quick return for one-column matrix
                    615: *
                    616:       IF( N.EQ.1 ) THEN
1.4       bertrand  617:          IF( LSVEC )CALL DLASCL( 'G', 0, 0, SVA( 1 ), SKL, M, 1,
1.6       bertrand  618:      $                           A( 1, 1 ), LDA, IERR )
1.4       bertrand  619:          WORK( 1 ) = ONE / SKL
1.1       bertrand  620:          IF( SVA( 1 ).GE.SFMIN ) THEN
                    621:             WORK( 2 ) = ONE
                    622:          ELSE
                    623:             WORK( 2 ) = ZERO
                    624:          END IF
                    625:          WORK( 3 ) = ZERO
                    626:          WORK( 4 ) = ZERO
                    627:          WORK( 5 ) = ZERO
                    628:          WORK( 6 ) = ZERO
                    629:          RETURN
                    630:       END IF
                    631: *
                    632: *     Protect small singular values from underflow, and try to
                    633: *     avoid underflows/overflows in computing Jacobi rotations.
                    634: *
1.4       bertrand  635:       SN = DSQRT( SFMIN / EPSLN )
1.1       bertrand  636:       TEMP1 = DSQRT( BIG / DBLE( N ) )
                    637:       IF( ( AAPP.LE.SN ) .OR. ( AAQQ.GE.TEMP1 ) .OR.
1.6       bertrand  638:      $    ( ( SN.LE.AAQQ ) .AND. ( AAPP.LE.TEMP1 ) ) ) THEN
1.1       bertrand  639:          TEMP1 = DMIN1( BIG, TEMP1 / AAPP )
                    640: *         AAQQ  = AAQQ*TEMP1
                    641: *         AAPP  = AAPP*TEMP1
                    642:       ELSE IF( ( AAQQ.LE.SN ) .AND. ( AAPP.LE.TEMP1 ) ) THEN
                    643:          TEMP1 = DMIN1( SN / AAQQ, BIG / ( AAPP*DSQRT( DBLE( N ) ) ) )
                    644: *         AAQQ  = AAQQ*TEMP1
                    645: *         AAPP  = AAPP*TEMP1
                    646:       ELSE IF( ( AAQQ.GE.SN ) .AND. ( AAPP.GE.TEMP1 ) ) THEN
                    647:          TEMP1 = DMAX1( SN / AAQQ, TEMP1 / AAPP )
                    648: *         AAQQ  = AAQQ*TEMP1
                    649: *         AAPP  = AAPP*TEMP1
                    650:       ELSE IF( ( AAQQ.LE.SN ) .AND. ( AAPP.GE.TEMP1 ) ) THEN
                    651:          TEMP1 = DMIN1( SN / AAQQ, BIG / ( DSQRT( DBLE( N ) )*AAPP ) )
                    652: *         AAQQ  = AAQQ*TEMP1
                    653: *         AAPP  = AAPP*TEMP1
                    654:       ELSE
                    655:          TEMP1 = ONE
                    656:       END IF
                    657: *
                    658: *     Scale, if necessary
                    659: *
                    660:       IF( TEMP1.NE.ONE ) THEN
                    661:          CALL DLASCL( 'G', 0, 0, ONE, TEMP1, N, 1, SVA, N, IERR )
                    662:       END IF
1.4       bertrand  663:       SKL= TEMP1*SKL
                    664:       IF( SKL.NE.ONE ) THEN
                    665:          CALL DLASCL( JOBA, 0, 0, ONE, SKL, M, N, A, LDA, IERR )
                    666:          SKL= ONE / SKL
1.1       bertrand  667:       END IF
                    668: *
                    669: *     Row-cyclic Jacobi SVD algorithm with column pivoting
                    670: *
                    671:       EMPTSW = ( N*( N-1 ) ) / 2
                    672:       NOTROT = 0
                    673:       FASTR( 1 ) = ZERO
                    674: *
                    675: *     A is represented in factored form A = A * diag(WORK), where diag(WORK)
                    676: *     is initialized to identity. WORK is updated during fast scaled
                    677: *     rotations.
                    678: *
                    679:       DO 1868 q = 1, N
                    680:          WORK( q ) = ONE
                    681:  1868 CONTINUE
                    682: *
                    683: *
                    684:       SWBAND = 3
                    685: *[TP] SWBAND is a tuning parameter [TP]. It is meaningful and effective
                    686: *     if DGESVJ is used as a computational routine in the preconditioned
                    687: *     Jacobi SVD algorithm DGESVJ. For sweeps i=1:SWBAND the procedure
                    688: *     works on pivots inside a band-like region around the diagonal.
                    689: *     The boundaries are determined dynamically, based on the number of
                    690: *     pivots above a threshold.
                    691: *
                    692:       KBL = MIN0( 8, N )
                    693: *[TP] KBL is a tuning parameter that defines the tile size in the
                    694: *     tiling of the p-q loops of pivot pairs. In general, an optimal
                    695: *     value of KBL depends on the matrix dimensions and on the
                    696: *     parameters of the computer's memory.
                    697: *
                    698:       NBL = N / KBL
                    699:       IF( ( NBL*KBL ).NE.N )NBL = NBL + 1
                    700: *
                    701:       BLSKIP = KBL**2
                    702: *[TP] BLKSKIP is a tuning parameter that depends on SWBAND and KBL.
                    703: *
                    704:       ROWSKIP = MIN0( 5, KBL )
                    705: *[TP] ROWSKIP is a tuning parameter.
                    706: *
                    707:       LKAHEAD = 1
                    708: *[TP] LKAHEAD is a tuning parameter.
                    709: *
                    710: *     Quasi block transformations, using the lower (upper) triangular
                    711: *     structure of the input matrix. The quasi-block-cycling usually
                    712: *     invokes cubic convergence. Big part of this cycle is done inside
                    713: *     canonical subspaces of dimensions less than M.
                    714: *
                    715:       IF( ( LOWER .OR. UPPER ) .AND. ( N.GT.MAX0( 64, 4*KBL ) ) ) THEN
                    716: *[TP] The number of partition levels and the actual partition are
                    717: *     tuning parameters.
                    718:          N4 = N / 4
                    719:          N2 = N / 2
                    720:          N34 = 3*N4
                    721:          IF( APPLV ) THEN
                    722:             q = 0
                    723:          ELSE
                    724:             q = 1
                    725:          END IF
                    726: *
                    727:          IF( LOWER ) THEN
                    728: *
                    729: *     This works very well on lower triangular matrices, in particular
                    730: *     in the framework of the preconditioned Jacobi SVD (xGEJSV).
                    731: *     The idea is simple:
                    732: *     [+ 0 0 0]   Note that Jacobi transformations of [0 0]
                    733: *     [+ + 0 0]                                       [0 0]
                    734: *     [+ + x 0]   actually work on [x 0]              [x 0]
                    735: *     [+ + x x]                    [x x].             [x x]
                    736: *
                    737:             CALL DGSVJ0( JOBV, M-N34, N-N34, A( N34+1, N34+1 ), LDA,
1.6       bertrand  738:      $                   WORK( N34+1 ), SVA( N34+1 ), MVL,
                    739:      $                   V( N34*q+1, N34+1 ), LDV, EPSLN, SFMIN, TOL,
                    740:      $                   2, WORK( N+1 ), LWORK-N, IERR )
1.1       bertrand  741: *
                    742:             CALL DGSVJ0( JOBV, M-N2, N34-N2, A( N2+1, N2+1 ), LDA,
1.6       bertrand  743:      $                   WORK( N2+1 ), SVA( N2+1 ), MVL,
                    744:      $                   V( N2*q+1, N2+1 ), LDV, EPSLN, SFMIN, TOL, 2,
                    745:      $                   WORK( N+1 ), LWORK-N, IERR )
1.1       bertrand  746: *
                    747:             CALL DGSVJ1( JOBV, M-N2, N-N2, N4, A( N2+1, N2+1 ), LDA,
1.6       bertrand  748:      $                   WORK( N2+1 ), SVA( N2+1 ), MVL,
                    749:      $                   V( N2*q+1, N2+1 ), LDV, EPSLN, SFMIN, TOL, 1,
                    750:      $                   WORK( N+1 ), LWORK-N, IERR )
1.1       bertrand  751: *
                    752:             CALL DGSVJ0( JOBV, M-N4, N2-N4, A( N4+1, N4+1 ), LDA,
1.6       bertrand  753:      $                   WORK( N4+1 ), SVA( N4+1 ), MVL,
                    754:      $                   V( N4*q+1, N4+1 ), LDV, EPSLN, SFMIN, TOL, 1,
                    755:      $                   WORK( N+1 ), LWORK-N, IERR )
1.1       bertrand  756: *
                    757:             CALL DGSVJ0( JOBV, M, N4, A, LDA, WORK, SVA, MVL, V, LDV,
1.6       bertrand  758:      $                   EPSLN, SFMIN, TOL, 1, WORK( N+1 ), LWORK-N,
                    759:      $                   IERR )
1.1       bertrand  760: *
                    761:             CALL DGSVJ1( JOBV, M, N2, N4, A, LDA, WORK, SVA, MVL, V,
1.6       bertrand  762:      $                   LDV, EPSLN, SFMIN, TOL, 1, WORK( N+1 ),
                    763:      $                   LWORK-N, IERR )
1.1       bertrand  764: *
                    765: *
                    766:          ELSE IF( UPPER ) THEN
                    767: *
                    768: *
                    769:             CALL DGSVJ0( JOBV, N4, N4, A, LDA, WORK, SVA, MVL, V, LDV,
1.6       bertrand  770:      $                   EPSLN, SFMIN, TOL, 2, WORK( N+1 ), LWORK-N,
                    771:      $                   IERR )
1.1       bertrand  772: *
                    773:             CALL DGSVJ0( JOBV, N2, N4, A( 1, N4+1 ), LDA, WORK( N4+1 ),
1.6       bertrand  774:      $                   SVA( N4+1 ), MVL, V( N4*q+1, N4+1 ), LDV,
                    775:      $                   EPSLN, SFMIN, TOL, 1, WORK( N+1 ), LWORK-N,
                    776:      $                   IERR )
1.1       bertrand  777: *
                    778:             CALL DGSVJ1( JOBV, N2, N2, N4, A, LDA, WORK, SVA, MVL, V,
1.6       bertrand  779:      $                   LDV, EPSLN, SFMIN, TOL, 1, WORK( N+1 ),
                    780:      $                   LWORK-N, IERR )
1.1       bertrand  781: *
                    782:             CALL DGSVJ0( JOBV, N2+N4, N4, A( 1, N2+1 ), LDA,
1.6       bertrand  783:      $                   WORK( N2+1 ), SVA( N2+1 ), MVL,
                    784:      $                   V( N2*q+1, N2+1 ), LDV, EPSLN, SFMIN, TOL, 1,
                    785:      $                   WORK( N+1 ), LWORK-N, IERR )
1.1       bertrand  786: 
                    787:          END IF
                    788: *
                    789:       END IF
                    790: *
                    791: *     .. Row-cyclic pivot strategy with de Rijk's pivoting ..
                    792: *
                    793:       DO 1993 i = 1, NSWEEP
                    794: *
                    795: *     .. go go go ...
                    796: *
                    797:          MXAAPQ = ZERO
                    798:          MXSINJ = ZERO
                    799:          ISWROT = 0
                    800: *
                    801:          NOTROT = 0
                    802:          PSKIPPED = 0
                    803: *
                    804: *     Each sweep is unrolled using KBL-by-KBL tiles over the pivot pairs
                    805: *     1 <= p < q <= N. This is the first step toward a blocked implementation
                    806: *     of the rotations. New implementation, based on block transformations,
                    807: *     is under development.
                    808: *
                    809:          DO 2000 ibr = 1, NBL
                    810: *
                    811:             igl = ( ibr-1 )*KBL + 1
                    812: *
                    813:             DO 1002 ir1 = 0, MIN0( LKAHEAD, NBL-ibr )
                    814: *
                    815:                igl = igl + ir1*KBL
                    816: *
                    817:                DO 2001 p = igl, MIN0( igl+KBL-1, N-1 )
                    818: *
                    819: *     .. de Rijk's pivoting
                    820: *
                    821:                   q = IDAMAX( N-p+1, SVA( p ), 1 ) + p - 1
                    822:                   IF( p.NE.q ) THEN
                    823:                      CALL DSWAP( M, A( 1, p ), 1, A( 1, q ), 1 )
                    824:                      IF( RSVEC )CALL DSWAP( MVL, V( 1, p ), 1,
1.6       bertrand  825:      $                                      V( 1, q ), 1 )
1.1       bertrand  826:                      TEMP1 = SVA( p )
                    827:                      SVA( p ) = SVA( q )
                    828:                      SVA( q ) = TEMP1
                    829:                      TEMP1 = WORK( p )
                    830:                      WORK( p ) = WORK( q )
                    831:                      WORK( q ) = TEMP1
                    832:                   END IF
                    833: *
                    834:                   IF( ir1.EQ.0 ) THEN
                    835: *
                    836: *        Column norms are periodically updated by explicit
                    837: *        norm computation.
                    838: *        Caveat:
                    839: *        Unfortunately, some BLAS implementations compute DNRM2(M,A(1,p),1)
                    840: *        as DSQRT(DDOT(M,A(1,p),1,A(1,p),1)), which may cause the result to
                    841: *        overflow for ||A(:,p)||_2 > DSQRT(overflow_threshold), and to
                    842: *        underflow for ||A(:,p)||_2 < DSQRT(underflow_threshold).
                    843: *        Hence, DNRM2 cannot be trusted, not even in the case when
                    844: *        the true norm is far from the under(over)flow boundaries.
                    845: *        If properly implemented DNRM2 is available, the IF-THEN-ELSE
                    846: *        below should read "AAPP = DNRM2( M, A(1,p), 1 ) * WORK(p)".
                    847: *
                    848:                      IF( ( SVA( p ).LT.ROOTBIG ) .AND.
1.6       bertrand  849:      $                   ( SVA( p ).GT.ROOTSFMIN ) ) THEN
1.1       bertrand  850:                         SVA( p ) = DNRM2( M, A( 1, p ), 1 )*WORK( p )
                    851:                      ELSE
                    852:                         TEMP1 = ZERO
1.4       bertrand  853:                         AAPP = ONE
1.1       bertrand  854:                         CALL DLASSQ( M, A( 1, p ), 1, TEMP1, AAPP )
                    855:                         SVA( p ) = TEMP1*DSQRT( AAPP )*WORK( p )
                    856:                      END IF
                    857:                      AAPP = SVA( p )
                    858:                   ELSE
                    859:                      AAPP = SVA( p )
                    860:                   END IF
                    861: *
                    862:                   IF( AAPP.GT.ZERO ) THEN
                    863: *
                    864:                      PSKIPPED = 0
                    865: *
                    866:                      DO 2002 q = p + 1, MIN0( igl+KBL-1, N )
                    867: *
                    868:                         AAQQ = SVA( q )
                    869: *
                    870:                         IF( AAQQ.GT.ZERO ) THEN
                    871: *
                    872:                            AAPP0 = AAPP
                    873:                            IF( AAQQ.GE.ONE ) THEN
                    874:                               ROTOK = ( SMALL*AAPP ).LE.AAQQ
                    875:                               IF( AAPP.LT.( BIG / AAQQ ) ) THEN
                    876:                                  AAPQ = ( DDOT( M, A( 1, p ), 1, A( 1,
1.6       bertrand  877:      $                                  q ), 1 )*WORK( p )*WORK( q ) /
                    878:      $                                  AAQQ ) / AAPP
1.1       bertrand  879:                               ELSE
                    880:                                  CALL DCOPY( M, A( 1, p ), 1,
1.6       bertrand  881:      $                                       WORK( N+1 ), 1 )
1.1       bertrand  882:                                  CALL DLASCL( 'G', 0, 0, AAPP,
1.6       bertrand  883:      $                                        WORK( p ), M, 1,
                    884:      $                                        WORK( N+1 ), LDA, IERR )
1.1       bertrand  885:                                  AAPQ = DDOT( M, WORK( N+1 ), 1,
1.6       bertrand  886:      $                                  A( 1, q ), 1 )*WORK( q ) / AAQQ
1.1       bertrand  887:                               END IF
                    888:                            ELSE
                    889:                               ROTOK = AAPP.LE.( AAQQ / SMALL )
                    890:                               IF( AAPP.GT.( SMALL / AAQQ ) ) THEN
                    891:                                  AAPQ = ( DDOT( M, A( 1, p ), 1, A( 1,
1.6       bertrand  892:      $                                  q ), 1 )*WORK( p )*WORK( q ) /
                    893:      $                                  AAQQ ) / AAPP
1.1       bertrand  894:                               ELSE
                    895:                                  CALL DCOPY( M, A( 1, q ), 1,
1.6       bertrand  896:      $                                       WORK( N+1 ), 1 )
1.1       bertrand  897:                                  CALL DLASCL( 'G', 0, 0, AAQQ,
1.6       bertrand  898:      $                                        WORK( q ), M, 1,
                    899:      $                                        WORK( N+1 ), LDA, IERR )
1.1       bertrand  900:                                  AAPQ = DDOT( M, WORK( N+1 ), 1,
1.6       bertrand  901:      $                                  A( 1, p ), 1 )*WORK( p ) / AAPP
1.1       bertrand  902:                               END IF
                    903:                            END IF
                    904: *
                    905:                            MXAAPQ = DMAX1( MXAAPQ, DABS( AAPQ ) )
                    906: *
                    907: *        TO rotate or NOT to rotate, THAT is the question ...
                    908: *
                    909:                            IF( DABS( AAPQ ).GT.TOL ) THEN
                    910: *
                    911: *           .. rotate
                    912: *[RTD]      ROTATED = ROTATED + ONE
                    913: *
                    914:                               IF( ir1.EQ.0 ) THEN
                    915:                                  NOTROT = 0
                    916:                                  PSKIPPED = 0
                    917:                                  ISWROT = ISWROT + 1
                    918:                               END IF
                    919: *
                    920:                               IF( ROTOK ) THEN
                    921: *
                    922:                                  AQOAP = AAQQ / AAPP
                    923:                                  APOAQ = AAPP / AAQQ
1.6       bertrand  924:                                  THETA = -HALF*DABS(AQOAP-APOAQ)/AAPQ
1.1       bertrand  925: *
                    926:                                  IF( DABS( THETA ).GT.BIGTHETA ) THEN
                    927: *
                    928:                                     T = HALF / THETA
                    929:                                     FASTR( 3 ) = T*WORK( p ) / WORK( q )
                    930:                                     FASTR( 4 ) = -T*WORK( q ) /
1.6       bertrand  931:      $                                           WORK( p )
1.1       bertrand  932:                                     CALL DROTM( M, A( 1, p ), 1,
1.6       bertrand  933:      $                                          A( 1, q ), 1, FASTR )
1.1       bertrand  934:                                     IF( RSVEC )CALL DROTM( MVL,
1.6       bertrand  935:      $                                              V( 1, p ), 1,
                    936:      $                                              V( 1, q ), 1,
                    937:      $                                              FASTR )
1.1       bertrand  938:                                     SVA( q ) = AAQQ*DSQRT( DMAX1( ZERO,
1.6       bertrand  939:      $                                         ONE+T*APOAQ*AAPQ ) )
1.4       bertrand  940:                                     AAPP = AAPP*DSQRT( DMAX1( ZERO,
1.6       bertrand  941:      $                                     ONE-T*AQOAP*AAPQ ) )
1.1       bertrand  942:                                     MXSINJ = DMAX1( MXSINJ, DABS( T ) )
                    943: *
                    944:                                  ELSE
                    945: *
                    946: *                 .. choose correct signum for THETA and rotate
                    947: *
                    948:                                     THSIGN = -DSIGN( ONE, AAPQ )
                    949:                                     T = ONE / ( THETA+THSIGN*
1.6       bertrand  950:      $                                  DSQRT( ONE+THETA*THETA ) )
1.1       bertrand  951:                                     CS = DSQRT( ONE / ( ONE+T*T ) )
                    952:                                     SN = T*CS
                    953: *
                    954:                                     MXSINJ = DMAX1( MXSINJ, DABS( SN ) )
                    955:                                     SVA( q ) = AAQQ*DSQRT( DMAX1( ZERO,
1.6       bertrand  956:      $                                         ONE+T*APOAQ*AAPQ ) )
1.1       bertrand  957:                                     AAPP = AAPP*DSQRT( DMAX1( ZERO,
1.6       bertrand  958:      $                                     ONE-T*AQOAP*AAPQ ) )
1.1       bertrand  959: *
                    960:                                     APOAQ = WORK( p ) / WORK( q )
                    961:                                     AQOAP = WORK( q ) / WORK( p )
                    962:                                     IF( WORK( p ).GE.ONE ) THEN
                    963:                                        IF( WORK( q ).GE.ONE ) THEN
                    964:                                           FASTR( 3 ) = T*APOAQ
                    965:                                           FASTR( 4 ) = -T*AQOAP
                    966:                                           WORK( p ) = WORK( p )*CS
                    967:                                           WORK( q ) = WORK( q )*CS
                    968:                                           CALL DROTM( M, A( 1, p ), 1,
1.6       bertrand  969:      $                                                A( 1, q ), 1,
                    970:      $                                                FASTR )
1.1       bertrand  971:                                           IF( RSVEC )CALL DROTM( MVL,
1.6       bertrand  972:      $                                        V( 1, p ), 1, V( 1, q ),
                    973:      $                                        1, FASTR )
1.1       bertrand  974:                                        ELSE
                    975:                                           CALL DAXPY( M, -T*AQOAP,
1.6       bertrand  976:      $                                                A( 1, q ), 1,
                    977:      $                                                A( 1, p ), 1 )
1.1       bertrand  978:                                           CALL DAXPY( M, CS*SN*APOAQ,
1.6       bertrand  979:      $                                                A( 1, p ), 1,
                    980:      $                                                A( 1, q ), 1 )
1.1       bertrand  981:                                           WORK( p ) = WORK( p )*CS
                    982:                                           WORK( q ) = WORK( q ) / CS
                    983:                                           IF( RSVEC ) THEN
                    984:                                              CALL DAXPY( MVL, -T*AQOAP,
1.6       bertrand  985:      $                                                   V( 1, q ), 1,
                    986:      $                                                   V( 1, p ), 1 )
1.1       bertrand  987:                                              CALL DAXPY( MVL,
1.6       bertrand  988:      $                                                   CS*SN*APOAQ,
                    989:      $                                                   V( 1, p ), 1,
                    990:      $                                                   V( 1, q ), 1 )
1.1       bertrand  991:                                           END IF
                    992:                                        END IF
                    993:                                     ELSE
                    994:                                        IF( WORK( q ).GE.ONE ) THEN
                    995:                                           CALL DAXPY( M, T*APOAQ,
1.6       bertrand  996:      $                                                A( 1, p ), 1,
                    997:      $                                                A( 1, q ), 1 )
1.1       bertrand  998:                                           CALL DAXPY( M, -CS*SN*AQOAP,
1.6       bertrand  999:      $                                                A( 1, q ), 1,
                   1000:      $                                                A( 1, p ), 1 )
1.1       bertrand 1001:                                           WORK( p ) = WORK( p ) / CS
                   1002:                                           WORK( q ) = WORK( q )*CS
                   1003:                                           IF( RSVEC ) THEN
                   1004:                                              CALL DAXPY( MVL, T*APOAQ,
1.6       bertrand 1005:      $                                                   V( 1, p ), 1,
                   1006:      $                                                   V( 1, q ), 1 )
1.1       bertrand 1007:                                              CALL DAXPY( MVL,
1.6       bertrand 1008:      $                                                   -CS*SN*AQOAP,
                   1009:      $                                                   V( 1, q ), 1,
                   1010:      $                                                   V( 1, p ), 1 )
1.1       bertrand 1011:                                           END IF
                   1012:                                        ELSE
                   1013:                                           IF( WORK( p ).GE.WORK( q ) )
1.6       bertrand 1014:      $                                        THEN
1.1       bertrand 1015:                                              CALL DAXPY( M, -T*AQOAP,
1.6       bertrand 1016:      $                                                   A( 1, q ), 1,
                   1017:      $                                                   A( 1, p ), 1 )
1.1       bertrand 1018:                                              CALL DAXPY( M, CS*SN*APOAQ,
1.6       bertrand 1019:      $                                                   A( 1, p ), 1,
                   1020:      $                                                   A( 1, q ), 1 )
1.1       bertrand 1021:                                              WORK( p ) = WORK( p )*CS
                   1022:                                              WORK( q ) = WORK( q ) / CS
                   1023:                                              IF( RSVEC ) THEN
                   1024:                                                 CALL DAXPY( MVL,
1.6       bertrand 1025:      $                                               -T*AQOAP,
                   1026:      $                                               V( 1, q ), 1,
                   1027:      $                                               V( 1, p ), 1 )
1.1       bertrand 1028:                                                 CALL DAXPY( MVL,
1.6       bertrand 1029:      $                                               CS*SN*APOAQ,
                   1030:      $                                               V( 1, p ), 1,
                   1031:      $                                               V( 1, q ), 1 )
1.1       bertrand 1032:                                              END IF
                   1033:                                           ELSE
                   1034:                                              CALL DAXPY( M, T*APOAQ,
1.6       bertrand 1035:      $                                                   A( 1, p ), 1,
                   1036:      $                                                   A( 1, q ), 1 )
1.1       bertrand 1037:                                              CALL DAXPY( M,
1.6       bertrand 1038:      $                                                   -CS*SN*AQOAP,
                   1039:      $                                                   A( 1, q ), 1,
                   1040:      $                                                   A( 1, p ), 1 )
1.1       bertrand 1041:                                              WORK( p ) = WORK( p ) / CS
                   1042:                                              WORK( q ) = WORK( q )*CS
                   1043:                                              IF( RSVEC ) THEN
                   1044:                                                 CALL DAXPY( MVL,
1.6       bertrand 1045:      $                                               T*APOAQ, V( 1, p ),
                   1046:      $                                               1, V( 1, q ), 1 )
1.1       bertrand 1047:                                                 CALL DAXPY( MVL,
1.6       bertrand 1048:      $                                               -CS*SN*AQOAP,
                   1049:      $                                               V( 1, q ), 1,
                   1050:      $                                               V( 1, p ), 1 )
1.1       bertrand 1051:                                              END IF
                   1052:                                           END IF
                   1053:                                        END IF
                   1054:                                     END IF
                   1055:                                  END IF
                   1056: *
                   1057:                               ELSE
                   1058: *              .. have to use modified Gram-Schmidt like transformation
                   1059:                                  CALL DCOPY( M, A( 1, p ), 1,
1.6       bertrand 1060:      $                                       WORK( N+1 ), 1 )
1.1       bertrand 1061:                                  CALL DLASCL( 'G', 0, 0, AAPP, ONE, M,
1.6       bertrand 1062:      $                                        1, WORK( N+1 ), LDA,
                   1063:      $                                        IERR )
1.1       bertrand 1064:                                  CALL DLASCL( 'G', 0, 0, AAQQ, ONE, M,
1.6       bertrand 1065:      $                                        1, A( 1, q ), LDA, IERR )
1.1       bertrand 1066:                                  TEMP1 = -AAPQ*WORK( p ) / WORK( q )
                   1067:                                  CALL DAXPY( M, TEMP1, WORK( N+1 ), 1,
1.6       bertrand 1068:      $                                       A( 1, q ), 1 )
1.1       bertrand 1069:                                  CALL DLASCL( 'G', 0, 0, ONE, AAQQ, M,
1.6       bertrand 1070:      $                                        1, A( 1, q ), LDA, IERR )
1.1       bertrand 1071:                                  SVA( q ) = AAQQ*DSQRT( DMAX1( ZERO,
1.6       bertrand 1072:      $                                      ONE-AAPQ*AAPQ ) )
1.1       bertrand 1073:                                  MXSINJ = DMAX1( MXSINJ, SFMIN )
                   1074:                               END IF
                   1075: *           END IF ROTOK THEN ... ELSE
                   1076: *
                   1077: *           In the case of cancellation in updating SVA(q), SVA(p)
                   1078: *           recompute SVA(q), SVA(p).
                   1079: *
                   1080:                               IF( ( SVA( q ) / AAQQ )**2.LE.ROOTEPS )
1.6       bertrand 1081:      $                            THEN
1.1       bertrand 1082:                                  IF( ( AAQQ.LT.ROOTBIG ) .AND.
1.6       bertrand 1083:      $                               ( AAQQ.GT.ROOTSFMIN ) ) THEN
1.1       bertrand 1084:                                     SVA( q ) = DNRM2( M, A( 1, q ), 1 )*
1.6       bertrand 1085:      $                                         WORK( q )
1.1       bertrand 1086:                                  ELSE
                   1087:                                     T = ZERO
1.4       bertrand 1088:                                     AAQQ = ONE
1.1       bertrand 1089:                                     CALL DLASSQ( M, A( 1, q ), 1, T,
1.6       bertrand 1090:      $                                           AAQQ )
1.1       bertrand 1091:                                     SVA( q ) = T*DSQRT( AAQQ )*WORK( q )
                   1092:                                  END IF
                   1093:                               END IF
                   1094:                               IF( ( AAPP / AAPP0 ).LE.ROOTEPS ) THEN
                   1095:                                  IF( ( AAPP.LT.ROOTBIG ) .AND.
1.6       bertrand 1096:      $                               ( AAPP.GT.ROOTSFMIN ) ) THEN
1.1       bertrand 1097:                                     AAPP = DNRM2( M, A( 1, p ), 1 )*
1.6       bertrand 1098:      $                                     WORK( p )
1.1       bertrand 1099:                                  ELSE
                   1100:                                     T = ZERO
1.4       bertrand 1101:                                     AAPP = ONE
1.1       bertrand 1102:                                     CALL DLASSQ( M, A( 1, p ), 1, T,
1.6       bertrand 1103:      $                                           AAPP )
1.1       bertrand 1104:                                     AAPP = T*DSQRT( AAPP )*WORK( p )
                   1105:                                  END IF
                   1106:                                  SVA( p ) = AAPP
                   1107:                               END IF
                   1108: *
                   1109:                            ELSE
                   1110: *        A(:,p) and A(:,q) already numerically orthogonal
                   1111:                               IF( ir1.EQ.0 )NOTROT = NOTROT + 1
                   1112: *[RTD]      SKIPPED  = SKIPPED  + 1
                   1113:                               PSKIPPED = PSKIPPED + 1
                   1114:                            END IF
                   1115:                         ELSE
                   1116: *        A(:,q) is zero column
                   1117:                            IF( ir1.EQ.0 )NOTROT = NOTROT + 1
                   1118:                            PSKIPPED = PSKIPPED + 1
                   1119:                         END IF
                   1120: *
                   1121:                         IF( ( i.LE.SWBAND ) .AND.
1.6       bertrand 1122:      $                      ( PSKIPPED.GT.ROWSKIP ) ) THEN
1.1       bertrand 1123:                            IF( ir1.EQ.0 )AAPP = -AAPP
                   1124:                            NOTROT = 0
                   1125:                            GO TO 2103
                   1126:                         END IF
                   1127: *
                   1128:  2002                CONTINUE
                   1129: *     END q-LOOP
                   1130: *
                   1131:  2103                CONTINUE
                   1132: *     bailed out of q-loop
                   1133: *
                   1134:                      SVA( p ) = AAPP
                   1135: *
                   1136:                   ELSE
                   1137:                      SVA( p ) = AAPP
                   1138:                      IF( ( ir1.EQ.0 ) .AND. ( AAPP.EQ.ZERO ) )
1.6       bertrand 1139:      $                   NOTROT = NOTROT + MIN0( igl+KBL-1, N ) - p
1.1       bertrand 1140:                   END IF
                   1141: *
                   1142:  2001          CONTINUE
                   1143: *     end of the p-loop
                   1144: *     end of doing the block ( ibr, ibr )
                   1145:  1002       CONTINUE
                   1146: *     end of ir1-loop
                   1147: *
                   1148: * ... go to the off diagonal blocks
                   1149: *
                   1150:             igl = ( ibr-1 )*KBL + 1
                   1151: *
                   1152:             DO 2010 jbc = ibr + 1, NBL
                   1153: *
                   1154:                jgl = ( jbc-1 )*KBL + 1
                   1155: *
                   1156: *        doing the block at ( ibr, jbc )
                   1157: *
                   1158:                IJBLSK = 0
                   1159:                DO 2100 p = igl, MIN0( igl+KBL-1, N )
                   1160: *
                   1161:                   AAPP = SVA( p )
                   1162:                   IF( AAPP.GT.ZERO ) THEN
                   1163: *
                   1164:                      PSKIPPED = 0
                   1165: *
                   1166:                      DO 2200 q = jgl, MIN0( jgl+KBL-1, N )
                   1167: *
                   1168:                         AAQQ = SVA( q )
                   1169:                         IF( AAQQ.GT.ZERO ) THEN
                   1170:                            AAPP0 = AAPP
                   1171: *
                   1172: *     .. M x 2 Jacobi SVD ..
                   1173: *
                   1174: *        Safe Gram matrix computation
                   1175: *
                   1176:                            IF( AAQQ.GE.ONE ) THEN
                   1177:                               IF( AAPP.GE.AAQQ ) THEN
                   1178:                                  ROTOK = ( SMALL*AAPP ).LE.AAQQ
                   1179:                               ELSE
                   1180:                                  ROTOK = ( SMALL*AAQQ ).LE.AAPP
                   1181:                               END IF
                   1182:                               IF( AAPP.LT.( BIG / AAQQ ) ) THEN
                   1183:                                  AAPQ = ( DDOT( M, A( 1, p ), 1, A( 1,
1.6       bertrand 1184:      $                                  q ), 1 )*WORK( p )*WORK( q ) /
                   1185:      $                                  AAQQ ) / AAPP
1.1       bertrand 1186:                               ELSE
                   1187:                                  CALL DCOPY( M, A( 1, p ), 1,
1.6       bertrand 1188:      $                                       WORK( N+1 ), 1 )
1.1       bertrand 1189:                                  CALL DLASCL( 'G', 0, 0, AAPP,
1.6       bertrand 1190:      $                                        WORK( p ), M, 1,
                   1191:      $                                        WORK( N+1 ), LDA, IERR )
1.1       bertrand 1192:                                  AAPQ = DDOT( M, WORK( N+1 ), 1,
1.6       bertrand 1193:      $                                  A( 1, q ), 1 )*WORK( q ) / AAQQ
1.1       bertrand 1194:                               END IF
                   1195:                            ELSE
                   1196:                               IF( AAPP.GE.AAQQ ) THEN
                   1197:                                  ROTOK = AAPP.LE.( AAQQ / SMALL )
                   1198:                               ELSE
                   1199:                                  ROTOK = AAQQ.LE.( AAPP / SMALL )
                   1200:                               END IF
                   1201:                               IF( AAPP.GT.( SMALL / AAQQ ) ) THEN
                   1202:                                  AAPQ = ( DDOT( M, A( 1, p ), 1, A( 1,
1.6       bertrand 1203:      $                                  q ), 1 )*WORK( p )*WORK( q ) /
                   1204:      $                                  AAQQ ) / AAPP
1.1       bertrand 1205:                               ELSE
                   1206:                                  CALL DCOPY( M, A( 1, q ), 1,
1.6       bertrand 1207:      $                                       WORK( N+1 ), 1 )
1.1       bertrand 1208:                                  CALL DLASCL( 'G', 0, 0, AAQQ,
1.6       bertrand 1209:      $                                        WORK( q ), M, 1,
                   1210:      $                                        WORK( N+1 ), LDA, IERR )
1.1       bertrand 1211:                                  AAPQ = DDOT( M, WORK( N+1 ), 1,
1.6       bertrand 1212:      $                                  A( 1, p ), 1 )*WORK( p ) / AAPP
1.1       bertrand 1213:                               END IF
                   1214:                            END IF
                   1215: *
                   1216:                            MXAAPQ = DMAX1( MXAAPQ, DABS( AAPQ ) )
                   1217: *
                   1218: *        TO rotate or NOT to rotate, THAT is the question ...
                   1219: *
                   1220:                            IF( DABS( AAPQ ).GT.TOL ) THEN
                   1221:                               NOTROT = 0
                   1222: *[RTD]      ROTATED  = ROTATED + 1
                   1223:                               PSKIPPED = 0
                   1224:                               ISWROT = ISWROT + 1
                   1225: *
                   1226:                               IF( ROTOK ) THEN
                   1227: *
                   1228:                                  AQOAP = AAQQ / AAPP
                   1229:                                  APOAQ = AAPP / AAQQ
1.6       bertrand 1230:                                  THETA = -HALF*DABS(AQOAP-APOAQ)/AAPQ
1.1       bertrand 1231:                                  IF( AAQQ.GT.AAPP0 )THETA = -THETA
                   1232: *
                   1233:                                  IF( DABS( THETA ).GT.BIGTHETA ) THEN
                   1234:                                     T = HALF / THETA
                   1235:                                     FASTR( 3 ) = T*WORK( p ) / WORK( q )
                   1236:                                     FASTR( 4 ) = -T*WORK( q ) /
1.6       bertrand 1237:      $                                           WORK( p )
1.1       bertrand 1238:                                     CALL DROTM( M, A( 1, p ), 1,
1.6       bertrand 1239:      $                                          A( 1, q ), 1, FASTR )
1.1       bertrand 1240:                                     IF( RSVEC )CALL DROTM( MVL,
1.6       bertrand 1241:      $                                              V( 1, p ), 1,
                   1242:      $                                              V( 1, q ), 1,
                   1243:      $                                              FASTR )
1.1       bertrand 1244:                                     SVA( q ) = AAQQ*DSQRT( DMAX1( ZERO,
1.6       bertrand 1245:      $                                         ONE+T*APOAQ*AAPQ ) )
1.1       bertrand 1246:                                     AAPP = AAPP*DSQRT( DMAX1( ZERO,
1.6       bertrand 1247:      $                                     ONE-T*AQOAP*AAPQ ) )
1.1       bertrand 1248:                                     MXSINJ = DMAX1( MXSINJ, DABS( T ) )
                   1249:                                  ELSE
                   1250: *
                   1251: *                 .. choose correct signum for THETA and rotate
                   1252: *
                   1253:                                     THSIGN = -DSIGN( ONE, AAPQ )
                   1254:                                     IF( AAQQ.GT.AAPP0 )THSIGN = -THSIGN
                   1255:                                     T = ONE / ( THETA+THSIGN*
1.6       bertrand 1256:      $                                  DSQRT( ONE+THETA*THETA ) )
1.1       bertrand 1257:                                     CS = DSQRT( ONE / ( ONE+T*T ) )
                   1258:                                     SN = T*CS
                   1259:                                     MXSINJ = DMAX1( MXSINJ, DABS( SN ) )
                   1260:                                     SVA( q ) = AAQQ*DSQRT( DMAX1( ZERO,
1.6       bertrand 1261:      $                                         ONE+T*APOAQ*AAPQ ) )
1.4       bertrand 1262:                                     AAPP = AAPP*DSQRT( DMAX1( ZERO, 
1.6       bertrand 1263:      $                                     ONE-T*AQOAP*AAPQ ) )
1.1       bertrand 1264: *
                   1265:                                     APOAQ = WORK( p ) / WORK( q )
                   1266:                                     AQOAP = WORK( q ) / WORK( p )
                   1267:                                     IF( WORK( p ).GE.ONE ) THEN
                   1268: *
                   1269:                                        IF( WORK( q ).GE.ONE ) THEN
                   1270:                                           FASTR( 3 ) = T*APOAQ
                   1271:                                           FASTR( 4 ) = -T*AQOAP
                   1272:                                           WORK( p ) = WORK( p )*CS
                   1273:                                           WORK( q ) = WORK( q )*CS
                   1274:                                           CALL DROTM( M, A( 1, p ), 1,
1.6       bertrand 1275:      $                                                A( 1, q ), 1,
                   1276:      $                                                FASTR )
1.1       bertrand 1277:                                           IF( RSVEC )CALL DROTM( MVL,
1.6       bertrand 1278:      $                                        V( 1, p ), 1, V( 1, q ),
                   1279:      $                                        1, FASTR )
1.1       bertrand 1280:                                        ELSE
                   1281:                                           CALL DAXPY( M, -T*AQOAP,
1.6       bertrand 1282:      $                                                A( 1, q ), 1,
                   1283:      $                                                A( 1, p ), 1 )
1.1       bertrand 1284:                                           CALL DAXPY( M, CS*SN*APOAQ,
1.6       bertrand 1285:      $                                                A( 1, p ), 1,
                   1286:      $                                                A( 1, q ), 1 )
1.1       bertrand 1287:                                           IF( RSVEC ) THEN
                   1288:                                              CALL DAXPY( MVL, -T*AQOAP,
1.6       bertrand 1289:      $                                                   V( 1, q ), 1,
                   1290:      $                                                   V( 1, p ), 1 )
1.1       bertrand 1291:                                              CALL DAXPY( MVL,
1.6       bertrand 1292:      $                                                   CS*SN*APOAQ,
                   1293:      $                                                   V( 1, p ), 1,
                   1294:      $                                                   V( 1, q ), 1 )
1.1       bertrand 1295:                                           END IF
                   1296:                                           WORK( p ) = WORK( p )*CS
                   1297:                                           WORK( q ) = WORK( q ) / CS
                   1298:                                        END IF
                   1299:                                     ELSE
                   1300:                                        IF( WORK( q ).GE.ONE ) THEN
                   1301:                                           CALL DAXPY( M, T*APOAQ,
1.6       bertrand 1302:      $                                                A( 1, p ), 1,
                   1303:      $                                                A( 1, q ), 1 )
1.1       bertrand 1304:                                           CALL DAXPY( M, -CS*SN*AQOAP,
1.6       bertrand 1305:      $                                                A( 1, q ), 1,
                   1306:      $                                                A( 1, p ), 1 )
1.1       bertrand 1307:                                           IF( RSVEC ) THEN
                   1308:                                              CALL DAXPY( MVL, T*APOAQ,
1.6       bertrand 1309:      $                                                   V( 1, p ), 1,
                   1310:      $                                                   V( 1, q ), 1 )
1.1       bertrand 1311:                                              CALL DAXPY( MVL,
1.6       bertrand 1312:      $                                                   -CS*SN*AQOAP,
                   1313:      $                                                   V( 1, q ), 1,
                   1314:      $                                                   V( 1, p ), 1 )
1.1       bertrand 1315:                                           END IF
                   1316:                                           WORK( p ) = WORK( p ) / CS
                   1317:                                           WORK( q ) = WORK( q )*CS
                   1318:                                        ELSE
                   1319:                                           IF( WORK( p ).GE.WORK( q ) )
1.6       bertrand 1320:      $                                        THEN
1.1       bertrand 1321:                                              CALL DAXPY( M, -T*AQOAP,
1.6       bertrand 1322:      $                                                   A( 1, q ), 1,
                   1323:      $                                                   A( 1, p ), 1 )
1.1       bertrand 1324:                                              CALL DAXPY( M, CS*SN*APOAQ,
1.6       bertrand 1325:      $                                                   A( 1, p ), 1,
                   1326:      $                                                   A( 1, q ), 1 )
1.1       bertrand 1327:                                              WORK( p ) = WORK( p )*CS
                   1328:                                              WORK( q ) = WORK( q ) / CS
                   1329:                                              IF( RSVEC ) THEN
                   1330:                                                 CALL DAXPY( MVL,
1.6       bertrand 1331:      $                                               -T*AQOAP,
                   1332:      $                                               V( 1, q ), 1,
                   1333:      $                                               V( 1, p ), 1 )
1.1       bertrand 1334:                                                 CALL DAXPY( MVL,
1.6       bertrand 1335:      $                                               CS*SN*APOAQ,
                   1336:      $                                               V( 1, p ), 1,
                   1337:      $                                               V( 1, q ), 1 )
1.1       bertrand 1338:                                              END IF
                   1339:                                           ELSE
                   1340:                                              CALL DAXPY( M, T*APOAQ,
1.6       bertrand 1341:      $                                                   A( 1, p ), 1,
                   1342:      $                                                   A( 1, q ), 1 )
1.1       bertrand 1343:                                              CALL DAXPY( M,
1.6       bertrand 1344:      $                                                   -CS*SN*AQOAP,
                   1345:      $                                                   A( 1, q ), 1,
                   1346:      $                                                   A( 1, p ), 1 )
1.1       bertrand 1347:                                              WORK( p ) = WORK( p ) / CS
                   1348:                                              WORK( q ) = WORK( q )*CS
                   1349:                                              IF( RSVEC ) THEN
                   1350:                                                 CALL DAXPY( MVL,
1.6       bertrand 1351:      $                                               T*APOAQ, V( 1, p ),
                   1352:      $                                               1, V( 1, q ), 1 )
1.1       bertrand 1353:                                                 CALL DAXPY( MVL,
1.6       bertrand 1354:      $                                               -CS*SN*AQOAP,
                   1355:      $                                               V( 1, q ), 1,
                   1356:      $                                               V( 1, p ), 1 )
1.1       bertrand 1357:                                              END IF
                   1358:                                           END IF
                   1359:                                        END IF
                   1360:                                     END IF
                   1361:                                  END IF
                   1362: *
                   1363:                               ELSE
                   1364:                                  IF( AAPP.GT.AAQQ ) THEN
                   1365:                                     CALL DCOPY( M, A( 1, p ), 1,
1.6       bertrand 1366:      $                                          WORK( N+1 ), 1 )
1.1       bertrand 1367:                                     CALL DLASCL( 'G', 0, 0, AAPP, ONE,
1.6       bertrand 1368:      $                                           M, 1, WORK( N+1 ), LDA,
                   1369:      $                                           IERR )
1.1       bertrand 1370:                                     CALL DLASCL( 'G', 0, 0, AAQQ, ONE,
1.6       bertrand 1371:      $                                           M, 1, A( 1, q ), LDA,
                   1372:      $                                           IERR )
1.1       bertrand 1373:                                     TEMP1 = -AAPQ*WORK( p ) / WORK( q )
                   1374:                                     CALL DAXPY( M, TEMP1, WORK( N+1 ),
1.6       bertrand 1375:      $                                          1, A( 1, q ), 1 )
1.1       bertrand 1376:                                     CALL DLASCL( 'G', 0, 0, ONE, AAQQ,
1.6       bertrand 1377:      $                                           M, 1, A( 1, q ), LDA,
                   1378:      $                                           IERR )
1.1       bertrand 1379:                                     SVA( q ) = AAQQ*DSQRT( DMAX1( ZERO,
1.6       bertrand 1380:      $                                         ONE-AAPQ*AAPQ ) )
1.1       bertrand 1381:                                     MXSINJ = DMAX1( MXSINJ, SFMIN )
                   1382:                                  ELSE
                   1383:                                     CALL DCOPY( M, A( 1, q ), 1,
1.6       bertrand 1384:      $                                          WORK( N+1 ), 1 )
1.1       bertrand 1385:                                     CALL DLASCL( 'G', 0, 0, AAQQ, ONE,
1.6       bertrand 1386:      $                                           M, 1, WORK( N+1 ), LDA,
                   1387:      $                                           IERR )
1.1       bertrand 1388:                                     CALL DLASCL( 'G', 0, 0, AAPP, ONE,
1.6       bertrand 1389:      $                                           M, 1, A( 1, p ), LDA,
                   1390:      $                                           IERR )
1.1       bertrand 1391:                                     TEMP1 = -AAPQ*WORK( q ) / WORK( p )
                   1392:                                     CALL DAXPY( M, TEMP1, WORK( N+1 ),
1.6       bertrand 1393:      $                                          1, A( 1, p ), 1 )
1.1       bertrand 1394:                                     CALL DLASCL( 'G', 0, 0, ONE, AAPP,
1.6       bertrand 1395:      $                                           M, 1, A( 1, p ), LDA,
                   1396:      $                                           IERR )
1.1       bertrand 1397:                                     SVA( p ) = AAPP*DSQRT( DMAX1( ZERO,
1.6       bertrand 1398:      $                                         ONE-AAPQ*AAPQ ) )
1.1       bertrand 1399:                                     MXSINJ = DMAX1( MXSINJ, SFMIN )
                   1400:                                  END IF
                   1401:                               END IF
                   1402: *           END IF ROTOK THEN ... ELSE
                   1403: *
                   1404: *           In the case of cancellation in updating SVA(q)
                   1405: *           .. recompute SVA(q)
                   1406:                               IF( ( SVA( q ) / AAQQ )**2.LE.ROOTEPS )
1.6       bertrand 1407:      $                            THEN
1.1       bertrand 1408:                                  IF( ( AAQQ.LT.ROOTBIG ) .AND.
1.6       bertrand 1409:      $                               ( AAQQ.GT.ROOTSFMIN ) ) THEN
1.1       bertrand 1410:                                     SVA( q ) = DNRM2( M, A( 1, q ), 1 )*
1.6       bertrand 1411:      $                                         WORK( q )
1.1       bertrand 1412:                                  ELSE
                   1413:                                     T = ZERO
1.4       bertrand 1414:                                     AAQQ = ONE
1.1       bertrand 1415:                                     CALL DLASSQ( M, A( 1, q ), 1, T,
1.6       bertrand 1416:      $                                           AAQQ )
1.1       bertrand 1417:                                     SVA( q ) = T*DSQRT( AAQQ )*WORK( q )
                   1418:                                  END IF
                   1419:                               END IF
                   1420:                               IF( ( AAPP / AAPP0 )**2.LE.ROOTEPS ) THEN
                   1421:                                  IF( ( AAPP.LT.ROOTBIG ) .AND.
1.6       bertrand 1422:      $                               ( AAPP.GT.ROOTSFMIN ) ) THEN
1.1       bertrand 1423:                                     AAPP = DNRM2( M, A( 1, p ), 1 )*
1.6       bertrand 1424:      $                                     WORK( p )
1.1       bertrand 1425:                                  ELSE
                   1426:                                     T = ZERO
1.4       bertrand 1427:                                     AAPP = ONE
1.1       bertrand 1428:                                     CALL DLASSQ( M, A( 1, p ), 1, T,
1.6       bertrand 1429:      $                                           AAPP )
1.1       bertrand 1430:                                     AAPP = T*DSQRT( AAPP )*WORK( p )
                   1431:                                  END IF
                   1432:                                  SVA( p ) = AAPP
                   1433:                               END IF
                   1434: *              end of OK rotation
                   1435:                            ELSE
                   1436:                               NOTROT = NOTROT + 1
                   1437: *[RTD]      SKIPPED  = SKIPPED  + 1
                   1438:                               PSKIPPED = PSKIPPED + 1
                   1439:                               IJBLSK = IJBLSK + 1
                   1440:                            END IF
                   1441:                         ELSE
                   1442:                            NOTROT = NOTROT + 1
                   1443:                            PSKIPPED = PSKIPPED + 1
                   1444:                            IJBLSK = IJBLSK + 1
                   1445:                         END IF
                   1446: *
                   1447:                         IF( ( i.LE.SWBAND ) .AND. ( IJBLSK.GE.BLSKIP ) )
1.6       bertrand 1448:      $                      THEN
1.1       bertrand 1449:                            SVA( p ) = AAPP
                   1450:                            NOTROT = 0
                   1451:                            GO TO 2011
                   1452:                         END IF
                   1453:                         IF( ( i.LE.SWBAND ) .AND.
1.6       bertrand 1454:      $                      ( PSKIPPED.GT.ROWSKIP ) ) THEN
1.1       bertrand 1455:                            AAPP = -AAPP
                   1456:                            NOTROT = 0
                   1457:                            GO TO 2203
                   1458:                         END IF
                   1459: *
                   1460:  2200                CONTINUE
                   1461: *        end of the q-loop
                   1462:  2203                CONTINUE
                   1463: *
                   1464:                      SVA( p ) = AAPP
                   1465: *
                   1466:                   ELSE
                   1467: *
                   1468:                      IF( AAPP.EQ.ZERO )NOTROT = NOTROT +
1.6       bertrand 1469:      $                   MIN0( jgl+KBL-1, N ) - jgl + 1
1.1       bertrand 1470:                      IF( AAPP.LT.ZERO )NOTROT = 0
                   1471: *
                   1472:                   END IF
                   1473: *
                   1474:  2100          CONTINUE
                   1475: *     end of the p-loop
                   1476:  2010       CONTINUE
                   1477: *     end of the jbc-loop
                   1478:  2011       CONTINUE
                   1479: *2011 bailed out of the jbc-loop
                   1480:             DO 2012 p = igl, MIN0( igl+KBL-1, N )
                   1481:                SVA( p ) = DABS( SVA( p ) )
                   1482:  2012       CONTINUE
                   1483: ***
                   1484:  2000    CONTINUE
                   1485: *2000 :: end of the ibr-loop
                   1486: *
                   1487: *     .. update SVA(N)
                   1488:          IF( ( SVA( N ).LT.ROOTBIG ) .AND. ( SVA( N ).GT.ROOTSFMIN ) )
1.6       bertrand 1489:      $       THEN
1.1       bertrand 1490:             SVA( N ) = DNRM2( M, A( 1, N ), 1 )*WORK( N )
                   1491:          ELSE
                   1492:             T = ZERO
1.4       bertrand 1493:             AAPP = ONE
1.1       bertrand 1494:             CALL DLASSQ( M, A( 1, N ), 1, T, AAPP )
                   1495:             SVA( N ) = T*DSQRT( AAPP )*WORK( N )
                   1496:          END IF
                   1497: *
                   1498: *     Additional steering devices
                   1499: *
                   1500:          IF( ( i.LT.SWBAND ) .AND. ( ( MXAAPQ.LE.ROOTTOL ) .OR.
1.6       bertrand 1501:      $       ( ISWROT.LE.N ) ) )SWBAND = i
1.1       bertrand 1502: *
                   1503:          IF( ( i.GT.SWBAND+1 ) .AND. ( MXAAPQ.LT.DSQRT( DBLE( N ) )*
1.6       bertrand 1504:      $       TOL ) .AND. ( DBLE( N )*MXAAPQ*MXSINJ.LT.TOL ) ) THEN
1.1       bertrand 1505:             GO TO 1994
                   1506:          END IF
                   1507: *
                   1508:          IF( NOTROT.GE.EMPTSW )GO TO 1994
                   1509: *
                   1510:  1993 CONTINUE
                   1511: *     end i=1:NSWEEP loop
                   1512: *
                   1513: * #:( Reaching this point means that the procedure has not converged.
                   1514:       INFO = NSWEEP - 1
                   1515:       GO TO 1995
                   1516: *
                   1517:  1994 CONTINUE
                   1518: * #:) Reaching this point means numerical convergence after the i-th
                   1519: *     sweep.
                   1520: *
                   1521:       INFO = 0
                   1522: * #:) INFO = 0 confirms successful iterations.
                   1523:  1995 CONTINUE
                   1524: *
                   1525: *     Sort the singular values and find how many are above
                   1526: *     the underflow threshold.
                   1527: *
                   1528:       N2 = 0
                   1529:       N4 = 0
                   1530:       DO 5991 p = 1, N - 1
                   1531:          q = IDAMAX( N-p+1, SVA( p ), 1 ) + p - 1
                   1532:          IF( p.NE.q ) THEN
                   1533:             TEMP1 = SVA( p )
                   1534:             SVA( p ) = SVA( q )
                   1535:             SVA( q ) = TEMP1
                   1536:             TEMP1 = WORK( p )
                   1537:             WORK( p ) = WORK( q )
                   1538:             WORK( q ) = TEMP1
                   1539:             CALL DSWAP( M, A( 1, p ), 1, A( 1, q ), 1 )
                   1540:             IF( RSVEC )CALL DSWAP( MVL, V( 1, p ), 1, V( 1, q ), 1 )
                   1541:          END IF
                   1542:          IF( SVA( p ).NE.ZERO ) THEN
                   1543:             N4 = N4 + 1
1.4       bertrand 1544:             IF( SVA( p )*SKL.GT.SFMIN )N2 = N2 + 1
1.1       bertrand 1545:          END IF
                   1546:  5991 CONTINUE
                   1547:       IF( SVA( N ).NE.ZERO ) THEN
                   1548:          N4 = N4 + 1
1.4       bertrand 1549:          IF( SVA( N )*SKL.GT.SFMIN )N2 = N2 + 1
1.1       bertrand 1550:       END IF
                   1551: *
                   1552: *     Normalize the left singular vectors.
                   1553: *
                   1554:       IF( LSVEC .OR. UCTOL ) THEN
                   1555:          DO 1998 p = 1, N2
                   1556:             CALL DSCAL( M, WORK( p ) / SVA( p ), A( 1, p ), 1 )
                   1557:  1998    CONTINUE
                   1558:       END IF
                   1559: *
                   1560: *     Scale the product of Jacobi rotations (assemble the fast rotations).
                   1561: *
                   1562:       IF( RSVEC ) THEN
                   1563:          IF( APPLV ) THEN
                   1564:             DO 2398 p = 1, N
                   1565:                CALL DSCAL( MVL, WORK( p ), V( 1, p ), 1 )
                   1566:  2398       CONTINUE
                   1567:          ELSE
                   1568:             DO 2399 p = 1, N
                   1569:                TEMP1 = ONE / DNRM2( MVL, V( 1, p ), 1 )
                   1570:                CALL DSCAL( MVL, TEMP1, V( 1, p ), 1 )
                   1571:  2399       CONTINUE
                   1572:          END IF
                   1573:       END IF
                   1574: *
                   1575: *     Undo scaling, if necessary (and possible).
1.11      bertrand 1576:       IF( ( ( SKL.GT.ONE ) .AND. ( SVA( 1 ).LT.( BIG / SKL) ) )
                   1577:      $    .OR. ( ( SKL.LT.ONE ) .AND. ( SVA( MAX( N2, 1 ) ) .GT.
1.6       bertrand 1578:      $    ( SFMIN / SKL) ) ) ) THEN
1.1       bertrand 1579:          DO 2400 p = 1, N
1.11      bertrand 1580:             SVA( P ) = SKL*SVA( P )
1.1       bertrand 1581:  2400    CONTINUE
1.4       bertrand 1582:          SKL= ONE
1.1       bertrand 1583:       END IF
                   1584: *
1.4       bertrand 1585:       WORK( 1 ) = SKL
                   1586: *     The singular values of A are SKL*SVA(1:N). If SKL.NE.ONE
1.1       bertrand 1587: *     then some of the singular values may overflow or underflow and
                   1588: *     the spectrum is given in this factored representation.
                   1589: *
                   1590:       WORK( 2 ) = DBLE( N4 )
                   1591: *     N4 is the number of computed nonzero singular values of A.
                   1592: *
                   1593:       WORK( 3 ) = DBLE( N2 )
                   1594: *     N2 is the number of singular values of A greater than SFMIN.
                   1595: *     If N2<N, SVA(N2:N) contains ZEROS and/or denormalized numbers
                   1596: *     that may carry some information.
                   1597: *
                   1598:       WORK( 4 ) = DBLE( i )
                   1599: *     i is the index of the last sweep before declaring convergence.
                   1600: *
                   1601:       WORK( 5 ) = MXAAPQ
                   1602: *     MXAAPQ is the largest absolute value of scaled pivots in the
                   1603: *     last sweep
                   1604: *
                   1605:       WORK( 6 ) = MXSINJ
                   1606: *     MXSINJ is the largest absolute value of the sines of Jacobi angles
                   1607: *     in the last sweep
                   1608: *
                   1609:       RETURN
                   1610: *     ..
                   1611: *     .. END OF DGESVJ
                   1612: *     ..
                   1613:       END

CVSweb interface <joel.bertrand@systella.fr>