LAPACK  3.10.0
LAPACK: Linear Algebra PACKage

◆ zgesdd()

subroutine zgesdd ( character  JOBZ,
integer  M,
integer  N,
complex*16, dimension( lda, * )  A,
integer  LDA,
double precision, dimension( * )  S,
complex*16, dimension( ldu, * )  U,
integer  LDU,
complex*16, dimension( ldvt, * )  VT,
integer  LDVT,
complex*16, dimension( * )  WORK,
integer  LWORK,
double precision, dimension( * )  RWORK,
integer, dimension( * )  IWORK,
integer  INFO 
)

ZGESDD

Download ZGESDD + dependencies [TGZ] [ZIP] [TXT]

Purpose:
 ZGESDD computes the singular value decomposition (SVD) of a complex
 M-by-N matrix A, optionally computing the left and/or right singular
 vectors, by using divide-and-conquer method. The SVD is written

      A = U * SIGMA * conjugate-transpose(V)

 where SIGMA is an M-by-N matrix which is zero except for its
 min(m,n) diagonal elements, U is an M-by-M unitary matrix, and
 V is an N-by-N unitary matrix.  The diagonal elements of SIGMA
 are the singular values of A; they are real and non-negative, and
 are returned in descending order.  The first min(m,n) columns of
 U and V are the left and right singular vectors of A.

 Note that the routine returns VT = V**H, not V.

 The divide and conquer algorithm makes very mild assumptions about
 floating point arithmetic. It will work on machines with a guard
 digit in add/subtract, or on those binary machines without guard
 digits which subtract like the Cray X-MP, Cray Y-MP, Cray C-90, or
 Cray-2. It could conceivably fail on hexadecimal or decimal machines
 without guard digits, but we know of none.
Parameters
[in]JOBZ
          JOBZ is CHARACTER*1
          Specifies options for computing all or part of the matrix U:
          = 'A':  all M columns of U and all N rows of V**H are
                  returned in the arrays U and VT;
          = 'S':  the first min(M,N) columns of U and the first
                  min(M,N) rows of V**H are returned in the arrays U
                  and VT;
          = 'O':  If M >= N, the first N columns of U are overwritten
                  in the array A and all rows of V**H are returned in
                  the array VT;
                  otherwise, all columns of U are returned in the
                  array U and the first M rows of V**H are overwritten
                  in the array A;
          = 'N':  no columns of U or rows of V**H are computed.
[in]M
          M is INTEGER
          The number of rows of the input matrix A.  M >= 0.
[in]N
          N is INTEGER
          The number of columns of the input matrix A.  N >= 0.
[in,out]A
          A is COMPLEX*16 array, dimension (LDA,N)
          On entry, the M-by-N matrix A.
          On exit,
          if JOBZ = 'O',  A is overwritten with the first N columns
                          of U (the left singular vectors, stored
                          columnwise) if M >= N;
                          A is overwritten with the first M rows
                          of V**H (the right singular vectors, stored
                          rowwise) otherwise.
          if JOBZ .ne. 'O', the contents of A are destroyed.
[in]LDA
          LDA is INTEGER
          The leading dimension of the array A.  LDA >= max(1,M).
[out]S
          S is DOUBLE PRECISION array, dimension (min(M,N))
          The singular values of A, sorted so that S(i) >= S(i+1).
[out]U
          U is COMPLEX*16 array, dimension (LDU,UCOL)
          UCOL = M if JOBZ = 'A' or JOBZ = 'O' and M < N;
          UCOL = min(M,N) if JOBZ = 'S'.
          If JOBZ = 'A' or JOBZ = 'O' and M < N, U contains the M-by-M
          unitary matrix U;
          if JOBZ = 'S', U contains the first min(M,N) columns of U
          (the left singular vectors, stored columnwise);
          if JOBZ = 'O' and M >= N, or JOBZ = 'N', U is not referenced.
[in]LDU
          LDU is INTEGER
          The leading dimension of the array U.  LDU >= 1;
          if JOBZ = 'S' or 'A' or JOBZ = 'O' and M < N, LDU >= M.
[out]VT
          VT is COMPLEX*16 array, dimension (LDVT,N)
          If JOBZ = 'A' or JOBZ = 'O' and M >= N, VT contains the
          N-by-N unitary matrix V**H;
          if JOBZ = 'S', VT contains the first min(M,N) rows of
          V**H (the right singular vectors, stored rowwise);
          if JOBZ = 'O' and M < N, or JOBZ = 'N', VT is not referenced.
[in]LDVT
          LDVT is INTEGER
          The leading dimension of the array VT.  LDVT >= 1;
          if JOBZ = 'A' or JOBZ = 'O' and M >= N, LDVT >= N;
          if JOBZ = 'S', LDVT >= min(M,N).
[out]WORK
          WORK is COMPLEX*16 array, dimension (MAX(1,LWORK))
          On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
[in]LWORK
          LWORK is INTEGER
          The dimension of the array WORK. LWORK >= 1.
          If LWORK = -1, a workspace query is assumed.  The optimal
          size for the WORK array is calculated and stored in WORK(1),
          and no other work except argument checking is performed.

          Let mx = max(M,N) and mn = min(M,N).
          If JOBZ = 'N', LWORK >= 2*mn + mx.
          If JOBZ = 'O', LWORK >= 2*mn*mn + 2*mn + mx.
          If JOBZ = 'S', LWORK >=   mn*mn + 3*mn.
          If JOBZ = 'A', LWORK >=   mn*mn + 2*mn + mx.
          These are not tight minimums in all cases; see comments inside code.
          For good performance, LWORK should generally be larger;
          a query is recommended.
[out]RWORK
          RWORK is DOUBLE PRECISION array, dimension (MAX(1,LRWORK))
          Let mx = max(M,N) and mn = min(M,N).
          If JOBZ = 'N',    LRWORK >= 5*mn (LAPACK <= 3.6 needs 7*mn);
          else if mx >> mn, LRWORK >= 5*mn*mn + 5*mn;
          else              LRWORK >= max( 5*mn*mn + 5*mn,
                                           2*mx*mn + 2*mn*mn + mn ).
[out]IWORK
          IWORK is INTEGER array, dimension (8*min(M,N))
[out]INFO
          INFO is INTEGER
          <  0:  if INFO = -i, the i-th argument had an illegal value.
          = -4:  if A had a NAN entry.
          >  0:  The updating process of DBDSDC did not converge.
          =  0:  successful exit.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Contributors:
Ming Gu and Huan Ren, Computer Science Division, University of California at Berkeley, USA

Definition at line 225 of file zgesdd.f.

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