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