LAPACK  3.4.2
LAPACK: Linear Algebra PACKage
 All Files Functions Groups
cgelsd.f
Go to the documentation of this file.
1 *> \brief <b> CGELSD 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 CGELSD + dependencies
10 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/cgelsd.f">
11 *> [TGZ]</a>
12 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/cgelsd.f">
13 *> [ZIP]</a>
14 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/cgelsd.f">
15 *> [TXT]</a>
16 *> \endhtmlonly
17 *
18 * Definition:
19 * ===========
20 *
21 * SUBROUTINE CGELSD( 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 * REAL RCOND
27 * ..
28 * .. Array Arguments ..
29 * INTEGER IWORK( * )
30 * REAL RWORK( * ), S( * )
31 * COMPLEX A( LDA, * ), B( LDB, * ), WORK( * )
32 * ..
33 *
34 *
35 *> \par Purpose:
36 * =============
37 *>
38 *> \verbatim
39 *>
40 *> CGELSD 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,out] A
94 *> \verbatim
95 *> A is COMPLEX 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 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 REAL 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 REAL
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 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 REAL 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 complexGEsolve
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 cgelsd( 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  REAL rcond
236 * ..
237 * .. Array Arguments ..
238  INTEGER iwork( * )
239  REAL rwork( * ), s( * )
240  COMPLEX a( lda, * ), b( ldb, * ), work( * )
241 * ..
242 *
243 * =====================================================================
244 *
245 * .. Parameters ..
246  REAL zero, one, two
247  parameter( zero = 0.0e+0, one = 1.0e+0, two = 2.0e+0 )
248  COMPLEX czero
249  parameter( czero = ( 0.0e+0, 0.0e+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  REAL anrm, bignum, bnrm, eps, sfmin, smlnum
257 * ..
258 * .. External Subroutines ..
259  EXTERNAL cgebrd, cgelqf, cgeqrf, clacpy,
260  $ clalsd, clascl, claset, cunmbr,
261  $ cunmlq, cunmqr, slabad, slascl,
262  $ slaset, xerbla
263 * ..
264 * .. External Functions ..
265  INTEGER ilaenv
266  REAL clange, slamch
267  EXTERNAL clange, slamch, ilaenv
268 * ..
269 * .. Intrinsic Functions ..
270  INTRINSIC int, log, max, min, real
271 * ..
272 * .. Executable Statements ..
273 *
274 * Test the input arguments.
275 *
276  info = 0
277  minmn = min( m, n )
278  maxmn = max( m, n )
279  lquery = ( lwork.EQ.-1 )
280  IF( m.LT.0 ) THEN
281  info = -1
282  ELSE IF( n.LT.0 ) THEN
283  info = -2
284  ELSE IF( nrhs.LT.0 ) THEN
285  info = -3
286  ELSE IF( lda.LT.max( 1, m ) ) THEN
287  info = -5
288  ELSE IF( ldb.LT.max( 1, maxmn ) ) THEN
289  info = -7
290  END IF
291 *
292 * Compute workspace.
293 * (Note: Comments in the code beginning "Workspace:" describe the
294 * minimal amount of workspace needed at that point in the code,
295 * as well as the preferred amount for good performance.
296 * NB refers to the optimal block size for the immediately
297 * following subroutine, as returned by ILAENV.)
298 *
299  IF( info.EQ.0 ) THEN
300  minwrk = 1
301  maxwrk = 1
302  liwork = 1
303  lrwork = 1
304  IF( minmn.GT.0 ) THEN
305  smlsiz = ilaenv( 9, 'CGELSD', ' ', 0, 0, 0, 0 )
306  mnthr = ilaenv( 6, 'CGELSD', ' ', m, n, nrhs, -1 )
307  nlvl = max( int( log( REAL( MINMN ) / REAL( SMLSIZ + 1 ) ) /
308  $ log( two ) ) + 1, 0 )
309  liwork = 3*minmn*nlvl + 11*minmn
310  mm = m
311  IF( m.GE.n .AND. m.GE.mnthr ) THEN
312 *
313 * Path 1a - overdetermined, with many more rows than
314 * columns.
315 *
316  mm = n
317  maxwrk = max( maxwrk, n*ilaenv( 1, 'CGEQRF', ' ', m, n,
318  $ -1, -1 ) )
319  maxwrk = max( maxwrk, nrhs*ilaenv( 1, 'CUNMQR', 'LC', m,
320  $ nrhs, n, -1 ) )
321  END IF
322  IF( m.GE.n ) THEN
323 *
324 * Path 1 - overdetermined or exactly determined.
325 *
326  lrwork = 10*n + 2*n*smlsiz + 8*n*nlvl + 3*smlsiz*nrhs +
327  $ max( (smlsiz+1)**2, n*(1+nrhs) + 2*nrhs )
328  maxwrk = max( maxwrk, 2*n + ( mm + n )*ilaenv( 1,
329  $ 'CGEBRD', ' ', mm, n, -1, -1 ) )
330  maxwrk = max( maxwrk, 2*n + nrhs*ilaenv( 1, 'CUNMBR',
331  $ 'QLC', mm, nrhs, n, -1 ) )
332  maxwrk = max( maxwrk, 2*n + ( n - 1 )*ilaenv( 1,
333  $ 'CUNMBR', 'PLN', n, nrhs, n, -1 ) )
334  maxwrk = max( maxwrk, 2*n + n*nrhs )
335  minwrk = max( 2*n + mm, 2*n + n*nrhs )
336  END IF
337  IF( n.GT.m ) THEN
338  lrwork = 10*m + 2*m*smlsiz + 8*m*nlvl + 3*smlsiz*nrhs +
339  $ max( (smlsiz+1)**2, n*(1+nrhs) + 2*nrhs )
340  IF( n.GE.mnthr ) THEN
341 *
342 * Path 2a - underdetermined, with many more columns
343 * than rows.
344 *
345  maxwrk = m + m*ilaenv( 1, 'CGELQF', ' ', m, n, -1,
346  $ -1 )
347  maxwrk = max( maxwrk, m*m + 4*m + 2*m*ilaenv( 1,
348  $ 'CGEBRD', ' ', m, m, -1, -1 ) )
349  maxwrk = max( maxwrk, m*m + 4*m + nrhs*ilaenv( 1,
350  $ 'CUNMBR', 'QLC', m, nrhs, m, -1 ) )
351  maxwrk = max( maxwrk, m*m + 4*m + ( m - 1 )*ilaenv( 1,
352  $ 'CUNMLQ', 'LC', n, nrhs, m, -1 ) )
353  IF( nrhs.GT.1 ) THEN
354  maxwrk = max( maxwrk, m*m + m + m*nrhs )
355  ELSE
356  maxwrk = max( maxwrk, m*m + 2*m )
357  END IF
358  maxwrk = max( maxwrk, m*m + 4*m + m*nrhs )
359 ! XXX: Ensure the Path 2a case below is triggered. The workspace
360 ! calculation should use queries for all routines eventually.
361  maxwrk = max( maxwrk,
362  $ 4*m+m*m+max( m, 2*m-4, nrhs, n-3*m ) )
363  ELSE
364 *
365 * Path 2 - underdetermined.
366 *
367  maxwrk = 2*m + ( n + m )*ilaenv( 1, 'CGEBRD', ' ', m,
368  $ n, -1, -1 )
369  maxwrk = max( maxwrk, 2*m + nrhs*ilaenv( 1, 'CUNMBR',
370  $ 'QLC', m, nrhs, m, -1 ) )
371  maxwrk = max( maxwrk, 2*m + m*ilaenv( 1, 'CUNMBR',
372  $ 'PLN', n, nrhs, m, -1 ) )
373  maxwrk = max( maxwrk, 2*m + m*nrhs )
374  END IF
375  minwrk = max( 2*m + n, 2*m + m*nrhs )
376  END IF
377  END IF
378  minwrk = min( minwrk, maxwrk )
379  work( 1 ) = maxwrk
380  iwork( 1 ) = liwork
381  rwork( 1 ) = lrwork
382 *
383  IF( lwork.LT.minwrk .AND. .NOT.lquery ) THEN
384  info = -12
385  END IF
386  END IF
387 *
388  IF( info.NE.0 ) THEN
389  CALL xerbla( 'CGELSD', -info )
390  return
391  ELSE IF( lquery ) THEN
392  return
393  END IF
394 *
395 * Quick return if possible.
396 *
397  IF( m.EQ.0 .OR. n.EQ.0 ) THEN
398  rank = 0
399  return
400  END IF
401 *
402 * Get machine parameters.
403 *
404  eps = slamch( 'P' )
405  sfmin = slamch( 'S' )
406  smlnum = sfmin / eps
407  bignum = one / smlnum
408  CALL slabad( smlnum, bignum )
409 *
410 * Scale A if max entry outside range [SMLNUM,BIGNUM].
411 *
412  anrm = clange( 'M', m, n, a, lda, rwork )
413  iascl = 0
414  IF( anrm.GT.zero .AND. anrm.LT.smlnum ) THEN
415 *
416 * Scale matrix norm up to SMLNUM
417 *
418  CALL clascl( 'G', 0, 0, anrm, smlnum, m, n, a, lda, info )
419  iascl = 1
420  ELSE IF( anrm.GT.bignum ) THEN
421 *
422 * Scale matrix norm down to BIGNUM.
423 *
424  CALL clascl( 'G', 0, 0, anrm, bignum, m, n, a, lda, info )
425  iascl = 2
426  ELSE IF( anrm.EQ.zero ) THEN
427 *
428 * Matrix all zero. Return zero solution.
429 *
430  CALL claset( 'F', max( m, n ), nrhs, czero, czero, b, ldb )
431  CALL slaset( 'F', minmn, 1, zero, zero, s, 1 )
432  rank = 0
433  go to 10
434  END IF
435 *
436 * Scale B if max entry outside range [SMLNUM,BIGNUM].
437 *
438  bnrm = clange( 'M', m, nrhs, b, ldb, rwork )
439  ibscl = 0
440  IF( bnrm.GT.zero .AND. bnrm.LT.smlnum ) THEN
441 *
442 * Scale matrix norm up to SMLNUM.
443 *
444  CALL clascl( 'G', 0, 0, bnrm, smlnum, m, nrhs, b, ldb, info )
445  ibscl = 1
446  ELSE IF( bnrm.GT.bignum ) THEN
447 *
448 * Scale matrix norm down to BIGNUM.
449 *
450  CALL clascl( 'G', 0, 0, bnrm, bignum, m, nrhs, b, ldb, info )
451  ibscl = 2
452  END IF
453 *
454 * If M < N make sure B(M+1:N,:) = 0
455 *
456  IF( m.LT.n )
457  $ CALL claset( 'F', n-m, nrhs, czero, czero, b( m+1, 1 ), ldb )
458 *
459 * Overdetermined case.
460 *
461  IF( m.GE.n ) THEN
462 *
463 * Path 1 - overdetermined or exactly determined.
464 *
465  mm = m
466  IF( m.GE.mnthr ) THEN
467 *
468 * Path 1a - overdetermined, with many more rows than columns
469 *
470  mm = n
471  itau = 1
472  nwork = itau + n
473 *
474 * Compute A=Q*R.
475 * (RWorkspace: need N)
476 * (CWorkspace: need N, prefer N*NB)
477 *
478  CALL cgeqrf( m, n, a, lda, work( itau ), work( nwork ),
479  $ lwork-nwork+1, info )
480 *
481 * Multiply B by transpose(Q).
482 * (RWorkspace: need N)
483 * (CWorkspace: need NRHS, prefer NRHS*NB)
484 *
485  CALL cunmqr( 'L', 'C', m, nrhs, n, a, lda, work( itau ), b,
486  $ ldb, work( nwork ), lwork-nwork+1, info )
487 *
488 * Zero out below R.
489 *
490  IF( n.GT.1 ) THEN
491  CALL claset( 'L', n-1, n-1, czero, czero, a( 2, 1 ),
492  $ lda )
493  END IF
494  END IF
495 *
496  itauq = 1
497  itaup = itauq + n
498  nwork = itaup + n
499  ie = 1
500  nrwork = ie + n
501 *
502 * Bidiagonalize R in A.
503 * (RWorkspace: need N)
504 * (CWorkspace: need 2*N+MM, prefer 2*N+(MM+N)*NB)
505 *
506  CALL cgebrd( mm, n, a, lda, s, rwork( ie ), work( itauq ),
507  $ work( itaup ), work( nwork ), lwork-nwork+1,
508  $ info )
509 *
510 * Multiply B by transpose of left bidiagonalizing vectors of R.
511 * (CWorkspace: need 2*N+NRHS, prefer 2*N+NRHS*NB)
512 *
513  CALL cunmbr( 'Q', 'L', 'C', mm, nrhs, n, a, lda, work( itauq ),
514  $ b, ldb, work( nwork ), lwork-nwork+1, info )
515 *
516 * Solve the bidiagonal least squares problem.
517 *
518  CALL clalsd( 'U', smlsiz, n, nrhs, s, rwork( ie ), b, ldb,
519  $ rcond, rank, work( nwork ), rwork( nrwork ),
520  $ iwork, info )
521  IF( info.NE.0 ) THEN
522  go to 10
523  END IF
524 *
525 * Multiply B by right bidiagonalizing vectors of R.
526 *
527  CALL cunmbr( 'P', 'L', 'N', n, nrhs, n, a, lda, work( itaup ),
528  $ b, ldb, work( nwork ), lwork-nwork+1, info )
529 *
530  ELSE IF( n.GE.mnthr .AND. lwork.GE.4*m+m*m+
531  $ max( m, 2*m-4, nrhs, n-3*m ) ) THEN
532 *
533 * Path 2a - underdetermined, with many more columns than rows
534 * and sufficient workspace for an efficient algorithm.
535 *
536  ldwork = m
537  IF( lwork.GE.max( 4*m+m*lda+max( m, 2*m-4, nrhs, n-3*m ),
538  $ m*lda+m+m*nrhs ) )ldwork = lda
539  itau = 1
540  nwork = m + 1
541 *
542 * Compute A=L*Q.
543 * (CWorkspace: need 2*M, prefer M+M*NB)
544 *
545  CALL cgelqf( m, n, a, lda, work( itau ), work( nwork ),
546  $ lwork-nwork+1, info )
547  il = nwork
548 *
549 * Copy L to WORK(IL), zeroing out above its diagonal.
550 *
551  CALL clacpy( 'L', m, m, a, lda, work( il ), ldwork )
552  CALL claset( 'U', m-1, m-1, czero, czero, work( il+ldwork ),
553  $ ldwork )
554  itauq = il + ldwork*m
555  itaup = itauq + m
556  nwork = itaup + m
557  ie = 1
558  nrwork = ie + m
559 *
560 * Bidiagonalize L in WORK(IL).
561 * (RWorkspace: need M)
562 * (CWorkspace: need M*M+4*M, prefer M*M+4*M+2*M*NB)
563 *
564  CALL cgebrd( m, m, work( il ), ldwork, s, rwork( ie ),
565  $ work( itauq ), work( itaup ), work( nwork ),
566  $ lwork-nwork+1, info )
567 *
568 * Multiply B by transpose of left bidiagonalizing vectors of L.
569 * (CWorkspace: need M*M+4*M+NRHS, prefer M*M+4*M+NRHS*NB)
570 *
571  CALL cunmbr( 'Q', 'L', 'C', m, nrhs, m, work( il ), ldwork,
572  $ work( itauq ), b, ldb, work( nwork ),
573  $ lwork-nwork+1, info )
574 *
575 * Solve the bidiagonal least squares problem.
576 *
577  CALL clalsd( 'U', smlsiz, m, nrhs, s, rwork( ie ), b, ldb,
578  $ rcond, rank, work( nwork ), rwork( nrwork ),
579  $ iwork, info )
580  IF( info.NE.0 ) THEN
581  go to 10
582  END IF
583 *
584 * Multiply B by right bidiagonalizing vectors of L.
585 *
586  CALL cunmbr( 'P', 'L', 'N', m, nrhs, m, work( il ), ldwork,
587  $ work( itaup ), b, ldb, work( nwork ),
588  $ lwork-nwork+1, info )
589 *
590 * Zero out below first M rows of B.
591 *
592  CALL claset( 'F', n-m, nrhs, czero, czero, b( m+1, 1 ), ldb )
593  nwork = itau + m
594 *
595 * Multiply transpose(Q) by B.
596 * (CWorkspace: need NRHS, prefer NRHS*NB)
597 *
598  CALL cunmlq( 'L', 'C', n, nrhs, m, a, lda, work( itau ), b,
599  $ ldb, work( nwork ), lwork-nwork+1, info )
600 *
601  ELSE
602 *
603 * Path 2 - remaining underdetermined cases.
604 *
605  itauq = 1
606  itaup = itauq + m
607  nwork = itaup + m
608  ie = 1
609  nrwork = ie + m
610 *
611 * Bidiagonalize A.
612 * (RWorkspace: need M)
613 * (CWorkspace: need 2*M+N, prefer 2*M+(M+N)*NB)
614 *
615  CALL cgebrd( m, n, a, lda, s, rwork( ie ), work( itauq ),
616  $ work( itaup ), work( nwork ), lwork-nwork+1,
617  $ info )
618 *
619 * Multiply B by transpose of left bidiagonalizing vectors.
620 * (CWorkspace: need 2*M+NRHS, prefer 2*M+NRHS*NB)
621 *
622  CALL cunmbr( 'Q', 'L', 'C', m, nrhs, n, a, lda, work( itauq ),
623  $ b, ldb, work( nwork ), lwork-nwork+1, info )
624 *
625 * Solve the bidiagonal least squares problem.
626 *
627  CALL clalsd( 'L', smlsiz, m, nrhs, s, rwork( ie ), b, ldb,
628  $ rcond, rank, work( nwork ), rwork( nrwork ),
629  $ iwork, info )
630  IF( info.NE.0 ) THEN
631  go to 10
632  END IF
633 *
634 * Multiply B by right bidiagonalizing vectors of A.
635 *
636  CALL cunmbr( 'P', 'L', 'N', n, nrhs, m, a, lda, work( itaup ),
637  $ b, ldb, work( nwork ), lwork-nwork+1, info )
638 *
639  END IF
640 *
641 * Undo scaling.
642 *
643  IF( iascl.EQ.1 ) THEN
644  CALL clascl( 'G', 0, 0, anrm, smlnum, n, nrhs, b, ldb, info )
645  CALL slascl( 'G', 0, 0, smlnum, anrm, minmn, 1, s, minmn,
646  $ info )
647  ELSE IF( iascl.EQ.2 ) THEN
648  CALL clascl( 'G', 0, 0, anrm, bignum, n, nrhs, b, ldb, info )
649  CALL slascl( 'G', 0, 0, bignum, anrm, minmn, 1, s, minmn,
650  $ info )
651  END IF
652  IF( ibscl.EQ.1 ) THEN
653  CALL clascl( 'G', 0, 0, smlnum, bnrm, n, nrhs, b, ldb, info )
654  ELSE IF( ibscl.EQ.2 ) THEN
655  CALL clascl( 'G', 0, 0, bignum, bnrm, n, nrhs, b, ldb, info )
656  END IF
657 *
658  10 continue
659  work( 1 ) = maxwrk
660  iwork( 1 ) = liwork
661  rwork( 1 ) = lrwork
662  return
663 *
664 * End of CGELSD
665 *
666  END