LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
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 *> \date November 2011
112 *
113 *> \ingroup doubleGEcomputational
114 *
115 *> \par Further Details:
116 * =====================
117 *>
118 *> \verbatim
119 *>
120 *> The matrix Q is represented as a product of elementary reflectors
121 *>
122 *> Q = H(1) H(2) . . . H(n)
123 *>
124 *> Each H(i) has the form
125 *>
126 *> H = I - tau * v * v**T
127 *>
128 *> where tau is a real scalar, and v is a real vector with
129 *> v(1:i-1) = 0 and v(i) = 1; v(i+1:m) is stored on exit in A(i+1:m,i).
130 *>
131 *> The matrix P is represented in jpvt as follows: If
132 *> jpvt(j) = i
133 *> then the jth column of P is the ith canonical unit vector.
134 *>
135 *> Partial column norm updating strategy modified by
136 *> Z. Drmac and Z. Bujanovic, Dept. of Mathematics,
137 *> University of Zagreb, Croatia.
138 *> -- April 2011 --
139 *> For more details see LAPACK Working Note 176.
140 *> \endverbatim
141 *>
142 * =====================================================================
143  SUBROUTINE dgeqpf( M, N, A, LDA, JPVT, TAU, WORK, INFO )
144 *
145 * -- LAPACK computational routine (version 3.4.0) --
146 * -- LAPACK is a software package provided by Univ. of Tennessee, --
147 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
148 * November 2011
149 *
150 * .. Scalar Arguments ..
151  INTEGER INFO, LDA, M, N
152 * ..
153 * .. Array Arguments ..
154  INTEGER JPVT( * )
155  DOUBLE PRECISION A( lda, * ), TAU( * ), WORK( * )
156 * ..
157 *
158 * =====================================================================
159 *
160 * .. Parameters ..
161  DOUBLE PRECISION ZERO, ONE
162  parameter ( zero = 0.0d+0, one = 1.0d+0 )
163 * ..
164 * .. Local Scalars ..
165  INTEGER I, ITEMP, J, MA, MN, PVT
166  DOUBLE PRECISION AII, TEMP, TEMP2, TOL3Z
167 * ..
168 * .. External Subroutines ..
169  EXTERNAL dgeqr2, dlarf, dlarfg, dorm2r, dswap, xerbla
170 * ..
171 * .. Intrinsic Functions ..
172  INTRINSIC abs, max, min, sqrt
173 * ..
174 * .. External Functions ..
175  INTEGER IDAMAX
176  DOUBLE PRECISION DLAMCH, DNRM2
177  EXTERNAL idamax, dlamch, dnrm2
178 * ..
179 * .. Executable Statements ..
180 *
181 * Test the input arguments
182 *
183  info = 0
184  IF( m.LT.0 ) THEN
185  info = -1
186  ELSE IF( n.LT.0 ) THEN
187  info = -2
188  ELSE IF( lda.LT.max( 1, m ) ) THEN
189  info = -4
190  END IF
191  IF( info.NE.0 ) THEN
192  CALL xerbla( 'DGEQPF', -info )
193  RETURN
194  END IF
195 *
196  mn = min( m, n )
197  tol3z = sqrt(dlamch('Epsilon'))
198 *
199 * Move initial columns up front
200 *
201  itemp = 1
202  DO 10 i = 1, n
203  IF( jpvt( i ).NE.0 ) THEN
204  IF( i.NE.itemp ) THEN
205  CALL dswap( m, a( 1, i ), 1, a( 1, itemp ), 1 )
206  jpvt( i ) = jpvt( itemp )
207  jpvt( itemp ) = i
208  ELSE
209  jpvt( i ) = i
210  END IF
211  itemp = itemp + 1
212  ELSE
213  jpvt( i ) = i
214  END IF
215  10 CONTINUE
216  itemp = itemp - 1
217 *
218 * Compute the QR factorization and update remaining columns
219 *
220  IF( itemp.GT.0 ) THEN
221  ma = min( itemp, m )
222  CALL dgeqr2( m, ma, a, lda, tau, work, info )
223  IF( ma.LT.n ) THEN
224  CALL dorm2r( 'Left', 'Transpose', m, n-ma, ma, a, lda, tau,
225  $ a( 1, ma+1 ), lda, work, info )
226  END IF
227  END IF
228 *
229  IF( itemp.LT.mn ) THEN
230 *
231 * Initialize partial column norms. The first n elements of
232 * work store the exact column norms.
233 *
234  DO 20 i = itemp + 1, n
235  work( i ) = dnrm2( m-itemp, a( itemp+1, i ), 1 )
236  work( n+i ) = work( i )
237  20 CONTINUE
238 *
239 * Compute factorization
240 *
241  DO 40 i = itemp + 1, mn
242 *
243 * Determine ith pivot column and swap if necessary
244 *
245  pvt = ( i-1 ) + idamax( n-i+1, work( i ), 1 )
246 *
247  IF( pvt.NE.i ) THEN
248  CALL dswap( m, a( 1, pvt ), 1, a( 1, i ), 1 )
249  itemp = jpvt( pvt )
250  jpvt( pvt ) = jpvt( i )
251  jpvt( i ) = itemp
252  work( pvt ) = work( i )
253  work( n+pvt ) = work( n+i )
254  END IF
255 *
256 * Generate elementary reflector H(i)
257 *
258  IF( i.LT.m ) THEN
259  CALL dlarfg( m-i+1, a( i, i ), a( i+1, i ), 1, tau( i ) )
260  ELSE
261  CALL dlarfg( 1, a( m, m ), a( m, m ), 1, tau( m ) )
262  END IF
263 *
264  IF( i.LT.n ) THEN
265 *
266 * Apply H(i) to A(i:m,i+1:n) from the left
267 *
268  aii = a( i, i )
269  a( i, i ) = one
270  CALL dlarf( 'LEFT', m-i+1, n-i, a( i, i ), 1, tau( i ),
271  $ a( i, i+1 ), lda, work( 2*n+1 ) )
272  a( i, i ) = aii
273  END IF
274 *
275 * Update partial column norms
276 *
277  DO 30 j = i + 1, n
278  IF( work( j ).NE.zero ) THEN
279 *
280 * NOTE: The following 4 lines follow from the analysis in
281 * Lapack Working Note 176.
282 *
283  temp = abs( a( i, j ) ) / work( j )
284  temp = max( zero, ( one+temp )*( one-temp ) )
285  temp2 = temp*( work( j ) / work( n+j ) )**2
286  IF( temp2 .LE. tol3z ) THEN
287  IF( m-i.GT.0 ) THEN
288  work( j ) = dnrm2( m-i, a( i+1, j ), 1 )
289  work( n+j ) = work( j )
290  ELSE
291  work( j ) = zero
292  work( n+j ) = zero
293  END IF
294  ELSE
295  work( j ) = work( j )*sqrt( temp )
296  END IF
297  END IF
298  30 CONTINUE
299 *
300  40 CONTINUE
301  END IF
302  RETURN
303 *
304 * End of DGEQPF
305 *
306  END
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:123
subroutine dlarf(SIDE, M, N, V, INCV, TAU, C, LDC, WORK)
DLARF applies an elementary reflector to a general rectangular matrix.
Definition: dlarf.f:126
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine dswap(N, DX, INCX, DY, INCY)
DSWAP
Definition: dswap.f:53
subroutine dgeqpf(M, N, A, LDA, JPVT, TAU, WORK, INFO)
DGEQPF
Definition: dgeqpf.f:144
subroutine dlarfg(N, ALPHA, X, INCX, TAU)
DLARFG generates an elementary reflector (Householder matrix).
Definition: dlarfg.f:108
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:161