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