001:       SUBROUTINE ZLAQPS( M, N, OFFSET, NB, KB, A, LDA, JPVT, TAU, VN1,
002:      $                   VN2, AUXV, F, LDF )
003: *
004: *  -- LAPACK auxiliary routine (version 3.2) --
005: *     Univ. of Tennessee, Univ. of California Berkeley and NAG Ltd..
006: *     November 2006
007: *
008: *     .. Scalar Arguments ..
009:       INTEGER            KB, LDA, LDF, M, N, NB, OFFSET
010: *     ..
011: *     .. Array Arguments ..
012:       INTEGER            JPVT( * )
013:       DOUBLE PRECISION   VN1( * ), VN2( * )
014:       COMPLEX*16         A( LDA, * ), AUXV( * ), F( LDF, * ), TAU( * )
015: *     ..
016: *
017: *  Purpose
018: *  =======
019: *
020: *  ZLAQPS computes a step of QR factorization with column pivoting
021: *  of a complex M-by-N matrix A by using Blas-3.  It tries to factorize
022: *  NB columns from A starting from the row OFFSET+1, and updates all
023: *  of the matrix with Blas-3 xGEMM.
024: *
025: *  In some cases, due to catastrophic cancellations, it cannot
026: *  factorize NB columns.  Hence, the actual number of factorized
027: *  columns is returned in KB.
028: *
029: *  Block A(1:OFFSET,1:N) is accordingly pivoted, but not factorized.
030: *
031: *  Arguments
032: *  =========
033: *
034: *  M       (input) INTEGER
035: *          The number of rows of the matrix A. M >= 0.
036: *
037: *  N       (input) INTEGER
038: *          The number of columns of the matrix A. N >= 0
039: *
040: *  OFFSET  (input) INTEGER
041: *          The number of rows of A that have been factorized in
042: *          previous steps.
043: *
044: *  NB      (input) INTEGER
045: *          The number of columns to factorize.
046: *
047: *  KB      (output) INTEGER
048: *          The number of columns actually factorized.
049: *
050: *  A       (input/output) COMPLEX*16 array, dimension (LDA,N)
051: *          On entry, the M-by-N matrix A.
052: *          On exit, block A(OFFSET+1:M,1:KB) is the triangular
053: *          factor obtained and block A(1:OFFSET,1:N) has been
054: *          accordingly pivoted, but no factorized.
055: *          The rest of the matrix, block A(OFFSET+1:M,KB+1:N) has
056: *          been updated.
057: *
058: *  LDA     (input) INTEGER
059: *          The leading dimension of the array A. LDA >= max(1,M).
060: *
061: *  JPVT    (input/output) INTEGER array, dimension (N)
062: *          JPVT(I) = K <==> Column K of the full matrix A has been
063: *          permuted into position I in AP.
064: *
065: *  TAU     (output) COMPLEX*16 array, dimension (KB)
066: *          The scalar factors of the elementary reflectors.
067: *
068: *  VN1     (input/output) DOUBLE PRECISION array, dimension (N)
069: *          The vector with the partial column norms.
070: *
071: *  VN2     (input/output) DOUBLE PRECISION array, dimension (N)
072: *          The vector with the exact column norms.
073: *
074: *  AUXV    (input/output) COMPLEX*16 array, dimension (NB)
075: *          Auxiliar vector.
076: *
077: *  F       (input/output) COMPLEX*16 array, dimension (LDF,NB)
078: *          Matrix F' = L*Y'*A.
079: *
080: *  LDF     (input) INTEGER
081: *          The leading dimension of the array F. LDF >= max(1,N).
082: *
083: *  Further Details
084: *  ===============
085: *
086: *  Based on contributions by
087: *    G. Quintana-Orti, Depto. de Informatica, Universidad Jaime I, Spain
088: *    X. Sun, Computer Science Dept., Duke University, USA
089: *
090: *  =====================================================================
091: *
092: *     .. Parameters ..
093:       DOUBLE PRECISION   ZERO, ONE
094:       COMPLEX*16         CZERO, CONE
095:       PARAMETER          ( ZERO = 0.0D+0, ONE = 1.0D+0,
096:      $                   CZERO = ( 0.0D+0, 0.0D+0 ),
097:      $                   CONE = ( 1.0D+0, 0.0D+0 ) )
098: *     ..
099: *     .. Local Scalars ..
100:       INTEGER            ITEMP, J, K, LASTRK, LSTICC, PVT, RK
101:       DOUBLE PRECISION   TEMP, TEMP2, TOL3Z
102:       COMPLEX*16         AKK
103: *     ..
104: *     .. External Subroutines ..
105:       EXTERNAL           ZGEMM, ZGEMV, ZLARFP, ZSWAP
106: *     ..
107: *     .. Intrinsic Functions ..
108:       INTRINSIC          ABS, DBLE, DCONJG, MAX, MIN, NINT, SQRT
109: *     ..
110: *     .. External Functions ..
111:       INTEGER            IDAMAX
112:       DOUBLE PRECISION   DLAMCH, DZNRM2
113:       EXTERNAL           IDAMAX, DLAMCH, DZNRM2
114: *     ..
115: *     .. Executable Statements ..
116: *
117:       LASTRK = MIN( M, N+OFFSET )
118:       LSTICC = 0
119:       K = 0
120:       TOL3Z = SQRT(DLAMCH('Epsilon'))
121: *
122: *     Beginning of while loop.
123: *
124:    10 CONTINUE
125:       IF( ( K.LT.NB ) .AND. ( LSTICC.EQ.0 ) ) THEN
126:          K = K + 1
127:          RK = OFFSET + K
128: *
129: *        Determine ith pivot column and swap if necessary
130: *
131:          PVT = ( K-1 ) + IDAMAX( N-K+1, VN1( K ), 1 )
132:          IF( PVT.NE.K ) THEN
133:             CALL ZSWAP( M, A( 1, PVT ), 1, A( 1, K ), 1 )
134:             CALL ZSWAP( K-1, F( PVT, 1 ), LDF, F( K, 1 ), LDF )
135:             ITEMP = JPVT( PVT )
136:             JPVT( PVT ) = JPVT( K )
137:             JPVT( K ) = ITEMP
138:             VN1( PVT ) = VN1( K )
139:             VN2( PVT ) = VN2( K )
140:          END IF
141: *
142: *        Apply previous Householder reflectors to column K:
143: *        A(RK:M,K) := A(RK:M,K) - A(RK:M,1:K-1)*F(K,1:K-1)'.
144: *
145:          IF( K.GT.1 ) THEN
146:             DO 20 J = 1, K - 1
147:                F( K, J ) = DCONJG( F( K, J ) )
148:    20       CONTINUE
149:             CALL ZGEMV( 'No transpose', M-RK+1, K-1, -CONE, A( RK, 1 ),
150:      $                  LDA, F( K, 1 ), LDF, CONE, A( RK, K ), 1 )
151:             DO 30 J = 1, K - 1
152:                F( K, J ) = DCONJG( F( K, J ) )
153:    30       CONTINUE
154:          END IF
155: *
156: *        Generate elementary reflector H(k).
157: *
158:          IF( RK.LT.M ) THEN
159:             CALL ZLARFP( M-RK+1, A( RK, K ), A( RK+1, K ), 1, TAU( K ) )
160:          ELSE
161:             CALL ZLARFP( 1, A( RK, K ), A( RK, K ), 1, TAU( K ) )
162:          END IF
163: *
164:          AKK = A( RK, K )
165:          A( RK, K ) = CONE
166: *
167: *        Compute Kth column of F:
168: *
169: *        Compute  F(K+1:N,K) := tau(K)*A(RK:M,K+1:N)'*A(RK:M,K).
170: *
171:          IF( K.LT.N ) THEN
172:             CALL ZGEMV( 'Conjugate transpose', M-RK+1, N-K, TAU( K ),
173:      $                  A( RK, K+1 ), LDA, A( RK, K ), 1, CZERO,
174:      $                  F( K+1, K ), 1 )
175:          END IF
176: *
177: *        Padding F(1:K,K) with zeros.
178: *
179:          DO 40 J = 1, K
180:             F( J, K ) = CZERO
181:    40    CONTINUE
182: *
183: *        Incremental updating of F:
184: *        F(1:N,K) := F(1:N,K) - tau(K)*F(1:N,1:K-1)*A(RK:M,1:K-1)'
185: *                    *A(RK:M,K).
186: *
187:          IF( K.GT.1 ) THEN
188:             CALL ZGEMV( 'Conjugate transpose', M-RK+1, K-1, -TAU( K ),
189:      $                  A( RK, 1 ), LDA, A( RK, K ), 1, CZERO,
190:      $                  AUXV( 1 ), 1 )
191: *
192:             CALL ZGEMV( 'No transpose', N, K-1, CONE, F( 1, 1 ), LDF,
193:      $                  AUXV( 1 ), 1, CONE, F( 1, K ), 1 )
194:          END IF
195: *
196: *        Update the current row of A:
197: *        A(RK,K+1:N) := A(RK,K+1:N) - A(RK,1:K)*F(K+1:N,1:K)'.
198: *
199:          IF( K.LT.N ) THEN
200:             CALL ZGEMM( 'No transpose', 'Conjugate transpose', 1, N-K,
201:      $                  K, -CONE, A( RK, 1 ), LDA, F( K+1, 1 ), LDF,
202:      $                  CONE, A( RK, K+1 ), LDA )
203:          END IF
204: *
205: *        Update partial column norms.
206: *
207:          IF( RK.LT.LASTRK ) THEN
208:             DO 50 J = K + 1, N
209:                IF( VN1( J ).NE.ZERO ) THEN
210: *
211: *                 NOTE: The following 4 lines follow from the analysis in
212: *                 Lapack Working Note 176.
213: *
214:                   TEMP = ABS( A( RK, J ) ) / VN1( J )
215:                   TEMP = MAX( ZERO, ( ONE+TEMP )*( ONE-TEMP ) )
216:                   TEMP2 = TEMP*( VN1( J ) / VN2( J ) )**2
217:                   IF( TEMP2 .LE. TOL3Z ) THEN
218:                      VN2( J ) = DBLE( LSTICC )
219:                      LSTICC = J
220:                   ELSE
221:                      VN1( J ) = VN1( J )*SQRT( TEMP )
222:                   END IF
223:                END IF
224:    50       CONTINUE
225:          END IF
226: *
227:          A( RK, K ) = AKK
228: *
229: *        End of while loop.
230: *
231:          GO TO 10
232:       END IF
233:       KB = K
234:       RK = OFFSET + KB
235: *
236: *     Apply the block reflector to the rest of the matrix:
237: *     A(OFFSET+KB+1:M,KB+1:N) := A(OFFSET+KB+1:M,KB+1:N) -
238: *                         A(OFFSET+KB+1:M,1:KB)*F(KB+1:N,1:KB)'.
239: *
240:       IF( KB.LT.MIN( N, M-OFFSET ) ) THEN
241:          CALL ZGEMM( 'No transpose', 'Conjugate transpose', M-RK, N-KB,
242:      $               KB, -CONE, A( RK+1, 1 ), LDA, F( KB+1, 1 ), LDF,
243:      $               CONE, A( RK+1, KB+1 ), LDA )
244:       END IF
245: *
246: *     Recomputation of difficult columns.
247: *
248:    60 CONTINUE
249:       IF( LSTICC.GT.0 ) THEN
250:          ITEMP = NINT( VN2( LSTICC ) )
251:          VN1( LSTICC ) = DZNRM2( M-RK, A( RK+1, LSTICC ), 1 )
252: *
253: *        NOTE: The computation of VN1( LSTICC ) relies on the fact that 
254: *        SNRM2 does not fail on vectors with norm below the value of
255: *        SQRT(DLAMCH('S')) 
256: *
257:          VN2( LSTICC ) = VN1( LSTICC )
258:          LSTICC = ITEMP
259:          GO TO 60
260:       END IF
261: *
262:       RETURN
263: *
264: *     End of ZLAQPS
265: *
266:       END
267: