3 * =========== DOCUMENTATION ===========
5 * Online html documentation available at
6 * http://www.netlib.org/lapack/explore-html/
9 *> Download SGESDD + dependencies
10 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/sgesdd.f">
12 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/sgesdd.f">
14 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/sgesdd.f">
21 * SUBROUTINE SGESDD( JOBZ, M, N, A, LDA, S, U, LDU, VT, LDVT,
22 * WORK, LWORK, IWORK, INFO )
24 * .. Scalar Arguments ..
26 * INTEGER INFO, LDA, LDU, LDVT, LWORK, M, N
28 * .. Array Arguments ..
30 * REAL A( LDA, * ), S( * ), U( LDU, * ),
31 * $ VT( LDVT, * ), WORK( * )
40 *> SGESDD computes the singular value decomposition (SVD) of a real
41 *> M-by-N matrix A, optionally computing the left and right singular
42 *> vectors. If singular vectors are desired, it uses a
43 *> divide-and-conquer algorithm.
47 *> A = U * SIGMA * transpose(V)
49 *> where SIGMA is an M-by-N matrix which is zero except for its
50 *> min(m,n) diagonal elements, U is an M-by-M orthogonal matrix, and
51 *> V is an N-by-N orthogonal matrix. The diagonal elements of SIGMA
52 *> are the singular values of A; they are real and non-negative, and
53 *> are returned in descending order. The first min(m,n) columns of
54 *> U and V are the left and right singular vectors of A.
56 *> Note that the routine returns VT = V**T, not V.
58 *> The divide and conquer algorithm makes very mild assumptions about
59 *> floating point arithmetic. It will work on machines with a guard
60 *> digit in add/subtract, or on those binary machines without guard
61 *> digits which subtract like the Cray X-MP, Cray Y-MP, Cray C-90, or
62 *> Cray-2. It could conceivably fail on hexadecimal or decimal machines
63 *> without guard digits, but we know of none.
71 *> JOBZ is CHARACTER*1
72 *> Specifies options for computing all or part of the matrix U:
73 *> = 'A': all M columns of U and all N rows of V**T are
74 *> returned in the arrays U and VT;
75 *> = 'S': the first min(M,N) columns of U and the first
76 *> min(M,N) rows of V**T are returned in the arrays U
78 *> = 'O': If M >= N, the first N columns of U are overwritten
79 *> on the array A and all rows of V**T are returned in
81 *> otherwise, all columns of U are returned in the
82 *> array U and the first M rows of V**T are overwritten
84 *> = 'N': no columns of U or rows of V**T are computed.
90 *> The number of rows of the input matrix A. M >= 0.
96 *> The number of columns of the input matrix A. N >= 0.
101 *> A is REAL array, dimension (LDA,N)
102 *> On entry, the M-by-N matrix A.
104 *> if JOBZ = 'O', A is overwritten with the first N columns
105 *> of U (the left singular vectors, stored
106 *> columnwise) if M >= N;
107 *> A is overwritten with the first M rows
108 *> of V**T (the right singular vectors, stored
109 *> rowwise) otherwise.
110 *> if JOBZ .ne. 'O', the contents of A are destroyed.
116 *> The leading dimension of the array A. LDA >= max(1,M).
121 *> S is REAL array, dimension (min(M,N))
122 *> The singular values of A, sorted so that S(i) >= S(i+1).
127 *> U is REAL array, dimension (LDU,UCOL)
128 *> UCOL = M if JOBZ = 'A' or JOBZ = 'O' and M < N;
129 *> UCOL = min(M,N) if JOBZ = 'S'.
130 *> If JOBZ = 'A' or JOBZ = 'O' and M < N, U contains the M-by-M
131 *> orthogonal matrix U;
132 *> if JOBZ = 'S', U contains the first min(M,N) columns of U
133 *> (the left singular vectors, stored columnwise);
134 *> if JOBZ = 'O' and M >= N, or JOBZ = 'N', U is not referenced.
140 *> The leading dimension of the array U. LDU >= 1; if
141 *> JOBZ = 'S' or 'A' or JOBZ = 'O' and M < N, LDU >= M.
146 *> VT is REAL array, dimension (LDVT,N)
147 *> If JOBZ = 'A' or JOBZ = 'O' and M >= N, VT contains the
148 *> N-by-N orthogonal matrix V**T;
149 *> if JOBZ = 'S', VT contains the first min(M,N) rows of
150 *> V**T (the right singular vectors, stored rowwise);
151 *> if JOBZ = 'O' and M < N, or JOBZ = 'N', VT is not referenced.
157 *> The leading dimension of the array VT. LDVT >= 1;
158 *> if JOBZ = 'A' or JOBZ = 'O' and M >= N, LDVT >= N;
159 *> if JOBZ = 'S', LDVT >= min(M,N).
164 *> WORK is REAL array, dimension (MAX(1,LWORK))
165 *> On exit, if INFO = 0, WORK(1) returns the optimal LWORK;
171 *> The dimension of the array WORK. LWORK >= 1.
172 *> If LWORK = -1, a workspace query is assumed. The optimal
173 *> size for the WORK array is calculated and stored in WORK(1),
174 *> and no other work except argument checking is performed.
176 *> Let mx = max(M,N) and mn = min(M,N).
177 *> If JOBZ = 'N', LWORK >= 3*mn + max( mx, 7*mn ).
178 *> If JOBZ = 'O', LWORK >= 3*mn + max( mx, 5*mn*mn + 4*mn ).
179 *> If JOBZ = 'S', LWORK >= 4*mn*mn + 7*mn.
180 *> If JOBZ = 'A', LWORK >= 4*mn*mn + 6*mn + mx.
181 *> These are not tight minimums in all cases; see comments inside code.
182 *> For good performance, LWORK should generally be larger;
183 *> a query is recommended.
188 *> IWORK is INTEGER array, dimension (8*min(M,N))
194 *> = 0: successful exit.
195 *> < 0: if INFO = -i, the i-th argument had an illegal value.
196 *> > 0: SBDSDC did not converge, updating process failed.
202 *> \author Univ. of Tennessee
203 *> \author Univ. of California Berkeley
204 *> \author Univ. of Colorado Denver
209 *> \ingroup realGEsing
211 *> \par Contributors:
214 *> Ming Gu and Huan Ren, Computer Science Division, University of
215 *> California at Berkeley, USA
217 * =====================================================================
218 SUBROUTINE SGESDD( JOBZ, M, N, A, LDA, S, U, LDU, VT, LDVT,
219 $ WORK, LWORK, IWORK, INFO )
222 * -- LAPACK driver routine (version 3.6.1) --
223 * -- LAPACK is a software package provided by Univ. of Tennessee, --
224 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
227 * .. Scalar Arguments ..
229 INTEGER INFO, LDA, LDU, LDVT, LWORK, M, N
231 * .. Array Arguments ..
233 REAL A( LDA, * ), S( * ), U( LDU, * ),
234 $ VT( LDVT, * ), WORK( * )
237 * =====================================================================
241 PARAMETER ( ZERO = 0.0E0, ONE = 1.0E0 )
243 * .. Local Scalars ..
244 LOGICAL LQUERY, WNTQA, WNTQAS, WNTQN, WNTQO, WNTQS
245 INTEGER BDSPAC, BLK, CHUNK, I, IE, IERR, IL,
246 $ IR, ISCL, ITAU, ITAUP, ITAUQ, IU, IVT, LDWKVT,
247 $ LDWRKL, LDWRKR, LDWRKU, MAXWRK, MINMN, MINWRK,
248 $ MNTHR, NWORK, WRKBL
249 INTEGER LWORK_SGEBRD_MN, LWORK_SGEBRD_MM,
250 $ LWORK_SGEBRD_NN, LWORK_SGELQF_MN,
252 $ LWORK_SORGBR_P_MM, LWORK_SORGBR_Q_NN,
253 $ LWORK_SORGLQ_MN, LWORK_SORGLQ_NN,
254 $ LWORK_SORGQR_MM, LWORK_SORGQR_MN,
255 $ LWORK_SORMBR_PRT_MM, LWORK_SORMBR_QLN_MM,
256 $ LWORK_SORMBR_PRT_MN, LWORK_SORMBR_QLN_MN,
257 $ LWORK_SORMBR_PRT_NN, LWORK_SORMBR_QLN_NN
258 REAL ANRM, BIGNUM, EPS, SMLNUM
264 * .. External Subroutines ..
265 EXTERNAL SBDSDC, SGEBRD, SGELQF, SGEMM, SGEQRF, SLACPY,
266 $ SLASCL, SLASET, SORGBR, SORGLQ, SORGQR, SORMBR,
269 * .. External Functions ..
272 EXTERNAL SLAMCH, SLANGE, LSAME
274 * .. Intrinsic Functions ..
275 INTRINSIC INT, MAX, MIN, SQRT
277 * .. Executable Statements ..
279 * Test the input arguments
283 WNTQA = LSAME( JOBZ, 'A' )
284 WNTQS = LSAME( JOBZ, 'S' )
285 WNTQAS = WNTQA .OR. WNTQS
286 WNTQO = LSAME( JOBZ, 'O' )
287 WNTQN = LSAME( JOBZ, 'N' )
288 LQUERY = ( LWORK.EQ.-1 )
290 IF( .NOT.( WNTQA .OR. WNTQS .OR. WNTQO .OR. WNTQN ) ) THEN
292 ELSE IF( M.LT.0 ) THEN
294 ELSE IF( N.LT.0 ) THEN
296 ELSE IF( LDA.LT.MAX( 1, M ) ) THEN
298 ELSE IF( LDU.LT.1 .OR. ( WNTQAS .AND. LDU.LT.M ) .OR.
299 $ ( WNTQO .AND. M.LT.N .AND. LDU.LT.M ) ) THEN
301 ELSE IF( LDVT.LT.1 .OR. ( WNTQA .AND. LDVT.LT.N ) .OR.
302 $ ( WNTQS .AND. LDVT.LT.MINMN ) .OR.
303 $ ( WNTQO .AND. M.GE.N .AND. LDVT.LT.N ) ) THEN
308 * Note: Comments in the code beginning "Workspace:" describe the
309 * minimal amount of workspace allocated at that point in the code,
310 * as well as the preferred amount for good performance.
311 * NB refers to the optimal block size for the immediately
312 * following subroutine, as returned by ILAENV.
318 MNTHR = INT( MINMN*11.0E0 / 6.0E0 )
319 IF( M.GE.N .AND. MINMN.GT.0 ) THEN
321 * Compute space needed for SBDSDC
324 * sbdsdc needs only 4*N (or 6*N for uplo=L for LAPACK <= 3.6)
325 * keep 7*N for backwards compatability.
331 * Compute space preferred for each routine
332 CALL SGEBRD( M, N, DUM(1), M, DUM(1), DUM(1), DUM(1),
333 $ DUM(1), DUM(1), -1, IERR )
334 LWORK_SGEBRD_MN = INT( DUM(1) )
336 CALL SGEBRD( N, N, DUM(1), N, DUM(1), DUM(1), DUM(1),
337 $ DUM(1), DUM(1), -1, IERR )
338 LWORK_SGEBRD_NN = INT( DUM(1) )
340 CALL SGEQRF( M, N, DUM(1), M, DUM(1), DUM(1), -1, IERR )
341 LWORK_SGEQRF_MN = INT( DUM(1) )
343 CALL SORGBR( 'Q', N, N, N, DUM(1), N, DUM(1), DUM(1), -1,
345 LWORK_SORGBR_Q_NN = INT( DUM(1) )
347 CALL SORGQR( M, M, N, DUM(1), M, DUM(1), DUM(1), -1, IERR )
348 LWORK_SORGQR_MM = INT( DUM(1) )
350 CALL SORGQR( M, N, N, DUM(1), M, DUM(1), DUM(1), -1, IERR )
351 LWORK_SORGQR_MN = INT( DUM(1) )
353 CALL SORMBR( 'P', 'R', 'T', N, N, N, DUM(1), N,
354 $ DUM(1), DUM(1), N, DUM(1), -1, IERR )
355 LWORK_SORMBR_PRT_NN = INT( DUM(1) )
357 CALL SORMBR( 'Q', 'L', 'N', N, N, N, DUM(1), N,
358 $ DUM(1), DUM(1), N, DUM(1), -1, IERR )
359 LWORK_SORMBR_QLN_NN = INT( DUM(1) )
361 CALL SORMBR( 'Q', 'L', 'N', M, N, N, DUM(1), M,
362 $ DUM(1), DUM(1), M, DUM(1), -1, IERR )
363 LWORK_SORMBR_QLN_MN = INT( DUM(1) )
365 CALL SORMBR( 'Q', 'L', 'N', M, M, N, DUM(1), M,
366 $ DUM(1), DUM(1), M, DUM(1), -1, IERR )
367 LWORK_SORMBR_QLN_MM = INT( DUM(1) )
369 IF( M.GE.MNTHR ) THEN
372 * Path 1 (M >> N, JOBZ='N')
374 WRKBL = N + LWORK_SGEQRF_MN
375 WRKBL = MAX( WRKBL, 3*N + LWORK_SGEBRD_NN )
376 MAXWRK = MAX( WRKBL, BDSPAC + N )
378 ELSE IF( WNTQO ) THEN
380 * Path 2 (M >> N, JOBZ='O')
382 WRKBL = N + LWORK_SGEQRF_MN
383 WRKBL = MAX( WRKBL, N + LWORK_SORGQR_MN )
384 WRKBL = MAX( WRKBL, 3*N + LWORK_SGEBRD_NN )
385 WRKBL = MAX( WRKBL, 3*N + LWORK_SORMBR_QLN_NN )
386 WRKBL = MAX( WRKBL, 3*N + LWORK_SORMBR_PRT_NN )
387 WRKBL = MAX( WRKBL, 3*N + BDSPAC )
388 MAXWRK = WRKBL + 2*N*N
389 MINWRK = BDSPAC + 2*N*N + 3*N
390 ELSE IF( WNTQS ) THEN
392 * Path 3 (M >> N, JOBZ='S')
394 WRKBL = N + LWORK_SGEQRF_MN
395 WRKBL = MAX( WRKBL, N + LWORK_SORGQR_MN )
396 WRKBL = MAX( WRKBL, 3*N + LWORK_SGEBRD_NN )
397 WRKBL = MAX( WRKBL, 3*N + LWORK_SORMBR_QLN_NN )
398 WRKBL = MAX( WRKBL, 3*N + LWORK_SORMBR_PRT_NN )
399 WRKBL = MAX( WRKBL, 3*N + BDSPAC )
401 MINWRK = BDSPAC + N*N + 3*N
402 ELSE IF( WNTQA ) THEN
404 * Path 4 (M >> N, JOBZ='A')
406 WRKBL = N + LWORK_SGEQRF_MN
407 WRKBL = MAX( WRKBL, N + LWORK_SORGQR_MM )
408 WRKBL = MAX( WRKBL, 3*N + LWORK_SGEBRD_NN )
409 WRKBL = MAX( WRKBL, 3*N + LWORK_SORMBR_QLN_NN )
410 WRKBL = MAX( WRKBL, 3*N + LWORK_SORMBR_PRT_NN )
411 WRKBL = MAX( WRKBL, 3*N + BDSPAC )
413 MINWRK = N*N + MAX( 3*N + BDSPAC, N + M )
417 * Path 5 (M >= N, but not much larger)
419 WRKBL = 3*N + LWORK_SGEBRD_MN
421 * Path 5n (M >= N, jobz='N')
422 MAXWRK = MAX( WRKBL, 3*N + BDSPAC )
423 MINWRK = 3*N + MAX( M, BDSPAC )
424 ELSE IF( WNTQO ) THEN
425 * Path 5o (M >= N, jobz='O')
426 WRKBL = MAX( WRKBL, 3*N + LWORK_SORMBR_PRT_NN )
427 WRKBL = MAX( WRKBL, 3*N + LWORK_SORMBR_QLN_MN )
428 WRKBL = MAX( WRKBL, 3*N + BDSPAC )
430 MINWRK = 3*N + MAX( M, N*N + BDSPAC )
431 ELSE IF( WNTQS ) THEN
432 * Path 5s (M >= N, jobz='S')
433 WRKBL = MAX( WRKBL, 3*N + LWORK_SORMBR_QLN_MN )
434 WRKBL = MAX( WRKBL, 3*N + LWORK_SORMBR_PRT_NN )
435 MAXWRK = MAX( WRKBL, 3*N + BDSPAC )
436 MINWRK = 3*N + MAX( M, BDSPAC )
437 ELSE IF( WNTQA ) THEN
438 * Path 5a (M >= N, jobz='A')
439 WRKBL = MAX( WRKBL, 3*N + LWORK_SORMBR_QLN_MM )
440 WRKBL = MAX( WRKBL, 3*N + LWORK_SORMBR_PRT_NN )
441 MAXWRK = MAX( WRKBL, 3*N + BDSPAC )
442 MINWRK = 3*N + MAX( M, BDSPAC )
445 ELSE IF( MINMN.GT.0 ) THEN
447 * Compute space needed for SBDSDC
450 * sbdsdc needs only 4*N (or 6*N for uplo=L for LAPACK <= 3.6)
451 * keep 7*N for backwards compatability.
457 * Compute space preferred for each routine
458 CALL SGEBRD( M, N, DUM(1), M, DUM(1), DUM(1), DUM(1),
459 $ DUM(1), DUM(1), -1, IERR )
460 LWORK_SGEBRD_MN = INT( DUM(1) )
462 CALL SGEBRD( M, M, A, M, S, DUM(1), DUM(1),
463 $ DUM(1), DUM(1), -1, IERR )
464 LWORK_SGEBRD_MM = INT( DUM(1) )
466 CALL SGELQF( M, N, A, M, DUM(1), DUM(1), -1, IERR )
467 LWORK_SGELQF_MN = INT( DUM(1) )
469 CALL SORGLQ( N, N, M, DUM(1), N, DUM(1), DUM(1), -1, IERR )
470 LWORK_SORGLQ_NN = INT( DUM(1) )
472 CALL SORGLQ( M, N, M, A, M, DUM(1), DUM(1), -1, IERR )
473 LWORK_SORGLQ_MN = INT( DUM(1) )
475 CALL SORGBR( 'P', M, M, M, A, N, DUM(1), DUM(1), -1, IERR )
476 LWORK_SORGBR_P_MM = INT( DUM(1) )
478 CALL SORMBR( 'P', 'R', 'T', M, M, M, DUM(1), M,
479 $ DUM(1), DUM(1), M, DUM(1), -1, IERR )
480 LWORK_SORMBR_PRT_MM = INT( DUM(1) )
482 CALL SORMBR( 'P', 'R', 'T', M, N, M, DUM(1), M,
483 $ DUM(1), DUM(1), M, DUM(1), -1, IERR )
484 LWORK_SORMBR_PRT_MN = INT( DUM(1) )
486 CALL SORMBR( 'P', 'R', 'T', N, N, M, DUM(1), N,
487 $ DUM(1), DUM(1), N, DUM(1), -1, IERR )
488 LWORK_SORMBR_PRT_NN = INT( DUM(1) )
490 CALL SORMBR( 'Q', 'L', 'N', M, M, M, DUM(1), M,
491 $ DUM(1), DUM(1), M, DUM(1), -1, IERR )
492 LWORK_SORMBR_QLN_MM = INT( DUM(1) )
494 IF( N.GE.MNTHR ) THEN
497 * Path 1t (N >> M, JOBZ='N')
499 WRKBL = M + LWORK_SGELQF_MN
500 WRKBL = MAX( WRKBL, 3*M + LWORK_SGEBRD_MM )
501 MAXWRK = MAX( WRKBL, BDSPAC + M )
503 ELSE IF( WNTQO ) THEN
505 * Path 2t (N >> M, JOBZ='O')
507 WRKBL = M + LWORK_SGELQF_MN
508 WRKBL = MAX( WRKBL, M + LWORK_SORGLQ_MN )
509 WRKBL = MAX( WRKBL, 3*M + LWORK_SGEBRD_MM )
510 WRKBL = MAX( WRKBL, 3*M + LWORK_SORMBR_QLN_MM )
511 WRKBL = MAX( WRKBL, 3*M + LWORK_SORMBR_PRT_MM )
512 WRKBL = MAX( WRKBL, 3*M + BDSPAC )
513 MAXWRK = WRKBL + 2*M*M
514 MINWRK = BDSPAC + 2*M*M + 3*M
515 ELSE IF( WNTQS ) THEN
517 * Path 3t (N >> M, JOBZ='S')
519 WRKBL = M + LWORK_SGELQF_MN
520 WRKBL = MAX( WRKBL, M + LWORK_SORGLQ_MN )
521 WRKBL = MAX( WRKBL, 3*M + LWORK_SGEBRD_MM )
522 WRKBL = MAX( WRKBL, 3*M + LWORK_SORMBR_QLN_MM )
523 WRKBL = MAX( WRKBL, 3*M + LWORK_SORMBR_PRT_MM )
524 WRKBL = MAX( WRKBL, 3*M + BDSPAC )
526 MINWRK = BDSPAC + M*M + 3*M
527 ELSE IF( WNTQA ) THEN
529 * Path 4t (N >> M, JOBZ='A')
531 WRKBL = M + LWORK_SGELQF_MN
532 WRKBL = MAX( WRKBL, M + LWORK_SORGLQ_NN )
533 WRKBL = MAX( WRKBL, 3*M + LWORK_SGEBRD_MM )
534 WRKBL = MAX( WRKBL, 3*M + LWORK_SORMBR_QLN_MM )
535 WRKBL = MAX( WRKBL, 3*M + LWORK_SORMBR_PRT_MM )
536 WRKBL = MAX( WRKBL, 3*M + BDSPAC )
538 MINWRK = M*M + MAX( 3*M + BDSPAC, M + N )
542 * Path 5t (N > M, but not much larger)
544 WRKBL = 3*M + LWORK_SGEBRD_MN
546 * Path 5tn (N > M, jobz='N')
547 MAXWRK = MAX( WRKBL, 3*M + BDSPAC )
548 MINWRK = 3*M + MAX( N, BDSPAC )
549 ELSE IF( WNTQO ) THEN
550 * Path 5to (N > M, jobz='O')
551 WRKBL = MAX( WRKBL, 3*M + LWORK_SORMBR_QLN_MM )
552 WRKBL = MAX( WRKBL, 3*M + LWORK_SORMBR_PRT_MN )
553 WRKBL = MAX( WRKBL, 3*M + BDSPAC )
555 MINWRK = 3*M + MAX( N, M*M + BDSPAC )
556 ELSE IF( WNTQS ) THEN
557 * Path 5ts (N > M, jobz='S')
558 WRKBL = MAX( WRKBL, 3*M + LWORK_SORMBR_QLN_MM )
559 WRKBL = MAX( WRKBL, 3*M + LWORK_SORMBR_PRT_MN )
560 MAXWRK = MAX( WRKBL, 3*M + BDSPAC )
561 MINWRK = 3*M + MAX( N, BDSPAC )
562 ELSE IF( WNTQA ) THEN
563 * Path 5ta (N > M, jobz='A')
564 WRKBL = MAX( WRKBL, 3*M + LWORK_SORMBR_QLN_MM )
565 WRKBL = MAX( WRKBL, 3*M + LWORK_SORMBR_PRT_NN )
566 MAXWRK = MAX( WRKBL, 3*M + BDSPAC )
567 MINWRK = 3*M + MAX( N, BDSPAC )
572 MAXWRK = MAX( MAXWRK, MINWRK )
575 IF( LWORK.LT.MINWRK .AND. .NOT.LQUERY ) THEN
581 CALL XERBLA( 'SGESDD', -INFO )
583 ELSE IF( LQUERY ) THEN
587 * Quick return if possible
589 IF( M.EQ.0 .OR. N.EQ.0 ) THEN
593 * Get machine constants
596 SMLNUM = SQRT( SLAMCH( 'S' ) ) / EPS
597 BIGNUM = ONE / SMLNUM
599 * Scale A if max element outside range [SMLNUM,BIGNUM]
601 ANRM = SLANGE( 'M', M, N, A, LDA, DUM )
603 IF( ANRM.GT.ZERO .AND. ANRM.LT.SMLNUM ) THEN
605 CALL SLASCL( 'G', 0, 0, ANRM, SMLNUM, M, N, A, LDA, IERR )
606 ELSE IF( ANRM.GT.BIGNUM ) THEN
608 CALL SLASCL( 'G', 0, 0, ANRM, BIGNUM, M, N, A, LDA, IERR )
613 * A has at least as many rows as columns. If A has sufficiently
614 * more rows than columns, first reduce using the QR
615 * decomposition (if sufficient workspace available)
617 IF( M.GE.MNTHR ) THEN
621 * Path 1 (M >> N, JOBZ='N')
622 * No singular vectors to be computed
628 * Workspace: need N [tau] + N [work]
629 * Workspace: prefer N [tau] + N*NB [work]
631 CALL SGEQRF( M, N, A, LDA, WORK( ITAU ), WORK( NWORK ),
632 $ LWORK - NWORK + 1, IERR )
636 CALL SLASET( 'L', N-1, N-1, ZERO, ZERO, A( 2, 1 ), LDA )
642 * Bidiagonalize R in A
643 * Workspace: need 3*N [e, tauq, taup] + N [work]
644 * Workspace: prefer 3*N [e, tauq, taup] + 2*N*NB [work]
646 CALL SGEBRD( N, N, A, LDA, S, WORK( IE ), WORK( ITAUQ ),
647 $ WORK( ITAUP ), WORK( NWORK ), LWORK-NWORK+1,
651 * Perform bidiagonal SVD, computing singular values only
652 * Workspace: need N [e] + BDSPAC
654 CALL SBDSDC( 'U', 'N', N, S, WORK( IE ), DUM, 1, DUM, 1,
655 $ DUM, IDUM, WORK( NWORK ), IWORK, INFO )
657 ELSE IF( WNTQO ) THEN
659 * Path 2 (M >> N, JOBZ = 'O')
660 * N left singular vectors to be overwritten on A and
661 * N right singular vectors to be computed in VT
665 * WORK(IR) is LDWRKR by N
667 IF( LWORK .GE. LDA*N + N*N + 3*N + BDSPAC ) THEN
670 LDWRKR = ( LWORK - N*N - 3*N - BDSPAC ) / N
676 * Workspace: need N*N [R] + N [tau] + N [work]
677 * Workspace: prefer N*N [R] + N [tau] + N*NB [work]
679 CALL SGEQRF( M, N, A, LDA, WORK( ITAU ), WORK( NWORK ),
680 $ LWORK - NWORK + 1, IERR )
682 * Copy R to WORK(IR), zeroing out below it
684 CALL SLACPY( 'U', N, N, A, LDA, WORK( IR ), LDWRKR )
685 CALL SLASET( 'L', N - 1, N - 1, ZERO, ZERO, WORK(IR+1),
689 * Workspace: need N*N [R] + N [tau] + N [work]
690 * Workspace: prefer N*N [R] + N [tau] + N*NB [work]
692 CALL SORGQR( M, N, N, A, LDA, WORK( ITAU ),
693 $ WORK( NWORK ), LWORK - NWORK + 1, IERR )
699 * Bidiagonalize R in WORK(IR)
700 * Workspace: need N*N [R] + 3*N [e, tauq, taup] + N [work]
701 * Workspace: prefer N*N [R] + 3*N [e, tauq, taup] + 2*N*NB [work]
703 CALL SGEBRD( N, N, WORK( IR ), LDWRKR, S, WORK( IE ),
704 $ WORK( ITAUQ ), WORK( ITAUP ), WORK( NWORK ),
705 $ LWORK - NWORK + 1, IERR )
712 * Perform bidiagonal SVD, computing left singular vectors
713 * of bidiagonal matrix in WORK(IU) and computing right
714 * singular vectors of bidiagonal matrix in VT
715 * Workspace: need N*N [R] + 3*N [e, tauq, taup] + N*N [U] + BDSPAC
717 CALL SBDSDC( 'U', 'I', N, S, WORK( IE ), WORK( IU ), N,
718 $ VT, LDVT, DUM, IDUM, WORK( NWORK ), IWORK,
721 * Overwrite WORK(IU) by left singular vectors of R
722 * and VT by right singular vectors of R
723 * Workspace: need N*N [R] + 3*N [e, tauq, taup] + N*N [U] + N [work]
724 * Workspace: prefer N*N [R] + 3*N [e, tauq, taup] + N*N [U] + N*NB [work]
726 CALL SORMBR( 'Q', 'L', 'N', N, N, N, WORK( IR ), LDWRKR,
727 $ WORK( ITAUQ ), WORK( IU ), N, WORK( NWORK ),
728 $ LWORK - NWORK + 1, IERR )
729 CALL SORMBR( 'P', 'R', 'T', N, N, N, WORK( IR ), LDWRKR,
730 $ WORK( ITAUP ), VT, LDVT, WORK( NWORK ),
731 $ LWORK - NWORK + 1, IERR )
733 * Multiply Q in A by left singular vectors of R in
734 * WORK(IU), storing result in WORK(IR) and copying to A
735 * Workspace: need N*N [R] + 3*N [e, tauq, taup] + N*N [U]
736 * Workspace: prefer M*N [R] + 3*N [e, tauq, taup] + N*N [U]
738 DO 10 I = 1, M, LDWRKR
739 CHUNK = MIN( M - I + 1, LDWRKR )
740 CALL SGEMM( 'N', 'N', CHUNK, N, N, ONE, A( I, 1 ),
741 $ LDA, WORK( IU ), N, ZERO, WORK( IR ),
743 CALL SLACPY( 'F', CHUNK, N, WORK( IR ), LDWRKR,
747 ELSE IF( WNTQS ) THEN
749 * Path 3 (M >> N, JOBZ='S')
750 * N left singular vectors to be computed in U and
751 * N right singular vectors to be computed in VT
762 * Workspace: need N*N [R] + N [tau] + N [work]
763 * Workspace: prefer N*N [R] + N [tau] + N*NB [work]
765 CALL SGEQRF( M, N, A, LDA, WORK( ITAU ), WORK( NWORK ),
766 $ LWORK - NWORK + 1, IERR )
768 * Copy R to WORK(IR), zeroing out below it
770 CALL SLACPY( 'U', N, N, A, LDA, WORK( IR ), LDWRKR )
771 CALL SLASET( 'L', N - 1, N - 1, ZERO, ZERO, WORK(IR+1),
775 * Workspace: need N*N [R] + N [tau] + N [work]
776 * Workspace: prefer N*N [R] + N [tau] + N*NB [work]
778 CALL SORGQR( M, N, N, A, LDA, WORK( ITAU ),
779 $ WORK( NWORK ), LWORK - NWORK + 1, IERR )
785 * Bidiagonalize R in WORK(IR)
786 * Workspace: need N*N [R] + 3*N [e, tauq, taup] + N [work]
787 * Workspace: prefer N*N [R] + 3*N [e, tauq, taup] + 2*N*NB [work]
789 CALL SGEBRD( N, N, WORK( IR ), LDWRKR, S, WORK( IE ),
790 $ WORK( ITAUQ ), WORK( ITAUP ), WORK( NWORK ),
791 $ LWORK - NWORK + 1, IERR )
793 * Perform bidiagonal SVD, computing left singular vectors
794 * of bidiagoal matrix in U and computing right singular
795 * vectors of bidiagonal matrix in VT
796 * Workspace: need N*N [R] + 3*N [e, tauq, taup] + BDSPAC
798 CALL SBDSDC( 'U', 'I', N, S, WORK( IE ), U, LDU, VT,
799 $ LDVT, DUM, IDUM, WORK( NWORK ), IWORK,
802 * Overwrite U by left singular vectors of R and VT
803 * by right singular vectors of R
804 * Workspace: need N*N [R] + 3*N [e, tauq, taup] + N [work]
805 * Workspace: prefer N*N [R] + 3*N [e, tauq, taup] + N*NB [work]
807 CALL SORMBR( 'Q', 'L', 'N', N, N, N, WORK( IR ), LDWRKR,
808 $ WORK( ITAUQ ), U, LDU, WORK( NWORK ),
809 $ LWORK - NWORK + 1, IERR )
811 CALL SORMBR( 'P', 'R', 'T', N, N, N, WORK( IR ), LDWRKR,
812 $ WORK( ITAUP ), VT, LDVT, WORK( NWORK ),
813 $ LWORK - NWORK + 1, IERR )
815 * Multiply Q in A by left singular vectors of R in
816 * WORK(IR), storing result in U
817 * Workspace: need N*N [R]
819 CALL SLACPY( 'F', N, N, U, LDU, WORK( IR ), LDWRKR )
820 CALL SGEMM( 'N', 'N', M, N, N, ONE, A, LDA, WORK( IR ),
821 $ LDWRKR, ZERO, U, LDU )
823 ELSE IF( WNTQA ) THEN
825 * Path 4 (M >> N, JOBZ='A')
826 * M left singular vectors to be computed in U and
827 * N right singular vectors to be computed in VT
837 * Compute A=Q*R, copying result to U
838 * Workspace: need N*N [U] + N [tau] + N [work]
839 * Workspace: prefer N*N [U] + N [tau] + N*NB [work]
841 CALL SGEQRF( M, N, A, LDA, WORK( ITAU ), WORK( NWORK ),
842 $ LWORK - NWORK + 1, IERR )
843 CALL SLACPY( 'L', M, N, A, LDA, U, LDU )
846 * Workspace: need N*N [U] + N [tau] + M [work]
847 * Workspace: prefer N*N [U] + N [tau] + M*NB [work]
848 CALL SORGQR( M, M, N, U, LDU, WORK( ITAU ),
849 $ WORK( NWORK ), LWORK - NWORK + 1, IERR )
851 * Produce R in A, zeroing out other entries
853 CALL SLASET( 'L', N-1, N-1, ZERO, ZERO, A( 2, 1 ), LDA )
859 * Bidiagonalize R in A
860 * Workspace: need N*N [U] + 3*N [e, tauq, taup] + N [work]
861 * Workspace: prefer N*N [U] + 3*N [e, tauq, taup] + 2*N*NB [work]
863 CALL SGEBRD( N, N, A, LDA, S, WORK( IE ), WORK( ITAUQ ),
864 $ WORK( ITAUP ), WORK( NWORK ), LWORK-NWORK+1,
867 * Perform bidiagonal SVD, computing left singular vectors
868 * of bidiagonal matrix in WORK(IU) and computing right
869 * singular vectors of bidiagonal matrix in VT
870 * Workspace: need N*N [U] + 3*N [e, tauq, taup] + BDSPAC
872 CALL SBDSDC( 'U', 'I', N, S, WORK( IE ), WORK( IU ), N,
873 $ VT, LDVT, DUM, IDUM, WORK( NWORK ), IWORK,
876 * Overwrite WORK(IU) by left singular vectors of R and VT
877 * by right singular vectors of R
878 * Workspace: need N*N [U] + 3*N [e, tauq, taup] + N [work]
879 * Workspace: prefer N*N [U] + 3*N [e, tauq, taup] + N*NB [work]
881 CALL SORMBR( 'Q', 'L', 'N', N, N, N, A, LDA,
882 $ WORK( ITAUQ ), WORK( IU ), LDWRKU,
883 $ WORK( NWORK ), LWORK - NWORK + 1, IERR )
884 CALL SORMBR( 'P', 'R', 'T', N, N, N, A, LDA,
885 $ WORK( ITAUP ), VT, LDVT, WORK( NWORK ),
886 $ LWORK - NWORK + 1, IERR )
888 * Multiply Q in U by left singular vectors of R in
889 * WORK(IU), storing result in A
890 * Workspace: need N*N [U]
892 CALL SGEMM( 'N', 'N', M, N, N, ONE, U, LDU, WORK( IU ),
893 $ LDWRKU, ZERO, A, LDA )
895 * Copy left singular vectors of A from A to U
897 CALL SLACPY( 'F', M, N, A, LDA, U, LDU )
905 * Path 5 (M >= N, but not much larger)
906 * Reduce to bidiagonal form without QR decomposition
914 * Workspace: need 3*N [e, tauq, taup] + M [work]
915 * Workspace: prefer 3*N [e, tauq, taup] + (M+N)*NB [work]
917 CALL SGEBRD( M, N, A, LDA, S, WORK( IE ), WORK( ITAUQ ),
918 $ WORK( ITAUP ), WORK( NWORK ), LWORK-NWORK+1,
922 * Path 5n (M >= N, JOBZ='N')
923 * Perform bidiagonal SVD, only computing singular values
924 * Workspace: need 3*N [e, tauq, taup] + BDSPAC
926 CALL SBDSDC( 'U', 'N', N, S, WORK( IE ), DUM, 1, DUM, 1,
927 $ DUM, IDUM, WORK( NWORK ), IWORK, INFO )
928 ELSE IF( WNTQO ) THEN
929 * Path 5o (M >= N, JOBZ='O')
931 IF( LWORK .GE. M*N + 3*N + BDSPAC ) THEN
933 * WORK( IU ) is M by N
936 NWORK = IU + LDWRKU*N
937 CALL SLASET( 'F', M, N, ZERO, ZERO, WORK( IU ),
939 * IR is unused; silence compile warnings
943 * WORK( IU ) is N by N
946 NWORK = IU + LDWRKU*N
948 * WORK(IR) is LDWRKR by N
951 LDWRKR = ( LWORK - N*N - 3*N ) / N
953 NWORK = IU + LDWRKU*N
955 * Perform bidiagonal SVD, computing left singular vectors
956 * of bidiagonal matrix in WORK(IU) and computing right
957 * singular vectors of bidiagonal matrix in VT
958 * Workspace: need 3*N [e, tauq, taup] + N*N [U] + BDSPAC
960 CALL SBDSDC( 'U', 'I', N, S, WORK( IE ), WORK( IU ),
961 $ LDWRKU, VT, LDVT, DUM, IDUM, WORK( NWORK ),
964 * Overwrite VT by right singular vectors of A
965 * Workspace: need 3*N [e, tauq, taup] + N*N [U] + N [work]
966 * Workspace: prefer 3*N [e, tauq, taup] + N*N [U] + N*NB [work]
968 CALL SORMBR( 'P', 'R', 'T', N, N, N, A, LDA,
969 $ WORK( ITAUP ), VT, LDVT, WORK( NWORK ),
970 $ LWORK - NWORK + 1, IERR )
972 IF( LWORK .GE. M*N + 3*N + BDSPAC ) THEN
975 * Overwrite WORK(IU) by left singular vectors of A
976 * Workspace: need 3*N [e, tauq, taup] + M*N [U] + N [work]
977 * Workspace: prefer 3*N [e, tauq, taup] + M*N [U] + N*NB [work]
979 CALL SORMBR( 'Q', 'L', 'N', M, N, N, A, LDA,
980 $ WORK( ITAUQ ), WORK( IU ), LDWRKU,
981 $ WORK( NWORK ), LWORK - NWORK + 1, IERR )
983 * Copy left singular vectors of A from WORK(IU) to A
985 CALL SLACPY( 'F', M, N, WORK( IU ), LDWRKU, A, LDA )
990 * Workspace: need 3*N [e, tauq, taup] + N*N [U] + N [work]
991 * Workspace: prefer 3*N [e, tauq, taup] + N*N [U] + N*NB [work]
993 CALL SORGBR( 'Q', M, N, N, A, LDA, WORK( ITAUQ ),
994 $ WORK( NWORK ), LWORK - NWORK + 1, IERR )
996 * Multiply Q in A by left singular vectors of
997 * bidiagonal matrix in WORK(IU), storing result in
998 * WORK(IR) and copying to A
999 * Workspace: need 3*N [e, tauq, taup] + N*N [U] + NB*N [R]
1000 * Workspace: prefer 3*N [e, tauq, taup] + N*N [U] + M*N [R]
1002 DO 20 I = 1, M, LDWRKR
1003 CHUNK = MIN( M - I + 1, LDWRKR )
1004 CALL SGEMM( 'N', 'N', CHUNK, N, N, ONE, A( I, 1 ),
1005 $ LDA, WORK( IU ), LDWRKU, ZERO,
1006 $ WORK( IR ), LDWRKR )
1007 CALL SLACPY( 'F', CHUNK, N, WORK( IR ), LDWRKR,
1012 ELSE IF( WNTQS ) THEN
1014 * Path 5s (M >= N, JOBZ='S')
1015 * Perform bidiagonal SVD, computing left singular vectors
1016 * of bidiagonal matrix in U and computing right singular
1017 * vectors of bidiagonal matrix in VT
1018 * Workspace: need 3*N [e, tauq, taup] + BDSPAC
1020 CALL SLASET( 'F', M, N, ZERO, ZERO, U, LDU )
1021 CALL SBDSDC( 'U', 'I', N, S, WORK( IE ), U, LDU, VT,
1022 $ LDVT, DUM, IDUM, WORK( NWORK ), IWORK,
1025 * Overwrite U by left singular vectors of A and VT
1026 * by right singular vectors of A
1027 * Workspace: need 3*N [e, tauq, taup] + N [work]
1028 * Workspace: prefer 3*N [e, tauq, taup] + N*NB [work]
1030 CALL SORMBR( 'Q', 'L', 'N', M, N, N, A, LDA,
1031 $ WORK( ITAUQ ), U, LDU, WORK( NWORK ),
1032 $ LWORK - NWORK + 1, IERR )
1033 CALL SORMBR( 'P', 'R', 'T', N, N, N, A, LDA,
1034 $ WORK( ITAUP ), VT, LDVT, WORK( NWORK ),
1035 $ LWORK - NWORK + 1, IERR )
1036 ELSE IF( WNTQA ) THEN
1038 * Path 5a (M >= N, JOBZ='A')
1039 * Perform bidiagonal SVD, computing left singular vectors
1040 * of bidiagonal matrix in U and computing right singular
1041 * vectors of bidiagonal matrix in VT
1042 * Workspace: need 3*N [e, tauq, taup] + BDSPAC
1044 CALL SLASET( 'F', M, M, ZERO, ZERO, U, LDU )
1045 CALL SBDSDC( 'U', 'I', N, S, WORK( IE ), U, LDU, VT,
1046 $ LDVT, DUM, IDUM, WORK( NWORK ), IWORK,
1049 * Set the right corner of U to identity matrix
1052 CALL SLASET( 'F', M - N, M - N, ZERO, ONE, U(N+1,N+1),
1056 * Overwrite U by left singular vectors of A and VT
1057 * by right singular vectors of A
1058 * Workspace: need 3*N [e, tauq, taup] + M [work]
1059 * Workspace: prefer 3*N [e, tauq, taup] + M*NB [work]
1061 CALL SORMBR( 'Q', 'L', 'N', M, M, N, A, LDA,
1062 $ WORK( ITAUQ ), U, LDU, WORK( NWORK ),
1063 $ LWORK - NWORK + 1, IERR )
1064 CALL SORMBR( 'P', 'R', 'T', N, N, M, A, LDA,
1065 $ WORK( ITAUP ), VT, LDVT, WORK( NWORK ),
1066 $ LWORK - NWORK + 1, IERR )
1073 * A has more columns than rows. If A has sufficiently more
1074 * columns than rows, first reduce using the LQ decomposition (if
1075 * sufficient workspace available)
1077 IF( N.GE.MNTHR ) THEN
1081 * Path 1t (N >> M, JOBZ='N')
1082 * No singular vectors to be computed
1088 * Workspace: need M [tau] + M [work]
1089 * Workspace: prefer M [tau] + M*NB [work]
1091 CALL SGELQF( M, N, A, LDA, WORK( ITAU ), WORK( NWORK ),
1092 $ LWORK - NWORK + 1, IERR )
1096 CALL SLASET( 'U', M-1, M-1, ZERO, ZERO, A( 1, 2 ), LDA )
1102 * Bidiagonalize L in A
1103 * Workspace: need 3*M [e, tauq, taup] + M [work]
1104 * Workspace: prefer 3*M [e, tauq, taup] + 2*M*NB [work]
1106 CALL SGEBRD( M, M, A, LDA, S, WORK( IE ), WORK( ITAUQ ),
1107 $ WORK( ITAUP ), WORK( NWORK ), LWORK-NWORK+1,
1111 * Perform bidiagonal SVD, computing singular values only
1112 * Workspace: need M [e] + BDSPAC
1114 CALL SBDSDC( 'U', 'N', M, S, WORK( IE ), DUM, 1, DUM, 1,
1115 $ DUM, IDUM, WORK( NWORK ), IWORK, INFO )
1117 ELSE IF( WNTQO ) THEN
1119 * Path 2t (N >> M, JOBZ='O')
1120 * M right singular vectors to be overwritten on A and
1121 * M left singular vectors to be computed in U
1125 * WORK(IVT) is M by M
1126 * WORK(IL) is M by M; it is later resized to M by chunk for gemm
1129 IF( LWORK .GE. M*N + M*M + 3*M + BDSPAC ) THEN
1134 CHUNK = ( LWORK - M*M ) / M
1136 ITAU = IL + LDWRKL*M
1140 * Workspace: need M*M [VT] + M*M [L] + M [tau] + M [work]
1141 * Workspace: prefer M*M [VT] + M*M [L] + M [tau] + M*NB [work]
1143 CALL SGELQF( M, N, A, LDA, WORK( ITAU ), WORK( NWORK ),
1144 $ LWORK - NWORK + 1, IERR )
1146 * Copy L to WORK(IL), zeroing about above it
1148 CALL SLACPY( 'L', M, M, A, LDA, WORK( IL ), LDWRKL )
1149 CALL SLASET( 'U', M - 1, M - 1, ZERO, ZERO,
1150 $ WORK( IL + LDWRKL ), LDWRKL )
1153 * Workspace: need M*M [VT] + M*M [L] + M [tau] + M [work]
1154 * Workspace: prefer M*M [VT] + M*M [L] + M [tau] + M*NB [work]
1156 CALL SORGLQ( M, N, M, A, LDA, WORK( ITAU ),
1157 $ WORK( NWORK ), LWORK - NWORK + 1, IERR )
1163 * Bidiagonalize L in WORK(IL)
1164 * Workspace: need M*M [VT] + M*M [L] + 3*M [e, tauq, taup] + M [work]
1165 * Workspace: prefer M*M [VT] + M*M [L] + 3*M [e, tauq, taup] + 2*M*NB [work]
1167 CALL SGEBRD( M, M, WORK( IL ), LDWRKL, S, WORK( IE ),
1168 $ WORK( ITAUQ ), WORK( ITAUP ), WORK( NWORK ),
1169 $ LWORK - NWORK + 1, IERR )
1171 * Perform bidiagonal SVD, computing left singular vectors
1172 * of bidiagonal matrix in U, and computing right singular
1173 * vectors of bidiagonal matrix in WORK(IVT)
1174 * Workspace: need M*M [VT] + M*M [L] + 3*M [e, tauq, taup] + BDSPAC
1176 CALL SBDSDC( 'U', 'I', M, S, WORK( IE ), U, LDU,
1177 $ WORK( IVT ), M, DUM, IDUM, WORK( NWORK ),
1180 * Overwrite U by left singular vectors of L and WORK(IVT)
1181 * by right singular vectors of L
1182 * Workspace: need M*M [VT] + M*M [L] + 3*M [e, tauq, taup] + M [work]
1183 * Workspace: prefer M*M [VT] + M*M [L] + 3*M [e, tauq, taup] + M*NB [work]
1185 CALL SORMBR( 'Q', 'L', 'N', M, M, M, WORK( IL ), LDWRKL,
1186 $ WORK( ITAUQ ), U, LDU, WORK( NWORK ),
1187 $ LWORK - NWORK + 1, IERR )
1188 CALL SORMBR( 'P', 'R', 'T', M, M, M, WORK( IL ), LDWRKL,
1189 $ WORK( ITAUP ), WORK( IVT ), M,
1190 $ WORK( NWORK ), LWORK - NWORK + 1, IERR )
1192 * Multiply right singular vectors of L in WORK(IVT) by Q
1193 * in A, storing result in WORK(IL) and copying to A
1194 * Workspace: need M*M [VT] + M*M [L]
1195 * Workspace: prefer M*M [VT] + M*N [L]
1196 * At this point, L is resized as M by chunk.
1198 DO 30 I = 1, N, CHUNK
1199 BLK = MIN( N - I + 1, CHUNK )
1200 CALL SGEMM( 'N', 'N', M, BLK, M, ONE, WORK( IVT ), M,
1201 $ A( 1, I ), LDA, ZERO, WORK( IL ), LDWRKL )
1202 CALL SLACPY( 'F', M, BLK, WORK( IL ), LDWRKL,
1206 ELSE IF( WNTQS ) THEN
1208 * Path 3t (N >> M, JOBZ='S')
1209 * M right singular vectors to be computed in VT and
1210 * M left singular vectors to be computed in U
1214 * WORK(IL) is M by M
1217 ITAU = IL + LDWRKL*M
1221 * Workspace: need M*M [L] + M [tau] + M [work]
1222 * Workspace: prefer M*M [L] + M [tau] + M*NB [work]
1224 CALL SGELQF( M, N, A, LDA, WORK( ITAU ), WORK( NWORK ),
1225 $ LWORK - NWORK + 1, IERR )
1227 * Copy L to WORK(IL), zeroing out above it
1229 CALL SLACPY( 'L', M, M, A, LDA, WORK( IL ), LDWRKL )
1230 CALL SLASET( 'U', M - 1, M - 1, ZERO, ZERO,
1231 $ WORK( IL + LDWRKL ), LDWRKL )
1234 * Workspace: need M*M [L] + M [tau] + M [work]
1235 * Workspace: prefer M*M [L] + M [tau] + M*NB [work]
1237 CALL SORGLQ( M, N, M, A, LDA, WORK( ITAU ),
1238 $ WORK( NWORK ), LWORK - NWORK + 1, IERR )
1244 * Bidiagonalize L in WORK(IU).
1245 * Workspace: need M*M [L] + 3*M [e, tauq, taup] + M [work]
1246 * Workspace: prefer M*M [L] + 3*M [e, tauq, taup] + 2*M*NB [work]
1248 CALL SGEBRD( M, M, WORK( IL ), LDWRKL, S, WORK( IE ),
1249 $ WORK( ITAUQ ), WORK( ITAUP ), WORK( NWORK ),
1250 $ LWORK - NWORK + 1, IERR )
1252 * Perform bidiagonal SVD, computing left singular vectors
1253 * of bidiagonal matrix in U and computing right singular
1254 * vectors of bidiagonal matrix in VT
1255 * Workspace: need M*M [L] + 3*M [e, tauq, taup] + BDSPAC
1257 CALL SBDSDC( 'U', 'I', M, S, WORK( IE ), U, LDU, VT,
1258 $ LDVT, DUM, IDUM, WORK( NWORK ), IWORK,
1261 * Overwrite U by left singular vectors of L and VT
1262 * by right singular vectors of L
1263 * Workspace: need M*M [L] + 3*M [e, tauq, taup] + M [work]
1264 * Workspace: prefer M*M [L] + 3*M [e, tauq, taup] + M*NB [work]
1266 CALL SORMBR( 'Q', 'L', 'N', M, M, M, WORK( IL ), LDWRKL,
1267 $ WORK( ITAUQ ), U, LDU, WORK( NWORK ),
1268 $ LWORK - NWORK + 1, IERR )
1269 CALL SORMBR( 'P', 'R', 'T', M, M, M, WORK( IL ), LDWRKL,
1270 $ WORK( ITAUP ), VT, LDVT, WORK( NWORK ),
1271 $ LWORK - NWORK + 1, IERR )
1273 * Multiply right singular vectors of L in WORK(IL) by
1274 * Q in A, storing result in VT
1275 * Workspace: need M*M [L]
1277 CALL SLACPY( 'F', M, M, VT, LDVT, WORK( IL ), LDWRKL )
1278 CALL SGEMM( 'N', 'N', M, N, M, ONE, WORK( IL ), LDWRKL,
1279 $ A, LDA, ZERO, VT, LDVT )
1281 ELSE IF( WNTQA ) THEN
1283 * Path 4t (N >> M, JOBZ='A')
1284 * N right singular vectors to be computed in VT and
1285 * M left singular vectors to be computed in U
1289 * WORK(IVT) is M by M
1292 ITAU = IVT + LDWKVT*M
1295 * Compute A=L*Q, copying result to VT
1296 * Workspace: need M*M [VT] + M [tau] + M [work]
1297 * Workspace: prefer M*M [VT] + M [tau] + M*NB [work]
1299 CALL SGELQF( M, N, A, LDA, WORK( ITAU ), WORK( NWORK ),
1300 $ LWORK - NWORK + 1, IERR )
1301 CALL SLACPY( 'U', M, N, A, LDA, VT, LDVT )
1304 * Workspace: need M*M [VT] + M [tau] + N [work]
1305 * Workspace: prefer M*M [VT] + M [tau] + N*NB [work]
1307 CALL SORGLQ( N, N, M, VT, LDVT, WORK( ITAU ),
1308 $ WORK( NWORK ), LWORK - NWORK + 1, IERR )
1310 * Produce L in A, zeroing out other entries
1312 CALL SLASET( 'U', M-1, M-1, ZERO, ZERO, A( 1, 2 ), LDA )
1318 * Bidiagonalize L in A
1319 * Workspace: need M*M [VT] + 3*M [e, tauq, taup] + M [work]
1320 * Workspace: prefer M*M [VT] + 3*M [e, tauq, taup] + 2*M*NB [work]
1322 CALL SGEBRD( M, M, A, LDA, S, WORK( IE ), WORK( ITAUQ ),
1323 $ WORK( ITAUP ), WORK( NWORK ), LWORK-NWORK+1,
1326 * Perform bidiagonal SVD, computing left singular vectors
1327 * of bidiagonal matrix in U and computing right singular
1328 * vectors of bidiagonal matrix in WORK(IVT)
1329 * Workspace: need M*M [VT] + 3*M [e, tauq, taup] + BDSPAC
1331 CALL SBDSDC( 'U', 'I', M, S, WORK( IE ), U, LDU,
1332 $ WORK( IVT ), LDWKVT, DUM, IDUM,
1333 $ WORK( NWORK ), IWORK, INFO )
1335 * Overwrite U by left singular vectors of L and WORK(IVT)
1336 * by right singular vectors of L
1337 * Workspace: need M*M [VT] + 3*M [e, tauq, taup]+ M [work]
1338 * Workspace: prefer M*M [VT] + 3*M [e, tauq, taup]+ M*NB [work]
1340 CALL SORMBR( 'Q', 'L', 'N', M, M, M, A, LDA,
1341 $ WORK( ITAUQ ), U, LDU, WORK( NWORK ),
1342 $ LWORK - NWORK + 1, IERR )
1343 CALL SORMBR( 'P', 'R', 'T', M, M, M, A, LDA,
1344 $ WORK( ITAUP ), WORK( IVT ), LDWKVT,
1345 $ WORK( NWORK ), LWORK - NWORK + 1, IERR )
1347 * Multiply right singular vectors of L in WORK(IVT) by
1348 * Q in VT, storing result in A
1349 * Workspace: need M*M [VT]
1351 CALL SGEMM( 'N', 'N', M, N, M, ONE, WORK( IVT ), LDWKVT,
1352 $ VT, LDVT, ZERO, A, LDA )
1354 * Copy right singular vectors of A from A to VT
1356 CALL SLACPY( 'F', M, N, A, LDA, VT, LDVT )
1364 * Path 5t (N > M, but not much larger)
1365 * Reduce to bidiagonal form without LQ decomposition
1373 * Workspace: need 3*M [e, tauq, taup] + N [work]
1374 * Workspace: prefer 3*M [e, tauq, taup] + (M+N)*NB [work]
1376 CALL SGEBRD( M, N, A, LDA, S, WORK( IE ), WORK( ITAUQ ),
1377 $ WORK( ITAUP ), WORK( NWORK ), LWORK-NWORK+1,
1381 * Path 5tn (N > M, JOBZ='N')
1382 * Perform bidiagonal SVD, only computing singular values
1383 * Workspace: need 3*M [e, tauq, taup] + BDSPAC
1385 CALL SBDSDC( 'L', 'N', M, S, WORK( IE ), DUM, 1, DUM, 1,
1386 $ DUM, IDUM, WORK( NWORK ), IWORK, INFO )
1387 ELSE IF( WNTQO ) THEN
1388 * Path 5to (N > M, JOBZ='O')
1391 IF( LWORK .GE. M*N + 3*M + BDSPAC ) THEN
1393 * WORK( IVT ) is M by N
1395 CALL SLASET( 'F', M, N, ZERO, ZERO, WORK( IVT ),
1397 NWORK = IVT + LDWKVT*N
1398 * IL is unused; silence compile warnings
1402 * WORK( IVT ) is M by M
1404 NWORK = IVT + LDWKVT*M
1407 * WORK(IL) is M by CHUNK
1409 CHUNK = ( LWORK - M*M - 3*M ) / M
1412 * Perform bidiagonal SVD, computing left singular vectors
1413 * of bidiagonal matrix in U and computing right singular
1414 * vectors of bidiagonal matrix in WORK(IVT)
1415 * Workspace: need 3*M [e, tauq, taup] + M*M [VT] + BDSPAC
1417 CALL SBDSDC( 'L', 'I', M, S, WORK( IE ), U, LDU,
1418 $ WORK( IVT ), LDWKVT, DUM, IDUM,
1419 $ WORK( NWORK ), IWORK, INFO )
1421 * Overwrite U by left singular vectors of A
1422 * Workspace: need 3*M [e, tauq, taup] + M*M [VT] + M [work]
1423 * Workspace: prefer 3*M [e, tauq, taup] + M*M [VT] + M*NB [work]
1425 CALL SORMBR( 'Q', 'L', 'N', M, M, N, A, LDA,
1426 $ WORK( ITAUQ ), U, LDU, WORK( NWORK ),
1427 $ LWORK - NWORK + 1, IERR )
1429 IF( LWORK .GE. M*N + 3*M + BDSPAC ) THEN
1432 * Overwrite WORK(IVT) by left singular vectors of A
1433 * Workspace: need 3*M [e, tauq, taup] + M*N [VT] + M [work]
1434 * Workspace: prefer 3*M [e, tauq, taup] + M*N [VT] + M*NB [work]
1436 CALL SORMBR( 'P', 'R', 'T', M, N, M, A, LDA,
1437 $ WORK( ITAUP ), WORK( IVT ), LDWKVT,
1438 $ WORK( NWORK ), LWORK - NWORK + 1, IERR )
1440 * Copy right singular vectors of A from WORK(IVT) to A
1442 CALL SLACPY( 'F', M, N, WORK( IVT ), LDWKVT, A, LDA )
1446 * Generate P**T in A
1447 * Workspace: need 3*M [e, tauq, taup] + M*M [VT] + M [work]
1448 * Workspace: prefer 3*M [e, tauq, taup] + M*M [VT] + M*NB [work]
1450 CALL SORGBR( 'P', M, N, M, A, LDA, WORK( ITAUP ),
1451 $ WORK( NWORK ), LWORK - NWORK + 1, IERR )
1453 * Multiply Q in A by right singular vectors of
1454 * bidiagonal matrix in WORK(IVT), storing result in
1455 * WORK(IL) and copying to A
1456 * Workspace: need 3*M [e, tauq, taup] + M*M [VT] + M*NB [L]
1457 * Workspace: prefer 3*M [e, tauq, taup] + M*M [VT] + M*N [L]
1459 DO 40 I = 1, N, CHUNK
1460 BLK = MIN( N - I + 1, CHUNK )
1461 CALL SGEMM( 'N', 'N', M, BLK, M, ONE, WORK( IVT ),
1462 $ LDWKVT, A( 1, I ), LDA, ZERO,
1464 CALL SLACPY( 'F', M, BLK, WORK( IL ), M, A( 1, I ),
1468 ELSE IF( WNTQS ) THEN
1470 * Path 5ts (N > M, JOBZ='S')
1471 * Perform bidiagonal SVD, computing left singular vectors
1472 * of bidiagonal matrix in U and computing right singular
1473 * vectors of bidiagonal matrix in VT
1474 * Workspace: need 3*M [e, tauq, taup] + BDSPAC
1476 CALL SLASET( 'F', M, N, ZERO, ZERO, VT, LDVT )
1477 CALL SBDSDC( 'L', 'I', M, S, WORK( IE ), U, LDU, VT,
1478 $ LDVT, DUM, IDUM, WORK( NWORK ), IWORK,
1481 * Overwrite U by left singular vectors of A and VT
1482 * by right singular vectors of A
1483 * Workspace: need 3*M [e, tauq, taup] + M [work]
1484 * Workspace: prefer 3*M [e, tauq, taup] + M*NB [work]
1486 CALL SORMBR( 'Q', 'L', 'N', M, M, N, A, LDA,
1487 $ WORK( ITAUQ ), U, LDU, WORK( NWORK ),
1488 $ LWORK - NWORK + 1, IERR )
1489 CALL SORMBR( 'P', 'R', 'T', M, N, M, A, LDA,
1490 $ WORK( ITAUP ), VT, LDVT, WORK( NWORK ),
1491 $ LWORK - NWORK + 1, IERR )
1492 ELSE IF( WNTQA ) THEN
1494 * Path 5ta (N > M, JOBZ='A')
1495 * Perform bidiagonal SVD, computing left singular vectors
1496 * of bidiagonal matrix in U and computing right singular
1497 * vectors of bidiagonal matrix in VT
1498 * Workspace: need 3*M [e, tauq, taup] + BDSPAC
1500 CALL SLASET( 'F', N, N, ZERO, ZERO, VT, LDVT )
1501 CALL SBDSDC( 'L', 'I', M, S, WORK( IE ), U, LDU, VT,
1502 $ LDVT, DUM, IDUM, WORK( NWORK ), IWORK,
1505 * Set the right corner of VT to identity matrix
1508 CALL SLASET( 'F', N-M, N-M, ZERO, ONE, VT(M+1,M+1),
1512 * Overwrite U by left singular vectors of A and VT
1513 * by right singular vectors of A
1514 * Workspace: need 3*M [e, tauq, taup] + N [work]
1515 * Workspace: prefer 3*M [e, tauq, taup] + N*NB [work]
1517 CALL SORMBR( 'Q', 'L', 'N', M, M, N, A, LDA,
1518 $ WORK( ITAUQ ), U, LDU, WORK( NWORK ),
1519 $ LWORK - NWORK + 1, IERR )
1520 CALL SORMBR( 'P', 'R', 'T', N, N, M, A, LDA,
1521 $ WORK( ITAUP ), VT, LDVT, WORK( NWORK ),
1522 $ LWORK - NWORK + 1, IERR )
1529 * Undo scaling if necessary
1531 IF( ISCL.EQ.1 ) THEN
1532 IF( ANRM.GT.BIGNUM )
1533 $ CALL SLASCL( 'G', 0, 0, BIGNUM, ANRM, MINMN, 1, S, MINMN,
1535 IF( ANRM.LT.SMLNUM )
1536 $ CALL SLASCL( 'G', 0, 0, SMLNUM, ANRM, MINMN, 1, S, MINMN,
1540 * Return optimal workspace in WORK(1)