LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
zgelsd.f
Go to the documentation of this file.
1 *> \brief <b> ZGELSD computes the minimum-norm solution to a linear least squares problem for GE matrices</b>
2 *
3 * =========== DOCUMENTATION ===========
4 *
5 * Online html documentation available at
6 * http://www.netlib.org/lapack/explore-html/
7 *
8 *> \htmlonly
9 *> Download ZGELSD + dependencies
10 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/zgelsd.f">
11 *> [TGZ]</a>
12 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/zgelsd.f">
13 *> [ZIP]</a>
14 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/zgelsd.f">
15 *> [TXT]</a>
16 *> \endhtmlonly
17 *
18 * Definition:
19 * ===========
20 *
21 * SUBROUTINE ZGELSD( M, N, NRHS, A, LDA, B, LDB, S, RCOND, RANK,
22 * WORK, LWORK, RWORK, IWORK, INFO )
23 *
24 * .. Scalar Arguments ..
25 * INTEGER INFO, LDA, LDB, LWORK, M, N, NRHS, RANK
26 * DOUBLE PRECISION RCOND
27 * ..
28 * .. Array Arguments ..
29 * INTEGER IWORK( * )
30 * DOUBLE PRECISION RWORK( * ), S( * )
31 * COMPLEX*16 A( LDA, * ), B( LDB, * ), WORK( * )
32 * ..
33 *
34 *
35 *> \par Purpose:
36 * =============
37 *>
38 *> \verbatim
39 *>
40 *> ZGELSD computes the minimum-norm solution to a real linear least
41 *> squares problem:
42 *> minimize 2-norm(| b - A*x |)
43 *> using the singular value decomposition (SVD) of A. A is an M-by-N
44 *> matrix which may be rank-deficient.
45 *>
46 *> Several right hand side vectors b and solution vectors x can be
47 *> handled in a single call; they are stored as the columns of the
48 *> M-by-NRHS right hand side matrix B and the N-by-NRHS solution
49 *> matrix X.
50 *>
51 *> The problem is solved in three steps:
52 *> (1) Reduce the coefficient matrix A to bidiagonal form with
53 *> Householder tranformations, reducing the original problem
54 *> into a "bidiagonal least squares problem" (BLS)
55 *> (2) Solve the BLS using a divide and conquer approach.
56 *> (3) Apply back all the Householder tranformations to solve
57 *> the original least squares problem.
58 *>
59 *> The effective rank of A is determined by treating as zero those
60 *> singular values which are less than RCOND times the largest singular
61 *> value.
62 *>
63 *> The divide and conquer algorithm makes very mild assumptions about
64 *> floating point arithmetic. It will work on machines with a guard
65 *> digit in add/subtract, or on those binary machines without guard
66 *> digits which subtract like the Cray X-MP, Cray Y-MP, Cray C-90, or
67 *> Cray-2. It could conceivably fail on hexadecimal or decimal machines
68 *> without guard digits, but we know of none.
69 *> \endverbatim
70 *
71 * Arguments:
72 * ==========
73 *
74 *> \param[in] M
75 *> \verbatim
76 *> M is INTEGER
77 *> The number of rows of the matrix A. M >= 0.
78 *> \endverbatim
79 *>
80 *> \param[in] N
81 *> \verbatim
82 *> N is INTEGER
83 *> The number of columns of the matrix A. N >= 0.
84 *> \endverbatim
85 *>
86 *> \param[in] NRHS
87 *> \verbatim
88 *> NRHS is INTEGER
89 *> The number of right hand sides, i.e., the number of columns
90 *> of the matrices B and X. NRHS >= 0.
91 *> \endverbatim
92 *>
93 *> \param[in] A
94 *> \verbatim
95 *> A is COMPLEX*16 array, dimension (LDA,N)
96 *> On entry, the M-by-N matrix A.
97 *> On exit, A has been destroyed.
98 *> \endverbatim
99 *>
100 *> \param[in] LDA
101 *> \verbatim
102 *> LDA is INTEGER
103 *> The leading dimension of the array A. LDA >= max(1,M).
104 *> \endverbatim
105 *>
106 *> \param[in,out] B
107 *> \verbatim
108 *> B is COMPLEX*16 array, dimension (LDB,NRHS)
109 *> On entry, the M-by-NRHS right hand side matrix B.
110 *> On exit, B is overwritten by the N-by-NRHS solution matrix X.
111 *> If m >= n and RANK = n, the residual sum-of-squares for
112 *> the solution in the i-th column is given by the sum of
113 *> squares of the modulus of elements n+1:m in that column.
114 *> \endverbatim
115 *>
116 *> \param[in] LDB
117 *> \verbatim
118 *> LDB is INTEGER
119 *> The leading dimension of the array B. LDB >= max(1,M,N).
120 *> \endverbatim
121 *>
122 *> \param[out] S
123 *> \verbatim
124 *> S is DOUBLE PRECISION array, dimension (min(M,N))
125 *> The singular values of A in decreasing order.
126 *> The condition number of A in the 2-norm = S(1)/S(min(m,n)).
127 *> \endverbatim
128 *>
129 *> \param[in] RCOND
130 *> \verbatim
131 *> RCOND is DOUBLE PRECISION
132 *> RCOND is used to determine the effective rank of A.
133 *> Singular values S(i) <= RCOND*S(1) are treated as zero.
134 *> If RCOND < 0, machine precision is used instead.
135 *> \endverbatim
136 *>
137 *> \param[out] RANK
138 *> \verbatim
139 *> RANK is INTEGER
140 *> The effective rank of A, i.e., the number of singular values
141 *> which are greater than RCOND*S(1).
142 *> \endverbatim
143 *>
144 *> \param[out] WORK
145 *> \verbatim
146 *> WORK is COMPLEX*16 array, dimension (MAX(1,LWORK))
147 *> On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
148 *> \endverbatim
149 *>
150 *> \param[in] LWORK
151 *> \verbatim
152 *> LWORK is INTEGER
153 *> The dimension of the array WORK. LWORK must be at least 1.
154 *> The exact minimum amount of workspace needed depends on M,
155 *> N and NRHS. As long as LWORK is at least
156 *> 2*N + N*NRHS
157 *> if M is greater than or equal to N or
158 *> 2*M + M*NRHS
159 *> if M is less than N, the code will execute correctly.
160 *> For good performance, LWORK should generally be larger.
161 *>
162 *> If LWORK = -1, then a workspace query is assumed; the routine
163 *> only calculates the optimal size of the array WORK and the
164 *> minimum sizes of the arrays RWORK and IWORK, and returns
165 *> these values as the first entries of the WORK, RWORK and
166 *> IWORK arrays, and no error message related to LWORK is issued
167 *> by XERBLA.
168 *> \endverbatim
169 *>
170 *> \param[out] RWORK
171 *> \verbatim
172 *> RWORK is DOUBLE PRECISION array, dimension (MAX(1,LRWORK))
173 *> LRWORK >=
174 *> 10*N + 2*N*SMLSIZ + 8*N*NLVL + 3*SMLSIZ*NRHS +
175 *> MAX( (SMLSIZ+1)**2, N*(1+NRHS) + 2*NRHS )
176 *> if M is greater than or equal to N or
177 *> 10*M + 2*M*SMLSIZ + 8*M*NLVL + 3*SMLSIZ*NRHS +
178 *> MAX( (SMLSIZ+1)**2, N*(1+NRHS) + 2*NRHS )
179 *> if M is less than N, the code will execute correctly.
180 *> SMLSIZ is returned by ILAENV and is equal to the maximum
181 *> size of the subproblems at the bottom of the computation
182 *> tree (usually about 25), and
183 *> NLVL = MAX( 0, INT( LOG_2( MIN( M,N )/(SMLSIZ+1) ) ) + 1 )
184 *> On exit, if INFO = 0, RWORK(1) returns the minimum LRWORK.
185 *> \endverbatim
186 *>
187 *> \param[out] IWORK
188 *> \verbatim
189 *> IWORK is INTEGER array, dimension (MAX(1,LIWORK))
190 *> LIWORK >= max(1, 3*MINMN*NLVL + 11*MINMN),
191 *> where MINMN = MIN( M,N ).
192 *> On exit, if INFO = 0, IWORK(1) returns the minimum LIWORK.
193 *> \endverbatim
194 *>
195 *> \param[out] INFO
196 *> \verbatim
197 *> INFO is INTEGER
198 *> = 0: successful exit
199 *> < 0: if INFO = -i, the i-th argument had an illegal value.
200 *> > 0: the algorithm for computing the SVD failed to converge;
201 *> if INFO = i, i off-diagonal elements of an intermediate
202 *> bidiagonal form did not converge to zero.
203 *> \endverbatim
204 *
205 * Authors:
206 * ========
207 *
208 *> \author Univ. of Tennessee
209 *> \author Univ. of California Berkeley
210 *> \author Univ. of Colorado Denver
211 *> \author NAG Ltd.
212 *
213 *> \date November 2011
214 *
215 *> \ingroup complex16GEsolve
216 *
217 *> \par Contributors:
218 * ==================
219 *>
220 *> Ming Gu and Ren-Cang Li, Computer Science Division, University of
221 *> California at Berkeley, USA \n
222 *> Osni Marques, LBNL/NERSC, USA \n
223 *
224 * =====================================================================
225  SUBROUTINE zgelsd( M, N, NRHS, A, LDA, B, LDB, S, RCOND, RANK,
226  $ work, lwork, rwork, iwork, info )
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 *
665  END
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
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 zgelsd(M, N, NRHS, A, LDA, B, LDB, S, RCOND, RANK, WORK, LWORK, RWORK, IWORK, INFO)
ZGELSD computes the minimum-norm solution to a linear least squares problem for GE matrices ...
Definition: zgelsd.f:227
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
subroutine zunmlq(SIDE, TRANS, M, N, K, A, LDA, TAU, C, LDC, WORK, LWORK, INFO)
ZUNMLQ
Definition: zunmlq.f:169
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