LAPACK 3.12.1
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
cgetsls.f
Go to the documentation of this file.
1*> \brief \b CGETSLS
2*
3* Definition:
4* ===========
5*
6* SUBROUTINE CGETSLS( TRANS, M, N, NRHS, A, LDA, B, LDB,
7* $ WORK, LWORK, INFO )
8*
9* .. Scalar Arguments ..
10* CHARACTER TRANS
11* INTEGER INFO, LDA, LDB, LWORK, M, N, NRHS
12* ..
13* .. Array Arguments ..
14* COMPLEX A( LDA, * ), B( LDB, * ), WORK( * )
15* ..
16*
17*
18*> \par Purpose:
19* =============
20*>
21*> \verbatim
22*>
23*> CGETSLS solves overdetermined or underdetermined complex linear systems
24*> involving an M-by-N matrix A, using a tall skinny QR or short wide LQ
25*> factorization of A.
26*>
27*> It is assumed that A has full rank, and only a rudimentary protection
28*> against rank-deficient matrices is provided. This subroutine only detects
29*> exact rank-deficiency, where a diagonal element of the triangular factor
30*> of A is exactly zero.
31*>
32*> It is conceivable for one (or more) of the diagonal elements of the triangular
33*> factor of A to be subnormally tiny numbers without this subroutine signalling
34*> an error. The solutions computed for such almost-rank-deficient matrices may
35*> be less accurate due to a loss of numerical precision.
36*>
37*>
38*> The following options are provided:
39*>
40*> 1. If TRANS = 'N' and m >= n: find the least squares solution of
41*> an overdetermined system, i.e., solve the least squares problem
42*> minimize || B - A*X ||.
43*>
44*> 2. If TRANS = 'N' and m < n: find the minimum norm solution of
45*> an underdetermined system A * X = B.
46*>
47*> 3. If TRANS = 'C' and m >= n: find the minimum norm solution of
48*> an undetermined system A**T * X = B.
49*>
50*> 4. If TRANS = 'C' and m < n: find the least squares solution of
51*> an overdetermined system, i.e., solve the least squares problem
52*> minimize || B - A**T * X ||.
53*>
54*> Several right hand side vectors b and solution vectors x can be
55*> handled in a single call; they are stored as the columns of the
56*> M-by-NRHS right hand side matrix B and the N-by-NRHS solution
57*> matrix X.
58*> \endverbatim
59*
60* Arguments:
61* ==========
62*
63*> \param[in] TRANS
64*> \verbatim
65*> TRANS is CHARACTER*1
66*> = 'N': the linear system involves A;
67*> = 'C': the linear system involves A**H.
68*> \endverbatim
69*>
70*> \param[in] M
71*> \verbatim
72*> M is INTEGER
73*> The number of rows of the matrix A. M >= 0.
74*> \endverbatim
75*>
76*> \param[in] N
77*> \verbatim
78*> N is INTEGER
79*> The number of columns of the matrix A. N >= 0.
80*> \endverbatim
81*>
82*> \param[in] NRHS
83*> \verbatim
84*> NRHS is INTEGER
85*> The number of right hand sides, i.e., the number of
86*> columns of the matrices B and X. NRHS >=0.
87*> \endverbatim
88*>
89*> \param[in,out] A
90*> \verbatim
91*> A is COMPLEX array, dimension (LDA,N)
92*> On entry, the M-by-N matrix A.
93*> On exit,
94*> A is overwritten by details of its QR or LQ
95*> factorization as returned by CGEQR or CGELQ.
96*> \endverbatim
97*>
98*> \param[in] LDA
99*> \verbatim
100*> LDA is INTEGER
101*> The leading dimension of the array A. LDA >= max(1,M).
102*> \endverbatim
103*>
104*> \param[in,out] B
105*> \verbatim
106*> B is COMPLEX array, dimension (LDB,NRHS)
107*> On entry, the matrix B of right hand side vectors, stored
108*> columnwise; B is M-by-NRHS if TRANS = 'N', or N-by-NRHS
109*> if TRANS = 'C'.
110*> On exit, if INFO = 0, B is overwritten by the solution
111*> vectors, stored columnwise:
112*> if TRANS = 'N' and m >= n, rows 1 to n of B contain the least
113*> squares solution vectors.
114*> if TRANS = 'N' and m < n, rows 1 to N of B contain the
115*> minimum norm solution vectors;
116*> if TRANS = 'C' and m >= n, rows 1 to M of B contain the
117*> minimum norm solution vectors;
118*> if TRANS = 'C' and m < n, rows 1 to M of B contain the
119*> least squares solution vectors.
120*> \endverbatim
121*>
122*> \param[in] LDB
123*> \verbatim
124*> LDB is INTEGER
125*> The leading dimension of the array B. LDB >= MAX(1,M,N).
126*> \endverbatim
127*>
128*> \param[out] WORK
129*> \verbatim
130*> (workspace) COMPLEX array, dimension (MAX(1,LWORK))
131*> On exit, if INFO = 0, WORK(1) contains optimal (or either minimal
132*> or optimal, if query was assumed) LWORK.
133*> See LWORK for details.
134*> \endverbatim
135*>
136*> \param[in] LWORK
137*> \verbatim
138*> LWORK is INTEGER
139*> The dimension of the array WORK. LWORK >= 1.
140*> If LWORK = -1 or -2, then a workspace query is assumed.
141*> If LWORK = -1, the routine calculates optimal size of WORK for the
142*> optimal performance and returns this value in WORK(1).
143*> If LWORK = -2, the routine calculates minimal size of WORK and
144*> returns this value in WORK(1).
145*> \endverbatim
146*>
147*> \param[out] INFO
148*> \verbatim
149*> INFO is INTEGER
150*> = 0: successful exit
151*> < 0: if INFO = -i, the i-th argument had an illegal value
152*> > 0: if INFO = i, the i-th diagonal element of the
153*> triangular factor of A is exactly zero, so that A does not have
154*> full rank; the least squares solution could not be
155*> computed.
156*> \endverbatim
157*
158* Authors:
159* ========
160*
161*> \author Univ. of Tennessee
162*> \author Univ. of California Berkeley
163*> \author Univ. of Colorado Denver
164*> \author NAG Ltd.
165*
166*> \ingroup getsls
167*
168* =====================================================================
169 SUBROUTINE cgetsls( TRANS, M, N, NRHS, A, LDA, B, LDB,
170 $ WORK, LWORK, INFO )
171*
172* -- LAPACK driver routine --
173* -- LAPACK is a software package provided by Univ. of Tennessee, --
174* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
175*
176* .. Scalar Arguments ..
177 CHARACTER TRANS
178 INTEGER INFO, LDA, LDB, LWORK, M, N, NRHS
179* ..
180* .. Array Arguments ..
181 COMPLEX A( LDA, * ), B( LDB, * ), WORK( * )
182*
183* ..
184*
185* =====================================================================
186*
187* .. Parameters ..
188 REAL ZERO, ONE
189 parameter( zero = 0.0e0, one = 1.0e0 )
190 COMPLEX CZERO
191 parameter( czero = ( 0.0e+0, 0.0e+0 ) )
192* ..
193* .. Local Scalars ..
194 LOGICAL LQUERY, TRAN
195 INTEGER I, IASCL, IBSCL, J, MAXMN, BROW,
196 $ scllen, tszo, tszm, lwo, lwm, lw1, lw2,
197 $ wsizeo, wsizem, info2
198 REAL ANRM, BIGNUM, BNRM, SMLNUM, DUM( 1 )
199 COMPLEX TQ( 5 ), WORKQ( 1 )
200* ..
201* .. External Functions ..
202 LOGICAL LSAME
203 REAL SLAMCH, CLANGE, SROUNDUP_LWORK
204 EXTERNAL lsame, slamch, clange,
205 $ sroundup_lwork
206* ..
207* .. External Subroutines ..
208 EXTERNAL cgeqr, cgemqr, clascl, claset,
210* ..
211* .. Intrinsic Functions ..
212 INTRINSIC max, min, int
213* ..
214* .. Executable Statements ..
215*
216* Test the input arguments.
217*
218 info = 0
219 maxmn = max( m, n )
220 tran = lsame( trans, 'C' )
221*
222 lquery = ( lwork.EQ.-1 .OR. lwork.EQ.-2 )
223 IF( .NOT.( lsame( trans, 'N' ) .OR.
224 $ lsame( trans, 'C' ) ) ) THEN
225 info = -1
226 ELSE IF( m.LT.0 ) THEN
227 info = -2
228 ELSE IF( n.LT.0 ) THEN
229 info = -3
230 ELSE IF( nrhs.LT.0 ) THEN
231 info = -4
232 ELSE IF( lda.LT.max( 1, m ) ) THEN
233 info = -6
234 ELSE IF( ldb.LT.max( 1, m, n ) ) THEN
235 info = -8
236 END IF
237*
238 IF( info.EQ.0 ) THEN
239*
240* Determine the optimum and minimum LWORK
241*
242 IF( min( m, n, nrhs ).EQ.0 ) THEN
243 wsizeo = 1
244 wsizem = 1
245 ELSE IF ( m.GE.n ) THEN
246 CALL cgeqr( m, n, a, lda, tq, -1, workq, -1, info2 )
247 tszo = int( tq( 1 ) )
248 lwo = int( workq( 1 ) )
249 CALL cgemqr( 'L', trans, m, nrhs, n, a, lda, tq,
250 $ tszo, b, ldb, workq, -1, info2 )
251 lwo = max( lwo, int( workq( 1 ) ) )
252 CALL cgeqr( m, n, a, lda, tq, -2, workq, -2, info2 )
253 tszm = int( tq( 1 ) )
254 lwm = int( workq( 1 ) )
255 CALL cgemqr( 'L', trans, m, nrhs, n, a, lda, tq,
256 $ tszm, b, ldb, workq, -1, info2 )
257 lwm = max( lwm, int( workq( 1 ) ) )
258 wsizeo = tszo + lwo
259 wsizem = tszm + lwm
260 ELSE
261 CALL cgelq( m, n, a, lda, tq, -1, workq, -1, info2 )
262 tszo = int( tq( 1 ) )
263 lwo = int( workq( 1 ) )
264 CALL cgemlq( 'L', trans, n, nrhs, m, a, lda, tq,
265 $ tszo, b, ldb, workq, -1, info2 )
266 lwo = max( lwo, int( workq( 1 ) ) )
267 CALL cgelq( m, n, a, lda, tq, -2, workq, -2, info2 )
268 tszm = int( tq( 1 ) )
269 lwm = int( workq( 1 ) )
270 CALL cgemlq( 'L', trans, n, nrhs, m, a, lda, tq,
271 $ tszm, b, ldb, workq, -1, info2 )
272 lwm = max( lwm, int( workq( 1 ) ) )
273 wsizeo = tszo + lwo
274 wsizem = tszm + lwm
275 END IF
276*
277 IF( ( lwork.LT.wsizem ).AND.( .NOT.lquery ) ) THEN
278 info = -10
279 END IF
280*
281 work( 1 ) = sroundup_lwork( wsizeo )
282*
283 END IF
284*
285 IF( info.NE.0 ) THEN
286 CALL xerbla( 'CGETSLS', -info )
287 RETURN
288 END IF
289 IF( lquery ) THEN
290 IF( lwork.EQ.-2 ) work( 1 ) = sroundup_lwork( wsizem )
291 RETURN
292 END IF
293 IF( lwork.LT.wsizeo ) THEN
294 lw1 = tszm
295 lw2 = lwm
296 ELSE
297 lw1 = tszo
298 lw2 = lwo
299 END IF
300*
301* Quick return if possible
302*
303 IF( min( m, n, nrhs ).EQ.0 ) THEN
304 CALL claset( 'FULL', max( m, n ), nrhs, czero, czero,
305 $ b, ldb )
306 RETURN
307 END IF
308*
309* Get machine parameters
310*
311 smlnum = slamch( 'S' ) / slamch( 'P' )
312 bignum = one / smlnum
313*
314* Scale A, B if max element outside range [SMLNUM,BIGNUM]
315*
316 anrm = clange( 'M', m, n, a, lda, dum )
317 iascl = 0
318 IF( anrm.GT.zero .AND. anrm.LT.smlnum ) THEN
319*
320* Scale matrix norm up to SMLNUM
321*
322 CALL clascl( 'G', 0, 0, anrm, smlnum, m, n, a, lda, info )
323 iascl = 1
324 ELSE IF( anrm.GT.bignum ) THEN
325*
326* Scale matrix norm down to BIGNUM
327*
328 CALL clascl( 'G', 0, 0, anrm, bignum, m, n, a, lda, info )
329 iascl = 2
330 ELSE IF( anrm.EQ.zero ) THEN
331*
332* Matrix all zero. Return zero solution.
333*
334 CALL claset( 'F', maxmn, nrhs, czero, czero, b, ldb )
335 GO TO 50
336 END IF
337*
338 brow = m
339 IF ( tran ) THEN
340 brow = n
341 END IF
342 bnrm = clange( 'M', brow, nrhs, b, ldb, dum )
343 ibscl = 0
344 IF( bnrm.GT.zero .AND. bnrm.LT.smlnum ) THEN
345*
346* Scale matrix norm up to SMLNUM
347*
348 CALL clascl( 'G', 0, 0, bnrm, smlnum, brow, nrhs, b, ldb,
349 $ info )
350 ibscl = 1
351 ELSE IF( bnrm.GT.bignum ) THEN
352*
353* Scale matrix norm down to BIGNUM
354*
355 CALL clascl( 'G', 0, 0, bnrm, bignum, brow, nrhs, b, ldb,
356 $ info )
357 ibscl = 2
358 END IF
359*
360 IF ( m.GE.n ) THEN
361*
362* compute QR factorization of A
363*
364 CALL cgeqr( m, n, a, lda, work( lw2+1 ), lw1,
365 $ work( 1 ), lw2, info )
366 IF ( .NOT.tran ) THEN
367*
368* Least-Squares Problem min || A * X - B ||
369*
370* B(1:M,1:NRHS) := Q**T * B(1:M,1:NRHS)
371*
372 CALL cgemqr( 'L' , 'C', m, nrhs, n, a, lda,
373 $ work( lw2+1 ), lw1, b, ldb, work( 1 ), lw2,
374 $ info )
375*
376* B(1:N,1:NRHS) := inv(R) * B(1:N,1:NRHS)
377*
378 CALL ctrtrs( 'U', 'N', 'N', n, nrhs,
379 $ a, lda, b, ldb, info )
380 IF( info.GT.0 ) THEN
381 RETURN
382 END IF
383 scllen = n
384 ELSE
385*
386* Overdetermined system of equations A**T * X = B
387*
388* B(1:N,1:NRHS) := inv(R**T) * B(1:N,1:NRHS)
389*
390 CALL ctrtrs( 'U', 'C', 'N', n, nrhs,
391 $ a, lda, b, ldb, info )
392*
393 IF( info.GT.0 ) THEN
394 RETURN
395 END IF
396*
397* B(N+1:M,1:NRHS) = CZERO
398*
399 DO 20 j = 1, nrhs
400 DO 10 i = n + 1, m
401 b( i, j ) = czero
402 10 CONTINUE
403 20 CONTINUE
404*
405* B(1:M,1:NRHS) := Q(1:N,:) * B(1:N,1:NRHS)
406*
407 CALL cgemqr( 'L', 'N', m, nrhs, n, a, lda,
408 $ work( lw2+1 ), lw1, b, ldb, work( 1 ), lw2,
409 $ info )
410*
411 scllen = m
412*
413 END IF
414*
415 ELSE
416*
417* Compute LQ factorization of A
418*
419 CALL cgelq( m, n, a, lda, work( lw2+1 ), lw1,
420 $ work( 1 ), lw2, info )
421*
422* workspace at least M, optimally M*NB.
423*
424 IF( .NOT.tran ) THEN
425*
426* underdetermined system of equations A * X = B
427*
428* B(1:M,1:NRHS) := inv(L) * B(1:M,1:NRHS)
429*
430 CALL ctrtrs( 'L', 'N', 'N', m, nrhs,
431 $ a, lda, b, ldb, info )
432*
433 IF( info.GT.0 ) THEN
434 RETURN
435 END IF
436*
437* B(M+1:N,1:NRHS) = 0
438*
439 DO 40 j = 1, nrhs
440 DO 30 i = m + 1, n
441 b( i, j ) = czero
442 30 CONTINUE
443 40 CONTINUE
444*
445* B(1:N,1:NRHS) := Q(1:N,:)**T * B(1:M,1:NRHS)
446*
447 CALL cgemlq( 'L', 'C', n, nrhs, m, a, lda,
448 $ work( lw2+1 ), lw1, b, ldb, work( 1 ), lw2,
449 $ info )
450*
451* workspace at least NRHS, optimally NRHS*NB
452*
453 scllen = n
454*
455 ELSE
456*
457* overdetermined system min || A**T * X - B ||
458*
459* B(1:N,1:NRHS) := Q * B(1:N,1:NRHS)
460*
461 CALL cgemlq( 'L', 'N', n, nrhs, m, a, lda,
462 $ work( lw2+1 ), lw1, b, ldb, work( 1 ), lw2,
463 $ info )
464*
465* workspace at least NRHS, optimally NRHS*NB
466*
467* B(1:M,1:NRHS) := inv(L**T) * B(1:M,1:NRHS)
468*
469 CALL ctrtrs( 'L', 'C', 'N', m, nrhs,
470 $ a, lda, b, ldb, info )
471*
472 IF( info.GT.0 ) THEN
473 RETURN
474 END IF
475*
476 scllen = m
477*
478 END IF
479*
480 END IF
481*
482* Undo scaling
483*
484 IF( iascl.EQ.1 ) THEN
485 CALL clascl( 'G', 0, 0, anrm, smlnum, scllen, nrhs, b, ldb,
486 $ info )
487 ELSE IF( iascl.EQ.2 ) THEN
488 CALL clascl( 'G', 0, 0, anrm, bignum, scllen, nrhs, b, ldb,
489 $ info )
490 END IF
491 IF( ibscl.EQ.1 ) THEN
492 CALL clascl( 'G', 0, 0, smlnum, bnrm, scllen, nrhs, b, ldb,
493 $ info )
494 ELSE IF( ibscl.EQ.2 ) THEN
495 CALL clascl( 'G', 0, 0, bignum, bnrm, scllen, nrhs, b, ldb,
496 $ info )
497 END IF
498*
499 50 CONTINUE
500 work( 1 ) = sroundup_lwork( tszo + lwo )
501 RETURN
502*
503* End of CGETSLS
504*
505 END
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine cgelq(m, n, a, lda, t, tsize, work, lwork, info)
CGELQ
Definition cgelq.f:174
subroutine cgemlq(side, trans, m, n, k, a, lda, t, tsize, c, ldc, work, lwork, info)
CGEMLQ
Definition cgemlq.f:173
subroutine cgemqr(side, trans, m, n, k, a, lda, t, tsize, c, ldc, work, lwork, info)
CGEMQR
Definition cgemqr.f:175
subroutine cgeqr(m, n, a, lda, t, tsize, work, lwork, info)
CGEQR
Definition cgeqr.f:176
subroutine cgetsls(trans, m, n, nrhs, a, lda, b, ldb, work, lwork, info)
CGETSLS
Definition cgetsls.f:171
subroutine clascl(type, kl, ku, cfrom, cto, m, n, a, lda, info)
CLASCL multiplies a general rectangular matrix by a real scalar defined as cto/cfrom.
Definition clascl.f:142
subroutine claset(uplo, m, n, alpha, beta, a, lda)
CLASET initializes the off-diagonal elements and the diagonal elements of a matrix to given values.
Definition claset.f:104
subroutine ctrtrs(uplo, trans, diag, n, nrhs, a, lda, b, ldb, info)
CTRTRS
Definition ctrtrs.f:144