178 SUBROUTINE slaqps( 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 A( lda, * ), AUXV( * ), F( ldf, * ), TAU( * ),
199 parameter ( zero = 0.0e+0, one = 1.0e+0 )
202 INTEGER ITEMP, J, K, LASTRK, LSTICC, PVT, RK
203 REAL AKK, TEMP, TEMP2, TOL3Z
209 INTRINSIC abs, max, min, nint,
REAL, SQRT
214 EXTERNAL isamax, slamch, snrm2
218 lastrk = min( m, n+offset )
221 tol3z = sqrt(slamch(
'Epsilon'))
226 IF( ( k.LT.nb ) .AND. ( lsticc.EQ.0 ) )
THEN
232 pvt = ( k-1 ) + isamax( n-k+1, vn1( k ), 1 )
234 CALL sswap( m, a( 1, pvt ), 1, a( 1, k ), 1 )
235 CALL sswap( k-1, f( pvt, 1 ), ldf, f( k, 1 ), ldf )
237 jpvt( pvt ) = jpvt( k )
239 vn1( pvt ) = vn1( k )
240 vn2( pvt ) = vn2( k )
247 CALL sgemv(
'No transpose', m-rk+1, k-1, -one, a( rk, 1 ),
248 $ lda, f( k, 1 ), ldf, one, a( rk, k ), 1 )
254 CALL slarfg( m-rk+1, a( rk, k ), a( rk+1, k ), 1, tau( k ) )
256 CALL slarfg( 1, a( rk, k ), a( rk, k ), 1, tau( k ) )
267 CALL sgemv(
'Transpose', m-rk+1, n-k, tau( k ),
268 $ a( rk, k+1 ), lda, a( rk, k ), 1, zero,
283 CALL sgemv(
'Transpose', m-rk+1, k-1, -tau( k ), a( rk, 1 ),
284 $ lda, a( rk, k ), 1, zero, auxv( 1 ), 1 )
286 CALL sgemv(
'No transpose', n, k-1, one, f( 1, 1 ), ldf,
287 $ auxv( 1 ), 1, one, f( 1, k ), 1 )
294 CALL sgemv(
'No transpose', n-k, k, -one, f( k+1, 1 ), ldf,
295 $ a( rk, 1 ), lda, one, a( rk, k+1 ), lda )
300 IF( rk.LT.lastrk )
THEN
302 IF( vn1( j ).NE.zero )
THEN
307 temp = abs( a( rk, j ) ) / vn1( j )
308 temp = max( zero, ( one+temp )*( one-temp ) )
309 temp2 = temp*( vn1( j ) / vn2( j ) )**2
310 IF( temp2 .LE. tol3z )
THEN
311 vn2( j ) =
REAL( lsticc )
314 vn1( j ) = vn1( j )*sqrt( temp )
333 IF( kb.LT.min( n, m-offset ) )
THEN
334 CALL sgemm(
'No transpose',
'Transpose', m-rk, n-kb, kb, -one,
335 $ a( rk+1, 1 ), lda, f( kb+1, 1 ), ldf, one,
336 $ a( rk+1, kb+1 ), lda )
342 IF( lsticc.GT.0 )
THEN
343 itemp = nint( vn2( lsticc ) )
344 vn1( lsticc ) = snrm2( m-rk, a( rk+1, lsticc ), 1 )
350 vn2( lsticc ) = vn1( lsticc )
subroutine slarfg(N, ALPHA, X, INCX, TAU)
SLARFG generates an elementary reflector (Householder matrix).
subroutine sgemm(TRANSA, TRANSB, M, N, K, ALPHA, A, LDA, B, LDB, BETA, C, LDC)
SGEMM
subroutine sgemv(TRANS, M, N, ALPHA, A, LDA, X, INCX, BETA, Y, INCY)
SGEMV
subroutine slaqps(M, N, OFFSET, NB, KB, A, LDA, JPVT, TAU, VN1, VN2, AUXV, F, LDF)
SLAQPS computes a step of QR factorization with column pivoting of a real m-by-n matrix A by using BL...
subroutine sswap(N, SX, INCX, SY, INCY)
SSWAP