398 SUBROUTINE dlaqp3rk( M, N, NRHS, IOFFSET, NB, ABSTOL,
399 $ RELTOL, KP1, MAXC2NRM, A, LDA, DONE, KB,
400 $ MAXC2NRMK, RELMAXC2NRMK, JPIV, TAU,
401 $ VN1, VN2, AUXV, F, LDF, IWORK, INFO )
410 INTEGER INFO, IOFFSET, KB, KP1, LDA, LDF, M, N,
412 DOUBLE PRECISION ABSTOL, MAXC2NRM, MAXC2NRMK, RELMAXC2NRMK,
416 INTEGER IWORK( * ), JPIV( * )
417 DOUBLE PRECISION A( LDA, * ), AUXV( * ), F( LDF, * ), TAU( * ),
424 DOUBLE PRECISION ZERO, ONE
425 PARAMETER ( ZERO = 0.0d+0, one = 1.0d+0 )
428 INTEGER ITEMP, J, K, MINMNFACT, MINMNUPDT,
430 DOUBLE PRECISION AIK, HUGEVAL, TEMP, TEMP2, TOL3Z
436 INTRINSIC abs, max, min, sqrt
441 DOUBLE PRECISION DLAMCH, DNRM2
442 EXTERNAL disnan, dlamch, idamax, dnrm2
453 minmnfact = min( m-ioffset, n )
454 minmnupdt = min( m-ioffset, n+nrhs )
455 nb = min( nb, minmnfact )
456 tol3z = sqrt( dlamch(
'Epsilon' ) )
457 hugeval = dlamch(
'Overflow' )
466 DO WHILE ( k.LT.nb .AND. lsticc.EQ.0 )
484 kp = ( k-1 ) + idamax( n-k+1, vn1( k ), 1 )
489 maxc2nrmk = vn1( kp )
501 IF( disnan( maxc2nrmk ) )
THEN
518 relmaxc2nrmk = maxc2nrmk
532 IF( nrhs.GT.0 .AND. kb.LT.(m-ioffset) )
THEN
533 CALL dgemm(
'No transpose',
'Transpose',
534 $ m-
IF, nrhs, kb, -one, a( if+1, 1 ), lda,
535 $ f( n+1, 1 ), ldf, one, a( if+1, n+1 ), lda )
555 IF( maxc2nrmk.EQ.zero )
THEN
581 IF( nrhs.GT.0 .AND. kb.LT.(m-ioffset) )
THEN
582 CALL dgemm(
'No transpose',
'Transpose',
583 $ m-
IF, nrhs, kb, -one, a( if+1, 1 ), lda,
584 $ f( n+1, 1 ), ldf, one, a( if+1, n+1 ), lda )
615 IF( info.EQ.0 .AND. maxc2nrmk.GT.hugeval )
THEN
616 info = n + k - 1 + kp
631 relmaxc2nrmk = maxc2nrmk / maxc2nrm
633 IF( maxc2nrmk.LE.abstol .OR. relmaxc2nrmk.LE.reltol )
THEN
657 IF( kb.LT.minmnupdt )
THEN
658 CALL dgemm(
'No transpose',
'Transpose',
659 $ m-
IF, n+nrhs-kb, kb,-one, a( if+1, 1 ), lda,
660 $ f( kb+1, 1 ), ldf, one, a( if+1, kb+1 ), lda )
701 CALL dswap( m, a( 1, kp ), 1, a( 1, k ), 1 )
702 CALL dswap( k-1, f( kp, 1 ), ldf, f( k, 1 ), ldf )
706 jpiv( kp ) = jpiv( k )
714 CALL dgemv(
'No transpose', m-i+1, k-1, -one, a( i, 1 ),
715 $ lda, f( k, 1 ), ldf, one, a( i, k ), 1 )
721 CALL dlarfg( m-i+1, a( i, k ), a( i+1, k ), 1, tau( k ) )
736 IF( disnan( tau(k) ) )
THEN
754 relmaxc2nrmk = tau( k )
768 IF( nrhs.GT.0 .AND. kb.LT.(m-ioffset) )
THEN
769 CALL dgemm(
'No transpose',
'Transpose',
770 $ m-
IF, nrhs, kb, -one, a( if+1, 1 ), lda,
771 $ f( n+1, 1 ), ldf, one, a( if+1, n+1 ), lda )
795 IF( k.LT.n+nrhs )
THEN
796 CALL dgemv(
'Transpose', m-i+1, n+nrhs-k,
797 $ tau( k ), a( i, k+1 ), lda, a( i, k ), 1,
798 $ zero, f( k+1, k ), 1 )
813 CALL dgemv(
'Transpose', m-i+1, k-1, -tau( k ),
814 $ a( i, 1 ), lda, a( i, k ), 1, zero,
817 CALL dgemv(
'No transpose', n+nrhs, k-1, one,
818 $ f( 1, 1 ), ldf, auxv( 1 ), 1, one,
828 IF( k.LT.n+nrhs )
THEN
829 CALL dgemv(
'No transpose', n+nrhs-k, k, -one,
830 $ f( k+1, 1 ), ldf, a( i, 1 ), lda, one,
840 IF( k.LT.minmnfact )
THEN
843 IF( vn1( j ).NE.zero )
THEN
848 temp = abs( a( i, j ) ) / vn1( j )
849 temp = max( zero, ( one+temp )*( one-temp ) )
850 temp2 = temp*( vn1( j ) / vn2( j ) )**2
851 IF( temp2.LE.tol3z )
THEN
860 iwork( j-1 ) = lsticc
867 vn1( j ) = vn1( j )*sqrt( temp )
896 IF( kb.LT.minmnupdt )
THEN
897 CALL dgemm(
'No transpose',
'Transpose',
898 $ m-
IF, n+nrhs-kb, kb, -one, a( if+1, 1 ), lda,
899 $ f( kb+1, 1 ), ldf, one, a( if+1, kb+1 ), lda )
906 DO WHILE( lsticc.GT.0 )
912 itemp = iwork( lsticc-1 )
921 vn1( lsticc ) = dnrm2( m-
IF, a( if+1, lsticc ), 1 )
922 vn2( lsticc ) = vn1( lsticc )
subroutine dlaqp3rk(m, n, nrhs, ioffset, nb, abstol, reltol, kp1, maxc2nrm, a, lda, done, kb, maxc2nrmk, relmaxc2nrmk, jpiv, tau, vn1, vn2, auxv, f, ldf, iwork, info)
DLAQP3RK computes a step of truncated QR factorization with column pivoting of a real m-by-n matrix A...
subroutine dgemm(transa, transb, m, n, k, alpha, a, lda, b, ldb, beta, c, ldc)
DGEMM
subroutine dgemv(trans, m, n, alpha, a, lda, x, incx, beta, y, incy)
DGEMV
subroutine dlarfg(n, alpha, x, incx, tau)
DLARFG generates an elementary reflector (Householder matrix).
subroutine dswap(n, dx, incx, dy, incy)
DSWAP