0001:       SUBROUTINE SGESVJ( JOBA, JOBU, JOBV, M, N, A, LDA, SVA,
0002:      &                     MV, V, LDV, WORK, LWORK, INFO )
0003: *
0004: *  -- LAPACK routine (version 3.2)                                    --
0005: *
0006: *  -- Contributed by Zlatko Drmac of the University of Zagreb and     --
0007: *  -- Kresimir Veselic of the Fernuniversitaet Hagen                  --
0008: *  -- November 2008                                                   --
0009: *
0010: *  -- LAPACK is a software package provided by Univ. of Tennessee,    --
0011: *  -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
0012: *
0013: * This routine is also part of SIGMA (version 1.23, October 23. 2008.)
0014: * SIGMA is a library of algorithms for highly accurate algorithms for
0015: * computation of SVD, PSVD, QSVD, (H,K)-SVD, and for solution of the
0016: * eigenvalue problems Hx = lambda M x, H M x = lambda x with H, M > 0.
0017: *
0018: *     -#- Scalar Arguments -#-
0019: *
0020:       IMPLICIT    NONE
0021:       INTEGER     INFO, LDA, LDV, LWORK, M, MV, N
0022:       CHARACTER*1 JOBA, JOBU, JOBV
0023: *
0024: *     -#- Array Arguments -#-
0025: *
0026:       REAL        A( LDA, * ), SVA( N ), V( LDV, * ), WORK( LWORK )
0027: *     ..
0028: *
0029: *  Purpose
0030: *  ~~~~~~~
0031: *  SGESVJ computes the singular value decomposition (SVD) of a real
0032: *  M-by-N matrix A, where M >= N. The SVD of A is written as
0033: *                                     [++]   [xx]   [x0]   [xx]
0034: *               A = U * SIGMA * V^t,  [++] = [xx] * [ox] * [xx]
0035: *                                     [++]   [xx]
0036: *  where SIGMA is an N-by-N diagonal matrix, U is an M-by-N orthonormal
0037: *  matrix, and V is an N-by-N orthogonal matrix. The diagonal elements
0038: *  of SIGMA are the singular values of A. The columns of U and V are the
0039: *  left and the right singular vectors of A, respectively.
0040: *
0041: *  Further Details
0042: *  ~~~~~~~~~~~~~~~
0043: *  The orthogonal N-by-N matrix V is obtained as a product of Jacobi plane
0044: *  rotations. The rotations are implemented as fast scaled rotations of
0045: *  Anda and Park [1]. In the case of underflow of the Jacobi angle, a
0046: *  modified Jacobi transformation of Drmac [4] is used. Pivot strategy uses
0047: *  column interchanges of de Rijk [2]. The relative accuracy of the computed
0048: *  singular values and the accuracy of the computed singular vectors (in
0049: *  angle metric) is as guaranteed by the theory of Demmel and Veselic [3].
0050: *  The condition number that determines the accuracy in the full rank case
0051: *  is essentially min_{D=diag} kappa(A*D), where kappa(.) is the
0052: *  spectral condition number. The best performance of this Jacobi SVD
0053: *  procedure is achieved if used in an  accelerated version of Drmac and
0054: *  Veselic [5,6], and it is the kernel routine in the SIGMA library [7].
0055: *  Some tunning parameters (marked with [TP]) are available for the
0056: *  implementer.
0057: *  The computational range for the nonzero singular values is the  machine
0058: *  number interval ( UNDERFLOW , OVERFLOW ). In extreme cases, even
0059: *  denormalized singular values can be computed with the corresponding
0060: *  gradual loss of accurate digits.
0061: *
0062: *  Contributors
0063: *  ~~~~~~~~~~~~
0064: *  Zlatko Drmac (Zagreb, Croatia) and Kresimir Veselic (Hagen, Germany)
0065: *
0066: *  References
0067: *  ~~~~~~~~~~
0068: * [1] A. A. Anda and H. Park: Fast plane rotations with dynamic scaling.
0069: *     SIAM J. matrix Anal. Appl., Vol. 15 (1994), pp. 162-174.
0070: * [2] P. P. M. De Rijk: A one-sided Jacobi algorithm for computing the
0071: *     singular value decomposition on a vector computer.
0072: *     SIAM J. Sci. Stat. Comp., Vol. 10 (1998), pp. 359-371.
0073: * [3] J. Demmel and K. Veselic: Jacobi method is more accurate than QR.
0074: * [4] Z. Drmac: Implementation of Jacobi rotations for accurate singular
0075: *     value computation in floating point arithmetic.
0076: *     SIAM J. Sci. Comp., Vol. 18 (1997), pp. 1200-1222.
0077: * [5] Z. Drmac and K. Veselic: New fast and accurate Jacobi SVD algorithm I.
0078: *     SIAM J. Matrix Anal. Appl. Vol. 35, No. 2 (2008), pp. 1322-1342.
0079: *     LAPACK Working note 169.
0080: * [6] Z. Drmac and K. Veselic: New fast and accurate Jacobi SVD algorithm II.
0081: *     SIAM J. Matrix Anal. Appl. Vol. 35, No. 2 (2008), pp. 1343-1362.
0082: *     LAPACK Working note 170.
0083: * [7] Z. Drmac: SIGMA - mathematical software library for accurate SVD, PSV,
0084: *     QSVD, (H,K)-SVD computations.
0085: *     Department of Mathematics, University of Zagreb, 2008.
0086: *
0087: *  Bugs, Examples and Comments
0088: *  ~~~~~~~~~~~~~~~~~~~~~~~~~~~
0089: *  Please report all bugs and send interesting test examples and comments to
0090: *  drmac@math.hr. Thank you.
0091: *
0092: *  Arguments
0093: *  ~~~~~~~~~
0094: *
0095: *  JOBA    (input) CHARACTER* 1
0096: *          Specifies the structure of A.
0097: *          = 'L': The input matrix A is lower triangular;
0098: *          = 'U': The input matrix A is upper triangular;
0099: *          = 'G': The input matrix A is general M-by-N matrix, M >= N.
0100: *
0101: *  JOBU    (input) CHARACTER*1
0102: *          Specifies whether to compute the left singular vectors
0103: *          (columns of U):
0104: *
0105: *          = 'U': The left singular vectors corresponding to the nonzero
0106: *                 singular values are computed and returned in the leading
0107: *                 columns of A. See more details in the description of A.
0108: *                 The default numerical orthogonality threshold is set to
0109: *                 approximately TOL=CTOL*EPS, CTOL=SQRT(M), EPS=SLAMCH('E').
0110: *          = 'C': Analogous to JOBU='U', except that user can control the
0111: *                 level of numerical orthogonality of the computed left
0112: *                 singular vectors. TOL can be set to TOL = CTOL*EPS, where
0113: *                 CTOL is given on input in the array WORK.
0114: *                 No CTOL smaller than ONE is allowed. CTOL greater
0115: *                 than 1 / EPS is meaningless. The option 'C'
0116: *                 can be used if M*EPS is satisfactory orthogonality
0117: *                 of the computed left singular vectors, so CTOL=M could
0118: *                 save few sweeps of Jacobi rotations.
0119: *                 See the descriptions of A and WORK(1).
0120: *          = 'N': The matrix U is not computed. However, see the
0121: *                 description of A.
0122: *
0123: *  JOBV    (input) CHARACTER*1
0124: *          Specifies whether to compute the right singular vectors, that
0125: *          is, the matrix V:
0126: *          = 'V' : the matrix V is computed and returned in the array V
0127: *          = 'A' : the Jacobi rotations are applied to the MV-by-N
0128: *                  array V. In other words, the right singular vector
0129: *                  matrix V is not computed explicitly; instead it is
0130: *                  applied to an MV-by-N matrix initially stored in the
0131: *                  first MV rows of V.
0132: *          = 'N' : the matrix V is not computed and the array V is not
0133: *                  referenced
0134: *
0135: *  M       (input) INTEGER
0136: *          The number of rows of the input matrix A.  M >= 0.
0137: *
0138: *  N       (input) INTEGER
0139: *          The number of columns of the input matrix A.
0140: *          M >= N >= 0.
0141: *
0142: *  A       (input/output) REAL array, dimension (LDA,N)
0143: *          On entry, the M-by-N matrix A.
0144: *          On exit,
0145: *          If JOBU .EQ. 'U' .OR. JOBU .EQ. 'C':
0146: *          ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
0147: *                 If INFO .EQ. 0,
0148: *                 ~~~~~~~~~~~~~~~
0149: *                 RANKA orthonormal columns of U are returned in the
0150: *                 leading RANKA columns of the array A. Here RANKA <= N
0151: *                 is the number of computed singular values of A that are
0152: *                 above the underflow threshold SLAMCH('S'). The singular
0153: *                 vectors corresponding to underflowed or zero singular
0154: *                 values are not computed. The value of RANKA is returned
0155: *                 in the array WORK as RANKA=NINT(WORK(2)). Also see the
0156: *                 descriptions of SVA and WORK. The computed columns of U
0157: *                 are mutually numerically orthogonal up to approximately
0158: *                 TOL=SQRT(M)*EPS (default); or TOL=CTOL*EPS (JOBU.EQ.'C'),
0159: *                 see the description of JOBU.
0160: *                 If INFO .GT. 0,
0161: *                 ~~~~~~~~~~~~~~~
0162: *                 the procedure SGESVJ did not converge in the given number
0163: *                 of iterations (sweeps). In that case, the computed
0164: *                 columns of U may not be orthogonal up to TOL. The output
0165: *                 U (stored in A), SIGMA (given by the computed singular
0166: *                 values in SVA(1:N)) and V is still a decomposition of the
0167: *                 input matrix A in the sense that the residual
0168: *                 ||A-SCALE*U*SIGMA*V^T||_2 / ||A||_2 is small.
0169: *
0170: *          If JOBU .EQ. 'N':
0171: *          ~~~~~~~~~~~~~~~~~
0172: *                 If INFO .EQ. 0
0173: *                 ~~~~~~~~~~~~~~
0174: *                 Note that the left singular vectors are 'for free' in the
0175: *                 one-sided Jacobi SVD algorithm. However, if only the
0176: *                 singular values are needed, the level of numerical
0177: *                 orthogonality of U is not an issue and iterations are
0178: *                 stopped when the columns of the iterated matrix are
0179: *                 numerically orthogonal up to approximately M*EPS. Thus,
0180: *                 on exit, A contains the columns of U scaled with the
0181: *                 corresponding singular values.
0182: *                 If INFO .GT. 0,
0183: *                 ~~~~~~~~~~~~~~~
0184: *                 the procedure SGESVJ did not converge in the given number
0185: *                 of iterations (sweeps).
0186: *
0187: *  LDA     (input) INTEGER
0188: *          The leading dimension of the array A.  LDA >= max(1,M).
0189: *
0190: *  SVA     (workspace/output) REAL array, dimension (N)
0191: *          On exit,
0192: *          If INFO .EQ. 0,
0193: *          ~~~~~~~~~~~~~~~
0194: *          depending on the value SCALE = WORK(1), we have:
0195: *                 If SCALE .EQ. ONE:
0196: *                 ~~~~~~~~~~~~~~~~~~
0197: *                 SVA(1:N) contains the computed singular values of A.
0198: *                 During the computation SVA contains the Euclidean column
0199: *                 norms of the iterated matrices in the array A.
0200: *                 If SCALE .NE. ONE:
0201: *                 ~~~~~~~~~~~~~~~~~~
0202: *                 The singular values of A are SCALE*SVA(1:N), and this
0203: *                 factored representation is due to the fact that some of the
0204: *                 singular values of A might underflow or overflow.
0205: *
0206: *          If INFO .GT. 0,
0207: *          ~~~~~~~~~~~~~~~
0208: *          the procedure SGESVJ did not converge in the given number of
0209: *          iterations (sweeps) and SCALE*SVA(1:N) may not be accurate.
0210: *
0211: *  MV      (input) INTEGER
0212: *          If JOBV .EQ. 'A', then the product of Jacobi rotations in SGESVJ
0213: *          is applied to the first MV rows of V. See the description of JOBV.
0214: *
0215: *  V       (input/output) REAL array, dimension (LDV,N)
0216: *          If JOBV = 'V', then V contains on exit the N-by-N matrix of
0217: *                         the right singular vectors;
0218: *          If JOBV = 'A', then V contains the product of the computed right
0219: *                         singular vector matrix and the initial matrix in
0220: *                         the array V.
0221: *          If JOBV = 'N', then V is not referenced.
0222: *
0223: *  LDV     (input) INTEGER
0224: *          The leading dimension of the array V, LDV .GE. 1.
0225: *          If JOBV .EQ. 'V', then LDV .GE. max(1,N).
0226: *          If JOBV .EQ. 'A', then LDV .GE. max(1,MV) .
0227: *
0228: *  WORK    (input/workspace/output) REAL array, dimension max(4,M+N).
0229: *          On entry,
0230: *          If JOBU .EQ. 'C',
0231: *          ~~~~~~~~~~~~~~~~~
0232: *          WORK(1) = CTOL, where CTOL defines the threshold for convergence.
0233: *                    The process stops if all columns of A are mutually
0234: *                    orthogonal up to CTOL*EPS, EPS=SLAMCH('E').
0235: *                    It is required that CTOL >= ONE, i.e. it is not
0236: *                    allowed to force the routine to obtain orthogonality
0237: *                    below EPSILON.
0238: *          On exit,
0239: *          WORK(1) = SCALE is the scaling factor such that SCALE*SVA(1:N)
0240: *                    are the computed singular vcalues of A.
0241: *                    (See description of SVA().)
0242: *          WORK(2) = NINT(WORK(2)) is the number of the computed nonzero
0243: *                    singular values.
0244: *          WORK(3) = NINT(WORK(3)) is the number of the computed singular
0245: *                    values that are larger than the underflow threshold.
0246: *          WORK(4) = NINT(WORK(4)) is the number of sweeps of Jacobi
0247: *                    rotations needed for numerical convergence.
0248: *          WORK(5) = max_{i.NE.j} |COS(A(:,i),A(:,j))| in the last sweep.
0249: *                    This is useful information in cases when SGESVJ did
0250: *                    not converge, as it can be used to estimate whether
0251: *                    the output is stil useful and for post festum analysis.
0252: *          WORK(6) = the largest absolute value over all sines of the
0253: *                    Jacobi rotation angles in the last sweep. It can be
0254: *                    useful for a post festum analysis.
0255: *
0256: *  LWORK   length of WORK, WORK >= MAX(6,M+N)
0257: *
0258: *  INFO    (output) INTEGER
0259: *          = 0 : successful exit.
0260: *          < 0 : if INFO = -i, then the i-th argument had an illegal value
0261: *          > 0 : SGESVJ did not converge in the maximal allowed number (30)
0262: *                of sweeps. The output may still be useful. See the
0263: *                description of WORK.
0264: *
0265: *     Local Parameters
0266: *
0267:       REAL        ZERO,         HALF,         ONE,         TWO
0268:       PARAMETER ( ZERO = 0.0E0, HALF = 0.5E0, ONE = 1.0E0, TWO = 2.0E0 )
0269:       INTEGER     NSWEEP
0270:       PARAMETER ( NSWEEP = 30 )
0271: *
0272: *     Local Scalars
0273: *
0274:       REAL    AAPP,    AAPP0,    AAPQ,     AAQQ,    APOAQ,     AQOAP,
0275:      &        BIG,     BIGTHETA, CS,       CTOL,    EPSILON,   LARGE,
0276:      &        MXAAPQ,  MXSINJ,   ROOTBIG,  ROOTEPS, ROOTSFMIN, ROOTTOL,
0277:      &        SCALE,   SFMIN,    SMALL,    SN,      T,         TEMP1,
0278:      &        THETA,   THSIGN,   TOL
0279:       INTEGER BLSKIP,  EMPTSW,   i,        ibr,     IERR,      igl,
0280:      &        IJBLSK,  ir1,      ISWROT,   jbc,     jgl,       KBL,
0281:      &        LKAHEAD, MVL,      N2,       N34,     N4,        NBL,
0282:      &        NOTROT,  p,        PSKIPPED, q,       ROWSKIP,   SWBAND
0283:       LOGICAL APPLV,   GOSCALE,  LOWER,    LSVEC,   NOSCALE,   ROTOK,
0284:      &        RSVEC,   UCTOL,    UPPER
0285: *
0286: *     Local Arrays
0287: *
0288:       REAL      FASTR(5)
0289: *
0290: *     Intrinsic Functions
0291: *
0292:       INTRINSIC ABS, AMAX1, AMIN1, FLOAT, MIN0, SIGN, SQRT
0293: *
0294: *     External Functions
0295: *     .. from BLAS
0296:       REAL             SDOT, SNRM2
0297:       EXTERNAL         SDOT, SNRM2
0298:       INTEGER          ISAMAX
0299:       EXTERNAL         ISAMAX
0300: *     .. from LAPACK
0301:       REAL             SLAMCH
0302:       EXTERNAL         SLAMCH
0303:       LOGICAL          LSAME
0304:       EXTERNAL         LSAME
0305: *
0306: *     External Subroutines
0307: *     .. from BLAS
0308:       EXTERNAL  SAXPY,  SCOPY, SROTM, SSCAL, SSWAP
0309: *     .. from LAPACK
0310:       EXTERNAL  SLASCL, SLASET, SLASSQ, XERBLA
0311: *
0312:       EXTERNAL  SGSVJ0, SGSVJ1
0313: *
0314: *     Test the input arguments
0315: *
0316:       LSVEC = LSAME( JOBU, 'U' )
0317:       UCTOL = LSAME( JOBU, 'C' )
0318:       RSVEC = LSAME( JOBV, 'V' )
0319:       APPLV = LSAME( JOBV, 'A' )
0320:       UPPER = LSAME( JOBA, 'U' )
0321:       LOWER = LSAME( JOBA, 'L' )
0322: *
0323:       IF ( .NOT.( UPPER .OR. LOWER .OR. LSAME(JOBA,'G') ) ) THEN
0324:          INFO = - 1
0325:       ELSE IF ( .NOT.( LSVEC .OR. UCTOL .OR. LSAME(JOBU,'N') ) ) THEN
0326:          INFO = - 2
0327:       ELSE IF ( .NOT.( RSVEC .OR. APPLV .OR. LSAME(JOBV,'N') ) ) THEN
0328:          INFO = - 3
0329:       ELSE IF ( M .LT. 0 ) THEN
0330:          INFO = - 4
0331:       ELSE IF ( ( N .LT. 0 ) .OR. ( N .GT. M ) ) THEN
0332:          INFO = - 5
0333:       ELSE IF ( LDA .LT. M ) THEN
0334:          INFO = - 7
0335:       ELSE IF ( MV .LT. 0 ) THEN
0336:          INFO = - 9
0337:       ELSE IF ( ( RSVEC .AND. (LDV .LT. N ) ) .OR.
0338:      &          ( APPLV .AND. (LDV .LT. MV) ) ) THEN
0339:          INFO = -11
0340:       ELSE IF ( UCTOL .AND. (WORK(1) .LE. ONE) ) THEN
0341:          INFO = - 12
0342:       ELSE IF ( LWORK .LT. MAX0( M + N , 6 ) ) THEN
0343:          INFO = - 13
0344:       ELSE
0345:          INFO =   0
0346:       END IF
0347: *
0348: *     #:(
0349:       IF ( INFO .NE. 0 ) THEN
0350:          CALL XERBLA( 'SGESVJ', -INFO )
0351:          RETURN
0352:       END IF
0353: *
0354: * #:) Quick return for void matrix
0355: *
0356:       IF ( ( M .EQ. 0 ) .OR. ( N .EQ. 0 ) ) RETURN
0357: *
0358: *     Set numerical parameters
0359: *     The stopping criterion for Jacobi rotations is
0360: *
0361: *     max_{i<>j}|A(:,i)^T * A(:,j)|/(||A(:,i)||*||A(:,j)||) < CTOL*EPS
0362: *
0363: *     where EPS is the round-off and CTOL is defined as follows:
0364: *
0365:       IF ( UCTOL ) THEN
0366: *        ... user controlled
0367:          CTOL = WORK(1)
0368:       ELSE
0369: *        ... default
0370:          IF ( LSVEC .OR. RSVEC .OR. APPLV ) THEN
0371:             CTOL = SQRT(FLOAT(M))
0372:          ELSE
0373:             CTOL = FLOAT(M)
0374:          END IF
0375:       END IF
0376: *     ... and the machine dependent parameters are
0377: *[!]  (Make sure that SLAMCH() works properly on the target machine.)
0378: *
0379:       EPSILON     = SLAMCH('Epsilon')
0380:       ROOTEPS     = SQRT(EPSILON)
0381:       SFMIN       = SLAMCH('SafeMinimum')
0382:       ROOTSFMIN   = SQRT(SFMIN)
0383:       SMALL       = SFMIN   / EPSILON
0384:       BIG         = SLAMCH('Overflow')
0385:       ROOTBIG     = ONE   / ROOTSFMIN
0386:       LARGE       = BIG  / SQRT(FLOAT(M*N))
0387:       BIGTHETA    = ONE / ROOTEPS
0388: *
0389:       TOL         = CTOL * EPSILON
0390:       ROOTTOL     = SQRT(TOL)
0391: *
0392:       IF ( FLOAT(M)*EPSILON .GE. ONE ) THEN
0393:          INFO = - 5
0394:          CALL XERBLA( 'SGESVJ', -INFO )
0395:          RETURN
0396:       END IF
0397: *
0398: *     Initialize the right singular vector matrix.
0399: *
0400:       IF ( RSVEC )  THEN
0401:          MVL = N
0402:          CALL SLASET( 'A', MVL, N, ZERO, ONE, V, LDV )
0403:       ELSE IF ( APPLV ) THEN
0404:          MVL = MV
0405:       END IF
0406:       RSVEC = RSVEC .OR. APPLV
0407: *
0408: *     Initialize SVA( 1:N ) = ( ||A e_i||_2, i = 1:N )
0409: *(!)  If necessary, scale A to protect the largest singular value
0410: *     from overflow. It is possible that saving the largest singular
0411: *     value destroys the information about the small ones.
0412: *     This initial scaling is almost minimal in the sense that the
0413: *     goal is to make sure that no column norm overflows, and that
0414: *     SQRT(N)*max_i SVA(i) does not overflow. If INFinite entries
0415: *     in A are detected, the procedure returns with INFO=-6.
0416: *
0417:       SCALE    = ONE / SQRT(FLOAT(M)*FLOAT(N))
0418:       NOSCALE  = .TRUE.
0419:       GOSCALE  = .TRUE.
0420: *
0421:       IF ( LOWER ) THEN
0422: *        the input matrix is M-by-N lower triangular (trapezoidal)
0423:          DO 1874 p = 1, N
0424:             AAPP = ZERO
0425:             AAQQ = ZERO
0426:             CALL SLASSQ( M-p+1, A(p,p), 1, AAPP, AAQQ )
0427:             IF ( AAPP .GT. BIG ) THEN
0428:                INFO = - 6
0429:                CALL XERBLA( 'SGESVJ', -INFO )
0430:                RETURN
0431:             END IF
0432:             AAQQ = SQRT(AAQQ)
0433:             IF ( ( AAPP .LT. (BIG / AAQQ) ) .AND. NOSCALE ) THEN
0434:                SVA(p)   = AAPP * AAQQ
0435:             ELSE
0436:                NOSCALE = .FALSE.
0437:                SVA(p)   = AAPP * ( AAQQ * SCALE )
0438:                IF ( GOSCALE ) THEN
0439:                   GOSCALE = .FALSE.
0440:                   DO 1873 q = 1, p - 1
0441:                      SVA(q) = SVA(q)*SCALE
0442:  1873             CONTINUE
0443:                END IF
0444:             END IF
0445:  1874    CONTINUE
0446:       ELSE IF ( UPPER ) THEN
0447: *        the input matrix is M-by-N upper triangular (trapezoidal)
0448:          DO 2874 p = 1, N
0449:             AAPP = ZERO
0450:             AAQQ = ZERO
0451:             CALL SLASSQ( p, A(1,p), 1, AAPP, AAQQ )
0452:             IF ( AAPP .GT. BIG ) THEN
0453:                INFO = - 6
0454:                CALL XERBLA( 'SGESVJ', -INFO )
0455:                RETURN
0456:             END IF
0457:             AAQQ = SQRT(AAQQ)
0458:             IF ( ( AAPP .LT. (BIG / AAQQ) ) .AND. NOSCALE ) THEN
0459:                SVA(p)   = AAPP * AAQQ
0460:             ELSE
0461:                NOSCALE = .FALSE.
0462:                SVA(p)   = AAPP * ( AAQQ * SCALE )
0463:                IF ( GOSCALE ) THEN
0464:                   GOSCALE = .FALSE.
0465:                   DO 2873 q = 1, p - 1
0466:                      SVA(q) = SVA(q)*SCALE
0467:  2873             CONTINUE
0468:                END IF
0469:             END IF
0470:  2874    CONTINUE
0471:       ELSE
0472: *        the input matrix is M-by-N general dense
0473:          DO 3874 p = 1, N
0474:             AAPP = ZERO
0475:             AAQQ = ZERO
0476:             CALL SLASSQ( M, A(1,p), 1, AAPP, AAQQ )
0477:             IF ( AAPP .GT. BIG ) THEN
0478:                INFO = - 6
0479:                CALL XERBLA( 'SGESVJ', -INFO )
0480:                RETURN
0481:             END IF
0482:             AAQQ = SQRT(AAQQ)
0483:             IF ( ( AAPP .LT. (BIG / AAQQ) ) .AND. NOSCALE ) THEN
0484:                SVA(p)  = AAPP * AAQQ
0485:             ELSE
0486:                NOSCALE = .FALSE.
0487:                SVA(p)  = AAPP * ( AAQQ * SCALE )
0488:                IF ( GOSCALE ) THEN
0489:                   GOSCALE = .FALSE.
0490:                   DO 3873 q = 1, p - 1
0491:                      SVA(q) = SVA(q)*SCALE
0492:  3873             CONTINUE
0493:                END IF
0494:             END IF
0495:  3874    CONTINUE
0496:       END IF
0497: *
0498:       IF ( NOSCALE ) SCALE = ONE
0499: *
0500: *     Move the smaller part of the spectrum from the underflow threshold
0501: *(!)  Start by determining the position of the nonzero entries of the
0502: *     array SVA() relative to ( SFMIN, BIG ).
0503: *
0504:       AAPP = ZERO
0505:       AAQQ = BIG
0506:       DO 4781 p = 1, N
0507:          IF ( SVA(p) .NE. ZERO ) AAQQ = AMIN1( AAQQ, SVA(p) )
0508:          AAPP = AMAX1( AAPP, SVA(p) )
0509:  4781 CONTINUE
0510: *
0511: * #:) Quick return for zero matrix
0512: *
0513:       IF ( AAPP .EQ. ZERO ) THEN
0514:          IF ( LSVEC ) CALL SLASET( 'G', M, N, ZERO, ONE, A, LDA )
0515:          WORK(1) = ONE
0516:          WORK(2) = ZERO
0517:          WORK(3) = ZERO
0518:          WORK(4) = ZERO
0519:          WORK(5) = ZERO
0520:          WORK(6) = ZERO
0521:          RETURN
0522:       END IF
0523: *
0524: * #:) Quick return for one-column matrix
0525: *
0526:       IF ( N .EQ. 1 ) THEN
0527:          IF ( LSVEC )
0528:      &      CALL SLASCL( 'G',0,0,SVA(1),SCALE,M,1,A(1,1),LDA,IERR )
0529:          WORK(1) = ONE / SCALE
0530:          IF ( SVA(1) .GE. SFMIN ) THEN
0531:             WORK(2) = ONE
0532:          ELSE
0533:             WORK(2) = ZERO
0534:          END IF
0535:          WORK(3) = ZERO
0536:          WORK(4) = ZERO
0537:          WORK(5) = ZERO
0538:          WORK(6) = ZERO
0539:          RETURN
0540:       END IF
0541: *
0542: *     Protect small singular values from underflow, and try to
0543: *     avoid underflows/overflows in computing Jacobi rotations.
0544: *
0545:       SN    = SQRT( SFMIN / EPSILON  )
0546:       TEMP1 = SQRT(  BIG /  FLOAT(N) )
0547:       IF ( (AAPP.LE.SN).OR.(AAQQ.GE.TEMP1)
0548:      &   .OR.((SN.LE.AAQQ).AND.(AAPP.LE.TEMP1)) ) THEN
0549:          TEMP1 = AMIN1(BIG,TEMP1/AAPP)
0550: *         AAQQ  = AAQQ*TEMP1
0551: *         AAPP  = AAPP*TEMP1
0552:       ELSE IF ( (AAQQ.LE.SN).AND.(AAPP.LE.TEMP1) ) THEN
0553:          TEMP1 = AMIN1( SN / AAQQ, BIG/(AAPP*SQRT(FLOAT(N))) )
0554: *         AAQQ  = AAQQ*TEMP1
0555: *         AAPP  = AAPP*TEMP1
0556:       ELSE IF ( (AAQQ.GE.SN).AND.(AAPP.GE.TEMP1) ) THEN
0557:          TEMP1 = AMAX1( SN / AAQQ, TEMP1 / AAPP )
0558: *         AAQQ  = AAQQ*TEMP1
0559: *         AAPP  = AAPP*TEMP1
0560:       ELSE IF ( (AAQQ.LE.SN).AND.(AAPP.GE.TEMP1) ) THEN
0561:          TEMP1 = AMIN1( SN / AAQQ, BIG / (SQRT(FLOAT(N))*AAPP))
0562: *         AAQQ  = AAQQ*TEMP1
0563: *         AAPP  = AAPP*TEMP1
0564:       ELSE
0565:          TEMP1 = ONE
0566:       END IF
0567: *
0568: *     Scale, if necessary
0569: *
0570:       IF ( TEMP1 .NE. ONE ) THEN
0571:          CALL SLASCL( 'G', 0, 0, ONE, TEMP1, N, 1, SVA, N, IERR )
0572:       END IF
0573:       SCALE = TEMP1 * SCALE
0574:       IF ( SCALE .NE. ONE ) THEN
0575:          CALL SLASCL( JOBA, 0, 0, ONE, SCALE, M, N, A, LDA, IERR )
0576:          SCALE = ONE / SCALE
0577:       END IF
0578: *
0579: *     Row-cyclic Jacobi SVD algorithm with column pivoting
0580: *
0581:       EMPTSW   = ( N * ( N - 1 ) ) / 2
0582:       NOTROT   = 0
0583:       FASTR(1) = ZERO
0584: *
0585: *     A is represented in factored form A = A * diag(WORK), where diag(WORK)
0586: *     is initialized to identity. WORK is updated during fast scaled
0587: *     rotations.
0588: *
0589:       DO 1868 q = 1, N
0590:          WORK(q) = ONE
0591:  1868 CONTINUE
0592: *
0593: *
0594:       SWBAND = 3
0595: *[TP] SWBAND is a tuning parameter [TP]. It is meaningful and effective
0596: *     if SGESVJ is used as a computational routine in the preconditioned
0597: *     Jacobi SVD algorithm SGESVJ. For sweeps i=1:SWBAND the procedure
0598: *     works on pivots inside a band-like region around the diagonal.
0599: *     The boundaries are determined dynamically, based on the number of
0600: *     pivots above a threshold.
0601: *
0602:       KBL = MIN0( 8, N )
0603: *[TP] KBL is a tuning parameter that defines the tile size in the
0604: *     tiling of the p-q loops of pivot pairs. In general, an optimal
0605: *     value of KBL depends on the matrix dimensions and on the
0606: *     parameters of the computer's memory.
0607: *
0608:       NBL = N / KBL
0609:       IF ( ( NBL * KBL ) .NE. N ) NBL = NBL + 1
0610: *
0611:       BLSKIP = KBL**2
0612: *[TP] BLKSKIP is a tuning parameter that depends on SWBAND and KBL.
0613: *
0614:       ROWSKIP = MIN0( 5, KBL )
0615: *[TP] ROWSKIP is a tuning parameter.
0616: *
0617:       LKAHEAD = 1
0618: *[TP] LKAHEAD is a tuning parameter.
0619: *
0620: *     Quasi block transformations, using the lower (upper) triangular
0621: *     structure of the input matrix. The quasi-block-cycling usually
0622: *     invokes cubic convergence. Big part of this cycle is done inside
0623: *     canonical subspaces of dimensions less than M.
0624: *
0625:       IF ( (LOWER .OR. UPPER) .AND. (N .GT. MAX0(64, 4*KBL)) ) THEN
0626: *[TP] The number of partition levels and the actual partition are
0627: *     tuning parameters.
0628:       N4     = N / 4
0629:       N2     = N / 2
0630:       N34    = 3 * N4
0631:       IF ( APPLV ) THEN
0632:          q = 0
0633:       ELSE
0634:          q = 1
0635:       END IF
0636: *
0637:       IF ( LOWER ) THEN
0638: *
0639: *     This works very well on lower triangular matrices, in particular
0640: *     in the framework of the preconditioned Jacobi SVD (xGEJSV).
0641: *     The idea is simple:
0642: *     [+ 0 0 0]   Note that Jacobi transformations of [0 0]
0643: *     [+ + 0 0]                                       [0 0]
0644: *     [+ + x 0]   actually work on [x 0]              [x 0]
0645: *     [+ + x x]                    [x x].             [x x]
0646: *
0647:       CALL SGSVJ0(JOBV,M-N34,N-N34,A(N34+1,N34+1),LDA,WORK(N34+1),
0648:      &     SVA(N34+1),MVL,V(N34*q+1,N34+1),LDV,EPSILON,SFMIN,TOL,2,
0649:      &     WORK(N+1),LWORK-N,IERR )
0650: *
0651:       CALL SGSVJ0( JOBV,M-N2,N34-N2,A(N2+1,N2+1),LDA,WORK(N2+1),
0652:      &     SVA(N2+1),MVL,V(N2*q+1,N2+1),LDV,EPSILON,SFMIN,TOL,2,
0653:      &     WORK(N+1),LWORK-N,IERR )
0654: *
0655:       CALL SGSVJ1( JOBV,M-N2,N-N2,N4,A(N2+1,N2+1),LDA,WORK(N2+1),
0656:      &     SVA(N2+1),MVL,V(N2*q+1,N2+1),LDV,EPSILON,SFMIN,TOL,1,
0657:      &     WORK(N+1),LWORK-N,IERR )
0658: *
0659:       CALL SGSVJ0( JOBV,M-N4,N2-N4,A(N4+1,N4+1),LDA,WORK(N4+1),
0660:      &     SVA(N4+1),MVL,V(N4*q+1,N4+1),LDV,EPSILON,SFMIN,TOL,1,
0661:      &     WORK(N+1),LWORK-N,IERR )
0662: *
0663:       CALL SGSVJ0( JOBV,M,N4,A,LDA,WORK,SVA,MVL,V,LDV,EPSILON,
0664:      &     SFMIN,TOL,1,WORK(N+1),LWORK-N,IERR )
0665: *
0666:       CALL SGSVJ1( JOBV,M,N2,N4,A,LDA,WORK,SVA,MVL,V,LDV,EPSILON,
0667:      &     SFMIN,TOL,1,WORK(N+1),LWORK-N,IERR )
0668: *
0669: *
0670:       ELSE IF ( UPPER ) THEN
0671: *
0672: *
0673:       CALL SGSVJ0( JOBV,N4,N4,A,LDA,WORK,SVA,MVL,V,LDV,EPSILON,
0674:      &     SFMIN,TOL,2,WORK(N+1),LWORK-N,IERR )
0675: *
0676:       CALL SGSVJ0(JOBV,N2,N4,A(1,N4+1),LDA,WORK(N4+1),SVA(N4+1),MVL,
0677:      &     V(N4*q+1,N4+1),LDV,EPSILON,SFMIN,TOL,1,WORK(N+1),LWORK-N,
0678:      &     IERR )
0679: *
0680:       CALL SGSVJ1( JOBV,N2,N2,N4,A,LDA,WORK,SVA,MVL,V,LDV,EPSILON,
0681:      &     SFMIN,TOL,1,WORK(N+1),LWORK-N,IERR )
0682: *
0683:       CALL SGSVJ0( JOBV,N2+N4,N4,A(1,N2+1),LDA,WORK(N2+1),SVA(N2+1),MVL,
0684:      &     V(N2*q+1,N2+1),LDV,EPSILON,SFMIN,TOL,1,
0685:      &     WORK(N+1),LWORK-N,IERR )
0686: 
0687:       END IF
0688: *
0689:       END IF
0690: *
0691: *     -#- Row-cyclic pivot strategy with de Rijk's pivoting -#-
0692: *
0693:       DO 1993 i = 1, NSWEEP
0694: *     .. go go go ...
0695: *
0696:       MXAAPQ = ZERO
0697:       MXSINJ = ZERO
0698:       ISWROT = 0
0699: *
0700:       NOTROT   = 0
0701:       PSKIPPED = 0
0702: *
0703: *     Each sweep is unrolled using KBL-by-KBL tiles over the pivot pairs
0704: *     1 <= p < q <= N. This is the first step toward a blocked implementation
0705: *     of the rotations. New implementation, based on block transformations,
0706: *     is under development.
0707: *
0708:       DO 2000 ibr = 1, NBL
0709: *
0710:       igl = ( ibr - 1 ) * KBL + 1
0711: *
0712:       DO 1002 ir1 = 0, MIN0( LKAHEAD, NBL - ibr )
0713: *
0714:       igl = igl + ir1 * KBL
0715: *
0716:       DO 2001 p = igl, MIN0( igl + KBL - 1, N - 1)
0717: *
0718: *     .. de Rijk's pivoting
0719: *
0720:       q   = ISAMAX( N-p+1, SVA(p), 1 ) + p - 1
0721:       IF ( p .NE. q ) THEN
0722:          CALL SSWAP( M, A(1,p), 1, A(1,q), 1 )
0723:          IF ( RSVEC ) CALL SSWAP( MVL, V(1,p), 1, V(1,q), 1 )
0724:          TEMP1   = SVA(p)
0725:          SVA(p)  = SVA(q)
0726:          SVA(q)  = TEMP1
0727:          TEMP1   = WORK(p)
0728:          WORK(p) = WORK(q)
0729:          WORK(q) = TEMP1
0730:       END IF
0731: *
0732:       IF ( ir1 .EQ. 0 ) THEN
0733: *
0734: *        Column norms are periodically updated by explicit
0735: *        norm computation.
0736: *        Caveat:
0737: *        Unfortunately, some BLAS implementations compute SNRM2(M,A(1,p),1)
0738: *        as SQRT(SDOT(M,A(1,p),1,A(1,p),1)), which may cause the result to
0739: *        overflow for ||A(:,p)||_2 > SQRT(overflow_threshold), and to
0740: *        underflow for ||A(:,p)||_2 < SQRT(underflow_threshold).
0741: *        Hence, SNRM2 cannot be trusted, not even in the case when
0742: *        the true norm is far from the under(over)flow boundaries.
0743: *        If properly implemented SNRM2 is available, the IF-THEN-ELSE
0744: *        below should read "AAPP = SNRM2( M, A(1,p), 1 ) * WORK(p)".
0745: *
0746:          IF ((SVA(p) .LT. ROOTBIG) .AND. (SVA(p) .GT. ROOTSFMIN)) THEN
0747:             SVA(p) = SNRM2( M, A(1,p), 1 ) * WORK(p)
0748:          ELSE
0749:             TEMP1 = ZERO
0750:             AAPP  = ZERO
0751:             CALL SLASSQ( M, A(1,p), 1, TEMP1, AAPP )
0752:             SVA(p) = TEMP1 * SQRT(AAPP) * WORK(p)
0753:          END IF
0754:          AAPP = SVA(p)
0755:       ELSE
0756:          AAPP = SVA(p)
0757:       END IF
0758: *
0759:       IF ( AAPP .GT. ZERO ) THEN
0760: *
0761:       PSKIPPED = 0
0762: *
0763:       DO 2002 q = p + 1, MIN0( igl + KBL - 1, N )
0764: *
0765:       AAQQ = SVA(q)
0766: *
0767:       IF ( AAQQ .GT. ZERO ) THEN
0768: *
0769:          AAPP0 = AAPP
0770:          IF ( AAQQ .GE. ONE ) THEN
0771:             ROTOK  = ( SMALL*AAPP ) .LE. AAQQ
0772:             IF ( AAPP .LT. ( BIG / AAQQ ) ) THEN
0773:                AAPQ = ( SDOT(M, A(1,p), 1, A(1,q), 1 ) *
0774:      &                  WORK(p) * WORK(q) / AAQQ ) / AAPP
0775:             ELSE
0776:                CALL SCOPY( M, A(1,p), 1, WORK(N+1), 1 )
0777:                CALL SLASCL( 'G', 0, 0, AAPP, WORK(p), M,
0778:      &              1, WORK(N+1), LDA, IERR )
0779:                AAPQ = SDOT( M, WORK(N+1),1, A(1,q),1 )*WORK(q) / AAQQ
0780:             END IF
0781:          ELSE
0782:             ROTOK  = AAPP .LE. ( AAQQ / SMALL )
0783:             IF ( AAPP .GT.  ( SMALL / AAQQ ) ) THEN
0784:                AAPQ = ( SDOT( M, A(1,p), 1, A(1,q), 1 ) *
0785:      &               WORK(p) * WORK(q) / AAQQ ) / AAPP
0786:             ELSE
0787:                CALL SCOPY( M, A(1,q), 1, WORK(N+1), 1 )
0788:                CALL SLASCL( 'G', 0, 0, AAQQ, WORK(q), M,
0789:      &              1, WORK(N+1), LDA, IERR )
0790:                AAPQ = SDOT( M, WORK(N+1),1, A(1,p),1 )*WORK(p) / AAPP
0791:             END IF
0792:          END IF
0793: *
0794:          MXAAPQ = AMAX1( MXAAPQ, ABS(AAPQ) )
0795: *
0796: *        TO rotate or NOT to rotate, THAT is the question ...
0797: *
0798:          IF ( ABS( AAPQ ) .GT. TOL ) THEN
0799: *
0800: *           .. rotate
0801: *[RTD]      ROTATED = ROTATED + ONE
0802: *
0803:             IF ( ir1 .EQ. 0 ) THEN
0804:                NOTROT   = 0
0805:                PSKIPPED = 0
0806:                ISWROT   = ISWROT  + 1
0807:             END IF
0808: *
0809:             IF ( ROTOK ) THEN
0810: *
0811:                AQOAP = AAQQ / AAPP
0812:                APOAQ = AAPP / AAQQ
0813:                THETA = - HALF * ABS( AQOAP - APOAQ ) / AAPQ
0814: *
0815:                IF ( ABS( THETA ) .GT. BIGTHETA ) THEN
0816: *
0817:                   T        = HALF / THETA
0818:                   FASTR(3) =   T * WORK(p) / WORK(q)
0819:                   FASTR(4) = - T * WORK(q) / WORK(p)
0820:                   CALL SROTM( M,   A(1,p), 1, A(1,q), 1, FASTR )
0821:                   IF ( RSVEC )
0822:      &            CALL SROTM( MVL, V(1,p), 1, V(1,q), 1, FASTR )
0823:                   SVA(q) = AAQQ*SQRT( AMAX1(ZERO,ONE + T*APOAQ*AAPQ) )
0824:                   AAPP   = AAPP*SQRT( ONE - T*AQOAP*AAPQ )
0825:                   MXSINJ = AMAX1( MXSINJ, ABS(T) )
0826: *
0827:                ELSE
0828: *
0829: *                 .. choose correct signum for THETA and rotate
0830: *
0831:                   THSIGN =  - SIGN(ONE,AAPQ)
0832:                   T  = ONE / ( THETA + THSIGN*SQRT(ONE+THETA*THETA) )
0833:                   CS = SQRT( ONE / ( ONE + T*T ) )
0834:                   SN = T * CS
0835: *
0836:                   MXSINJ = AMAX1( MXSINJ, ABS(SN) )
0837:                   SVA(q) = AAQQ*SQRT( AMAX1(ZERO, ONE+T*APOAQ*AAPQ) )
0838:                   AAPP   = AAPP*SQRT( AMAX1(ZERO, ONE-T*AQOAP*AAPQ) )
0839: *
0840:                   APOAQ = WORK(p) / WORK(q)
0841:                   AQOAP = WORK(q) / WORK(p)
0842:                   IF ( WORK(p) .GE. ONE ) THEN
0843:                      IF ( WORK(q) .GE.  ONE ) THEN
0844:                         FASTR(3) =   T * APOAQ
0845:                         FASTR(4) = - T * AQOAP
0846:                         WORK(p)  = WORK(p) * CS
0847:                         WORK(q)  = WORK(q) * CS
0848:                         CALL SROTM( M,   A(1,p),1, A(1,q),1, FASTR )
0849:                         IF ( RSVEC )
0850:      &                  CALL SROTM( MVL, V(1,p),1, V(1,q),1, FASTR )
0851:                      ELSE
0852:                         CALL SAXPY( M,    -T*AQOAP, A(1,q),1, A(1,p),1 )
0853:                         CALL SAXPY( M, CS*SN*APOAQ, A(1,p),1, A(1,q),1 )
0854:                         WORK(p) = WORK(p) * CS
0855:                         WORK(q) = WORK(q) / CS
0856:                         IF ( RSVEC ) THEN
0857:                         CALL SAXPY(MVL,   -T*AQOAP, V(1,q),1,V(1,p),1)
0858:                         CALL SAXPY(MVL,CS*SN*APOAQ, V(1,p),1,V(1,q),1)
0859:                         END IF
0860:                      END IF
0861:                   ELSE
0862:                      IF ( WORK(q) .GE. ONE ) THEN
0863:                         CALL SAXPY( M,     T*APOAQ, A(1,p),1, A(1,q),1 )
0864:                         CALL SAXPY( M,-CS*SN*AQOAP, A(1,q),1, A(1,p),1 )
0865:                         WORK(p) = WORK(p) / CS
0866:                         WORK(q) = WORK(q) * CS
0867:                         IF ( RSVEC ) THEN
0868:                         CALL SAXPY(MVL, T*APOAQ,    V(1,p),1,V(1,q),1)
0869:                         CALL SAXPY(MVL,-CS*SN*AQOAP,V(1,q),1,V(1,p),1)
0870:                         END IF
0871:                      ELSE
0872:                         IF ( WORK(p) .GE. WORK(q) ) THEN
0873:                            CALL SAXPY( M,-T*AQOAP,   A(1,q),1,A(1,p),1 )
0874:                            CALL SAXPY( M,CS*SN*APOAQ,A(1,p),1,A(1,q),1 )
0875:                            WORK(p) = WORK(p) * CS
0876:                            WORK(q) = WORK(q) / CS
0877:                            IF ( RSVEC ) THEN
0878:                            CALL SAXPY(MVL, -T*AQOAP,  V(1,q),1,V(1,p),1)
0879:                            CALL SAXPY(MVL,CS*SN*APOAQ,V(1,p),1,V(1,q),1)
0880:                            END IF
0881:                         ELSE
0882:                            CALL SAXPY( M, T*APOAQ,    A(1,p),1,A(1,q),1)
0883:                            CALL SAXPY( M,-CS*SN*AQOAP,A(1,q),1,A(1,p),1)
0884:                            WORK(p) = WORK(p) / CS
0885:                            WORK(q) = WORK(q) * CS
0886:                           IF ( RSVEC ) THEN
0887:                           CALL SAXPY(MVL, T*APOAQ,    V(1,p),1,V(1,q),1)
0888:                           CALL SAXPY(MVL,-CS*SN*AQOAP,V(1,q),1,V(1,p),1)
0889:                           END IF
0890:                         END IF
0891:                      END IF
0892:                   ENDIF
0893:                END IF
0894: *
0895:             ELSE
0896: *              .. have to use modified Gram-Schmidt like transformation
0897:                CALL SCOPY( M, A(1,p), 1, WORK(N+1), 1 )
0898:                CALL SLASCL( 'G',0,0,AAPP,ONE,M,1,WORK(N+1),LDA,IERR )
0899:                CALL SLASCL( 'G',0,0,AAQQ,ONE,M,1,   A(1,q),LDA,IERR )
0900:                TEMP1 = -AAPQ * WORK(p) / WORK(q)
0901:                CALL SAXPY ( M, TEMP1, WORK(N+1), 1, A(1,q), 1 )
0902:                CALL SLASCL( 'G',0,0,ONE,AAQQ,M,1,   A(1,q),LDA,IERR )
0903:                SVA(q) = AAQQ*SQRT( AMAX1( ZERO, ONE - AAPQ*AAPQ ) )
0904:                MXSINJ = AMAX1( MXSINJ, SFMIN )
0905:             END IF
0906: *           END IF ROTOK THEN ... ELSE
0907: *
0908: *           In the case of cancellation in updating SVA(q), SVA(p)
0909: *           recompute SVA(q), SVA(p).
0910: *
0911:             IF ( (SVA(q) / AAQQ )**2 .LE. ROOTEPS  ) THEN
0912:                IF ((AAQQ .LT. ROOTBIG).AND.(AAQQ .GT. ROOTSFMIN)) THEN
0913:                   SVA(q) = SNRM2( M, A(1,q), 1 ) * WORK(q)
0914:                ELSE
0915:                   T    = ZERO
0916:                   AAQQ = ZERO
0917:                   CALL SLASSQ( M, A(1,q), 1, T, AAQQ )
0918:                   SVA(q) = T * SQRT(AAQQ) * WORK(q)
0919:                END IF
0920:             END IF
0921:             IF ( ( AAPP / AAPP0) .LE. ROOTEPS  ) THEN
0922:                IF ((AAPP .LT. ROOTBIG).AND.(AAPP .GT. ROOTSFMIN)) THEN
0923:                   AAPP = SNRM2( M, A(1,p), 1 ) * WORK(p)
0924:                ELSE
0925:                   T    = ZERO
0926:                   AAPP = ZERO
0927:                   CALL SLASSQ( M, A(1,p), 1, T, AAPP )
0928:                   AAPP = T * SQRT(AAPP) * WORK(p)
0929:                END IF
0930:                SVA(p) = AAPP
0931:             END IF
0932: *
0933:          ELSE
0934: *        A(:,p) and A(:,q) already numerically orthogonal
0935:             IF ( ir1 .EQ. 0 ) NOTROT   = NOTROT + 1
0936: *[RTD]      SKIPPED  = SKIPPED  + 1
0937:             PSKIPPED = PSKIPPED + 1
0938:          END IF
0939:       ELSE
0940: *        A(:,q) is zero column
0941:          IF ( ir1. EQ. 0 ) NOTROT = NOTROT + 1
0942:          PSKIPPED = PSKIPPED + 1
0943:       END IF
0944: *
0945:       IF ( ( i .LE. SWBAND ) .AND. ( PSKIPPED .GT. ROWSKIP ) ) THEN
0946:          IF ( ir1 .EQ. 0 ) AAPP = - AAPP
0947:          NOTROT = 0
0948:          GO TO 2103
0949:       END IF
0950: *
0951:  2002 CONTINUE
0952: *     END q-LOOP
0953: *
0954:  2103 CONTINUE
0955: *     bailed out of q-loop
0956: *
0957:       SVA(p) = AAPP
0958: *
0959:       ELSE
0960:          SVA(p) = AAPP
0961:          IF ( ( ir1 .EQ. 0 ) .AND. (AAPP .EQ. ZERO) )
0962:      &        NOTROT=NOTROT+MIN0(igl+KBL-1,N)-p
0963:       END IF
0964: *
0965:  2001 CONTINUE
0966: *     end of the p-loop
0967: *     end of doing the block ( ibr, ibr )
0968:  1002 CONTINUE
0969: *     end of ir1-loop
0970: *
0971: * ... go to the off diagonal blocks
0972: *
0973:       igl = ( ibr - 1 ) * KBL + 1
0974: *
0975:       DO 2010 jbc = ibr + 1, NBL
0976: *
0977:          jgl = ( jbc - 1 ) * KBL + 1
0978: *
0979: *        doing the block at ( ibr, jbc )
0980: *
0981:          IJBLSK = 0
0982:          DO 2100 p = igl, MIN0( igl + KBL - 1, N )
0983: *
0984:          AAPP = SVA(p)
0985:          IF ( AAPP .GT. ZERO ) THEN
0986: *
0987:          PSKIPPED = 0
0988: *
0989:          DO 2200 q = jgl, MIN0( jgl + KBL - 1, N )
0990: *
0991:          AAQQ = SVA(q)
0992:          IF ( AAQQ .GT. ZERO ) THEN
0993:             AAPP0 = AAPP
0994: *
0995: *     -#- M x 2 Jacobi SVD -#-
0996: *
0997: *        Safe Gram matrix computation
0998: *
0999:          IF ( AAQQ .GE. ONE ) THEN
1000:             IF ( AAPP .GE. AAQQ ) THEN
1001:                ROTOK = ( SMALL*AAPP ) .LE. AAQQ
1002:             ELSE
1003:                ROTOK = ( SMALL*AAQQ ) .LE. AAPP
1004:             END IF
1005:             IF ( AAPP .LT. ( BIG / AAQQ ) ) THEN
1006:                AAPQ = ( SDOT(M, A(1,p), 1, A(1,q), 1 ) *
1007:      &                  WORK(p) * WORK(q) / AAQQ ) / AAPP
1008:             ELSE
1009:                CALL SCOPY( M, A(1,p), 1, WORK(N+1), 1 )
1010:                CALL SLASCL( 'G', 0, 0, AAPP, WORK(p), M,
1011:      &              1, WORK(N+1), LDA, IERR )
1012:                AAPQ = SDOT( M, WORK(N+1), 1, A(1,q), 1 ) *
1013:      &                WORK(q) / AAQQ
1014:             END IF
1015:          ELSE
1016:             IF ( AAPP .GE. AAQQ ) THEN
1017:                ROTOK = AAPP .LE. ( AAQQ / SMALL )
1018:             ELSE
1019:                ROTOK = AAQQ .LE. ( AAPP / SMALL )
1020:             END IF
1021:             IF ( AAPP .GT.  ( SMALL / AAQQ ) ) THEN
1022:                AAPQ = ( SDOT( M, A(1,p), 1, A(1,q), 1 ) *
1023:      &               WORK(p) * WORK(q) / AAQQ ) / AAPP
1024:             ELSE
1025:                CALL SCOPY( M, A(1,q), 1, WORK(N+1), 1 )
1026:                CALL SLASCL( 'G', 0, 0, AAQQ, WORK(q), M, 1,
1027:      &              WORK(N+1), LDA, IERR )
1028:                AAPQ = SDOT(M,WORK(N+1),1,A(1,p),1) * WORK(p) / AAPP
1029:             END IF
1030:          END IF
1031: *
1032:          MXAAPQ = AMAX1( MXAAPQ, ABS(AAPQ) )
1033: *
1034: *        TO rotate or NOT to rotate, THAT is the question ...
1035: *
1036:          IF ( ABS( AAPQ ) .GT. TOL ) THEN
1037:             NOTROT   = 0
1038: *[RTD]      ROTATED  = ROTATED + 1
1039:             PSKIPPED = 0
1040:             ISWROT   = ISWROT  + 1
1041: *
1042:             IF ( ROTOK ) THEN
1043: *
1044:                AQOAP = AAQQ / AAPP
1045:                APOAQ = AAPP / AAQQ
1046:                THETA = - HALF * ABS( AQOAP - APOAQ ) / AAPQ
1047:                IF ( AAQQ .GT. AAPP0 ) THETA = - THETA
1048: *
1049:                IF ( ABS( THETA ) .GT. BIGTHETA ) THEN
1050:                   T = HALF / THETA
1051:                   FASTR(3) =  T * WORK(p) / WORK(q)
1052:                   FASTR(4) = -T * WORK(q) / WORK(p)
1053:                   CALL SROTM( M,   A(1,p), 1, A(1,q), 1, FASTR )
1054:                   IF ( RSVEC )
1055:      &            CALL SROTM( MVL, V(1,p), 1, V(1,q), 1, FASTR )
1056:                   SVA(q) = AAQQ*SQRT( AMAX1(ZERO,ONE + T*APOAQ*AAPQ) )
1057:                   AAPP   = AAPP*SQRT( AMAX1(ZERO,ONE - T*AQOAP*AAPQ) )
1058:                   MXSINJ = AMAX1( MXSINJ, ABS(T) )
1059:                ELSE
1060: *
1061: *                 .. choose correct signum for THETA and rotate
1062: *
1063:                   THSIGN = - SIGN(ONE,AAPQ)
1064:                   IF ( AAQQ .GT. AAPP0 ) THSIGN = - THSIGN
1065:                   T  = ONE / ( THETA + THSIGN*SQRT(ONE+THETA*THETA) )
1066:                   CS = SQRT( ONE / ( ONE + T*T ) )
1067:                   SN = T * CS
1068:                   MXSINJ = AMAX1( MXSINJ, ABS(SN) )
1069:                   SVA(q) = AAQQ*SQRT( AMAX1(ZERO, ONE+T*APOAQ*AAPQ) )
1070:                   AAPP   = AAPP*SQRT( ONE - T*AQOAP*AAPQ)
1071: *
1072:                   APOAQ = WORK(p) / WORK(q)
1073:                   AQOAP = WORK(q) / WORK(p)
1074:                   IF ( WORK(p) .GE. ONE ) THEN
1075: *
1076:                      IF ( WORK(q) .GE.  ONE ) THEN
1077:                         FASTR(3) =   T * APOAQ
1078:                         FASTR(4) = - T * AQOAP
1079:                         WORK(p)  = WORK(p) * CS
1080:                         WORK(q)  = WORK(q) * CS
1081:                         CALL SROTM( M,   A(1,p),1, A(1,q),1, FASTR )
1082:                         IF ( RSVEC )
1083:      &                  CALL SROTM( MVL, V(1,p),1, V(1,q),1, FASTR )
1084:                      ELSE
1085:                         CALL SAXPY( M,    -T*AQOAP, A(1,q),1, A(1,p),1 )
1086:                         CALL SAXPY( M, CS*SN*APOAQ, A(1,p),1, A(1,q),1 )
1087:                         IF ( RSVEC ) THEN
1088:                         CALL SAXPY( MVL, -T*AQOAP,  V(1,q),1, V(1,p),1 )
1089:                         CALL SAXPY( MVL,CS*SN*APOAQ,V(1,p),1, V(1,q),1 )
1090:                         END IF
1091:                         WORK(p) = WORK(p) * CS
1092:                         WORK(q) = WORK(q) / CS
1093:                      END IF
1094:                   ELSE
1095:                      IF ( WORK(q) .GE. ONE ) THEN
1096:                         CALL SAXPY( M,     T*APOAQ, A(1,p),1, A(1,q),1 )
1097:                         CALL SAXPY( M,-CS*SN*AQOAP, A(1,q),1, A(1,p),1 )
1098:                         IF ( RSVEC ) THEN
1099:                         CALL SAXPY(MVL,T*APOAQ,     V(1,p),1, V(1,q),1 )
1100:                         CALL SAXPY(MVL,-CS*SN*AQOAP,V(1,q),1, V(1,p),1 )
1101:                         END IF
1102:                         WORK(p) = WORK(p) / CS
1103:                         WORK(q) = WORK(q) * CS
1104:                      ELSE
1105:                         IF ( WORK(p) .GE. WORK(q) ) THEN
1106:                            CALL SAXPY( M,-T*AQOAP,   A(1,q),1,A(1,p),1 )
1107:                            CALL SAXPY( M,CS*SN*APOAQ,A(1,p),1,A(1,q),1 )
1108:                            WORK(p) = WORK(p) * CS
1109:                            WORK(q) = WORK(q) / CS
1110:                            IF ( RSVEC ) THEN
1111:                            CALL SAXPY( MVL, -T*AQOAP, V(1,q),1,V(1,p),1)
1112:                            CALL SAXPY(MVL,CS*SN*APOAQ,V(1,p),1,V(1,q),1)
1113:                            END IF
1114:                         ELSE
1115:                            CALL SAXPY(M, T*APOAQ,    A(1,p),1,A(1,q),1)
1116:                            CALL SAXPY(M,-CS*SN*AQOAP,A(1,q),1,A(1,p),1)
1117:                            WORK(p) = WORK(p) / CS
1118:                            WORK(q) = WORK(q) * CS
1119:                           IF ( RSVEC ) THEN
1120:                           CALL SAXPY(MVL, T*APOAQ,    V(1,p),1,V(1,q),1)
1121:                           CALL SAXPY(MVL,-CS*SN*AQOAP,V(1,q),1,V(1,p),1)
1122:                           END IF
1123:                         END IF
1124:                      END IF
1125:                   ENDIF
1126:                END IF
1127: *
1128:             ELSE
1129:                IF ( AAPP .GT. AAQQ ) THEN
1130:                   CALL SCOPY( M, A(1,p), 1, WORK(N+1), 1 )
1131:                   CALL SLASCL('G',0,0,AAPP,ONE,M,1,WORK(N+1),LDA,IERR)
1132:                   CALL SLASCL('G',0,0,AAQQ,ONE,M,1,   A(1,q),LDA,IERR)
1133:                   TEMP1 = -AAPQ * WORK(p) / WORK(q)
1134:                   CALL SAXPY(M,TEMP1,WORK(N+1),1,A(1,q),1)
1135:                   CALL SLASCL('G',0,0,ONE,AAQQ,M,1,A(1,q),LDA,IERR)
1136:                   SVA(q) = AAQQ*SQRT(AMAX1(ZERO, ONE - AAPQ*AAPQ))
1137:                   MXSINJ = AMAX1( MXSINJ, SFMIN )
1138:                ELSE
1139:                   CALL SCOPY( M, A(1,q), 1, WORK(N+1), 1 )
1140:                   CALL SLASCL('G',0,0,AAQQ,ONE,M,1,WORK(N+1),LDA,IERR)
1141:                   CALL SLASCL('G',0,0,AAPP,ONE,M,1,   A(1,p),LDA,IERR)
1142:                   TEMP1 = -AAPQ * WORK(q) / WORK(p)
1143:                   CALL SAXPY(M,TEMP1,WORK(N+1),1,A(1,p),1)
1144:                   CALL SLASCL('G',0,0,ONE,AAPP,M,1,A(1,p),LDA,IERR)
1145:                   SVA(p) = AAPP*SQRT(AMAX1(ZERO, ONE - AAPQ*AAPQ))
1146:                   MXSINJ = AMAX1( MXSINJ, SFMIN )
1147:                END IF
1148:             END IF
1149: *           END IF ROTOK THEN ... ELSE
1150: *
1151: *           In the case of cancellation in updating SVA(q)
1152: *           .. recompute SVA(q)
1153:             IF ( (SVA(q) / AAQQ )**2 .LE. ROOTEPS  ) THEN
1154:                IF ((AAQQ .LT. ROOTBIG).AND.(AAQQ .GT. ROOTSFMIN)) THEN
1155:                   SVA(q) = SNRM2( M, A(1,q), 1 ) * WORK(q)
1156:                ELSE
1157:                   T    = ZERO
1158:                   AAQQ = ZERO
1159:                   CALL SLASSQ( M, A(1,q), 1, T, AAQQ )
1160:                   SVA(q) = T * SQRT(AAQQ) * WORK(q)
1161:                END IF
1162:             END IF
1163:             IF ( (AAPP / AAPP0 )**2 .LE. ROOTEPS  ) THEN
1164:                IF ((AAPP .LT. ROOTBIG).AND.(AAPP .GT. ROOTSFMIN)) THEN
1165:                   AAPP = SNRM2( M, A(1,p), 1 ) * WORK(p)
1166:                ELSE
1167:                   T    = ZERO
1168:                   AAPP = ZERO
1169:                   CALL SLASSQ( M, A(1,p), 1, T, AAPP )
1170:                   AAPP = T * SQRT(AAPP) * WORK(p)
1171:                END IF
1172:                SVA(p) = AAPP
1173:             END IF
1174: *              end of OK rotation
1175:          ELSE
1176:             NOTROT   = NOTROT   + 1
1177: *[RTD]      SKIPPED  = SKIPPED  + 1
1178:             PSKIPPED = PSKIPPED + 1
1179:             IJBLSK   = IJBLSK   + 1
1180:          END IF
1181:       ELSE
1182:          NOTROT   = NOTROT   + 1
1183:          PSKIPPED = PSKIPPED + 1
1184:          IJBLSK   = IJBLSK   + 1
1185:       END IF
1186: *
1187:       IF ( ( i .LE. SWBAND ) .AND. ( IJBLSK .GE. BLSKIP ) ) THEN
1188:          SVA(p) = AAPP
1189:          NOTROT = 0
1190:          GO TO 2011
1191:       END IF
1192:       IF ( ( i .LE. SWBAND ) .AND. ( PSKIPPED .GT. ROWSKIP ) ) THEN
1193:          AAPP = -AAPP
1194:          NOTROT = 0
1195:          GO TO 2203
1196:       END IF
1197: *
1198:  2200    CONTINUE
1199: *        end of the q-loop
1200:  2203    CONTINUE
1201: *
1202:          SVA(p) = AAPP
1203: *
1204:       ELSE
1205: *
1206:          IF ( AAPP .EQ. ZERO ) NOTROT=NOTROT+MIN0(jgl+KBL-1,N)-jgl+1
1207:          IF ( AAPP .LT. ZERO ) NOTROT = 0
1208: *
1209:       END IF
1210: *
1211:  2100 CONTINUE
1212: *     end of the p-loop
1213:  2010 CONTINUE
1214: *     end of the jbc-loop
1215:  2011 CONTINUE
1216: *2011 bailed out of the jbc-loop
1217:       DO 2012 p = igl, MIN0( igl + KBL - 1, N )
1218:          SVA(p) = ABS(SVA(p))
1219:  2012 CONTINUE
1220: ***
1221:  2000 CONTINUE
1222: *2000 :: end of the ibr-loop
1223: *
1224: *     .. update SVA(N)
1225:       IF ((SVA(N) .LT. ROOTBIG).AND.(SVA(N) .GT. ROOTSFMIN)) THEN
1226:          SVA(N) = SNRM2( M, A(1,N), 1 ) * WORK(N)
1227:       ELSE
1228:          T    = ZERO
1229:          AAPP = ZERO
1230:          CALL SLASSQ( M, A(1,N), 1, T, AAPP )
1231:          SVA(N) = T * SQRT(AAPP) * WORK(N)
1232:       END IF
1233: *
1234: *     Additional steering devices
1235: *
1236:       IF ( ( i.LT.SWBAND ) .AND. ( ( MXAAPQ.LE.ROOTTOL ) .OR.
1237:      &     ( ISWROT .LE. N ) ) )
1238:      &   SWBAND = i
1239: *
1240:       IF ( (i .GT. SWBAND+1) .AND. (MXAAPQ .LT. SQRT(FLOAT(N))*TOL)
1241:      &   .AND. (FLOAT(N)*MXAAPQ*MXSINJ .LT. TOL) ) THEN
1242:         GO TO 1994
1243:       END IF
1244: *
1245:       IF ( NOTROT .GE. EMPTSW ) GO TO 1994
1246: *
1247:  1993 CONTINUE
1248: *     end i=1:NSWEEP loop
1249: *
1250: * #:( Reaching this point means that the procedure has not converged.
1251:       INFO = NSWEEP - 1
1252:       GO TO 1995
1253: *
1254:  1994 CONTINUE
1255: * #:) Reaching this point means numerical convergence after the i-th
1256: *     sweep.
1257: *
1258:       INFO = 0
1259: * #:) INFO = 0 confirms successful iterations.
1260:  1995 CONTINUE
1261: *
1262: *     Sort the singular values and find how many are above
1263: *     the underflow threshold.
1264: *
1265:       N2 = 0
1266:       N4 = 0
1267:       DO 5991 p = 1, N - 1
1268:          q = ISAMAX( N-p+1, SVA(p), 1 ) + p - 1
1269:          IF ( p .NE. q ) THEN
1270:             TEMP1  = SVA(p)
1271:             SVA(p) = SVA(q)
1272:             SVA(q) = TEMP1
1273:             TEMP1   = WORK(p)
1274:             WORK(p) = WORK(q)
1275:             WORK(q) = TEMP1
1276:             CALL SSWAP( M, A(1,p), 1, A(1,q), 1 )
1277:             IF ( RSVEC ) CALL SSWAP( MVL, V(1,p), 1, V(1,q), 1 )
1278:          END IF
1279:          IF ( SVA(p) .NE. ZERO ) THEN
1280:             N4 = N4 + 1
1281:             IF ( SVA(p)*SCALE .GT. SFMIN ) N2 = N2 + 1
1282:          END IF
1283:  5991 CONTINUE
1284:       IF ( SVA(N) .NE. ZERO ) THEN
1285:          N4 = N4 + 1
1286:          IF ( SVA(N)*SCALE .GT. SFMIN ) N2 = N2 + 1
1287:       END IF
1288: *
1289: *     Normalize the left singular vectors.
1290: *
1291:       IF ( LSVEC .OR. UCTOL ) THEN
1292:          DO 1998 p = 1, N2
1293:             CALL SSCAL( M, WORK(p) / SVA(p), A(1,p), 1 )
1294:  1998    CONTINUE
1295:       END IF
1296: *
1297: *     Scale the product of Jacobi rotations (assemble the fast rotations).
1298: *
1299:       IF ( RSVEC ) THEN
1300:          IF ( APPLV ) THEN
1301:             DO 2398 p = 1, N
1302:                CALL SSCAL( MVL, WORK(p), V(1,p), 1 )
1303:  2398       CONTINUE
1304:          ELSE
1305:             DO 2399 p = 1, N
1306:                TEMP1 = ONE / SNRM2(MVL, V(1,p), 1 )
1307:                CALL SSCAL( MVL, TEMP1, V(1,p), 1 )
1308:  2399       CONTINUE
1309:          END IF
1310:       END IF
1311: *
1312: *     Undo scaling, if necessary (and possible).
1313:       IF ( ((SCALE.GT.ONE).AND.(SVA(1).LT.(BIG/SCALE)))
1314:      & .OR.((SCALE.LT.ONE).AND.(SVA(N2).GT.(SFMIN/SCALE))) ) THEN
1315:          DO 2400 p = 1, N
1316:             SVA(p) = SCALE*SVA(p)
1317:  2400    CONTINUE
1318:          SCALE = ONE
1319:       END IF
1320: *
1321:       WORK(1) = SCALE
1322: *     The singular values of A are SCALE*SVA(1:N). If SCALE.NE.ONE
1323: *     then some of the singular values may overflow or underflow and
1324: *     the spectrum is given in this factored representation.
1325: *
1326:       WORK(2) = FLOAT(N4)
1327: *     N4 is the number of computed nonzero singular values of A.
1328: *
1329:       WORK(3) = FLOAT(N2)
1330: *     N2 is the number of singular values of A greater than SFMIN.
1331: *     If N2<N, SVA(N2:N) contains ZEROS and/or denormalized numbers
1332: *     that may carry some information.
1333: *
1334:       WORK(4) = FLOAT(i)
1335: *     i is the index of the last sweep before declaring convergence.
1336: *
1337:       WORK(5) = MXAAPQ
1338: *     MXAAPQ is the largest absolute value of scaled pivots in the
1339: *     last sweep
1340: *
1341:       WORK(6) = MXSINJ
1342: *     MXSINJ is the largest absolute value of the sines of Jacobi angles
1343: *     in the last sweep
1344: *
1345:       RETURN
1346: *     ..
1347: *     .. END OF SGESVJ
1348: *     ..
1349:       END
1350: *
1351: