175 SUBROUTINE zlaqps( M, N, OFFSET, NB, KB, A, LDA, JPVT, TAU, VN1,
176 $ VN2, AUXV, F, LDF )
183 INTEGER KB, LDA, LDF, M, N, NB, OFFSET
187 DOUBLE PRECISION VN1( * ), VN2( * )
188 COMPLEX*16 A( LDA, * ), AUXV( * ), F( LDF, * ), TAU( * )
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 ) )
201 INTEGER ITEMP, J, K, LASTRK, LSTICC, PVT, RK
202 DOUBLE PRECISION TEMP, TEMP2, TOL3Z
209 INTRINSIC abs, dble, dconjg, max, min, nint, sqrt
213 DOUBLE PRECISION DLAMCH, DZNRM2
214 EXTERNAL idamax, dlamch, dznrm2
218 lastrk = min( m, n+offset )
221 tol3z = sqrt(dlamch(
'Epsilon'))
226 IF( ( k.LT.nb ) .AND. ( lsticc.EQ.0 ) )
THEN
232 pvt = ( k-1 ) + idamax( n-k+1, vn1( k ), 1 )
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 )
237 jpvt( pvt ) = jpvt( k )
239 vn1( pvt ) = vn1( k )
240 vn2( pvt ) = vn2( k )
248 f( k, j ) = dconjg( f( k, j ) )
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 )
253 f( k, j ) = dconjg( f( k, j ) )
260 CALL zlarfg( m-rk+1, a( rk, k ), a( rk+1, k ), 1, tau( k ) )
262 CALL zlarfg( 1, a( rk, k ), a( rk, k ), 1, tau( k ) )
273 CALL zgemv(
'Conjugate transpose', m-rk+1, n-k, tau( k ),
274 $ a( rk, k+1 ), lda, a( rk, k ), 1, czero,
289 CALL zgemv(
'Conjugate transpose', m-rk+1, k-1, -tau( k ),
290 $ a( rk, 1 ), lda, a( rk, k ), 1, czero,
293 CALL zgemv(
'No transpose', n, k-1, cone, f( 1, 1 ), ldf,
294 $ auxv( 1 ), 1, cone, f( 1, k ), 1 )
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 )
308 IF( rk.LT.lastrk )
THEN
310 IF( vn1( j ).NE.zero )
THEN
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 )
322 vn1( j ) = vn1( j )*sqrt( temp )
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 )
350 IF( lsticc.GT.0 )
THEN
351 itemp = nint( vn2( lsticc ) )
352 vn1( lsticc ) = dznrm2( m-rk, a( rk+1, lsticc ), 1 )
358 vn2( lsticc ) = vn1( lsticc )
subroutine zgemm(transa, transb, m, n, k, alpha, a, lda, b, ldb, beta, c, ldc)
ZGEMM
subroutine zgemv(trans, m, n, alpha, a, lda, x, incx, beta, y, incy)
ZGEMV
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...
subroutine zlarfg(n, alpha, x, incx, tau)
ZLARFG generates an elementary reflector (Householder matrix).
subroutine zswap(n, zx, incx, zy, incy)
ZSWAP