LAPACK 3.12.1
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.
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 = 'T' and m >= n: find the minimum norm solution of
48*> an undetermined system A**T * X = B.
49*>
50*> 4. If TRANS = 'T' 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*> = 'T': the linear system involves A**T.
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 DOUBLE PRECISION 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 DGEQR or DGELQ.
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 DOUBLE PRECISION 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 = 'T'.
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 = 'T' and m >= n, rows 1 to M of B contain the
117*> minimum norm solution vectors;
118*> if TRANS = 'T' 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) DOUBLE PRECISION 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 dgetsls( 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 DOUBLE PRECISION A( LDA, * ), B( LDB, * ), WORK( * )
182*
183* ..
184*
185* =====================================================================
186*
187* .. Parameters ..
188 DOUBLE PRECISION ZERO, ONE
189 parameter( zero = 0.0d0, one = 1.0d0 )
190* ..
191* .. Local Scalars ..
192 LOGICAL LQUERY, TRAN
193 INTEGER I, IASCL, IBSCL, J, MAXMN, BROW,
194 $ scllen, tszo, tszm, lwo, lwm, lw1, lw2,
195 $ wsizeo, wsizem, info2
196 DOUBLE PRECISION ANRM, BIGNUM, BNRM, SMLNUM, TQ( 5 ), WORKQ( 1 )
197* ..
198* .. External Functions ..
199 LOGICAL LSAME
200 DOUBLE PRECISION DLAMCH, DLANGE
201 EXTERNAL lsame, dlamch, dlange
202* ..
203* .. External Subroutines ..
204 EXTERNAL dgeqr, dgemqr, dlascl, dlaset,
206* ..
207* .. Intrinsic Functions ..
208 INTRINSIC dble, max, min, int
209* ..
210* .. Executable Statements ..
211*
212* Test the input arguments.
213*
214 info = 0
215 maxmn = max( m, n )
216 tran = lsame( trans, 'T' )
217*
218 lquery = ( lwork.EQ.-1 .OR. lwork.EQ.-2 )
219 IF( .NOT.( lsame( trans, 'N' ) .OR.
220 $ lsame( trans, 'T' ) ) ) THEN
221 info = -1
222 ELSE IF( m.LT.0 ) THEN
223 info = -2
224 ELSE IF( n.LT.0 ) THEN
225 info = -3
226 ELSE IF( nrhs.LT.0 ) THEN
227 info = -4
228 ELSE IF( lda.LT.max( 1, m ) ) THEN
229 info = -6
230 ELSE IF( ldb.LT.max( 1, m, n ) ) THEN
231 info = -8
232 END IF
233*
234 IF( info.EQ.0 ) THEN
235*
236* Determine the optimum and minimum LWORK
237*
238 IF( min( m, n, nrhs ).EQ.0 ) THEN
239 wsizem = 1
240 wsizeo = 1
241 ELSE IF( m.GE.n ) THEN
242 CALL dgeqr( m, n, a, lda, tq, -1, workq, -1, info2 )
243 tszo = int( tq( 1 ) )
244 lwo = int( workq( 1 ) )
245 CALL dgemqr( 'L', trans, m, nrhs, n, a, lda, tq,
246 $ tszo, b, ldb, workq, -1, info2 )
247 lwo = max( lwo, int( workq( 1 ) ) )
248 CALL dgeqr( m, n, a, lda, tq, -2, workq, -2, info2 )
249 tszm = int( tq( 1 ) )
250 lwm = int( workq( 1 ) )
251 CALL dgemqr( 'L', trans, m, nrhs, n, a, lda, tq,
252 $ tszm, b, ldb, workq, -1, info2 )
253 lwm = max( lwm, int( workq( 1 ) ) )
254 wsizeo = tszo + lwo
255 wsizem = tszm + lwm
256 ELSE
257 CALL dgelq( m, n, a, lda, tq, -1, workq, -1, info2 )
258 tszo = int( tq( 1 ) )
259 lwo = int( workq( 1 ) )
260 CALL dgemlq( 'L', trans, n, nrhs, m, a, lda, tq,
261 $ tszo, b, ldb, workq, -1, info2 )
262 lwo = max( lwo, int( workq( 1 ) ) )
263 CALL dgelq( m, n, a, lda, tq, -2, workq, -2, info2 )
264 tszm = int( tq( 1 ) )
265 lwm = int( workq( 1 ) )
266 CALL dgemlq( 'L', trans, n, nrhs, m, a, lda, tq,
267 $ tszm, b, ldb, workq, -1, info2 )
268 lwm = max( lwm, int( workq( 1 ) ) )
269 wsizeo = tszo + lwo
270 wsizem = tszm + lwm
271 END IF
272*
273 IF( ( lwork.LT.wsizem ).AND.( .NOT.lquery ) ) THEN
274 info = -10
275 END IF
276*
277 work( 1 ) = dble( wsizeo )
278*
279 END IF
280*
281 IF( info.NE.0 ) THEN
282 CALL xerbla( 'DGETSLS', -info )
283 RETURN
284 END IF
285 IF( lquery ) THEN
286 IF( lwork.EQ.-2 ) work( 1 ) = dble( wsizem )
287 RETURN
288 END IF
289 IF( lwork.LT.wsizeo ) THEN
290 lw1 = tszm
291 lw2 = lwm
292 ELSE
293 lw1 = tszo
294 lw2 = lwo
295 END IF
296*
297* Quick return if possible
298*
299 IF( min( m, n, nrhs ).EQ.0 ) THEN
300 CALL dlaset( 'FULL', max( m, n ), nrhs, zero, zero,
301 $ b, ldb )
302 RETURN
303 END IF
304*
305* Get machine parameters
306*
307 smlnum = dlamch( 'S' ) / dlamch( 'P' )
308 bignum = one / smlnum
309*
310* Scale A, B if max element outside range [SMLNUM,BIGNUM]
311*
312 anrm = dlange( 'M', m, n, a, lda, work )
313 iascl = 0
314 IF( anrm.GT.zero .AND. anrm.LT.smlnum ) THEN
315*
316* Scale matrix norm up to SMLNUM
317*
318 CALL dlascl( 'G', 0, 0, anrm, smlnum, m, n, a, lda, info )
319 iascl = 1
320 ELSE IF( anrm.GT.bignum ) THEN
321*
322* Scale matrix norm down to BIGNUM
323*
324 CALL dlascl( 'G', 0, 0, anrm, bignum, m, n, a, lda, info )
325 iascl = 2
326 ELSE IF( anrm.EQ.zero ) THEN
327*
328* Matrix all zero. Return zero solution.
329*
330 CALL dlaset( 'F', maxmn, nrhs, zero, zero, b, ldb )
331 GO TO 50
332 END IF
333*
334 brow = m
335 IF ( tran ) THEN
336 brow = n
337 END IF
338 bnrm = dlange( 'M', brow, nrhs, b, ldb, work )
339 ibscl = 0
340 IF( bnrm.GT.zero .AND. bnrm.LT.smlnum ) THEN
341*
342* Scale matrix norm up to SMLNUM
343*
344 CALL dlascl( 'G', 0, 0, bnrm, smlnum, brow, nrhs, b, ldb,
345 $ info )
346 ibscl = 1
347 ELSE IF( bnrm.GT.bignum ) THEN
348*
349* Scale matrix norm down to BIGNUM
350*
351 CALL dlascl( 'G', 0, 0, bnrm, bignum, brow, nrhs, b, ldb,
352 $ info )
353 ibscl = 2
354 END IF
355*
356 IF ( m.GE.n ) THEN
357*
358* compute QR factorization of A
359*
360 CALL dgeqr( m, n, a, lda, work( lw2+1 ), lw1,
361 $ work( 1 ), lw2, info )
362 IF ( .NOT.tran ) THEN
363*
364* Least-Squares Problem min || A * X - B ||
365*
366* B(1:M,1:NRHS) := Q**T * B(1:M,1:NRHS)
367*
368 CALL dgemqr( 'L' , 'T', m, nrhs, n, a, lda,
369 $ work( lw2+1 ), lw1, b, ldb, work( 1 ), lw2,
370 $ info )
371*
372* B(1:N,1:NRHS) := inv(R) * B(1:N,1:NRHS)
373*
374 CALL dtrtrs( 'U', 'N', 'N', n, nrhs,
375 $ a, lda, b, ldb, info )
376 IF( info.GT.0 ) THEN
377 RETURN
378 END IF
379 scllen = n
380 ELSE
381*
382* Overdetermined system of equations A**T * X = B
383*
384* B(1:N,1:NRHS) := inv(R**T) * B(1:N,1:NRHS)
385*
386 CALL dtrtrs( 'U', 'T', 'N', n, nrhs,
387 $ a, lda, b, ldb, info )
388*
389 IF( info.GT.0 ) THEN
390 RETURN
391 END IF
392*
393* B(N+1:M,1:NRHS) = ZERO
394*
395 DO 20 j = 1, nrhs
396 DO 10 i = n + 1, m
397 b( i, j ) = zero
398 10 CONTINUE
399 20 CONTINUE
400*
401* B(1:M,1:NRHS) := Q(1:N,:) * B(1:N,1:NRHS)
402*
403 CALL dgemqr( 'L', 'N', m, nrhs, n, a, lda,
404 $ work( lw2+1 ), lw1, b, ldb, work( 1 ), lw2,
405 $ info )
406*
407 scllen = m
408*
409 END IF
410*
411 ELSE
412*
413* Compute LQ factorization of A
414*
415 CALL dgelq( m, n, a, lda, work( lw2+1 ), lw1,
416 $ work( 1 ), lw2, info )
417*
418* workspace at least M, optimally M*NB.
419*
420 IF( .NOT.tran ) THEN
421*
422* underdetermined system of equations A * X = B
423*
424* B(1:M,1:NRHS) := inv(L) * B(1:M,1:NRHS)
425*
426 CALL dtrtrs( 'L', 'N', 'N', m, nrhs,
427 $ a, lda, b, ldb, info )
428*
429 IF( info.GT.0 ) THEN
430 RETURN
431 END IF
432*
433* B(M+1:N,1:NRHS) = 0
434*
435 DO 40 j = 1, nrhs
436 DO 30 i = m + 1, n
437 b( i, j ) = zero
438 30 CONTINUE
439 40 CONTINUE
440*
441* B(1:N,1:NRHS) := Q(1:N,:)**T * B(1:M,1:NRHS)
442*
443 CALL dgemlq( 'L', 'T', n, nrhs, m, a, lda,
444 $ work( lw2+1 ), lw1, b, ldb, work( 1 ), lw2,
445 $ info )
446*
447* workspace at least NRHS, optimally NRHS*NB
448*
449 scllen = n
450*
451 ELSE
452*
453* overdetermined system min || A**T * X - B ||
454*
455* B(1:N,1:NRHS) := Q * B(1:N,1:NRHS)
456*
457 CALL dgemlq( 'L', 'N', n, nrhs, m, a, lda,
458 $ work( lw2+1 ), lw1, b, ldb, work( 1 ), lw2,
459 $ info )
460*
461* workspace at least NRHS, optimally NRHS*NB
462*
463* B(1:M,1:NRHS) := inv(L**T) * B(1:M,1:NRHS)
464*
465 CALL dtrtrs( 'Lower', 'Transpose', 'Non-unit', m, nrhs,
466 $ a, lda, b, ldb, info )
467*
468 IF( info.GT.0 ) THEN
469 RETURN
470 END IF
471*
472 scllen = m
473*
474 END IF
475*
476 END IF
477*
478* Undo scaling
479*
480 IF( iascl.EQ.1 ) THEN
481 CALL dlascl( 'G', 0, 0, anrm, smlnum, scllen, nrhs, b, ldb,
482 $ info )
483 ELSE IF( iascl.EQ.2 ) THEN
484 CALL dlascl( 'G', 0, 0, anrm, bignum, scllen, nrhs, b, ldb,
485 $ info )
486 END IF
487 IF( ibscl.EQ.1 ) THEN
488 CALL dlascl( 'G', 0, 0, smlnum, bnrm, scllen, nrhs, b, ldb,
489 $ info )
490 ELSE IF( ibscl.EQ.2 ) THEN
491 CALL dlascl( 'G', 0, 0, bignum, bnrm, scllen, nrhs, b, ldb,
492 $ info )
493 END IF
494*
495 50 CONTINUE
496 work( 1 ) = dble( tszo + lwo )
497 RETURN
498*
499* End of DGETSLS
500*
501 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:174
subroutine dgemqr(side, trans, m, n, k, a, lda, t, tsize, c, ldc, work, lwork, info)
DGEMQR
Definition dgemqr.f:175
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:171
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:142
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:108
subroutine dtrtrs(uplo, trans, diag, n, nrhs, a, lda, b, ldb, info)
DTRTRS
Definition dtrtrs.f:144