178 SUBROUTINE claqps( M, N, OFFSET, NB, KB, A, LDA, JPVT, TAU, VN1,
179 $ vn2, auxv, f, ldf )
187 INTEGER KB, LDA, LDF, M, N, NB, OFFSET
191 REAL VN1( * ), VN2( * )
192 COMPLEX A( lda, * ), AUXV( * ), F( ldf, * ), TAU( * )
200 parameter ( zero = 0.0e+0, one = 1.0e+0,
201 $ czero = ( 0.0e+0, 0.0e+0 ),
202 $ cone = ( 1.0e+0, 0.0e+0 ) )
205 INTEGER ITEMP, J, K, LASTRK, LSTICC, PVT, RK
206 REAL TEMP, TEMP2, TOL3Z
213 INTRINSIC abs, conjg, max, min, nint,
REAL, SQRT
218 EXTERNAL isamax, scnrm2, slamch
222 lastrk = min( m, n+offset )
225 tol3z = sqrt(slamch(
'Epsilon'))
230 IF( ( k.LT.nb ) .AND. ( lsticc.EQ.0 ) )
THEN
236 pvt = ( k-1 ) + isamax( n-k+1, vn1( k ), 1 )
238 CALL cswap( m, a( 1, pvt ), 1, a( 1, k ), 1 )
239 CALL cswap( k-1, f( pvt, 1 ), ldf, f( k, 1 ), ldf )
241 jpvt( pvt ) = jpvt( k )
243 vn1( pvt ) = vn1( k )
244 vn2( pvt ) = vn2( k )
252 f( k, j ) = conjg( f( k, j ) )
254 CALL cgemv(
'No transpose', m-rk+1, k-1, -cone, a( rk, 1 ),
255 $ lda, f( k, 1 ), ldf, cone, a( rk, k ), 1 )
257 f( k, j ) = conjg( f( k, j ) )
264 CALL clarfg( m-rk+1, a( rk, k ), a( rk+1, k ), 1, tau( k ) )
266 CALL clarfg( 1, a( rk, k ), a( rk, k ), 1, tau( k ) )
277 CALL cgemv(
'Conjugate transpose', m-rk+1, n-k, tau( k ),
278 $ a( rk, k+1 ), lda, a( rk, k ), 1, czero,
293 CALL cgemv(
'Conjugate transpose', m-rk+1, k-1, -tau( k ),
294 $ a( rk, 1 ), lda, a( rk, k ), 1, czero,
297 CALL cgemv(
'No transpose', n, k-1, cone, f( 1, 1 ), ldf,
298 $ auxv( 1 ), 1, cone, f( 1, k ), 1 )
305 CALL cgemm(
'No transpose',
'Conjugate transpose', 1, n-k,
306 $ k, -cone, a( rk, 1 ), lda, f( k+1, 1 ), ldf,
307 $ cone, a( rk, k+1 ), lda )
312 IF( rk.LT.lastrk )
THEN
314 IF( vn1( j ).NE.zero )
THEN
319 temp = abs( a( rk, j ) ) / vn1( j )
320 temp = max( zero, ( one+temp )*( one-temp ) )
321 temp2 = temp*( vn1( j ) / vn2( j ) )**2
322 IF( temp2 .LE. tol3z )
THEN
323 vn2( j ) =
REAL( lsticc )
326 vn1( j ) = vn1( j )*sqrt( temp )
345 IF( kb.LT.min( n, m-offset ) )
THEN
346 CALL cgemm(
'No transpose',
'Conjugate transpose', m-rk, n-kb,
347 $ kb, -cone, a( rk+1, 1 ), lda, f( kb+1, 1 ), ldf,
348 $ cone, a( rk+1, kb+1 ), lda )
354 IF( lsticc.GT.0 )
THEN
355 itemp = nint( vn2( lsticc ) )
356 vn1( lsticc ) = scnrm2( m-rk, a( rk+1, lsticc ), 1 )
362 vn2( lsticc ) = vn1( lsticc )
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...
subroutine cgemv(TRANS, M, N, ALPHA, A, LDA, X, INCX, BETA, Y, INCY)
CGEMV
subroutine cswap(N, CX, INCX, CY, INCY)
CSWAP
subroutine cgemm(TRANSA, TRANSB, M, N, K, ALPHA, A, LDA, B, LDB, BETA, C, LDC)
CGEMM
subroutine clarfg(N, ALPHA, X, INCX, TAU)
CLARFG generates an elementary reflector (Householder matrix).