LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
dgelqt.f
Go to the documentation of this file.
1*> \brief \b DGELQT
2*
3* =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6* http://www.netlib.org/lapack/explore-html/
7*
8*> \htmlonly
9*> Download DGELQT + dependencies
10*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dgelqt.f">
11*> [TGZ]</a>
12*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dgelqt.f">
13*> [ZIP]</a>
14*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dgelqt.f">
15*> [TXT]</a>
16*> \endhtmlonly
17*
18* Definition:
19* ===========
20*
21* SUBROUTINE DGELQT( M, N, MB, A, LDA, T, LDT, WORK, INFO )
22*
23* .. Scalar Arguments ..
24* INTEGER INFO, LDA, LDT, M, N, MB
25* ..
26* .. Array Arguments ..
27* DOUBLE PRECISION A( LDA, * ), T( LDT, * ), WORK( * )
28* ..
29*
30*
31*> \par Purpose:
32* =============
33*>
34*> \verbatim
35*>
36*> DGELQT computes a blocked LQ factorization of a real M-by-N matrix A
37*> using the compact WY representation of Q.
38*> \endverbatim
39*
40* Arguments:
41* ==========
42*
43*> \param[in] M
44*> \verbatim
45*> M is INTEGER
46*> The number of rows of the matrix A. M >= 0.
47*> \endverbatim
48*>
49*> \param[in] N
50*> \verbatim
51*> N is INTEGER
52*> The number of columns of the matrix A. N >= 0.
53*> \endverbatim
54*>
55*> \param[in] MB
56*> \verbatim
57*> MB is INTEGER
58*> The block size to be used in the blocked QR. MIN(M,N) >= MB >= 1.
59*> \endverbatim
60*>
61*> \param[in,out] A
62*> \verbatim
63*> A is DOUBLE PRECISION array, dimension (LDA,N)
64*> On entry, the M-by-N matrix A.
65*> On exit, the elements on and below the diagonal of the array
66*> contain the M-by-MIN(M,N) lower trapezoidal matrix L (L is
67*> lower triangular if M <= N); the elements above the diagonal
68*> are the rows of V.
69*> \endverbatim
70*>
71*> \param[in] LDA
72*> \verbatim
73*> LDA is INTEGER
74*> The leading dimension of the array A. LDA >= max(1,M).
75*> \endverbatim
76*>
77*> \param[out] T
78*> \verbatim
79*> T is DOUBLE PRECISION array, dimension (LDT,MIN(M,N))
80*> The upper triangular block reflectors stored in compact form
81*> as a sequence of upper triangular blocks. See below
82*> for further details.
83*> \endverbatim
84*>
85*> \param[in] LDT
86*> \verbatim
87*> LDT is INTEGER
88*> The leading dimension of the array T. LDT >= MB.
89*> \endverbatim
90*>
91*> \param[out] WORK
92*> \verbatim
93*> WORK is DOUBLE PRECISION array, dimension (MB*N)
94*> \endverbatim
95*>
96*> \param[out] INFO
97*> \verbatim
98*> INFO is INTEGER
99*> = 0: successful exit
100*> < 0: if INFO = -i, the i-th argument had an illegal value
101*> \endverbatim
102*
103* Authors:
104* ========
105*
106*> \author Univ. of Tennessee
107*> \author Univ. of California Berkeley
108*> \author Univ. of Colorado Denver
109*> \author NAG Ltd.
110*
111*> \ingroup gelqt
112*
113*> \par Further Details:
114* =====================
115*>
116*> \verbatim
117*>
118*> The matrix V stores the elementary reflectors H(i) in the i-th row
119*> above the diagonal. For example, if M=5 and N=3, the matrix V is
120*>
121*> V = ( 1 v1 v1 v1 v1 )
122*> ( 1 v2 v2 v2 )
123*> ( 1 v3 v3 )
124*>
125*>
126*> where the vi's represent the vectors which define H(i), which are returned
127*> in the matrix A. The 1's along the diagonal of V are not stored in A.
128*> Let K=MIN(M,N). The number of blocks is B = ceiling(K/MB), where each
129*> block is of order MB except for the last block, which is of order
130*> IB = K - (B-1)*MB. For each of the B blocks, a upper triangular block
131*> reflector factor is computed: T1, T2, ..., TB. The MB-by-MB (and IB-by-IB
132*> for the last block) T's are stored in the MB-by-K matrix T as
133*>
134*> T = (T1 T2 ... TB).
135*> \endverbatim
136*>
137* =====================================================================
138 SUBROUTINE dgelqt( M, N, MB, A, LDA, T, LDT, WORK, INFO )
139*
140* -- LAPACK computational routine --
141* -- LAPACK is a software package provided by Univ. of Tennessee, --
142* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
143*
144* .. Scalar Arguments ..
145 INTEGER INFO, LDA, LDT, M, N, MB
146* ..
147* .. Array Arguments ..
148 DOUBLE PRECISION A( LDA, * ), T( LDT, * ), WORK( * )
149* ..
150*
151* =====================================================================
152*
153* ..
154* .. Local Scalars ..
155 INTEGER I, IB, IINFO, K
156* ..
157* .. External Subroutines ..
158 EXTERNAL dgelqt3, dlarfb, xerbla
159* ..
160* .. Executable Statements ..
161*
162* Test the input arguments
163*
164 info = 0
165 IF( m.LT.0 ) THEN
166 info = -1
167 ELSE IF( n.LT.0 ) THEN
168 info = -2
169 ELSE IF( mb.LT.1 .OR. ( mb.GT.min(m,n) .AND. min(m,n).GT.0 ) )THEN
170 info = -3
171 ELSE IF( lda.LT.max( 1, m ) ) THEN
172 info = -5
173 ELSE IF( ldt.LT.mb ) THEN
174 info = -7
175 END IF
176 IF( info.NE.0 ) THEN
177 CALL xerbla( 'DGELQT', -info )
178 RETURN
179 END IF
180*
181* Quick return if possible
182*
183 k = min( m, n )
184 IF( k.EQ.0 ) RETURN
185*
186* Blocked loop of length K
187*
188 DO i = 1, k, mb
189 ib = min( k-i+1, mb )
190*
191* Compute the LQ factorization of the current block A(I:M,I:I+IB-1)
192*
193 CALL dgelqt3( ib, n-i+1, a(i,i), lda, t(1,i), ldt, iinfo )
194 IF( i+ib.LE.m ) THEN
195*
196* Update by applying H**T to A(I:M,I+IB:N) from the right
197*
198 CALL dlarfb( 'R', 'N', 'F', 'R', m-i-ib+1, n-i+1, ib,
199 $ a( i, i ), lda, t( 1, i ), ldt,
200 $ a( i+ib, i ), lda, work , m-i-ib+1 )
201 END IF
202 END DO
203 RETURN
204*
205* End of DGELQT
206*
207 END
subroutine xerbla(srname, info)
Definition cblat2.f:3285
recursive subroutine dgelqt3(m, n, a, lda, t, ldt, info)
DGELQT3 recursively computes a LQ factorization of a general real or complex matrix using the compact...
Definition dgelqt3.f:131
subroutine dgelqt(m, n, mb, a, lda, t, ldt, work, info)
DGELQT
Definition dgelqt.f:139
subroutine dlarfb(side, trans, direct, storev, m, n, k, v, ldv, t, ldt, c, ldc, work, ldwork)
DLARFB applies a block reflector or its transpose to a general rectangular matrix.
Definition dlarfb.f:197