OpenRadioss 2025.1.11
OpenRadioss project
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 doubleGEsolve
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, dlabad, 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 CALL dlabad( smlnum, bignum )
298*
299* Scale A, B if max element outside range [SMLNUM,BIGNUM]
300*
301 anrm = dlange( 'M', m, n, a, lda, work )
302 iascl = 0
303 IF( anrm.GT.zero .AND. anrm.LT.smlnum ) THEN
304*
305* Scale matrix norm up to SMLNUM
306*
307 CALL dlascl( 'g', 0, 0, ANRM, SMLNUM, M, N, A, LDA, INFO )
308 IASCL = 1
309.GT. ELSE IF( ANRMBIGNUM ) THEN
310*
311* Scale matrix norm down to BIGNUM
312*
313 CALL DLASCL( 'g', 0, 0, ANRM, BIGNUM, M, N, A, LDA, INFO )
314 IASCL = 2
315.EQ. ELSE IF( ANRMZERO ) THEN
316*
317* Matrix all zero. Return zero solution.
318*
319 CALL DLASET( 'f', MAXMN, NRHS, ZERO, ZERO, B, LDB )
320 GO TO 50
321 END IF
322*
323 BROW = M
324 IF ( TRAN ) THEN
325 BROW = N
326 END IF
327 BNRM = DLANGE( 'm', BROW, NRHS, B, LDB, WORK )
328 IBSCL = 0
329.GT..AND..LT. IF( BNRMZERO BNRMSMLNUM ) THEN
330*
331* Scale matrix norm up to SMLNUM
332*
333 CALL DLASCL( 'g', 0, 0, BNRM, SMLNUM, BROW, NRHS, B, LDB,
334 $ INFO )
335 IBSCL = 1
336.GT. ELSE IF( BNRMBIGNUM ) THEN
337*
338* Scale matrix norm down to BIGNUM
339*
340 CALL DLASCL( 'g', 0, 0, BNRM, BIGNUM, BROW, NRHS, B, LDB,
341 $ INFO )
342 IBSCL = 2
343 END IF
344*
345.GE. IF ( MN ) THEN
346*
347* compute QR factorization of A
348*
349 CALL DGEQR( M, N, A, LDA, WORK( LW2+1 ), LW1,
350 $ WORK( 1 ), LW2, INFO )
351.NOT. IF ( TRAN ) THEN
352*
353* Least-Squares Problem min || A * X - B ||
354*
355* B(1:M,1:NRHS) := Q**T * B(1:M,1:NRHS)
356*
357 CALL DGEMQR( 'l' , 't', M, NRHS, N, A, LDA,
358 $ WORK( LW2+1 ), LW1, B, LDB, WORK( 1 ), LW2,
359 $ INFO )
360*
361* B(1:N,1:NRHS) := inv(R) * B(1:N,1:NRHS)
362*
363 CALL DTRTRS( 'u', 'n', 'n', N, NRHS,
364 $ A, LDA, B, LDB, INFO )
365.GT. IF( INFO0 ) THEN
366 RETURN
367 END IF
368 SCLLEN = N
369 ELSE
370*
371* Overdetermined system of equations A**T * X = B
372*
373* B(1:N,1:NRHS) := inv(R**T) * B(1:N,1:NRHS)
374*
375 CALL DTRTRS( 'u', 't', 'n', N, NRHS,
376 $ A, LDA, B, LDB, INFO )
377*
378.GT. IF( INFO0 ) THEN
379 RETURN
380 END IF
381*
382* B(N+1:M,1:NRHS) = ZERO
383*
384 DO 20 J = 1, NRHS
385 DO 10 I = N + 1, M
386 B( I, J ) = ZERO
387 10 CONTINUE
388 20 CONTINUE
389*
390* B(1:M,1:NRHS) := Q(1:N,:) * B(1:N,1:NRHS)
391*
392 CALL DGEMQR( 'l', 'n', M, NRHS, N, A, LDA,
393 $ WORK( LW2+1 ), LW1, B, LDB, WORK( 1 ), LW2,
394 $ INFO )
395*
396 SCLLEN = M
397*
398 END IF
399*
400 ELSE
401*
402* Compute LQ factorization of A
403*
404 CALL DGELQ( M, N, A, LDA, WORK( LW2+1 ), LW1,
405 $ WORK( 1 ), LW2, INFO )
406*
407* workspace at least M, optimally M*NB.
408*
409.NOT. IF( TRAN ) THEN
410*
411* underdetermined system of equations A * X = B
412*
413* B(1:M,1:NRHS) := inv(L) * B(1:M,1:NRHS)
414*
415 CALL DTRTRS( 'l', 'n', 'n', M, NRHS,
416 $ A, LDA, B, LDB, INFO )
417*
418.GT. IF( INFO0 ) THEN
419 RETURN
420 END IF
421*
422* B(M+1:N,1:NRHS) = 0
423*
424 DO 40 J = 1, NRHS
425 DO 30 I = M + 1, N
426 B( I, J ) = ZERO
427 30 CONTINUE
428 40 CONTINUE
429*
430* B(1:N,1:NRHS) := Q(1:N,:)**T * B(1:M,1:NRHS)
431*
432 CALL DGEMLQ( 'l', 't', N, NRHS, M, A, LDA,
433 $ WORK( LW2+1 ), LW1, B, LDB, WORK( 1 ), LW2,
434 $ INFO )
435*
436* workspace at least NRHS, optimally NRHS*NB
437*
438 SCLLEN = N
439*
440 ELSE
441*
442* overdetermined system min || A**T * X - B ||
443*
444* B(1:N,1:NRHS) := Q * B(1:N,1:NRHS)
445*
446 CALL DGEMLQ( 'l', 'n', N, NRHS, M, A, LDA,
447 $ WORK( LW2+1 ), LW1, B, LDB, WORK( 1 ), LW2,
448 $ INFO )
449*
450* workspace at least NRHS, optimally NRHS*NB
451*
452* B(1:M,1:NRHS) := inv(L**T) * B(1:M,1:NRHS)
453*
454 CALL DTRTRS( 'lower', 'transpose', 'non-unit', M, NRHS,
455 $ A, LDA, B, LDB, INFO )
456*
457.GT. IF( INFO0 ) THEN
458 RETURN
459 END IF
460*
461 SCLLEN = M
462*
463 END IF
464*
465 END IF
466*
467* Undo scaling
468*
469.EQ. IF( IASCL1 ) THEN
470 CALL DLASCL( 'g', 0, 0, ANRM, SMLNUM, SCLLEN, NRHS, B, LDB,
471 $ INFO )
472.EQ. ELSE IF( IASCL2 ) THEN
473 CALL DLASCL( 'g', 0, 0, ANRM, BIGNUM, SCLLEN, NRHS, B, LDB,
474 $ INFO )
475 END IF
476.EQ. IF( IBSCL1 ) THEN
477 CALL DLASCL( 'g', 0, 0, SMLNUM, BNRM, SCLLEN, NRHS, B, LDB,
478 $ INFO )
479.EQ. ELSE IF( IBSCL2 ) THEN
480 CALL DLASCL( 'g', 0, 0, BIGNUM, BNRM, SCLLEN, NRHS, B, LDB,
481 $ INFO )
482 END IF
483*
484 50 CONTINUE
485 WORK( 1 ) = DBLE( TSZO + LWO )
486 RETURN
487*
488* End of DGETSLS
489*
490 END
subroutine dgelq(m, n, a, lda, t, tsize, work, lwork, info)
DGELQ
Definition dgelq.f:172
subroutine dgemlq(side, trans, m, n, k, a, lda, t, tsize, c, ldc, work, lwork, info)
DGEMLQ
Definition dgemlq.f:171
subroutine dgemqr(side, trans, m, n, k, a, lda, t, tsize, c, ldc, work, lwork, info)
DGEMQR
Definition dgemqr.f:172
subroutine dgeqr(m, n, a, lda, t, tsize, work, lwork, info)
DGEQR
Definition dgeqr.f:174
subroutine dlabad(small, large)
DLABAD
Definition dlabad.f:74
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 xerbla(srname, info)
XERBLA
Definition xerbla.f:60
subroutine dgetsls(trans, m, n, nrhs, a, lda, b, ldb, work, lwork, info)
DGETSLS
Definition dgetsls.f:162
subroutine dtrtrs(uplo, trans, diag, n, nrhs, a, lda, b, ldb, info)
DTRTRS
Definition dtrtrs.f:140
#define min(a, b)
Definition macros.h:20
#define max(a, b)
Definition macros.h:21