LAPACK 3.12.1
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
zgelqf.f
Go to the documentation of this file.
1*> \brief \b ZGELQF
2*
3* =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6* http://www.netlib.org/lapack/explore-html/
7*
8*> Download ZGELQF + dependencies
9*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/zgelqf.f">
10*> [TGZ]</a>
11*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/zgelqf.f">
12*> [ZIP]</a>
13*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/zgelqf.f">
14*> [TXT]</a>
15*
16* Definition:
17* ===========
18*
19* SUBROUTINE ZGELQF( M, N, A, LDA, TAU, WORK, LWORK, INFO )
20*
21* .. Scalar Arguments ..
22* INTEGER INFO, LDA, LWORK, M, N
23* ..
24* .. Array Arguments ..
25* COMPLEX*16 A( LDA, * ), TAU( * ), WORK( * )
26* ..
27*
28*
29*> \par Purpose:
30* =============
31*>
32*> \verbatim
33*>
34*> ZGELQF computes an LQ factorization of a complex M-by-N matrix A:
35*>
36*> A = ( L 0 ) * Q
37*>
38*> where:
39*>
40*> Q is a N-by-N orthogonal matrix;
41*> L is a lower-triangular M-by-M matrix;
42*> 0 is a M-by-(N-M) zero matrix, if M < N.
43*>
44*> \endverbatim
45*
46* Arguments:
47* ==========
48*
49*> \param[in] M
50*> \verbatim
51*> M is INTEGER
52*> The number of rows of the matrix A. M >= 0.
53*> \endverbatim
54*>
55*> \param[in] N
56*> \verbatim
57*> N is INTEGER
58*> The number of columns of the matrix A. N >= 0.
59*> \endverbatim
60*>
61*> \param[in,out] A
62*> \verbatim
63*> A is COMPLEX*16 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*> with the array TAU, represent the unitary matrix Q as a
69*> product of elementary reflectors (see Further Details).
70*> \endverbatim
71*>
72*> \param[in] LDA
73*> \verbatim
74*> LDA is INTEGER
75*> The leading dimension of the array A. LDA >= max(1,M).
76*> \endverbatim
77*>
78*> \param[out] TAU
79*> \verbatim
80*> TAU is COMPLEX*16 array, dimension (min(M,N))
81*> The scalar factors of the elementary reflectors (see Further
82*> Details).
83*> \endverbatim
84*>
85*> \param[out] WORK
86*> \verbatim
87*> WORK is COMPLEX*16 array, dimension (MAX(1,LWORK))
88*> On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
89*> \endverbatim
90*>
91*> \param[in] LWORK
92*> \verbatim
93*> LWORK is INTEGER
94*> The dimension of the array WORK.
95*> LWORK >= 1, if MIN(M,N) = 0, and LWORK >= M, otherwise.
96*> For optimum performance LWORK >= M*NB, where NB is the
97*> optimal blocksize.
98*>
99*> If LWORK = -1, then a workspace query is assumed; the routine
100*> only calculates the optimal size of the WORK array, returns
101*> this value as the first entry of the WORK array, and no error
102*> message related to LWORK is issued by XERBLA.
103*> \endverbatim
104*>
105*> \param[out] INFO
106*> \verbatim
107*> INFO is INTEGER
108*> = 0: successful exit
109*> < 0: if INFO = -i, the i-th argument had an illegal value
110*> \endverbatim
111*
112* Authors:
113* ========
114*
115*> \author Univ. of Tennessee
116*> \author Univ. of California Berkeley
117*> \author Univ. of Colorado Denver
118*> \author NAG Ltd.
119*
120*> \ingroup gelqf
121*
122*> \par Further Details:
123* =====================
124*>
125*> \verbatim
126*>
127*> The matrix Q is represented as a product of elementary reflectors
128*>
129*> Q = H(k)**H . . . H(2)**H H(1)**H, where k = min(m,n).
130*>
131*> Each H(i) has the form
132*>
133*> H(i) = I - tau * v * v**H
134*>
135*> where tau is a complex scalar, and v is a complex vector with
136*> v(1:i-1) = 0 and v(i) = 1; conjg(v(i+1:n)) is stored on exit in
137*> A(i,i+1:n), and tau in TAU(i).
138*> \endverbatim
139*>
140* =====================================================================
141 SUBROUTINE zgelqf( M, N, A, LDA, TAU, WORK, LWORK, INFO )
142*
143* -- LAPACK computational routine --
144* -- LAPACK is a software package provided by Univ. of Tennessee, --
145* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
146*
147* .. Scalar Arguments ..
148 INTEGER INFO, LDA, LWORK, M, N
149* ..
150* .. Array Arguments ..
151 COMPLEX*16 A( LDA, * ), TAU( * ), WORK( * )
152* ..
153*
154* =====================================================================
155*
156* .. Local Scalars ..
157 LOGICAL LQUERY
158 INTEGER I, IB, IINFO, IWS, K, LDWORK, LWKOPT, NB,
159 $ NBMIN, NX
160* ..
161* .. External Subroutines ..
162 EXTERNAL xerbla, zgelq2, zlarfb, zlarft
163* ..
164* .. Intrinsic Functions ..
165 INTRINSIC max, min
166* ..
167* .. External Functions ..
168 INTEGER ILAENV
169 EXTERNAL ilaenv
170* ..
171* .. Executable Statements ..
172*
173* Test the input arguments
174*
175 info = 0
176 k = min( m, n )
177 nb = ilaenv( 1, 'ZGELQF', ' ', m, n, -1, -1 )
178 lquery = ( lwork.EQ.-1 )
179 IF( m.LT.0 ) THEN
180 info = -1
181 ELSE IF( n.LT.0 ) THEN
182 info = -2
183 ELSE IF( lda.LT.max( 1, m ) ) THEN
184 info = -4
185 ELSE IF( .NOT.lquery ) THEN
186 IF( lwork.LE.0 .OR. ( n.GT.0 .AND. lwork.LT.max( 1, m ) ) )
187 $ info = -7
188 END IF
189 IF( info.NE.0 ) THEN
190 CALL xerbla( 'ZGELQF', -info )
191 RETURN
192 ELSE IF( lquery ) THEN
193 IF( k.EQ.0 ) THEN
194 lwkopt = 1
195 ELSE
196 lwkopt = m*nb
197 END IF
198 work( 1 ) = lwkopt
199 RETURN
200 END IF
201*
202* Quick return if possible
203*
204 IF( k.EQ.0 ) THEN
205 work( 1 ) = 1
206 RETURN
207 END IF
208*
209 nbmin = 2
210 nx = 0
211 iws = m
212 IF( nb.GT.1 .AND. nb.LT.k ) THEN
213*
214* Determine when to cross over from blocked to unblocked code.
215*
216 nx = max( 0, ilaenv( 3, 'ZGELQF', ' ', m, n, -1, -1 ) )
217 IF( nx.LT.k ) THEN
218*
219* Determine if workspace is large enough for blocked code.
220*
221 ldwork = m
222 iws = ldwork*nb
223 IF( lwork.LT.iws ) THEN
224*
225* Not enough workspace to use optimal NB: reduce NB and
226* determine the minimum value of NB.
227*
228 nb = lwork / ldwork
229 nbmin = max( 2, ilaenv( 2, 'ZGELQF', ' ', m, n, -1,
230 $ -1 ) )
231 END IF
232 END IF
233 END IF
234*
235 IF( nb.GE.nbmin .AND. nb.LT.k .AND. nx.LT.k ) THEN
236*
237* Use blocked code initially
238*
239 DO 10 i = 1, k - nx, nb
240 ib = min( k-i+1, nb )
241*
242* Compute the LQ factorization of the current block
243* A(i:i+ib-1,i:n)
244*
245 CALL zgelq2( ib, n-i+1, a( i, i ), lda, tau( i ), work,
246 $ iinfo )
247 IF( i+ib.LE.m ) THEN
248*
249* Form the triangular factor of the block reflector
250* H = H(i) H(i+1) . . . H(i+ib-1)
251*
252 CALL zlarft( 'Forward', 'Rowwise', n-i+1, ib, a( i,
253 $ i ),
254 $ lda, tau( i ), work, ldwork )
255*
256* Apply H to A(i+ib:m,i:n) from the right
257*
258 CALL zlarfb( 'Right', 'No transpose', 'Forward',
259 $ 'Rowwise', m-i-ib+1, n-i+1, ib, a( i, i ),
260 $ lda, work, ldwork, a( i+ib, i ), lda,
261 $ work( ib+1 ), ldwork )
262 END IF
263 10 CONTINUE
264 ELSE
265 i = 1
266 END IF
267*
268* Use unblocked code to factor the last or only block.
269*
270 IF( i.LE.k )
271 $ CALL zgelq2( m-i+1, n-i+1, a( i, i ), lda, tau( i ), work,
272 $ iinfo )
273*
274 work( 1 ) = iws
275 RETURN
276*
277* End of ZGELQF
278*
279 END
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine zgelq2(m, n, a, lda, tau, work, info)
ZGELQ2 computes the LQ factorization of a general rectangular matrix using an unblocked algorithm.
Definition zgelq2.f:127
subroutine zgelqf(m, n, a, lda, tau, work, lwork, info)
ZGELQF
Definition zgelqf.f:142
subroutine zlarfb(side, trans, direct, storev, m, n, k, v, ldv, t, ldt, c, ldc, work, ldwork)
ZLARFB applies a block reflector or its conjugate-transpose to a general rectangular matrix.
Definition zlarfb.f:195
recursive subroutine zlarft(direct, storev, n, k, v, ldv, tau, t, ldt)
ZLARFT forms the triangular factor T of a block reflector H = I - vtvH
Definition zlarft.f:162