1 | SUBROUTINE DLAQPS( M, N, OFFSET, NB, KB, A, LDA, JPVT, TAU, VN1,
|
---|
2 | $ VN2, AUXV, F, LDF )
|
---|
3 | *
|
---|
4 | * -- LAPACK auxiliary routine (version 3.1) --
|
---|
5 | * Univ. of Tennessee, Univ. of California Berkeley and NAG Ltd..
|
---|
6 | * November 2006
|
---|
7 | *
|
---|
8 | * .. Scalar Arguments ..
|
---|
9 | INTEGER KB, LDA, LDF, M, N, NB, OFFSET
|
---|
10 | * ..
|
---|
11 | * .. Array Arguments ..
|
---|
12 | INTEGER JPVT( * )
|
---|
13 | DOUBLE PRECISION A( LDA, * ), AUXV( * ), F( LDF, * ), TAU( * ),
|
---|
14 | $ VN1( * ), VN2( * )
|
---|
15 | * ..
|
---|
16 | *
|
---|
17 | * Purpose
|
---|
18 | * =======
|
---|
19 | *
|
---|
20 | * DLAQPS computes a step of QR factorization with column pivoting
|
---|
21 | * of a real M-by-N matrix A by using Blas-3. It tries to factorize
|
---|
22 | * NB columns from A starting from the row OFFSET+1, and updates all
|
---|
23 | * of the matrix with Blas-3 xGEMM.
|
---|
24 | *
|
---|
25 | * In some cases, due to catastrophic cancellations, it cannot
|
---|
26 | * factorize NB columns. Hence, the actual number of factorized
|
---|
27 | * columns is returned in KB.
|
---|
28 | *
|
---|
29 | * Block A(1:OFFSET,1:N) is accordingly pivoted, but not factorized.
|
---|
30 | *
|
---|
31 | * Arguments
|
---|
32 | * =========
|
---|
33 | *
|
---|
34 | * M (input) INTEGER
|
---|
35 | * The number of rows of the matrix A. M >= 0.
|
---|
36 | *
|
---|
37 | * N (input) INTEGER
|
---|
38 | * The number of columns of the matrix A. N >= 0
|
---|
39 | *
|
---|
40 | * OFFSET (input) INTEGER
|
---|
41 | * The number of rows of A that have been factorized in
|
---|
42 | * previous steps.
|
---|
43 | *
|
---|
44 | * NB (input) INTEGER
|
---|
45 | * The number of columns to factorize.
|
---|
46 | *
|
---|
47 | * KB (output) INTEGER
|
---|
48 | * The number of columns actually factorized.
|
---|
49 | *
|
---|
50 | * A (input/output) DOUBLE PRECISION array, dimension (LDA,N)
|
---|
51 | * On entry, the M-by-N matrix A.
|
---|
52 | * On exit, block A(OFFSET+1:M,1:KB) is the triangular
|
---|
53 | * factor obtained and block A(1:OFFSET,1:N) has been
|
---|
54 | * accordingly pivoted, but no factorized.
|
---|
55 | * The rest of the matrix, block A(OFFSET+1:M,KB+1:N) has
|
---|
56 | * been updated.
|
---|
57 | *
|
---|
58 | * LDA (input) INTEGER
|
---|
59 | * The leading dimension of the array A. LDA >= max(1,M).
|
---|
60 | *
|
---|
61 | * JPVT (input/output) INTEGER array, dimension (N)
|
---|
62 | * JPVT(I) = K <==> Column K of the full matrix A has been
|
---|
63 | * permuted into position I in AP.
|
---|
64 | *
|
---|
65 | * TAU (output) DOUBLE PRECISION array, dimension (KB)
|
---|
66 | * The scalar factors of the elementary reflectors.
|
---|
67 | *
|
---|
68 | * VN1 (input/output) DOUBLE PRECISION array, dimension (N)
|
---|
69 | * The vector with the partial column norms.
|
---|
70 | *
|
---|
71 | * VN2 (input/output) DOUBLE PRECISION array, dimension (N)
|
---|
72 | * The vector with the exact column norms.
|
---|
73 | *
|
---|
74 | * AUXV (input/output) DOUBLE PRECISION array, dimension (NB)
|
---|
75 | * Auxiliar vector.
|
---|
76 | *
|
---|
77 | * F (input/output) DOUBLE PRECISION array, dimension (LDF,NB)
|
---|
78 | * Matrix F' = L*Y'*A.
|
---|
79 | *
|
---|
80 | * LDF (input) INTEGER
|
---|
81 | * The leading dimension of the array F. LDF >= max(1,N).
|
---|
82 | *
|
---|
83 | * Further Details
|
---|
84 | * ===============
|
---|
85 | *
|
---|
86 | * Based on contributions by
|
---|
87 | * G. Quintana-Orti, Depto. de Informatica, Universidad Jaime I, Spain
|
---|
88 | * X. Sun, Computer Science Dept., Duke University, USA
|
---|
89 | *
|
---|
90 | * Partial column norm updating strategy modified by
|
---|
91 | * Z. Drmac and Z. Bujanovic, Dept. of Mathematics,
|
---|
92 | * University of Zagreb, Croatia.
|
---|
93 | * June 2006.
|
---|
94 | * For more details see LAPACK Working Note 176.
|
---|
95 | * =====================================================================
|
---|
96 | *
|
---|
97 | * .. Parameters ..
|
---|
98 | DOUBLE PRECISION ZERO, ONE
|
---|
99 | PARAMETER ( ZERO = 0.0D+0, ONE = 1.0D+0 )
|
---|
100 | * ..
|
---|
101 | * .. Local Scalars ..
|
---|
102 | INTEGER ITEMP, J, K, LASTRK, LSTICC, PVT, RK
|
---|
103 | DOUBLE PRECISION AKK, TEMP, TEMP2, TOL3Z
|
---|
104 | * ..
|
---|
105 | * .. External Subroutines ..
|
---|
106 | EXTERNAL DGEMM, DGEMV, DLARFG, DSWAP
|
---|
107 | * ..
|
---|
108 | * .. Intrinsic Functions ..
|
---|
109 | INTRINSIC ABS, DBLE, MAX, MIN, NINT, SQRT
|
---|
110 | * ..
|
---|
111 | * .. External Functions ..
|
---|
112 | INTEGER IDAMAX
|
---|
113 | DOUBLE PRECISION DLAMCH, DNRM2
|
---|
114 | EXTERNAL IDAMAX, DLAMCH, DNRM2
|
---|
115 | * ..
|
---|
116 | * .. Executable Statements ..
|
---|
117 | *
|
---|
118 | LASTRK = MIN( M, N+OFFSET )
|
---|
119 | LSTICC = 0
|
---|
120 | K = 0
|
---|
121 | TOL3Z = SQRT(DLAMCH('Epsilon'))
|
---|
122 | *
|
---|
123 | * Beginning of while loop.
|
---|
124 | *
|
---|
125 | 10 CONTINUE
|
---|
126 | IF( ( K.LT.NB ) .AND. ( LSTICC.EQ.0 ) ) THEN
|
---|
127 | K = K + 1
|
---|
128 | RK = OFFSET + K
|
---|
129 | *
|
---|
130 | * Determine ith pivot column and swap if necessary
|
---|
131 | *
|
---|
132 | PVT = ( K-1 ) + IDAMAX( N-K+1, VN1( K ), 1 )
|
---|
133 | IF( PVT.NE.K ) THEN
|
---|
134 | CALL DSWAP( M, A( 1, PVT ), 1, A( 1, K ), 1 )
|
---|
135 | CALL DSWAP( K-1, F( PVT, 1 ), LDF, F( K, 1 ), LDF )
|
---|
136 | ITEMP = JPVT( PVT )
|
---|
137 | JPVT( PVT ) = JPVT( K )
|
---|
138 | JPVT( K ) = ITEMP
|
---|
139 | VN1( PVT ) = VN1( K )
|
---|
140 | VN2( PVT ) = VN2( K )
|
---|
141 | END IF
|
---|
142 | *
|
---|
143 | * Apply previous Householder reflectors to column K:
|
---|
144 | * A(RK:M,K) := A(RK:M,K) - A(RK:M,1:K-1)*F(K,1:K-1)'.
|
---|
145 | *
|
---|
146 | IF( K.GT.1 ) THEN
|
---|
147 | CALL DGEMV( 'No transpose', M-RK+1, K-1, -ONE, A( RK, 1 ),
|
---|
148 | $ LDA, F( K, 1 ), LDF, ONE, A( RK, K ), 1 )
|
---|
149 | END IF
|
---|
150 | *
|
---|
151 | * Generate elementary reflector H(k).
|
---|
152 | *
|
---|
153 | IF( RK.LT.M ) THEN
|
---|
154 | CALL DLARFG( M-RK+1, A( RK, K ), A( RK+1, K ), 1, TAU( K ) )
|
---|
155 | ELSE
|
---|
156 | CALL DLARFG( 1, A( RK, K ), A( RK, K ), 1, TAU( K ) )
|
---|
157 | END IF
|
---|
158 | *
|
---|
159 | AKK = A( RK, K )
|
---|
160 | A( RK, K ) = ONE
|
---|
161 | *
|
---|
162 | * Compute Kth column of F:
|
---|
163 | *
|
---|
164 | * Compute F(K+1:N,K) := tau(K)*A(RK:M,K+1:N)'*A(RK:M,K).
|
---|
165 | *
|
---|
166 | IF( K.LT.N ) THEN
|
---|
167 | CALL DGEMV( 'Transpose', M-RK+1, N-K, TAU( K ),
|
---|
168 | $ A( RK, K+1 ), LDA, A( RK, K ), 1, ZERO,
|
---|
169 | $ F( K+1, K ), 1 )
|
---|
170 | END IF
|
---|
171 | *
|
---|
172 | * Padding F(1:K,K) with zeros.
|
---|
173 | *
|
---|
174 | DO 20 J = 1, K
|
---|
175 | F( J, K ) = ZERO
|
---|
176 | 20 CONTINUE
|
---|
177 | *
|
---|
178 | * Incremental updating of F:
|
---|
179 | * F(1:N,K) := F(1:N,K) - tau(K)*F(1:N,1:K-1)*A(RK:M,1:K-1)'
|
---|
180 | * *A(RK:M,K).
|
---|
181 | *
|
---|
182 | IF( K.GT.1 ) THEN
|
---|
183 | CALL DGEMV( 'Transpose', M-RK+1, K-1, -TAU( K ), A( RK, 1 ),
|
---|
184 | $ LDA, A( RK, K ), 1, ZERO, AUXV( 1 ), 1 )
|
---|
185 | *
|
---|
186 | CALL DGEMV( 'No transpose', N, K-1, ONE, F( 1, 1 ), LDF,
|
---|
187 | $ AUXV( 1 ), 1, ONE, F( 1, K ), 1 )
|
---|
188 | END IF
|
---|
189 | *
|
---|
190 | * Update the current row of A:
|
---|
191 | * A(RK,K+1:N) := A(RK,K+1:N) - A(RK,1:K)*F(K+1:N,1:K)'.
|
---|
192 | *
|
---|
193 | IF( K.LT.N ) THEN
|
---|
194 | CALL DGEMV( 'No transpose', N-K, K, -ONE, F( K+1, 1 ), LDF,
|
---|
195 | $ A( RK, 1 ), LDA, ONE, A( RK, K+1 ), LDA )
|
---|
196 | END IF
|
---|
197 | *
|
---|
198 | * Update partial column norms.
|
---|
199 | *
|
---|
200 | IF( RK.LT.LASTRK ) THEN
|
---|
201 | DO 30 J = K + 1, N
|
---|
202 | IF( VN1( J ).NE.ZERO ) THEN
|
---|
203 | *
|
---|
204 | * NOTE: The following 4 lines follow from the analysis in
|
---|
205 | * Lapack Working Note 176.
|
---|
206 | *
|
---|
207 | TEMP = ABS( A( RK, J ) ) / VN1( J )
|
---|
208 | TEMP = MAX( ZERO, ( ONE+TEMP )*( ONE-TEMP ) )
|
---|
209 | TEMP2 = TEMP*( VN1( J ) / VN2( J ) )**2
|
---|
210 | IF( TEMP2 .LE. TOL3Z ) THEN
|
---|
211 | VN2( J ) = DBLE( LSTICC )
|
---|
212 | LSTICC = J
|
---|
213 | ELSE
|
---|
214 | VN1( J ) = VN1( J )*SQRT( TEMP )
|
---|
215 | END IF
|
---|
216 | END IF
|
---|
217 | 30 CONTINUE
|
---|
218 | END IF
|
---|
219 | *
|
---|
220 | A( RK, K ) = AKK
|
---|
221 | *
|
---|
222 | * End of while loop.
|
---|
223 | *
|
---|
224 | GO TO 10
|
---|
225 | END IF
|
---|
226 | KB = K
|
---|
227 | RK = OFFSET + KB
|
---|
228 | *
|
---|
229 | * Apply the block reflector to the rest of the matrix:
|
---|
230 | * A(OFFSET+KB+1:M,KB+1:N) := A(OFFSET+KB+1:M,KB+1:N) -
|
---|
231 | * A(OFFSET+KB+1:M,1:KB)*F(KB+1:N,1:KB)'.
|
---|
232 | *
|
---|
233 | IF( KB.LT.MIN( N, M-OFFSET ) ) THEN
|
---|
234 | CALL DGEMM( 'No transpose', 'Transpose', M-RK, N-KB, KB, -ONE,
|
---|
235 | $ A( RK+1, 1 ), LDA, F( KB+1, 1 ), LDF, ONE,
|
---|
236 | $ A( RK+1, KB+1 ), LDA )
|
---|
237 | END IF
|
---|
238 | *
|
---|
239 | * Recomputation of difficult columns.
|
---|
240 | *
|
---|
241 | 40 CONTINUE
|
---|
242 | IF( LSTICC.GT.0 ) THEN
|
---|
243 | ITEMP = NINT( VN2( LSTICC ) )
|
---|
244 | VN1( LSTICC ) = DNRM2( M-RK, A( RK+1, LSTICC ), 1 )
|
---|
245 | *
|
---|
246 | * NOTE: The computation of VN1( LSTICC ) relies on the fact that
|
---|
247 | * SNRM2 does not fail on vectors with norm below the value of
|
---|
248 | * SQRT(DLAMCH('S'))
|
---|
249 | *
|
---|
250 | VN2( LSTICC ) = VN1( LSTICC )
|
---|
251 | LSTICC = ITEMP
|
---|
252 | GO TO 40
|
---|
253 | END IF
|
---|
254 | *
|
---|
255 | RETURN
|
---|
256 | *
|
---|
257 | * End of DLAQPS
|
---|
258 | *
|
---|
259 | END
|
---|