LAPACK  3.10.0
LAPACK: Linear Algebra PACKage
zlaqps.f
Go to the documentation of this file.
1 *> \brief \b ZLAQPS computes a step of QR factorization with column pivoting of a real m-by-n matrix A by using BLAS level 3.
2 *
3 * =========== DOCUMENTATION ===========
4 *
5 * Online html documentation available at
6 * http://www.netlib.org/lapack/explore-html/
7 *
8 *> \htmlonly
9 *> Download ZLAQPS + dependencies
10 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/zlaqps.f">
11 *> [TGZ]</a>
12 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/zlaqps.f">
13 *> [ZIP]</a>
14 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/zlaqps.f">
15 *> [TXT]</a>
16 *> \endhtmlonly
17 *
18 * Definition:
19 * ===========
20 *
21 * SUBROUTINE ZLAQPS( M, N, OFFSET, NB, KB, A, LDA, JPVT, TAU, VN1,
22 * VN2, AUXV, F, LDF )
23 *
24 * .. Scalar Arguments ..
25 * INTEGER KB, LDA, LDF, M, N, NB, OFFSET
26 * ..
27 * .. Array Arguments ..
28 * INTEGER JPVT( * )
29 * DOUBLE PRECISION VN1( * ), VN2( * )
30 * COMPLEX*16 A( LDA, * ), AUXV( * ), F( LDF, * ), TAU( * )
31 * ..
32 *
33 *
34 *> \par Purpose:
35 * =============
36 *>
37 *> \verbatim
38 *>
39 *> ZLAQPS computes a step of QR factorization with column pivoting
40 *> of a complex M-by-N matrix A by using Blas-3. It tries to factorize
41 *> NB columns from A starting from the row OFFSET+1, and updates all
42 *> of the matrix with Blas-3 xGEMM.
43 *>
44 *> In some cases, due to catastrophic cancellations, it cannot
45 *> factorize NB columns. Hence, the actual number of factorized
46 *> columns is returned in KB.
47 *>
48 *> Block A(1:OFFSET,1:N) is accordingly pivoted, but not factorized.
49 *> \endverbatim
50 *
51 * Arguments:
52 * ==========
53 *
54 *> \param[in] M
55 *> \verbatim
56 *> M is INTEGER
57 *> The number of rows of the matrix A. M >= 0.
58 *> \endverbatim
59 *>
60 *> \param[in] N
61 *> \verbatim
62 *> N is INTEGER
63 *> The number of columns of the matrix A. N >= 0
64 *> \endverbatim
65 *>
66 *> \param[in] OFFSET
67 *> \verbatim
68 *> OFFSET is INTEGER
69 *> The number of rows of A that have been factorized in
70 *> previous steps.
71 *> \endverbatim
72 *>
73 *> \param[in] NB
74 *> \verbatim
75 *> NB is INTEGER
76 *> The number of columns to factorize.
77 *> \endverbatim
78 *>
79 *> \param[out] KB
80 *> \verbatim
81 *> KB is INTEGER
82 *> The number of columns actually factorized.
83 *> \endverbatim
84 *>
85 *> \param[in,out] A
86 *> \verbatim
87 *> A is COMPLEX*16 array, dimension (LDA,N)
88 *> On entry, the M-by-N matrix A.
89 *> On exit, block A(OFFSET+1:M,1:KB) is the triangular
90 *> factor obtained and block A(1:OFFSET,1:N) has been
91 *> accordingly pivoted, but no factorized.
92 *> The rest of the matrix, block A(OFFSET+1:M,KB+1:N) has
93 *> been updated.
94 *> \endverbatim
95 *>
96 *> \param[in] LDA
97 *> \verbatim
98 *> LDA is INTEGER
99 *> The leading dimension of the array A. LDA >= max(1,M).
100 *> \endverbatim
101 *>
102 *> \param[in,out] JPVT
103 *> \verbatim
104 *> JPVT is INTEGER array, dimension (N)
105 *> JPVT(I) = K <==> Column K of the full matrix A has been
106 *> permuted into position I in AP.
107 *> \endverbatim
108 *>
109 *> \param[out] TAU
110 *> \verbatim
111 *> TAU is COMPLEX*16 array, dimension (KB)
112 *> The scalar factors of the elementary reflectors.
113 *> \endverbatim
114 *>
115 *> \param[in,out] VN1
116 *> \verbatim
117 *> VN1 is DOUBLE PRECISION array, dimension (N)
118 *> The vector with the partial column norms.
119 *> \endverbatim
120 *>
121 *> \param[in,out] VN2
122 *> \verbatim
123 *> VN2 is DOUBLE PRECISION array, dimension (N)
124 *> The vector with the exact column norms.
125 *> \endverbatim
126 *>
127 *> \param[in,out] AUXV
128 *> \verbatim
129 *> AUXV is COMPLEX*16 array, dimension (NB)
130 *> Auxiliary vector.
131 *> \endverbatim
132 *>
133 *> \param[in,out] F
134 *> \verbatim
135 *> F is COMPLEX*16 array, dimension (LDF,NB)
136 *> Matrix F**H = L * Y**H * A.
137 *> \endverbatim
138 *>
139 *> \param[in] LDF
140 *> \verbatim
141 *> LDF is INTEGER
142 *> The leading dimension of the array F. LDF >= max(1,N).
143 *> \endverbatim
144 *
145 * Authors:
146 * ========
147 *
148 *> \author Univ. of Tennessee
149 *> \author Univ. of California Berkeley
150 *> \author Univ. of Colorado Denver
151 *> \author NAG Ltd.
152 *
153 *> \ingroup complex16OTHERauxiliary
154 *
155 *> \par Contributors:
156 * ==================
157 *>
158 *> G. Quintana-Orti, Depto. de Informatica, Universidad Jaime I, Spain
159 *> X. Sun, Computer Science Dept., Duke University, USA
160 *> \n
161 *> Partial column norm updating strategy modified on April 2011
162 *> Z. Drmac and Z. Bujanovic, Dept. of Mathematics,
163 *> University of Zagreb, Croatia.
164 *
165 *> \par References:
166 * ================
167 *>
168 *> LAPACK Working Note 176
169 *
170 *> \htmlonly
171 *> <a href="http://www.netlib.org/lapack/lawnspdf/lawn176.pdf">[PDF]</a>
172 *> \endhtmlonly
173 *
174 * =====================================================================
175  SUBROUTINE zlaqps( M, N, OFFSET, NB, KB, A, LDA, JPVT, TAU, VN1,
176  $ VN2, AUXV, F, LDF )
177 *
178 * -- LAPACK auxiliary routine --
179 * -- LAPACK is a software package provided by Univ. of Tennessee, --
180 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
181 *
182 * .. Scalar Arguments ..
183  INTEGER KB, LDA, LDF, M, N, NB, OFFSET
184 * ..
185 * .. Array Arguments ..
186  INTEGER JPVT( * )
187  DOUBLE PRECISION VN1( * ), VN2( * )
188  COMPLEX*16 A( LDA, * ), AUXV( * ), F( LDF, * ), TAU( * )
189 * ..
190 *
191 * =====================================================================
192 *
193 * .. Parameters ..
194  DOUBLE PRECISION ZERO, ONE
195  COMPLEX*16 CZERO, CONE
196  parameter( zero = 0.0d+0, one = 1.0d+0,
197  $ czero = ( 0.0d+0, 0.0d+0 ),
198  $ cone = ( 1.0d+0, 0.0d+0 ) )
199 * ..
200 * .. Local Scalars ..
201  INTEGER ITEMP, J, K, LASTRK, LSTICC, PVT, RK
202  DOUBLE PRECISION TEMP, TEMP2, TOL3Z
203  COMPLEX*16 AKK
204 * ..
205 * .. External Subroutines ..
206  EXTERNAL zgemm, zgemv, zlarfg, zswap
207 * ..
208 * .. Intrinsic Functions ..
209  INTRINSIC abs, dble, dconjg, max, min, nint, sqrt
210 * ..
211 * .. External Functions ..
212  INTEGER IDAMAX
213  DOUBLE PRECISION DLAMCH, DZNRM2
214  EXTERNAL idamax, dlamch, dznrm2
215 * ..
216 * .. Executable Statements ..
217 *
218  lastrk = min( m, n+offset )
219  lsticc = 0
220  k = 0
221  tol3z = sqrt(dlamch('Epsilon'))
222 *
223 * Beginning of while loop.
224 *
225  10 CONTINUE
226  IF( ( k.LT.nb ) .AND. ( lsticc.EQ.0 ) ) THEN
227  k = k + 1
228  rk = offset + k
229 *
230 * Determine ith pivot column and swap if necessary
231 *
232  pvt = ( k-1 ) + idamax( n-k+1, vn1( k ), 1 )
233  IF( pvt.NE.k ) THEN
234  CALL zswap( m, a( 1, pvt ), 1, a( 1, k ), 1 )
235  CALL zswap( k-1, f( pvt, 1 ), ldf, f( k, 1 ), ldf )
236  itemp = jpvt( pvt )
237  jpvt( pvt ) = jpvt( k )
238  jpvt( k ) = itemp
239  vn1( pvt ) = vn1( k )
240  vn2( pvt ) = vn2( k )
241  END IF
242 *
243 * Apply previous Householder reflectors to column K:
244 * A(RK:M,K) := A(RK:M,K) - A(RK:M,1:K-1)*F(K,1:K-1)**H.
245 *
246  IF( k.GT.1 ) THEN
247  DO 20 j = 1, k - 1
248  f( k, j ) = dconjg( f( k, j ) )
249  20 CONTINUE
250  CALL zgemv( 'No transpose', m-rk+1, k-1, -cone, a( rk, 1 ),
251  $ lda, f( k, 1 ), ldf, cone, a( rk, k ), 1 )
252  DO 30 j = 1, k - 1
253  f( k, j ) = dconjg( f( k, j ) )
254  30 CONTINUE
255  END IF
256 *
257 * Generate elementary reflector H(k).
258 *
259  IF( rk.LT.m ) THEN
260  CALL zlarfg( m-rk+1, a( rk, k ), a( rk+1, k ), 1, tau( k ) )
261  ELSE
262  CALL zlarfg( 1, a( rk, k ), a( rk, k ), 1, tau( k ) )
263  END IF
264 *
265  akk = a( rk, k )
266  a( rk, k ) = cone
267 *
268 * Compute Kth column of F:
269 *
270 * Compute F(K+1:N,K) := tau(K)*A(RK:M,K+1:N)**H*A(RK:M,K).
271 *
272  IF( k.LT.n ) THEN
273  CALL zgemv( 'Conjugate transpose', m-rk+1, n-k, tau( k ),
274  $ a( rk, k+1 ), lda, a( rk, k ), 1, czero,
275  $ f( k+1, k ), 1 )
276  END IF
277 *
278 * Padding F(1:K,K) with zeros.
279 *
280  DO 40 j = 1, k
281  f( j, k ) = czero
282  40 CONTINUE
283 *
284 * Incremental updating of F:
285 * F(1:N,K) := F(1:N,K) - tau(K)*F(1:N,1:K-1)*A(RK:M,1:K-1)**H
286 * *A(RK:M,K).
287 *
288  IF( k.GT.1 ) THEN
289  CALL zgemv( 'Conjugate transpose', m-rk+1, k-1, -tau( k ),
290  $ a( rk, 1 ), lda, a( rk, k ), 1, czero,
291  $ auxv( 1 ), 1 )
292 *
293  CALL zgemv( 'No transpose', n, k-1, cone, f( 1, 1 ), ldf,
294  $ auxv( 1 ), 1, cone, f( 1, k ), 1 )
295  END IF
296 *
297 * Update the current row of A:
298 * A(RK,K+1:N) := A(RK,K+1:N) - A(RK,1:K)*F(K+1:N,1:K)**H.
299 *
300  IF( k.LT.n ) THEN
301  CALL zgemm( 'No transpose', 'Conjugate transpose', 1, n-k,
302  $ k, -cone, a( rk, 1 ), lda, f( k+1, 1 ), ldf,
303  $ cone, a( rk, k+1 ), lda )
304  END IF
305 *
306 * Update partial column norms.
307 *
308  IF( rk.LT.lastrk ) THEN
309  DO 50 j = k + 1, n
310  IF( vn1( j ).NE.zero ) THEN
311 *
312 * NOTE: The following 4 lines follow from the analysis in
313 * Lapack Working Note 176.
314 *
315  temp = abs( a( rk, j ) ) / vn1( j )
316  temp = max( zero, ( one+temp )*( one-temp ) )
317  temp2 = temp*( vn1( j ) / vn2( j ) )**2
318  IF( temp2 .LE. tol3z ) THEN
319  vn2( j ) = dble( lsticc )
320  lsticc = j
321  ELSE
322  vn1( j ) = vn1( j )*sqrt( temp )
323  END IF
324  END IF
325  50 CONTINUE
326  END IF
327 *
328  a( rk, k ) = akk
329 *
330 * End of while loop.
331 *
332  GO TO 10
333  END IF
334  kb = k
335  rk = offset + kb
336 *
337 * Apply the block reflector to the rest of the matrix:
338 * A(OFFSET+KB+1:M,KB+1:N) := A(OFFSET+KB+1:M,KB+1:N) -
339 * A(OFFSET+KB+1:M,1:KB)*F(KB+1:N,1:KB)**H.
340 *
341  IF( kb.LT.min( n, m-offset ) ) THEN
342  CALL zgemm( 'No transpose', 'Conjugate transpose', m-rk, n-kb,
343  $ kb, -cone, a( rk+1, 1 ), lda, f( kb+1, 1 ), ldf,
344  $ cone, a( rk+1, kb+1 ), lda )
345  END IF
346 *
347 * Recomputation of difficult columns.
348 *
349  60 CONTINUE
350  IF( lsticc.GT.0 ) THEN
351  itemp = nint( vn2( lsticc ) )
352  vn1( lsticc ) = dznrm2( m-rk, a( rk+1, lsticc ), 1 )
353 *
354 * NOTE: The computation of VN1( LSTICC ) relies on the fact that
355 * SNRM2 does not fail on vectors with norm below the value of
356 * SQRT(DLAMCH('S'))
357 *
358  vn2( lsticc ) = vn1( lsticc )
359  lsticc = itemp
360  GO TO 60
361  END IF
362 *
363  RETURN
364 *
365 * End of ZLAQPS
366 *
367  END
subroutine zswap(N, ZX, INCX, ZY, INCY)
ZSWAP
Definition: zswap.f:81
subroutine zgemv(TRANS, M, N, ALPHA, A, LDA, X, INCX, BETA, Y, INCY)
ZGEMV
Definition: zgemv.f:158
subroutine zgemm(TRANSA, TRANSB, M, N, K, ALPHA, A, LDA, B, LDB, BETA, C, LDC)
ZGEMM
Definition: zgemm.f:187
subroutine zlaqps(M, N, OFFSET, NB, KB, A, LDA, JPVT, TAU, VN1, VN2, AUXV, F, LDF)
ZLAQPS computes a step of QR factorization with column pivoting of a real m-by-n matrix A by using BL...
Definition: zlaqps.f:177
subroutine zlarfg(N, ALPHA, X, INCX, TAU)
ZLARFG generates an elementary reflector (Householder matrix).
Definition: zlarfg.f:106