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