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