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