OpenRadioss 2025.1.11
OpenRadioss project
Loading...
Searching...
No Matches
cgels.f
Go to the documentation of this file.
1*> \brief <b> CGELS 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*> \htmlonly
9*> Download CGELS + dependencies
10*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/cgels.f">
11*> [TGZ]</a>
12*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/cgels.f">
13*> [ZIP]</a>
14*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/cgels.f">
15*> [TXT]</a>
16*> \endhtmlonly
17*
18* Definition:
19* ===========
20*
21* SUBROUTINE CGELS( TRANS, M, N, NRHS, A, LDA, B, LDB, WORK, LWORK,
22* INFO )
23*
24* .. Scalar Arguments ..
25* CHARACTER TRANS
26* INTEGER INFO, LDA, LDB, LWORK, M, N, NRHS
27* ..
28* .. Array Arguments ..
29* COMPLEX A( LDA, * ), B( LDB, * ), WORK( * )
30* ..
31*
32*
33*> \par Purpose:
34* =============
35*>
36*> \verbatim
37*>
38*> CGELS solves overdetermined or underdetermined complex linear systems
39*> involving an M-by-N matrix A, or its conjugate-transpose, using a QR
40*> or LQ factorization of A. It is assumed that A has full rank.
41*>
42*> The following options are provided:
43*>
44*> 1. If TRANS = 'N' and m >= n: find the least squares solution of
45*> an overdetermined system, i.e., solve the least squares problem
46*> minimize || B - A*X ||.
47*>
48*> 2. If TRANS = 'N' and m < n: find the minimum norm solution of
49*> an underdetermined system A * X = B.
50*>
51*> 3. If TRANS = 'C' and m >= n: find the minimum norm solution of
52*> an underdetermined system A**H * X = B.
53*>
54*> 4. If TRANS = 'C' and m < n: find the least squares solution of
55*> an overdetermined system, i.e., solve the least squares problem
56*> minimize || B - A**H * X ||.
57*>
58*> Several right hand side vectors b and solution vectors x can be
59*> handled in a single call; they are stored as the columns of the
60*> M-by-NRHS right hand side matrix B and the N-by-NRHS solution
61*> matrix X.
62*> \endverbatim
63*
64* Arguments:
65* ==========
66*
67*> \param[in] TRANS
68*> \verbatim
69*> TRANS is CHARACTER*1
70*> = 'N': the linear system involves A;
71*> = 'C': the linear system involves A**H.
72*> \endverbatim
73*>
74*> \param[in] M
75*> \verbatim
76*> M is INTEGER
77*> The number of rows of the matrix A. M >= 0.
78*> \endverbatim
79*>
80*> \param[in] N
81*> \verbatim
82*> N is INTEGER
83*> The number of columns of the matrix A. N >= 0.
84*> \endverbatim
85*>
86*> \param[in] NRHS
87*> \verbatim
88*> NRHS is INTEGER
89*> The number of right hand sides, i.e., the number of
90*> columns of the matrices B and X. NRHS >= 0.
91*> \endverbatim
92*>
93*> \param[in,out] A
94*> \verbatim
95*> A is COMPLEX array, dimension (LDA,N)
96*> On entry, the M-by-N matrix A.
97*> if M >= N, A is overwritten by details of its QR
98*> factorization as returned by CGEQRF;
99*> if M < N, A is overwritten by details of its LQ
100*> factorization as returned by CGELQF.
101*> \endverbatim
102*>
103*> \param[in] LDA
104*> \verbatim
105*> LDA is INTEGER
106*> The leading dimension of the array A. LDA >= max(1,M).
107*> \endverbatim
108*>
109*> \param[in,out] B
110*> \verbatim
111*> B is COMPLEX array, dimension (LDB,NRHS)
112*> On entry, the matrix B of right hand side vectors, stored
113*> columnwise; B is M-by-NRHS if TRANS = 'N', or N-by-NRHS
114*> if TRANS = 'C'.
115*> On exit, if INFO = 0, B is overwritten by the solution
116*> vectors, stored columnwise:
117*> if TRANS = 'N' and m >= n, rows 1 to n of B contain the least
118*> squares solution vectors; the residual sum of squares for the
119*> solution in each column is given by the sum of squares of the
120*> modulus of elements N+1 to M in that column;
121*> if TRANS = 'N' and m < n, rows 1 to N of B contain the
122*> minimum norm solution vectors;
123*> if TRANS = 'C' and m >= n, rows 1 to M of B contain the
124*> minimum norm solution vectors;
125*> if TRANS = 'C' and m < n, rows 1 to M of B contain the
126*> least squares solution vectors; the residual sum of squares
127*> for the solution in each column is given by the sum of
128*> squares of the modulus of elements M+1 to N in that column.
129*> \endverbatim
130*>
131*> \param[in] LDB
132*> \verbatim
133*> LDB is INTEGER
134*> The leading dimension of the array B. LDB >= MAX(1,M,N).
135*> \endverbatim
136*>
137*> \param[out] WORK
138*> \verbatim
139*> WORK is COMPLEX array, dimension (MAX(1,LWORK))
140*> On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
141*> \endverbatim
142*>
143*> \param[in] LWORK
144*> \verbatim
145*> LWORK is INTEGER
146*> The dimension of the array WORK.
147*> LWORK >= max( 1, MN + max( MN, NRHS ) ).
148*> For optimal performance,
149*> LWORK >= max( 1, MN + max( MN, NRHS )*NB ).
150*> where MN = min(M,N) and NB is the optimum block size.
151*>
152*> If LWORK = -1, then a workspace query is assumed; the routine
153*> only calculates the optimal size of the WORK array, returns
154*> this value as the first entry of the WORK array, and no error
155*> message related to LWORK is issued by XERBLA.
156*> \endverbatim
157*>
158*> \param[out] INFO
159*> \verbatim
160*> INFO is INTEGER
161*> = 0: successful exit
162*> < 0: if INFO = -i, the i-th argument had an illegal value
163*> > 0: if INFO = i, the i-th diagonal element of the
164*> triangular factor of A is zero, so that A does not have
165*> full rank; the least squares solution could not be
166*> computed.
167*> \endverbatim
168*
169* Authors:
170* ========
171*
172*> \author Univ. of Tennessee
173*> \author Univ. of California Berkeley
174*> \author Univ. of Colorado Denver
175*> \author NAG Ltd.
176*
177*> \ingroup complexGEsolve
178*
179* =====================================================================
180 SUBROUTINE cgels( TRANS, M, N, NRHS, A, LDA, B, LDB, WORK, LWORK,
181 $ INFO )
182*
183* -- LAPACK driver routine --
184* -- LAPACK is a software package provided by Univ. of Tennessee, --
185* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
186*
187* .. Scalar Arguments ..
188 CHARACTER TRANS
189 INTEGER INFO, LDA, LDB, LWORK, M, N, NRHS
190* ..
191* .. Array Arguments ..
192 COMPLEX A( LDA, * ), B( LDB, * ), WORK( * )
193* ..
194*
195* =====================================================================
196*
197* .. Parameters ..
198 REAL ZERO, ONE
199 parameter( zero = 0.0e+0, one = 1.0e+0 )
200 COMPLEX CZERO
201 parameter( czero = ( 0.0e+0, 0.0e+0 ) )
202* ..
203* .. Local Scalars ..
204 LOGICAL LQUERY, TPSD
205 INTEGER BROW, I, IASCL, IBSCL, J, MN, NB, SCLLEN, WSIZE
206 REAL ANRM, BIGNUM, BNRM, SMLNUM
207* ..
208* .. Local Arrays ..
209 REAL RWORK( 1 )
210* ..
211* .. External Functions ..
212 LOGICAL LSAME
213 INTEGER ILAENV
214 REAL CLANGE, SLAMCH
215 EXTERNAL lsame, ilaenv, clange, slamch
216* ..
217* .. External Subroutines ..
218 EXTERNAL cgelqf, cgeqrf, clascl, claset, ctrtrs, cunmlq,
220* ..
221* .. Intrinsic Functions ..
222 INTRINSIC max, min, real
223* ..
224* .. Executable Statements ..
225*
226* Test the input arguments.
227*
228 info = 0
229 mn = min( m, n )
230 lquery = ( lwork.EQ.-1 )
231 IF( .NOT.( lsame( trans, 'N' ) .OR. lsame( trans, 'C' ) ) ) THEN
232 info = -1
233 ELSE IF( m.LT.0 ) THEN
234 info = -2
235 ELSE IF( n.LT.0 ) THEN
236 info = -3
237 ELSE IF( nrhs.LT.0 ) THEN
238 info = -4
239 ELSE IF( lda.LT.max( 1, m ) ) THEN
240 info = -6
241 ELSE IF( ldb.LT.max( 1, m, n ) ) THEN
242 info = -8
243 ELSE IF( lwork.LT.max( 1, mn+max( mn, nrhs ) ) .AND.
244 $ .NOT.lquery ) THEN
245 info = -10
246 END IF
247*
248* Figure out optimal block size
249*
250 IF( info.EQ.0 .OR. info.EQ.-10 ) THEN
251*
252 tpsd = .true.
253 IF( lsame( trans, 'N' ) )
254 $ tpsd = .false.
255*
256 IF( m.GE.n ) THEN
257 nb = ilaenv( 1, 'CGEQRF', ' ', m, n, -1, -1 )
258 IF( tpsd ) THEN
259 nb = max( nb, ilaenv( 1, 'CUNMQR', 'ln', M, NRHS, N,
260 $ -1 ) )
261 ELSE
262 NB = MAX( NB, ILAENV( 1, 'cunmqr', 'lc', M, NRHS, N,
263 $ -1 ) )
264 END IF
265 ELSE
266 NB = ILAENV( 1, 'cgelqf', ' ', M, N, -1, -1 )
267 IF( TPSD ) THEN
268 NB = MAX( NB, ILAENV( 1, 'cunmlq', 'lc', N, NRHS, M,
269 $ -1 ) )
270 ELSE
271 NB = MAX( NB, ILAENV( 1, 'cunmlq', 'ln', N, NRHS, M,
272 $ -1 ) )
273 END IF
274 END IF
275*
276 WSIZE = MAX( 1, MN + MAX( MN, NRHS )*NB )
277 WORK( 1 ) = REAL( WSIZE )
278*
279 END IF
280*
281.NE. IF( INFO0 ) THEN
282 CALL XERBLA( 'cgels ', -INFO )
283 RETURN
284 ELSE IF( LQUERY ) THEN
285 RETURN
286 END IF
287*
288* Quick return if possible
289*
290.EQ. IF( MIN( M, N, NRHS )0 ) THEN
291 CALL CLASET( 'full', MAX( M, N ), NRHS, CZERO, CZERO, B, LDB )
292 RETURN
293 END IF
294*
295* Get machine parameters
296*
297 SMLNUM = SLAMCH( 's' ) / SLAMCH( 'p' )
298 BIGNUM = ONE / SMLNUM
299 CALL SLABAD( SMLNUM, BIGNUM )
300*
301* Scale A, B if max element outside range [SMLNUM,BIGNUM]
302*
303 ANRM = CLANGE( 'm', M, N, A, LDA, RWORK )
304 IASCL = 0
305.GT..AND..LT. IF( ANRMZERO ANRMSMLNUM ) THEN
306*
307* Scale matrix norm up to SMLNUM
308*
309 CALL CLASCL( 'g', 0, 0, ANRM, SMLNUM, M, N, A, LDA, INFO )
310 IASCL = 1
311.GT. ELSE IF( ANRMBIGNUM ) THEN
312*
313* Scale matrix norm down to BIGNUM
314*
315 CALL CLASCL( 'g', 0, 0, ANRM, BIGNUM, M, N, A, LDA, INFO )
316 IASCL = 2
317.EQ. ELSE IF( ANRMZERO ) THEN
318*
319* Matrix all zero. Return zero solution.
320*
321 CALL CLASET( 'f', MAX( M, N ), NRHS, CZERO, CZERO, B, LDB )
322 GO TO 50
323 END IF
324*
325 BROW = M
326 IF( TPSD )
327 $ BROW = N
328 BNRM = CLANGE( 'm', BROW, NRHS, B, LDB, RWORK )
329 IBSCL = 0
330.GT..AND..LT. IF( BNRMZERO BNRMSMLNUM ) THEN
331*
332* Scale matrix norm up to SMLNUM
333*
334 CALL CLASCL( 'g', 0, 0, BNRM, SMLNUM, BROW, NRHS, B, LDB,
335 $ INFO )
336 IBSCL = 1
337.GT. ELSE IF( BNRMBIGNUM ) THEN
338*
339* Scale matrix norm down to BIGNUM
340*
341 CALL CLASCL( 'g', 0, 0, BNRM, BIGNUM, BROW, NRHS, B, LDB,
342 $ INFO )
343 IBSCL = 2
344 END IF
345*
346.GE. IF( MN ) THEN
347*
348* compute QR factorization of A
349*
350 CALL CGEQRF( M, N, A, LDA, WORK( 1 ), WORK( MN+1 ), LWORK-MN,
351 $ INFO )
352*
353* workspace at least N, optimally N*NB
354*
355.NOT. IF( TPSD ) THEN
356*
357* Least-Squares Problem min || A * X - B ||
358*
359* B(1:M,1:NRHS) := Q**H * B(1:M,1:NRHS)
360*
361 CALL CUNMQR( 'left', 'conjugate transpose', M, NRHS, N, A,
362 $ LDA, WORK( 1 ), B, LDB, WORK( MN+1 ), LWORK-MN,
363 $ INFO )
364*
365* workspace at least NRHS, optimally NRHS*NB
366*
367* B(1:N,1:NRHS) := inv(R) * B(1:N,1:NRHS)
368*
369 CALL CTRTRS( 'upper', 'no transpose', 'non-unit', N, NRHS,
370 $ A, LDA, B, LDB, INFO )
371*
372.GT. IF( INFO0 ) THEN
373 RETURN
374 END IF
375*
376 SCLLEN = N
377*
378 ELSE
379*
380* Underdetermined system of equations A**T * X = B
381*
382* B(1:N,1:NRHS) := inv(R**H) * B(1:N,1:NRHS)
383*
384 CALL CTRTRS( 'upper', 'conjugate transpose','non-unit',
385 $ N, NRHS, A, LDA, B, LDB, INFO )
386*
387.GT. IF( INFO0 ) THEN
388 RETURN
389 END IF
390*
391* B(N+1:M,1:NRHS) = ZERO
392*
393 DO 20 J = 1, NRHS
394 DO 10 I = N + 1, M
395 B( I, J ) = CZERO
396 10 CONTINUE
397 20 CONTINUE
398*
399* B(1:M,1:NRHS) := Q(1:N,:) * B(1:N,1:NRHS)
400*
401 CALL CUNMQR( 'left', 'no transpose', M, NRHS, N, A, LDA,
402 $ WORK( 1 ), B, LDB, WORK( MN+1 ), LWORK-MN,
403 $ INFO )
404*
405* workspace at least NRHS, optimally NRHS*NB
406*
407 SCLLEN = M
408*
409 END IF
410*
411 ELSE
412*
413* Compute LQ factorization of A
414*
415 CALL CGELQF( M, N, A, LDA, WORK( 1 ), WORK( MN+1 ), LWORK-MN,
416 $ INFO )
417*
418* workspace at least M, optimally M*NB.
419*
420.NOT. IF( TPSD ) 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 CTRTRS( 'lower', 'no transpose', 'non-unit', M, NRHS,
427 $ A, LDA, B, LDB, INFO )
428*
429.GT. IF( INFO0 ) 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 ) = CZERO
438 30 CONTINUE
439 40 CONTINUE
440*
441* B(1:N,1:NRHS) := Q(1:N,:)**H * B(1:M,1:NRHS)
442*
443 CALL CUNMLQ( 'left', 'conjugate transpose', N, NRHS, M, A,
444 $ LDA, WORK( 1 ), B, LDB, WORK( MN+1 ), LWORK-MN,
445 $ INFO )
446*
447* workspace at least NRHS, optimally NRHS*NB
448*
449 SCLLEN = N
450*
451 ELSE
452*
453* overdetermined system min || A**H * X - B ||
454*
455* B(1:N,1:NRHS) := Q * B(1:N,1:NRHS)
456*
457 CALL CUNMLQ( 'left', 'no transpose', N, NRHS, M, A, LDA,
458 $ WORK( 1 ), B, LDB, WORK( MN+1 ), LWORK-MN,
459 $ INFO )
460*
461* workspace at least NRHS, optimally NRHS*NB
462*
463* B(1:M,1:NRHS) := inv(L**H) * B(1:M,1:NRHS)
464*
465 CALL CTRTRS( 'lower', 'conjugate transpose', 'non-unit',
466 $ M, NRHS, A, LDA, B, LDB, INFO )
467*
468.GT. IF( INFO0 ) THEN
469 RETURN
470 END IF
471*
472 SCLLEN = M
473*
474 END IF
475*
476 END IF
477*
478* Undo scaling
479*
480.EQ. IF( IASCL1 ) THEN
481 CALL CLASCL( 'g', 0, 0, ANRM, SMLNUM, SCLLEN, NRHS, B, LDB,
482 $ INFO )
483.EQ. ELSE IF( IASCL2 ) THEN
484 CALL CLASCL( 'g', 0, 0, ANRM, BIGNUM, SCLLEN, NRHS, B, LDB,
485 $ INFO )
486 END IF
487.EQ. IF( IBSCL1 ) THEN
488 CALL CLASCL( 'g', 0, 0, SMLNUM, BNRM, SCLLEN, NRHS, B, LDB,
489 $ INFO )
490.EQ. ELSE IF( IBSCL2 ) THEN
491 CALL CLASCL( 'g', 0, 0, BIGNUM, BNRM, SCLLEN, NRHS, B, LDB,
492 $ INFO )
493 END IF
494*
495 50 CONTINUE
496 WORK( 1 ) = REAL( WSIZE )
497*
498 RETURN
499*
500* End of CGELS
501*
502 END
subroutine slabad(small, large)
SLABAD
Definition slabad.f:74
subroutine xerbla(srname, info)
XERBLA
Definition xerbla.f:60
logical function lsame(ca, cb)
LSAME
Definition lsame.f:53
subroutine cgeqrf(m, n, a, lda, tau, work, lwork, info)
CGEQRF
Definition cgeqrf.f:146
subroutine cgelqf(m, n, a, lda, tau, work, lwork, info)
CGELQF
Definition cgelqf.f:143
subroutine cgels(trans, m, n, nrhs, a, lda, b, ldb, work, lwork, info)
CGELS solves overdetermined or underdetermined systems for GE matrices
Definition cgels.f:182
subroutine clascl(type, kl, ku, cfrom, cto, m, n, a, lda, info)
CLASCL multiplies a general rectangular matrix by a real scalar defined as cto/cfrom.
Definition clascl.f:143
subroutine claset(uplo, m, n, alpha, beta, a, lda)
CLASET initializes the off-diagonal elements and the diagonal elements of a matrix to given values.
Definition claset.f:106
subroutine cunmqr(side, trans, m, n, k, a, lda, tau, c, ldc, work, lwork, info)
CUNMQR
Definition cunmqr.f:168
subroutine cunmlq(side, trans, m, n, k, a, lda, tau, c, ldc, work, lwork, info)
CUNMLQ
Definition cunmlq.f:168
subroutine ctrtrs(uplo, trans, diag, n, nrhs, a, lda, b, ldb, info)
CTRTRS
Definition ctrtrs.f:140
#define min(a, b)
Definition macros.h:20
#define max(a, b)
Definition macros.h:21