LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
zlaqp2.f
Go to the documentation of this file.
1 *> \brief \b ZLAQP2 computes a QR factorization with column pivoting of the matrix block.
2 *
3 * =========== DOCUMENTATION ===========
4 *
5 * Online html documentation available at
6 * http://www.netlib.org/lapack/explore-html/
7 *
8 *> \htmlonly
9 *> Download ZLAQP2 + dependencies
10 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/zlaqp2.f">
11 *> [TGZ]</a>
12 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/zlaqp2.f">
13 *> [ZIP]</a>
14 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/zlaqp2.f">
15 *> [TXT]</a>
16 *> \endhtmlonly
17 *
18 * Definition:
19 * ===========
20 *
21 * SUBROUTINE ZLAQP2( M, N, OFFSET, A, LDA, JPVT, TAU, VN1, VN2,
22 * WORK )
23 *
24 * .. Scalar Arguments ..
25 * INTEGER LDA, M, N, OFFSET
26 * ..
27 * .. Array Arguments ..
28 * INTEGER JPVT( * )
29 * DOUBLE PRECISION VN1( * ), VN2( * )
30 * COMPLEX*16 A( LDA, * ), TAU( * ), WORK( * )
31 * ..
32 *
33 *
34 *> \par Purpose:
35 * =============
36 *>
37 *> \verbatim
38 *>
39 *> ZLAQP2 computes a QR factorization with column pivoting of
40 *> the block A(OFFSET+1:M,1:N).
41 *> The block A(1:OFFSET,1:N) is accordingly pivoted, but not factorized.
42 *> \endverbatim
43 *
44 * Arguments:
45 * ==========
46 *
47 *> \param[in] M
48 *> \verbatim
49 *> M is INTEGER
50 *> The number of rows of the matrix A. M >= 0.
51 *> \endverbatim
52 *>
53 *> \param[in] N
54 *> \verbatim
55 *> N is INTEGER
56 *> The number of columns of the matrix A. N >= 0.
57 *> \endverbatim
58 *>
59 *> \param[in] OFFSET
60 *> \verbatim
61 *> OFFSET is INTEGER
62 *> The number of rows of the matrix A that must be pivoted
63 *> but no factorized. OFFSET >= 0.
64 *> \endverbatim
65 *>
66 *> \param[in,out] A
67 *> \verbatim
68 *> A is COMPLEX*16 array, dimension (LDA,N)
69 *> On entry, the M-by-N matrix A.
70 *> On exit, the upper triangle of block A(OFFSET+1:M,1:N) is
71 *> the triangular factor obtained; the elements in block
72 *> A(OFFSET+1:M,1:N) below the diagonal, together with the
73 *> array TAU, represent the orthogonal matrix Q as a product of
74 *> elementary reflectors. Block A(1:OFFSET,1:N) has been
75 *> accordingly pivoted, but no factorized.
76 *> \endverbatim
77 *>
78 *> \param[in] LDA
79 *> \verbatim
80 *> LDA is INTEGER
81 *> The leading dimension of the array A. LDA >= max(1,M).
82 *> \endverbatim
83 *>
84 *> \param[in,out] JPVT
85 *> \verbatim
86 *> JPVT is INTEGER array, dimension (N)
87 *> On entry, if JPVT(i) .ne. 0, the i-th column of A is permuted
88 *> to the front of A*P (a leading column); if JPVT(i) = 0,
89 *> the i-th column of A is a free column.
90 *> On exit, if JPVT(i) = k, then the i-th column of A*P
91 *> was the k-th column of A.
92 *> \endverbatim
93 *>
94 *> \param[out] TAU
95 *> \verbatim
96 *> TAU is COMPLEX*16 array, dimension (min(M,N))
97 *> The scalar factors of the elementary reflectors.
98 *> \endverbatim
99 *>
100 *> \param[in,out] VN1
101 *> \verbatim
102 *> VN1 is DOUBLE PRECISION array, dimension (N)
103 *> The vector with the partial column norms.
104 *> \endverbatim
105 *>
106 *> \param[in,out] VN2
107 *> \verbatim
108 *> VN2 is DOUBLE PRECISION array, dimension (N)
109 *> The vector with the exact column norms.
110 *> \endverbatim
111 *>
112 *> \param[out] WORK
113 *> \verbatim
114 *> WORK is COMPLEX*16 array, dimension (N)
115 *> \endverbatim
116 *
117 * Authors:
118 * ========
119 *
120 *> \author Univ. of Tennessee
121 *> \author Univ. of California Berkeley
122 *> \author Univ. of Colorado Denver
123 *> \author NAG Ltd.
124 *
125 *> \date September 2012
126 *
127 *> \ingroup complex16OTHERauxiliary
128 *
129 *> \par Contributors:
130 * ==================
131 *>
132 *> G. Quintana-Orti, Depto. de Informatica, Universidad Jaime I, Spain
133 *> X. Sun, Computer Science Dept., Duke University, USA
134 *> \n
135 *> Partial column norm updating strategy modified on April 2011
136 *> Z. Drmac and Z. Bujanovic, Dept. of Mathematics,
137 *> University of Zagreb, Croatia.
138 *
139 *> \par References:
140 * ================
141 *>
142 *> LAPACK Working Note 176
143 *
144 *> \htmlonly
145 *> <a href="http://www.netlib.org/lapack/lawnspdf/lawn176.pdf">[PDF]</a>
146 *> \endhtmlonly
147 *
148 * =====================================================================
149  SUBROUTINE zlaqp2( M, N, OFFSET, A, LDA, JPVT, TAU, VN1, VN2,
150  $ work )
151 *
152 * -- LAPACK auxiliary routine (version 3.4.2) --
153 * -- LAPACK is a software package provided by Univ. of Tennessee, --
154 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
155 * September 2012
156 *
157 * .. Scalar Arguments ..
158  INTEGER LDA, M, N, OFFSET
159 * ..
160 * .. Array Arguments ..
161  INTEGER JPVT( * )
162  DOUBLE PRECISION VN1( * ), VN2( * )
163  COMPLEX*16 A( lda, * ), TAU( * ), WORK( * )
164 * ..
165 *
166 * =====================================================================
167 *
168 * .. Parameters ..
169  DOUBLE PRECISION ZERO, ONE
170  COMPLEX*16 CONE
171  parameter ( zero = 0.0d+0, one = 1.0d+0,
172  $ cone = ( 1.0d+0, 0.0d+0 ) )
173 * ..
174 * .. Local Scalars ..
175  INTEGER I, ITEMP, J, MN, OFFPI, PVT
176  DOUBLE PRECISION TEMP, TEMP2, TOL3Z
177  COMPLEX*16 AII
178 * ..
179 * .. External Subroutines ..
180  EXTERNAL zlarf, zlarfg, zswap
181 * ..
182 * .. Intrinsic Functions ..
183  INTRINSIC abs, dconjg, max, min, sqrt
184 * ..
185 * .. External Functions ..
186  INTEGER IDAMAX
187  DOUBLE PRECISION DLAMCH, DZNRM2
188  EXTERNAL idamax, dlamch, dznrm2
189 * ..
190 * .. Executable Statements ..
191 *
192  mn = min( m-offset, n )
193  tol3z = sqrt(dlamch('Epsilon'))
194 *
195 * Compute factorization.
196 *
197  DO 20 i = 1, mn
198 *
199  offpi = offset + i
200 *
201 * Determine ith pivot column and swap if necessary.
202 *
203  pvt = ( i-1 ) + idamax( n-i+1, vn1( i ), 1 )
204 *
205  IF( pvt.NE.i ) THEN
206  CALL zswap( m, a( 1, pvt ), 1, a( 1, i ), 1 )
207  itemp = jpvt( pvt )
208  jpvt( pvt ) = jpvt( i )
209  jpvt( i ) = itemp
210  vn1( pvt ) = vn1( i )
211  vn2( pvt ) = vn2( i )
212  END IF
213 *
214 * Generate elementary reflector H(i).
215 *
216  IF( offpi.LT.m ) THEN
217  CALL zlarfg( m-offpi+1, a( offpi, i ), a( offpi+1, i ), 1,
218  $ tau( i ) )
219  ELSE
220  CALL zlarfg( 1, a( m, i ), a( m, i ), 1, tau( i ) )
221  END IF
222 *
223  IF( i.LT.n ) THEN
224 *
225 * Apply H(i)**H to A(offset+i:m,i+1:n) from the left.
226 *
227  aii = a( offpi, i )
228  a( offpi, i ) = cone
229  CALL zlarf( 'Left', m-offpi+1, n-i, a( offpi, i ), 1,
230  $ dconjg( tau( i ) ), a( offpi, i+1 ), lda,
231  $ work( 1 ) )
232  a( offpi, i ) = aii
233  END IF
234 *
235 * Update partial column norms.
236 *
237  DO 10 j = i + 1, n
238  IF( vn1( j ).NE.zero ) THEN
239 *
240 * NOTE: The following 4 lines follow from the analysis in
241 * Lapack Working Note 176.
242 *
243  temp = one - ( abs( a( offpi, j ) ) / vn1( j ) )**2
244  temp = max( temp, zero )
245  temp2 = temp*( vn1( j ) / vn2( j ) )**2
246  IF( temp2 .LE. tol3z ) THEN
247  IF( offpi.LT.m ) THEN
248  vn1( j ) = dznrm2( m-offpi, a( offpi+1, j ), 1 )
249  vn2( j ) = vn1( j )
250  ELSE
251  vn1( j ) = zero
252  vn2( j ) = zero
253  END IF
254  ELSE
255  vn1( j ) = vn1( j )*sqrt( temp )
256  END IF
257  END IF
258  10 CONTINUE
259 *
260  20 CONTINUE
261 *
262  RETURN
263 *
264 * End of ZLAQP2
265 *
266  END
subroutine zlaqp2(M, N, OFFSET, A, LDA, JPVT, TAU, VN1, VN2, WORK)
ZLAQP2 computes a QR factorization with column pivoting of the matrix block.
Definition: zlaqp2.f:151
subroutine zlarfg(N, ALPHA, X, INCX, TAU)
ZLARFG generates an elementary reflector (Householder matrix).
Definition: zlarfg.f:108
subroutine zswap(N, ZX, INCX, ZY, INCY)
ZSWAP
Definition: zswap.f:52
subroutine zlarf(SIDE, M, N, V, INCV, TAU, C, LDC, WORK)
ZLARF applies an elementary reflector to a general rectangular matrix.
Definition: zlarf.f:130