LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine zgelsd ( integer  M,
integer  N,
integer  NRHS,
complex*16, dimension( lda, * )  A,
integer  LDA,
complex*16, dimension( ldb, * )  B,
integer  LDB,
double precision, dimension( * )  S,
double precision  RCOND,
integer  RANK,
complex*16, dimension( * )  WORK,
integer  LWORK,
double precision, dimension( * )  RWORK,
integer, dimension( * )  IWORK,
integer  INFO 
)

ZGELSD computes the minimum-norm solution to a linear least squares problem for GE matrices

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

Purpose:
 ZGELSD computes the minimum-norm solution to a real linear least
 squares problem:
     minimize 2-norm(| b - A*x |)
 using the singular value decomposition (SVD) of A. A is an M-by-N
 matrix which may be rank-deficient.

 Several right hand side vectors b and solution vectors x can be
 handled in a single call; they are stored as the columns of the
 M-by-NRHS right hand side matrix B and the N-by-NRHS solution
 matrix X.

 The problem is solved in three steps:
 (1) Reduce the coefficient matrix A to bidiagonal form with
     Householder tranformations, reducing the original problem
     into a "bidiagonal least squares problem" (BLS)
 (2) Solve the BLS using a divide and conquer approach.
 (3) Apply back all the Householder tranformations to solve
     the original least squares problem.

 The effective rank of A is determined by treating as zero those
 singular values which are less than RCOND times the largest singular
 value.

 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]M
          M is INTEGER
          The number of rows of the matrix A. M >= 0.
[in]N
          N is INTEGER
          The number of columns of the matrix A. N >= 0.
[in]NRHS
          NRHS is INTEGER
          The number of right hand sides, i.e., the number of columns
          of the matrices B and X. NRHS >= 0.
[in]A
          A is COMPLEX*16 array, dimension (LDA,N)
          On entry, the M-by-N matrix A.
          On exit, A has been destroyed.
[in]LDA
          LDA is INTEGER
          The leading dimension of the array A. LDA >= max(1,M).
[in,out]B
          B is COMPLEX*16 array, dimension (LDB,NRHS)
          On entry, the M-by-NRHS right hand side matrix B.
          On exit, B is overwritten by the N-by-NRHS solution matrix X.
          If m >= n and RANK = n, the residual sum-of-squares for
          the solution in the i-th column is given by the sum of
          squares of the modulus of elements n+1:m in that column.
[in]LDB
          LDB is INTEGER
          The leading dimension of the array B.  LDB >= max(1,M,N).
[out]S
          S is DOUBLE PRECISION array, dimension (min(M,N))
          The singular values of A in decreasing order.
          The condition number of A in the 2-norm = S(1)/S(min(m,n)).
[in]RCOND
          RCOND is DOUBLE PRECISION
          RCOND is used to determine the effective rank of A.
          Singular values S(i) <= RCOND*S(1) are treated as zero.
          If RCOND < 0, machine precision is used instead.
[out]RANK
          RANK is INTEGER
          The effective rank of A, i.e., the number of singular values
          which are greater than RCOND*S(1).
[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 must be at least 1.
          The exact minimum amount of workspace needed depends on M,
          N and NRHS. As long as LWORK is at least
              2*N + N*NRHS
          if M is greater than or equal to N or
              2*M + M*NRHS
          if M is less than N, the code will execute correctly.
          For good performance, LWORK should generally be larger.

          If LWORK = -1, then a workspace query is assumed; the routine
          only calculates the optimal size of the array WORK and the
          minimum sizes of the arrays RWORK and IWORK, and returns
          these values as the first entries of the WORK, RWORK and
          IWORK arrays, and no error message related to LWORK is issued
          by XERBLA.
[out]RWORK
          RWORK is DOUBLE PRECISION array, dimension (MAX(1,LRWORK))
          LRWORK >=
             10*N + 2*N*SMLSIZ + 8*N*NLVL + 3*SMLSIZ*NRHS +
             MAX( (SMLSIZ+1)**2, N*(1+NRHS) + 2*NRHS )
          if M is greater than or equal to N or
             10*M + 2*M*SMLSIZ + 8*M*NLVL + 3*SMLSIZ*NRHS +
             MAX( (SMLSIZ+1)**2, N*(1+NRHS) + 2*NRHS )
          if M is less than N, the code will execute correctly.
          SMLSIZ is returned by ILAENV and is equal to the maximum
          size of the subproblems at the bottom of the computation
          tree (usually about 25), and
             NLVL = MAX( 0, INT( LOG_2( MIN( M,N )/(SMLSIZ+1) ) ) + 1 )
          On exit, if INFO = 0, RWORK(1) returns the minimum LRWORK.
[out]IWORK
          IWORK is INTEGER array, dimension (MAX(1,LIWORK))
          LIWORK >= max(1, 3*MINMN*NLVL + 11*MINMN),
          where MINMN = MIN( M,N ).
          On exit, if INFO = 0, IWORK(1) returns the minimum LIWORK.
[out]INFO
          INFO is INTEGER
          = 0: successful exit
          < 0: if INFO = -i, the i-th argument had an illegal value.
          > 0:  the algorithm for computing the SVD failed to converge;
                if INFO = i, i off-diagonal elements of an intermediate
                bidiagonal form did not converge to zero.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
November 2011
Contributors:
Ming Gu and Ren-Cang Li, Computer Science Division, University of California at Berkeley, USA
Osni Marques, LBNL/NERSC, USA

Definition at line 227 of file zgelsd.f.

227 *
228 * -- LAPACK driver routine (version 3.4.0) --
229 * -- LAPACK is a software package provided by Univ. of Tennessee, --
230 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
231 * November 2011
232 *
233 * .. Scalar Arguments ..
234  INTEGER info, lda, ldb, lwork, m, n, nrhs, rank
235  DOUBLE PRECISION rcond
236 * ..
237 * .. Array Arguments ..
238  INTEGER iwork( * )
239  DOUBLE PRECISION rwork( * ), s( * )
240  COMPLEX*16 a( lda, * ), b( ldb, * ), work( * )
241 * ..
242 *
243 * =====================================================================
244 *
245 * .. Parameters ..
246  DOUBLE PRECISION zero, one, two
247  parameter ( zero = 0.0d+0, one = 1.0d+0, two = 2.0d+0 )
248  COMPLEX*16 czero
249  parameter ( czero = ( 0.0d+0, 0.0d+0 ) )
250 * ..
251 * .. Local Scalars ..
252  LOGICAL lquery
253  INTEGER iascl, ibscl, ie, il, itau, itaup, itauq,
254  $ ldwork, liwork, lrwork, maxmn, maxwrk, minmn,
255  $ minwrk, mm, mnthr, nlvl, nrwork, nwork, smlsiz
256  DOUBLE PRECISION anrm, bignum, bnrm, eps, sfmin, smlnum
257 * ..
258 * .. External Subroutines ..
259  EXTERNAL dlabad, dlascl, dlaset, xerbla, zgebrd, zgelqf,
261  $ zunmlq, zunmqr
262 * ..
263 * .. External Functions ..
264  INTEGER ilaenv
265  DOUBLE PRECISION dlamch, zlange
266  EXTERNAL ilaenv, dlamch, zlange
267 * ..
268 * .. Intrinsic Functions ..
269  INTRINSIC int, log, max, min, dble
270 * ..
271 * .. Executable Statements ..
272 *
273 * Test the input arguments.
274 *
275  info = 0
276  minmn = min( m, n )
277  maxmn = max( m, n )
278  lquery = ( lwork.EQ.-1 )
279  IF( m.LT.0 ) THEN
280  info = -1
281  ELSE IF( n.LT.0 ) THEN
282  info = -2
283  ELSE IF( nrhs.LT.0 ) THEN
284  info = -3
285  ELSE IF( lda.LT.max( 1, m ) ) THEN
286  info = -5
287  ELSE IF( ldb.LT.max( 1, maxmn ) ) THEN
288  info = -7
289  END IF
290 *
291 * Compute workspace.
292 * (Note: Comments in the code beginning "Workspace:" describe the
293 * minimal amount of workspace needed at that point in the code,
294 * as well as the preferred amount for good performance.
295 * NB refers to the optimal block size for the immediately
296 * following subroutine, as returned by ILAENV.)
297 *
298  IF( info.EQ.0 ) THEN
299  minwrk = 1
300  maxwrk = 1
301  liwork = 1
302  lrwork = 1
303  IF( minmn.GT.0 ) THEN
304  smlsiz = ilaenv( 9, 'ZGELSD', ' ', 0, 0, 0, 0 )
305  mnthr = ilaenv( 6, 'ZGELSD', ' ', m, n, nrhs, -1 )
306  nlvl = max( int( log( dble( minmn ) / dble( smlsiz + 1 ) ) /
307  $ log( two ) ) + 1, 0 )
308  liwork = 3*minmn*nlvl + 11*minmn
309  mm = m
310  IF( m.GE.n .AND. m.GE.mnthr ) THEN
311 *
312 * Path 1a - overdetermined, with many more rows than
313 * columns.
314 *
315  mm = n
316  maxwrk = max( maxwrk, n*ilaenv( 1, 'ZGEQRF', ' ', m, n,
317  $ -1, -1 ) )
318  maxwrk = max( maxwrk, nrhs*ilaenv( 1, 'ZUNMQR', 'LC', m,
319  $ nrhs, n, -1 ) )
320  END IF
321  IF( m.GE.n ) THEN
322 *
323 * Path 1 - overdetermined or exactly determined.
324 *
325  lrwork = 10*n + 2*n*smlsiz + 8*n*nlvl + 3*smlsiz*nrhs +
326  $ max( (smlsiz+1)**2, n*(1+nrhs) + 2*nrhs )
327  maxwrk = max( maxwrk, 2*n + ( mm + n )*ilaenv( 1,
328  $ 'ZGEBRD', ' ', mm, n, -1, -1 ) )
329  maxwrk = max( maxwrk, 2*n + nrhs*ilaenv( 1, 'ZUNMBR',
330  $ 'QLC', mm, nrhs, n, -1 ) )
331  maxwrk = max( maxwrk, 2*n + ( n - 1 )*ilaenv( 1,
332  $ 'ZUNMBR', 'PLN', n, nrhs, n, -1 ) )
333  maxwrk = max( maxwrk, 2*n + n*nrhs )
334  minwrk = max( 2*n + mm, 2*n + n*nrhs )
335  END IF
336  IF( n.GT.m ) THEN
337  lrwork = 10*m + 2*m*smlsiz + 8*m*nlvl + 3*smlsiz*nrhs +
338  $ max( (smlsiz+1)**2, n*(1+nrhs) + 2*nrhs )
339  IF( n.GE.mnthr ) THEN
340 *
341 * Path 2a - underdetermined, with many more columns
342 * than rows.
343 *
344  maxwrk = m + m*ilaenv( 1, 'ZGELQF', ' ', m, n, -1,
345  $ -1 )
346  maxwrk = max( maxwrk, m*m + 4*m + 2*m*ilaenv( 1,
347  $ 'ZGEBRD', ' ', m, m, -1, -1 ) )
348  maxwrk = max( maxwrk, m*m + 4*m + nrhs*ilaenv( 1,
349  $ 'ZUNMBR', 'QLC', m, nrhs, m, -1 ) )
350  maxwrk = max( maxwrk, m*m + 4*m + ( m - 1 )*ilaenv( 1,
351  $ 'ZUNMLQ', 'LC', n, nrhs, m, -1 ) )
352  IF( nrhs.GT.1 ) THEN
353  maxwrk = max( maxwrk, m*m + m + m*nrhs )
354  ELSE
355  maxwrk = max( maxwrk, m*m + 2*m )
356  END IF
357  maxwrk = max( maxwrk, m*m + 4*m + m*nrhs )
358 ! XXX: Ensure the Path 2a case below is triggered. The workspace
359 ! calculation should use queries for all routines eventually.
360  maxwrk = max( maxwrk,
361  $ 4*m+m*m+max( m, 2*m-4, nrhs, n-3*m ) )
362  ELSE
363 *
364 * Path 2 - underdetermined.
365 *
366  maxwrk = 2*m + ( n + m )*ilaenv( 1, 'ZGEBRD', ' ', m,
367  $ n, -1, -1 )
368  maxwrk = max( maxwrk, 2*m + nrhs*ilaenv( 1, 'ZUNMBR',
369  $ 'QLC', m, nrhs, m, -1 ) )
370  maxwrk = max( maxwrk, 2*m + m*ilaenv( 1, 'ZUNMBR',
371  $ 'PLN', n, nrhs, m, -1 ) )
372  maxwrk = max( maxwrk, 2*m + m*nrhs )
373  END IF
374  minwrk = max( 2*m + n, 2*m + m*nrhs )
375  END IF
376  END IF
377  minwrk = min( minwrk, maxwrk )
378  work( 1 ) = maxwrk
379  iwork( 1 ) = liwork
380  rwork( 1 ) = lrwork
381 *
382  IF( lwork.LT.minwrk .AND. .NOT.lquery ) THEN
383  info = -12
384  END IF
385  END IF
386 *
387  IF( info.NE.0 ) THEN
388  CALL xerbla( 'ZGELSD', -info )
389  RETURN
390  ELSE IF( lquery ) THEN
391  RETURN
392  END IF
393 *
394 * Quick return if possible.
395 *
396  IF( m.EQ.0 .OR. n.EQ.0 ) THEN
397  rank = 0
398  RETURN
399  END IF
400 *
401 * Get machine parameters.
402 *
403  eps = dlamch( 'P' )
404  sfmin = dlamch( 'S' )
405  smlnum = sfmin / eps
406  bignum = one / smlnum
407  CALL dlabad( smlnum, bignum )
408 *
409 * Scale A if max entry outside range [SMLNUM,BIGNUM].
410 *
411  anrm = zlange( 'M', m, n, a, lda, rwork )
412  iascl = 0
413  IF( anrm.GT.zero .AND. anrm.LT.smlnum ) THEN
414 *
415 * Scale matrix norm up to SMLNUM
416 *
417  CALL zlascl( 'G', 0, 0, anrm, smlnum, m, n, a, lda, info )
418  iascl = 1
419  ELSE IF( anrm.GT.bignum ) THEN
420 *
421 * Scale matrix norm down to BIGNUM.
422 *
423  CALL zlascl( 'G', 0, 0, anrm, bignum, m, n, a, lda, info )
424  iascl = 2
425  ELSE IF( anrm.EQ.zero ) THEN
426 *
427 * Matrix all zero. Return zero solution.
428 *
429  CALL zlaset( 'F', max( m, n ), nrhs, czero, czero, b, ldb )
430  CALL dlaset( 'F', minmn, 1, zero, zero, s, 1 )
431  rank = 0
432  GO TO 10
433  END IF
434 *
435 * Scale B if max entry outside range [SMLNUM,BIGNUM].
436 *
437  bnrm = zlange( 'M', m, nrhs, b, ldb, rwork )
438  ibscl = 0
439  IF( bnrm.GT.zero .AND. bnrm.LT.smlnum ) THEN
440 *
441 * Scale matrix norm up to SMLNUM.
442 *
443  CALL zlascl( 'G', 0, 0, bnrm, smlnum, m, nrhs, b, ldb, info )
444  ibscl = 1
445  ELSE IF( bnrm.GT.bignum ) THEN
446 *
447 * Scale matrix norm down to BIGNUM.
448 *
449  CALL zlascl( 'G', 0, 0, bnrm, bignum, m, nrhs, b, ldb, info )
450  ibscl = 2
451  END IF
452 *
453 * If M < N make sure B(M+1:N,:) = 0
454 *
455  IF( m.LT.n )
456  $ CALL zlaset( 'F', n-m, nrhs, czero, czero, b( m+1, 1 ), ldb )
457 *
458 * Overdetermined case.
459 *
460  IF( m.GE.n ) THEN
461 *
462 * Path 1 - overdetermined or exactly determined.
463 *
464  mm = m
465  IF( m.GE.mnthr ) THEN
466 *
467 * Path 1a - overdetermined, with many more rows than columns
468 *
469  mm = n
470  itau = 1
471  nwork = itau + n
472 *
473 * Compute A=Q*R.
474 * (RWorkspace: need N)
475 * (CWorkspace: need N, prefer N*NB)
476 *
477  CALL zgeqrf( m, n, a, lda, work( itau ), work( nwork ),
478  $ lwork-nwork+1, info )
479 *
480 * Multiply B by transpose(Q).
481 * (RWorkspace: need N)
482 * (CWorkspace: need NRHS, prefer NRHS*NB)
483 *
484  CALL zunmqr( 'L', 'C', m, nrhs, n, a, lda, work( itau ), b,
485  $ ldb, work( nwork ), lwork-nwork+1, info )
486 *
487 * Zero out below R.
488 *
489  IF( n.GT.1 ) THEN
490  CALL zlaset( 'L', n-1, n-1, czero, czero, a( 2, 1 ),
491  $ lda )
492  END IF
493  END IF
494 *
495  itauq = 1
496  itaup = itauq + n
497  nwork = itaup + n
498  ie = 1
499  nrwork = ie + n
500 *
501 * Bidiagonalize R in A.
502 * (RWorkspace: need N)
503 * (CWorkspace: need 2*N+MM, prefer 2*N+(MM+N)*NB)
504 *
505  CALL zgebrd( mm, n, a, lda, s, rwork( ie ), work( itauq ),
506  $ work( itaup ), work( nwork ), lwork-nwork+1,
507  $ info )
508 *
509 * Multiply B by transpose of left bidiagonalizing vectors of R.
510 * (CWorkspace: need 2*N+NRHS, prefer 2*N+NRHS*NB)
511 *
512  CALL zunmbr( 'Q', 'L', 'C', mm, nrhs, n, a, lda, work( itauq ),
513  $ b, ldb, work( nwork ), lwork-nwork+1, info )
514 *
515 * Solve the bidiagonal least squares problem.
516 *
517  CALL zlalsd( 'U', smlsiz, n, nrhs, s, rwork( ie ), b, ldb,
518  $ rcond, rank, work( nwork ), rwork( nrwork ),
519  $ iwork, info )
520  IF( info.NE.0 ) THEN
521  GO TO 10
522  END IF
523 *
524 * Multiply B by right bidiagonalizing vectors of R.
525 *
526  CALL zunmbr( 'P', 'L', 'N', n, nrhs, n, a, lda, work( itaup ),
527  $ b, ldb, work( nwork ), lwork-nwork+1, info )
528 *
529  ELSE IF( n.GE.mnthr .AND. lwork.GE.4*m+m*m+
530  $ max( m, 2*m-4, nrhs, n-3*m ) ) THEN
531 *
532 * Path 2a - underdetermined, with many more columns than rows
533 * and sufficient workspace for an efficient algorithm.
534 *
535  ldwork = m
536  IF( lwork.GE.max( 4*m+m*lda+max( m, 2*m-4, nrhs, n-3*m ),
537  $ m*lda+m+m*nrhs ) )ldwork = lda
538  itau = 1
539  nwork = m + 1
540 *
541 * Compute A=L*Q.
542 * (CWorkspace: need 2*M, prefer M+M*NB)
543 *
544  CALL zgelqf( m, n, a, lda, work( itau ), work( nwork ),
545  $ lwork-nwork+1, info )
546  il = nwork
547 *
548 * Copy L to WORK(IL), zeroing out above its diagonal.
549 *
550  CALL zlacpy( 'L', m, m, a, lda, work( il ), ldwork )
551  CALL zlaset( 'U', m-1, m-1, czero, czero, work( il+ldwork ),
552  $ ldwork )
553  itauq = il + ldwork*m
554  itaup = itauq + m
555  nwork = itaup + m
556  ie = 1
557  nrwork = ie + m
558 *
559 * Bidiagonalize L in WORK(IL).
560 * (RWorkspace: need M)
561 * (CWorkspace: need M*M+4*M, prefer M*M+4*M+2*M*NB)
562 *
563  CALL zgebrd( m, m, work( il ), ldwork, s, rwork( ie ),
564  $ work( itauq ), work( itaup ), work( nwork ),
565  $ lwork-nwork+1, info )
566 *
567 * Multiply B by transpose of left bidiagonalizing vectors of L.
568 * (CWorkspace: need M*M+4*M+NRHS, prefer M*M+4*M+NRHS*NB)
569 *
570  CALL zunmbr( 'Q', 'L', 'C', m, nrhs, m, work( il ), ldwork,
571  $ work( itauq ), b, ldb, work( nwork ),
572  $ lwork-nwork+1, info )
573 *
574 * Solve the bidiagonal least squares problem.
575 *
576  CALL zlalsd( 'U', smlsiz, m, nrhs, s, rwork( ie ), b, ldb,
577  $ rcond, rank, work( nwork ), rwork( nrwork ),
578  $ iwork, info )
579  IF( info.NE.0 ) THEN
580  GO TO 10
581  END IF
582 *
583 * Multiply B by right bidiagonalizing vectors of L.
584 *
585  CALL zunmbr( 'P', 'L', 'N', m, nrhs, m, work( il ), ldwork,
586  $ work( itaup ), b, ldb, work( nwork ),
587  $ lwork-nwork+1, info )
588 *
589 * Zero out below first M rows of B.
590 *
591  CALL zlaset( 'F', n-m, nrhs, czero, czero, b( m+1, 1 ), ldb )
592  nwork = itau + m
593 *
594 * Multiply transpose(Q) by B.
595 * (CWorkspace: need NRHS, prefer NRHS*NB)
596 *
597  CALL zunmlq( 'L', 'C', n, nrhs, m, a, lda, work( itau ), b,
598  $ ldb, work( nwork ), lwork-nwork+1, info )
599 *
600  ELSE
601 *
602 * Path 2 - remaining underdetermined cases.
603 *
604  itauq = 1
605  itaup = itauq + m
606  nwork = itaup + m
607  ie = 1
608  nrwork = ie + m
609 *
610 * Bidiagonalize A.
611 * (RWorkspace: need M)
612 * (CWorkspace: need 2*M+N, prefer 2*M+(M+N)*NB)
613 *
614  CALL zgebrd( m, n, a, lda, s, rwork( ie ), work( itauq ),
615  $ work( itaup ), work( nwork ), lwork-nwork+1,
616  $ info )
617 *
618 * Multiply B by transpose of left bidiagonalizing vectors.
619 * (CWorkspace: need 2*M+NRHS, prefer 2*M+NRHS*NB)
620 *
621  CALL zunmbr( 'Q', 'L', 'C', m, nrhs, n, a, lda, work( itauq ),
622  $ b, ldb, work( nwork ), lwork-nwork+1, info )
623 *
624 * Solve the bidiagonal least squares problem.
625 *
626  CALL zlalsd( 'L', smlsiz, m, nrhs, s, rwork( ie ), b, ldb,
627  $ rcond, rank, work( nwork ), rwork( nrwork ),
628  $ iwork, info )
629  IF( info.NE.0 ) THEN
630  GO TO 10
631  END IF
632 *
633 * Multiply B by right bidiagonalizing vectors of A.
634 *
635  CALL zunmbr( 'P', 'L', 'N', n, nrhs, m, a, lda, work( itaup ),
636  $ b, ldb, work( nwork ), lwork-nwork+1, info )
637 *
638  END IF
639 *
640 * Undo scaling.
641 *
642  IF( iascl.EQ.1 ) THEN
643  CALL zlascl( 'G', 0, 0, anrm, smlnum, n, nrhs, b, ldb, info )
644  CALL dlascl( 'G', 0, 0, smlnum, anrm, minmn, 1, s, minmn,
645  $ info )
646  ELSE IF( iascl.EQ.2 ) THEN
647  CALL zlascl( 'G', 0, 0, anrm, bignum, n, nrhs, b, ldb, info )
648  CALL dlascl( 'G', 0, 0, bignum, anrm, minmn, 1, s, minmn,
649  $ info )
650  END IF
651  IF( ibscl.EQ.1 ) THEN
652  CALL zlascl( 'G', 0, 0, smlnum, bnrm, n, nrhs, b, ldb, info )
653  ELSE IF( ibscl.EQ.2 ) THEN
654  CALL zlascl( 'G', 0, 0, bignum, bnrm, n, nrhs, b, ldb, info )
655  END IF
656 *
657  10 CONTINUE
658  work( 1 ) = maxwrk
659  iwork( 1 ) = liwork
660  rwork( 1 ) = lrwork
661  RETURN
662 *
663 * End of ZGELSD
664 *
subroutine dlaset(UPLO, M, N, ALPHA, BETA, A, LDA)
DLASET initializes the off-diagonal elements and the diagonal elements of a matrix to given values...
Definition: dlaset.f:112
subroutine zlacpy(UPLO, M, N, A, LDA, B, LDB)
ZLACPY copies all or part of one two-dimensional array to another.
Definition: zlacpy.f:105
double precision function dlamch(CMACH)
DLAMCH
Definition: dlamch.f:65
subroutine zgeqrf(M, N, A, LDA, TAU, WORK, LWORK, INFO)
ZGEQRF VARIANT: left-looking Level 3 BLAS of the algorithm.
Definition: zgeqrf.f:151
subroutine zunmbr(VECT, SIDE, TRANS, M, N, K, A, LDA, TAU, C, LDC, WORK, LWORK, INFO)
ZUNMBR
Definition: zunmbr.f:198
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:145
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:108
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine dlabad(SMALL, LARGE)
DLABAD
Definition: dlabad.f:76
subroutine zgebrd(M, N, A, LDA, D, E, TAUQ, TAUP, WORK, LWORK, INFO)
ZGEBRD
Definition: zgebrd.f:207
subroutine zunmqr(SIDE, TRANS, M, N, K, A, LDA, TAU, C, LDC, WORK, LWORK, INFO)
ZUNMQR
Definition: zunmqr.f:169
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:117
subroutine zunmlq(SIDE, TRANS, M, N, K, A, LDA, TAU, C, LDC, WORK, LWORK, INFO)
ZUNMLQ
Definition: zunmlq.f:169
integer function ilaenv(ISPEC, NAME, OPTS, N1, N2, N3, N4)
Definition: tstiee.f:83
subroutine zgelqf(M, N, A, LDA, TAU, WORK, LWORK, INFO)
ZGELQF
Definition: zgelqf.f:137
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:145
subroutine zlalsd(UPLO, SMLSIZ, N, NRHS, D, E, B, LDB, RCOND, RANK, WORK, RWORK, IWORK, INFO)
ZLALSD uses the singular value decomposition of A to solve the least squares problem.
Definition: zlalsd.f:190

Here is the call graph for this function:

Here is the caller graph for this function: