Annotation of rpl/lapack/lapack/dgesvdq.f, revision 1.1
1.1 ! bertrand 1: *> \brief <b> DGESVDQ computes the singular value decomposition (SVD) with a QR-Preconditioned QR SVD Method for GE matrices</b>
! 2: *
! 3: * =========== DOCUMENTATION ===========
! 4: *
! 5: * Online html documentation available at
! 6: * http://www.netlib.org/lapack/explore-html/
! 7: *
! 8: *> \htmlonly
! 9: *> Download DGESVDQ + dependencies
! 10: *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dgesvdq.f">
! 11: *> [TGZ]</a>
! 12: *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dgesvdq.f">
! 13: *> [ZIP]</a>
! 14: *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dgesvdq.f">
! 15: *> [TXT]</a>
! 16: *> \endhtmlonly
! 17: *
! 18: * Definition:
! 19: * ===========
! 20: *
! 21: * SUBROUTINE DGESVDQ( JOBA, JOBP, JOBR, JOBU, JOBV, M, N, A, LDA,
! 22: * S, U, LDU, V, LDV, NUMRANK, IWORK, LIWORK,
! 23: * WORK, LWORK, RWORK, LRWORK, INFO )
! 24: *
! 25: * .. Scalar Arguments ..
! 26: * IMPLICIT NONE
! 27: * CHARACTER JOBA, JOBP, JOBR, JOBU, JOBV
! 28: * INTEGER M, N, LDA, LDU, LDV, NUMRANK, LIWORK, LWORK, LRWORK,
! 29: * INFO
! 30: * ..
! 31: * .. Array Arguments ..
! 32: * DOUBLE PRECISION A( LDA, * ), U( LDU, * ), V( LDV, * ), WORK( * )
! 33: * DOUBLE PRECISION S( * ), RWORK( * )
! 34: * INTEGER IWORK( * )
! 35: * ..
! 36: *
! 37: *
! 38: *> \par Purpose:
! 39: * =============
! 40: *>
! 41: *> \verbatim
! 42: *>
! 43: *> DGESVDQ computes the singular value decomposition (SVD) of a real
! 44: *> M-by-N matrix A, where M >= N. The SVD of A is written as
! 45: *> [++] [xx] [x0] [xx]
! 46: *> A = U * SIGMA * V^*, [++] = [xx] * [ox] * [xx]
! 47: *> [++] [xx]
! 48: *> where SIGMA is an N-by-N diagonal matrix, U is an M-by-N orthonormal
! 49: *> matrix, and V is an N-by-N orthogonal matrix. The diagonal elements
! 50: *> of SIGMA are the singular values of A. The columns of U and V are the
! 51: *> left and the right singular vectors of A, respectively.
! 52: *> \endverbatim
! 53: *
! 54: * Arguments:
! 55: * ==========
! 56: *
! 57: *> \param[in] JOBA
! 58: *> \verbatim
! 59: *> JOBA is CHARACTER*1
! 60: *> Specifies the level of accuracy in the computed SVD
! 61: *> = 'A' The requested accuracy corresponds to having the backward
! 62: *> error bounded by || delta A ||_F <= f(m,n) * EPS * || A ||_F,
! 63: *> where EPS = DLAMCH('Epsilon'). This authorises DGESVDQ to
! 64: *> truncate the computed triangular factor in a rank revealing
! 65: *> QR factorization whenever the truncated part is below the
! 66: *> threshold of the order of EPS * ||A||_F. This is aggressive
! 67: *> truncation level.
! 68: *> = 'M' Similarly as with 'A', but the truncation is more gentle: it
! 69: *> is allowed only when there is a drop on the diagonal of the
! 70: *> triangular factor in the QR factorization. This is medium
! 71: *> truncation level.
! 72: *> = 'H' High accuracy requested. No numerical rank determination based
! 73: *> on the rank revealing QR factorization is attempted.
! 74: *> = 'E' Same as 'H', and in addition the condition number of column
! 75: *> scaled A is estimated and returned in RWORK(1).
! 76: *> N^(-1/4)*RWORK(1) <= ||pinv(A_scaled)||_2 <= N^(1/4)*RWORK(1)
! 77: *> \endverbatim
! 78: *>
! 79: *> \param[in] JOBP
! 80: *> \verbatim
! 81: *> JOBP is CHARACTER*1
! 82: *> = 'P' The rows of A are ordered in decreasing order with respect to
! 83: *> ||A(i,:)||_\infty. This enhances numerical accuracy at the cost
! 84: *> of extra data movement. Recommended for numerical robustness.
! 85: *> = 'N' No row pivoting.
! 86: *> \endverbatim
! 87: *>
! 88: *> \param[in] JOBR
! 89: *> \verbatim
! 90: *> JOBR is CHARACTER*1
! 91: *> = 'T' After the initial pivoted QR factorization, DGESVD is applied to
! 92: *> the transposed R**T of the computed triangular factor R. This involves
! 93: *> some extra data movement (matrix transpositions). Useful for
! 94: *> experiments, research and development.
! 95: *> = 'N' The triangular factor R is given as input to DGESVD. This may be
! 96: *> preferred as it involves less data movement.
! 97: *> \endverbatim
! 98: *>
! 99: *> \param[in] JOBU
! 100: *> \verbatim
! 101: *> JOBU is CHARACTER*1
! 102: *> = 'A' All M left singular vectors are computed and returned in the
! 103: *> matrix U. See the description of U.
! 104: *> = 'S' or 'U' N = min(M,N) left singular vectors are computed and returned
! 105: *> in the matrix U. See the description of U.
! 106: *> = 'R' Numerical rank NUMRANK is determined and only NUMRANK left singular
! 107: *> vectors are computed and returned in the matrix U.
! 108: *> = 'F' The N left singular vectors are returned in factored form as the
! 109: *> product of the Q factor from the initial QR factorization and the
! 110: *> N left singular vectors of (R**T , 0)**T. If row pivoting is used,
! 111: *> then the necessary information on the row pivoting is stored in
! 112: *> IWORK(N+1:N+M-1).
! 113: *> = 'N' The left singular vectors are not computed.
! 114: *> \endverbatim
! 115: *>
! 116: *> \param[in] JOBV
! 117: *> \verbatim
! 118: *> JOBV is CHARACTER*1
! 119: *> = 'A', 'V' All N right singular vectors are computed and returned in
! 120: *> the matrix V.
! 121: *> = 'R' Numerical rank NUMRANK is determined and only NUMRANK right singular
! 122: *> vectors are computed and returned in the matrix V. This option is
! 123: *> allowed only if JOBU = 'R' or JOBU = 'N'; otherwise it is illegal.
! 124: *> = 'N' The right singular vectors are not computed.
! 125: *> \endverbatim
! 126: *>
! 127: *> \param[in] M
! 128: *> \verbatim
! 129: *> M is INTEGER
! 130: *> The number of rows of the input matrix A. M >= 0.
! 131: *> \endverbatim
! 132: *>
! 133: *> \param[in] N
! 134: *> \verbatim
! 135: *> N is INTEGER
! 136: *> The number of columns of the input matrix A. M >= N >= 0.
! 137: *> \endverbatim
! 138: *>
! 139: *> \param[in,out] A
! 140: *> \verbatim
! 141: *> A is DOUBLE PRECISION array of dimensions LDA x N
! 142: *> On entry, the input matrix A.
! 143: *> On exit, if JOBU .NE. 'N' or JOBV .NE. 'N', the lower triangle of A contains
! 144: *> the Householder vectors as stored by DGEQP3. If JOBU = 'F', these Householder
! 145: *> vectors together with WORK(1:N) can be used to restore the Q factors from
! 146: *> the initial pivoted QR factorization of A. See the description of U.
! 147: *> \endverbatim
! 148: *>
! 149: *> \param[in] LDA
! 150: *> \verbatim
! 151: *> LDA is INTEGER.
! 152: *> The leading dimension of the array A. LDA >= max(1,M).
! 153: *> \endverbatim
! 154: *>
! 155: *> \param[out] S
! 156: *> \verbatim
! 157: *> S is DOUBLE PRECISION array of dimension N.
! 158: *> The singular values of A, ordered so that S(i) >= S(i+1).
! 159: *> \endverbatim
! 160: *>
! 161: *> \param[out] U
! 162: *> \verbatim
! 163: *> U is DOUBLE PRECISION array, dimension
! 164: *> LDU x M if JOBU = 'A'; see the description of LDU. In this case,
! 165: *> on exit, U contains the M left singular vectors.
! 166: *> LDU x N if JOBU = 'S', 'U', 'R' ; see the description of LDU. In this
! 167: *> case, U contains the leading N or the leading NUMRANK left singular vectors.
! 168: *> LDU x N if JOBU = 'F' ; see the description of LDU. In this case U
! 169: *> contains N x N orthogonal matrix that can be used to form the left
! 170: *> singular vectors.
! 171: *> If JOBU = 'N', U is not referenced.
! 172: *> \endverbatim
! 173: *>
! 174: *> \param[in] LDU
! 175: *> \verbatim
! 176: *> LDU is INTEGER.
! 177: *> The leading dimension of the array U.
! 178: *> If JOBU = 'A', 'S', 'U', 'R', LDU >= max(1,M).
! 179: *> If JOBU = 'F', LDU >= max(1,N).
! 180: *> Otherwise, LDU >= 1.
! 181: *> \endverbatim
! 182: *>
! 183: *> \param[out] V
! 184: *> \verbatim
! 185: *> V is DOUBLE PRECISION array, dimension
! 186: *> LDV x N if JOBV = 'A', 'V', 'R' or if JOBA = 'E' .
! 187: *> If JOBV = 'A', or 'V', V contains the N-by-N orthogonal matrix V**T;
! 188: *> If JOBV = 'R', V contains the first NUMRANK rows of V**T (the right
! 189: *> singular vectors, stored rowwise, of the NUMRANK largest singular values).
! 190: *> If JOBV = 'N' and JOBA = 'E', V is used as a workspace.
! 191: *> If JOBV = 'N', and JOBA.NE.'E', V is not referenced.
! 192: *> \endverbatim
! 193: *>
! 194: *> \param[in] LDV
! 195: *> \verbatim
! 196: *> LDV is INTEGER
! 197: *> The leading dimension of the array V.
! 198: *> If JOBV = 'A', 'V', 'R', or JOBA = 'E', LDV >= max(1,N).
! 199: *> Otherwise, LDV >= 1.
! 200: *> \endverbatim
! 201: *>
! 202: *> \param[out] NUMRANK
! 203: *> \verbatim
! 204: *> NUMRANK is INTEGER
! 205: *> NUMRANK is the numerical rank first determined after the rank
! 206: *> revealing QR factorization, following the strategy specified by the
! 207: *> value of JOBA. If JOBV = 'R' and JOBU = 'R', only NUMRANK
! 208: *> leading singular values and vectors are then requested in the call
! 209: *> of DGESVD. The final value of NUMRANK might be further reduced if
! 210: *> some singular values are computed as zeros.
! 211: *> \endverbatim
! 212: *>
! 213: *> \param[out] IWORK
! 214: *> \verbatim
! 215: *> IWORK is INTEGER array, dimension (max(1, LIWORK)).
! 216: *> On exit, IWORK(1:N) contains column pivoting permutation of the
! 217: *> rank revealing QR factorization.
! 218: *> If JOBP = 'P', IWORK(N+1:N+M-1) contains the indices of the sequence
! 219: *> of row swaps used in row pivoting. These can be used to restore the
! 220: *> left singular vectors in the case JOBU = 'F'.
! 221: *>
! 222: *> If LIWORK, LWORK, or LRWORK = -1, then on exit, if INFO = 0,
! 223: *> LIWORK(1) returns the minimal LIWORK.
! 224: *> \endverbatim
! 225: *>
! 226: *> \param[in] LIWORK
! 227: *> \verbatim
! 228: *> LIWORK is INTEGER
! 229: *> The dimension of the array IWORK.
! 230: *> LIWORK >= N + M - 1, if JOBP = 'P' and JOBA .NE. 'E';
! 231: *> LIWORK >= N if JOBP = 'N' and JOBA .NE. 'E';
! 232: *> LIWORK >= N + M - 1 + N, if JOBP = 'P' and JOBA = 'E';
! 233: *> LIWORK >= N + N if JOBP = 'N' and JOBA = 'E'.
! 234: *
! 235: *> If LIWORK = -1, then a workspace query is assumed; the routine
! 236: *> only calculates and returns the optimal and minimal sizes
! 237: *> for the WORK, IWORK, and RWORK arrays, and no error
! 238: *> message related to LWORK is issued by XERBLA.
! 239: *> \endverbatim
! 240: *>
! 241: *> \param[out] WORK
! 242: *> \verbatim
! 243: *> WORK is DOUBLE PRECISION array, dimension (max(2, LWORK)), used as a workspace.
! 244: *> On exit, if, on entry, LWORK.NE.-1, WORK(1:N) contains parameters
! 245: *> needed to recover the Q factor from the QR factorization computed by
! 246: *> DGEQP3.
! 247: *>
! 248: *> If LIWORK, LWORK, or LRWORK = -1, then on exit, if INFO = 0,
! 249: *> WORK(1) returns the optimal LWORK, and
! 250: *> WORK(2) returns the minimal LWORK.
! 251: *> \endverbatim
! 252: *>
! 253: *> \param[in,out] LWORK
! 254: *> \verbatim
! 255: *> LWORK is INTEGER
! 256: *> The dimension of the array WORK. It is determined as follows:
! 257: *> Let LWQP3 = 3*N+1, LWCON = 3*N, and let
! 258: *> LWORQ = { MAX( N, 1 ), if JOBU = 'R', 'S', or 'U'
! 259: *> { MAX( M, 1 ), if JOBU = 'A'
! 260: *> LWSVD = MAX( 5*N, 1 )
! 261: *> LWLQF = MAX( N/2, 1 ), LWSVD2 = MAX( 5*(N/2), 1 ), LWORLQ = MAX( N, 1 ),
! 262: *> LWQRF = MAX( N/2, 1 ), LWORQ2 = MAX( N, 1 )
! 263: *> Then the minimal value of LWORK is:
! 264: *> = MAX( N + LWQP3, LWSVD ) if only the singular values are needed;
! 265: *> = MAX( N + LWQP3, LWCON, LWSVD ) if only the singular values are needed,
! 266: *> and a scaled condition estimate requested;
! 267: *>
! 268: *> = N + MAX( LWQP3, LWSVD, LWORQ ) if the singular values and the left
! 269: *> singular vectors are requested;
! 270: *> = N + MAX( LWQP3, LWCON, LWSVD, LWORQ ) if the singular values and the left
! 271: *> singular vectors are requested, and also
! 272: *> a scaled condition estimate requested;
! 273: *>
! 274: *> = N + MAX( LWQP3, LWSVD ) if the singular values and the right
! 275: *> singular vectors are requested;
! 276: *> = N + MAX( LWQP3, LWCON, LWSVD ) if the singular values and the right
! 277: *> singular vectors are requested, and also
! 278: *> a scaled condition etimate requested;
! 279: *>
! 280: *> = N + MAX( LWQP3, LWSVD, LWORQ ) if the full SVD is requested with JOBV = 'R';
! 281: *> independent of JOBR;
! 282: *> = N + MAX( LWQP3, LWCON, LWSVD, LWORQ ) if the full SVD is requested,
! 283: *> JOBV = 'R' and, also a scaled condition
! 284: *> estimate requested; independent of JOBR;
! 285: *> = MAX( N + MAX( LWQP3, LWSVD, LWORQ ),
! 286: *> N + MAX( LWQP3, N/2+LWLQF, N/2+LWSVD2, N/2+LWORLQ, LWORQ) ) if the
! 287: *> full SVD is requested with JOBV = 'A' or 'V', and
! 288: *> JOBR ='N'
! 289: *> = MAX( N + MAX( LWQP3, LWCON, LWSVD, LWORQ ),
! 290: *> N + MAX( LWQP3, LWCON, N/2+LWLQF, N/2+LWSVD2, N/2+LWORLQ, LWORQ ) )
! 291: *> if the full SVD is requested with JOBV = 'A' or 'V', and
! 292: *> JOBR ='N', and also a scaled condition number estimate
! 293: *> requested.
! 294: *> = MAX( N + MAX( LWQP3, LWSVD, LWORQ ),
! 295: *> N + MAX( LWQP3, N/2+LWQRF, N/2+LWSVD2, N/2+LWORQ2, LWORQ ) ) if the
! 296: *> full SVD is requested with JOBV = 'A', 'V', and JOBR ='T'
! 297: *> = MAX( N + MAX( LWQP3, LWCON, LWSVD, LWORQ ),
! 298: *> N + MAX( LWQP3, LWCON, N/2+LWQRF, N/2+LWSVD2, N/2+LWORQ2, LWORQ ) )
! 299: *> if the full SVD is requested with JOBV = 'A' or 'V', and
! 300: *> JOBR ='T', and also a scaled condition number estimate
! 301: *> requested.
! 302: *> Finally, LWORK must be at least two: LWORK = MAX( 2, LWORK ).
! 303: *>
! 304: *> If LWORK = -1, then a workspace query is assumed; the routine
! 305: *> only calculates and returns the optimal and minimal sizes
! 306: *> for the WORK, IWORK, and RWORK arrays, and no error
! 307: *> message related to LWORK is issued by XERBLA.
! 308: *> \endverbatim
! 309: *>
! 310: *> \param[out] RWORK
! 311: *> \verbatim
! 312: *> RWORK is DOUBLE PRECISION array, dimension (max(1, LRWORK)).
! 313: *> On exit,
! 314: *> 1. If JOBA = 'E', RWORK(1) contains an estimate of the condition
! 315: *> number of column scaled A. If A = C * D where D is diagonal and C
! 316: *> has unit columns in the Euclidean norm, then, assuming full column rank,
! 317: *> N^(-1/4) * RWORK(1) <= ||pinv(C)||_2 <= N^(1/4) * RWORK(1).
! 318: *> Otherwise, RWORK(1) = -1.
! 319: *> 2. RWORK(2) contains the number of singular values computed as
! 320: *> exact zeros in DGESVD applied to the upper triangular or trapeziodal
! 321: *> R (from the initial QR factorization). In case of early exit (no call to
! 322: *> DGESVD, such as in the case of zero matrix) RWORK(2) = -1.
! 323: *>
! 324: *> If LIWORK, LWORK, or LRWORK = -1, then on exit, if INFO = 0,
! 325: *> RWORK(1) returns the minimal LRWORK.
! 326: *> \endverbatim
! 327: *>
! 328: *> \param[in] LRWORK
! 329: *> \verbatim
! 330: *> LRWORK is INTEGER.
! 331: *> The dimension of the array RWORK.
! 332: *> If JOBP ='P', then LRWORK >= MAX(2, M).
! 333: *> Otherwise, LRWORK >= 2
! 334: *
! 335: *> If LRWORK = -1, then a workspace query is assumed; the routine
! 336: *> only calculates and returns the optimal and minimal sizes
! 337: *> for the WORK, IWORK, and RWORK arrays, and no error
! 338: *> message related to LWORK is issued by XERBLA.
! 339: *> \endverbatim
! 340: *>
! 341: *> \param[out] INFO
! 342: *> \verbatim
! 343: *> INFO is INTEGER
! 344: *> = 0: successful exit.
! 345: *> < 0: if INFO = -i, the i-th argument had an illegal value.
! 346: *> > 0: if DBDSQR did not converge, INFO specifies how many superdiagonals
! 347: *> of an intermediate bidiagonal form B (computed in DGESVD) did not
! 348: *> converge to zero.
! 349: *> \endverbatim
! 350: *
! 351: *> \par Further Details:
! 352: * ========================
! 353: *>
! 354: *> \verbatim
! 355: *>
! 356: *> 1. The data movement (matrix transpose) is coded using simple nested
! 357: *> DO-loops because BLAS and LAPACK do not provide corresponding subroutines.
! 358: *> Those DO-loops are easily identified in this source code - by the CONTINUE
! 359: *> statements labeled with 11**. In an optimized version of this code, the
! 360: *> nested DO loops should be replaced with calls to an optimized subroutine.
! 361: *> 2. This code scales A by 1/SQRT(M) if the largest ABS(A(i,j)) could cause
! 362: *> column norm overflow. This is the minial precaution and it is left to the
! 363: *> SVD routine (CGESVD) to do its own preemptive scaling if potential over-
! 364: *> or underflows are detected. To avoid repeated scanning of the array A,
! 365: *> an optimal implementation would do all necessary scaling before calling
! 366: *> CGESVD and the scaling in CGESVD can be switched off.
! 367: *> 3. Other comments related to code optimization are given in comments in the
! 368: *> code, enlosed in [[double brackets]].
! 369: *> \endverbatim
! 370: *
! 371: *> \par Bugs, examples and comments
! 372: * ===========================
! 373: *
! 374: *> \verbatim
! 375: *> Please report all bugs and send interesting examples and/or comments to
! 376: *> drmac@math.hr. Thank you.
! 377: *> \endverbatim
! 378: *
! 379: *> \par References
! 380: * ===============
! 381: *
! 382: *> \verbatim
! 383: *> [1] Zlatko Drmac, Algorithm 977: A QR-Preconditioned QR SVD Method for
! 384: *> Computing the SVD with High Accuracy. ACM Trans. Math. Softw.
! 385: *> 44(1): 11:1-11:30 (2017)
! 386: *>
! 387: *> SIGMA library, xGESVDQ section updated February 2016.
! 388: *> Developed and coded by Zlatko Drmac, Department of Mathematics
! 389: *> University of Zagreb, Croatia, drmac@math.hr
! 390: *> \endverbatim
! 391: *
! 392: *
! 393: *> \par Contributors:
! 394: * ==================
! 395: *>
! 396: *> \verbatim
! 397: *> Developed and coded by Zlatko Drmac, Department of Mathematics
! 398: *> University of Zagreb, Croatia, drmac@math.hr
! 399: *> \endverbatim
! 400: *
! 401: * Authors:
! 402: * ========
! 403: *
! 404: *> \author Univ. of Tennessee
! 405: *> \author Univ. of California Berkeley
! 406: *> \author Univ. of Colorado Denver
! 407: *> \author NAG Ltd.
! 408: *
! 409: *> \date November 2018
! 410: *
! 411: *> \ingroup doubleGEsing
! 412: *
! 413: * =====================================================================
! 414: SUBROUTINE DGESVDQ( JOBA, JOBP, JOBR, JOBU, JOBV, M, N, A, LDA,
! 415: $ S, U, LDU, V, LDV, NUMRANK, IWORK, LIWORK,
! 416: $ WORK, LWORK, RWORK, LRWORK, INFO )
! 417: * .. Scalar Arguments ..
! 418: IMPLICIT NONE
! 419: CHARACTER JOBA, JOBP, JOBR, JOBU, JOBV
! 420: INTEGER M, N, LDA, LDU, LDV, NUMRANK, LIWORK, LWORK, LRWORK,
! 421: $ INFO
! 422: * ..
! 423: * .. Array Arguments ..
! 424: DOUBLE PRECISION A( LDA, * ), U( LDU, * ), V( LDV, * ), WORK( * )
! 425: DOUBLE PRECISION S( * ), RWORK( * )
! 426: INTEGER IWORK( * )
! 427: *
! 428: * =====================================================================
! 429: *
! 430: * .. Parameters ..
! 431: DOUBLE PRECISION ZERO, ONE
! 432: PARAMETER ( ZERO = 0.0D0, ONE = 1.0D0 )
! 433: * .. Local Scalars ..
! 434: INTEGER IERR, IWOFF, NR, N1, OPTRATIO, p, q
! 435: INTEGER LWCON, LWQP3, LWRK_DGELQF, LWRK_DGESVD, LWRK_DGESVD2,
! 436: $ LWRK_DGEQP3, LWRK_DGEQRF, LWRK_DORMLQ, LWRK_DORMQR,
! 437: $ LWRK_DORMQR2, LWLQF, LWQRF, LWSVD, LWSVD2, LWORQ,
! 438: $ LWORQ2, LWORLQ, MINWRK, MINWRK2, OPTWRK, OPTWRK2,
! 439: $ IMINWRK, RMINWRK
! 440: LOGICAL ACCLA, ACCLM, ACCLH, ASCALED, CONDA, DNTWU, DNTWV,
! 441: $ LQUERY, LSVC0, LSVEC, ROWPRM, RSVEC, RTRANS, WNTUA,
! 442: $ WNTUF, WNTUR, WNTUS, WNTVA, WNTVR
! 443: DOUBLE PRECISION BIG, EPSLN, RTMP, SCONDA, SFMIN
! 444: * .. Local Arrays
! 445: DOUBLE PRECISION RDUMMY(1)
! 446: * ..
! 447: * .. External Subroutines (BLAS, LAPACK)
! 448: EXTERNAL DGELQF, DGEQP3, DGEQRF, DGESVD, DLACPY, DLAPMT,
! 449: $ DLASCL, DLASET, DLASWP, DSCAL, DPOCON, DORMLQ,
! 450: $ DORMQR, XERBLA
! 451: * ..
! 452: * .. External Functions (BLAS, LAPACK)
! 453: LOGICAL LSAME
! 454: INTEGER IDAMAX
! 455: DOUBLE PRECISION DLANGE, DNRM2, DLAMCH
! 456: EXTERNAL DLANGE, LSAME, IDAMAX, DNRM2, DLAMCH
! 457: * ..
! 458: * .. Intrinsic Functions ..
! 459: *
! 460: INTRINSIC ABS, MAX, MIN, DBLE, SQRT
! 461: *
! 462: * Test the input arguments
! 463: *
! 464: WNTUS = LSAME( JOBU, 'S' ) .OR. LSAME( JOBU, 'U' )
! 465: WNTUR = LSAME( JOBU, 'R' )
! 466: WNTUA = LSAME( JOBU, 'A' )
! 467: WNTUF = LSAME( JOBU, 'F' )
! 468: LSVC0 = WNTUS .OR. WNTUR .OR. WNTUA
! 469: LSVEC = LSVC0 .OR. WNTUF
! 470: DNTWU = LSAME( JOBU, 'N' )
! 471: *
! 472: WNTVR = LSAME( JOBV, 'R' )
! 473: WNTVA = LSAME( JOBV, 'A' ) .OR. LSAME( JOBV, 'V' )
! 474: RSVEC = WNTVR .OR. WNTVA
! 475: DNTWV = LSAME( JOBV, 'N' )
! 476: *
! 477: ACCLA = LSAME( JOBA, 'A' )
! 478: ACCLM = LSAME( JOBA, 'M' )
! 479: CONDA = LSAME( JOBA, 'E' )
! 480: ACCLH = LSAME( JOBA, 'H' ) .OR. CONDA
! 481: *
! 482: ROWPRM = LSAME( JOBP, 'P' )
! 483: RTRANS = LSAME( JOBR, 'T' )
! 484: *
! 485: IF ( ROWPRM ) THEN
! 486: IF ( CONDA ) THEN
! 487: IMINWRK = MAX( 1, N + M - 1 + N )
! 488: ELSE
! 489: IMINWRK = MAX( 1, N + M - 1 )
! 490: END IF
! 491: RMINWRK = MAX( 2, M )
! 492: ELSE
! 493: IF ( CONDA ) THEN
! 494: IMINWRK = MAX( 1, N + N )
! 495: ELSE
! 496: IMINWRK = MAX( 1, N )
! 497: END IF
! 498: RMINWRK = 2
! 499: END IF
! 500: LQUERY = (LIWORK .EQ. -1 .OR. LWORK .EQ. -1 .OR. LRWORK .EQ. -1)
! 501: INFO = 0
! 502: IF ( .NOT. ( ACCLA .OR. ACCLM .OR. ACCLH ) ) THEN
! 503: INFO = -1
! 504: ELSE IF ( .NOT.( ROWPRM .OR. LSAME( JOBP, 'N' ) ) ) THEN
! 505: INFO = -2
! 506: ELSE IF ( .NOT.( RTRANS .OR. LSAME( JOBR, 'N' ) ) ) THEN
! 507: INFO = -3
! 508: ELSE IF ( .NOT.( LSVEC .OR. DNTWU ) ) THEN
! 509: INFO = -4
! 510: ELSE IF ( WNTUR .AND. WNTVA ) THEN
! 511: INFO = -5
! 512: ELSE IF ( .NOT.( RSVEC .OR. DNTWV )) THEN
! 513: INFO = -5
! 514: ELSE IF ( M.LT.0 ) THEN
! 515: INFO = -6
! 516: ELSE IF ( ( N.LT.0 ) .OR. ( N.GT.M ) ) THEN
! 517: INFO = -7
! 518: ELSE IF ( LDA.LT.MAX( 1, M ) ) THEN
! 519: INFO = -9
! 520: ELSE IF ( LDU.LT.1 .OR. ( LSVC0 .AND. LDU.LT.M ) .OR.
! 521: $ ( WNTUF .AND. LDU.LT.N ) ) THEN
! 522: INFO = -12
! 523: ELSE IF ( LDV.LT.1 .OR. ( RSVEC .AND. LDV.LT.N ) .OR.
! 524: $ ( CONDA .AND. LDV.LT.N ) ) THEN
! 525: INFO = -14
! 526: ELSE IF ( LIWORK .LT. IMINWRK .AND. .NOT. LQUERY ) THEN
! 527: INFO = -17
! 528: END IF
! 529: *
! 530: *
! 531: IF ( INFO .EQ. 0 ) THEN
! 532: * .. compute the minimal and the optimal workspace lengths
! 533: * [[The expressions for computing the minimal and the optimal
! 534: * values of LWORK are written with a lot of redundancy and
! 535: * can be simplified. However, this detailed form is easier for
! 536: * maintenance and modifications of the code.]]
! 537: *
! 538: * .. minimal workspace length for DGEQP3 of an M x N matrix
! 539: LWQP3 = 3 * N + 1
! 540: * .. minimal workspace length for DORMQR to build left singular vectors
! 541: IF ( WNTUS .OR. WNTUR ) THEN
! 542: LWORQ = MAX( N , 1 )
! 543: ELSE IF ( WNTUA ) THEN
! 544: LWORQ = MAX( M , 1 )
! 545: END IF
! 546: * .. minimal workspace length for DPOCON of an N x N matrix
! 547: LWCON = 3 * N
! 548: * .. DGESVD of an N x N matrix
! 549: LWSVD = MAX( 5 * N, 1 )
! 550: IF ( LQUERY ) THEN
! 551: CALL DGEQP3( M, N, A, LDA, IWORK, RDUMMY, RDUMMY, -1,
! 552: $ IERR )
! 553: LWRK_DGEQP3 = INT( RDUMMY(1) )
! 554: IF ( WNTUS .OR. WNTUR ) THEN
! 555: CALL DORMQR( 'L', 'N', M, N, N, A, LDA, RDUMMY, U,
! 556: $ LDU, RDUMMY, -1, IERR )
! 557: LWRK_DORMQR = INT( RDUMMY(1) )
! 558: ELSE IF ( WNTUA ) THEN
! 559: CALL DORMQR( 'L', 'N', M, M, N, A, LDA, RDUMMY, U,
! 560: $ LDU, RDUMMY, -1, IERR )
! 561: LWRK_DORMQR = INT( RDUMMY(1) )
! 562: ELSE
! 563: LWRK_DORMQR = 0
! 564: END IF
! 565: END IF
! 566: MINWRK = 2
! 567: OPTWRK = 2
! 568: IF ( .NOT. (LSVEC .OR. RSVEC )) THEN
! 569: * .. minimal and optimal sizes of the workspace if
! 570: * only the singular values are requested
! 571: IF ( CONDA ) THEN
! 572: MINWRK = MAX( N+LWQP3, LWCON, LWSVD )
! 573: ELSE
! 574: MINWRK = MAX( N+LWQP3, LWSVD )
! 575: END IF
! 576: IF ( LQUERY ) THEN
! 577: CALL DGESVD( 'N', 'N', N, N, A, LDA, S, U, LDU,
! 578: $ V, LDV, RDUMMY, -1, IERR )
! 579: LWRK_DGESVD = INT( RDUMMY(1) )
! 580: IF ( CONDA ) THEN
! 581: OPTWRK = MAX( N+LWRK_DGEQP3, N+LWCON, LWRK_DGESVD )
! 582: ELSE
! 583: OPTWRK = MAX( N+LWRK_DGEQP3, LWRK_DGESVD )
! 584: END IF
! 585: END IF
! 586: ELSE IF ( LSVEC .AND. (.NOT.RSVEC) ) THEN
! 587: * .. minimal and optimal sizes of the workspace if the
! 588: * singular values and the left singular vectors are requested
! 589: IF ( CONDA ) THEN
! 590: MINWRK = N + MAX( LWQP3, LWCON, LWSVD, LWORQ )
! 591: ELSE
! 592: MINWRK = N + MAX( LWQP3, LWSVD, LWORQ )
! 593: END IF
! 594: IF ( LQUERY ) THEN
! 595: IF ( RTRANS ) THEN
! 596: CALL DGESVD( 'N', 'O', N, N, A, LDA, S, U, LDU,
! 597: $ V, LDV, RDUMMY, -1, IERR )
! 598: ELSE
! 599: CALL DGESVD( 'O', 'N', N, N, A, LDA, S, U, LDU,
! 600: $ V, LDV, RDUMMY, -1, IERR )
! 601: END IF
! 602: LWRK_DGESVD = INT( RDUMMY(1) )
! 603: IF ( CONDA ) THEN
! 604: OPTWRK = N + MAX( LWRK_DGEQP3, LWCON, LWRK_DGESVD,
! 605: $ LWRK_DORMQR )
! 606: ELSE
! 607: OPTWRK = N + MAX( LWRK_DGEQP3, LWRK_DGESVD,
! 608: $ LWRK_DORMQR )
! 609: END IF
! 610: END IF
! 611: ELSE IF ( RSVEC .AND. (.NOT.LSVEC) ) THEN
! 612: * .. minimal and optimal sizes of the workspace if the
! 613: * singular values and the right singular vectors are requested
! 614: IF ( CONDA ) THEN
! 615: MINWRK = N + MAX( LWQP3, LWCON, LWSVD )
! 616: ELSE
! 617: MINWRK = N + MAX( LWQP3, LWSVD )
! 618: END IF
! 619: IF ( LQUERY ) THEN
! 620: IF ( RTRANS ) THEN
! 621: CALL DGESVD( 'O', 'N', N, N, A, LDA, S, U, LDU,
! 622: $ V, LDV, RDUMMY, -1, IERR )
! 623: ELSE
! 624: CALL DGESVD( 'N', 'O', N, N, A, LDA, S, U, LDU,
! 625: $ V, LDV, RDUMMY, -1, IERR )
! 626: END IF
! 627: LWRK_DGESVD = INT( RDUMMY(1) )
! 628: IF ( CONDA ) THEN
! 629: OPTWRK = N + MAX( LWRK_DGEQP3, LWCON, LWRK_DGESVD )
! 630: ELSE
! 631: OPTWRK = N + MAX( LWRK_DGEQP3, LWRK_DGESVD )
! 632: END IF
! 633: END IF
! 634: ELSE
! 635: * .. minimal and optimal sizes of the workspace if the
! 636: * full SVD is requested
! 637: IF ( RTRANS ) THEN
! 638: MINWRK = MAX( LWQP3, LWSVD, LWORQ )
! 639: IF ( CONDA ) MINWRK = MAX( MINWRK, LWCON )
! 640: MINWRK = MINWRK + N
! 641: IF ( WNTVA ) THEN
! 642: * .. minimal workspace length for N x N/2 DGEQRF
! 643: LWQRF = MAX( N/2, 1 )
! 644: * .. minimal workspace lengt for N/2 x N/2 DGESVD
! 645: LWSVD2 = MAX( 5 * (N/2), 1 )
! 646: LWORQ2 = MAX( N, 1 )
! 647: MINWRK2 = MAX( LWQP3, N/2+LWQRF, N/2+LWSVD2,
! 648: $ N/2+LWORQ2, LWORQ )
! 649: IF ( CONDA ) MINWRK2 = MAX( MINWRK2, LWCON )
! 650: MINWRK2 = N + MINWRK2
! 651: MINWRK = MAX( MINWRK, MINWRK2 )
! 652: END IF
! 653: ELSE
! 654: MINWRK = MAX( LWQP3, LWSVD, LWORQ )
! 655: IF ( CONDA ) MINWRK = MAX( MINWRK, LWCON )
! 656: MINWRK = MINWRK + N
! 657: IF ( WNTVA ) THEN
! 658: * .. minimal workspace length for N/2 x N DGELQF
! 659: LWLQF = MAX( N/2, 1 )
! 660: LWSVD2 = MAX( 5 * (N/2), 1 )
! 661: LWORLQ = MAX( N , 1 )
! 662: MINWRK2 = MAX( LWQP3, N/2+LWLQF, N/2+LWSVD2,
! 663: $ N/2+LWORLQ, LWORQ )
! 664: IF ( CONDA ) MINWRK2 = MAX( MINWRK2, LWCON )
! 665: MINWRK2 = N + MINWRK2
! 666: MINWRK = MAX( MINWRK, MINWRK2 )
! 667: END IF
! 668: END IF
! 669: IF ( LQUERY ) THEN
! 670: IF ( RTRANS ) THEN
! 671: CALL DGESVD( 'O', 'A', N, N, A, LDA, S, U, LDU,
! 672: $ V, LDV, RDUMMY, -1, IERR )
! 673: LWRK_DGESVD = INT( RDUMMY(1) )
! 674: OPTWRK = MAX(LWRK_DGEQP3,LWRK_DGESVD,LWRK_DORMQR)
! 675: IF ( CONDA ) OPTWRK = MAX( OPTWRK, LWCON )
! 676: OPTWRK = N + OPTWRK
! 677: IF ( WNTVA ) THEN
! 678: CALL DGEQRF(N,N/2,U,LDU,RDUMMY,RDUMMY,-1,IERR)
! 679: LWRK_DGEQRF = INT( RDUMMY(1) )
! 680: CALL DGESVD( 'S', 'O', N/2,N/2, V,LDV, S, U,LDU,
! 681: $ V, LDV, RDUMMY, -1, IERR )
! 682: LWRK_DGESVD2 = INT( RDUMMY(1) )
! 683: CALL DORMQR( 'R', 'C', N, N, N/2, U, LDU, RDUMMY,
! 684: $ V, LDV, RDUMMY, -1, IERR )
! 685: LWRK_DORMQR2 = INT( RDUMMY(1) )
! 686: OPTWRK2 = MAX( LWRK_DGEQP3, N/2+LWRK_DGEQRF,
! 687: $ N/2+LWRK_DGESVD2, N/2+LWRK_DORMQR2 )
! 688: IF ( CONDA ) OPTWRK2 = MAX( OPTWRK2, LWCON )
! 689: OPTWRK2 = N + OPTWRK2
! 690: OPTWRK = MAX( OPTWRK, OPTWRK2 )
! 691: END IF
! 692: ELSE
! 693: CALL DGESVD( 'S', 'O', N, N, A, LDA, S, U, LDU,
! 694: $ V, LDV, RDUMMY, -1, IERR )
! 695: LWRK_DGESVD = INT( RDUMMY(1) )
! 696: OPTWRK = MAX(LWRK_DGEQP3,LWRK_DGESVD,LWRK_DORMQR)
! 697: IF ( CONDA ) OPTWRK = MAX( OPTWRK, LWCON )
! 698: OPTWRK = N + OPTWRK
! 699: IF ( WNTVA ) THEN
! 700: CALL DGELQF(N/2,N,U,LDU,RDUMMY,RDUMMY,-1,IERR)
! 701: LWRK_DGELQF = INT( RDUMMY(1) )
! 702: CALL DGESVD( 'S','O', N/2,N/2, V, LDV, S, U, LDU,
! 703: $ V, LDV, RDUMMY, -1, IERR )
! 704: LWRK_DGESVD2 = INT( RDUMMY(1) )
! 705: CALL DORMLQ( 'R', 'N', N, N, N/2, U, LDU, RDUMMY,
! 706: $ V, LDV, RDUMMY,-1,IERR )
! 707: LWRK_DORMLQ = INT( RDUMMY(1) )
! 708: OPTWRK2 = MAX( LWRK_DGEQP3, N/2+LWRK_DGELQF,
! 709: $ N/2+LWRK_DGESVD2, N/2+LWRK_DORMLQ )
! 710: IF ( CONDA ) OPTWRK2 = MAX( OPTWRK2, LWCON )
! 711: OPTWRK2 = N + OPTWRK2
! 712: OPTWRK = MAX( OPTWRK, OPTWRK2 )
! 713: END IF
! 714: END IF
! 715: END IF
! 716: END IF
! 717: *
! 718: MINWRK = MAX( 2, MINWRK )
! 719: OPTWRK = MAX( 2, OPTWRK )
! 720: IF ( LWORK .LT. MINWRK .AND. (.NOT.LQUERY) ) INFO = -19
! 721: *
! 722: END IF
! 723: *
! 724: IF (INFO .EQ. 0 .AND. LRWORK .LT. RMINWRK .AND. .NOT. LQUERY) THEN
! 725: INFO = -21
! 726: END IF
! 727: IF( INFO.NE.0 ) THEN
! 728: CALL XERBLA( 'DGESVDQ', -INFO )
! 729: RETURN
! 730: ELSE IF ( LQUERY ) THEN
! 731: *
! 732: * Return optimal workspace
! 733: *
! 734: IWORK(1) = IMINWRK
! 735: WORK(1) = OPTWRK
! 736: WORK(2) = MINWRK
! 737: RWORK(1) = RMINWRK
! 738: RETURN
! 739: END IF
! 740: *
! 741: * Quick return if the matrix is void.
! 742: *
! 743: IF( ( M.EQ.0 ) .OR. ( N.EQ.0 ) ) THEN
! 744: * .. all output is void.
! 745: RETURN
! 746: END IF
! 747: *
! 748: BIG = DLAMCH('O')
! 749: ASCALED = .FALSE.
! 750: IWOFF = 1
! 751: IF ( ROWPRM ) THEN
! 752: IWOFF = M
! 753: * .. reordering the rows in decreasing sequence in the
! 754: * ell-infinity norm - this enhances numerical robustness in
! 755: * the case of differently scaled rows.
! 756: DO 1904 p = 1, M
! 757: * RWORK(p) = ABS( A(p,ICAMAX(N,A(p,1),LDA)) )
! 758: * [[DLANGE will return NaN if an entry of the p-th row is Nan]]
! 759: RWORK(p) = DLANGE( 'M', 1, N, A(p,1), LDA, RDUMMY )
! 760: * .. check for NaN's and Inf's
! 761: IF ( ( RWORK(p) .NE. RWORK(p) ) .OR.
! 762: $ ( (RWORK(p)*ZERO) .NE. ZERO ) ) THEN
! 763: INFO = -8
! 764: CALL XERBLA( 'DGESVDQ', -INFO )
! 765: RETURN
! 766: END IF
! 767: 1904 CONTINUE
! 768: DO 1952 p = 1, M - 1
! 769: q = IDAMAX( M-p+1, RWORK(p), 1 ) + p - 1
! 770: IWORK(N+p) = q
! 771: IF ( p .NE. q ) THEN
! 772: RTMP = RWORK(p)
! 773: RWORK(p) = RWORK(q)
! 774: RWORK(q) = RTMP
! 775: END IF
! 776: 1952 CONTINUE
! 777: *
! 778: IF ( RWORK(1) .EQ. ZERO ) THEN
! 779: * Quick return: A is the M x N zero matrix.
! 780: NUMRANK = 0
! 781: CALL DLASET( 'G', N, 1, ZERO, ZERO, S, N )
! 782: IF ( WNTUS ) CALL DLASET('G', M, N, ZERO, ONE, U, LDU)
! 783: IF ( WNTUA ) CALL DLASET('G', M, M, ZERO, ONE, U, LDU)
! 784: IF ( WNTVA ) CALL DLASET('G', N, N, ZERO, ONE, V, LDV)
! 785: IF ( WNTUF ) THEN
! 786: CALL DLASET( 'G', N, 1, ZERO, ZERO, WORK, N )
! 787: CALL DLASET( 'G', M, N, ZERO, ONE, U, LDU )
! 788: END IF
! 789: DO 5001 p = 1, N
! 790: IWORK(p) = p
! 791: 5001 CONTINUE
! 792: IF ( ROWPRM ) THEN
! 793: DO 5002 p = N + 1, N + M - 1
! 794: IWORK(p) = p - N
! 795: 5002 CONTINUE
! 796: END IF
! 797: IF ( CONDA ) RWORK(1) = -1
! 798: RWORK(2) = -1
! 799: RETURN
! 800: END IF
! 801: *
! 802: IF ( RWORK(1) .GT. BIG / SQRT(DBLE(M)) ) THEN
! 803: * .. to prevent overflow in the QR factorization, scale the
! 804: * matrix by 1/sqrt(M) if too large entry detected
! 805: CALL DLASCL('G',0,0,SQRT(DBLE(M)),ONE, M,N, A,LDA, IERR)
! 806: ASCALED = .TRUE.
! 807: END IF
! 808: CALL DLASWP( N, A, LDA, 1, M-1, IWORK(N+1), 1 )
! 809: END IF
! 810: *
! 811: * .. At this stage, preemptive scaling is done only to avoid column
! 812: * norms overflows during the QR factorization. The SVD procedure should
! 813: * have its own scaling to save the singular values from overflows and
! 814: * underflows. That depends on the SVD procedure.
! 815: *
! 816: IF ( .NOT.ROWPRM ) THEN
! 817: RTMP = DLANGE( 'M', M, N, A, LDA, RDUMMY )
! 818: IF ( ( RTMP .NE. RTMP ) .OR.
! 819: $ ( (RTMP*ZERO) .NE. ZERO ) ) THEN
! 820: INFO = -8
! 821: CALL XERBLA( 'DGESVDQ', -INFO )
! 822: RETURN
! 823: END IF
! 824: IF ( RTMP .GT. BIG / SQRT(DBLE(M)) ) THEN
! 825: * .. to prevent overflow in the QR factorization, scale the
! 826: * matrix by 1/sqrt(M) if too large entry detected
! 827: CALL DLASCL('G',0,0, SQRT(DBLE(M)),ONE, M,N, A,LDA, IERR)
! 828: ASCALED = .TRUE.
! 829: END IF
! 830: END IF
! 831: *
! 832: * .. QR factorization with column pivoting
! 833: *
! 834: * A * P = Q * [ R ]
! 835: * [ 0 ]
! 836: *
! 837: DO 1963 p = 1, N
! 838: * .. all columns are free columns
! 839: IWORK(p) = 0
! 840: 1963 CONTINUE
! 841: CALL DGEQP3( M, N, A, LDA, IWORK, WORK, WORK(N+1), LWORK-N,
! 842: $ IERR )
! 843: *
! 844: * If the user requested accuracy level allows truncation in the
! 845: * computed upper triangular factor, the matrix R is examined and,
! 846: * if possible, replaced with its leading upper trapezoidal part.
! 847: *
! 848: EPSLN = DLAMCH('E')
! 849: SFMIN = DLAMCH('S')
! 850: * SMALL = SFMIN / EPSLN
! 851: NR = N
! 852: *
! 853: IF ( ACCLA ) THEN
! 854: *
! 855: * Standard absolute error bound suffices. All sigma_i with
! 856: * sigma_i < N*EPS*||A||_F are flushed to zero. This is an
! 857: * aggressive enforcement of lower numerical rank by introducing a
! 858: * backward error of the order of N*EPS*||A||_F.
! 859: NR = 1
! 860: RTMP = SQRT(DBLE(N))*EPSLN
! 861: DO 3001 p = 2, N
! 862: IF ( ABS(A(p,p)) .LT. (RTMP*ABS(A(1,1))) ) GO TO 3002
! 863: NR = NR + 1
! 864: 3001 CONTINUE
! 865: 3002 CONTINUE
! 866: *
! 867: ELSEIF ( ACCLM ) THEN
! 868: * .. similarly as above, only slightly more gentle (less aggressive).
! 869: * Sudden drop on the diagonal of R is used as the criterion for being
! 870: * close-to-rank-deficient. The threshold is set to EPSLN=DLAMCH('E').
! 871: * [[This can be made more flexible by replacing this hard-coded value
! 872: * with a user specified threshold.]] Also, the values that underflow
! 873: * will be truncated.
! 874: NR = 1
! 875: DO 3401 p = 2, N
! 876: IF ( ( ABS(A(p,p)) .LT. (EPSLN*ABS(A(p-1,p-1))) ) .OR.
! 877: $ ( ABS(A(p,p)) .LT. SFMIN ) ) GO TO 3402
! 878: NR = NR + 1
! 879: 3401 CONTINUE
! 880: 3402 CONTINUE
! 881: *
! 882: ELSE
! 883: * .. RRQR not authorized to determine numerical rank except in the
! 884: * obvious case of zero pivots.
! 885: * .. inspect R for exact zeros on the diagonal;
! 886: * R(i,i)=0 => R(i:N,i:N)=0.
! 887: NR = 1
! 888: DO 3501 p = 2, N
! 889: IF ( ABS(A(p,p)) .EQ. ZERO ) GO TO 3502
! 890: NR = NR + 1
! 891: 3501 CONTINUE
! 892: 3502 CONTINUE
! 893: *
! 894: IF ( CONDA ) THEN
! 895: * Estimate the scaled condition number of A. Use the fact that it is
! 896: * the same as the scaled condition number of R.
! 897: * .. V is used as workspace
! 898: CALL DLACPY( 'U', N, N, A, LDA, V, LDV )
! 899: * Only the leading NR x NR submatrix of the triangular factor
! 900: * is considered. Only if NR=N will this give a reliable error
! 901: * bound. However, even for NR < N, this can be used on an
! 902: * expert level and obtain useful information in the sense of
! 903: * perturbation theory.
! 904: DO 3053 p = 1, NR
! 905: RTMP = DNRM2( p, V(1,p), 1 )
! 906: CALL DSCAL( p, ONE/RTMP, V(1,p), 1 )
! 907: 3053 CONTINUE
! 908: IF ( .NOT. ( LSVEC .OR. RSVEC ) ) THEN
! 909: CALL DPOCON( 'U', NR, V, LDV, ONE, RTMP,
! 910: $ WORK, IWORK(N+IWOFF), IERR )
! 911: ELSE
! 912: CALL DPOCON( 'U', NR, V, LDV, ONE, RTMP,
! 913: $ WORK(N+1), IWORK(N+IWOFF), IERR )
! 914: END IF
! 915: SCONDA = ONE / SQRT(RTMP)
! 916: * For NR=N, SCONDA is an estimate of SQRT(||(R^* * R)^(-1)||_1),
! 917: * N^(-1/4) * SCONDA <= ||R^(-1)||_2 <= N^(1/4) * SCONDA
! 918: * See the reference [1] for more details.
! 919: END IF
! 920: *
! 921: ENDIF
! 922: *
! 923: IF ( WNTUR ) THEN
! 924: N1 = NR
! 925: ELSE IF ( WNTUS .OR. WNTUF) THEN
! 926: N1 = N
! 927: ELSE IF ( WNTUA ) THEN
! 928: N1 = M
! 929: END IF
! 930: *
! 931: IF ( .NOT. ( RSVEC .OR. LSVEC ) ) THEN
! 932: *.......................................................................
! 933: * .. only the singular values are requested
! 934: *.......................................................................
! 935: IF ( RTRANS ) THEN
! 936: *
! 937: * .. compute the singular values of R**T = [A](1:NR,1:N)**T
! 938: * .. set the lower triangle of [A] to [A](1:NR,1:N)**T and
! 939: * the upper triangle of [A] to zero.
! 940: DO 1146 p = 1, MIN( N, NR )
! 941: DO 1147 q = p + 1, N
! 942: A(q,p) = A(p,q)
! 943: IF ( q .LE. NR ) A(p,q) = ZERO
! 944: 1147 CONTINUE
! 945: 1146 CONTINUE
! 946: *
! 947: CALL DGESVD( 'N', 'N', N, NR, A, LDA, S, U, LDU,
! 948: $ V, LDV, WORK, LWORK, INFO )
! 949: *
! 950: ELSE
! 951: *
! 952: * .. compute the singular values of R = [A](1:NR,1:N)
! 953: *
! 954: IF ( NR .GT. 1 )
! 955: $ CALL DLASET( 'L', NR-1,NR-1, ZERO,ZERO, A(2,1), LDA )
! 956: CALL DGESVD( 'N', 'N', NR, N, A, LDA, S, U, LDU,
! 957: $ V, LDV, WORK, LWORK, INFO )
! 958: *
! 959: END IF
! 960: *
! 961: ELSE IF ( LSVEC .AND. ( .NOT. RSVEC) ) THEN
! 962: *.......................................................................
! 963: * .. the singular values and the left singular vectors requested
! 964: *.......................................................................""""""""
! 965: IF ( RTRANS ) THEN
! 966: * .. apply DGESVD to R**T
! 967: * .. copy R**T into [U] and overwrite [U] with the right singular
! 968: * vectors of R
! 969: DO 1192 p = 1, NR
! 970: DO 1193 q = p, N
! 971: U(q,p) = A(p,q)
! 972: 1193 CONTINUE
! 973: 1192 CONTINUE
! 974: IF ( NR .GT. 1 )
! 975: $ CALL DLASET( 'U', NR-1,NR-1, ZERO,ZERO, U(1,2), LDU )
! 976: * .. the left singular vectors not computed, the NR right singular
! 977: * vectors overwrite [U](1:NR,1:NR) as transposed. These
! 978: * will be pre-multiplied by Q to build the left singular vectors of A.
! 979: CALL DGESVD( 'N', 'O', N, NR, U, LDU, S, U, LDU,
! 980: $ U, LDU, WORK(N+1), LWORK-N, INFO )
! 981: *
! 982: DO 1119 p = 1, NR
! 983: DO 1120 q = p + 1, NR
! 984: RTMP = U(q,p)
! 985: U(q,p) = U(p,q)
! 986: U(p,q) = RTMP
! 987: 1120 CONTINUE
! 988: 1119 CONTINUE
! 989: *
! 990: ELSE
! 991: * .. apply DGESVD to R
! 992: * .. copy R into [U] and overwrite [U] with the left singular vectors
! 993: CALL DLACPY( 'U', NR, N, A, LDA, U, LDU )
! 994: IF ( NR .GT. 1 )
! 995: $ CALL DLASET( 'L', NR-1, NR-1, ZERO, ZERO, U(2,1), LDU )
! 996: * .. the right singular vectors not computed, the NR left singular
! 997: * vectors overwrite [U](1:NR,1:NR)
! 998: CALL DGESVD( 'O', 'N', NR, N, U, LDU, S, U, LDU,
! 999: $ V, LDV, WORK(N+1), LWORK-N, INFO )
! 1000: * .. now [U](1:NR,1:NR) contains the NR left singular vectors of
! 1001: * R. These will be pre-multiplied by Q to build the left singular
! 1002: * vectors of A.
! 1003: END IF
! 1004: *
! 1005: * .. assemble the left singular vector matrix U of dimensions
! 1006: * (M x NR) or (M x N) or (M x M).
! 1007: IF ( ( NR .LT. M ) .AND. ( .NOT.WNTUF ) ) THEN
! 1008: CALL DLASET('A', M-NR, NR, ZERO, ZERO, U(NR+1,1), LDU)
! 1009: IF ( NR .LT. N1 ) THEN
! 1010: CALL DLASET( 'A',NR,N1-NR,ZERO,ZERO,U(1,NR+1), LDU )
! 1011: CALL DLASET( 'A',M-NR,N1-NR,ZERO,ONE,
! 1012: $ U(NR+1,NR+1), LDU )
! 1013: END IF
! 1014: END IF
! 1015: *
! 1016: * The Q matrix from the first QRF is built into the left singular
! 1017: * vectors matrix U.
! 1018: *
! 1019: IF ( .NOT.WNTUF )
! 1020: $ CALL DORMQR( 'L', 'N', M, N1, N, A, LDA, WORK, U,
! 1021: $ LDU, WORK(N+1), LWORK-N, IERR )
! 1022: IF ( ROWPRM .AND. .NOT.WNTUF )
! 1023: $ CALL DLASWP( N1, U, LDU, 1, M-1, IWORK(N+1), -1 )
! 1024: *
! 1025: ELSE IF ( RSVEC .AND. ( .NOT. LSVEC ) ) THEN
! 1026: *.......................................................................
! 1027: * .. the singular values and the right singular vectors requested
! 1028: *.......................................................................
! 1029: IF ( RTRANS ) THEN
! 1030: * .. apply DGESVD to R**T
! 1031: * .. copy R**T into V and overwrite V with the left singular vectors
! 1032: DO 1165 p = 1, NR
! 1033: DO 1166 q = p, N
! 1034: V(q,p) = (A(p,q))
! 1035: 1166 CONTINUE
! 1036: 1165 CONTINUE
! 1037: IF ( NR .GT. 1 )
! 1038: $ CALL DLASET( 'U', NR-1,NR-1, ZERO,ZERO, V(1,2), LDV )
! 1039: * .. the left singular vectors of R**T overwrite V, the right singular
! 1040: * vectors not computed
! 1041: IF ( WNTVR .OR. ( NR .EQ. N ) ) THEN
! 1042: CALL DGESVD( 'O', 'N', N, NR, V, LDV, S, U, LDU,
! 1043: $ U, LDU, WORK(N+1), LWORK-N, INFO )
! 1044: *
! 1045: DO 1121 p = 1, NR
! 1046: DO 1122 q = p + 1, NR
! 1047: RTMP = V(q,p)
! 1048: V(q,p) = V(p,q)
! 1049: V(p,q) = RTMP
! 1050: 1122 CONTINUE
! 1051: 1121 CONTINUE
! 1052: *
! 1053: IF ( NR .LT. N ) THEN
! 1054: DO 1103 p = 1, NR
! 1055: DO 1104 q = NR + 1, N
! 1056: V(p,q) = V(q,p)
! 1057: 1104 CONTINUE
! 1058: 1103 CONTINUE
! 1059: END IF
! 1060: CALL DLAPMT( .FALSE., NR, N, V, LDV, IWORK )
! 1061: ELSE
! 1062: * .. need all N right singular vectors and NR < N
! 1063: * [!] This is simple implementation that augments [V](1:N,1:NR)
! 1064: * by padding a zero block. In the case NR << N, a more efficient
! 1065: * way is to first use the QR factorization. For more details
! 1066: * how to implement this, see the " FULL SVD " branch.
! 1067: CALL DLASET('G', N, N-NR, ZERO, ZERO, V(1,NR+1), LDV)
! 1068: CALL DGESVD( 'O', 'N', N, N, V, LDV, S, U, LDU,
! 1069: $ U, LDU, WORK(N+1), LWORK-N, INFO )
! 1070: *
! 1071: DO 1123 p = 1, N
! 1072: DO 1124 q = p + 1, N
! 1073: RTMP = V(q,p)
! 1074: V(q,p) = V(p,q)
! 1075: V(p,q) = RTMP
! 1076: 1124 CONTINUE
! 1077: 1123 CONTINUE
! 1078: CALL DLAPMT( .FALSE., N, N, V, LDV, IWORK )
! 1079: END IF
! 1080: *
! 1081: ELSE
! 1082: * .. aply DGESVD to R
! 1083: * .. copy R into V and overwrite V with the right singular vectors
! 1084: CALL DLACPY( 'U', NR, N, A, LDA, V, LDV )
! 1085: IF ( NR .GT. 1 )
! 1086: $ CALL DLASET( 'L', NR-1, NR-1, ZERO, ZERO, V(2,1), LDV )
! 1087: * .. the right singular vectors overwrite V, the NR left singular
! 1088: * vectors stored in U(1:NR,1:NR)
! 1089: IF ( WNTVR .OR. ( NR .EQ. N ) ) THEN
! 1090: CALL DGESVD( 'N', 'O', NR, N, V, LDV, S, U, LDU,
! 1091: $ V, LDV, WORK(N+1), LWORK-N, INFO )
! 1092: CALL DLAPMT( .FALSE., NR, N, V, LDV, IWORK )
! 1093: * .. now [V](1:NR,1:N) contains V(1:N,1:NR)**T
! 1094: ELSE
! 1095: * .. need all N right singular vectors and NR < N
! 1096: * [!] This is simple implementation that augments [V](1:NR,1:N)
! 1097: * by padding a zero block. In the case NR << N, a more efficient
! 1098: * way is to first use the LQ factorization. For more details
! 1099: * how to implement this, see the " FULL SVD " branch.
! 1100: CALL DLASET('G', N-NR, N, ZERO,ZERO, V(NR+1,1), LDV)
! 1101: CALL DGESVD( 'N', 'O', N, N, V, LDV, S, U, LDU,
! 1102: $ V, LDV, WORK(N+1), LWORK-N, INFO )
! 1103: CALL DLAPMT( .FALSE., N, N, V, LDV, IWORK )
! 1104: END IF
! 1105: * .. now [V] contains the transposed matrix of the right singular
! 1106: * vectors of A.
! 1107: END IF
! 1108: *
! 1109: ELSE
! 1110: *.......................................................................
! 1111: * .. FULL SVD requested
! 1112: *.......................................................................
! 1113: IF ( RTRANS ) THEN
! 1114: *
! 1115: * .. apply DGESVD to R**T [[this option is left for R&D&T]]
! 1116: *
! 1117: IF ( WNTVR .OR. ( NR .EQ. N ) ) THEN
! 1118: * .. copy R**T into [V] and overwrite [V] with the left singular
! 1119: * vectors of R**T
! 1120: DO 1168 p = 1, NR
! 1121: DO 1169 q = p, N
! 1122: V(q,p) = A(p,q)
! 1123: 1169 CONTINUE
! 1124: 1168 CONTINUE
! 1125: IF ( NR .GT. 1 )
! 1126: $ CALL DLASET( 'U', NR-1,NR-1, ZERO,ZERO, V(1,2), LDV )
! 1127: *
! 1128: * .. the left singular vectors of R**T overwrite [V], the NR right
! 1129: * singular vectors of R**T stored in [U](1:NR,1:NR) as transposed
! 1130: CALL DGESVD( 'O', 'A', N, NR, V, LDV, S, V, LDV,
! 1131: $ U, LDU, WORK(N+1), LWORK-N, INFO )
! 1132: * .. assemble V
! 1133: DO 1115 p = 1, NR
! 1134: DO 1116 q = p + 1, NR
! 1135: RTMP = V(q,p)
! 1136: V(q,p) = V(p,q)
! 1137: V(p,q) = RTMP
! 1138: 1116 CONTINUE
! 1139: 1115 CONTINUE
! 1140: IF ( NR .LT. N ) THEN
! 1141: DO 1101 p = 1, NR
! 1142: DO 1102 q = NR+1, N
! 1143: V(p,q) = V(q,p)
! 1144: 1102 CONTINUE
! 1145: 1101 CONTINUE
! 1146: END IF
! 1147: CALL DLAPMT( .FALSE., NR, N, V, LDV, IWORK )
! 1148: *
! 1149: DO 1117 p = 1, NR
! 1150: DO 1118 q = p + 1, NR
! 1151: RTMP = U(q,p)
! 1152: U(q,p) = U(p,q)
! 1153: U(p,q) = RTMP
! 1154: 1118 CONTINUE
! 1155: 1117 CONTINUE
! 1156: *
! 1157: IF ( ( NR .LT. M ) .AND. .NOT.(WNTUF)) THEN
! 1158: CALL DLASET('A', M-NR,NR, ZERO,ZERO, U(NR+1,1), LDU)
! 1159: IF ( NR .LT. N1 ) THEN
! 1160: CALL DLASET('A',NR,N1-NR,ZERO,ZERO,U(1,NR+1),LDU)
! 1161: CALL DLASET( 'A',M-NR,N1-NR,ZERO,ONE,
! 1162: $ U(NR+1,NR+1), LDU )
! 1163: END IF
! 1164: END IF
! 1165: *
! 1166: ELSE
! 1167: * .. need all N right singular vectors and NR < N
! 1168: * .. copy R**T into [V] and overwrite [V] with the left singular
! 1169: * vectors of R**T
! 1170: * [[The optimal ratio N/NR for using QRF instead of padding
! 1171: * with zeros. Here hard coded to 2; it must be at least
! 1172: * two due to work space constraints.]]
! 1173: * OPTRATIO = ILAENV(6, 'DGESVD', 'S' // 'O', NR,N,0,0)
! 1174: * OPTRATIO = MAX( OPTRATIO, 2 )
! 1175: OPTRATIO = 2
! 1176: IF ( OPTRATIO*NR .GT. N ) THEN
! 1177: DO 1198 p = 1, NR
! 1178: DO 1199 q = p, N
! 1179: V(q,p) = A(p,q)
! 1180: 1199 CONTINUE
! 1181: 1198 CONTINUE
! 1182: IF ( NR .GT. 1 )
! 1183: $ CALL DLASET('U',NR-1,NR-1, ZERO,ZERO, V(1,2),LDV)
! 1184: *
! 1185: CALL DLASET('A',N,N-NR,ZERO,ZERO,V(1,NR+1),LDV)
! 1186: CALL DGESVD( 'O', 'A', N, N, V, LDV, S, V, LDV,
! 1187: $ U, LDU, WORK(N+1), LWORK-N, INFO )
! 1188: *
! 1189: DO 1113 p = 1, N
! 1190: DO 1114 q = p + 1, N
! 1191: RTMP = V(q,p)
! 1192: V(q,p) = V(p,q)
! 1193: V(p,q) = RTMP
! 1194: 1114 CONTINUE
! 1195: 1113 CONTINUE
! 1196: CALL DLAPMT( .FALSE., N, N, V, LDV, IWORK )
! 1197: * .. assemble the left singular vector matrix U of dimensions
! 1198: * (M x N1), i.e. (M x N) or (M x M).
! 1199: *
! 1200: DO 1111 p = 1, N
! 1201: DO 1112 q = p + 1, N
! 1202: RTMP = U(q,p)
! 1203: U(q,p) = U(p,q)
! 1204: U(p,q) = RTMP
! 1205: 1112 CONTINUE
! 1206: 1111 CONTINUE
! 1207: *
! 1208: IF ( ( N .LT. M ) .AND. .NOT.(WNTUF)) THEN
! 1209: CALL DLASET('A',M-N,N,ZERO,ZERO,U(N+1,1),LDU)
! 1210: IF ( N .LT. N1 ) THEN
! 1211: CALL DLASET('A',N,N1-N,ZERO,ZERO,U(1,N+1),LDU)
! 1212: CALL DLASET('A',M-N,N1-N,ZERO,ONE,
! 1213: $ U(N+1,N+1), LDU )
! 1214: END IF
! 1215: END IF
! 1216: ELSE
! 1217: * .. copy R**T into [U] and overwrite [U] with the right
! 1218: * singular vectors of R
! 1219: DO 1196 p = 1, NR
! 1220: DO 1197 q = p, N
! 1221: U(q,NR+p) = A(p,q)
! 1222: 1197 CONTINUE
! 1223: 1196 CONTINUE
! 1224: IF ( NR .GT. 1 )
! 1225: $ CALL DLASET('U',NR-1,NR-1,ZERO,ZERO,U(1,NR+2),LDU)
! 1226: CALL DGEQRF( N, NR, U(1,NR+1), LDU, WORK(N+1),
! 1227: $ WORK(N+NR+1), LWORK-N-NR, IERR )
! 1228: DO 1143 p = 1, NR
! 1229: DO 1144 q = 1, N
! 1230: V(q,p) = U(p,NR+q)
! 1231: 1144 CONTINUE
! 1232: 1143 CONTINUE
! 1233: CALL DLASET('U',NR-1,NR-1,ZERO,ZERO,V(1,2),LDV)
! 1234: CALL DGESVD( 'S', 'O', NR, NR, V, LDV, S, U, LDU,
! 1235: $ V,LDV, WORK(N+NR+1),LWORK-N-NR, INFO )
! 1236: CALL DLASET('A',N-NR,NR,ZERO,ZERO,V(NR+1,1),LDV)
! 1237: CALL DLASET('A',NR,N-NR,ZERO,ZERO,V(1,NR+1),LDV)
! 1238: CALL DLASET('A',N-NR,N-NR,ZERO,ONE,V(NR+1,NR+1),LDV)
! 1239: CALL DORMQR('R','C', N, N, NR, U(1,NR+1), LDU,
! 1240: $ WORK(N+1),V,LDV,WORK(N+NR+1),LWORK-N-NR,IERR)
! 1241: CALL DLAPMT( .FALSE., N, N, V, LDV, IWORK )
! 1242: * .. assemble the left singular vector matrix U of dimensions
! 1243: * (M x NR) or (M x N) or (M x M).
! 1244: IF ( ( NR .LT. M ) .AND. .NOT.(WNTUF)) THEN
! 1245: CALL DLASET('A',M-NR,NR,ZERO,ZERO,U(NR+1,1),LDU)
! 1246: IF ( NR .LT. N1 ) THEN
! 1247: CALL DLASET('A',NR,N1-NR,ZERO,ZERO,U(1,NR+1),LDU)
! 1248: CALL DLASET( 'A',M-NR,N1-NR,ZERO,ONE,
! 1249: $ U(NR+1,NR+1),LDU)
! 1250: END IF
! 1251: END IF
! 1252: END IF
! 1253: END IF
! 1254: *
! 1255: ELSE
! 1256: *
! 1257: * .. apply DGESVD to R [[this is the recommended option]]
! 1258: *
! 1259: IF ( WNTVR .OR. ( NR .EQ. N ) ) THEN
! 1260: * .. copy R into [V] and overwrite V with the right singular vectors
! 1261: CALL DLACPY( 'U', NR, N, A, LDA, V, LDV )
! 1262: IF ( NR .GT. 1 )
! 1263: $ CALL DLASET( 'L', NR-1,NR-1, ZERO,ZERO, V(2,1), LDV )
! 1264: * .. the right singular vectors of R overwrite [V], the NR left
! 1265: * singular vectors of R stored in [U](1:NR,1:NR)
! 1266: CALL DGESVD( 'S', 'O', NR, N, V, LDV, S, U, LDU,
! 1267: $ V, LDV, WORK(N+1), LWORK-N, INFO )
! 1268: CALL DLAPMT( .FALSE., NR, N, V, LDV, IWORK )
! 1269: * .. now [V](1:NR,1:N) contains V(1:N,1:NR)**T
! 1270: * .. assemble the left singular vector matrix U of dimensions
! 1271: * (M x NR) or (M x N) or (M x M).
! 1272: IF ( ( NR .LT. M ) .AND. .NOT.(WNTUF)) THEN
! 1273: CALL DLASET('A', M-NR,NR, ZERO,ZERO, U(NR+1,1), LDU)
! 1274: IF ( NR .LT. N1 ) THEN
! 1275: CALL DLASET('A',NR,N1-NR,ZERO,ZERO,U(1,NR+1),LDU)
! 1276: CALL DLASET( 'A',M-NR,N1-NR,ZERO,ONE,
! 1277: $ U(NR+1,NR+1), LDU )
! 1278: END IF
! 1279: END IF
! 1280: *
! 1281: ELSE
! 1282: * .. need all N right singular vectors and NR < N
! 1283: * .. the requested number of the left singular vectors
! 1284: * is then N1 (N or M)
! 1285: * [[The optimal ratio N/NR for using LQ instead of padding
! 1286: * with zeros. Here hard coded to 2; it must be at least
! 1287: * two due to work space constraints.]]
! 1288: * OPTRATIO = ILAENV(6, 'DGESVD', 'S' // 'O', NR,N,0,0)
! 1289: * OPTRATIO = MAX( OPTRATIO, 2 )
! 1290: OPTRATIO = 2
! 1291: IF ( OPTRATIO * NR .GT. N ) THEN
! 1292: CALL DLACPY( 'U', NR, N, A, LDA, V, LDV )
! 1293: IF ( NR .GT. 1 )
! 1294: $ CALL DLASET('L', NR-1,NR-1, ZERO,ZERO, V(2,1),LDV)
! 1295: * .. the right singular vectors of R overwrite [V], the NR left
! 1296: * singular vectors of R stored in [U](1:NR,1:NR)
! 1297: CALL DLASET('A', N-NR,N, ZERO,ZERO, V(NR+1,1),LDV)
! 1298: CALL DGESVD( 'S', 'O', N, N, V, LDV, S, U, LDU,
! 1299: $ V, LDV, WORK(N+1), LWORK-N, INFO )
! 1300: CALL DLAPMT( .FALSE., N, N, V, LDV, IWORK )
! 1301: * .. now [V] contains the transposed matrix of the right
! 1302: * singular vectors of A. The leading N left singular vectors
! 1303: * are in [U](1:N,1:N)
! 1304: * .. assemble the left singular vector matrix U of dimensions
! 1305: * (M x N1), i.e. (M x N) or (M x M).
! 1306: IF ( ( N .LT. M ) .AND. .NOT.(WNTUF)) THEN
! 1307: CALL DLASET('A',M-N,N,ZERO,ZERO,U(N+1,1),LDU)
! 1308: IF ( N .LT. N1 ) THEN
! 1309: CALL DLASET('A',N,N1-N,ZERO,ZERO,U(1,N+1),LDU)
! 1310: CALL DLASET( 'A',M-N,N1-N,ZERO,ONE,
! 1311: $ U(N+1,N+1), LDU )
! 1312: END IF
! 1313: END IF
! 1314: ELSE
! 1315: CALL DLACPY( 'U', NR, N, A, LDA, U(NR+1,1), LDU )
! 1316: IF ( NR .GT. 1 )
! 1317: $ CALL DLASET('L',NR-1,NR-1,ZERO,ZERO,U(NR+2,1),LDU)
! 1318: CALL DGELQF( NR, N, U(NR+1,1), LDU, WORK(N+1),
! 1319: $ WORK(N+NR+1), LWORK-N-NR, IERR )
! 1320: CALL DLACPY('L',NR,NR,U(NR+1,1),LDU,V,LDV)
! 1321: IF ( NR .GT. 1 )
! 1322: $ CALL DLASET('U',NR-1,NR-1,ZERO,ZERO,V(1,2),LDV)
! 1323: CALL DGESVD( 'S', 'O', NR, NR, V, LDV, S, U, LDU,
! 1324: $ V, LDV, WORK(N+NR+1), LWORK-N-NR, INFO )
! 1325: CALL DLASET('A',N-NR,NR,ZERO,ZERO,V(NR+1,1),LDV)
! 1326: CALL DLASET('A',NR,N-NR,ZERO,ZERO,V(1,NR+1),LDV)
! 1327: CALL DLASET('A',N-NR,N-NR,ZERO,ONE,V(NR+1,NR+1),LDV)
! 1328: CALL DORMLQ('R','N',N,N,NR,U(NR+1,1),LDU,WORK(N+1),
! 1329: $ V, LDV, WORK(N+NR+1),LWORK-N-NR,IERR)
! 1330: CALL DLAPMT( .FALSE., N, N, V, LDV, IWORK )
! 1331: * .. assemble the left singular vector matrix U of dimensions
! 1332: * (M x NR) or (M x N) or (M x M).
! 1333: IF ( ( NR .LT. M ) .AND. .NOT.(WNTUF)) THEN
! 1334: CALL DLASET('A',M-NR,NR,ZERO,ZERO,U(NR+1,1),LDU)
! 1335: IF ( NR .LT. N1 ) THEN
! 1336: CALL DLASET('A',NR,N1-NR,ZERO,ZERO,U(1,NR+1),LDU)
! 1337: CALL DLASET( 'A',M-NR,N1-NR,ZERO,ONE,
! 1338: $ U(NR+1,NR+1), LDU )
! 1339: END IF
! 1340: END IF
! 1341: END IF
! 1342: END IF
! 1343: * .. end of the "R**T or R" branch
! 1344: END IF
! 1345: *
! 1346: * The Q matrix from the first QRF is built into the left singular
! 1347: * vectors matrix U.
! 1348: *
! 1349: IF ( .NOT. WNTUF )
! 1350: $ CALL DORMQR( 'L', 'N', M, N1, N, A, LDA, WORK, U,
! 1351: $ LDU, WORK(N+1), LWORK-N, IERR )
! 1352: IF ( ROWPRM .AND. .NOT.WNTUF )
! 1353: $ CALL DLASWP( N1, U, LDU, 1, M-1, IWORK(N+1), -1 )
! 1354: *
! 1355: * ... end of the "full SVD" branch
! 1356: END IF
! 1357: *
! 1358: * Check whether some singular values are returned as zeros, e.g.
! 1359: * due to underflow, and update the numerical rank.
! 1360: p = NR
! 1361: DO 4001 q = p, 1, -1
! 1362: IF ( S(q) .GT. ZERO ) GO TO 4002
! 1363: NR = NR - 1
! 1364: 4001 CONTINUE
! 1365: 4002 CONTINUE
! 1366: *
! 1367: * .. if numerical rank deficiency is detected, the truncated
! 1368: * singular values are set to zero.
! 1369: IF ( NR .LT. N ) CALL DLASET( 'G', N-NR,1, ZERO,ZERO, S(NR+1), N )
! 1370: * .. undo scaling; this may cause overflow in the largest singular
! 1371: * values.
! 1372: IF ( ASCALED )
! 1373: $ CALL DLASCL( 'G',0,0, ONE,SQRT(DBLE(M)), NR,1, S, N, IERR )
! 1374: IF ( CONDA ) RWORK(1) = SCONDA
! 1375: RWORK(2) = p - NR
! 1376: * .. p-NR is the number of singular values that are computed as
! 1377: * exact zeros in DGESVD() applied to the (possibly truncated)
! 1378: * full row rank triangular (trapezoidal) factor of A.
! 1379: NUMRANK = NR
! 1380: *
! 1381: RETURN
! 1382: *
! 1383: * End of DGESVDQ
! 1384: *
! 1385: END
CVSweb interface <joel.bertrand@systella.fr>