LAPACK  3.9.1
LAPACK: Linear Algebra PACKage
zgesdd.f
Go to the documentation of this file.
1 *> \brief \b ZGESDD
2 *
3 * =========== DOCUMENTATION ===========
4 *
5 * Online html documentation available at
6 * http://www.netlib.org/lapack/explore-html/
7 *
8 *> \htmlonly
9 *> Download ZGESDD + dependencies
10 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/zgesdd.f">
11 *> [TGZ]</a>
12 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/zgesdd.f">
13 *> [ZIP]</a>
14 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/zgesdd.f">
15 *> [TXT]</a>
16 *> \endhtmlonly
17 *
18 * Definition:
19 * ===========
20 *
21 * SUBROUTINE ZGESDD( JOBZ, M, N, A, LDA, S, U, LDU, VT, LDVT,
22 * WORK, LWORK, RWORK, IWORK, INFO )
23 *
24 * .. Scalar Arguments ..
25 * CHARACTER JOBZ
26 * INTEGER INFO, LDA, LDU, LDVT, LWORK, M, N
27 * ..
28 * .. Array Arguments ..
29 * INTEGER IWORK( * )
30 * DOUBLE PRECISION RWORK( * ), S( * )
31 * COMPLEX*16 A( LDA, * ), U( LDU, * ), VT( LDVT, * ),
32 * $ WORK( * )
33 * ..
34 *
35 *
36 *> \par Purpose:
37 * =============
38 *>
39 *> \verbatim
40 *>
41 *> ZGESDD computes the singular value decomposition (SVD) of a complex
42 *> M-by-N matrix A, optionally computing the left and/or right singular
43 *> vectors, by using divide-and-conquer method. The SVD is written
44 *>
45 *> A = U * SIGMA * conjugate-transpose(V)
46 *>
47 *> where SIGMA is an M-by-N matrix which is zero except for its
48 *> min(m,n) diagonal elements, U is an M-by-M unitary matrix, and
49 *> V is an N-by-N unitary matrix. The diagonal elements of SIGMA
50 *> are the singular values of A; they are real and non-negative, and
51 *> are returned in descending order. The first min(m,n) columns of
52 *> U and V are the left and right singular vectors of A.
53 *>
54 *> Note that the routine returns VT = V**H, not V.
55 *>
56 *> The divide and conquer algorithm makes very mild assumptions about
57 *> floating point arithmetic. It will work on machines with a guard
58 *> digit in add/subtract, or on those binary machines without guard
59 *> digits which subtract like the Cray X-MP, Cray Y-MP, Cray C-90, or
60 *> Cray-2. It could conceivably fail on hexadecimal or decimal machines
61 *> without guard digits, but we know of none.
62 *> \endverbatim
63 *
64 * Arguments:
65 * ==========
66 *
67 *> \param[in] JOBZ
68 *> \verbatim
69 *> JOBZ is CHARACTER*1
70 *> Specifies options for computing all or part of the matrix U:
71 *> = 'A': all M columns of U and all N rows of V**H are
72 *> returned in the arrays U and VT;
73 *> = 'S': the first min(M,N) columns of U and the first
74 *> min(M,N) rows of V**H are returned in the arrays U
75 *> and VT;
76 *> = 'O': If M >= N, the first N columns of U are overwritten
77 *> in the array A and all rows of V**H are returned in
78 *> the array VT;
79 *> otherwise, all columns of U are returned in the
80 *> array U and the first M rows of V**H are overwritten
81 *> in the array A;
82 *> = 'N': no columns of U or rows of V**H are computed.
83 *> \endverbatim
84 *>
85 *> \param[in] M
86 *> \verbatim
87 *> M is INTEGER
88 *> The number of rows of the input matrix A. M >= 0.
89 *> \endverbatim
90 *>
91 *> \param[in] N
92 *> \verbatim
93 *> N is INTEGER
94 *> The number of columns of the input matrix A. N >= 0.
95 *> \endverbatim
96 *>
97 *> \param[in,out] A
98 *> \verbatim
99 *> A is COMPLEX*16 array, dimension (LDA,N)
100 *> On entry, the M-by-N matrix A.
101 *> On exit,
102 *> if JOBZ = 'O', A is overwritten with the first N columns
103 *> of U (the left singular vectors, stored
104 *> columnwise) if M >= N;
105 *> A is overwritten with the first M rows
106 *> of V**H (the right singular vectors, stored
107 *> rowwise) otherwise.
108 *> if JOBZ .ne. 'O', the contents of A are destroyed.
109 *> \endverbatim
110 *>
111 *> \param[in] LDA
112 *> \verbatim
113 *> LDA is INTEGER
114 *> The leading dimension of the array A. LDA >= max(1,M).
115 *> \endverbatim
116 *>
117 *> \param[out] S
118 *> \verbatim
119 *> S is DOUBLE PRECISION array, dimension (min(M,N))
120 *> The singular values of A, sorted so that S(i) >= S(i+1).
121 *> \endverbatim
122 *>
123 *> \param[out] U
124 *> \verbatim
125 *> U is COMPLEX*16 array, dimension (LDU,UCOL)
126 *> UCOL = M if JOBZ = 'A' or JOBZ = 'O' and M < N;
127 *> UCOL = min(M,N) if JOBZ = 'S'.
128 *> If JOBZ = 'A' or JOBZ = 'O' and M < N, U contains the M-by-M
129 *> unitary matrix U;
130 *> if JOBZ = 'S', U contains the first min(M,N) columns of U
131 *> (the left singular vectors, stored columnwise);
132 *> if JOBZ = 'O' and M >= N, or JOBZ = 'N', U is not referenced.
133 *> \endverbatim
134 *>
135 *> \param[in] LDU
136 *> \verbatim
137 *> LDU is INTEGER
138 *> The leading dimension of the array U. LDU >= 1;
139 *> if JOBZ = 'S' or 'A' or JOBZ = 'O' and M < N, LDU >= M.
140 *> \endverbatim
141 *>
142 *> \param[out] VT
143 *> \verbatim
144 *> VT is COMPLEX*16 array, dimension (LDVT,N)
145 *> If JOBZ = 'A' or JOBZ = 'O' and M >= N, VT contains the
146 *> N-by-N unitary matrix V**H;
147 *> if JOBZ = 'S', VT contains the first min(M,N) rows of
148 *> V**H (the right singular vectors, stored rowwise);
149 *> if JOBZ = 'O' and M < N, or JOBZ = 'N', VT is not referenced.
150 *> \endverbatim
151 *>
152 *> \param[in] LDVT
153 *> \verbatim
154 *> LDVT is INTEGER
155 *> The leading dimension of the array VT. LDVT >= 1;
156 *> if JOBZ = 'A' or JOBZ = 'O' and M >= N, LDVT >= N;
157 *> if JOBZ = 'S', LDVT >= min(M,N).
158 *> \endverbatim
159 *>
160 *> \param[out] WORK
161 *> \verbatim
162 *> WORK is COMPLEX*16 array, dimension (MAX(1,LWORK))
163 *> On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
164 *> \endverbatim
165 *>
166 *> \param[in] LWORK
167 *> \verbatim
168 *> LWORK is INTEGER
169 *> The dimension of the array WORK. LWORK >= 1.
170 *> If LWORK = -1, a workspace query is assumed. The optimal
171 *> size for the WORK array is calculated and stored in WORK(1),
172 *> and no other work except argument checking is performed.
173 *>
174 *> Let mx = max(M,N) and mn = min(M,N).
175 *> If JOBZ = 'N', LWORK >= 2*mn + mx.
176 *> If JOBZ = 'O', LWORK >= 2*mn*mn + 2*mn + mx.
177 *> If JOBZ = 'S', LWORK >= mn*mn + 3*mn.
178 *> If JOBZ = 'A', LWORK >= mn*mn + 2*mn + mx.
179 *> These are not tight minimums in all cases; see comments inside code.
180 *> For good performance, LWORK should generally be larger;
181 *> a query is recommended.
182 *> \endverbatim
183 *>
184 *> \param[out] RWORK
185 *> \verbatim
186 *> RWORK is DOUBLE PRECISION array, dimension (MAX(1,LRWORK))
187 *> Let mx = max(M,N) and mn = min(M,N).
188 *> If JOBZ = 'N', LRWORK >= 5*mn (LAPACK <= 3.6 needs 7*mn);
189 *> else if mx >> mn, LRWORK >= 5*mn*mn + 5*mn;
190 *> else LRWORK >= max( 5*mn*mn + 5*mn,
191 *> 2*mx*mn + 2*mn*mn + mn ).
192 *> \endverbatim
193 *>
194 *> \param[out] IWORK
195 *> \verbatim
196 *> IWORK is INTEGER array, dimension (8*min(M,N))
197 *> \endverbatim
198 *>
199 *> \param[out] INFO
200 *> \verbatim
201 *> INFO is INTEGER
202 *> = 0: successful exit.
203 *> < 0: if INFO = -i, the i-th argument had an illegal value.
204 *> > 0: The updating process of DBDSDC did not converge.
205 *> \endverbatim
206 *
207 * Authors:
208 * ========
209 *
210 *> \author Univ. of Tennessee
211 *> \author Univ. of California Berkeley
212 *> \author Univ. of Colorado Denver
213 *> \author NAG Ltd.
214 *
215 *> \ingroup complex16GEsing
216 *
217 *> \par Contributors:
218 * ==================
219 *>
220 *> Ming Gu and Huan Ren, Computer Science Division, University of
221 *> California at Berkeley, USA
222 *>
223 * =====================================================================
224  SUBROUTINE zgesdd( JOBZ, M, N, A, LDA, S, U, LDU, VT, LDVT,
225  $ WORK, LWORK, RWORK, IWORK, INFO )
226  implicit none
227 *
228 * -- LAPACK driver routine --
229 * -- LAPACK is a software package provided by Univ. of Tennessee, --
230 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
231 *
232 * .. Scalar Arguments ..
233  CHARACTER JOBZ
234  INTEGER INFO, LDA, LDU, LDVT, LWORK, M, N
235 * ..
236 * .. Array Arguments ..
237  INTEGER IWORK( * )
238  DOUBLE PRECISION RWORK( * ), S( * )
239  COMPLEX*16 A( LDA, * ), U( LDU, * ), VT( LDVT, * ),
240  $ work( * )
241 * ..
242 *
243 * =====================================================================
244 *
245 * .. Parameters ..
246  COMPLEX*16 CZERO, CONE
247  parameter( czero = ( 0.0d+0, 0.0d+0 ),
248  $ cone = ( 1.0d+0, 0.0d+0 ) )
249  DOUBLE PRECISION ZERO, ONE
250  parameter( zero = 0.0d+0, one = 1.0d+0 )
251 * ..
252 * .. Local Scalars ..
253  LOGICAL LQUERY, WNTQA, WNTQAS, WNTQN, WNTQO, WNTQS
254  INTEGER BLK, CHUNK, I, IE, IERR, IL, IR, IRU, IRVT,
255  $ iscl, itau, itaup, itauq, iu, ivt, ldwkvt,
256  $ ldwrkl, ldwrkr, ldwrku, maxwrk, minmn, minwrk,
257  $ mnthr1, mnthr2, nrwork, nwork, wrkbl
258  INTEGER LWORK_ZGEBRD_MN, LWORK_ZGEBRD_MM,
259  $ lwork_zgebrd_nn, lwork_zgelqf_mn,
260  $ lwork_zgeqrf_mn,
261  $ lwork_zungbr_p_mn, lwork_zungbr_p_nn,
262  $ lwork_zungbr_q_mn, lwork_zungbr_q_mm,
263  $ lwork_zunglq_mn, lwork_zunglq_nn,
264  $ lwork_zungqr_mm, lwork_zungqr_mn,
265  $ lwork_zunmbr_prc_mm, lwork_zunmbr_qln_mm,
266  $ lwork_zunmbr_prc_mn, lwork_zunmbr_qln_mn,
267  $ lwork_zunmbr_prc_nn, lwork_zunmbr_qln_nn
268  DOUBLE PRECISION ANRM, BIGNUM, EPS, SMLNUM
269 * ..
270 * .. Local Arrays ..
271  INTEGER IDUM( 1 )
272  DOUBLE PRECISION DUM( 1 )
273  COMPLEX*16 CDUM( 1 )
274 * ..
275 * .. External Subroutines ..
276  EXTERNAL dbdsdc, dlascl, xerbla, zgebrd, zgelqf, zgemm,
279 * ..
280 * .. External Functions ..
281  LOGICAL LSAME, DISNAN
282  DOUBLE PRECISION DLAMCH, ZLANGE
283  EXTERNAL lsame, dlamch, zlange, disnan
284 * ..
285 * .. Intrinsic Functions ..
286  INTRINSIC int, max, min, sqrt
287 * ..
288 * .. Executable Statements ..
289 *
290 * Test the input arguments
291 *
292  info = 0
293  minmn = min( m, n )
294  mnthr1 = int( minmn*17.0d0 / 9.0d0 )
295  mnthr2 = int( minmn*5.0d0 / 3.0d0 )
296  wntqa = lsame( jobz, 'A' )
297  wntqs = lsame( jobz, 'S' )
298  wntqas = wntqa .OR. wntqs
299  wntqo = lsame( jobz, 'O' )
300  wntqn = lsame( jobz, 'N' )
301  lquery = ( lwork.EQ.-1 )
302  minwrk = 1
303  maxwrk = 1
304 *
305  IF( .NOT.( wntqa .OR. wntqs .OR. wntqo .OR. wntqn ) ) THEN
306  info = -1
307  ELSE IF( m.LT.0 ) THEN
308  info = -2
309  ELSE IF( n.LT.0 ) THEN
310  info = -3
311  ELSE IF( lda.LT.max( 1, m ) ) THEN
312  info = -5
313  ELSE IF( ldu.LT.1 .OR. ( wntqas .AND. ldu.LT.m ) .OR.
314  $ ( wntqo .AND. m.LT.n .AND. ldu.LT.m ) ) THEN
315  info = -8
316  ELSE IF( ldvt.LT.1 .OR. ( wntqa .AND. ldvt.LT.n ) .OR.
317  $ ( wntqs .AND. ldvt.LT.minmn ) .OR.
318  $ ( wntqo .AND. m.GE.n .AND. ldvt.LT.n ) ) THEN
319  info = -10
320  END IF
321 *
322 * Compute workspace
323 * Note: Comments in the code beginning "Workspace:" describe the
324 * minimal amount of workspace allocated at that point in the code,
325 * as well as the preferred amount for good performance.
326 * CWorkspace refers to complex workspace, and RWorkspace to
327 * real workspace. NB refers to the optimal block size for the
328 * immediately following subroutine, as returned by ILAENV.)
329 *
330  IF( info.EQ.0 ) THEN
331  minwrk = 1
332  maxwrk = 1
333  IF( m.GE.n .AND. minmn.GT.0 ) THEN
334 *
335 * There is no complex work space needed for bidiagonal SVD
336 * The real work space needed for bidiagonal SVD (dbdsdc) is
337 * BDSPAC = 3*N*N + 4*N for singular values and vectors;
338 * BDSPAC = 4*N for singular values only;
339 * not including e, RU, and RVT matrices.
340 *
341 * Compute space preferred for each routine
342  CALL zgebrd( m, n, cdum(1), m, dum(1), dum(1), cdum(1),
343  $ cdum(1), cdum(1), -1, ierr )
344  lwork_zgebrd_mn = int( cdum(1) )
345 *
346  CALL zgebrd( n, n, cdum(1), n, dum(1), dum(1), cdum(1),
347  $ cdum(1), cdum(1), -1, ierr )
348  lwork_zgebrd_nn = int( cdum(1) )
349 *
350  CALL zgeqrf( m, n, cdum(1), m, cdum(1), cdum(1), -1, ierr )
351  lwork_zgeqrf_mn = int( cdum(1) )
352 *
353  CALL zungbr( 'P', n, n, n, cdum(1), n, cdum(1), cdum(1),
354  $ -1, ierr )
355  lwork_zungbr_p_nn = int( cdum(1) )
356 *
357  CALL zungbr( 'Q', m, m, n, cdum(1), m, cdum(1), cdum(1),
358  $ -1, ierr )
359  lwork_zungbr_q_mm = int( cdum(1) )
360 *
361  CALL zungbr( 'Q', m, n, n, cdum(1), m, cdum(1), cdum(1),
362  $ -1, ierr )
363  lwork_zungbr_q_mn = int( cdum(1) )
364 *
365  CALL zungqr( m, m, n, cdum(1), m, cdum(1), cdum(1),
366  $ -1, ierr )
367  lwork_zungqr_mm = int( cdum(1) )
368 *
369  CALL zungqr( m, n, n, cdum(1), m, cdum(1), cdum(1),
370  $ -1, ierr )
371  lwork_zungqr_mn = int( cdum(1) )
372 *
373  CALL zunmbr( 'P', 'R', 'C', n, n, n, cdum(1), n, cdum(1),
374  $ cdum(1), n, cdum(1), -1, ierr )
375  lwork_zunmbr_prc_nn = int( cdum(1) )
376 *
377  CALL zunmbr( 'Q', 'L', 'N', m, m, n, cdum(1), m, cdum(1),
378  $ cdum(1), m, cdum(1), -1, ierr )
379  lwork_zunmbr_qln_mm = int( cdum(1) )
380 *
381  CALL zunmbr( 'Q', 'L', 'N', m, n, n, cdum(1), m, cdum(1),
382  $ cdum(1), m, cdum(1), -1, ierr )
383  lwork_zunmbr_qln_mn = int( cdum(1) )
384 *
385  CALL zunmbr( 'Q', 'L', 'N', n, n, n, cdum(1), n, cdum(1),
386  $ cdum(1), n, cdum(1), -1, ierr )
387  lwork_zunmbr_qln_nn = int( cdum(1) )
388 *
389  IF( m.GE.mnthr1 ) THEN
390  IF( wntqn ) THEN
391 *
392 * Path 1 (M >> N, JOBZ='N')
393 *
394  maxwrk = n + lwork_zgeqrf_mn
395  maxwrk = max( maxwrk, 2*n + lwork_zgebrd_nn )
396  minwrk = 3*n
397  ELSE IF( wntqo ) THEN
398 *
399 * Path 2 (M >> N, JOBZ='O')
400 *
401  wrkbl = n + lwork_zgeqrf_mn
402  wrkbl = max( wrkbl, n + lwork_zungqr_mn )
403  wrkbl = max( wrkbl, 2*n + lwork_zgebrd_nn )
404  wrkbl = max( wrkbl, 2*n + lwork_zunmbr_qln_nn )
405  wrkbl = max( wrkbl, 2*n + lwork_zunmbr_prc_nn )
406  maxwrk = m*n + n*n + wrkbl
407  minwrk = 2*n*n + 3*n
408  ELSE IF( wntqs ) THEN
409 *
410 * Path 3 (M >> N, JOBZ='S')
411 *
412  wrkbl = n + lwork_zgeqrf_mn
413  wrkbl = max( wrkbl, n + lwork_zungqr_mn )
414  wrkbl = max( wrkbl, 2*n + lwork_zgebrd_nn )
415  wrkbl = max( wrkbl, 2*n + lwork_zunmbr_qln_nn )
416  wrkbl = max( wrkbl, 2*n + lwork_zunmbr_prc_nn )
417  maxwrk = n*n + wrkbl
418  minwrk = n*n + 3*n
419  ELSE IF( wntqa ) THEN
420 *
421 * Path 4 (M >> N, JOBZ='A')
422 *
423  wrkbl = n + lwork_zgeqrf_mn
424  wrkbl = max( wrkbl, n + lwork_zungqr_mm )
425  wrkbl = max( wrkbl, 2*n + lwork_zgebrd_nn )
426  wrkbl = max( wrkbl, 2*n + lwork_zunmbr_qln_nn )
427  wrkbl = max( wrkbl, 2*n + lwork_zunmbr_prc_nn )
428  maxwrk = n*n + wrkbl
429  minwrk = n*n + max( 3*n, n + m )
430  END IF
431  ELSE IF( m.GE.mnthr2 ) THEN
432 *
433 * Path 5 (M >> N, but not as much as MNTHR1)
434 *
435  maxwrk = 2*n + lwork_zgebrd_mn
436  minwrk = 2*n + m
437  IF( wntqo ) THEN
438 * Path 5o (M >> N, JOBZ='O')
439  maxwrk = max( maxwrk, 2*n + lwork_zungbr_p_nn )
440  maxwrk = max( maxwrk, 2*n + lwork_zungbr_q_mn )
441  maxwrk = maxwrk + m*n
442  minwrk = minwrk + n*n
443  ELSE IF( wntqs ) THEN
444 * Path 5s (M >> N, JOBZ='S')
445  maxwrk = max( maxwrk, 2*n + lwork_zungbr_p_nn )
446  maxwrk = max( maxwrk, 2*n + lwork_zungbr_q_mn )
447  ELSE IF( wntqa ) THEN
448 * Path 5a (M >> N, JOBZ='A')
449  maxwrk = max( maxwrk, 2*n + lwork_zungbr_p_nn )
450  maxwrk = max( maxwrk, 2*n + lwork_zungbr_q_mm )
451  END IF
452  ELSE
453 *
454 * Path 6 (M >= N, but not much larger)
455 *
456  maxwrk = 2*n + lwork_zgebrd_mn
457  minwrk = 2*n + m
458  IF( wntqo ) THEN
459 * Path 6o (M >= N, JOBZ='O')
460  maxwrk = max( maxwrk, 2*n + lwork_zunmbr_prc_nn )
461  maxwrk = max( maxwrk, 2*n + lwork_zunmbr_qln_mn )
462  maxwrk = maxwrk + m*n
463  minwrk = minwrk + n*n
464  ELSE IF( wntqs ) THEN
465 * Path 6s (M >= N, JOBZ='S')
466  maxwrk = max( maxwrk, 2*n + lwork_zunmbr_qln_mn )
467  maxwrk = max( maxwrk, 2*n + lwork_zunmbr_prc_nn )
468  ELSE IF( wntqa ) THEN
469 * Path 6a (M >= N, JOBZ='A')
470  maxwrk = max( maxwrk, 2*n + lwork_zunmbr_qln_mm )
471  maxwrk = max( maxwrk, 2*n + lwork_zunmbr_prc_nn )
472  END IF
473  END IF
474  ELSE IF( minmn.GT.0 ) THEN
475 *
476 * There is no complex work space needed for bidiagonal SVD
477 * The real work space needed for bidiagonal SVD (dbdsdc) is
478 * BDSPAC = 3*M*M + 4*M for singular values and vectors;
479 * BDSPAC = 4*M for singular values only;
480 * not including e, RU, and RVT matrices.
481 *
482 * Compute space preferred for each routine
483  CALL zgebrd( m, n, cdum(1), m, dum(1), dum(1), cdum(1),
484  $ cdum(1), cdum(1), -1, ierr )
485  lwork_zgebrd_mn = int( cdum(1) )
486 *
487  CALL zgebrd( m, m, cdum(1), m, dum(1), dum(1), cdum(1),
488  $ cdum(1), cdum(1), -1, ierr )
489  lwork_zgebrd_mm = int( cdum(1) )
490 *
491  CALL zgelqf( m, n, cdum(1), m, cdum(1), cdum(1), -1, ierr )
492  lwork_zgelqf_mn = int( cdum(1) )
493 *
494  CALL zungbr( 'P', m, n, m, cdum(1), m, cdum(1), cdum(1),
495  $ -1, ierr )
496  lwork_zungbr_p_mn = int( cdum(1) )
497 *
498  CALL zungbr( 'P', n, n, m, cdum(1), n, cdum(1), cdum(1),
499  $ -1, ierr )
500  lwork_zungbr_p_nn = int( cdum(1) )
501 *
502  CALL zungbr( 'Q', m, m, n, cdum(1), m, cdum(1), cdum(1),
503  $ -1, ierr )
504  lwork_zungbr_q_mm = int( cdum(1) )
505 *
506  CALL zunglq( m, n, m, cdum(1), m, cdum(1), cdum(1),
507  $ -1, ierr )
508  lwork_zunglq_mn = int( cdum(1) )
509 *
510  CALL zunglq( n, n, m, cdum(1), n, cdum(1), cdum(1),
511  $ -1, ierr )
512  lwork_zunglq_nn = int( cdum(1) )
513 *
514  CALL zunmbr( 'P', 'R', 'C', m, m, m, cdum(1), m, cdum(1),
515  $ cdum(1), m, cdum(1), -1, ierr )
516  lwork_zunmbr_prc_mm = int( cdum(1) )
517 *
518  CALL zunmbr( 'P', 'R', 'C', m, n, m, cdum(1), m, cdum(1),
519  $ cdum(1), m, cdum(1), -1, ierr )
520  lwork_zunmbr_prc_mn = int( cdum(1) )
521 *
522  CALL zunmbr( 'P', 'R', 'C', n, n, m, cdum(1), n, cdum(1),
523  $ cdum(1), n, cdum(1), -1, ierr )
524  lwork_zunmbr_prc_nn = int( cdum(1) )
525 *
526  CALL zunmbr( 'Q', 'L', 'N', m, m, m, cdum(1), m, cdum(1),
527  $ cdum(1), m, cdum(1), -1, ierr )
528  lwork_zunmbr_qln_mm = int( cdum(1) )
529 *
530  IF( n.GE.mnthr1 ) THEN
531  IF( wntqn ) THEN
532 *
533 * Path 1t (N >> M, JOBZ='N')
534 *
535  maxwrk = m + lwork_zgelqf_mn
536  maxwrk = max( maxwrk, 2*m + lwork_zgebrd_mm )
537  minwrk = 3*m
538  ELSE IF( wntqo ) THEN
539 *
540 * Path 2t (N >> M, JOBZ='O')
541 *
542  wrkbl = m + lwork_zgelqf_mn
543  wrkbl = max( wrkbl, m + lwork_zunglq_mn )
544  wrkbl = max( wrkbl, 2*m + lwork_zgebrd_mm )
545  wrkbl = max( wrkbl, 2*m + lwork_zunmbr_qln_mm )
546  wrkbl = max( wrkbl, 2*m + lwork_zunmbr_prc_mm )
547  maxwrk = m*n + m*m + wrkbl
548  minwrk = 2*m*m + 3*m
549  ELSE IF( wntqs ) THEN
550 *
551 * Path 3t (N >> M, JOBZ='S')
552 *
553  wrkbl = m + lwork_zgelqf_mn
554  wrkbl = max( wrkbl, m + lwork_zunglq_mn )
555  wrkbl = max( wrkbl, 2*m + lwork_zgebrd_mm )
556  wrkbl = max( wrkbl, 2*m + lwork_zunmbr_qln_mm )
557  wrkbl = max( wrkbl, 2*m + lwork_zunmbr_prc_mm )
558  maxwrk = m*m + wrkbl
559  minwrk = m*m + 3*m
560  ELSE IF( wntqa ) THEN
561 *
562 * Path 4t (N >> M, JOBZ='A')
563 *
564  wrkbl = m + lwork_zgelqf_mn
565  wrkbl = max( wrkbl, m + lwork_zunglq_nn )
566  wrkbl = max( wrkbl, 2*m + lwork_zgebrd_mm )
567  wrkbl = max( wrkbl, 2*m + lwork_zunmbr_qln_mm )
568  wrkbl = max( wrkbl, 2*m + lwork_zunmbr_prc_mm )
569  maxwrk = m*m + wrkbl
570  minwrk = m*m + max( 3*m, m + n )
571  END IF
572  ELSE IF( n.GE.mnthr2 ) THEN
573 *
574 * Path 5t (N >> M, but not as much as MNTHR1)
575 *
576  maxwrk = 2*m + lwork_zgebrd_mn
577  minwrk = 2*m + n
578  IF( wntqo ) THEN
579 * Path 5to (N >> M, JOBZ='O')
580  maxwrk = max( maxwrk, 2*m + lwork_zungbr_q_mm )
581  maxwrk = max( maxwrk, 2*m + lwork_zungbr_p_mn )
582  maxwrk = maxwrk + m*n
583  minwrk = minwrk + m*m
584  ELSE IF( wntqs ) THEN
585 * Path 5ts (N >> M, JOBZ='S')
586  maxwrk = max( maxwrk, 2*m + lwork_zungbr_q_mm )
587  maxwrk = max( maxwrk, 2*m + lwork_zungbr_p_mn )
588  ELSE IF( wntqa ) THEN
589 * Path 5ta (N >> M, JOBZ='A')
590  maxwrk = max( maxwrk, 2*m + lwork_zungbr_q_mm )
591  maxwrk = max( maxwrk, 2*m + lwork_zungbr_p_nn )
592  END IF
593  ELSE
594 *
595 * Path 6t (N > M, but not much larger)
596 *
597  maxwrk = 2*m + lwork_zgebrd_mn
598  minwrk = 2*m + n
599  IF( wntqo ) THEN
600 * Path 6to (N > M, JOBZ='O')
601  maxwrk = max( maxwrk, 2*m + lwork_zunmbr_qln_mm )
602  maxwrk = max( maxwrk, 2*m + lwork_zunmbr_prc_mn )
603  maxwrk = maxwrk + m*n
604  minwrk = minwrk + m*m
605  ELSE IF( wntqs ) THEN
606 * Path 6ts (N > M, JOBZ='S')
607  maxwrk = max( maxwrk, 2*m + lwork_zunmbr_qln_mm )
608  maxwrk = max( maxwrk, 2*m + lwork_zunmbr_prc_mn )
609  ELSE IF( wntqa ) THEN
610 * Path 6ta (N > M, JOBZ='A')
611  maxwrk = max( maxwrk, 2*m + lwork_zunmbr_qln_mm )
612  maxwrk = max( maxwrk, 2*m + lwork_zunmbr_prc_nn )
613  END IF
614  END IF
615  END IF
616  maxwrk = max( maxwrk, minwrk )
617  END IF
618  IF( info.EQ.0 ) THEN
619  work( 1 ) = maxwrk
620  IF( lwork.LT.minwrk .AND. .NOT. lquery ) THEN
621  info = -12
622  END IF
623  END IF
624 *
625  IF( info.NE.0 ) THEN
626  CALL xerbla( 'ZGESDD', -info )
627  RETURN
628  ELSE IF( lquery ) THEN
629  RETURN
630  END IF
631 *
632 * Quick return if possible
633 *
634  IF( m.EQ.0 .OR. n.EQ.0 ) THEN
635  RETURN
636  END IF
637 *
638 * Get machine constants
639 *
640  eps = dlamch( 'P' )
641  smlnum = sqrt( dlamch( 'S' ) ) / eps
642  bignum = one / smlnum
643 *
644 * Scale A if max element outside range [SMLNUM,BIGNUM]
645 *
646  anrm = zlange( 'M', m, n, a, lda, dum )
647  IF( disnan( anrm ) ) THEN
648  info = -4
649  RETURN
650  END IF
651  iscl = 0
652  IF( anrm.GT.zero .AND. anrm.LT.smlnum ) THEN
653  iscl = 1
654  CALL zlascl( 'G', 0, 0, anrm, smlnum, m, n, a, lda, ierr )
655  ELSE IF( anrm.GT.bignum ) THEN
656  iscl = 1
657  CALL zlascl( 'G', 0, 0, anrm, bignum, m, n, a, lda, ierr )
658  END IF
659 *
660  IF( m.GE.n ) THEN
661 *
662 * A has at least as many rows as columns. If A has sufficiently
663 * more rows than columns, first reduce using the QR
664 * decomposition (if sufficient workspace available)
665 *
666  IF( m.GE.mnthr1 ) THEN
667 *
668  IF( wntqn ) THEN
669 *
670 * Path 1 (M >> N, JOBZ='N')
671 * No singular vectors to be computed
672 *
673  itau = 1
674  nwork = itau + n
675 *
676 * Compute A=Q*R
677 * CWorkspace: need N [tau] + N [work]
678 * CWorkspace: prefer N [tau] + N*NB [work]
679 * RWorkspace: need 0
680 *
681  CALL zgeqrf( m, n, a, lda, work( itau ), work( nwork ),
682  $ lwork-nwork+1, ierr )
683 *
684 * Zero out below R
685 *
686  CALL zlaset( 'L', n-1, n-1, czero, czero, a( 2, 1 ),
687  $ lda )
688  ie = 1
689  itauq = 1
690  itaup = itauq + n
691  nwork = itaup + n
692 *
693 * Bidiagonalize R in A
694 * CWorkspace: need 2*N [tauq, taup] + N [work]
695 * CWorkspace: prefer 2*N [tauq, taup] + 2*N*NB [work]
696 * RWorkspace: need N [e]
697 *
698  CALL zgebrd( n, n, a, lda, s, rwork( ie ), work( itauq ),
699  $ work( itaup ), work( nwork ), lwork-nwork+1,
700  $ ierr )
701  nrwork = ie + n
702 *
703 * Perform bidiagonal SVD, compute singular values only
704 * CWorkspace: need 0
705 * RWorkspace: need N [e] + BDSPAC
706 *
707  CALL dbdsdc( 'U', 'N', n, s, rwork( ie ), dum,1,dum,1,
708  $ dum, idum, rwork( nrwork ), iwork, info )
709 *
710  ELSE IF( wntqo ) THEN
711 *
712 * Path 2 (M >> N, JOBZ='O')
713 * N left singular vectors to be overwritten on A and
714 * N right singular vectors to be computed in VT
715 *
716  iu = 1
717 *
718 * WORK(IU) is N by N
719 *
720  ldwrku = n
721  ir = iu + ldwrku*n
722  IF( lwork .GE. m*n + n*n + 3*n ) THEN
723 *
724 * WORK(IR) is M by N
725 *
726  ldwrkr = m
727  ELSE
728  ldwrkr = ( lwork - n*n - 3*n ) / n
729  END IF
730  itau = ir + ldwrkr*n
731  nwork = itau + n
732 *
733 * Compute A=Q*R
734 * CWorkspace: need N*N [U] + N*N [R] + N [tau] + N [work]
735 * CWorkspace: prefer N*N [U] + N*N [R] + N [tau] + N*NB [work]
736 * RWorkspace: need 0
737 *
738  CALL zgeqrf( m, n, a, lda, work( itau ), work( nwork ),
739  $ lwork-nwork+1, ierr )
740 *
741 * Copy R to WORK( IR ), zeroing out below it
742 *
743  CALL zlacpy( 'U', n, n, a, lda, work( ir ), ldwrkr )
744  CALL zlaset( 'L', n-1, n-1, czero, czero, work( ir+1 ),
745  $ ldwrkr )
746 *
747 * Generate Q in A
748 * CWorkspace: need N*N [U] + N*N [R] + N [tau] + N [work]
749 * CWorkspace: prefer N*N [U] + N*N [R] + N [tau] + N*NB [work]
750 * RWorkspace: need 0
751 *
752  CALL zungqr( m, n, n, a, lda, work( itau ),
753  $ work( nwork ), lwork-nwork+1, ierr )
754  ie = 1
755  itauq = itau
756  itaup = itauq + n
757  nwork = itaup + n
758 *
759 * Bidiagonalize R in WORK(IR)
760 * CWorkspace: need N*N [U] + N*N [R] + 2*N [tauq, taup] + N [work]
761 * CWorkspace: prefer N*N [U] + N*N [R] + 2*N [tauq, taup] + 2*N*NB [work]
762 * RWorkspace: need N [e]
763 *
764  CALL zgebrd( n, n, work( ir ), ldwrkr, s, rwork( ie ),
765  $ work( itauq ), work( itaup ), work( nwork ),
766  $ lwork-nwork+1, ierr )
767 *
768 * Perform bidiagonal SVD, computing left singular vectors
769 * of R in WORK(IRU) and computing right singular vectors
770 * of R in WORK(IRVT)
771 * CWorkspace: need 0
772 * RWorkspace: need N [e] + N*N [RU] + N*N [RVT] + BDSPAC
773 *
774  iru = ie + n
775  irvt = iru + n*n
776  nrwork = irvt + n*n
777  CALL dbdsdc( 'U', 'I', n, s, rwork( ie ), rwork( iru ),
778  $ n, rwork( irvt ), n, dum, idum,
779  $ rwork( nrwork ), iwork, info )
780 *
781 * Copy real matrix RWORK(IRU) to complex matrix WORK(IU)
782 * Overwrite WORK(IU) by the left singular vectors of R
783 * CWorkspace: need N*N [U] + N*N [R] + 2*N [tauq, taup] + N [work]
784 * CWorkspace: prefer N*N [U] + N*N [R] + 2*N [tauq, taup] + N*NB [work]
785 * RWorkspace: need 0
786 *
787  CALL zlacp2( 'F', n, n, rwork( iru ), n, work( iu ),
788  $ ldwrku )
789  CALL zunmbr( 'Q', 'L', 'N', n, n, n, work( ir ), ldwrkr,
790  $ work( itauq ), work( iu ), ldwrku,
791  $ work( nwork ), lwork-nwork+1, ierr )
792 *
793 * Copy real matrix RWORK(IRVT) to complex matrix VT
794 * Overwrite VT by the right singular vectors of R
795 * CWorkspace: need N*N [U] + N*N [R] + 2*N [tauq, taup] + N [work]
796 * CWorkspace: prefer N*N [U] + N*N [R] + 2*N [tauq, taup] + N*NB [work]
797 * RWorkspace: need 0
798 *
799  CALL zlacp2( 'F', n, n, rwork( irvt ), n, vt, ldvt )
800  CALL zunmbr( 'P', 'R', 'C', n, n, n, work( ir ), ldwrkr,
801  $ work( itaup ), vt, ldvt, work( nwork ),
802  $ lwork-nwork+1, ierr )
803 *
804 * Multiply Q in A by left singular vectors of R in
805 * WORK(IU), storing result in WORK(IR) and copying to A
806 * CWorkspace: need N*N [U] + N*N [R]
807 * CWorkspace: prefer N*N [U] + M*N [R]
808 * RWorkspace: need 0
809 *
810  DO 10 i = 1, m, ldwrkr
811  chunk = min( m-i+1, ldwrkr )
812  CALL zgemm( 'N', 'N', chunk, n, n, cone, a( i, 1 ),
813  $ lda, work( iu ), ldwrku, czero,
814  $ work( ir ), ldwrkr )
815  CALL zlacpy( 'F', chunk, n, work( ir ), ldwrkr,
816  $ a( i, 1 ), lda )
817  10 CONTINUE
818 *
819  ELSE IF( wntqs ) THEN
820 *
821 * Path 3 (M >> N, JOBZ='S')
822 * N left singular vectors to be computed in U and
823 * N right singular vectors to be computed in VT
824 *
825  ir = 1
826 *
827 * WORK(IR) is N by N
828 *
829  ldwrkr = n
830  itau = ir + ldwrkr*n
831  nwork = itau + n
832 *
833 * Compute A=Q*R
834 * CWorkspace: need N*N [R] + N [tau] + N [work]
835 * CWorkspace: prefer N*N [R] + N [tau] + N*NB [work]
836 * RWorkspace: need 0
837 *
838  CALL zgeqrf( m, n, a, lda, work( itau ), work( nwork ),
839  $ lwork-nwork+1, ierr )
840 *
841 * Copy R to WORK(IR), zeroing out below it
842 *
843  CALL zlacpy( 'U', n, n, a, lda, work( ir ), ldwrkr )
844  CALL zlaset( 'L', n-1, n-1, czero, czero, work( ir+1 ),
845  $ ldwrkr )
846 *
847 * Generate Q in A
848 * CWorkspace: need N*N [R] + N [tau] + N [work]
849 * CWorkspace: prefer N*N [R] + N [tau] + N*NB [work]
850 * RWorkspace: need 0
851 *
852  CALL zungqr( m, n, n, a, lda, work( itau ),
853  $ work( nwork ), lwork-nwork+1, ierr )
854  ie = 1
855  itauq = itau
856  itaup = itauq + n
857  nwork = itaup + n
858 *
859 * Bidiagonalize R in WORK(IR)
860 * CWorkspace: need N*N [R] + 2*N [tauq, taup] + N [work]
861 * CWorkspace: prefer N*N [R] + 2*N [tauq, taup] + 2*N*NB [work]
862 * RWorkspace: need N [e]
863 *
864  CALL zgebrd( n, n, work( ir ), ldwrkr, s, rwork( ie ),
865  $ work( itauq ), work( itaup ), work( nwork ),
866  $ lwork-nwork+1, ierr )
867 *
868 * Perform bidiagonal SVD, computing left singular vectors
869 * of bidiagonal matrix in RWORK(IRU) and computing right
870 * singular vectors of bidiagonal matrix in RWORK(IRVT)
871 * CWorkspace: need 0
872 * RWorkspace: need N [e] + N*N [RU] + N*N [RVT] + BDSPAC
873 *
874  iru = ie + n
875  irvt = iru + n*n
876  nrwork = irvt + n*n
877  CALL dbdsdc( 'U', 'I', n, s, rwork( ie ), rwork( iru ),
878  $ n, rwork( irvt ), n, dum, idum,
879  $ rwork( nrwork ), iwork, info )
880 *
881 * Copy real matrix RWORK(IRU) to complex matrix U
882 * Overwrite U by left singular vectors of R
883 * CWorkspace: need N*N [R] + 2*N [tauq, taup] + N [work]
884 * CWorkspace: prefer N*N [R] + 2*N [tauq, taup] + N*NB [work]
885 * RWorkspace: need 0
886 *
887  CALL zlacp2( 'F', n, n, rwork( iru ), n, u, ldu )
888  CALL zunmbr( 'Q', 'L', 'N', n, n, n, work( ir ), ldwrkr,
889  $ work( itauq ), u, ldu, work( nwork ),
890  $ lwork-nwork+1, ierr )
891 *
892 * Copy real matrix RWORK(IRVT) to complex matrix VT
893 * Overwrite VT by right singular vectors of R
894 * CWorkspace: need N*N [R] + 2*N [tauq, taup] + N [work]
895 * CWorkspace: prefer N*N [R] + 2*N [tauq, taup] + N*NB [work]
896 * RWorkspace: need 0
897 *
898  CALL zlacp2( 'F', n, n, rwork( irvt ), n, vt, ldvt )
899  CALL zunmbr( 'P', 'R', 'C', n, n, n, work( ir ), ldwrkr,
900  $ work( itaup ), vt, ldvt, work( nwork ),
901  $ lwork-nwork+1, ierr )
902 *
903 * Multiply Q in A by left singular vectors of R in
904 * WORK(IR), storing result in U
905 * CWorkspace: need N*N [R]
906 * RWorkspace: need 0
907 *
908  CALL zlacpy( 'F', n, n, u, ldu, work( ir ), ldwrkr )
909  CALL zgemm( 'N', 'N', m, n, n, cone, a, lda, work( ir ),
910  $ ldwrkr, czero, u, ldu )
911 *
912  ELSE IF( wntqa ) THEN
913 *
914 * Path 4 (M >> N, JOBZ='A')
915 * M left singular vectors to be computed in U and
916 * N right singular vectors to be computed in VT
917 *
918  iu = 1
919 *
920 * WORK(IU) is N by N
921 *
922  ldwrku = n
923  itau = iu + ldwrku*n
924  nwork = itau + n
925 *
926 * Compute A=Q*R, copying result to U
927 * CWorkspace: need N*N [U] + N [tau] + N [work]
928 * CWorkspace: prefer N*N [U] + N [tau] + N*NB [work]
929 * RWorkspace: need 0
930 *
931  CALL zgeqrf( m, n, a, lda, work( itau ), work( nwork ),
932  $ lwork-nwork+1, ierr )
933  CALL zlacpy( 'L', m, n, a, lda, u, ldu )
934 *
935 * Generate Q in U
936 * CWorkspace: need N*N [U] + N [tau] + M [work]
937 * CWorkspace: prefer N*N [U] + N [tau] + M*NB [work]
938 * RWorkspace: need 0
939 *
940  CALL zungqr( m, m, n, u, ldu, work( itau ),
941  $ work( nwork ), lwork-nwork+1, ierr )
942 *
943 * Produce R in A, zeroing out below it
944 *
945  CALL zlaset( 'L', n-1, n-1, czero, czero, a( 2, 1 ),
946  $ lda )
947  ie = 1
948  itauq = itau
949  itaup = itauq + n
950  nwork = itaup + n
951 *
952 * Bidiagonalize R in A
953 * CWorkspace: need N*N [U] + 2*N [tauq, taup] + N [work]
954 * CWorkspace: prefer N*N [U] + 2*N [tauq, taup] + 2*N*NB [work]
955 * RWorkspace: need N [e]
956 *
957  CALL zgebrd( n, n, a, lda, s, rwork( ie ), work( itauq ),
958  $ work( itaup ), work( nwork ), lwork-nwork+1,
959  $ ierr )
960  iru = ie + n
961  irvt = iru + n*n
962  nrwork = irvt + n*n
963 *
964 * Perform bidiagonal SVD, computing left singular vectors
965 * of bidiagonal matrix in RWORK(IRU) and computing right
966 * singular vectors of bidiagonal matrix in RWORK(IRVT)
967 * CWorkspace: need 0
968 * RWorkspace: need N [e] + N*N [RU] + N*N [RVT] + BDSPAC
969 *
970  CALL dbdsdc( 'U', 'I', n, s, rwork( ie ), rwork( iru ),
971  $ n, rwork( irvt ), n, dum, idum,
972  $ rwork( nrwork ), iwork, info )
973 *
974 * Copy real matrix RWORK(IRU) to complex matrix WORK(IU)
975 * Overwrite WORK(IU) by left singular vectors of R
976 * CWorkspace: need N*N [U] + 2*N [tauq, taup] + N [work]
977 * CWorkspace: prefer N*N [U] + 2*N [tauq, taup] + N*NB [work]
978 * RWorkspace: need 0
979 *
980  CALL zlacp2( 'F', n, n, rwork( iru ), n, work( iu ),
981  $ ldwrku )
982  CALL zunmbr( 'Q', 'L', 'N', n, n, n, a, lda,
983  $ work( itauq ), work( iu ), ldwrku,
984  $ work( nwork ), lwork-nwork+1, ierr )
985 *
986 * Copy real matrix RWORK(IRVT) to complex matrix VT
987 * Overwrite VT by right singular vectors of R
988 * CWorkspace: need N*N [U] + 2*N [tauq, taup] + N [work]
989 * CWorkspace: prefer N*N [U] + 2*N [tauq, taup] + N*NB [work]
990 * RWorkspace: need 0
991 *
992  CALL zlacp2( 'F', n, n, rwork( irvt ), n, vt, ldvt )
993  CALL zunmbr( 'P', 'R', 'C', n, n, n, a, lda,
994  $ work( itaup ), vt, ldvt, work( nwork ),
995  $ lwork-nwork+1, ierr )
996 *
997 * Multiply Q in U by left singular vectors of R in
998 * WORK(IU), storing result in A
999 * CWorkspace: need N*N [U]
1000 * RWorkspace: need 0
1001 *
1002  CALL zgemm( 'N', 'N', m, n, n, cone, u, ldu, work( iu ),
1003  $ ldwrku, czero, a, lda )
1004 *
1005 * Copy left singular vectors of A from A to U
1006 *
1007  CALL zlacpy( 'F', m, n, a, lda, u, ldu )
1008 *
1009  END IF
1010 *
1011  ELSE IF( m.GE.mnthr2 ) THEN
1012 *
1013 * MNTHR2 <= M < MNTHR1
1014 *
1015 * Path 5 (M >> N, but not as much as MNTHR1)
1016 * Reduce to bidiagonal form without QR decomposition, use
1017 * ZUNGBR and matrix multiplication to compute singular vectors
1018 *
1019  ie = 1
1020  nrwork = ie + n
1021  itauq = 1
1022  itaup = itauq + n
1023  nwork = itaup + n
1024 *
1025 * Bidiagonalize A
1026 * CWorkspace: need 2*N [tauq, taup] + M [work]
1027 * CWorkspace: prefer 2*N [tauq, taup] + (M+N)*NB [work]
1028 * RWorkspace: need N [e]
1029 *
1030  CALL zgebrd( m, n, a, lda, s, rwork( ie ), work( itauq ),
1031  $ work( itaup ), work( nwork ), lwork-nwork+1,
1032  $ ierr )
1033  IF( wntqn ) THEN
1034 *
1035 * Path 5n (M >> N, JOBZ='N')
1036 * Compute singular values only
1037 * CWorkspace: need 0
1038 * RWorkspace: need N [e] + BDSPAC
1039 *
1040  CALL dbdsdc( 'U', 'N', n, s, rwork( ie ), dum, 1,dum,1,
1041  $ dum, idum, rwork( nrwork ), iwork, info )
1042  ELSE IF( wntqo ) THEN
1043  iu = nwork
1044  iru = nrwork
1045  irvt = iru + n*n
1046  nrwork = irvt + n*n
1047 *
1048 * Path 5o (M >> N, JOBZ='O')
1049 * Copy A to VT, generate P**H
1050 * CWorkspace: need 2*N [tauq, taup] + N [work]
1051 * CWorkspace: prefer 2*N [tauq, taup] + N*NB [work]
1052 * RWorkspace: need 0
1053 *
1054  CALL zlacpy( 'U', n, n, a, lda, vt, ldvt )
1055  CALL zungbr( 'P', n, n, n, vt, ldvt, work( itaup ),
1056  $ work( nwork ), lwork-nwork+1, ierr )
1057 *
1058 * Generate Q in A
1059 * CWorkspace: need 2*N [tauq, taup] + N [work]
1060 * CWorkspace: prefer 2*N [tauq, taup] + N*NB [work]
1061 * RWorkspace: need 0
1062 *
1063  CALL zungbr( 'Q', m, n, n, a, lda, work( itauq ),
1064  $ work( nwork ), lwork-nwork+1, ierr )
1065 *
1066  IF( lwork .GE. m*n + 3*n ) THEN
1067 *
1068 * WORK( IU ) is M by N
1069 *
1070  ldwrku = m
1071  ELSE
1072 *
1073 * WORK(IU) is LDWRKU by N
1074 *
1075  ldwrku = ( lwork - 3*n ) / n
1076  END IF
1077  nwork = iu + ldwrku*n
1078 *
1079 * Perform bidiagonal SVD, computing left singular vectors
1080 * of bidiagonal matrix in RWORK(IRU) and computing right
1081 * singular vectors of bidiagonal matrix in RWORK(IRVT)
1082 * CWorkspace: need 0
1083 * RWorkspace: need N [e] + N*N [RU] + N*N [RVT] + BDSPAC
1084 *
1085  CALL dbdsdc( 'U', 'I', n, s, rwork( ie ), rwork( iru ),
1086  $ n, rwork( irvt ), n, dum, idum,
1087  $ rwork( nrwork ), iwork, info )
1088 *
1089 * Multiply real matrix RWORK(IRVT) by P**H in VT,
1090 * storing the result in WORK(IU), copying to VT
1091 * CWorkspace: need 2*N [tauq, taup] + N*N [U]
1092 * RWorkspace: need N [e] + N*N [RU] + N*N [RVT] + 2*N*N [rwork]
1093 *
1094  CALL zlarcm( n, n, rwork( irvt ), n, vt, ldvt,
1095  $ work( iu ), ldwrku, rwork( nrwork ) )
1096  CALL zlacpy( 'F', n, n, work( iu ), ldwrku, vt, ldvt )
1097 *
1098 * Multiply Q in A by real matrix RWORK(IRU), storing the
1099 * result in WORK(IU), copying to A
1100 * CWorkspace: need 2*N [tauq, taup] + N*N [U]
1101 * CWorkspace: prefer 2*N [tauq, taup] + M*N [U]
1102 * RWorkspace: need N [e] + N*N [RU] + 2*N*N [rwork]
1103 * RWorkspace: prefer N [e] + N*N [RU] + 2*M*N [rwork] < N + 5*N*N since M < 2*N here
1104 *
1105  nrwork = irvt
1106  DO 20 i = 1, m, ldwrku
1107  chunk = min( m-i+1, ldwrku )
1108  CALL zlacrm( chunk, n, a( i, 1 ), lda, rwork( iru ),
1109  $ n, work( iu ), ldwrku, rwork( nrwork ) )
1110  CALL zlacpy( 'F', chunk, n, work( iu ), ldwrku,
1111  $ a( i, 1 ), lda )
1112  20 CONTINUE
1113 *
1114  ELSE IF( wntqs ) THEN
1115 *
1116 * Path 5s (M >> N, JOBZ='S')
1117 * Copy A to VT, generate P**H
1118 * CWorkspace: need 2*N [tauq, taup] + N [work]
1119 * CWorkspace: prefer 2*N [tauq, taup] + N*NB [work]
1120 * RWorkspace: need 0
1121 *
1122  CALL zlacpy( 'U', n, n, a, lda, vt, ldvt )
1123  CALL zungbr( 'P', n, n, n, vt, ldvt, work( itaup ),
1124  $ work( nwork ), lwork-nwork+1, ierr )
1125 *
1126 * Copy A to U, generate Q
1127 * CWorkspace: need 2*N [tauq, taup] + N [work]
1128 * CWorkspace: prefer 2*N [tauq, taup] + N*NB [work]
1129 * RWorkspace: need 0
1130 *
1131  CALL zlacpy( 'L', m, n, a, lda, u, ldu )
1132  CALL zungbr( 'Q', m, n, n, u, ldu, work( itauq ),
1133  $ work( nwork ), lwork-nwork+1, ierr )
1134 *
1135 * Perform bidiagonal SVD, computing left singular vectors
1136 * of bidiagonal matrix in RWORK(IRU) and computing right
1137 * singular vectors of bidiagonal matrix in RWORK(IRVT)
1138 * CWorkspace: need 0
1139 * RWorkspace: need N [e] + N*N [RU] + N*N [RVT] + BDSPAC
1140 *
1141  iru = nrwork
1142  irvt = iru + n*n
1143  nrwork = irvt + n*n
1144  CALL dbdsdc( 'U', 'I', n, s, rwork( ie ), rwork( iru ),
1145  $ n, rwork( irvt ), n, dum, idum,
1146  $ rwork( nrwork ), iwork, info )
1147 *
1148 * Multiply real matrix RWORK(IRVT) by P**H in VT,
1149 * storing the result in A, copying to VT
1150 * CWorkspace: need 0
1151 * RWorkspace: need N [e] + N*N [RU] + N*N [RVT] + 2*N*N [rwork]
1152 *
1153  CALL zlarcm( n, n, rwork( irvt ), n, vt, ldvt, a, lda,
1154  $ rwork( nrwork ) )
1155  CALL zlacpy( 'F', n, n, a, lda, vt, ldvt )
1156 *
1157 * Multiply Q in U by real matrix RWORK(IRU), storing the
1158 * result in A, copying to U
1159 * CWorkspace: need 0
1160 * RWorkspace: need N [e] + N*N [RU] + 2*M*N [rwork] < N + 5*N*N since M < 2*N here
1161 *
1162  nrwork = irvt
1163  CALL zlacrm( m, n, u, ldu, rwork( iru ), n, a, lda,
1164  $ rwork( nrwork ) )
1165  CALL zlacpy( 'F', m, n, a, lda, u, ldu )
1166  ELSE
1167 *
1168 * Path 5a (M >> N, JOBZ='A')
1169 * Copy A to VT, generate P**H
1170 * CWorkspace: need 2*N [tauq, taup] + N [work]
1171 * CWorkspace: prefer 2*N [tauq, taup] + N*NB [work]
1172 * RWorkspace: need 0
1173 *
1174  CALL zlacpy( 'U', n, n, a, lda, vt, ldvt )
1175  CALL zungbr( 'P', n, n, n, vt, ldvt, work( itaup ),
1176  $ work( nwork ), lwork-nwork+1, ierr )
1177 *
1178 * Copy A to U, generate Q
1179 * CWorkspace: need 2*N [tauq, taup] + M [work]
1180 * CWorkspace: prefer 2*N [tauq, taup] + M*NB [work]
1181 * RWorkspace: need 0
1182 *
1183  CALL zlacpy( 'L', m, n, a, lda, u, ldu )
1184  CALL zungbr( 'Q', m, m, n, u, ldu, work( itauq ),
1185  $ work( nwork ), lwork-nwork+1, ierr )
1186 *
1187 * Perform bidiagonal SVD, computing left singular vectors
1188 * of bidiagonal matrix in RWORK(IRU) and computing right
1189 * singular vectors of bidiagonal matrix in RWORK(IRVT)
1190 * CWorkspace: need 0
1191 * RWorkspace: need N [e] + N*N [RU] + N*N [RVT] + BDSPAC
1192 *
1193  iru = nrwork
1194  irvt = iru + n*n
1195  nrwork = irvt + n*n
1196  CALL dbdsdc( 'U', 'I', n, s, rwork( ie ), rwork( iru ),
1197  $ n, rwork( irvt ), n, dum, idum,
1198  $ rwork( nrwork ), iwork, info )
1199 *
1200 * Multiply real matrix RWORK(IRVT) by P**H in VT,
1201 * storing the result in A, copying to VT
1202 * CWorkspace: need 0
1203 * RWorkspace: need N [e] + N*N [RU] + N*N [RVT] + 2*N*N [rwork]
1204 *
1205  CALL zlarcm( n, n, rwork( irvt ), n, vt, ldvt, a, lda,
1206  $ rwork( nrwork ) )
1207  CALL zlacpy( 'F', n, n, a, lda, vt, ldvt )
1208 *
1209 * Multiply Q in U by real matrix RWORK(IRU), storing the
1210 * result in A, copying to U
1211 * CWorkspace: need 0
1212 * RWorkspace: need N [e] + N*N [RU] + 2*M*N [rwork] < N + 5*N*N since M < 2*N here
1213 *
1214  nrwork = irvt
1215  CALL zlacrm( m, n, u, ldu, rwork( iru ), n, a, lda,
1216  $ rwork( nrwork ) )
1217  CALL zlacpy( 'F', m, n, a, lda, u, ldu )
1218  END IF
1219 *
1220  ELSE
1221 *
1222 * M .LT. MNTHR2
1223 *
1224 * Path 6 (M >= N, but not much larger)
1225 * Reduce to bidiagonal form without QR decomposition
1226 * Use ZUNMBR to compute singular vectors
1227 *
1228  ie = 1
1229  nrwork = ie + n
1230  itauq = 1
1231  itaup = itauq + n
1232  nwork = itaup + n
1233 *
1234 * Bidiagonalize A
1235 * CWorkspace: need 2*N [tauq, taup] + M [work]
1236 * CWorkspace: prefer 2*N [tauq, taup] + (M+N)*NB [work]
1237 * RWorkspace: need N [e]
1238 *
1239  CALL zgebrd( m, n, a, lda, s, rwork( ie ), work( itauq ),
1240  $ work( itaup ), work( nwork ), lwork-nwork+1,
1241  $ ierr )
1242  IF( wntqn ) THEN
1243 *
1244 * Path 6n (M >= N, JOBZ='N')
1245 * Compute singular values only
1246 * CWorkspace: need 0
1247 * RWorkspace: need N [e] + BDSPAC
1248 *
1249  CALL dbdsdc( 'U', 'N', n, s, rwork( ie ), dum,1,dum,1,
1250  $ dum, idum, rwork( nrwork ), iwork, info )
1251  ELSE IF( wntqo ) THEN
1252  iu = nwork
1253  iru = nrwork
1254  irvt = iru + n*n
1255  nrwork = irvt + n*n
1256  IF( lwork .GE. m*n + 3*n ) THEN
1257 *
1258 * WORK( IU ) is M by N
1259 *
1260  ldwrku = m
1261  ELSE
1262 *
1263 * WORK( IU ) is LDWRKU by N
1264 *
1265  ldwrku = ( lwork - 3*n ) / n
1266  END IF
1267  nwork = iu + ldwrku*n
1268 *
1269 * Path 6o (M >= N, JOBZ='O')
1270 * Perform bidiagonal SVD, computing left singular vectors
1271 * of bidiagonal matrix in RWORK(IRU) and computing right
1272 * singular vectors of bidiagonal matrix in RWORK(IRVT)
1273 * CWorkspace: need 0
1274 * RWorkspace: need N [e] + N*N [RU] + N*N [RVT] + BDSPAC
1275 *
1276  CALL dbdsdc( 'U', 'I', n, s, rwork( ie ), rwork( iru ),
1277  $ n, rwork( irvt ), n, dum, idum,
1278  $ rwork( nrwork ), iwork, info )
1279 *
1280 * Copy real matrix RWORK(IRVT) to complex matrix VT
1281 * Overwrite VT by right singular vectors of A
1282 * CWorkspace: need 2*N [tauq, taup] + N*N [U] + N [work]
1283 * CWorkspace: prefer 2*N [tauq, taup] + N*N [U] + N*NB [work]
1284 * RWorkspace: need N [e] + N*N [RU] + N*N [RVT]
1285 *
1286  CALL zlacp2( 'F', n, n, rwork( irvt ), n, vt, ldvt )
1287  CALL zunmbr( 'P', 'R', 'C', n, n, n, a, lda,
1288  $ work( itaup ), vt, ldvt, work( nwork ),
1289  $ lwork-nwork+1, ierr )
1290 *
1291  IF( lwork .GE. m*n + 3*n ) THEN
1292 *
1293 * Path 6o-fast
1294 * Copy real matrix RWORK(IRU) to complex matrix WORK(IU)
1295 * Overwrite WORK(IU) by left singular vectors of A, copying
1296 * to A
1297 * CWorkspace: need 2*N [tauq, taup] + M*N [U] + N [work]
1298 * CWorkspace: prefer 2*N [tauq, taup] + M*N [U] + N*NB [work]
1299 * RWorkspace: need N [e] + N*N [RU]
1300 *
1301  CALL zlaset( 'F', m, n, czero, czero, work( iu ),
1302  $ ldwrku )
1303  CALL zlacp2( 'F', n, n, rwork( iru ), n, work( iu ),
1304  $ ldwrku )
1305  CALL zunmbr( 'Q', 'L', 'N', m, n, n, a, lda,
1306  $ work( itauq ), work( iu ), ldwrku,
1307  $ work( nwork ), lwork-nwork+1, ierr )
1308  CALL zlacpy( 'F', m, n, work( iu ), ldwrku, a, lda )
1309  ELSE
1310 *
1311 * Path 6o-slow
1312 * Generate Q in A
1313 * CWorkspace: need 2*N [tauq, taup] + N*N [U] + N [work]
1314 * CWorkspace: prefer 2*N [tauq, taup] + N*N [U] + N*NB [work]
1315 * RWorkspace: need 0
1316 *
1317  CALL zungbr( 'Q', m, n, n, a, lda, work( itauq ),
1318  $ work( nwork ), lwork-nwork+1, ierr )
1319 *
1320 * Multiply Q in A by real matrix RWORK(IRU), storing the
1321 * result in WORK(IU), copying to A
1322 * CWorkspace: need 2*N [tauq, taup] + N*N [U]
1323 * CWorkspace: prefer 2*N [tauq, taup] + M*N [U]
1324 * RWorkspace: need N [e] + N*N [RU] + 2*N*N [rwork]
1325 * RWorkspace: prefer N [e] + N*N [RU] + 2*M*N [rwork] < N + 5*N*N since M < 2*N here
1326 *
1327  nrwork = irvt
1328  DO 30 i = 1, m, ldwrku
1329  chunk = min( m-i+1, ldwrku )
1330  CALL zlacrm( chunk, n, a( i, 1 ), lda,
1331  $ rwork( iru ), n, work( iu ), ldwrku,
1332  $ rwork( nrwork ) )
1333  CALL zlacpy( 'F', chunk, n, work( iu ), ldwrku,
1334  $ a( i, 1 ), lda )
1335  30 CONTINUE
1336  END IF
1337 *
1338  ELSE IF( wntqs ) THEN
1339 *
1340 * Path 6s (M >= N, JOBZ='S')
1341 * Perform bidiagonal SVD, computing left singular vectors
1342 * of bidiagonal matrix in RWORK(IRU) and computing right
1343 * singular vectors of bidiagonal matrix in RWORK(IRVT)
1344 * CWorkspace: need 0
1345 * RWorkspace: need N [e] + N*N [RU] + N*N [RVT] + BDSPAC
1346 *
1347  iru = nrwork
1348  irvt = iru + n*n
1349  nrwork = irvt + n*n
1350  CALL dbdsdc( 'U', 'I', n, s, rwork( ie ), rwork( iru ),
1351  $ n, rwork( irvt ), n, dum, idum,
1352  $ rwork( nrwork ), iwork, info )
1353 *
1354 * Copy real matrix RWORK(IRU) to complex matrix U
1355 * Overwrite U by left singular vectors of A
1356 * CWorkspace: need 2*N [tauq, taup] + N [work]
1357 * CWorkspace: prefer 2*N [tauq, taup] + N*NB [work]
1358 * RWorkspace: need N [e] + N*N [RU] + N*N [RVT]
1359 *
1360  CALL zlaset( 'F', m, n, czero, czero, u, ldu )
1361  CALL zlacp2( 'F', n, n, rwork( iru ), n, u, ldu )
1362  CALL zunmbr( 'Q', 'L', 'N', m, n, n, a, lda,
1363  $ work( itauq ), u, ldu, work( nwork ),
1364  $ lwork-nwork+1, ierr )
1365 *
1366 * Copy real matrix RWORK(IRVT) to complex matrix VT
1367 * Overwrite VT by right singular vectors of A
1368 * CWorkspace: need 2*N [tauq, taup] + N [work]
1369 * CWorkspace: prefer 2*N [tauq, taup] + N*NB [work]
1370 * RWorkspace: need N [e] + N*N [RU] + N*N [RVT]
1371 *
1372  CALL zlacp2( 'F', n, n, rwork( irvt ), n, vt, ldvt )
1373  CALL zunmbr( 'P', 'R', 'C', n, n, n, a, lda,
1374  $ work( itaup ), vt, ldvt, work( nwork ),
1375  $ lwork-nwork+1, ierr )
1376  ELSE
1377 *
1378 * Path 6a (M >= N, JOBZ='A')
1379 * Perform bidiagonal SVD, computing left singular vectors
1380 * of bidiagonal matrix in RWORK(IRU) and computing right
1381 * singular vectors of bidiagonal matrix in RWORK(IRVT)
1382 * CWorkspace: need 0
1383 * RWorkspace: need N [e] + N*N [RU] + N*N [RVT] + BDSPAC
1384 *
1385  iru = nrwork
1386  irvt = iru + n*n
1387  nrwork = irvt + n*n
1388  CALL dbdsdc( 'U', 'I', n, s, rwork( ie ), rwork( iru ),
1389  $ n, rwork( irvt ), n, dum, idum,
1390  $ rwork( nrwork ), iwork, info )
1391 *
1392 * Set the right corner of U to identity matrix
1393 *
1394  CALL zlaset( 'F', m, m, czero, czero, u, ldu )
1395  IF( m.GT.n ) THEN
1396  CALL zlaset( 'F', m-n, m-n, czero, cone,
1397  $ u( n+1, n+1 ), ldu )
1398  END IF
1399 *
1400 * Copy real matrix RWORK(IRU) to complex matrix U
1401 * Overwrite U by left singular vectors of A
1402 * CWorkspace: need 2*N [tauq, taup] + M [work]
1403 * CWorkspace: prefer 2*N [tauq, taup] + M*NB [work]
1404 * RWorkspace: need N [e] + N*N [RU] + N*N [RVT]
1405 *
1406  CALL zlacp2( 'F', n, n, rwork( iru ), n, u, ldu )
1407  CALL zunmbr( 'Q', 'L', 'N', m, m, n, a, lda,
1408  $ work( itauq ), u, ldu, work( nwork ),
1409  $ lwork-nwork+1, ierr )
1410 *
1411 * Copy real matrix RWORK(IRVT) to complex matrix VT
1412 * Overwrite VT by right singular vectors of A
1413 * CWorkspace: need 2*N [tauq, taup] + N [work]
1414 * CWorkspace: prefer 2*N [tauq, taup] + N*NB [work]
1415 * RWorkspace: need N [e] + N*N [RU] + N*N [RVT]
1416 *
1417  CALL zlacp2( 'F', n, n, rwork( irvt ), n, vt, ldvt )
1418  CALL zunmbr( 'P', 'R', 'C', n, n, n, a, lda,
1419  $ work( itaup ), vt, ldvt, work( nwork ),
1420  $ lwork-nwork+1, ierr )
1421  END IF
1422 *
1423  END IF
1424 *
1425  ELSE
1426 *
1427 * A has more columns than rows. If A has sufficiently more
1428 * columns than rows, first reduce using the LQ decomposition (if
1429 * sufficient workspace available)
1430 *
1431  IF( n.GE.mnthr1 ) THEN
1432 *
1433  IF( wntqn ) THEN
1434 *
1435 * Path 1t (N >> M, JOBZ='N')
1436 * No singular vectors to be computed
1437 *
1438  itau = 1
1439  nwork = itau + m
1440 *
1441 * Compute A=L*Q
1442 * CWorkspace: need M [tau] + M [work]
1443 * CWorkspace: prefer M [tau] + M*NB [work]
1444 * RWorkspace: need 0
1445 *
1446  CALL zgelqf( m, n, a, lda, work( itau ), work( nwork ),
1447  $ lwork-nwork+1, ierr )
1448 *
1449 * Zero out above L
1450 *
1451  CALL zlaset( 'U', m-1, m-1, czero, czero, a( 1, 2 ),
1452  $ lda )
1453  ie = 1
1454  itauq = 1
1455  itaup = itauq + m
1456  nwork = itaup + m
1457 *
1458 * Bidiagonalize L in A
1459 * CWorkspace: need 2*M [tauq, taup] + M [work]
1460 * CWorkspace: prefer 2*M [tauq, taup] + 2*M*NB [work]
1461 * RWorkspace: need M [e]
1462 *
1463  CALL zgebrd( m, m, a, lda, s, rwork( ie ), work( itauq ),
1464  $ work( itaup ), work( nwork ), lwork-nwork+1,
1465  $ ierr )
1466  nrwork = ie + m
1467 *
1468 * Perform bidiagonal SVD, compute singular values only
1469 * CWorkspace: need 0
1470 * RWorkspace: need M [e] + BDSPAC
1471 *
1472  CALL dbdsdc( 'U', 'N', m, s, rwork( ie ), dum,1,dum,1,
1473  $ dum, idum, rwork( nrwork ), iwork, info )
1474 *
1475  ELSE IF( wntqo ) THEN
1476 *
1477 * Path 2t (N >> M, JOBZ='O')
1478 * M right singular vectors to be overwritten on A and
1479 * M left singular vectors to be computed in U
1480 *
1481  ivt = 1
1482  ldwkvt = m
1483 *
1484 * WORK(IVT) is M by M
1485 *
1486  il = ivt + ldwkvt*m
1487  IF( lwork .GE. m*n + m*m + 3*m ) THEN
1488 *
1489 * WORK(IL) M by N
1490 *
1491  ldwrkl = m
1492  chunk = n
1493  ELSE
1494 *
1495 * WORK(IL) is M by CHUNK
1496 *
1497  ldwrkl = m
1498  chunk = ( lwork - m*m - 3*m ) / m
1499  END IF
1500  itau = il + ldwrkl*chunk
1501  nwork = itau + m
1502 *
1503 * Compute A=L*Q
1504 * CWorkspace: need M*M [VT] + M*M [L] + M [tau] + M [work]
1505 * CWorkspace: prefer M*M [VT] + M*M [L] + M [tau] + M*NB [work]
1506 * RWorkspace: need 0
1507 *
1508  CALL zgelqf( m, n, a, lda, work( itau ), work( nwork ),
1509  $ lwork-nwork+1, ierr )
1510 *
1511 * Copy L to WORK(IL), zeroing about above it
1512 *
1513  CALL zlacpy( 'L', m, m, a, lda, work( il ), ldwrkl )
1514  CALL zlaset( 'U', m-1, m-1, czero, czero,
1515  $ work( il+ldwrkl ), ldwrkl )
1516 *
1517 * Generate Q in A
1518 * CWorkspace: need M*M [VT] + M*M [L] + M [tau] + M [work]
1519 * CWorkspace: prefer M*M [VT] + M*M [L] + M [tau] + M*NB [work]
1520 * RWorkspace: need 0
1521 *
1522  CALL zunglq( m, n, m, a, lda, work( itau ),
1523  $ work( nwork ), lwork-nwork+1, ierr )
1524  ie = 1
1525  itauq = itau
1526  itaup = itauq + m
1527  nwork = itaup + m
1528 *
1529 * Bidiagonalize L in WORK(IL)
1530 * CWorkspace: need M*M [VT] + M*M [L] + 2*M [tauq, taup] + M [work]
1531 * CWorkspace: prefer M*M [VT] + M*M [L] + 2*M [tauq, taup] + 2*M*NB [work]
1532 * RWorkspace: need M [e]
1533 *
1534  CALL zgebrd( m, m, work( il ), ldwrkl, s, rwork( ie ),
1535  $ work( itauq ), work( itaup ), work( nwork ),
1536  $ lwork-nwork+1, ierr )
1537 *
1538 * Perform bidiagonal SVD, computing left singular vectors
1539 * of bidiagonal matrix in RWORK(IRU) and computing right
1540 * singular vectors of bidiagonal matrix in RWORK(IRVT)
1541 * CWorkspace: need 0
1542 * RWorkspace: need M [e] + M*M [RU] + M*M [RVT] + BDSPAC
1543 *
1544  iru = ie + m
1545  irvt = iru + m*m
1546  nrwork = irvt + m*m
1547  CALL dbdsdc( 'U', 'I', m, s, rwork( ie ), rwork( iru ),
1548  $ m, rwork( irvt ), m, dum, idum,
1549  $ rwork( nrwork ), iwork, info )
1550 *
1551 * Copy real matrix RWORK(IRU) to complex matrix WORK(IU)
1552 * Overwrite WORK(IU) by the left singular vectors of L
1553 * CWorkspace: need M*M [VT] + M*M [L] + 2*M [tauq, taup] + M [work]
1554 * CWorkspace: prefer M*M [VT] + M*M [L] + 2*M [tauq, taup] + M*NB [work]
1555 * RWorkspace: need 0
1556 *
1557  CALL zlacp2( 'F', m, m, rwork( iru ), m, u, ldu )
1558  CALL zunmbr( 'Q', 'L', 'N', m, m, m, work( il ), ldwrkl,
1559  $ work( itauq ), u, ldu, work( nwork ),
1560  $ lwork-nwork+1, ierr )
1561 *
1562 * Copy real matrix RWORK(IRVT) to complex matrix WORK(IVT)
1563 * Overwrite WORK(IVT) by the right singular vectors of L
1564 * CWorkspace: need M*M [VT] + M*M [L] + 2*M [tauq, taup] + M [work]
1565 * CWorkspace: prefer M*M [VT] + M*M [L] + 2*M [tauq, taup] + M*NB [work]
1566 * RWorkspace: need 0
1567 *
1568  CALL zlacp2( 'F', m, m, rwork( irvt ), m, work( ivt ),
1569  $ ldwkvt )
1570  CALL zunmbr( 'P', 'R', 'C', m, m, m, work( il ), ldwrkl,
1571  $ work( itaup ), work( ivt ), ldwkvt,
1572  $ work( nwork ), lwork-nwork+1, ierr )
1573 *
1574 * Multiply right singular vectors of L in WORK(IL) by Q
1575 * in A, storing result in WORK(IL) and copying to A
1576 * CWorkspace: need M*M [VT] + M*M [L]
1577 * CWorkspace: prefer M*M [VT] + M*N [L]
1578 * RWorkspace: need 0
1579 *
1580  DO 40 i = 1, n, chunk
1581  blk = min( n-i+1, chunk )
1582  CALL zgemm( 'N', 'N', m, blk, m, cone, work( ivt ), m,
1583  $ a( 1, i ), lda, czero, work( il ),
1584  $ ldwrkl )
1585  CALL zlacpy( 'F', m, blk, work( il ), ldwrkl,
1586  $ a( 1, i ), lda )
1587  40 CONTINUE
1588 *
1589  ELSE IF( wntqs ) THEN
1590 *
1591 * Path 3t (N >> M, JOBZ='S')
1592 * M right singular vectors to be computed in VT and
1593 * M left singular vectors to be computed in U
1594 *
1595  il = 1
1596 *
1597 * WORK(IL) is M by M
1598 *
1599  ldwrkl = m
1600  itau = il + ldwrkl*m
1601  nwork = itau + m
1602 *
1603 * Compute A=L*Q
1604 * CWorkspace: need M*M [L] + M [tau] + M [work]
1605 * CWorkspace: prefer M*M [L] + M [tau] + M*NB [work]
1606 * RWorkspace: need 0
1607 *
1608  CALL zgelqf( m, n, a, lda, work( itau ), work( nwork ),
1609  $ lwork-nwork+1, ierr )
1610 *
1611 * Copy L to WORK(IL), zeroing out above it
1612 *
1613  CALL zlacpy( 'L', m, m, a, lda, work( il ), ldwrkl )
1614  CALL zlaset( 'U', m-1, m-1, czero, czero,
1615  $ work( il+ldwrkl ), ldwrkl )
1616 *
1617 * Generate Q in A
1618 * CWorkspace: need M*M [L] + M [tau] + M [work]
1619 * CWorkspace: prefer M*M [L] + M [tau] + M*NB [work]
1620 * RWorkspace: need 0
1621 *
1622  CALL zunglq( m, n, m, a, lda, work( itau ),
1623  $ work( nwork ), lwork-nwork+1, ierr )
1624  ie = 1
1625  itauq = itau
1626  itaup = itauq + m
1627  nwork = itaup + m
1628 *
1629 * Bidiagonalize L in WORK(IL)
1630 * CWorkspace: need M*M [L] + 2*M [tauq, taup] + M [work]
1631 * CWorkspace: prefer M*M [L] + 2*M [tauq, taup] + 2*M*NB [work]
1632 * RWorkspace: need M [e]
1633 *
1634  CALL zgebrd( m, m, work( il ), ldwrkl, s, rwork( ie ),
1635  $ work( itauq ), work( itaup ), work( nwork ),
1636  $ lwork-nwork+1, ierr )
1637 *
1638 * Perform bidiagonal SVD, computing left singular vectors
1639 * of bidiagonal matrix in RWORK(IRU) and computing right
1640 * singular vectors of bidiagonal matrix in RWORK(IRVT)
1641 * CWorkspace: need 0
1642 * RWorkspace: need M [e] + M*M [RU] + M*M [RVT] + BDSPAC
1643 *
1644  iru = ie + m
1645  irvt = iru + m*m
1646  nrwork = irvt + m*m
1647  CALL dbdsdc( 'U', 'I', m, s, rwork( ie ), rwork( iru ),
1648  $ m, rwork( irvt ), m, dum, idum,
1649  $ rwork( nrwork ), iwork, info )
1650 *
1651 * Copy real matrix RWORK(IRU) to complex matrix U
1652 * Overwrite U by left singular vectors of L
1653 * CWorkspace: need M*M [L] + 2*M [tauq, taup] + M [work]
1654 * CWorkspace: prefer M*M [L] + 2*M [tauq, taup] + M*NB [work]
1655 * RWorkspace: need 0
1656 *
1657  CALL zlacp2( 'F', m, m, rwork( iru ), m, u, ldu )
1658  CALL zunmbr( 'Q', 'L', 'N', m, m, m, work( il ), ldwrkl,
1659  $ work( itauq ), u, ldu, work( nwork ),
1660  $ lwork-nwork+1, ierr )
1661 *
1662 * Copy real matrix RWORK(IRVT) to complex matrix VT
1663 * Overwrite VT by left singular vectors of L
1664 * CWorkspace: need M*M [L] + 2*M [tauq, taup] + M [work]
1665 * CWorkspace: prefer M*M [L] + 2*M [tauq, taup] + M*NB [work]
1666 * RWorkspace: need 0
1667 *
1668  CALL zlacp2( 'F', m, m, rwork( irvt ), m, vt, ldvt )
1669  CALL zunmbr( 'P', 'R', 'C', m, m, m, work( il ), ldwrkl,
1670  $ work( itaup ), vt, ldvt, work( nwork ),
1671  $ lwork-nwork+1, ierr )
1672 *
1673 * Copy VT to WORK(IL), multiply right singular vectors of L
1674 * in WORK(IL) by Q in A, storing result in VT
1675 * CWorkspace: need M*M [L]
1676 * RWorkspace: need 0
1677 *
1678  CALL zlacpy( 'F', m, m, vt, ldvt, work( il ), ldwrkl )
1679  CALL zgemm( 'N', 'N', m, n, m, cone, work( il ), ldwrkl,
1680  $ a, lda, czero, vt, ldvt )
1681 *
1682  ELSE IF( wntqa ) THEN
1683 *
1684 * Path 4t (N >> M, JOBZ='A')
1685 * N right singular vectors to be computed in VT and
1686 * M left singular vectors to be computed in U
1687 *
1688  ivt = 1
1689 *
1690 * WORK(IVT) is M by M
1691 *
1692  ldwkvt = m
1693  itau = ivt + ldwkvt*m
1694  nwork = itau + m
1695 *
1696 * Compute A=L*Q, copying result to VT
1697 * CWorkspace: need M*M [VT] + M [tau] + M [work]
1698 * CWorkspace: prefer M*M [VT] + M [tau] + M*NB [work]
1699 * RWorkspace: need 0
1700 *
1701  CALL zgelqf( m, n, a, lda, work( itau ), work( nwork ),
1702  $ lwork-nwork+1, ierr )
1703  CALL zlacpy( 'U', m, n, a, lda, vt, ldvt )
1704 *
1705 * Generate Q in VT
1706 * CWorkspace: need M*M [VT] + M [tau] + N [work]
1707 * CWorkspace: prefer M*M [VT] + M [tau] + N*NB [work]
1708 * RWorkspace: need 0
1709 *
1710  CALL zunglq( n, n, m, vt, ldvt, work( itau ),
1711  $ work( nwork ), lwork-nwork+1, ierr )
1712 *
1713 * Produce L in A, zeroing out above it
1714 *
1715  CALL zlaset( 'U', m-1, m-1, czero, czero, a( 1, 2 ),
1716  $ lda )
1717  ie = 1
1718  itauq = itau
1719  itaup = itauq + m
1720  nwork = itaup + m
1721 *
1722 * Bidiagonalize L in A
1723 * CWorkspace: need M*M [VT] + 2*M [tauq, taup] + M [work]
1724 * CWorkspace: prefer M*M [VT] + 2*M [tauq, taup] + 2*M*NB [work]
1725 * RWorkspace: need M [e]
1726 *
1727  CALL zgebrd( m, m, a, lda, s, rwork( ie ), work( itauq ),
1728  $ work( itaup ), work( nwork ), lwork-nwork+1,
1729  $ ierr )
1730 *
1731 * Perform bidiagonal SVD, computing left singular vectors
1732 * of bidiagonal matrix in RWORK(IRU) and computing right
1733 * singular vectors of bidiagonal matrix in RWORK(IRVT)
1734 * CWorkspace: need 0
1735 * RWorkspace: need M [e] + M*M [RU] + M*M [RVT] + BDSPAC
1736 *
1737  iru = ie + m
1738  irvt = iru + m*m
1739  nrwork = irvt + m*m
1740  CALL dbdsdc( 'U', 'I', m, s, rwork( ie ), rwork( iru ),
1741  $ m, rwork( irvt ), m, dum, idum,
1742  $ rwork( nrwork ), iwork, info )
1743 *
1744 * Copy real matrix RWORK(IRU) to complex matrix U
1745 * Overwrite U by left singular vectors of L
1746 * CWorkspace: need M*M [VT] + 2*M [tauq, taup] + M [work]
1747 * CWorkspace: prefer M*M [VT] + 2*M [tauq, taup] + M*NB [work]
1748 * RWorkspace: need 0
1749 *
1750  CALL zlacp2( 'F', m, m, rwork( iru ), m, u, ldu )
1751  CALL zunmbr( 'Q', 'L', 'N', m, m, m, a, lda,
1752  $ work( itauq ), u, ldu, work( nwork ),
1753  $ lwork-nwork+1, ierr )
1754 *
1755 * Copy real matrix RWORK(IRVT) to complex matrix WORK(IVT)
1756 * Overwrite WORK(IVT) by right singular vectors of L
1757 * CWorkspace: need M*M [VT] + 2*M [tauq, taup] + M [work]
1758 * CWorkspace: prefer M*M [VT] + 2*M [tauq, taup] + M*NB [work]
1759 * RWorkspace: need 0
1760 *
1761  CALL zlacp2( 'F', m, m, rwork( irvt ), m, work( ivt ),
1762  $ ldwkvt )
1763  CALL zunmbr( 'P', 'R', 'C', m, m, m, a, lda,
1764  $ work( itaup ), work( ivt ), ldwkvt,
1765  $ work( nwork ), lwork-nwork+1, ierr )
1766 *
1767 * Multiply right singular vectors of L in WORK(IVT) by
1768 * Q in VT, storing result in A
1769 * CWorkspace: need M*M [VT]
1770 * RWorkspace: need 0
1771 *
1772  CALL zgemm( 'N', 'N', m, n, m, cone, work( ivt ), ldwkvt,
1773  $ vt, ldvt, czero, a, lda )
1774 *
1775 * Copy right singular vectors of A from A to VT
1776 *
1777  CALL zlacpy( 'F', m, n, a, lda, vt, ldvt )
1778 *
1779  END IF
1780 *
1781  ELSE IF( n.GE.mnthr2 ) THEN
1782 *
1783 * MNTHR2 <= N < MNTHR1
1784 *
1785 * Path 5t (N >> M, but not as much as MNTHR1)
1786 * Reduce to bidiagonal form without QR decomposition, use
1787 * ZUNGBR and matrix multiplication to compute singular vectors
1788 *
1789  ie = 1
1790  nrwork = ie + m
1791  itauq = 1
1792  itaup = itauq + m
1793  nwork = itaup + m
1794 *
1795 * Bidiagonalize A
1796 * CWorkspace: need 2*M [tauq, taup] + N [work]
1797 * CWorkspace: prefer 2*M [tauq, taup] + (M+N)*NB [work]
1798 * RWorkspace: need M [e]
1799 *
1800  CALL zgebrd( m, n, a, lda, s, rwork( ie ), work( itauq ),
1801  $ work( itaup ), work( nwork ), lwork-nwork+1,
1802  $ ierr )
1803 *
1804  IF( wntqn ) THEN
1805 *
1806 * Path 5tn (N >> M, JOBZ='N')
1807 * Compute singular values only
1808 * CWorkspace: need 0
1809 * RWorkspace: need M [e] + BDSPAC
1810 *
1811  CALL dbdsdc( 'L', 'N', m, s, rwork( ie ), dum,1,dum,1,
1812  $ dum, idum, rwork( nrwork ), iwork, info )
1813  ELSE IF( wntqo ) THEN
1814  irvt = nrwork
1815  iru = irvt + m*m
1816  nrwork = iru + m*m
1817  ivt = nwork
1818 *
1819 * Path 5to (N >> M, JOBZ='O')
1820 * Copy A to U, generate Q
1821 * CWorkspace: need 2*M [tauq, taup] + M [work]
1822 * CWorkspace: prefer 2*M [tauq, taup] + M*NB [work]
1823 * RWorkspace: need 0
1824 *
1825  CALL zlacpy( 'L', m, m, a, lda, u, ldu )
1826  CALL zungbr( 'Q', m, m, n, u, ldu, work( itauq ),
1827  $ work( nwork ), lwork-nwork+1, ierr )
1828 *
1829 * Generate P**H in A
1830 * CWorkspace: need 2*M [tauq, taup] + M [work]
1831 * CWorkspace: prefer 2*M [tauq, taup] + M*NB [work]
1832 * RWorkspace: need 0
1833 *
1834  CALL zungbr( 'P', m, n, m, a, lda, work( itaup ),
1835  $ work( nwork ), lwork-nwork+1, ierr )
1836 *
1837  ldwkvt = m
1838  IF( lwork .GE. m*n + 3*m ) THEN
1839 *
1840 * WORK( IVT ) is M by N
1841 *
1842  nwork = ivt + ldwkvt*n
1843  chunk = n
1844  ELSE
1845 *
1846 * WORK( IVT ) is M by CHUNK
1847 *
1848  chunk = ( lwork - 3*m ) / m
1849  nwork = ivt + ldwkvt*chunk
1850  END IF
1851 *
1852 * Perform bidiagonal SVD, computing left singular vectors
1853 * of bidiagonal matrix in RWORK(IRU) and computing right
1854 * singular vectors of bidiagonal matrix in RWORK(IRVT)
1855 * CWorkspace: need 0
1856 * RWorkspace: need M [e] + M*M [RVT] + M*M [RU] + BDSPAC
1857 *
1858  CALL dbdsdc( 'L', 'I', m, s, rwork( ie ), rwork( iru ),
1859  $ m, rwork( irvt ), m, dum, idum,
1860  $ rwork( nrwork ), iwork, info )
1861 *
1862 * Multiply Q in U by real matrix RWORK(IRVT)
1863 * storing the result in WORK(IVT), copying to U
1864 * CWorkspace: need 2*M [tauq, taup] + M*M [VT]
1865 * RWorkspace: need M [e] + M*M [RVT] + M*M [RU] + 2*M*M [rwork]
1866 *
1867  CALL zlacrm( m, m, u, ldu, rwork( iru ), m, work( ivt ),
1868  $ ldwkvt, rwork( nrwork ) )
1869  CALL zlacpy( 'F', m, m, work( ivt ), ldwkvt, u, ldu )
1870 *
1871 * Multiply RWORK(IRVT) by P**H in A, storing the
1872 * result in WORK(IVT), copying to A
1873 * CWorkspace: need 2*M [tauq, taup] + M*M [VT]
1874 * CWorkspace: prefer 2*M [tauq, taup] + M*N [VT]
1875 * RWorkspace: need M [e] + M*M [RVT] + 2*M*M [rwork]
1876 * RWorkspace: prefer M [e] + M*M [RVT] + 2*M*N [rwork] < M + 5*M*M since N < 2*M here
1877 *
1878  nrwork = iru
1879  DO 50 i = 1, n, chunk
1880  blk = min( n-i+1, chunk )
1881  CALL zlarcm( m, blk, rwork( irvt ), m, a( 1, i ), lda,
1882  $ work( ivt ), ldwkvt, rwork( nrwork ) )
1883  CALL zlacpy( 'F', m, blk, work( ivt ), ldwkvt,
1884  $ a( 1, i ), lda )
1885  50 CONTINUE
1886  ELSE IF( wntqs ) THEN
1887 *
1888 * Path 5ts (N >> M, JOBZ='S')
1889 * Copy A to U, generate Q
1890 * CWorkspace: need 2*M [tauq, taup] + M [work]
1891 * CWorkspace: prefer 2*M [tauq, taup] + M*NB [work]
1892 * RWorkspace: need 0
1893 *
1894  CALL zlacpy( 'L', m, m, a, lda, u, ldu )
1895  CALL zungbr( 'Q', m, m, n, u, ldu, work( itauq ),
1896  $ work( nwork ), lwork-nwork+1, ierr )
1897 *
1898 * Copy A to VT, generate P**H
1899 * CWorkspace: need 2*M [tauq, taup] + M [work]
1900 * CWorkspace: prefer 2*M [tauq, taup] + M*NB [work]
1901 * RWorkspace: need 0
1902 *
1903  CALL zlacpy( 'U', m, n, a, lda, vt, ldvt )
1904  CALL zungbr( 'P', m, n, m, vt, ldvt, work( itaup ),
1905  $ work( nwork ), lwork-nwork+1, ierr )
1906 *
1907 * Perform bidiagonal SVD, computing left singular vectors
1908 * of bidiagonal matrix in RWORK(IRU) and computing right
1909 * singular vectors of bidiagonal matrix in RWORK(IRVT)
1910 * CWorkspace: need 0
1911 * RWorkspace: need M [e] + M*M [RVT] + M*M [RU] + BDSPAC
1912 *
1913  irvt = nrwork
1914  iru = irvt + m*m
1915  nrwork = iru + m*m
1916  CALL dbdsdc( 'L', 'I', m, s, rwork( ie ), rwork( iru ),
1917  $ m, rwork( irvt ), m, dum, idum,
1918  $ rwork( nrwork ), iwork, info )
1919 *
1920 * Multiply Q in U by real matrix RWORK(IRU), storing the
1921 * result in A, copying to U
1922 * CWorkspace: need 0
1923 * RWorkspace: need M [e] + M*M [RVT] + M*M [RU] + 2*M*M [rwork]
1924 *
1925  CALL zlacrm( m, m, u, ldu, rwork( iru ), m, a, lda,
1926  $ rwork( nrwork ) )
1927  CALL zlacpy( 'F', m, m, a, lda, u, ldu )
1928 *
1929 * Multiply real matrix RWORK(IRVT) by P**H in VT,
1930 * storing the result in A, copying to VT
1931 * CWorkspace: need 0
1932 * RWorkspace: need M [e] + M*M [RVT] + 2*M*N [rwork] < M + 5*M*M since N < 2*M here
1933 *
1934  nrwork = iru
1935  CALL zlarcm( m, n, rwork( irvt ), m, vt, ldvt, a, lda,
1936  $ rwork( nrwork ) )
1937  CALL zlacpy( 'F', m, n, a, lda, vt, ldvt )
1938  ELSE
1939 *
1940 * Path 5ta (N >> M, JOBZ='A')
1941 * Copy A to U, generate Q
1942 * CWorkspace: need 2*M [tauq, taup] + M [work]
1943 * CWorkspace: prefer 2*M [tauq, taup] + M*NB [work]
1944 * RWorkspace: need 0
1945 *
1946  CALL zlacpy( 'L', m, m, a, lda, u, ldu )
1947  CALL zungbr( 'Q', m, m, n, u, ldu, work( itauq ),
1948  $ work( nwork ), lwork-nwork+1, ierr )
1949 *
1950 * Copy A to VT, generate P**H
1951 * CWorkspace: need 2*M [tauq, taup] + N [work]
1952 * CWorkspace: prefer 2*M [tauq, taup] + N*NB [work]
1953 * RWorkspace: need 0
1954 *
1955  CALL zlacpy( 'U', m, n, a, lda, vt, ldvt )
1956  CALL zungbr( 'P', n, n, m, vt, ldvt, work( itaup ),
1957  $ work( nwork ), lwork-nwork+1, ierr )
1958 *
1959 * Perform bidiagonal SVD, computing left singular vectors
1960 * of bidiagonal matrix in RWORK(IRU) and computing right
1961 * singular vectors of bidiagonal matrix in RWORK(IRVT)
1962 * CWorkspace: need 0
1963 * RWorkspace: need M [e] + M*M [RVT] + M*M [RU] + BDSPAC
1964 *
1965  irvt = nrwork
1966  iru = irvt + m*m
1967  nrwork = iru + m*m
1968  CALL dbdsdc( 'L', 'I', m, s, rwork( ie ), rwork( iru ),
1969  $ m, rwork( irvt ), m, dum, idum,
1970  $ rwork( nrwork ), iwork, info )
1971 *
1972 * Multiply Q in U by real matrix RWORK(IRU), storing the
1973 * result in A, copying to U
1974 * CWorkspace: need 0
1975 * RWorkspace: need M [e] + M*M [RVT] + M*M [RU] + 2*M*M [rwork]
1976 *
1977  CALL zlacrm( m, m, u, ldu, rwork( iru ), m, a, lda,
1978  $ rwork( nrwork ) )
1979  CALL zlacpy( 'F', m, m, a, lda, u, ldu )
1980 *
1981 * Multiply real matrix RWORK(IRVT) by P**H in VT,
1982 * storing the result in A, copying to VT
1983 * CWorkspace: need 0
1984 * RWorkspace: need M [e] + M*M [RVT] + 2*M*N [rwork] < M + 5*M*M since N < 2*M here
1985 *
1986  nrwork = iru
1987  CALL zlarcm( m, n, rwork( irvt ), m, vt, ldvt, a, lda,
1988  $ rwork( nrwork ) )
1989  CALL zlacpy( 'F', m, n, a, lda, vt, ldvt )
1990  END IF
1991 *
1992  ELSE
1993 *
1994 * N .LT. MNTHR2
1995 *
1996 * Path 6t (N > M, but not much larger)
1997 * Reduce to bidiagonal form without LQ decomposition
1998 * Use ZUNMBR to compute singular vectors
1999 *
2000  ie = 1
2001  nrwork = ie + m
2002  itauq = 1
2003  itaup = itauq + m
2004  nwork = itaup + m
2005 *
2006 * Bidiagonalize A
2007 * CWorkspace: need 2*M [tauq, taup] + N [work]
2008 * CWorkspace: prefer 2*M [tauq, taup] + (M+N)*NB [work]
2009 * RWorkspace: need M [e]
2010 *
2011  CALL zgebrd( m, n, a, lda, s, rwork( ie ), work( itauq ),
2012  $ work( itaup ), work( nwork ), lwork-nwork+1,
2013  $ ierr )
2014  IF( wntqn ) THEN
2015 *
2016 * Path 6tn (N > M, JOBZ='N')
2017 * Compute singular values only
2018 * CWorkspace: need 0
2019 * RWorkspace: need M [e] + BDSPAC
2020 *
2021  CALL dbdsdc( 'L', 'N', m, s, rwork( ie ), dum,1,dum,1,
2022  $ dum, idum, rwork( nrwork ), iwork, info )
2023  ELSE IF( wntqo ) THEN
2024 * Path 6to (N > M, JOBZ='O')
2025  ldwkvt = m
2026  ivt = nwork
2027  IF( lwork .GE. m*n + 3*m ) THEN
2028 *
2029 * WORK( IVT ) is M by N
2030 *
2031  CALL zlaset( 'F', m, n, czero, czero, work( ivt ),
2032  $ ldwkvt )
2033  nwork = ivt + ldwkvt*n
2034  ELSE
2035 *
2036 * WORK( IVT ) is M by CHUNK
2037 *
2038  chunk = ( lwork - 3*m ) / m
2039  nwork = ivt + ldwkvt*chunk
2040  END IF
2041 *
2042 * Perform bidiagonal SVD, computing left singular vectors
2043 * of bidiagonal matrix in RWORK(IRU) and computing right
2044 * singular vectors of bidiagonal matrix in RWORK(IRVT)
2045 * CWorkspace: need 0
2046 * RWorkspace: need M [e] + M*M [RVT] + M*M [RU] + BDSPAC
2047 *
2048  irvt = nrwork
2049  iru = irvt + m*m
2050  nrwork = iru + m*m
2051  CALL dbdsdc( 'L', 'I', m, s, rwork( ie ), rwork( iru ),
2052  $ m, rwork( irvt ), m, dum, idum,
2053  $ rwork( nrwork ), iwork, info )
2054 *
2055 * Copy real matrix RWORK(IRU) to complex matrix U
2056 * Overwrite U by left singular vectors of A
2057 * CWorkspace: need 2*M [tauq, taup] + M*M [VT] + M [work]
2058 * CWorkspace: prefer 2*M [tauq, taup] + M*M [VT] + M*NB [work]
2059 * RWorkspace: need M [e] + M*M [RVT] + M*M [RU]
2060 *
2061  CALL zlacp2( 'F', m, m, rwork( iru ), m, u, ldu )
2062  CALL zunmbr( 'Q', 'L', 'N', m, m, n, a, lda,
2063  $ work( itauq ), u, ldu, work( nwork ),
2064  $ lwork-nwork+1, ierr )
2065 *
2066  IF( lwork .GE. m*n + 3*m ) THEN
2067 *
2068 * Path 6to-fast
2069 * Copy real matrix RWORK(IRVT) to complex matrix WORK(IVT)
2070 * Overwrite WORK(IVT) by right singular vectors of A,
2071 * copying to A
2072 * CWorkspace: need 2*M [tauq, taup] + M*N [VT] + M [work]
2073 * CWorkspace: prefer 2*M [tauq, taup] + M*N [VT] + M*NB [work]
2074 * RWorkspace: need M [e] + M*M [RVT]
2075 *
2076  CALL zlacp2( 'F', m, m, rwork( irvt ), m, work( ivt ),
2077  $ ldwkvt )
2078  CALL zunmbr( 'P', 'R', 'C', m, n, m, a, lda,
2079  $ work( itaup ), work( ivt ), ldwkvt,
2080  $ work( nwork ), lwork-nwork+1, ierr )
2081  CALL zlacpy( 'F', m, n, work( ivt ), ldwkvt, a, lda )
2082  ELSE
2083 *
2084 * Path 6to-slow
2085 * Generate P**H in A
2086 * CWorkspace: need 2*M [tauq, taup] + M*M [VT] + M [work]
2087 * CWorkspace: prefer 2*M [tauq, taup] + M*M [VT] + M*NB [work]
2088 * RWorkspace: need 0
2089 *
2090  CALL zungbr( 'P', m, n, m, a, lda, work( itaup ),
2091  $ work( nwork ), lwork-nwork+1, ierr )
2092 *
2093 * Multiply Q in A by real matrix RWORK(IRU), storing the
2094 * result in WORK(IU), copying to A
2095 * CWorkspace: need 2*M [tauq, taup] + M*M [VT]
2096 * CWorkspace: prefer 2*M [tauq, taup] + M*N [VT]
2097 * RWorkspace: need M [e] + M*M [RVT] + 2*M*M [rwork]
2098 * RWorkspace: prefer M [e] + M*M [RVT] + 2*M*N [rwork] < M + 5*M*M since N < 2*M here
2099 *
2100  nrwork = iru
2101  DO 60 i = 1, n, chunk
2102  blk = min( n-i+1, chunk )
2103  CALL zlarcm( m, blk, rwork( irvt ), m, a( 1, i ),
2104  $ lda, work( ivt ), ldwkvt,
2105  $ rwork( nrwork ) )
2106  CALL zlacpy( 'F', m, blk, work( ivt ), ldwkvt,
2107  $ a( 1, i ), lda )
2108  60 CONTINUE
2109  END IF
2110  ELSE IF( wntqs ) THEN
2111 *
2112 * Path 6ts (N > M, JOBZ='S')
2113 * Perform bidiagonal SVD, computing left singular vectors
2114 * of bidiagonal matrix in RWORK(IRU) and computing right
2115 * singular vectors of bidiagonal matrix in RWORK(IRVT)
2116 * CWorkspace: need 0
2117 * RWorkspace: need M [e] + M*M [RVT] + M*M [RU] + BDSPAC
2118 *
2119  irvt = nrwork
2120  iru = irvt + m*m
2121  nrwork = iru + m*m
2122  CALL dbdsdc( 'L', 'I', m, s, rwork( ie ), rwork( iru ),
2123  $ m, rwork( irvt ), m, dum, idum,
2124  $ rwork( nrwork ), iwork, info )
2125 *
2126 * Copy real matrix RWORK(IRU) to complex matrix U
2127 * Overwrite U by left singular vectors of A
2128 * CWorkspace: need 2*M [tauq, taup] + M [work]
2129 * CWorkspace: prefer 2*M [tauq, taup] + M*NB [work]
2130 * RWorkspace: need M [e] + M*M [RVT] + M*M [RU]
2131 *
2132  CALL zlacp2( 'F', m, m, rwork( iru ), m, u, ldu )
2133  CALL zunmbr( 'Q', 'L', 'N', m, m, n, a, lda,
2134  $ work( itauq ), u, ldu, work( nwork ),
2135  $ lwork-nwork+1, ierr )
2136 *
2137 * Copy real matrix RWORK(IRVT) to complex matrix VT
2138 * Overwrite VT by right singular vectors of A
2139 * CWorkspace: need 2*M [tauq, taup] + M [work]
2140 * CWorkspace: prefer 2*M [tauq, taup] + M*NB [work]
2141 * RWorkspace: need M [e] + M*M [RVT]
2142 *
2143  CALL zlaset( 'F', m, n, czero, czero, vt, ldvt )
2144  CALL zlacp2( 'F', m, m, rwork( irvt ), m, vt, ldvt )
2145  CALL zunmbr( 'P', 'R', 'C', m, n, m, a, lda,
2146  $ work( itaup ), vt, ldvt, work( nwork ),
2147  $ lwork-nwork+1, ierr )
2148  ELSE
2149 *
2150 * Path 6ta (N > M, JOBZ='A')
2151 * Perform bidiagonal SVD, computing left singular vectors
2152 * of bidiagonal matrix in RWORK(IRU) and computing right
2153 * singular vectors of bidiagonal matrix in RWORK(IRVT)
2154 * CWorkspace: need 0
2155 * RWorkspace: need M [e] + M*M [RVT] + M*M [RU] + BDSPAC
2156 *
2157  irvt = nrwork
2158  iru = irvt + m*m
2159  nrwork = iru + m*m
2160 *
2161  CALL dbdsdc( 'L', 'I', m, s, rwork( ie ), rwork( iru ),
2162  $ m, rwork( irvt ), m, dum, idum,
2163  $ rwork( nrwork ), iwork, info )
2164 *
2165 * Copy real matrix RWORK(IRU) to complex matrix U
2166 * Overwrite U by left singular vectors of A
2167 * CWorkspace: need 2*M [tauq, taup] + M [work]
2168 * CWorkspace: prefer 2*M [tauq, taup] + M*NB [work]
2169 * RWorkspace: need M [e] + M*M [RVT] + M*M [RU]
2170 *
2171  CALL zlacp2( 'F', m, m, rwork( iru ), m, u, ldu )
2172  CALL zunmbr( 'Q', 'L', 'N', m, m, n, a, lda,
2173  $ work( itauq ), u, ldu, work( nwork ),
2174  $ lwork-nwork+1, ierr )
2175 *
2176 * Set all of VT to identity matrix
2177 *
2178  CALL zlaset( 'F', n, n, czero, cone, vt, ldvt )
2179 *
2180 * Copy real matrix RWORK(IRVT) to complex matrix VT
2181 * Overwrite VT by right singular vectors of A
2182 * CWorkspace: need 2*M [tauq, taup] + N [work]
2183 * CWorkspace: prefer 2*M [tauq, taup] + N*NB [work]
2184 * RWorkspace: need M [e] + M*M [RVT]
2185 *
2186  CALL zlacp2( 'F', m, m, rwork( irvt ), m, vt, ldvt )
2187  CALL zunmbr( 'P', 'R', 'C', n, n, m, a, lda,
2188  $ work( itaup ), vt, ldvt, work( nwork ),
2189  $ lwork-nwork+1, ierr )
2190  END IF
2191 *
2192  END IF
2193 *
2194  END IF
2195 *
2196 * Undo scaling if necessary
2197 *
2198  IF( iscl.EQ.1 ) THEN
2199  IF( anrm.GT.bignum )
2200  $ CALL dlascl( 'G', 0, 0, bignum, anrm, minmn, 1, s, minmn,
2201  $ ierr )
2202  IF( info.NE.0 .AND. anrm.GT.bignum )
2203  $ CALL dlascl( 'G', 0, 0, bignum, anrm, minmn-1, 1,
2204  $ rwork( ie ), minmn, ierr )
2205  IF( anrm.LT.smlnum )
2206  $ CALL dlascl( 'G', 0, 0, smlnum, anrm, minmn, 1, s, minmn,
2207  $ ierr )
2208  IF( info.NE.0 .AND. anrm.LT.smlnum )
2209  $ CALL dlascl( 'G', 0, 0, smlnum, anrm, minmn-1, 1,
2210  $ rwork( ie ), minmn, ierr )
2211  END IF
2212 *
2213 * Return optimal workspace in WORK(1)
2214 *
2215  work( 1 ) = maxwrk
2216 *
2217  RETURN
2218 *
2219 * End of ZGESDD
2220 *
2221  END
subroutine dlascl(TYPE, KL, KU, CFROM, CTO, M, N, A, LDA, INFO)
DLASCL multiplies a general rectangular matrix by a real scalar defined as cto/cfrom.
Definition: dlascl.f:143
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:60
subroutine dbdsdc(UPLO, COMPQ, N, D, E, U, LDU, VT, LDVT, Q, IQ, WORK, IWORK, INFO)
DBDSDC
Definition: dbdsdc.f:205
subroutine zgemm(TRANSA, TRANSB, M, N, K, ALPHA, A, LDA, B, LDB, BETA, C, LDC)
ZGEMM
Definition: zgemm.f:187
subroutine zungbr(VECT, M, N, K, A, LDA, TAU, WORK, LWORK, INFO)
ZUNGBR
Definition: zungbr.f:157
subroutine zgelqf(M, N, A, LDA, TAU, WORK, LWORK, INFO)
ZGELQF
Definition: zgelqf.f:143
subroutine zgebrd(M, N, A, LDA, D, E, TAUQ, TAUP, WORK, LWORK, INFO)
ZGEBRD
Definition: zgebrd.f:205
subroutine zgesdd(JOBZ, M, N, A, LDA, S, U, LDU, VT, LDVT, WORK, LWORK, RWORK, IWORK, INFO)
ZGESDD
Definition: zgesdd.f:226
subroutine zlarcm(M, N, A, LDA, B, LDB, C, LDC, RWORK)
ZLARCM copies all or part of a real two-dimensional array to a complex array.
Definition: zlarcm.f:114
subroutine zlascl(TYPE, KL, KU, CFROM, CTO, M, N, A, LDA, INFO)
ZLASCL multiplies a general rectangular matrix by a real scalar defined as cto/cfrom.
Definition: zlascl.f:143
subroutine zlacpy(UPLO, M, N, A, LDA, B, LDB)
ZLACPY copies all or part of one two-dimensional array to another.
Definition: zlacpy.f:103
subroutine zlacp2(UPLO, M, N, A, LDA, B, LDB)
ZLACP2 copies all or part of a real two-dimensional array to a complex array.
Definition: zlacp2.f:104
subroutine zlaset(UPLO, M, N, ALPHA, BETA, A, LDA)
ZLASET initializes the off-diagonal elements and the diagonal elements of a matrix to given values.
Definition: zlaset.f:106
subroutine zlacrm(M, N, A, LDA, B, LDB, C, LDC, RWORK)
ZLACRM multiplies a complex matrix by a square real matrix.
Definition: zlacrm.f:114
subroutine zunglq(M, N, K, A, LDA, TAU, WORK, LWORK, INFO)
ZUNGLQ
Definition: zunglq.f:127
subroutine zungqr(M, N, K, A, LDA, TAU, WORK, LWORK, INFO)
ZUNGQR
Definition: zungqr.f:128
subroutine zunmbr(VECT, SIDE, TRANS, M, N, K, A, LDA, TAU, C, LDC, WORK, LWORK, INFO)
ZUNMBR
Definition: zunmbr.f:196
subroutine zgeqrf(M, N, A, LDA, TAU, WORK, LWORK, INFO)
ZGEQRF VARIANT: left-looking Level 3 BLAS of the algorithm.
Definition: zgeqrf.f:151