LAPACK 3.12.1
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
zgeqpf.f
Go to the documentation of this file.
1*> \brief \b ZGEQPF
2*
3* =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6* http://www.netlib.org/lapack/explore-html/
7*
8*> Download ZGEQPF + dependencies
9*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/zgeqpf.f">
10*> [TGZ]</a>
11*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/zgeqpf.f">
12*> [ZIP]</a>
13*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/zgeqpf.f">
14*> [TXT]</a>
15*
16* Definition:
17* ===========
18*
19* SUBROUTINE ZGEQPF( M, N, A, LDA, JPVT, TAU, WORK, RWORK, INFO )
20*
21* .. Scalar Arguments ..
22* INTEGER INFO, LDA, M, N
23* ..
24* .. Array Arguments ..
25* INTEGER JPVT( * )
26* DOUBLE PRECISION RWORK( * )
27* COMPLEX*16 A( LDA, * ), TAU( * ), WORK( * )
28* ..
29*
30*
31*> \par Purpose:
32* =============
33*>
34*> \verbatim
35*>
36*> This routine is deprecated and has been replaced by routine ZGEQP3.
37*>
38*> ZGEQPF computes a QR factorization with column pivoting of a
39*> complex M-by-N matrix A: A*P = Q*R.
40*> \endverbatim
41*
42* Arguments:
43* ==========
44*
45*> \param[in] M
46*> \verbatim
47*> M is INTEGER
48*> The number of rows of the matrix A. M >= 0.
49*> \endverbatim
50*>
51*> \param[in] N
52*> \verbatim
53*> N is INTEGER
54*> The number of columns of the matrix A. N >= 0
55*> \endverbatim
56*>
57*> \param[in,out] A
58*> \verbatim
59*> A is COMPLEX*16 array, dimension (LDA,N)
60*> On entry, the M-by-N matrix A.
61*> On exit, the upper triangle of the array contains the
62*> min(M,N)-by-N upper triangular matrix R; the elements
63*> below the diagonal, together with the array TAU,
64*> represent the unitary matrix Q as a product of
65*> min(m,n) elementary reflectors.
66*> \endverbatim
67*>
68*> \param[in] LDA
69*> \verbatim
70*> LDA is INTEGER
71*> The leading dimension of the array A. LDA >= max(1,M).
72*> \endverbatim
73*>
74*> \param[in,out] JPVT
75*> \verbatim
76*> JPVT is INTEGER array, dimension (N)
77*> On entry, if JPVT(i) .ne. 0, the i-th column of A is permuted
78*> to the front of A*P (a leading column); if JPVT(i) = 0,
79*> the i-th column of A is a free column.
80*> On exit, if JPVT(i) = k, then the i-th column of A*P
81*> was the k-th column of A.
82*> \endverbatim
83*>
84*> \param[out] TAU
85*> \verbatim
86*> TAU is COMPLEX*16 array, dimension (min(M,N))
87*> The scalar factors of the elementary reflectors.
88*> \endverbatim
89*>
90*> \param[out] WORK
91*> \verbatim
92*> WORK is COMPLEX*16 array, dimension (N)
93*> \endverbatim
94*>
95*> \param[out] RWORK
96*> \verbatim
97*> RWORK is DOUBLE PRECISION array, dimension (2*N)
98*> \endverbatim
99*>
100*> \param[out] INFO
101*> \verbatim
102*> INFO is INTEGER
103*> = 0: successful exit
104*> < 0: if INFO = -i, the i-th argument had an illegal value
105*> \endverbatim
106*
107* Authors:
108* ========
109*
110*> \author Univ. of Tennessee
111*> \author Univ. of California Berkeley
112*> \author Univ. of Colorado Denver
113*> \author NAG Ltd.
114*
115*> \ingroup complex16GEcomputational
116*
117*> \par Further Details:
118* =====================
119*>
120*> \verbatim
121*>
122*> The matrix Q is represented as a product of elementary reflectors
123*>
124*> Q = H(1) H(2) . . . H(n)
125*>
126*> Each H(i) has the form
127*>
128*> H = I - tau * v * v**H
129*>
130*> where tau is a complex scalar, and v is a complex vector with
131*> v(1:i-1) = 0 and v(i) = 1; v(i+1:m) is stored on exit in A(i+1:m,i).
132*>
133*> The matrix P is represented in jpvt as follows: If
134*> jpvt(j) = i
135*> then the jth column of P is the ith canonical unit vector.
136*>
137*> Partial column norm updating strategy modified by
138*> Z. Drmac and Z. Bujanovic, Dept. of Mathematics,
139*> University of Zagreb, Croatia.
140*> -- April 2011 --
141*> For more details see LAPACK Working Note 176.
142*> \endverbatim
143*>
144* =====================================================================
145 SUBROUTINE zgeqpf( M, N, A, LDA, JPVT, TAU, WORK, RWORK, INFO )
146*
147* -- LAPACK computational routine --
148* -- LAPACK is a software package provided by Univ. of Tennessee, --
149* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
150*
151* .. Scalar Arguments ..
152 INTEGER INFO, LDA, M, N
153* ..
154* .. Array Arguments ..
155 INTEGER JPVT( * )
156 DOUBLE PRECISION RWORK( * )
157 COMPLEX*16 A( LDA, * ), TAU( * ), WORK( * )
158* ..
159*
160* =====================================================================
161*
162* .. Parameters ..
163 DOUBLE PRECISION ZERO, ONE
164 parameter( zero = 0.0d+0, one = 1.0d+0 )
165* ..
166* .. Local Scalars ..
167 INTEGER I, ITEMP, J, MA, MN, PVT
168 DOUBLE PRECISION TEMP, TEMP2, TOL3Z
169 COMPLEX*16 AII
170* ..
171* .. External Subroutines ..
172 EXTERNAL xerbla, zgeqr2, zlarf, zlarfg, zswap, zunm2r
173* ..
174* .. Intrinsic Functions ..
175 INTRINSIC abs, dcmplx, dconjg, max, min, sqrt
176* ..
177* .. External Functions ..
178 INTEGER IDAMAX
179 DOUBLE PRECISION DLAMCH, DZNRM2
180 EXTERNAL idamax, dlamch, dznrm2
181* ..
182* .. Executable Statements ..
183*
184* Test the input arguments
185*
186 info = 0
187 IF( m.LT.0 ) THEN
188 info = -1
189 ELSE IF( n.LT.0 ) THEN
190 info = -2
191 ELSE IF( lda.LT.max( 1, m ) ) THEN
192 info = -4
193 END IF
194 IF( info.NE.0 ) THEN
195 CALL xerbla( 'ZGEQPF', -info )
196 RETURN
197 END IF
198*
199 mn = min( m, n )
200 tol3z = sqrt(dlamch('Epsilon'))
201*
202* Move initial columns up front
203*
204 itemp = 1
205 DO 10 i = 1, n
206 IF( jpvt( i ).NE.0 ) THEN
207 IF( i.NE.itemp ) THEN
208 CALL zswap( m, a( 1, i ), 1, a( 1, itemp ), 1 )
209 jpvt( i ) = jpvt( itemp )
210 jpvt( itemp ) = i
211 ELSE
212 jpvt( i ) = i
213 END IF
214 itemp = itemp + 1
215 ELSE
216 jpvt( i ) = i
217 END IF
218 10 CONTINUE
219 itemp = itemp - 1
220*
221* Compute the QR factorization and update remaining columns
222*
223 IF( itemp.GT.0 ) THEN
224 ma = min( itemp, m )
225 CALL zgeqr2( m, ma, a, lda, tau, work, info )
226 IF( ma.LT.n ) THEN
227 CALL zunm2r( 'Left', 'Conjugate transpose', m, n-ma, ma, a,
228 $ lda, tau, a( 1, ma+1 ), lda, work, info )
229 END IF
230 END IF
231*
232 IF( itemp.LT.mn ) THEN
233*
234* Initialize partial column norms. The first n elements of
235* work store the exact column norms.
236*
237 DO 20 i = itemp + 1, n
238 rwork( i ) = dznrm2( m-itemp, a( itemp+1, i ), 1 )
239 rwork( n+i ) = rwork( i )
240 20 CONTINUE
241*
242* Compute factorization
243*
244 DO 40 i = itemp + 1, mn
245*
246* Determine ith pivot column and swap if necessary
247*
248 pvt = ( i-1 ) + idamax( n-i+1, rwork( i ), 1 )
249*
250 IF( pvt.NE.i ) THEN
251 CALL zswap( m, a( 1, pvt ), 1, a( 1, i ), 1 )
252 itemp = jpvt( pvt )
253 jpvt( pvt ) = jpvt( i )
254 jpvt( i ) = itemp
255 rwork( pvt ) = rwork( i )
256 rwork( n+pvt ) = rwork( n+i )
257 END IF
258*
259* Generate elementary reflector H(i)
260*
261 aii = a( i, i )
262 CALL zlarfg( m-i+1, aii, a( min( i+1, m ), i ), 1,
263 $ tau( i ) )
264 a( i, i ) = aii
265*
266 IF( i.LT.n ) THEN
267*
268* Apply H(i) to A(i:m,i+1:n) from the left
269*
270 aii = a( i, i )
271 a( i, i ) = dcmplx( one )
272 CALL zlarf( 'Left', m-i+1, n-i, a( i, i ), 1,
273 $ dconjg( tau( i ) ), a( i, i+1 ), lda, work )
274 a( i, i ) = aii
275 END IF
276*
277* Update partial column norms
278*
279 DO 30 j = i + 1, n
280 IF( rwork( j ).NE.zero ) THEN
281*
282* NOTE: The following 4 lines follow from the analysis in
283* Lapack Working Note 176.
284*
285 temp = abs( a( i, j ) ) / rwork( j )
286 temp = max( zero, ( one+temp )*( one-temp ) )
287 temp2 = temp*( rwork( j ) / rwork( n+j ) )**2
288 IF( temp2 .LE. tol3z ) THEN
289 IF( m-i.GT.0 ) THEN
290 rwork( j ) = dznrm2( m-i, a( i+1, j ), 1 )
291 rwork( n+j ) = rwork( j )
292 ELSE
293 rwork( j ) = zero
294 rwork( n+j ) = zero
295 END IF
296 ELSE
297 rwork( j ) = rwork( j )*sqrt( temp )
298 END IF
299 END IF
300 30 CONTINUE
301*
302 40 CONTINUE
303 END IF
304 RETURN
305*
306* End of ZGEQPF
307*
308 END
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine zgeqr2(m, n, a, lda, tau, work, info)
ZGEQR2 computes the QR factorization of a general rectangular matrix using an unblocked algorithm.
Definition zgeqr2.f:128
subroutine zlarf(side, m, n, v, incv, tau, c, ldc, work)
ZLARF applies an elementary reflector to a general rectangular matrix.
Definition zlarf.f:126
subroutine zlarfg(n, alpha, x, incx, tau)
ZLARFG generates an elementary reflector (Householder matrix).
Definition zlarfg.f:104
subroutine zswap(n, zx, incx, zy, incy)
ZSWAP
Definition zswap.f:81
subroutine zunm2r(side, trans, m, n, k, a, lda, tau, c, ldc, work, info)
ZUNM2R multiplies a general matrix by the unitary matrix from a QR factorization determined by cgeqrf...
Definition zunm2r.f:157
subroutine zgeqpf(m, n, a, lda, jpvt, tau, work, rwork, info)
ZGEQPF
Definition zgeqpf.f:146