0001:       SUBROUTINE SGEJSV( JOBA, JOBU, JOBV, JOBR, JOBT, JOBP,
0002:      &                   M, N, A, LDA, SVA, U, LDU, V, LDV,
0003:      &                   WORK, LWORK, IWORK, INFO )
0004: *
0005: *  -- LAPACK routine (version 3.2.1)                                    --
0006: *
0007: *  -- Contributed by Zlatko Drmac of the University of Zagreb and     --
0008: *  -- Kresimir Veselic of the Fernuniversitaet Hagen                  --
0009: *  -- April 2009                                                      --
0010: *
0011: *  -- LAPACK is a software package provided by Univ. of Tennessee,    --
0012: *  -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
0013: *
0014: * This routine is also part of SIGMA (version 1.23, October 23. 2008.)
0015: * SIGMA is a library of algorithms for highly accurate algorithms for
0016: * computation of SVD, PSVD, QSVD, (H,K)-SVD, and for solution of the
0017: * eigenvalue problems Hx = lambda M x, H M x = lambda x with H, M > 0.
0018: *
0019: *     .. Scalar Arguments ..
0020: *
0021:       IMPLICIT    NONE
0022:       INTEGER     INFO, LDA, LDU, LDV, LWORK, M, N
0023: *
0024: *     .. Array Arguments ..
0025: *
0026:       REAL        A( LDA, * ), SVA( N ), U( LDU, * ), V( LDV, * ),
0027:      &            WORK( LWORK )
0028:       INTEGER     IWORK( * )
0029:       CHARACTER*1 JOBA, JOBP, JOBR, JOBT, JOBU, JOBV
0030: *     ..
0031: *
0032: *  Purpose
0033: *  =======
0034: *  SGEJSV computes the singular value decomposition (SVD) of a real M-by-N
0035: *  matrix [A], where M >= N. The SVD of [A] is written as
0036: *
0037: *               [A] = [U] * [SIGMA] * [V]^t,
0038: *
0039: *  where [SIGMA] is an N-by-N (M-by-N) matrix which is zero except for its N
0040: *  diagonal elements, [U] is an M-by-N (or M-by-M) orthonormal matrix, and
0041: *  [V] is an N-by-N orthogonal matrix. The diagonal elements of [SIGMA] are
0042: *  the singular values of [A]. The columns of [U] and [V] are the left and
0043: *  the right singular vectors of [A], respectively. The matrices [U] and [V]
0044: *  are computed and stored in the arrays U and V, respectively. The diagonal
0045: *  of [SIGMA] is computed and stored in the array SVA.
0046: *
0047: *  Arguments
0048: *  =========
0049: *
0050: *  JOBA   (input) CHARACTER*1
0051: *         Specifies the level of accuracy:
0052: *       = 'C': This option works well (high relative accuracy) if A = B * D,
0053: *              with well-conditioned B and arbitrary diagonal matrix D.
0054: *              The accuracy cannot be spoiled by COLUMN scaling. The
0055: *              accuracy of the computed output depends on the condition of
0056: *              B, and the procedure aims at the best theoretical accuracy.
0057: *              The relative error max_{i=1:N}|d sigma_i| / sigma_i is
0058: *              bounded by f(M,N)*epsilon* cond(B), independent of D.
0059: *              The input matrix is preprocessed with the QRF with column
0060: *              pivoting. This initial preprocessing and preconditioning by
0061: *              a rank revealing QR factorization is common for all values of
0062: *              JOBA. Additional actions are specified as follows:
0063: *       = 'E': Computation as with 'C' with an additional estimate of the
0064: *              condition number of B. It provides a realistic error bound.
0065: *       = 'F': If A = D1 * C * D2 with ill-conditioned diagonal scalings
0066: *              D1, D2, and well-conditioned matrix C, this option gives
0067: *              higher accuracy than the 'C' option. If the structure of the
0068: *              input matrix is not known, and relative accuracy is
0069: *              desirable, then this option is advisable. The input matrix A
0070: *              is preprocessed with QR factorization with FULL (row and
0071: *              column) pivoting.
0072: *       = 'G'  Computation as with 'F' with an additional estimate of the
0073: *              condition number of B, where A=D*B. If A has heavily weighted
0074: *              rows, then using this condition number gives too pessimistic
0075: *              error bound.
0076: *       = 'A': Small singular values are the noise and the matrix is treated
0077: *              as numerically rank defficient. The error in the computed
0078: *              singular values is bounded by f(m,n)*epsilon*||A||.
0079: *              The computed SVD A = U * S * V^t restores A up to
0080: *              f(m,n)*epsilon*||A||.
0081: *              This gives the procedure the licence to discard (set to zero)
0082: *              all singular values below N*epsilon*||A||.
0083: *       = 'R': Similar as in 'A'. Rank revealing property of the initial
0084: *              QR factorization is used do reveal (using triangular factor)
0085: *              a gap sigma_{r+1} < epsilon * sigma_r in which case the
0086: *              numerical RANK is declared to be r. The SVD is computed with
0087: *              absolute error bounds, but more accurately than with 'A'.
0088: * 
0089: *  JOBU   (input) CHARACTER*1
0090: *         Specifies whether to compute the columns of U:
0091: *       = 'U': N columns of U are returned in the array U.
0092: *       = 'F': full set of M left sing. vectors is returned in the array U.
0093: *       = 'W': U may be used as workspace of length M*N. See the description
0094: *              of U.
0095: *       = 'N': U is not computed.
0096: * 
0097: *  JOBV   (input) CHARACTER*1
0098: *         Specifies whether to compute the matrix V:
0099: *       = 'V': N columns of V are returned in the array V; Jacobi rotations
0100: *              are not explicitly accumulated.
0101: *       = 'J': N columns of V are returned in the array V, but they are
0102: *              computed as the product of Jacobi rotations. This option is
0103: *              allowed only if JOBU .NE. 'N', i.e. in computing the full SVD.
0104: *       = 'W': V may be used as workspace of length N*N. See the description
0105: *              of V.
0106: *       = 'N': V is not computed.
0107: * 
0108: *  JOBR   (input) CHARACTER*1
0109: *         Specifies the RANGE for the singular values. Issues the licence to
0110: *         set to zero small positive singular values if they are outside
0111: *         specified range. If A .NE. 0 is scaled so that the largest singular
0112: *         value of c*A is around SQRT(BIG), BIG=SLAMCH('O'), then JOBR issues
0113: *         the licence to kill columns of A whose norm in c*A is less than
0114: *         SQRT(SFMIN) (for JOBR.EQ.'R'), or less than SMALL=SFMIN/EPSLN,
0115: *         where SFMIN=SLAMCH('S'), EPSLN=SLAMCH('E').
0116: *       = 'N': Do not kill small columns of c*A. This option assumes that
0117: *              BLAS and QR factorizations and triangular solvers are
0118: *              implemented to work in that range. If the condition of A
0119: *              is greater than BIG, use SGESVJ.
0120: *       = 'R': RESTRICTED range for sigma(c*A) is [SQRT(SFMIN), SQRT(BIG)]
0121: *              (roughly, as described above). This option is recommended.
0122: *                                             ===========================
0123: *         For computing the singular values in the FULL range [SFMIN,BIG]
0124: *         use SGESVJ.
0125: * 
0126: *  JOBT   (input) CHARACTER*1
0127: *         If the matrix is square then the procedure may determine to use
0128: *         transposed A if A^t seems to be better with respect to convergence.
0129: *         If the matrix is not square, JOBT is ignored. This is subject to
0130: *         changes in the future.
0131: *         The decision is based on two values of entropy over the adjoint
0132: *         orbit of A^t * A. See the descriptions of WORK(6) and WORK(7).
0133: *       = 'T': transpose if entropy test indicates possibly faster
0134: *         convergence of Jacobi process if A^t is taken as input. If A is
0135: *         replaced with A^t, then the row pivoting is included automatically.
0136: *       = 'N': do not speculate.
0137: *         This option can be used to compute only the singular values, or the
0138: *         full SVD (U, SIGMA and V). For only one set of singular vectors
0139: *         (U or V), the caller should provide both U and V, as one of the
0140: *         matrices is used as workspace if the matrix A is transposed.
0141: *         The implementer can easily remove this constraint and make the
0142: *         code more complicated. See the descriptions of U and V.
0143: * 
0144: *  JOBP   (input) CHARACTER*1
0145: *         Issues the licence to introduce structured perturbations to drown
0146: *         denormalized numbers. This licence should be active if the
0147: *         denormals are poorly implemented, causing slow computation,
0148: *         especially in cases of fast convergence (!). For details see [1,2].
0149: *         For the sake of simplicity, this perturbations are included only
0150: *         when the full SVD or only the singular values are requested. The
0151: *         implementer/user can easily add the perturbation for the cases of
0152: *         computing one set of singular vectors.
0153: *       = 'P': introduce perturbation
0154: *       = 'N': do not perturb
0155: *
0156: *  M      (input) INTEGER
0157: *         The number of rows of the input matrix A.  M >= 0.
0158: *
0159: *  N      (input) INTEGER
0160: *         The number of columns of the input matrix A. M >= N >= 0.
0161: *
0162: *  A       (input/workspace) REAL array, dimension (LDA,N)
0163: *          On entry, the M-by-N matrix A.
0164: *
0165: *  LDA     (input) INTEGER
0166: *          The leading dimension of the array A.  LDA >= max(1,M).
0167: *
0168: *  SVA     (workspace/output) REAL array, dimension (N)
0169: *          On exit,
0170: *          - For WORK(1)/WORK(2) = ONE: The singular values of A. During the
0171: *            computation SVA contains Euclidean column norms of the
0172: *            iterated matrices in the array A.
0173: *          - For WORK(1) .NE. WORK(2): The singular values of A are
0174: *            (WORK(1)/WORK(2)) * SVA(1:N). This factored form is used if
0175: *            sigma_max(A) overflows or if small singular values have been
0176: *            saved from underflow by scaling the input matrix A.
0177: *          - If JOBR='R' then some of the singular values may be returned
0178: *            as exact zeros obtained by "set to zero" because they are
0179: *            below the numerical rank threshold or are denormalized numbers.
0180: *
0181: *  U       (workspace/output) REAL array, dimension ( LDU, N )
0182: *          If JOBU = 'U', then U contains on exit the M-by-N matrix of
0183: *                         the left singular vectors.
0184: *          If JOBU = 'F', then U contains on exit the M-by-M matrix of
0185: *                         the left singular vectors, including an ONB
0186: *                         of the orthogonal complement of the Range(A).
0187: *          If JOBU = 'W'  .AND. (JOBV.EQ.'V' .AND. JOBT.EQ.'T' .AND. M.EQ.N),
0188: *                         then U is used as workspace if the procedure
0189: *                         replaces A with A^t. In that case, [V] is computed
0190: *                         in U as left singular vectors of A^t and then
0191: *                         copied back to the V array. This 'W' option is just
0192: *                         a reminder to the caller that in this case U is
0193: *                         reserved as workspace of length N*N.
0194: *          If JOBU = 'N'  U is not referenced.
0195: *
0196: * LDU      (input) INTEGER
0197: *          The leading dimension of the array U,  LDU >= 1.
0198: *          IF  JOBU = 'U' or 'F' or 'W',  then LDU >= M.
0199: *
0200: *  V       (workspace/output) REAL array, dimension ( LDV, N )
0201: *          If JOBV = 'V', 'J' then V contains on exit the N-by-N matrix of
0202: *                         the right singular vectors;
0203: *          If JOBV = 'W', AND (JOBU.EQ.'U' AND JOBT.EQ.'T' AND M.EQ.N),
0204: *                         then V is used as workspace if the pprocedure
0205: *                         replaces A with A^t. In that case, [U] is computed
0206: *                         in V as right singular vectors of A^t and then
0207: *                         copied back to the U array. This 'W' option is just
0208: *                         a reminder to the caller that in this case V is
0209: *                         reserved as workspace of length N*N.
0210: *          If JOBV = 'N'  V is not referenced.
0211: *
0212: *  LDV     (input) INTEGER
0213: *          The leading dimension of the array V,  LDV >= 1.
0214: *          If JOBV = 'V' or 'J' or 'W', then LDV >= N.
0215: *
0216: *  WORK    (workspace/output) REAL array, dimension at least LWORK.
0217: *          On exit,
0218: *          WORK(1) = SCALE = WORK(2) / WORK(1) is the scaling factor such
0219: *                    that SCALE*SVA(1:N) are the computed singular values
0220: *                    of A. (See the description of SVA().)
0221: *          WORK(2) = See the description of WORK(1).
0222: *          WORK(3) = SCONDA is an estimate for the condition number of
0223: *                    column equilibrated A. (If JOBA .EQ. 'E' or 'G')
0224: *                    SCONDA is an estimate of SQRT(||(R^t * R)^(-1)||_1).
0225: *                    It is computed using SPOCON. It holds
0226: *                    N^(-1/4) * SCONDA <= ||R^(-1)||_2 <= N^(1/4) * SCONDA
0227: *                    where R is the triangular factor from the QRF of A.
0228: *                    However, if R is truncated and the numerical rank is
0229: *                    determined to be strictly smaller than N, SCONDA is
0230: *                    returned as -1, thus indicating that the smallest
0231: *                    singular values might be lost.
0232: *
0233: *          If full SVD is needed, the following two condition numbers are
0234: *          useful for the analysis of the algorithm. They are provied for
0235: *          a developer/implementer who is familiar with the details of
0236: *          the method.
0237: *
0238: *          WORK(4) = an estimate of the scaled condition number of the
0239: *                    triangular factor in the first QR factorization.
0240: *          WORK(5) = an estimate of the scaled condition number of the
0241: *                    triangular factor in the second QR factorization.
0242: *          The following two parameters are computed if JOBT .EQ. 'T'.
0243: *          They are provided for a developer/implementer who is familiar
0244: *          with the details of the method.
0245: *
0246: *          WORK(6) = the entropy of A^t*A :: this is the Shannon entropy
0247: *                    of diag(A^t*A) / Trace(A^t*A) taken as point in the
0248: *                    probability simplex.
0249: *          WORK(7) = the entropy of A*A^t.
0250: *
0251: *  LWORK   (input) INTEGER
0252: *          Length of WORK to confirm proper allocation of work space.
0253: *          LWORK depends on the job:
0254: *
0255: *          If only SIGMA is needed ( JOBU.EQ.'N', JOBV.EQ.'N' ) and
0256: *            -> .. no scaled condition estimate required ( JOBE.EQ.'N'):
0257: *               LWORK >= max(2*M+N,4*N+1,7). This is the minimal requirement.
0258: *               For optimal performance (blocked code) the optimal value
0259: *               is LWORK >= max(2*M+N,3*N+(N+1)*NB,7). Here NB is the optimal
0260: *               block size for xGEQP3/xGEQRF.
0261: *            -> .. an estimate of the scaled condition number of A is
0262: *               required (JOBA='E', 'G'). In this case, LWORK is the maximum
0263: *               of the above and N*N+4*N, i.e. LWORK >= max(2*M+N,N*N+4N,7).
0264: *
0265: *          If SIGMA and the right singular vectors are needed (JOBV.EQ.'V'),
0266: *            -> the minimal requirement is LWORK >= max(2*N+M,7).
0267: *            -> For optimal performance, LWORK >= max(2*N+M,2*N+N*NB,7),
0268: *               where NB is the optimal block size.
0269: *
0270: *          If SIGMA and the left singular vectors are needed
0271: *            -> the minimal requirement is LWORK >= max(2*N+M,7).
0272: *            -> For optimal performance, LWORK >= max(2*N+M,2*N+N*NB,7),
0273: *               where NB is the optimal block size.
0274: *
0275: *          If full SVD is needed ( JOBU.EQ.'U' or 'F', JOBV.EQ.'V' ) and
0276: *            -> .. the singular vectors are computed without explicit
0277: *               accumulation of the Jacobi rotations, LWORK >= 6*N+2*N*N
0278: *            -> .. in the iterative part, the Jacobi rotations are
0279: *               explicitly accumulated (option, see the description of JOBV),
0280: *               then the minimal requirement is LWORK >= max(M+3*N+N*N,7).
0281: *               For better performance, if NB is the optimal block size,
0282: *               LWORK >= max(3*N+N*N+M,3*N+N*N+N*NB,7).
0283: *
0284: *  IWORK   (workspace/output) INTEGER array, dimension M+3*N.
0285: *          On exit,
0286: *          IWORK(1) = the numerical rank determined after the initial
0287: *                     QR factorization with pivoting. See the descriptions
0288: *                     of JOBA and JOBR.
0289: *          IWORK(2) = the number of the computed nonzero singular values
0290: *          IWORK(3) = if nonzero, a warning message:
0291: *                     If IWORK(3).EQ.1 then some of the column norms of A
0292: *                     were denormalized floats. The requested high accuracy
0293: *                     is not warranted by the data.
0294: *
0295: *  INFO    (output) INTEGER
0296: *           < 0  : if INFO = -i, then the i-th argument had an illegal value.
0297: *           = 0 :  successfull exit;
0298: *           > 0 :  SGEJSV  did not converge in the maximal allowed number
0299: *                  of sweeps. The computed values may be inaccurate.
0300: *
0301: *  Further Details
0302: *  ===============
0303: *
0304: *  SGEJSV implements a preconditioned Jacobi SVD algorithm. It uses SGEQP3,
0305: *  SGEQRF, and SGELQF as preprocessors and preconditioners. Optionally, an
0306: *  additional row pivoting can be used as a preprocessor, which in some
0307: *  cases results in much higher accuracy. An example is matrix A with the
0308: *  structure A = D1 * C * D2, where D1, D2 are arbitrarily ill-conditioned
0309: *  diagonal matrices and C is well-conditioned matrix. In that case, complete
0310: *  pivoting in the first QR factorizations provides accuracy dependent on the
0311: *  condition number of C, and independent of D1, D2. Such higher accuracy is
0312: *  not completely understood theoretically, but it works well in practice.
0313: *  Further, if A can be written as A = B*D, with well-conditioned B and some
0314: *  diagonal D, then the high accuracy is guaranteed, both theoretically and
0315: *  in software, independent of D. For more details see [1], [2].
0316: *     The computational range for the singular values can be the full range
0317: *  ( UNDERFLOW,OVERFLOW ), provided that the machine arithmetic and the BLAS
0318: *  & LAPACK routines called by SGEJSV are implemented to work in that range.
0319: *  If that is not the case, then the restriction for safe computation with
0320: *  the singular values in the range of normalized IEEE numbers is that the
0321: *  spectral condition number kappa(A)=sigma_max(A)/sigma_min(A) does not
0322: *  overflow. This code (SGEJSV) is best used in this restricted range,
0323: *  meaning that singular values of magnitude below ||A||_2 / SLAMCH('O') are
0324: *  returned as zeros. See JOBR for details on this.
0325: *     Further, this implementation is somewhat slower than the one described
0326: *  in [1,2] due to replacement of some non-LAPACK components, and because
0327: *  the choice of some tuning parameters in the iterative part (SGESVJ) is
0328: *  left to the implementer on a particular machine.
0329: *     The rank revealing QR factorization (in this code: SGEQP3) should be
0330: *  implemented as in [3]. We have a new version of SGEQP3 under development
0331: *  that is more robust than the current one in LAPACK, with a cleaner cut in
0332: *  rank defficient cases. It will be available in the SIGMA library [4].
0333: *  If M is much larger than N, it is obvious that the inital QRF with
0334: *  column pivoting can be preprocessed by the QRF without pivoting. That
0335: *  well known trick is not used in SGEJSV because in some cases heavy row
0336: *  weighting can be treated with complete pivoting. The overhead in cases
0337: *  M much larger than N is then only due to pivoting, but the benefits in
0338: *  terms of accuracy have prevailed. The implementer/user can incorporate
0339: *  this extra QRF step easily. The implementer can also improve data movement
0340: *  (matrix transpose, matrix copy, matrix transposed copy) - this
0341: *  implementation of SGEJSV uses only the simplest, naive data movement.
0342: *
0343: *  Contributors
0344: *
0345: *  Zlatko Drmac (Zagreb, Croatia) and Kresimir Veselic (Hagen, Germany)
0346: *
0347: *  References
0348: *
0349: * [1] Z. Drmac and K. Veselic: New fast and accurate Jacobi SVD algorithm I.
0350: *     SIAM J. Matrix Anal. Appl. Vol. 35, No. 2 (2008), pp. 1322-1342.
0351: *     LAPACK Working note 169.
0352: * [2] Z. Drmac and K. Veselic: New fast and accurate Jacobi SVD algorithm II.
0353: *     SIAM J. Matrix Anal. Appl. Vol. 35, No. 2 (2008), pp. 1343-1362.
0354: *     LAPACK Working note 170.
0355: * [3] Z. Drmac and Z. Bujanovic: On the failure of rank-revealing QR
0356: *     factorization software - a case study.
0357: *     ACM Trans. math. Softw. Vol. 35, No 2 (2008), pp. 1-28.
0358: *     LAPACK Working note 176.
0359: * [4] Z. Drmac: SIGMA - mathematical software library for accurate SVD, PSV,
0360: *     QSVD, (H,K)-SVD computations.
0361: *     Department of Mathematics, University of Zagreb, 2008.
0362: *
0363: *  Bugs, examples and comments
0364: *
0365: *  Please report all bugs and send interesting examples and/or comments to
0366: *  drmac@math.hr. Thank you.
0367: *
0368: *  ===========================================================================
0369: *
0370: *     .. Local Parameters ..
0371:       REAL        ZERO,         ONE
0372:       PARAMETER ( ZERO = 0.0E0, ONE = 1.0E0 )
0373: *     ..
0374: *     .. Local Scalars ..
0375:       REAL    AAPP,   AAQQ,   AATMAX, AATMIN, BIG,    BIG1,   COND_OK,
0376:      &        CONDR1, CONDR2, ENTRA,  ENTRAT, EPSLN,  MAXPRJ, SCALEM,
0377:      &        SCONDA, SFMIN,  SMALL,  TEMP1,  USCAL1, USCAL2, XSC
0378:       INTEGER IERR,   N1,     NR,     NUMRANK,        p, q,   WARNING
0379:       LOGICAL ALMORT, DEFR,   ERREST, GOSCAL, JRACC,  KILL,   LSVEC,
0380:      &        L2ABER, L2KILL, L2PERT, L2RANK, L2TRAN,
0381:      &        NOSCAL, ROWPIV, RSVEC,  TRANSP
0382: *     ..
0383: *     .. Intrinsic Functions ..
0384:       INTRINSIC ABS,  ALOG, AMAX1, AMIN1, FLOAT,
0385:      &          MAX0, MIN0, NINT,  SIGN,  SQRT
0386: *     ..
0387: *     .. External Functions ..
0388:       REAL      SLAMCH, SNRM2
0389:       INTEGER   ISAMAX
0390:       LOGICAL   LSAME
0391:       EXTERNAL  ISAMAX, LSAME, SLAMCH, SNRM2
0392: *     ..
0393: *     .. External Subroutines ..
0394:       EXTERNAL  SCOPY,  SGELQF, SGEQP3, SGEQRF, SLACPY, SLASCL,
0395:      &          SLASET, SLASSQ, SLASWP, SORGQR, SORMLQ,
0396:      &          SORMQR, SPOCON, SSCAL,  SSWAP,  STRSM,  XERBLA
0397: *
0398:       EXTERNAL  SGESVJ
0399: *     ..
0400: *
0401: *     Test the input arguments
0402: *
0403:       LSVEC  = LSAME( JOBU, 'U' ) .OR. LSAME( JOBU, 'F' )
0404:       JRACC  = LSAME( JOBV, 'J' )
0405:       RSVEC  = LSAME( JOBV, 'V' ) .OR. JRACC
0406:       ROWPIV = LSAME( JOBA, 'F' ) .OR. LSAME( JOBA, 'G' )
0407:       L2RANK = LSAME( JOBA, 'R' )
0408:       L2ABER = LSAME( JOBA, 'A' )
0409:       ERREST = LSAME( JOBA, 'E' ) .OR. LSAME( JOBA, 'G' )
0410:       L2TRAN = LSAME( JOBT, 'T' )
0411:       L2KILL = LSAME( JOBR, 'R' )
0412:       DEFR   = LSAME( JOBR, 'N' )
0413:       L2PERT = LSAME( JOBP, 'P' )
0414: *
0415:       IF ( .NOT.(ROWPIV .OR. L2RANK .OR. L2ABER .OR.
0416:      &     ERREST .OR. LSAME( JOBA, 'C' ) )) THEN
0417:          INFO = - 1
0418:       ELSE IF ( .NOT.( LSVEC  .OR. LSAME( JOBU, 'N' ) .OR.
0419:      &                             LSAME( JOBU, 'W' )) ) THEN
0420:          INFO = - 2
0421:       ELSE IF ( .NOT.( RSVEC .OR. LSAME( JOBV, 'N' ) .OR.
0422:      &   LSAME( JOBV, 'W' )) .OR. ( JRACC .AND. (.NOT.LSVEC) ) ) THEN
0423:          INFO = - 3
0424:       ELSE IF ( .NOT. ( L2KILL .OR. DEFR ) )    THEN
0425:          INFO = - 4
0426:       ELSE IF ( .NOT. ( L2TRAN .OR. LSAME( JOBT, 'N' ) ) ) THEN
0427:          INFO = - 5
0428:       ELSE IF ( .NOT. ( L2PERT .OR. LSAME( JOBP, 'N' ) ) ) THEN
0429:          INFO = - 6
0430:       ELSE IF ( M .LT. 0 ) THEN
0431:          INFO = - 7
0432:       ELSE IF ( ( N .LT. 0 ) .OR. ( N .GT. M ) ) THEN
0433:          INFO = - 8
0434:       ELSE IF ( LDA .LT. M ) THEN
0435:          INFO = - 10
0436:       ELSE IF ( LSVEC .AND. ( LDU .LT. M ) ) THEN
0437:          INFO = - 13
0438:       ELSE IF ( RSVEC .AND. ( LDV .LT. N ) ) THEN
0439:          INFO = - 14
0440:       ELSE IF ( (.NOT.(LSVEC .OR. RSVEC .OR. ERREST).AND.
0441:      &                           (LWORK .LT. MAX0(7,4*N+1,2*M+N))) .OR.
0442:      & (.NOT.(LSVEC .OR. LSVEC) .AND. ERREST .AND.
0443:      &                         (LWORK .LT. MAX0(7,4*N+N*N,2*M+N))) .OR.
0444:      & (LSVEC .AND. (.NOT.RSVEC) .AND. (LWORK .LT. MAX0(7,2*N+M))) .OR.
0445:      & (RSVEC .AND. (.NOT.LSVEC) .AND. (LWORK .LT. MAX0(7,2*N+M))) .OR.
0446:      & (LSVEC .AND. RSVEC .AND. .NOT.JRACC .AND. (LWORK.LT.6*N+2*N*N))
0447:      & .OR. (LSVEC.AND.RSVEC.AND.JRACC.AND.LWORK.LT.MAX0(7,M+3*N+N*N)))
0448:      &   THEN
0449:          INFO = - 17
0450:       ELSE
0451: *        #:)
0452:          INFO = 0
0453:       END IF
0454: *
0455:       IF ( INFO .NE. 0 ) THEN
0456: *       #:(
0457:          CALL XERBLA( 'SGEJSV', - INFO )
0458:       END IF
0459: *
0460: *     Quick return for void matrix (Y3K safe)
0461: * #:)
0462:       IF ( ( M .EQ. 0 ) .OR. ( N .EQ. 0 ) ) RETURN
0463: *
0464: *     Determine whether the matrix U should be M x N or M x M
0465: *
0466:       IF ( LSVEC ) THEN
0467:          N1 = N
0468:          IF ( LSAME( JOBU, 'F' ) ) N1 = M
0469:       END IF
0470: *
0471: *     Set numerical parameters
0472: *
0473: *!    NOTE: Make sure SLAMCH() does not fail on the target architecture.
0474: *
0475:       EPSLN = SLAMCH('Epsilon')
0476:       SFMIN = SLAMCH('SafeMinimum')
0477:       SMALL = SFMIN / EPSLN
0478:       BIG   = SLAMCH('O')
0479: *
0480: *     Initialize SVA(1:N) = diag( ||A e_i||_2 )_1^N
0481: *
0482: *(!)  If necessary, scale SVA() to protect the largest norm from
0483: *     overflow. It is possible that this scaling pushes the smallest
0484: *     column norm left from the underflow threshold (extreme case).
0485: *
0486:       SCALEM  = ONE / SQRT(FLOAT(M)*FLOAT(N))
0487:       NOSCAL  = .TRUE.
0488:       GOSCAL  = .TRUE.
0489:       DO 1874 p = 1, N
0490:          AAPP = ZERO
0491:          AAQQ = ZERO
0492:          CALL SLASSQ( M, A(1,p), 1, AAPP, AAQQ )
0493:          IF ( AAPP .GT. BIG ) THEN
0494:             INFO = - 9
0495:             CALL XERBLA( 'SGEJSV', -INFO )
0496:             RETURN
0497:          END IF
0498:          AAQQ = SQRT(AAQQ)
0499:          IF ( ( AAPP .LT. (BIG / AAQQ) ) .AND. NOSCAL  ) THEN
0500:             SVA(p)  = AAPP * AAQQ
0501:          ELSE
0502:             NOSCAL  = .FALSE.
0503:             SVA(p)  = AAPP * ( AAQQ * SCALEM )
0504:             IF ( GOSCAL ) THEN
0505:                GOSCAL = .FALSE.
0506:                CALL SSCAL( p-1, SCALEM, SVA, 1 )
0507:             END IF
0508:          END IF
0509:  1874 CONTINUE
0510: *
0511:       IF ( NOSCAL ) SCALEM = ONE
0512: *
0513:       AAPP = ZERO
0514:       AAQQ = BIG
0515:       DO 4781 p = 1, N
0516:          AAPP = AMAX1( AAPP, SVA(p) )
0517:          IF ( SVA(p) .NE. ZERO ) AAQQ = AMIN1( AAQQ, SVA(p) )
0518:  4781 CONTINUE
0519: *
0520: *     Quick return for zero M x N matrix
0521: * #:)
0522:       IF ( AAPP .EQ. ZERO ) THEN
0523:          IF ( LSVEC ) CALL SLASET( 'G', M, N1, ZERO, ONE, U, LDU )
0524:          IF ( RSVEC ) CALL SLASET( 'G', N, N,  ZERO, ONE, V, LDV )
0525:          WORK(1) = ONE
0526:          WORK(2) = ONE
0527:          IF ( ERREST ) WORK(3) = ONE
0528:          IF ( LSVEC .AND. RSVEC ) THEN
0529:             WORK(4) = ONE
0530:             WORK(5) = ONE
0531:          END IF
0532:          IF ( L2TRAN ) THEN
0533:             WORK(6) = ZERO
0534:             WORK(7) = ZERO
0535:          END IF
0536:          IWORK(1) = 0
0537:          IWORK(2) = 0
0538:          RETURN
0539:       END IF
0540: *
0541: *     Issue warning if denormalized column norms detected. Override the
0542: *     high relative accuracy request. Issue licence to kill columns
0543: *     (set them to zero) whose norm is less than sigma_max / BIG (roughly).
0544: * #:(
0545:       WARNING = 0
0546:       IF ( AAQQ .LE. SFMIN ) THEN
0547:          L2RANK = .TRUE.
0548:          L2KILL = .TRUE.
0549:          WARNING = 1
0550:       END IF
0551: *
0552: *     Quick return for one-column matrix
0553: * #:)
0554:       IF ( N .EQ. 1 ) THEN
0555: *
0556:          IF ( LSVEC ) THEN
0557:             CALL SLASCL( 'G',0,0,SVA(1),SCALEM, M,1,A(1,1),LDA,IERR )
0558:             CALL SLACPY( 'A', M, 1, A, LDA, U, LDU )
0559: *           computing all M left singular vectors of the M x 1 matrix
0560:             IF ( N1 .NE. N  ) THEN
0561:                CALL SGEQRF( M, N, U,LDU, WORK, WORK(N+1),LWORK-N,IERR )
0562:                CALL SORGQR( M,N1,1, U,LDU,WORK,WORK(N+1),LWORK-N,IERR )
0563:                CALL SCOPY( M, A(1,1), 1, U(1,1), 1 )
0564:             END IF
0565:          END IF
0566:          IF ( RSVEC ) THEN
0567:              V(1,1) = ONE
0568:          END IF
0569:          IF ( SVA(1) .LT. (BIG*SCALEM) ) THEN
0570:             SVA(1)  = SVA(1) / SCALEM
0571:             SCALEM  = ONE
0572:          END IF
0573:          WORK(1) = ONE / SCALEM
0574:          WORK(2) = ONE
0575:          IF ( SVA(1) .NE. ZERO ) THEN
0576:             IWORK(1) = 1
0577:             IF ( ( SVA(1) / SCALEM) .GE. SFMIN ) THEN
0578:                IWORK(2) = 1
0579:             ELSE
0580:                IWORK(2) = 0
0581:             END IF
0582:          ELSE
0583:             IWORK(1) = 0
0584:             IWORK(2) = 0
0585:          END IF
0586:          IF ( ERREST ) WORK(3) = ONE
0587:          IF ( LSVEC .AND. RSVEC ) THEN
0588:             WORK(4) = ONE
0589:             WORK(5) = ONE
0590:          END IF
0591:          IF ( L2TRAN ) THEN
0592:             WORK(6) = ZERO
0593:             WORK(7) = ZERO
0594:          END IF
0595:          RETURN
0596: *
0597:       END IF
0598: *
0599:       TRANSP = .FALSE.
0600:       L2TRAN = L2TRAN .AND. ( M .EQ. N )
0601: *
0602:       AATMAX = -ONE
0603:       AATMIN =  BIG
0604:       IF ( ROWPIV .OR. L2TRAN ) THEN
0605: *
0606: *     Compute the row norms, needed to determine row pivoting sequence
0607: *     (in the case of heavily row weighted A, row pivoting is strongly
0608: *     advised) and to collect information needed to compare the
0609: *     structures of A * A^t and A^t * A (in the case L2TRAN.EQ..TRUE.).
0610: *
0611:          IF ( L2TRAN ) THEN
0612:             DO 1950 p = 1, M
0613:                XSC   = ZERO
0614:                TEMP1 = ZERO
0615:                CALL SLASSQ( N, A(p,1), LDA, XSC, TEMP1 )
0616: *              SLASSQ gets both the ell_2 and the ell_infinity norm
0617: *              in one pass through the vector
0618:                WORK(M+N+p)  = XSC * SCALEM
0619:                WORK(N+p)    = XSC * (SCALEM*SQRT(TEMP1))
0620:                AATMAX = AMAX1( AATMAX, WORK(N+p) )
0621:                IF (WORK(N+p) .NE. ZERO) AATMIN = AMIN1(AATMIN,WORK(N+p))
0622:  1950       CONTINUE
0623:          ELSE
0624:             DO 1904 p = 1, M
0625:                WORK(M+N+p) = SCALEM*ABS( A(p,ISAMAX(N,A(p,1),LDA)) )
0626:                AATMAX = AMAX1( AATMAX, WORK(M+N+p) )
0627:                AATMIN = AMIN1( AATMIN, WORK(M+N+p) )
0628:  1904       CONTINUE
0629:          END IF
0630: *
0631:       END IF
0632: *
0633: *     For square matrix A try to determine whether A^t  would be  better
0634: *     input for the preconditioned Jacobi SVD, with faster convergence.
0635: *     The decision is based on an O(N) function of the vector of column
0636: *     and row norms of A, based on the Shannon entropy. This should give
0637: *     the right choice in most cases when the difference actually matters.
0638: *     It may fail and pick the slower converging side.
0639: *
0640:       ENTRA  = ZERO
0641:       ENTRAT = ZERO
0642:       IF ( L2TRAN ) THEN
0643: *
0644:          XSC   = ZERO
0645:          TEMP1 = ZERO
0646:          CALL SLASSQ( N, SVA, 1, XSC, TEMP1 )
0647:          TEMP1 = ONE / TEMP1
0648: *
0649:          ENTRA = ZERO
0650:          DO 1113 p = 1, N
0651:             BIG1  = ( ( SVA(p) / XSC )**2 ) * TEMP1
0652:             IF ( BIG1 .NE. ZERO ) ENTRA = ENTRA + BIG1 * ALOG(BIG1)
0653:  1113    CONTINUE
0654:          ENTRA = - ENTRA / ALOG(FLOAT(N))
0655: *
0656: *        Now, SVA().^2/Trace(A^t * A) is a point in the probability simplex.
0657: *        It is derived from the diagonal of  A^t * A.  Do the same with the
0658: *        diagonal of A * A^t, compute the entropy of the corresponding
0659: *        probability distribution. Note that A * A^t and A^t * A have the
0660: *        same trace.
0661: *
0662:          ENTRAT = ZERO
0663:          DO 1114 p = N+1, N+M
0664:             BIG1 = ( ( WORK(p) / XSC )**2 ) * TEMP1
0665:             IF ( BIG1 .NE. ZERO ) ENTRAT = ENTRAT + BIG1 * ALOG(BIG1)
0666:  1114    CONTINUE
0667:          ENTRAT = - ENTRAT / ALOG(FLOAT(M))
0668: *
0669: *        Analyze the entropies and decide A or A^t. Smaller entropy
0670: *        usually means better input for the algorithm.
0671: *
0672:          TRANSP = ( ENTRAT .LT. ENTRA )
0673: *
0674: *        If A^t is better than A, transpose A.
0675: *
0676:          IF ( TRANSP ) THEN
0677: *           In an optimal implementation, this trivial transpose
0678: *           should be replaced with faster transpose.
0679:             DO 1115 p = 1, N - 1
0680:                DO 1116 q = p + 1, N
0681:                    TEMP1 = A(q,p)
0682:                   A(q,p) = A(p,q)
0683:                   A(p,q) = TEMP1
0684:  1116          CONTINUE
0685:  1115       CONTINUE
0686:             DO 1117 p = 1, N
0687:                WORK(M+N+p) = SVA(p)
0688:                SVA(p)      = WORK(N+p)
0689:  1117       CONTINUE
0690:             TEMP1  = AAPP
0691:             AAPP   = AATMAX
0692:             AATMAX = TEMP1
0693:             TEMP1  = AAQQ
0694:             AAQQ   = AATMIN
0695:             AATMIN = TEMP1
0696:             KILL   = LSVEC
0697:             LSVEC  = RSVEC
0698:             RSVEC  = KILL
0699: *
0700:             ROWPIV = .TRUE.
0701:          END IF
0702: *
0703:       END IF
0704: *     END IF L2TRAN
0705: *
0706: *     Scale the matrix so that its maximal singular value remains less
0707: *     than SQRT(BIG) -- the matrix is scaled so that its maximal column
0708: *     has Euclidean norm equal to SQRT(BIG/N). The only reason to keep
0709: *     SQRT(BIG) instead of BIG is the fact that SGEJSV uses LAPACK and
0710: *     BLAS routines that, in some implementations, are not capable of
0711: *     working in the full interval [SFMIN,BIG] and that they may provoke
0712: *     overflows in the intermediate results. If the singular values spread
0713: *     from SFMIN to BIG, then SGESVJ will compute them. So, in that case,
0714: *     one should use SGESVJ instead of SGEJSV.
0715: *
0716:       BIG1   = SQRT( BIG )
0717:       TEMP1  = SQRT( BIG / FLOAT(N) )
0718: *
0719:       CALL SLASCL( 'G', 0, 0, AAPP, TEMP1, N, 1, SVA, N, IERR )
0720:       IF ( AAQQ .GT. (AAPP * SFMIN) ) THEN
0721:           AAQQ = ( AAQQ / AAPP ) * TEMP1
0722:       ELSE
0723:           AAQQ = ( AAQQ * TEMP1 ) / AAPP
0724:       END IF
0725:       TEMP1 = TEMP1 * SCALEM
0726:       CALL SLASCL( 'G', 0, 0, AAPP, TEMP1, M, N, A, LDA, IERR )
0727: *
0728: *     To undo scaling at the end of this procedure, multiply the
0729: *     computed singular values with USCAL2 / USCAL1.
0730: *
0731:       USCAL1 = TEMP1
0732:       USCAL2 = AAPP
0733: *
0734:       IF ( L2KILL ) THEN
0735: *        L2KILL enforces computation of nonzero singular values in
0736: *        the restricted range of condition number of the initial A,
0737: *        sigma_max(A) / sigma_min(A) approx. SQRT(BIG)/SQRT(SFMIN).
0738:          XSC = SQRT( SFMIN )
0739:       ELSE
0740:          XSC = SMALL
0741: *
0742: *        Now, if the condition number of A is too big,
0743: *        sigma_max(A) / sigma_min(A) .GT. SQRT(BIG/N) * EPSLN / SFMIN,
0744: *        as a precaution measure, the full SVD is computed using SGESVJ
0745: *        with accumulated Jacobi rotations. This provides numerically
0746: *        more robust computation, at the cost of slightly increased run
0747: *        time. Depending on the concrete implementation of BLAS and LAPACK
0748: *        (i.e. how they behave in presence of extreme ill-conditioning) the
0749: *        implementor may decide to remove this switch.
0750:          IF ( ( AAQQ.LT.SQRT(SFMIN) ) .AND. LSVEC .AND. RSVEC ) THEN
0751:             JRACC = .TRUE.
0752:          END IF
0753: *
0754:       END IF
0755:       IF ( AAQQ .LT. XSC ) THEN
0756:          DO 700 p = 1, N
0757:             IF ( SVA(p) .LT. XSC ) THEN
0758:                CALL SLASET( 'A', M, 1, ZERO, ZERO, A(1,p), LDA )
0759:                SVA(p) = ZERO
0760:             END IF
0761:  700     CONTINUE
0762:       END IF
0763: *
0764: *     Preconditioning using QR factorization with pivoting
0765: *
0766:       IF ( ROWPIV ) THEN
0767: *        Optional row permutation (Bjoerck row pivoting):
0768: *        A result by Cox and Higham shows that the Bjoerck's
0769: *        row pivoting combined with standard column pivoting
0770: *        has similar effect as Powell-Reid complete pivoting.
0771: *        The ell-infinity norms of A are made nonincreasing.
0772:          DO 1952 p = 1, M - 1
0773:             q = ISAMAX( M-p+1, WORK(M+N+p), 1 ) + p - 1
0774:             IWORK(2*N+p) = q
0775:             IF ( p .NE. q ) THEN
0776:                TEMP1       = WORK(M+N+p)
0777:                WORK(M+N+p) = WORK(M+N+q)
0778:                WORK(M+N+q) = TEMP1
0779:             END IF
0780:  1952    CONTINUE
0781:          CALL SLASWP( N, A, LDA, 1, M-1, IWORK(2*N+1), 1 )
0782:       END IF
0783: *
0784: *     End of the preparation phase (scaling, optional sorting and
0785: *     transposing, optional flushing of small columns).
0786: *
0787: *     Preconditioning
0788: *
0789: *     If the full SVD is needed, the right singular vectors are computed
0790: *     from a matrix equation, and for that we need theoretical analysis
0791: *     of the Businger-Golub pivoting. So we use SGEQP3 as the first RR QRF.
0792: *     In all other cases the first RR QRF can be chosen by other criteria
0793: *     (eg speed by replacing global with restricted window pivoting, such
0794: *     as in SGEQPX from TOMS # 782). Good results will be obtained using
0795: *     SGEQPX with properly (!) chosen numerical parameters.
0796: *     Any improvement of SGEQP3 improves overal performance of SGEJSV.
0797: *
0798: *     A * P1 = Q1 * [ R1^t 0]^t:
0799:       DO 1963 p = 1, N
0800: *        .. all columns are free columns
0801:          IWORK(p) = 0
0802:  1963 CONTINUE
0803:       CALL SGEQP3( M,N,A,LDA, IWORK,WORK, WORK(N+1),LWORK-N, IERR )
0804: *
0805: *     The upper triangular matrix R1 from the first QRF is inspected for
0806: *     rank deficiency and possibilities for deflation, or possible
0807: *     ill-conditioning. Depending on the user specified flag L2RANK,
0808: *     the procedure explores possibilities to reduce the numerical
0809: *     rank by inspecting the computed upper triangular factor. If
0810: *     L2RANK or L2ABER are up, then SGEJSV will compute the SVD of
0811: *     A + dA, where ||dA|| <= f(M,N)*EPSLN.
0812: *
0813:       NR = 1
0814:       IF ( L2ABER ) THEN
0815: *        Standard absolute error bound suffices. All sigma_i with
0816: *        sigma_i < N*EPSLN*||A|| are flushed to zero. This is an
0817: *        agressive enforcement of lower numerical rank by introducing a
0818: *        backward error of the order of N*EPSLN*||A||.
0819:          TEMP1 = SQRT(FLOAT(N))*EPSLN
0820:          DO 3001 p = 2, N
0821:             IF ( ABS(A(p,p)) .GE. (TEMP1*ABS(A(1,1))) ) THEN
0822:                NR = NR + 1
0823:             ELSE
0824:                GO TO 3002
0825:             END IF
0826:  3001    CONTINUE
0827:  3002    CONTINUE
0828:       ELSE IF ( L2RANK ) THEN
0829: *        .. similarly as above, only slightly more gentle (less agressive).
0830: *        Sudden drop on the diagonal of R1 is used as the criterion for
0831: *        close-to-rank-defficient.
0832:          TEMP1 = SQRT(SFMIN)
0833:          DO 3401 p = 2, N
0834:             IF ( ( ABS(A(p,p)) .LT. (EPSLN*ABS(A(p-1,p-1))) ) .OR.
0835:      &           ( ABS(A(p,p)) .LT. SMALL ) .OR.
0836:      &           ( L2KILL .AND. (ABS(A(p,p)) .LT. TEMP1) ) ) GO TO 3402
0837:             NR = NR + 1
0838:  3401    CONTINUE
0839:  3402    CONTINUE
0840: *
0841:       ELSE
0842: *        The goal is high relative accuracy. However, if the matrix
0843: *        has high scaled condition number the relative accuracy is in
0844: *        general not feasible. Later on, a condition number estimator
0845: *        will be deployed to estimate the scaled condition number.
0846: *        Here we just remove the underflowed part of the triangular
0847: *        factor. This prevents the situation in which the code is
0848: *        working hard to get the accuracy not warranted by the data.
0849:          TEMP1  = SQRT(SFMIN)
0850:          DO 3301 p = 2, N
0851:             IF ( ( ABS(A(p,p)) .LT. SMALL ) .OR.
0852:      &           ( L2KILL .AND. (ABS(A(p,p)) .LT. TEMP1) ) ) GO TO 3302
0853:             NR = NR + 1
0854:  3301    CONTINUE
0855:  3302    CONTINUE
0856: *
0857:       END IF
0858: *
0859:       ALMORT = .FALSE.
0860:       IF ( NR .EQ. N ) THEN
0861:          MAXPRJ = ONE
0862:          DO 3051 p = 2, N
0863:             TEMP1  = ABS(A(p,p)) / SVA(IWORK(p))
0864:             MAXPRJ = AMIN1( MAXPRJ, TEMP1 )
0865:  3051    CONTINUE
0866:          IF ( MAXPRJ**2 .GE. ONE - FLOAT(N)*EPSLN ) ALMORT = .TRUE.
0867:       END IF
0868: *
0869: *
0870:       SCONDA = - ONE
0871:       CONDR1 = - ONE
0872:       CONDR2 = - ONE
0873: *
0874:       IF ( ERREST ) THEN
0875:          IF ( N .EQ. NR ) THEN
0876:             IF ( RSVEC ) THEN
0877: *              .. V is available as workspace
0878:                CALL SLACPY( 'U', N, N, A, LDA, V, LDV )
0879:                DO 3053 p = 1, N
0880:                   TEMP1 = SVA(IWORK(p))
0881:                   CALL SSCAL( p, ONE/TEMP1, V(1,p), 1 )
0882:  3053          CONTINUE
0883:                CALL SPOCON( 'U', N, V, LDV, ONE, TEMP1,
0884:      &              WORK(N+1), IWORK(2*N+M+1), IERR )
0885:             ELSE IF ( LSVEC ) THEN
0886: *              .. U is available as workspace
0887:                CALL SLACPY( 'U', N, N, A, LDA, U, LDU )
0888:                DO 3054 p = 1, N
0889:                   TEMP1 = SVA(IWORK(p))
0890:                   CALL SSCAL( p, ONE/TEMP1, U(1,p), 1 )
0891:  3054          CONTINUE
0892:                CALL SPOCON( 'U', N, U, LDU, ONE, TEMP1,
0893:      &              WORK(N+1), IWORK(2*N+M+1), IERR )
0894:             ELSE
0895:                CALL SLACPY( 'U', N, N, A, LDA, WORK(N+1), N )
0896:                DO 3052 p = 1, N
0897:                   TEMP1 = SVA(IWORK(p))
0898:                   CALL SSCAL( p, ONE/TEMP1, WORK(N+(p-1)*N+1), 1 )
0899:  3052          CONTINUE
0900: *           .. the columns of R are scaled to have unit Euclidean lengths.
0901:                CALL SPOCON( 'U', N, WORK(N+1), N, ONE, TEMP1,
0902:      &              WORK(N+N*N+1), IWORK(2*N+M+1), IERR )
0903:             END IF
0904:             SCONDA = ONE / SQRT(TEMP1)
0905: *           SCONDA is an estimate of SQRT(||(R^t * R)^(-1)||_1).
0906: *           N^(-1/4) * SCONDA <= ||R^(-1)||_2 <= N^(1/4) * SCONDA
0907:          ELSE
0908:             SCONDA = - ONE
0909:          END IF
0910:       END IF
0911: *
0912:       L2PERT = L2PERT .AND. ( ABS( A(1,1)/A(NR,NR) ) .GT. SQRT(BIG1) )
0913: *     If there is no violent scaling, artificial perturbation is not needed.
0914: *
0915: *     Phase 3:
0916: *
0917:       IF ( .NOT. ( RSVEC .OR. LSVEC ) ) THEN
0918: *
0919: *         Singular Values only
0920: *
0921: *         .. transpose A(1:NR,1:N)
0922:          DO 1946 p = 1, MIN0( N-1, NR )
0923:             CALL SCOPY( N-p, A(p,p+1), LDA, A(p+1,p), 1 )
0924:  1946    CONTINUE
0925: *
0926: *        The following two DO-loops introduce small relative perturbation
0927: *        into the strict upper triangle of the lower triangular matrix.
0928: *        Small entries below the main diagonal are also changed.
0929: *        This modification is useful if the computing environment does not
0930: *        provide/allow FLUSH TO ZERO underflow, for it prevents many
0931: *        annoying denormalized numbers in case of strongly scaled matrices.
0932: *        The perturbation is structured so that it does not introduce any
0933: *        new perturbation of the singular values, and it does not destroy
0934: *        the job done by the preconditioner.
0935: *        The licence for this perturbation is in the variable L2PERT, which
0936: *        should be .FALSE. if FLUSH TO ZERO underflow is active.
0937: *
0938:          IF ( .NOT. ALMORT ) THEN
0939: *
0940:             IF ( L2PERT ) THEN
0941: *              XSC = SQRT(SMALL)
0942:                XSC = EPSLN / FLOAT(N)
0943:                DO 4947 q = 1, NR
0944:                   TEMP1 = XSC*ABS(A(q,q))
0945:                   DO 4949 p = 1, N
0946:                      IF ( ( (p.GT.q) .AND. (ABS(A(p,q)).LE.TEMP1) )
0947:      &                    .OR. ( p .LT. q ) )
0948:      &                     A(p,q) = SIGN( TEMP1, A(p,q) )
0949:  4949             CONTINUE
0950:  4947          CONTINUE
0951:             ELSE
0952:                CALL SLASET( 'U', NR-1,NR-1, ZERO,ZERO, A(1,2),LDA )
0953:             END IF
0954: *
0955: *            .. second preconditioning using the QR factorization
0956: *
0957:             CALL SGEQRF( N,NR, A,LDA, WORK, WORK(N+1),LWORK-N, IERR )
0958: *
0959: *           .. and transpose upper to lower triangular
0960:             DO 1948 p = 1, NR - 1
0961:                CALL SCOPY( NR-p, A(p,p+1), LDA, A(p+1,p), 1 )
0962:  1948       CONTINUE
0963: *
0964:          END IF
0965: *
0966: *           Row-cyclic Jacobi SVD algorithm with column pivoting
0967: *
0968: *           .. again some perturbation (a "background noise") is added
0969: *           to drown denormals
0970:             IF ( L2PERT ) THEN
0971: *              XSC = SQRT(SMALL)
0972:                XSC = EPSLN / FLOAT(N)
0973:                DO 1947 q = 1, NR
0974:                   TEMP1 = XSC*ABS(A(q,q))
0975:                   DO 1949 p = 1, NR
0976:                      IF ( ( (p.GT.q) .AND. (ABS(A(p,q)).LE.TEMP1) )
0977:      &                       .OR. ( p .LT. q ) )
0978:      &                   A(p,q) = SIGN( TEMP1, A(p,q) )
0979:  1949             CONTINUE
0980:  1947          CONTINUE
0981:             ELSE
0982:                CALL SLASET( 'U', NR-1, NR-1, ZERO, ZERO, A(1,2), LDA )
0983:             END IF
0984: *
0985: *           .. and one-sided Jacobi rotations are started on a lower
0986: *           triangular matrix (plus perturbation which is ignored in
0987: *           the part which destroys triangular form (confusing?!))
0988: *
0989:             CALL SGESVJ( 'L', 'NoU', 'NoV', NR, NR, A, LDA, SVA,
0990:      &                      N, V, LDV, WORK, LWORK, INFO )
0991: *
0992:             SCALEM  = WORK(1)
0993:             NUMRANK = NINT(WORK(2))
0994: *
0995: *
0996:       ELSE IF ( RSVEC .AND. ( .NOT. LSVEC ) ) THEN
0997: *
0998: *        -> Singular Values and Right Singular Vectors <-
0999: *
1000:          IF ( ALMORT ) THEN
1001: *
1002: *           .. in this case NR equals N
1003:             DO 1998 p = 1, NR
1004:                CALL SCOPY( N-p+1, A(p,p), LDA, V(p,p), 1 )
1005:  1998       CONTINUE
1006:             CALL SLASET( 'Upper', NR-1, NR-1, ZERO, ZERO, V(1,2), LDV )
1007: *
1008:             CALL SGESVJ( 'L','U','N', N, NR, V,LDV, SVA, NR, A,LDA,
1009:      &                  WORK, LWORK, INFO )
1010:             SCALEM  = WORK(1)
1011:             NUMRANK = NINT(WORK(2))
1012: 
1013:          ELSE
1014: *
1015: *        .. two more QR factorizations ( one QRF is not enough, two require
1016: *        accumulated product of Jacobi rotations, three are perfect )
1017: *
1018:             CALL SLASET( 'Lower', NR-1, NR-1, ZERO, ZERO, A(2,1), LDA )
1019:             CALL SGELQF( NR, N, A, LDA, WORK, WORK(N+1), LWORK-N, IERR)
1020:             CALL SLACPY( 'Lower', NR, NR, A, LDA, V, LDV )
1021:             CALL SLASET( 'Upper', NR-1, NR-1, ZERO, ZERO, V(1,2), LDV )
1022:             CALL SGEQRF( NR, NR, V, LDV, WORK(N+1), WORK(2*N+1),
1023:      &                   LWORK-2*N, IERR )
1024:             DO 8998 p = 1, NR
1025:                CALL SCOPY( NR-p+1, V(p,p), LDV, V(p,p), 1 )
1026:  8998       CONTINUE
1027:             CALL SLASET( 'Upper', NR-1, NR-1, ZERO, ZERO, V(1,2), LDV )
1028: *
1029:             CALL SGESVJ( 'Lower', 'U','N', NR, NR, V,LDV, SVA, NR, U,
1030:      &                  LDU, WORK(N+1), LWORK, INFO )
1031:             SCALEM  = WORK(N+1)
1032:             NUMRANK = NINT(WORK(N+2))
1033:             IF ( NR .LT. N ) THEN
1034:                CALL SLASET( 'A',N-NR, NR, ZERO,ZERO, V(NR+1,1),   LDV )
1035:                CALL SLASET( 'A',NR, N-NR, ZERO,ZERO, V(1,NR+1),   LDV )
1036:                CALL SLASET( 'A',N-NR,N-NR,ZERO,ONE, V(NR+1,NR+1), LDV )
1037:             END IF
1038: *
1039:          CALL SORMLQ( 'Left', 'Transpose', N, N, NR, A, LDA, WORK,
1040:      &               V, LDV, WORK(N+1), LWORK-N, IERR )
1041: *
1042:          END IF
1043: *
1044:          DO 8991 p = 1, N
1045:             CALL SCOPY( N, V(p,1), LDV, A(IWORK(p),1), LDA )
1046:  8991    CONTINUE
1047:          CALL SLACPY( 'All', N, N, A, LDA, V, LDV )
1048: *
1049:          IF ( TRANSP ) THEN
1050:             CALL SLACPY( 'All', N, N, V, LDV, U, LDU )
1051:          END IF
1052: *
1053:       ELSE IF ( LSVEC .AND. ( .NOT. RSVEC ) ) THEN
1054: *
1055: *        .. Singular Values and Left Singular Vectors                 ..
1056: *
1057: *        .. second preconditioning step to avoid need to accumulate
1058: *        Jacobi rotations in the Jacobi iterations.
1059:          DO 1965 p = 1, NR
1060:             CALL SCOPY( N-p+1, A(p,p), LDA, U(p,p), 1 )
1061:  1965    CONTINUE
1062:          CALL SLASET( 'Upper', NR-1, NR-1, ZERO, ZERO, U(1,2), LDU )
1063: *
1064:          CALL SGEQRF( N, NR, U, LDU, WORK(N+1), WORK(2*N+1),
1065:      &              LWORK-2*N, IERR )
1066: *
1067:          DO 1967 p = 1, NR - 1
1068:             CALL SCOPY( NR-p, U(p,p+1), LDU, U(p+1,p), 1 )
1069:  1967    CONTINUE
1070:          CALL SLASET( 'Upper', NR-1, NR-1, ZERO, ZERO, U(1,2), LDU )
1071: *
1072:          CALL SGESVJ( 'Lower', 'U', 'N', NR,NR, U, LDU, SVA, NR, A,
1073:      &        LDA, WORK(N+1), LWORK-N, INFO )
1074:          SCALEM  = WORK(N+1)
1075:          NUMRANK = NINT(WORK(N+2))
1076: *
1077:          IF ( NR .LT. M ) THEN
1078:             CALL SLASET( 'A',  M-NR, NR,ZERO, ZERO, U(NR+1,1), LDU )
1079:             IF ( NR .LT. N1 ) THEN
1080:                CALL SLASET( 'A',NR, N1-NR, ZERO, ZERO, U(1,NR+1), LDU )
1081:                CALL SLASET( 'A',M-NR,N1-NR,ZERO,ONE,U(NR+1,NR+1), LDU )
1082:             END IF
1083:          END IF
1084: *
1085:          CALL SORMQR( 'Left', 'No Tr', M, N1, N, A, LDA, WORK, U,
1086:      &               LDU, WORK(N+1), LWORK-N, IERR )
1087: *
1088:          IF ( ROWPIV )
1089:      &       CALL SLASWP( N1, U, LDU, 1, M-1, IWORK(2*N+1), -1 )
1090: *
1091:          DO 1974 p = 1, N1
1092:             XSC = ONE / SNRM2( M, U(1,p), 1 )
1093:             CALL SSCAL( M, XSC, U(1,p), 1 )
1094:  1974    CONTINUE
1095: *
1096:          IF ( TRANSP ) THEN
1097:             CALL SLACPY( 'All', N, N, U, LDU, V, LDV )
1098:          END IF
1099: *
1100:       ELSE
1101: *
1102: *        .. Full SVD ..
1103: *
1104:          IF ( .NOT. JRACC ) THEN
1105: *
1106:          IF ( .NOT. ALMORT ) THEN
1107: *
1108: *           Second Preconditioning Step (QRF [with pivoting])
1109: *           Note that the composition of TRANSPOSE, QRF and TRANSPOSE is
1110: *           equivalent to an LQF CALL. Since in many libraries the QRF
1111: *           seems to be better optimized than the LQF, we do explicit
1112: *           transpose and use the QRF. This is subject to changes in an
1113: *           optimized implementation of SGEJSV.
1114: *
1115:             DO 1968 p = 1, NR
1116:                CALL SCOPY( N-p+1, A(p,p), LDA, V(p,p), 1 )
1117:  1968       CONTINUE
1118: *
1119: *           .. the following two loops perturb small entries to avoid
1120: *           denormals in the second QR factorization, where they are
1121: *           as good as zeros. This is done to avoid painfully slow
1122: *           computation with denormals. The relative size of the perturbation
1123: *           is a parameter that can be changed by the implementer.
1124: *           This perturbation device will be obsolete on machines with
1125: *           properly implemented arithmetic.
1126: *           To switch it off, set L2PERT=.FALSE. To remove it from  the
1127: *           code, remove the action under L2PERT=.TRUE., leave the ELSE part.
1128: *           The following two loops should be blocked and fused with the
1129: *           transposed copy above.
1130: *
1131:             IF ( L2PERT ) THEN
1132:                XSC = SQRT(SMALL)
1133:                DO 2969 q = 1, NR
1134:                   TEMP1 = XSC*ABS( V(q,q) )
1135:                   DO 2968 p = 1, N
1136:                      IF ( ( p .GT. q ) .AND. ( ABS(V(p,q)) .LE. TEMP1 )
1137:      &                   .OR. ( p .LT. q ) )
1138:      &                   V(p,q) = SIGN( TEMP1, V(p,q) )
1139:                      IF ( p. LT. q ) V(p,q) = - V(p,q)
1140:  2968             CONTINUE
1141:  2969          CONTINUE
1142:             ELSE
1143:                CALL SLASET( 'U', NR-1, NR-1, ZERO, ZERO, V(1,2), LDV )
1144:             END IF
1145: *
1146: *           Estimate the row scaled condition number of R1
1147: *           (If R1 is rectangular, N > NR, then the condition number
1148: *           of the leading NR x NR submatrix is estimated.)
1149: *
1150:             CALL SLACPY( 'L', NR, NR, V, LDV, WORK(2*N+1), NR )
1151:             DO 3950 p = 1, NR
1152:                TEMP1 = SNRM2(NR-p+1,WORK(2*N+(p-1)*NR+p),1)
1153:                CALL SSCAL(NR-p+1,ONE/TEMP1,WORK(2*N+(p-1)*NR+p),1)
1154:  3950       CONTINUE
1155:             CALL SPOCON('Lower',NR,WORK(2*N+1),NR,ONE,TEMP1,
1156:      &                   WORK(2*N+NR*NR+1),IWORK(M+2*N+1),IERR)
1157:             CONDR1 = ONE / SQRT(TEMP1)
1158: *           .. here need a second oppinion on the condition number
1159: *           .. then assume worst case scenario
1160: *           R1 is OK for inverse <=> CONDR1 .LT. FLOAT(N)
1161: *           more conservative    <=> CONDR1 .LT. SQRT(FLOAT(N))
1162: *
1163:             COND_OK = SQRT(FLOAT(NR))
1164: *[TP]       COND_OK is a tuning parameter.
1165: 
1166:             IF ( CONDR1 .LT. COND_OK ) THEN
1167: *              .. the second QRF without pivoting. Note: in an optimized
1168: *              implementation, this QRF should be implemented as the QRF
1169: *              of a lower triangular matrix.
1170: *              R1^t = Q2 * R2
1171:                CALL SGEQRF( N, NR, V, LDV, WORK(N+1), WORK(2*N+1),
1172:      &              LWORK-2*N, IERR )
1173: *
1174:                IF ( L2PERT ) THEN
1175:                   XSC = SQRT(SMALL)/EPSLN
1176:                   DO 3959 p = 2, NR
1177:                      DO 3958 q = 1, p - 1
1178:                         TEMP1 = XSC * AMIN1(ABS(V(p,p)),ABS(V(q,q)))
1179:                         IF ( ABS(V(q,p)) .LE. TEMP1 )
1180:      &                     V(q,p) = SIGN( TEMP1, V(q,p) )
1181:  3958                CONTINUE
1182:  3959             CONTINUE
1183:                END IF
1184: *
1185:                IF ( NR .NE. N )
1186: *              .. save ...
1187:      &         CALL SLACPY( 'A', N, NR, V, LDV, WORK(2*N+1), N )
1188: *
1189: *           .. this transposed copy should be better than naive
1190:                DO 1969 p = 1, NR - 1
1191:                   CALL SCOPY( NR-p, V(p,p+1), LDV, V(p+1,p), 1 )
1192:  1969          CONTINUE
1193: *
1194:                CONDR2 = CONDR1
1195: *
1196:             ELSE
1197: *
1198: *              .. ill-conditioned case: second QRF with pivoting
1199: *              Note that windowed pivoting would be equaly good
1200: *              numerically, and more run-time efficient. So, in
1201: *              an optimal implementation, the next call to SGEQP3
1202: *              should be replaced with eg. CALL SGEQPX (ACM TOMS #782)
1203: *              with properly (carefully) chosen parameters.
1204: *
1205: *              R1^t * P2 = Q2 * R2
1206:                DO 3003 p = 1, NR
1207:                   IWORK(N+p) = 0
1208:  3003          CONTINUE
1209:                CALL SGEQP3( N, NR, V, LDV, IWORK(N+1), WORK(N+1),
1210:      &                  WORK(2*N+1), LWORK-2*N, IERR )
1211: **               CALL SGEQRF( N, NR, V, LDV, WORK(N+1), WORK(2*N+1),
1212: **     &              LWORK-2*N, IERR )
1213:                IF ( L2PERT ) THEN
1214:                   XSC = SQRT(SMALL)
1215:                   DO 3969 p = 2, NR
1216:                      DO 3968 q = 1, p - 1
1217:                         TEMP1 = XSC * AMIN1(ABS(V(p,p)),ABS(V(q,q)))
1218:                         IF ( ABS(V(q,p)) .LE. TEMP1 )
1219:      &                     V(q,p) = SIGN( TEMP1, V(q,p) )
1220:  3968                CONTINUE
1221:  3969             CONTINUE
1222:                END IF
1223: *
1224:                CALL SLACPY( 'A', N, NR, V, LDV, WORK(2*N+1), N )
1225: *
1226:                IF ( L2PERT ) THEN
1227:                   XSC = SQRT(SMALL)
1228:                   DO 8970 p = 2, NR
1229:                      DO 8971 q = 1, p - 1
1230:                         TEMP1 = XSC * AMIN1(ABS(V(p,p)),ABS(V(q,q)))
1231:                         V(p,q) = - SIGN( TEMP1, V(q,p) )
1232:  8971                CONTINUE
1233:  8970             CONTINUE
1234:                ELSE
1235:                   CALL SLASET( 'L',NR-1,NR-1,ZERO,ZERO,V(2,1),LDV )
1236:                END IF
1237: *              Now, compute R2 = L3 * Q3, the LQ factorization.
1238:                CALL SGELQF( NR, NR, V, LDV, WORK(2*N+N*NR+1),
1239:      &               WORK(2*N+N*NR+NR+1), LWORK-2*N-N*NR-NR, IERR )
1240: *              .. and estimate the condition number
1241:                CALL SLACPY( 'L',NR,NR,V,LDV,WORK(2*N+N*NR+NR+1),NR )
1242:                DO 4950 p = 1, NR
1243:                   TEMP1 = SNRM2( p, WORK(2*N+N*NR+NR+p), NR )
1244:                   CALL SSCAL( p, ONE/TEMP1, WORK(2*N+N*NR+NR+p), NR )
1245:  4950          CONTINUE
1246:                CALL SPOCON( 'L',NR,WORK(2*N+N*NR+NR+1),NR,ONE,TEMP1,
1247:      &              WORK(2*N+N*NR+NR+NR*NR+1),IWORK(M+2*N+1),IERR )
1248:                CONDR2 = ONE / SQRT(TEMP1)
1249: *
1250:                IF ( CONDR2 .GE. COND_OK ) THEN
1251: *                 .. save the Householder vectors used for Q3
1252: *                 (this overwrittes the copy of R2, as it will not be
1253: *                 needed in this branch, but it does not overwritte the
1254: *                 Huseholder vectors of Q2.).
1255:                   CALL SLACPY( 'U', NR, NR, V, LDV, WORK(2*N+1), N )
1256: *                 .. and the rest of the information on Q3 is in
1257: *                 WORK(2*N+N*NR+1:2*N+N*NR+N)
1258:                END IF
1259: *
1260:             END IF
1261: *
1262:             IF ( L2PERT ) THEN
1263:                XSC = SQRT(SMALL)
1264:                DO 4968 q = 2, NR
1265:                   TEMP1 = XSC * V(q,q)
1266:                   DO 4969 p = 1, q - 1
1267: *                    V(p,q) = - SIGN( TEMP1, V(q,p) )
1268:                      V(p,q) = - SIGN( TEMP1, V(p,q) )
1269:  4969             CONTINUE
1270:  4968          CONTINUE
1271:             ELSE
1272:                CALL SLASET( 'U', NR-1,NR-1, ZERO,ZERO, V(1,2), LDV )
1273:             END IF
1274: *
1275: *        Second preconditioning finished; continue with Jacobi SVD
1276: *        The input matrix is lower trinagular.
1277: *
1278: *        Recover the right singular vectors as solution of a well
1279: *        conditioned triangular matrix equation.
1280: *
1281:             IF ( CONDR1 .LT. COND_OK ) THEN
1282: *
1283:                CALL SGESVJ( 'L','U','N',NR,NR,V,LDV,SVA,NR,U,
1284:      &              LDU,WORK(2*N+N*NR+NR+1),LWORK-2*N-N*NR-NR,INFO )
1285:                SCALEM  = WORK(2*N+N*NR+NR+1)
1286:                NUMRANK = NINT(WORK(2*N+N*NR+NR+2))
1287:                DO 3970 p = 1, NR
1288:                   CALL SCOPY( NR, V(1,p), 1, U(1,p), 1 )
1289:                   CALL SSCAL( NR, SVA(p),    V(1,p), 1 )
1290:  3970          CONTINUE
1291: 
1292: *        .. pick the right matrix equation and solve it
1293: *
1294:                IF ( NR. EQ. N ) THEN
1295: * :))             .. best case, R1 is inverted. The solution of this matrix
1296: *                 equation is Q2*V2 = the product of the Jacobi rotations
1297: *                 used in SGESVJ, premultiplied with the orthogonal matrix
1298: *                 from the second QR factorization.
1299:                   CALL STRSM( 'L','U','N','N', NR,NR,ONE, A,LDA, V,LDV )
1300:                ELSE
1301: *                 .. R1 is well conditioned, but non-square. Transpose(R2)
1302: *                 is inverted to get the product of the Jacobi rotations
1303: *                 used in SGESVJ. The Q-factor from the second QR
1304: *                 factorization is then built in explicitly.
1305:                   CALL STRSM('L','U','T','N',NR,NR,ONE,WORK(2*N+1),
1306:      &                 N,V,LDV)
1307:                   IF ( NR .LT. N ) THEN
1308:                     CALL SLASET('A',N-NR,NR,ZERO,ZERO,V(NR+1,1),LDV)
1309:                     CALL SLASET('A',NR,N-NR,ZERO,ZERO,V(1,NR+1),LDV)
1310:                     CALL SLASET('A',N-NR,N-NR,ZERO,ONE,V(NR+1,NR+1),LDV)
1311:                   END IF
1312:                   CALL SORMQR('L','N',N,N,NR,WORK(2*N+1),N,WORK(N+1),
1313:      &                 V,LDV,WORK(2*N+N*NR+NR+1),LWORK-2*N-N*NR-NR,IERR)
1314:                END IF
1315: *
1316:             ELSE IF ( CONDR2 .LT. COND_OK ) THEN
1317: *
1318: * :)           .. the input matrix A is very likely a relative of
1319: *              the Kahan matrix :)
1320: *              The matrix R2 is inverted. The solution of the matrix equation
1321: *              is Q3^T*V3 = the product of the Jacobi rotations (appplied to
1322: *              the lower triangular L3 from the LQ factorization of
1323: *              R2=L3*Q3), pre-multiplied with the transposed Q3.
1324:                CALL SGESVJ( 'L', 'U', 'N', NR, NR, V, LDV, SVA, NR, U,
1325:      &              LDU, WORK(2*N+N*NR+NR+1), LWORK-2*N-N*NR-NR, INFO )
1326:                SCALEM  = WORK(2*N+N*NR+NR+1)
1327:                NUMRANK = NINT(WORK(2*N+N*NR+NR+2))
1328:                DO 3870 p = 1, NR
1329:                   CALL SCOPY( NR, V(1,p), 1, U(1,p), 1 )
1330:                   CALL SSCAL( NR, SVA(p),    U(1,p), 1 )
1331:  3870          CONTINUE
1332:                CALL STRSM('L','U','N','N',NR,NR,ONE,WORK(2*N+1),N,U,LDU)
1333: *              .. apply the permutation from the second QR factorization
1334:                DO 873 q = 1, NR
1335:                   DO 872 p = 1, NR
1336:                      WORK(2*N+N*NR+NR+IWORK(N+p)) = U(p,q)
1337:  872              CONTINUE
1338:                   DO 874 p = 1, NR
1339:                      U(p,q) = WORK(2*N+N*NR+NR+p)
1340:  874              CONTINUE
1341:  873           CONTINUE
1342:                IF ( NR .LT. N ) THEN
1343:                   CALL SLASET( 'A',N-NR,NR,ZERO,ZERO,V(NR+1,1),LDV )
1344:                   CALL SLASET( 'A',NR,N-NR,ZERO,ZERO,V(1,NR+1),LDV )
1345:                   CALL SLASET( 'A',N-NR,N-NR,ZERO,ONE,V(NR+1,NR+1),LDV )
1346:                END IF
1347:                CALL SORMQR( 'L','N',N,N,NR,WORK(2*N+1),N,WORK(N+1),
1348:      &              V,LDV,WORK(2*N+N*NR+NR+1),LWORK-2*N-N*NR-NR,IERR )
1349:             ELSE
1350: *              Last line of defense.
1351: * #:(          This is a rather pathological case: no scaled condition
1352: *              improvement after two pivoted QR factorizations. Other
1353: *              possibility is that the rank revealing QR factorization
1354: *              or the condition estimator has failed, or the COND_OK
1355: *              is set very close to ONE (which is unnecessary). Normally,
1356: *              this branch should never be executed, but in rare cases of
1357: *              failure of the RRQR or condition estimator, the last line of
1358: *              defense ensures that SGEJSV completes the task.
1359: *              Compute the full SVD of L3 using SGESVJ with explicit
1360: *              accumulation of Jacobi rotations.
1361:                CALL SGESVJ( 'L', 'U', 'V', NR, NR, V, LDV, SVA, NR, U,
1362:      &              LDU, WORK(2*N+N*NR+NR+1), LWORK-2*N-N*NR-NR, INFO )
1363:                SCALEM  = WORK(2*N+N*NR+NR+1)
1364:                NUMRANK = NINT(WORK(2*N+N*NR+NR+2))
1365:                IF ( NR .LT. N ) THEN
1366:                   CALL SLASET( 'A',N-NR,NR,ZERO,ZERO,V(NR+1,1),LDV )
1367:                   CALL SLASET( 'A',NR,N-NR,ZERO,ZERO,V(1,NR+1),LDV )
1368:                   CALL SLASET( 'A',N-NR,N-NR,ZERO,ONE,V(NR+1,NR+1),LDV )
1369:                END IF
1370:                CALL SORMQR( 'L','N',N,N,NR,WORK(2*N+1),N,WORK(N+1),
1371:      &              V,LDV,WORK(2*N+N*NR+NR+1),LWORK-2*N-N*NR-NR,IERR )
1372: *
1373:                CALL SORMLQ( 'L', 'T', NR, NR, NR, WORK(2*N+1), N,
1374:      &              WORK(2*N+N*NR+1), U, LDU, WORK(2*N+N*NR+NR+1),
1375:      &              LWORK-2*N-N*NR-NR, IERR )
1376:                DO 773 q = 1, NR
1377:                   DO 772 p = 1, NR
1378:                      WORK(2*N+N*NR+NR+IWORK(N+p)) = U(p,q)
1379:  772              CONTINUE
1380:                   DO 774 p = 1, NR
1381:                      U(p,q) = WORK(2*N+N*NR+NR+p)
1382:  774              CONTINUE
1383:  773           CONTINUE
1384: *
1385:             END IF
1386: *
1387: *           Permute the rows of V using the (column) permutation from the
1388: *           first QRF. Also, scale the columns to make them unit in
1389: *           Euclidean norm. This applies to all cases.
1390: *
1391:             TEMP1 = SQRT(FLOAT(N)) * EPSLN
1392:             DO 1972 q = 1, N
1393:                DO 972 p = 1, N
1394:                   WORK(2*N+N*NR+NR+IWORK(p)) = V(p,q)
1395:   972          CONTINUE
1396:                DO 973 p = 1, N
1397:                   V(p,q) = WORK(2*N+N*NR+NR+p)
1398:   973          CONTINUE
1399:                XSC = ONE / SNRM2( N, V(1,q), 1 )
1400:                IF ( (XSC .LT. (ONE-TEMP1)) .OR. (XSC .GT. (ONE+TEMP1)) )
1401:      &           CALL SSCAL( N, XSC, V(1,q), 1 )
1402:  1972       CONTINUE
1403: *           At this moment, V contains the right singular vectors of A.
1404: *           Next, assemble the left singular vector matrix U (M x N).
1405:             IF ( NR .LT. M ) THEN
1406:                CALL SLASET( 'A', M-NR, NR, ZERO, ZERO, U(NR+1,1), LDU )
1407:                IF ( NR .LT. N1 ) THEN
1408:                   CALL SLASET('A',NR,N1-NR,ZERO,ZERO,U(1,NR+1),LDU)
1409:                   CALL SLASET('A',M-NR,N1-NR,ZERO,ONE,U(NR+1,NR+1),LDU)
1410:                END IF
1411:             END IF
1412: *
1413: *           The Q matrix from the first QRF is built into the left singular
1414: *           matrix U. This applies to all cases.
1415: *
1416:             CALL SORMQR( 'Left', 'No_Tr', M, N1, N, A, LDA, WORK, U,
1417:      &           LDU, WORK(N+1), LWORK-N, IERR )
1418: 
1419: *           The columns of U are normalized. The cost is O(M*N) flops.
1420:             TEMP1 = SQRT(FLOAT(M)) * EPSLN
1421:             DO 1973 p = 1, NR
1422:                XSC = ONE / SNRM2( M, U(1,p), 1 )
1423:                IF ( (XSC .LT. (ONE-TEMP1)) .OR. (XSC .GT. (ONE+TEMP1)) )
1424:      &          CALL SSCAL( M, XSC, U(1,p), 1 )
1425:  1973       CONTINUE
1426: *
1427: *           If the initial QRF is computed with row pivoting, the left
1428: *           singular vectors must be adjusted.
1429: *
1430:             IF ( ROWPIV )
1431:      &          CALL SLASWP( N1, U, LDU, 1, M-1, IWORK(2*N+1), -1 )
1432: *
1433:          ELSE
1434: *
1435: *        .. the initial matrix A has almost orthogonal columns and
1436: *        the second QRF is not needed
1437: *
1438:             CALL SLACPY( 'Upper', N, N, A, LDA, WORK(N+1), N )
1439:             IF ( L2PERT ) THEN
1440:                XSC = SQRT(SMALL)
1441:                DO 5970 p = 2, N
1442:                   TEMP1 = XSC * WORK( N + (p-1)*N + p )
1443:                   DO 5971 q = 1, p - 1
1444:                      WORK(N+(q-1)*N+p)=-SIGN(TEMP1,WORK(N+(p-1)*N+q))
1445:  5971             CONTINUE
1446:  5970          CONTINUE
1447:             ELSE
1448:                CALL SLASET( 'Lower',N-1,N-1,ZERO,ZERO,WORK(N+2),N )
1449:             END IF
1450: *
1451:             CALL SGESVJ( 'Upper', 'U', 'N', N, N, WORK(N+1), N, SVA,
1452:      &           N, U, LDU, WORK(N+N*N+1), LWORK-N-N*N, INFO )
1453: *
1454:             SCALEM  = WORK(N+N*N+1)
1455:             NUMRANK = NINT(WORK(N+N*N+2))
1456:             DO 6970 p = 1, N
1457:                CALL SCOPY( N, WORK(N+(p-1)*N+1), 1, U(1,p), 1 )
1458:                CALL SSCAL( N, SVA(p), WORK(N+(p-1)*N+1), 1 )
1459:  6970       CONTINUE
1460: *
1461:             CALL STRSM( 'Left', 'Upper', 'NoTrans', 'No UD', N, N,
1462:      &           ONE, A, LDA, WORK(N+1), N )
1463:             DO 6972 p = 1, N
1464:                CALL SCOPY( N, WORK(N+p), N, V(IWORK(p),1), LDV )
1465:  6972       CONTINUE
1466:             TEMP1 = SQRT(FLOAT(N))*EPSLN
1467:             DO 6971 p = 1, N
1468:                XSC = ONE / SNRM2( N, V(1,p), 1 )
1469:                IF ( (XSC .LT. (ONE-TEMP1)) .OR. (XSC .GT. (ONE+TEMP1)) )
1470:      &            CALL SSCAL( N, XSC, V(1,p), 1 )
1471:  6971       CONTINUE
1472: *
1473: *           Assemble the left singular vector matrix U (M x N).
1474: *
1475:             IF ( N .LT. M ) THEN
1476:                CALL SLASET( 'A',  M-N, N, ZERO, ZERO, U(NR+1,1), LDU )
1477:                IF ( N .LT. N1 ) THEN
1478:                   CALL SLASET( 'A',N,  N1-N, ZERO, ZERO,  U(1,N+1),LDU )
1479:                   CALL SLASET( 'A',M-N,N1-N, ZERO, ONE,U(NR+1,N+1),LDU )
1480:                END IF
1481:             END IF
1482:             CALL SORMQR( 'Left', 'No Tr', M, N1, N, A, LDA, WORK, U,
1483:      &           LDU, WORK(N+1), LWORK-N, IERR )
1484:             TEMP1 = SQRT(FLOAT(M))*EPSLN
1485:             DO 6973 p = 1, N1
1486:                XSC = ONE / SNRM2( M, U(1,p), 1 )
1487:                IF ( (XSC .LT. (ONE-TEMP1)) .OR. (XSC .GT. (ONE+TEMP1)) )
1488:      &            CALL SSCAL( M, XSC, U(1,p), 1 )
1489:  6973       CONTINUE
1490: *
1491:             IF ( ROWPIV )
1492:      &         CALL SLASWP( N1, U, LDU, 1, M-1, IWORK(2*N+1), -1 )
1493: *
1494:          END IF
1495: *
1496: *        end of the  >> almost orthogonal case <<  in the full SVD
1497: *
1498:          ELSE
1499: *
1500: *        This branch deploys a preconditioned Jacobi SVD with explicitly
1501: *        accumulated rotations. It is included as optional, mainly for
1502: *        experimental purposes. It does perfom well, and can also be used.
1503: *        In this implementation, this branch will be automatically activated
1504: *        if the  condition number sigma_max(A) / sigma_min(A) is predicted
1505: *        to be greater than the overflow threshold. This is because the
1506: *        a posteriori computation of the singular vectors assumes robust
1507: *        implementation of BLAS and some LAPACK procedures, capable of working
1508: *        in presence of extreme values. Since that is not always the case, ...
1509: *
1510:          DO 7968 p = 1, NR
1511:             CALL SCOPY( N-p+1, A(p,p), LDA, V(p,p), 1 )
1512:  7968    CONTINUE
1513: *
1514:          IF ( L2PERT ) THEN
1515:             XSC = SQRT(SMALL/EPSLN)
1516:             DO 5969 q = 1, NR
1517:                TEMP1 = XSC*ABS( V(q,q) )
1518:                DO 5968 p = 1, N
1519:                   IF ( ( p .GT. q ) .AND. ( ABS(V(p,q)) .LE. TEMP1 )
1520:      &                .OR. ( p .LT. q ) )
1521:      &                V(p,q) = SIGN( TEMP1, V(p,q) )
1522:                   IF ( p. LT. q ) V(p,q) = - V(p,q)
1523:  5968          CONTINUE
1524:  5969       CONTINUE
1525:          ELSE
1526:             CALL SLASET( 'U', NR-1, NR-1, ZERO, ZERO, V(1,2), LDV )
1527:          END IF
1528: 
1529:          CALL SGEQRF( N, NR, V, LDV, WORK(N+1), WORK(2*N+1),
1530:      &        LWORK-2*N, IERR )
1531:          CALL SLACPY( 'L', N, NR, V, LDV, WORK(2*N+1), N )
1532: *
1533:          DO 7969 p = 1, NR
1534:             CALL SCOPY( NR-p+1, V(p,p), LDV, U(p,p), 1 )
1535:  7969    CONTINUE
1536: 
1537:          IF ( L2PERT ) THEN
1538:             XSC = SQRT(SMALL/EPSLN)
1539:             DO 9970 q = 2, NR
1540:                DO 9971 p = 1, q - 1
1541:                   TEMP1 = XSC * AMIN1(ABS(U(p,p)),ABS(U(q,q)))
1542:                   U(p,q) = - SIGN( TEMP1, U(q,p) )
1543:  9971          CONTINUE
1544:  9970       CONTINUE
1545:          ELSE
1546:             CALL SLASET('U', NR-1, NR-1, ZERO, ZERO, U(1,2), LDU )
1547:          END IF
1548: 
1549:          CALL SGESVJ( 'L', 'U', 'V', NR, NR, U, LDU, SVA,
1550:      &        N, V, LDV, WORK(2*N+N*NR+1), LWORK-2*N-N*NR, INFO )
1551:          SCALEM  = WORK(2*N+N*NR+1)
1552:          NUMRANK = NINT(WORK(2*N+N*NR+2))
1553: 
1554:          IF ( NR .LT. N ) THEN
1555:             CALL SLASET( 'A',N-NR,NR,ZERO,ZERO,V(NR+1,1),LDV )
1556:             CALL SLASET( 'A',NR,N-NR,ZERO,ZERO,V(1,NR+1),LDV )
1557:             CALL SLASET( 'A',N-NR,N-NR,ZERO,ONE,V(NR+1,NR+1),LDV )
1558:          END IF
1559: 
1560:          CALL SORMQR( 'L','N',N,N,NR,WORK(2*N+1),N,WORK(N+1),
1561:      &        V,LDV,WORK(2*N+N*NR+NR+1),LWORK-2*N-N*NR-NR,IERR )
1562: *
1563: *           Permute the rows of V using the (column) permutation from the
1564: *           first QRF. Also, scale the columns to make them unit in
1565: *           Euclidean norm. This applies to all cases.
1566: *
1567:             TEMP1 = SQRT(FLOAT(N)) * EPSLN
1568:             DO 7972 q = 1, N
1569:                DO 8972 p = 1, N
1570:                   WORK(2*N+N*NR+NR+IWORK(p)) = V(p,q)
1571:  8972          CONTINUE
1572:                DO 8973 p = 1, N
1573:                   V(p,q) = WORK(2*N+N*NR+NR+p)
1574:  8973          CONTINUE
1575:                XSC = ONE / SNRM2( N, V(1,q), 1 )
1576:                IF ( (XSC .LT. (ONE-TEMP1)) .OR. (XSC .GT. (ONE+TEMP1)) )
1577:      &           CALL SSCAL( N, XSC, V(1,q), 1 )
1578:  7972       CONTINUE
1579: *
1580: *           At this moment, V contains the right singular vectors of A.
1581: *           Next, assemble the left singular vector matrix U (M x N).
1582: *
1583:          IF ( N .LT. M ) THEN
1584:             CALL SLASET( 'A',  M-N, N, ZERO, ZERO, U(NR+1,1), LDU )
1585:             IF ( N .LT. N1 ) THEN
1586:                CALL SLASET( 'A',N,  N1-N, ZERO, ZERO,  U(1,N+1),LDU )
1587:                CALL SLASET( 'A',M-N,N1-N, ZERO, ONE,U(NR+1,N+1),LDU )
1588:             END IF
1589:          END IF
1590: *
1591:          CALL SORMQR( 'Left', 'No Tr', M, N1, N, A, LDA, WORK, U,
1592:      &        LDU, WORK(N+1), LWORK-N, IERR )
1593: *
1594:             IF ( ROWPIV )
1595:      &         CALL SLASWP( N1, U, LDU, 1, M-1, IWORK(2*N+1), -1 )
1596: *
1597: *
1598:          END IF
1599:          IF ( TRANSP ) THEN
1600: *           .. swap U and V because the procedure worked on A^t
1601:             DO 6974 p = 1, N
1602:                CALL SSWAP( N, U(1,p), 1, V(1,p), 1 )
1603:  6974       CONTINUE
1604:          END IF
1605: *
1606:       END IF
1607: *     end of the full SVD
1608: *
1609: *     Undo scaling, if necessary (and possible)
1610: *
1611:       IF ( USCAL2 .LE. (BIG/SVA(1))*USCAL1 ) THEN
1612:          CALL SLASCL( 'G', 0, 0, USCAL1, USCAL2, NR, 1, SVA, N, IERR )
1613:          USCAL1 = ONE
1614:          USCAL2 = ONE
1615:       END IF
1616: *
1617:       IF ( NR .LT. N ) THEN
1618:          DO 3004 p = NR+1, N
1619:             SVA(p) = ZERO
1620:  3004    CONTINUE
1621:       END IF
1622: *
1623:       WORK(1) = USCAL2 * SCALEM
1624:       WORK(2) = USCAL1
1625:       IF ( ERREST ) WORK(3) = SCONDA
1626:       IF ( LSVEC .AND. RSVEC ) THEN
1627:          WORK(4) = CONDR1
1628:          WORK(5) = CONDR2
1629:       END IF
1630:       IF ( L2TRAN ) THEN
1631:          WORK(6) = ENTRA
1632:          WORK(7) = ENTRAT
1633:       END IF
1634: *
1635:       IWORK(1) = NR
1636:       IWORK(2) = NUMRANK
1637:       IWORK(3) = WARNING
1638: *
1639:       RETURN
1640: *     ..
1641: *     .. END OF SGEJSV
1642: *     ..
1643:       END
1644: *
1645: