LAPACK 3.12.1
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches

◆ dgetsls()

subroutine dgetsls ( character trans,
integer m,
integer n,
integer nrhs,
double precision, dimension( lda, * ) a,
integer lda,
double precision, dimension( ldb, * ) b,
integer ldb,
double precision, dimension( * ) work,
integer lwork,
integer info )

DGETSLS

Purpose:
!>
!> DGETSLS solves overdetermined or underdetermined real linear systems
!> involving an M-by-N matrix A, using a tall skinny QR or short wide LQ
!> factorization of A.
!>
!> It is assumed that A has full rank, and only a rudimentary protection
!> against rank-deficient matrices is provided. This subroutine only detects
!> exact rank-deficiency, where a diagonal element of the triangular factor
!> of A is exactly zero.
!>
!> It is conceivable for one (or more) of the diagonal elements of the triangular
!> factor of A to be subnormally tiny numbers without this subroutine signalling
!> an error. The solutions computed for such almost-rank-deficient matrices may
!> be less accurate due to a loss of numerical precision.
!>
!>
!> The following options are provided:
!>
!> 1. If TRANS = 'N' and m >= n:  find the least squares solution of
!>    an overdetermined system, i.e., solve the least squares problem
!>                 minimize || B - A*X ||.
!>
!> 2. If TRANS = 'N' and m < n:  find the minimum norm solution of
!>    an underdetermined system A * X = B.
!>
!> 3. If TRANS = 'T' and m >= n:  find the minimum norm solution of
!>    an undetermined system A**T * X = B.
!>
!> 4. If TRANS = 'T' and m < n:  find the least squares solution of
!>    an overdetermined system, i.e., solve the least squares problem
!>                 minimize || B - A**T * X ||.
!>
!> Several right hand side vectors b and solution vectors x can be
!> handled in a single call; they are stored as the columns of the
!> M-by-NRHS right hand side matrix B and the N-by-NRHS solution
!> matrix X.
!> 
Parameters
[in]TRANS
!>          TRANS is CHARACTER*1
!>          = 'N': the linear system involves A;
!>          = 'T': the linear system involves A**T.
!> 
[in]M
!>          M is INTEGER
!>          The number of rows of the matrix A.  M >= 0.
!> 
[in]N
!>          N is INTEGER
!>          The number of columns of the matrix A.  N >= 0.
!> 
[in]NRHS
!>          NRHS is INTEGER
!>          The number of right hand sides, i.e., the number of
!>          columns of the matrices B and X. NRHS >=0.
!> 
[in,out]A
!>          A is DOUBLE PRECISION array, dimension (LDA,N)
!>          On entry, the M-by-N matrix A.
!>          On exit,
!>          A is overwritten by details of its QR or LQ
!>          factorization as returned by DGEQR or DGELQ.
!> 
[in]LDA
!>          LDA is INTEGER
!>          The leading dimension of the array A.  LDA >= max(1,M).
!> 
[in,out]B
!>          B is DOUBLE PRECISION array, dimension (LDB,NRHS)
!>          On entry, the matrix B of right hand side vectors, stored
!>          columnwise; B is M-by-NRHS if TRANS = 'N', or N-by-NRHS
!>          if TRANS = 'T'.
!>          On exit, if INFO = 0, B is overwritten by the solution
!>          vectors, stored columnwise:
!>          if TRANS = 'N' and m >= n, rows 1 to n of B contain the least
!>          squares solution vectors.
!>          if TRANS = 'N' and m < n, rows 1 to N of B contain the
!>          minimum norm solution vectors;
!>          if TRANS = 'T' and m >= n, rows 1 to M of B contain the
!>          minimum norm solution vectors;
!>          if TRANS = 'T' and m < n, rows 1 to M of B contain the
!>          least squares solution vectors.
!> 
[in]LDB
!>          LDB is INTEGER
!>          The leading dimension of the array B. LDB >= MAX(1,M,N).
!> 
[out]WORK
!>          (workspace) DOUBLE PRECISION array, dimension (MAX(1,LWORK))
!>          On exit, if INFO = 0, WORK(1) contains optimal (or either minimal
!>          or optimal, if query was assumed) LWORK.
!>          See LWORK for details.
!> 
[in]LWORK
!>          LWORK is INTEGER
!>          The dimension of the array WORK. LWORK >= 1.
!>          If LWORK = -1 or -2, then a workspace query is assumed.
!>          If LWORK = -1, the routine calculates optimal size of WORK for the
!>          optimal performance and returns this value in WORK(1).
!>          If LWORK = -2, the routine calculates minimal size of WORK and 
!>          returns this value in WORK(1).
!> 
[out]INFO
!>          INFO is INTEGER
!>          = 0:  successful exit
!>          < 0:  if INFO = -i, the i-th argument had an illegal value
!>          > 0:  if INFO =  i, the i-th diagonal element of the
!>                triangular factor of A is exactly zero, so that A does not have
!>                full rank; the least squares solution could not be
!>                computed.
!> 
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.

Definition at line 169 of file dgetsls.f.

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*
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
double precision function dlamch(cmach)
DLAMCH
Definition dlamch.f:69
double precision function dlange(norm, m, n, a, lda, work)
DLANGE returns the value of the 1-norm, Frobenius norm, infinity-norm, or the largest absolute value ...
Definition dlange.f:112
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
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
subroutine dtrtrs(uplo, trans, diag, n, nrhs, a, lda, b, ldb, info)
DTRTRS
Definition dtrtrs.f:144
Here is the call graph for this function:
Here is the caller graph for this function: