LAPACK 3.12.1
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
dgels.f
Go to the documentation of this file.
1*> \brief <b> DGELS solves overdetermined or underdetermined systems for GE matrices</b>
2*
3* =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6* http://www.netlib.org/lapack/explore-html/
7*
8*> Download DGELS + dependencies
9*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dgels.f">
10*> [TGZ]</a>
11*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dgels.f">
12*> [ZIP]</a>
13*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dgels.f">
14*> [TXT]</a>
15*
16* Definition:
17* ===========
18*
19* SUBROUTINE DGELS( TRANS, M, N, NRHS, A, LDA, B, LDB, WORK, LWORK,
20* INFO )
21*
22* .. Scalar Arguments ..
23* CHARACTER TRANS
24* INTEGER INFO, LDA, LDB, LWORK, M, N, NRHS
25* ..
26* .. Array Arguments ..
27* DOUBLE PRECISION A( LDA, * ), B( LDB, * ), WORK( * )
28* ..
29*
30*
31*> \par Purpose:
32* =============
33*>
34*> \verbatim
35*>
36*> DGELS solves overdetermined or underdetermined real linear systems
37*> involving an M-by-N matrix A, or its transpose, using a QR or LQ
38*> factorization of A.
39*>
40*> It is assumed that A has full rank, and only a rudimentary protection
41*> against rank-deficient matrices is provided. This subroutine only detects
42*> exact rank-deficiency, where a diagonal element of the triangular factor
43*> of A is exactly zero.
44*>
45*> It is conceivable for one (or more) of the diagonal elements of the triangular
46*> factor of A to be subnormally tiny numbers without this subroutine signalling
47*> an error. The solutions computed for such almost-rank-deficient matrices may
48*> be less accurate due to a loss of numerical precision.
49*>
50*> The following options are provided:
51*>
52*> 1. If TRANS = 'N' and m >= n: find the least squares solution of
53*> an overdetermined system, i.e., solve the least squares problem
54*> minimize || B - A*X ||.
55*>
56*> 2. If TRANS = 'N' and m < n: find the minimum norm solution of
57*> an underdetermined system A * X = B.
58*>
59*> 3. If TRANS = 'T' and m >= n: find the minimum norm solution of
60*> an underdetermined system A**T * X = B.
61*>
62*> 4. If TRANS = 'T' and m < n: find the least squares solution of
63*> an overdetermined system, i.e., solve the least squares problem
64*> minimize || B - A**T * X ||.
65*>
66*> Several right hand side vectors b and solution vectors x can be
67*> handled in a single call; they are stored as the columns of the
68*> M-by-NRHS right hand side matrix B and the N-by-NRHS solution
69*> matrix X.
70*> \endverbatim
71*
72* Arguments:
73* ==========
74*
75*> \param[in] TRANS
76*> \verbatim
77*> TRANS is CHARACTER*1
78*> = 'N': the linear system involves A;
79*> = 'T': the linear system involves A**T.
80*> \endverbatim
81*>
82*> \param[in] M
83*> \verbatim
84*> M is INTEGER
85*> The number of rows of the matrix A. M >= 0.
86*> \endverbatim
87*>
88*> \param[in] N
89*> \verbatim
90*> N is INTEGER
91*> The number of columns of the matrix A. N >= 0.
92*> \endverbatim
93*>
94*> \param[in] NRHS
95*> \verbatim
96*> NRHS is INTEGER
97*> The number of right hand sides, i.e., the number of
98*> columns of the matrices B and X. NRHS >=0.
99*> \endverbatim
100*>
101*> \param[in,out] A
102*> \verbatim
103*> A is DOUBLE PRECISION array, dimension (LDA,N)
104*> On entry, the M-by-N matrix A.
105*> On exit,
106*> if M >= N, A is overwritten by details of its QR
107*> factorization as returned by DGEQRF;
108*> if M < N, A is overwritten by details of its LQ
109*> factorization as returned by DGELQF.
110*> \endverbatim
111*>
112*> \param[in] LDA
113*> \verbatim
114*> LDA is INTEGER
115*> The leading dimension of the array A. LDA >= max(1,M).
116*> \endverbatim
117*>
118*> \param[in,out] B
119*> \verbatim
120*> B is DOUBLE PRECISION array, dimension (LDB,NRHS)
121*> On entry, the matrix B of right hand side vectors, stored
122*> columnwise; B is M-by-NRHS if TRANS = 'N', or N-by-NRHS
123*> if TRANS = 'T'.
124*> On exit, if INFO = 0, B is overwritten by the solution
125*> vectors, stored columnwise:
126*> if TRANS = 'N' and m >= n, rows 1 to n of B contain the least
127*> squares solution vectors; the residual sum of squares for the
128*> solution in each column is given by the sum of squares of
129*> elements N+1 to M in that column;
130*> if TRANS = 'N' and m < n, rows 1 to N of B contain the
131*> minimum norm solution vectors;
132*> if TRANS = 'T' and m >= n, rows 1 to M of B contain the
133*> minimum norm solution vectors;
134*> if TRANS = 'T' and m < n, rows 1 to M of B contain the
135*> least squares solution vectors; the residual sum of squares
136*> for the solution in each column is given by the sum of
137*> squares of elements M+1 to N in that column.
138*> \endverbatim
139*>
140*> \param[in] LDB
141*> \verbatim
142*> LDB is INTEGER
143*> The leading dimension of the array B. LDB >= MAX(1,M,N).
144*> \endverbatim
145*>
146*> \param[out] WORK
147*> \verbatim
148*> WORK is DOUBLE PRECISION array, dimension (MAX(1,LWORK))
149*> On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
150*> \endverbatim
151*>
152*> \param[in] LWORK
153*> \verbatim
154*> LWORK is INTEGER
155*> The dimension of the array WORK.
156*> LWORK >= max( 1, MN + max( MN, NRHS ) ).
157*> For optimal performance,
158*> LWORK >= max( 1, MN + max( MN, NRHS )*NB ).
159*> where MN = min(M,N) and NB is the optimum block size.
160*>
161*> If LWORK = -1, then a workspace query is assumed; the routine
162*> only calculates the optimal size of the WORK array, returns
163*> this value as the first entry of the WORK array, and no error
164*> message related to LWORK is issued by XERBLA.
165*> \endverbatim
166*>
167*> \param[out] INFO
168*> \verbatim
169*> INFO is INTEGER
170*> = 0: successful exit
171*> < 0: if INFO = -i, the i-th argument had an illegal value
172*> > 0: if INFO = i, the i-th diagonal element of the
173*> triangular factor of A is exactly zero, so that A does not have
174*> full rank; the least squares solution could not be
175*> computed.
176*> \endverbatim
177*
178* Authors:
179* ========
180*
181*> \author Univ. of Tennessee
182*> \author Univ. of California Berkeley
183*> \author Univ. of Colorado Denver
184*> \author NAG Ltd.
185*
186*> \ingroup gels
187*
188* =====================================================================
189 SUBROUTINE dgels( TRANS, M, N, NRHS, A, LDA, B, LDB, WORK,
190 $ LWORK,
191 $ INFO )
192*
193* -- LAPACK driver routine --
194* -- LAPACK is a software package provided by Univ. of Tennessee, --
195* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
196*
197* .. Scalar Arguments ..
198 CHARACTER TRANS
199 INTEGER INFO, LDA, LDB, LWORK, M, N, NRHS
200* ..
201* .. Array Arguments ..
202 DOUBLE PRECISION A( LDA, * ), B( LDB, * ), WORK( * )
203* ..
204*
205* =====================================================================
206*
207* .. Parameters ..
208 DOUBLE PRECISION ZERO, ONE
209 PARAMETER ( ZERO = 0.0d0, one = 1.0d0 )
210* ..
211* .. Local Scalars ..
212 LOGICAL LQUERY, TPSD
213 INTEGER BROW, I, IASCL, IBSCL, J, MN, NB, SCLLEN, WSIZE
214 DOUBLE PRECISION ANRM, BIGNUM, BNRM, SMLNUM
215* ..
216* .. Local Arrays ..
217 DOUBLE PRECISION RWORK( 1 )
218* ..
219* .. External Functions ..
220 LOGICAL LSAME
221 INTEGER ILAENV
222 DOUBLE PRECISION DLAMCH, DLANGE
223 EXTERNAL lsame, ilaenv, dlamch, dlange
224* ..
225* .. External Subroutines ..
226 EXTERNAL dgelqf, dgeqrf, dlascl, dlaset, dormlq,
227 $ dormqr,
228 $ dtrtrs, xerbla
229* ..
230* .. Intrinsic Functions ..
231 INTRINSIC dble, max, min
232* ..
233* .. Executable Statements ..
234*
235* Test the input arguments.
236*
237 info = 0
238 mn = min( m, n )
239 lquery = ( lwork.EQ.-1 )
240 IF( .NOT.( lsame( trans, 'N' ) .OR.
241 $ lsame( trans, 'T' ) ) ) THEN
242 info = -1
243 ELSE IF( m.LT.0 ) THEN
244 info = -2
245 ELSE IF( n.LT.0 ) THEN
246 info = -3
247 ELSE IF( nrhs.LT.0 ) THEN
248 info = -4
249 ELSE IF( lda.LT.max( 1, m ) ) THEN
250 info = -6
251 ELSE IF( ldb.LT.max( 1, m, n ) ) THEN
252 info = -8
253 ELSE IF( lwork.LT.max( 1, mn+max( mn, nrhs ) ) .AND. .NOT.lquery )
254 $ THEN
255 info = -10
256 END IF
257*
258* Figure out optimal block size
259*
260 IF( info.EQ.0 .OR. info.EQ.-10 ) THEN
261*
262 tpsd = .true.
263 IF( lsame( trans, 'N' ) )
264 $ tpsd = .false.
265*
266 IF( m.GE.n ) THEN
267 nb = ilaenv( 1, 'DGEQRF', ' ', m, n, -1, -1 )
268 IF( tpsd ) THEN
269 nb = max( nb, ilaenv( 1, 'DORMQR', 'LN', m, nrhs, n,
270 $ -1 ) )
271 ELSE
272 nb = max( nb, ilaenv( 1, 'DORMQR', 'LT', m, nrhs, n,
273 $ -1 ) )
274 END IF
275 ELSE
276 nb = ilaenv( 1, 'DGELQF', ' ', m, n, -1, -1 )
277 IF( tpsd ) THEN
278 nb = max( nb, ilaenv( 1, 'DORMLQ', 'LT', n, nrhs, m,
279 $ -1 ) )
280 ELSE
281 nb = max( nb, ilaenv( 1, 'DORMLQ', 'LN', n, nrhs, m,
282 $ -1 ) )
283 END IF
284 END IF
285*
286 wsize = max( 1, mn+max( mn, nrhs )*nb )
287 work( 1 ) = dble( wsize )
288*
289 END IF
290*
291 IF( info.NE.0 ) THEN
292 CALL xerbla( 'DGELS ', -info )
293 RETURN
294 ELSE IF( lquery ) THEN
295 RETURN
296 END IF
297*
298* Quick return if possible
299*
300 IF( min( m, n, nrhs ).EQ.0 ) THEN
301 CALL dlaset( 'Full', max( m, n ), nrhs, zero, zero, 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, rwork )
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', max( m, n ), nrhs, zero, zero, b, ldb )
331 GO TO 50
332 END IF
333*
334 brow = m
335 IF( tpsd )
336 $ brow = n
337 bnrm = dlange( 'M', brow, nrhs, b, ldb, rwork )
338 ibscl = 0
339 IF( bnrm.GT.zero .AND. bnrm.LT.smlnum ) THEN
340*
341* Scale matrix norm up to SMLNUM
342*
343 CALL dlascl( 'G', 0, 0, bnrm, smlnum, brow, nrhs, b, ldb,
344 $ info )
345 ibscl = 1
346 ELSE IF( bnrm.GT.bignum ) THEN
347*
348* Scale matrix norm down to BIGNUM
349*
350 CALL dlascl( 'G', 0, 0, bnrm, bignum, brow, nrhs, b, ldb,
351 $ info )
352 ibscl = 2
353 END IF
354*
355 IF( m.GE.n ) THEN
356*
357* compute QR factorization of A
358*
359 CALL dgeqrf( m, n, a, lda, work( 1 ), work( mn+1 ),
360 $ lwork-mn,
361 $ info )
362*
363* workspace at least N, optimally N*NB
364*
365 IF( .NOT.tpsd ) THEN
366*
367* Least-Squares Problem min || A * X - B ||
368*
369* B(1:M,1:NRHS) := Q**T * B(1:M,1:NRHS)
370*
371 CALL dormqr( 'Left', 'Transpose', m, nrhs, n, a, lda,
372 $ work( 1 ), b, ldb, work( mn+1 ), lwork-mn,
373 $ info )
374*
375* workspace at least NRHS, optimally NRHS*NB
376*
377* B(1:N,1:NRHS) := inv(R) * B(1:N,1:NRHS)
378*
379 CALL dtrtrs( 'Upper', 'No transpose', 'Non-unit', n,
380 $ nrhs,
381 $ a, lda, b, ldb, info )
382*
383 IF( info.GT.0 ) THEN
384 RETURN
385 END IF
386*
387 scllen = n
388*
389 ELSE
390*
391* Underdetermined system of equations A**T * X = B
392*
393* B(1:N,1:NRHS) := inv(R**T) * B(1:N,1:NRHS)
394*
395 CALL dtrtrs( 'Upper', 'Transpose', 'Non-unit', n, nrhs,
396 $ a, lda, b, ldb, info )
397*
398 IF( info.GT.0 ) THEN
399 RETURN
400 END IF
401*
402* B(N+1:M,1:NRHS) = ZERO
403*
404 DO 20 j = 1, nrhs
405 DO 10 i = n + 1, m
406 b( i, j ) = zero
407 10 CONTINUE
408 20 CONTINUE
409*
410* B(1:M,1:NRHS) := Q(1:N,:) * B(1:N,1:NRHS)
411*
412 CALL dormqr( 'Left', 'No transpose', m, nrhs, n, a, lda,
413 $ work( 1 ), b, ldb, work( mn+1 ), lwork-mn,
414 $ info )
415*
416* workspace at least NRHS, optimally NRHS*NB
417*
418 scllen = m
419*
420 END IF
421*
422 ELSE
423*
424* Compute LQ factorization of A
425*
426 CALL dgelqf( m, n, a, lda, work( 1 ), work( mn+1 ),
427 $ lwork-mn,
428 $ info )
429*
430* workspace at least M, optimally M*NB.
431*
432 IF( .NOT.tpsd ) THEN
433*
434* underdetermined system of equations A * X = B
435*
436* B(1:M,1:NRHS) := inv(L) * B(1:M,1:NRHS)
437*
438 CALL dtrtrs( 'Lower', 'No transpose', 'Non-unit', m,
439 $ nrhs,
440 $ a, lda, b, ldb, info )
441*
442 IF( info.GT.0 ) THEN
443 RETURN
444 END IF
445*
446* B(M+1:N,1:NRHS) = 0
447*
448 DO 40 j = 1, nrhs
449 DO 30 i = m + 1, n
450 b( i, j ) = zero
451 30 CONTINUE
452 40 CONTINUE
453*
454* B(1:N,1:NRHS) := Q(1:N,:)**T * B(1:M,1:NRHS)
455*
456 CALL dormlq( 'Left', 'Transpose', n, nrhs, m, a, lda,
457 $ work( 1 ), b, ldb, work( mn+1 ), lwork-mn,
458 $ info )
459*
460* workspace at least NRHS, optimally NRHS*NB
461*
462 scllen = n
463*
464 ELSE
465*
466* overdetermined system min || A**T * X - B ||
467*
468* B(1:N,1:NRHS) := Q * B(1:N,1:NRHS)
469*
470 CALL dormlq( 'Left', 'No transpose', n, nrhs, m, a, lda,
471 $ work( 1 ), b, ldb, work( mn+1 ), lwork-mn,
472 $ info )
473*
474* workspace at least NRHS, optimally NRHS*NB
475*
476* B(1:M,1:NRHS) := inv(L**T) * B(1:M,1:NRHS)
477*
478 CALL dtrtrs( 'Lower', 'Transpose', 'Non-unit', m, nrhs,
479 $ a, lda, b, ldb, info )
480*
481 IF( info.GT.0 ) THEN
482 RETURN
483 END IF
484*
485 scllen = m
486*
487 END IF
488*
489 END IF
490*
491* Undo scaling
492*
493 IF( iascl.EQ.1 ) THEN
494 CALL dlascl( 'G', 0, 0, anrm, smlnum, scllen, nrhs, b, ldb,
495 $ info )
496 ELSE IF( iascl.EQ.2 ) THEN
497 CALL dlascl( 'G', 0, 0, anrm, bignum, scllen, nrhs, b, ldb,
498 $ info )
499 END IF
500 IF( ibscl.EQ.1 ) THEN
501 CALL dlascl( 'G', 0, 0, smlnum, bnrm, scllen, nrhs, b, ldb,
502 $ info )
503 ELSE IF( ibscl.EQ.2 ) THEN
504 CALL dlascl( 'G', 0, 0, bignum, bnrm, scllen, nrhs, b, ldb,
505 $ info )
506 END IF
507*
508 50 CONTINUE
509 work( 1 ) = dble( wsize )
510*
511 RETURN
512*
513* End of DGELS
514*
515 END
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine dgelqf(m, n, a, lda, tau, work, lwork, info)
DGELQF
Definition dgelqf.f:142
subroutine dgels(trans, m, n, nrhs, a, lda, b, ldb, work, lwork, info)
DGELS solves overdetermined or underdetermined systems for GE matrices
Definition dgels.f:192
subroutine dgeqrf(m, n, a, lda, tau, work, lwork, info)
DGEQRF
Definition dgeqrf.f:144
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
subroutine dormlq(side, trans, m, n, k, a, lda, tau, c, ldc, work, lwork, info)
DORMLQ
Definition dormlq.f:165
subroutine dormqr(side, trans, m, n, k, a, lda, tau, c, ldc, work, lwork, info)
DORMQR
Definition dormqr.f:165