Annotation of rpl/lapack/lapack/dgelss.f, revision 1.9
1.9 ! bertrand 1: *> \brief <b> DGELSS solves overdetermined or underdetermined systems 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 DGELSS + dependencies
! 10: *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dgelss.f">
! 11: *> [TGZ]</a>
! 12: *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dgelss.f">
! 13: *> [ZIP]</a>
! 14: *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dgelss.f">
! 15: *> [TXT]</a>
! 16: *> \endhtmlonly
! 17: *
! 18: * Definition:
! 19: * ===========
! 20: *
! 21: * SUBROUTINE DGELSS( M, N, NRHS, A, LDA, B, LDB, S, RCOND, RANK,
! 22: * WORK, LWORK, INFO )
! 23: *
! 24: * .. Scalar Arguments ..
! 25: * INTEGER INFO, LDA, LDB, LWORK, M, N, NRHS, RANK
! 26: * DOUBLE PRECISION RCOND
! 27: * ..
! 28: * .. Array Arguments ..
! 29: * DOUBLE PRECISION A( LDA, * ), B( LDB, * ), S( * ), WORK( * )
! 30: * ..
! 31: *
! 32: *
! 33: *> \par Purpose:
! 34: * =============
! 35: *>
! 36: *> \verbatim
! 37: *>
! 38: *> DGELSS computes the minimum norm solution to a real linear least
! 39: *> squares problem:
! 40: *>
! 41: *> Minimize 2-norm(| b - A*x |).
! 42: *>
! 43: *> using the singular value decomposition (SVD) of A. A is an M-by-N
! 44: *> matrix which may be rank-deficient.
! 45: *>
! 46: *> Several right hand side vectors b and solution vectors x can be
! 47: *> handled in a single call; they are stored as the columns of the
! 48: *> M-by-NRHS right hand side matrix B and the N-by-NRHS solution matrix
! 49: *> X.
! 50: *>
! 51: *> The effective rank of A is determined by treating as zero those
! 52: *> singular values which are less than RCOND times the largest singular
! 53: *> value.
! 54: *> \endverbatim
! 55: *
! 56: * Arguments:
! 57: * ==========
! 58: *
! 59: *> \param[in] M
! 60: *> \verbatim
! 61: *> M is INTEGER
! 62: *> The number of rows of the matrix A. M >= 0.
! 63: *> \endverbatim
! 64: *>
! 65: *> \param[in] N
! 66: *> \verbatim
! 67: *> N is INTEGER
! 68: *> The number of columns of the matrix A. N >= 0.
! 69: *> \endverbatim
! 70: *>
! 71: *> \param[in] NRHS
! 72: *> \verbatim
! 73: *> NRHS is INTEGER
! 74: *> The number of right hand sides, i.e., the number of columns
! 75: *> of the matrices B and X. NRHS >= 0.
! 76: *> \endverbatim
! 77: *>
! 78: *> \param[in,out] A
! 79: *> \verbatim
! 80: *> A is DOUBLE PRECISION array, dimension (LDA,N)
! 81: *> On entry, the M-by-N matrix A.
! 82: *> On exit, the first min(m,n) rows of A are overwritten with
! 83: *> its right singular vectors, stored rowwise.
! 84: *> \endverbatim
! 85: *>
! 86: *> \param[in] LDA
! 87: *> \verbatim
! 88: *> LDA is INTEGER
! 89: *> The leading dimension of the array A. LDA >= max(1,M).
! 90: *> \endverbatim
! 91: *>
! 92: *> \param[in,out] B
! 93: *> \verbatim
! 94: *> B is DOUBLE PRECISION array, dimension (LDB,NRHS)
! 95: *> On entry, the M-by-NRHS right hand side matrix B.
! 96: *> On exit, B is overwritten by the N-by-NRHS solution
! 97: *> matrix X. If m >= n and RANK = n, the residual
! 98: *> sum-of-squares for the solution in the i-th column is given
! 99: *> by the sum of squares of elements n+1:m in that column.
! 100: *> \endverbatim
! 101: *>
! 102: *> \param[in] LDB
! 103: *> \verbatim
! 104: *> LDB is INTEGER
! 105: *> The leading dimension of the array B. LDB >= max(1,max(M,N)).
! 106: *> \endverbatim
! 107: *>
! 108: *> \param[out] S
! 109: *> \verbatim
! 110: *> S is DOUBLE PRECISION array, dimension (min(M,N))
! 111: *> The singular values of A in decreasing order.
! 112: *> The condition number of A in the 2-norm = S(1)/S(min(m,n)).
! 113: *> \endverbatim
! 114: *>
! 115: *> \param[in] RCOND
! 116: *> \verbatim
! 117: *> RCOND is DOUBLE PRECISION
! 118: *> RCOND is used to determine the effective rank of A.
! 119: *> Singular values S(i) <= RCOND*S(1) are treated as zero.
! 120: *> If RCOND < 0, machine precision is used instead.
! 121: *> \endverbatim
! 122: *>
! 123: *> \param[out] RANK
! 124: *> \verbatim
! 125: *> RANK is INTEGER
! 126: *> The effective rank of A, i.e., the number of singular values
! 127: *> which are greater than RCOND*S(1).
! 128: *> \endverbatim
! 129: *>
! 130: *> \param[out] WORK
! 131: *> \verbatim
! 132: *> WORK is DOUBLE PRECISION array, dimension (MAX(1,LWORK))
! 133: *> On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
! 134: *> \endverbatim
! 135: *>
! 136: *> \param[in] LWORK
! 137: *> \verbatim
! 138: *> LWORK is INTEGER
! 139: *> The dimension of the array WORK. LWORK >= 1, and also:
! 140: *> LWORK >= 3*min(M,N) + max( 2*min(M,N), max(M,N), NRHS )
! 141: *> For good performance, LWORK should generally be larger.
! 142: *>
! 143: *> If LWORK = -1, then a workspace query is assumed; the routine
! 144: *> only calculates the optimal size of the WORK array, returns
! 145: *> this value as the first entry of the WORK array, and no error
! 146: *> message related to LWORK is issued by XERBLA.
! 147: *> \endverbatim
! 148: *>
! 149: *> \param[out] INFO
! 150: *> \verbatim
! 151: *> INFO is INTEGER
! 152: *> = 0: successful exit
! 153: *> < 0: if INFO = -i, the i-th argument had an illegal value.
! 154: *> > 0: the algorithm for computing the SVD failed to converge;
! 155: *> if INFO = i, i off-diagonal elements of an intermediate
! 156: *> bidiagonal form did not converge to zero.
! 157: *> \endverbatim
! 158: *
! 159: * Authors:
! 160: * ========
! 161: *
! 162: *> \author Univ. of Tennessee
! 163: *> \author Univ. of California Berkeley
! 164: *> \author Univ. of Colorado Denver
! 165: *> \author NAG Ltd.
! 166: *
! 167: *> \date November 2011
! 168: *
! 169: *> \ingroup doubleGEsolve
! 170: *
! 171: * =====================================================================
1.1 bertrand 172: SUBROUTINE DGELSS( M, N, NRHS, A, LDA, B, LDB, S, RCOND, RANK,
173: $ WORK, LWORK, INFO )
174: *
1.9 ! bertrand 175: * -- LAPACK driver routine (version 3.4.0) --
1.1 bertrand 176: * -- LAPACK is a software package provided by Univ. of Tennessee, --
177: * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
1.9 ! bertrand 178: * November 2011
1.1 bertrand 179: *
180: * .. Scalar Arguments ..
181: INTEGER INFO, LDA, LDB, LWORK, M, N, NRHS, RANK
182: DOUBLE PRECISION RCOND
183: * ..
184: * .. Array Arguments ..
185: DOUBLE PRECISION A( LDA, * ), B( LDB, * ), S( * ), WORK( * )
186: * ..
187: *
188: * =====================================================================
189: *
190: * .. Parameters ..
191: DOUBLE PRECISION ZERO, ONE
192: PARAMETER ( ZERO = 0.0D+0, ONE = 1.0D+0 )
193: * ..
194: * .. Local Scalars ..
195: LOGICAL LQUERY
196: INTEGER BDSPAC, BL, CHUNK, I, IASCL, IBSCL, IE, IL,
197: $ ITAU, ITAUP, ITAUQ, IWORK, LDWORK, MAXMN,
198: $ MAXWRK, MINMN, MINWRK, MM, MNTHR
1.9 ! bertrand 199: INTEGER LWORK_DGEQRF, LWORK_DORMQR, LWORK_DGEBRD,
! 200: $ LWORK_DORMBR, LWORK_DORGBR, LWORK_DORMLQ,
! 201: $ LWORK_DGELQF
1.1 bertrand 202: DOUBLE PRECISION ANRM, BIGNUM, BNRM, EPS, SFMIN, SMLNUM, THR
203: * ..
204: * .. Local Arrays ..
1.9 ! bertrand 205: DOUBLE PRECISION DUM( 1 )
1.1 bertrand 206: * ..
207: * .. External Subroutines ..
208: EXTERNAL DBDSQR, DCOPY, DGEBRD, DGELQF, DGEMM, DGEMV,
209: $ DGEQRF, DLABAD, DLACPY, DLASCL, DLASET, DORGBR,
210: $ DORMBR, DORMLQ, DORMQR, DRSCL, XERBLA
211: * ..
212: * .. External Functions ..
213: INTEGER ILAENV
214: DOUBLE PRECISION DLAMCH, DLANGE
215: EXTERNAL ILAENV, DLAMCH, DLANGE
216: * ..
217: * .. Intrinsic Functions ..
218: INTRINSIC MAX, MIN
219: * ..
220: * .. Executable Statements ..
221: *
222: * Test the input arguments
223: *
224: INFO = 0
225: MINMN = MIN( M, N )
226: MAXMN = MAX( M, N )
227: LQUERY = ( LWORK.EQ.-1 )
228: IF( M.LT.0 ) THEN
229: INFO = -1
230: ELSE IF( N.LT.0 ) THEN
231: INFO = -2
232: ELSE IF( NRHS.LT.0 ) THEN
233: INFO = -3
234: ELSE IF( LDA.LT.MAX( 1, M ) ) THEN
235: INFO = -5
236: ELSE IF( LDB.LT.MAX( 1, MAXMN ) ) THEN
237: INFO = -7
238: END IF
239: *
240: * Compute workspace
241: * (Note: Comments in the code beginning "Workspace:" describe the
242: * minimal amount of workspace needed at that point in the code,
243: * as well as the preferred amount for good performance.
244: * NB refers to the optimal block size for the immediately
245: * following subroutine, as returned by ILAENV.)
246: *
247: IF( INFO.EQ.0 ) THEN
248: MINWRK = 1
249: MAXWRK = 1
250: IF( MINMN.GT.0 ) THEN
251: MM = M
252: MNTHR = ILAENV( 6, 'DGELSS', ' ', M, N, NRHS, -1 )
253: IF( M.GE.N .AND. M.GE.MNTHR ) THEN
254: *
255: * Path 1a - overdetermined, with many more rows than
256: * columns
257: *
1.9 ! bertrand 258: * Compute space needed for DGEQRF
! 259: CALL DGEQRF( M, N, A, LDA, DUM(1), DUM(1), -1, INFO )
! 260: LWORK_DGEQRF=DUM(1)
! 261: * Compute space needed for DORMQR
! 262: CALL DORMQR( 'L', 'T', M, NRHS, N, A, LDA, DUM(1), B,
! 263: $ LDB, DUM(1), -1, INFO )
! 264: LWORK_DORMQR=DUM(1)
1.1 bertrand 265: MM = N
1.9 ! bertrand 266: MAXWRK = MAX( MAXWRK, N + LWORK_DGEQRF )
! 267: MAXWRK = MAX( MAXWRK, N + LWORK_DORMQR )
1.1 bertrand 268: END IF
269: IF( M.GE.N ) THEN
270: *
271: * Path 1 - overdetermined or exactly determined
272: *
273: * Compute workspace needed for DBDSQR
274: *
275: BDSPAC = MAX( 1, 5*N )
1.9 ! bertrand 276: * Compute space needed for DGEBRD
! 277: CALL DGEBRD( MM, N, A, LDA, S, DUM(1), DUM(1),
! 278: $ DUM(1), DUM(1), -1, INFO )
! 279: LWORK_DGEBRD=DUM(1)
! 280: * Compute space needed for DORMBR
! 281: CALL DORMBR( 'Q', 'L', 'T', MM, NRHS, N, A, LDA, DUM(1),
! 282: $ B, LDB, DUM(1), -1, INFO )
! 283: LWORK_DORMBR=DUM(1)
! 284: * Compute space needed for DORGBR
! 285: CALL DORGBR( 'P', N, N, N, A, LDA, DUM(1),
! 286: $ DUM(1), -1, INFO )
! 287: LWORK_DORGBR=DUM(1)
! 288: * Compute total workspace needed
! 289: MAXWRK = MAX( MAXWRK, 3*N + LWORK_DGEBRD )
! 290: MAXWRK = MAX( MAXWRK, 3*N + LWORK_DORMBR )
! 291: MAXWRK = MAX( MAXWRK, 3*N + LWORK_DORGBR )
1.1 bertrand 292: MAXWRK = MAX( MAXWRK, BDSPAC )
293: MAXWRK = MAX( MAXWRK, N*NRHS )
294: MINWRK = MAX( 3*N + MM, 3*N + NRHS, BDSPAC )
295: MAXWRK = MAX( MINWRK, MAXWRK )
296: END IF
297: IF( N.GT.M ) THEN
298: *
299: * Compute workspace needed for DBDSQR
300: *
301: BDSPAC = MAX( 1, 5*M )
302: MINWRK = MAX( 3*M+NRHS, 3*M+N, BDSPAC )
303: IF( N.GE.MNTHR ) THEN
304: *
305: * Path 2a - underdetermined, with many more columns
306: * than rows
307: *
1.9 ! bertrand 308: * Compute space needed for DGELQF
! 309: CALL DGELQF( M, N, A, LDA, DUM(1), DUM(1),
! 310: $ -1, INFO )
! 311: LWORK_DGELQF=DUM(1)
! 312: * Compute space needed for DGEBRD
! 313: CALL DGEBRD( M, M, A, LDA, S, DUM(1), DUM(1),
! 314: $ DUM(1), DUM(1), -1, INFO )
! 315: LWORK_DGEBRD=DUM(1)
! 316: * Compute space needed for DORMBR
! 317: CALL DORMBR( 'Q', 'L', 'T', M, NRHS, N, A, LDA,
! 318: $ DUM(1), B, LDB, DUM(1), -1, INFO )
! 319: LWORK_DORMBR=DUM(1)
! 320: * Compute space needed for DORGBR
! 321: CALL DORGBR( 'P', M, M, M, A, LDA, DUM(1),
! 322: $ DUM(1), -1, INFO )
! 323: LWORK_DORGBR=DUM(1)
! 324: * Compute space needed for DORMLQ
! 325: CALL DORMLQ( 'L', 'T', N, NRHS, M, A, LDA, DUM(1),
! 326: $ B, LDB, DUM(1), -1, INFO )
! 327: LWORK_DORMLQ=DUM(1)
! 328: * Compute total workspace needed
! 329: MAXWRK = M + LWORK_DGELQF
! 330: MAXWRK = MAX( MAXWRK, M*M + 4*M + LWORK_DGEBRD )
! 331: MAXWRK = MAX( MAXWRK, M*M + 4*M + LWORK_DORMBR )
! 332: MAXWRK = MAX( MAXWRK, M*M + 4*M + LWORK_DORGBR )
1.1 bertrand 333: MAXWRK = MAX( MAXWRK, M*M + M + BDSPAC )
334: IF( NRHS.GT.1 ) THEN
335: MAXWRK = MAX( MAXWRK, M*M + M + M*NRHS )
336: ELSE
337: MAXWRK = MAX( MAXWRK, M*M + 2*M )
338: END IF
1.9 ! bertrand 339: MAXWRK = MAX( MAXWRK, M + LWORK_DORMLQ )
1.1 bertrand 340: ELSE
341: *
342: * Path 2 - underdetermined
343: *
1.9 ! bertrand 344: * Compute space needed for DGEBRD
! 345: CALL DGEBRD( M, N, A, LDA, S, DUM(1), DUM(1),
! 346: $ DUM(1), DUM(1), -1, INFO )
! 347: LWORK_DGEBRD=DUM(1)
! 348: * Compute space needed for DORMBR
! 349: CALL DORMBR( 'Q', 'L', 'T', M, NRHS, M, A, LDA,
! 350: $ DUM(1), B, LDB, DUM(1), -1, INFO )
! 351: LWORK_DORMBR=DUM(1)
! 352: * Compute space needed for DORGBR
! 353: CALL DORGBR( 'P', M, N, M, A, LDA, DUM(1),
! 354: $ DUM(1), -1, INFO )
! 355: LWORK_DORGBR=DUM(1)
! 356: MAXWRK = 3*M + LWORK_DGEBRD
! 357: MAXWRK = MAX( MAXWRK, 3*M + LWORK_DORMBR )
! 358: MAXWRK = MAX( MAXWRK, 3*M + LWORK_DORGBR )
1.1 bertrand 359: MAXWRK = MAX( MAXWRK, BDSPAC )
360: MAXWRK = MAX( MAXWRK, N*NRHS )
361: END IF
362: END IF
363: MAXWRK = MAX( MINWRK, MAXWRK )
364: END IF
365: WORK( 1 ) = MAXWRK
366: *
367: IF( LWORK.LT.MINWRK .AND. .NOT.LQUERY )
368: $ INFO = -12
369: END IF
370: *
371: IF( INFO.NE.0 ) THEN
372: CALL XERBLA( 'DGELSS', -INFO )
373: RETURN
374: ELSE IF( LQUERY ) THEN
375: RETURN
376: END IF
377: *
378: * Quick return if possible
379: *
380: IF( M.EQ.0 .OR. N.EQ.0 ) THEN
381: RANK = 0
382: RETURN
383: END IF
384: *
385: * Get machine parameters
386: *
387: EPS = DLAMCH( 'P' )
388: SFMIN = DLAMCH( 'S' )
389: SMLNUM = SFMIN / EPS
390: BIGNUM = ONE / SMLNUM
391: CALL DLABAD( SMLNUM, BIGNUM )
392: *
393: * Scale A if max element outside range [SMLNUM,BIGNUM]
394: *
395: ANRM = DLANGE( 'M', M, N, A, LDA, WORK )
396: IASCL = 0
397: IF( ANRM.GT.ZERO .AND. ANRM.LT.SMLNUM ) THEN
398: *
399: * Scale matrix norm up to SMLNUM
400: *
401: CALL DLASCL( 'G', 0, 0, ANRM, SMLNUM, M, N, A, LDA, INFO )
402: IASCL = 1
403: ELSE IF( ANRM.GT.BIGNUM ) THEN
404: *
405: * Scale matrix norm down to BIGNUM
406: *
407: CALL DLASCL( 'G', 0, 0, ANRM, BIGNUM, M, N, A, LDA, INFO )
408: IASCL = 2
409: ELSE IF( ANRM.EQ.ZERO ) THEN
410: *
411: * Matrix all zero. Return zero solution.
412: *
413: CALL DLASET( 'F', MAX( M, N ), NRHS, ZERO, ZERO, B, LDB )
1.5 bertrand 414: CALL DLASET( 'F', MINMN, 1, ZERO, ZERO, S, MINMN )
1.1 bertrand 415: RANK = 0
416: GO TO 70
417: END IF
418: *
419: * Scale B if max element outside range [SMLNUM,BIGNUM]
420: *
421: BNRM = DLANGE( 'M', M, NRHS, B, LDB, WORK )
422: IBSCL = 0
423: IF( BNRM.GT.ZERO .AND. BNRM.LT.SMLNUM ) THEN
424: *
425: * Scale matrix norm up to SMLNUM
426: *
427: CALL DLASCL( 'G', 0, 0, BNRM, SMLNUM, M, NRHS, B, LDB, INFO )
428: IBSCL = 1
429: ELSE IF( BNRM.GT.BIGNUM ) THEN
430: *
431: * Scale matrix norm down to BIGNUM
432: *
433: CALL DLASCL( 'G', 0, 0, BNRM, BIGNUM, M, NRHS, B, LDB, INFO )
434: IBSCL = 2
435: END IF
436: *
437: * Overdetermined case
438: *
439: IF( M.GE.N ) THEN
440: *
441: * Path 1 - overdetermined or exactly determined
442: *
443: MM = M
444: IF( M.GE.MNTHR ) THEN
445: *
446: * Path 1a - overdetermined, with many more rows than columns
447: *
448: MM = N
449: ITAU = 1
450: IWORK = ITAU + N
451: *
452: * Compute A=Q*R
453: * (Workspace: need 2*N, prefer N+N*NB)
454: *
455: CALL DGEQRF( M, N, A, LDA, WORK( ITAU ), WORK( IWORK ),
456: $ LWORK-IWORK+1, INFO )
457: *
458: * Multiply B by transpose(Q)
459: * (Workspace: need N+NRHS, prefer N+NRHS*NB)
460: *
461: CALL DORMQR( 'L', 'T', M, NRHS, N, A, LDA, WORK( ITAU ), B,
462: $ LDB, WORK( IWORK ), LWORK-IWORK+1, INFO )
463: *
464: * Zero out below R
465: *
466: IF( N.GT.1 )
467: $ CALL DLASET( 'L', N-1, N-1, ZERO, ZERO, A( 2, 1 ), LDA )
468: END IF
469: *
470: IE = 1
471: ITAUQ = IE + N
472: ITAUP = ITAUQ + N
473: IWORK = ITAUP + N
474: *
475: * Bidiagonalize R in A
476: * (Workspace: need 3*N+MM, prefer 3*N+(MM+N)*NB)
477: *
478: CALL DGEBRD( MM, N, A, LDA, S, WORK( IE ), WORK( ITAUQ ),
479: $ WORK( ITAUP ), WORK( IWORK ), LWORK-IWORK+1,
480: $ INFO )
481: *
482: * Multiply B by transpose of left bidiagonalizing vectors of R
483: * (Workspace: need 3*N+NRHS, prefer 3*N+NRHS*NB)
484: *
485: CALL DORMBR( 'Q', 'L', 'T', MM, NRHS, N, A, LDA, WORK( ITAUQ ),
486: $ B, LDB, WORK( IWORK ), LWORK-IWORK+1, INFO )
487: *
488: * Generate right bidiagonalizing vectors of R in A
489: * (Workspace: need 4*N-1, prefer 3*N+(N-1)*NB)
490: *
491: CALL DORGBR( 'P', N, N, N, A, LDA, WORK( ITAUP ),
492: $ WORK( IWORK ), LWORK-IWORK+1, INFO )
493: IWORK = IE + N
494: *
495: * Perform bidiagonal QR iteration
496: * multiply B by transpose of left singular vectors
497: * compute right singular vectors in A
498: * (Workspace: need BDSPAC)
499: *
1.9 ! bertrand 500: CALL DBDSQR( 'U', N, N, 0, NRHS, S, WORK( IE ), A, LDA, DUM,
1.1 bertrand 501: $ 1, B, LDB, WORK( IWORK ), INFO )
502: IF( INFO.NE.0 )
503: $ GO TO 70
504: *
505: * Multiply B by reciprocals of singular values
506: *
507: THR = MAX( RCOND*S( 1 ), SFMIN )
508: IF( RCOND.LT.ZERO )
509: $ THR = MAX( EPS*S( 1 ), SFMIN )
510: RANK = 0
511: DO 10 I = 1, N
512: IF( S( I ).GT.THR ) THEN
513: CALL DRSCL( NRHS, S( I ), B( I, 1 ), LDB )
514: RANK = RANK + 1
515: ELSE
516: CALL DLASET( 'F', 1, NRHS, ZERO, ZERO, B( I, 1 ), LDB )
517: END IF
518: 10 CONTINUE
519: *
520: * Multiply B by right singular vectors
521: * (Workspace: need N, prefer N*NRHS)
522: *
523: IF( LWORK.GE.LDB*NRHS .AND. NRHS.GT.1 ) THEN
524: CALL DGEMM( 'T', 'N', N, NRHS, N, ONE, A, LDA, B, LDB, ZERO,
525: $ WORK, LDB )
526: CALL DLACPY( 'G', N, NRHS, WORK, LDB, B, LDB )
527: ELSE IF( NRHS.GT.1 ) THEN
528: CHUNK = LWORK / N
529: DO 20 I = 1, NRHS, CHUNK
530: BL = MIN( NRHS-I+1, CHUNK )
531: CALL DGEMM( 'T', 'N', N, BL, N, ONE, A, LDA, B( 1, I ),
532: $ LDB, ZERO, WORK, N )
533: CALL DLACPY( 'G', N, BL, WORK, N, B( 1, I ), LDB )
534: 20 CONTINUE
535: ELSE
536: CALL DGEMV( 'T', N, N, ONE, A, LDA, B, 1, ZERO, WORK, 1 )
537: CALL DCOPY( N, WORK, 1, B, 1 )
538: END IF
539: *
540: ELSE IF( N.GE.MNTHR .AND. LWORK.GE.4*M+M*M+
541: $ MAX( M, 2*M-4, NRHS, N-3*M ) ) THEN
542: *
543: * Path 2a - underdetermined, with many more columns than rows
544: * and sufficient workspace for an efficient algorithm
545: *
546: LDWORK = M
547: IF( LWORK.GE.MAX( 4*M+M*LDA+MAX( M, 2*M-4, NRHS, N-3*M ),
548: $ M*LDA+M+M*NRHS ) )LDWORK = LDA
549: ITAU = 1
550: IWORK = M + 1
551: *
552: * Compute A=L*Q
553: * (Workspace: need 2*M, prefer M+M*NB)
554: *
555: CALL DGELQF( M, N, A, LDA, WORK( ITAU ), WORK( IWORK ),
556: $ LWORK-IWORK+1, INFO )
557: IL = IWORK
558: *
559: * Copy L to WORK(IL), zeroing out above it
560: *
561: CALL DLACPY( 'L', M, M, A, LDA, WORK( IL ), LDWORK )
562: CALL DLASET( 'U', M-1, M-1, ZERO, ZERO, WORK( IL+LDWORK ),
563: $ LDWORK )
564: IE = IL + LDWORK*M
565: ITAUQ = IE + M
566: ITAUP = ITAUQ + M
567: IWORK = ITAUP + M
568: *
569: * Bidiagonalize L in WORK(IL)
570: * (Workspace: need M*M+5*M, prefer M*M+4*M+2*M*NB)
571: *
572: CALL DGEBRD( M, M, WORK( IL ), LDWORK, S, WORK( IE ),
573: $ WORK( ITAUQ ), WORK( ITAUP ), WORK( IWORK ),
574: $ LWORK-IWORK+1, INFO )
575: *
576: * Multiply B by transpose of left bidiagonalizing vectors of L
577: * (Workspace: need M*M+4*M+NRHS, prefer M*M+4*M+NRHS*NB)
578: *
579: CALL DORMBR( 'Q', 'L', 'T', M, NRHS, M, WORK( IL ), LDWORK,
580: $ WORK( ITAUQ ), B, LDB, WORK( IWORK ),
581: $ LWORK-IWORK+1, INFO )
582: *
583: * Generate right bidiagonalizing vectors of R in WORK(IL)
584: * (Workspace: need M*M+5*M-1, prefer M*M+4*M+(M-1)*NB)
585: *
586: CALL DORGBR( 'P', M, M, M, WORK( IL ), LDWORK, WORK( ITAUP ),
587: $ WORK( IWORK ), LWORK-IWORK+1, INFO )
588: IWORK = IE + M
589: *
590: * Perform bidiagonal QR iteration,
591: * computing right singular vectors of L in WORK(IL) and
592: * multiplying B by transpose of left singular vectors
593: * (Workspace: need M*M+M+BDSPAC)
594: *
595: CALL DBDSQR( 'U', M, M, 0, NRHS, S, WORK( IE ), WORK( IL ),
596: $ LDWORK, A, LDA, B, LDB, WORK( IWORK ), INFO )
597: IF( INFO.NE.0 )
598: $ GO TO 70
599: *
600: * Multiply B by reciprocals of singular values
601: *
602: THR = MAX( RCOND*S( 1 ), SFMIN )
603: IF( RCOND.LT.ZERO )
604: $ THR = MAX( EPS*S( 1 ), SFMIN )
605: RANK = 0
606: DO 30 I = 1, M
607: IF( S( I ).GT.THR ) THEN
608: CALL DRSCL( NRHS, S( I ), B( I, 1 ), LDB )
609: RANK = RANK + 1
610: ELSE
611: CALL DLASET( 'F', 1, NRHS, ZERO, ZERO, B( I, 1 ), LDB )
612: END IF
613: 30 CONTINUE
614: IWORK = IE
615: *
616: * Multiply B by right singular vectors of L in WORK(IL)
617: * (Workspace: need M*M+2*M, prefer M*M+M+M*NRHS)
618: *
619: IF( LWORK.GE.LDB*NRHS+IWORK-1 .AND. NRHS.GT.1 ) THEN
620: CALL DGEMM( 'T', 'N', M, NRHS, M, ONE, WORK( IL ), LDWORK,
621: $ B, LDB, ZERO, WORK( IWORK ), LDB )
622: CALL DLACPY( 'G', M, NRHS, WORK( IWORK ), LDB, B, LDB )
623: ELSE IF( NRHS.GT.1 ) THEN
624: CHUNK = ( LWORK-IWORK+1 ) / M
625: DO 40 I = 1, NRHS, CHUNK
626: BL = MIN( NRHS-I+1, CHUNK )
627: CALL DGEMM( 'T', 'N', M, BL, M, ONE, WORK( IL ), LDWORK,
628: $ B( 1, I ), LDB, ZERO, WORK( IWORK ), M )
629: CALL DLACPY( 'G', M, BL, WORK( IWORK ), M, B( 1, I ),
630: $ LDB )
631: 40 CONTINUE
632: ELSE
633: CALL DGEMV( 'T', M, M, ONE, WORK( IL ), LDWORK, B( 1, 1 ),
634: $ 1, ZERO, WORK( IWORK ), 1 )
635: CALL DCOPY( M, WORK( IWORK ), 1, B( 1, 1 ), 1 )
636: END IF
637: *
638: * Zero out below first M rows of B
639: *
640: CALL DLASET( 'F', N-M, NRHS, ZERO, ZERO, B( M+1, 1 ), LDB )
641: IWORK = ITAU + M
642: *
643: * Multiply transpose(Q) by B
644: * (Workspace: need M+NRHS, prefer M+NRHS*NB)
645: *
646: CALL DORMLQ( 'L', 'T', N, NRHS, M, A, LDA, WORK( ITAU ), B,
647: $ LDB, WORK( IWORK ), LWORK-IWORK+1, INFO )
648: *
649: ELSE
650: *
651: * Path 2 - remaining underdetermined cases
652: *
653: IE = 1
654: ITAUQ = IE + M
655: ITAUP = ITAUQ + M
656: IWORK = ITAUP + M
657: *
658: * Bidiagonalize A
659: * (Workspace: need 3*M+N, prefer 3*M+(M+N)*NB)
660: *
661: CALL DGEBRD( M, N, A, LDA, S, WORK( IE ), WORK( ITAUQ ),
662: $ WORK( ITAUP ), WORK( IWORK ), LWORK-IWORK+1,
663: $ INFO )
664: *
665: * Multiply B by transpose of left bidiagonalizing vectors
666: * (Workspace: need 3*M+NRHS, prefer 3*M+NRHS*NB)
667: *
668: CALL DORMBR( 'Q', 'L', 'T', M, NRHS, N, A, LDA, WORK( ITAUQ ),
669: $ B, LDB, WORK( IWORK ), LWORK-IWORK+1, INFO )
670: *
671: * Generate right bidiagonalizing vectors in A
672: * (Workspace: need 4*M, prefer 3*M+M*NB)
673: *
674: CALL DORGBR( 'P', M, N, M, A, LDA, WORK( ITAUP ),
675: $ WORK( IWORK ), LWORK-IWORK+1, INFO )
676: IWORK = IE + M
677: *
678: * Perform bidiagonal QR iteration,
679: * computing right singular vectors of A in A and
680: * multiplying B by transpose of left singular vectors
681: * (Workspace: need BDSPAC)
682: *
1.9 ! bertrand 683: CALL DBDSQR( 'L', M, N, 0, NRHS, S, WORK( IE ), A, LDA, DUM,
1.1 bertrand 684: $ 1, B, LDB, WORK( IWORK ), INFO )
685: IF( INFO.NE.0 )
686: $ GO TO 70
687: *
688: * Multiply B by reciprocals of singular values
689: *
690: THR = MAX( RCOND*S( 1 ), SFMIN )
691: IF( RCOND.LT.ZERO )
692: $ THR = MAX( EPS*S( 1 ), SFMIN )
693: RANK = 0
694: DO 50 I = 1, M
695: IF( S( I ).GT.THR ) THEN
696: CALL DRSCL( NRHS, S( I ), B( I, 1 ), LDB )
697: RANK = RANK + 1
698: ELSE
699: CALL DLASET( 'F', 1, NRHS, ZERO, ZERO, B( I, 1 ), LDB )
700: END IF
701: 50 CONTINUE
702: *
703: * Multiply B by right singular vectors of A
704: * (Workspace: need N, prefer N*NRHS)
705: *
706: IF( LWORK.GE.LDB*NRHS .AND. NRHS.GT.1 ) THEN
707: CALL DGEMM( 'T', 'N', N, NRHS, M, ONE, A, LDA, B, LDB, ZERO,
708: $ WORK, LDB )
709: CALL DLACPY( 'F', N, NRHS, WORK, LDB, B, LDB )
710: ELSE IF( NRHS.GT.1 ) THEN
711: CHUNK = LWORK / N
712: DO 60 I = 1, NRHS, CHUNK
713: BL = MIN( NRHS-I+1, CHUNK )
714: CALL DGEMM( 'T', 'N', N, BL, M, ONE, A, LDA, B( 1, I ),
715: $ LDB, ZERO, WORK, N )
716: CALL DLACPY( 'F', N, BL, WORK, N, B( 1, I ), LDB )
717: 60 CONTINUE
718: ELSE
719: CALL DGEMV( 'T', M, N, ONE, A, LDA, B, 1, ZERO, WORK, 1 )
720: CALL DCOPY( N, WORK, 1, B, 1 )
721: END IF
722: END IF
723: *
724: * Undo scaling
725: *
726: IF( IASCL.EQ.1 ) THEN
727: CALL DLASCL( 'G', 0, 0, ANRM, SMLNUM, N, NRHS, B, LDB, INFO )
728: CALL DLASCL( 'G', 0, 0, SMLNUM, ANRM, MINMN, 1, S, MINMN,
729: $ INFO )
730: ELSE IF( IASCL.EQ.2 ) THEN
731: CALL DLASCL( 'G', 0, 0, ANRM, BIGNUM, N, NRHS, B, LDB, INFO )
732: CALL DLASCL( 'G', 0, 0, BIGNUM, ANRM, MINMN, 1, S, MINMN,
733: $ INFO )
734: END IF
735: IF( IBSCL.EQ.1 ) THEN
736: CALL DLASCL( 'G', 0, 0, SMLNUM, BNRM, N, NRHS, B, LDB, INFO )
737: ELSE IF( IBSCL.EQ.2 ) THEN
738: CALL DLASCL( 'G', 0, 0, BIGNUM, BNRM, N, NRHS, B, LDB, INFO )
739: END IF
740: *
741: 70 CONTINUE
742: WORK( 1 ) = MAXWRK
743: RETURN
744: *
745: * End of DGELSS
746: *
747: END
CVSweb interface <joel.bertrand@systella.fr>