b18fa9ebfd6d37d25bdde624a9bdb0f91ac54759
[platform/upstream/lapack.git] / SRC / zlaqps.f
1 *> \brief \b ZLAQPS computes a step of QR factorization with column pivoting of a real m-by-n matrix A by using BLAS level 3.
2 *
3 *  =========== DOCUMENTATION ===========
4 *
5 * Online html documentation available at 
6 *            http://www.netlib.org/lapack/explore-html/ 
7 *
8 *> \htmlonly
9 *> Download ZLAQPS + dependencies 
10 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/zlaqps.f"> 
11 *> [TGZ]</a> 
12 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/zlaqps.f"> 
13 *> [ZIP]</a> 
14 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/zlaqps.f"> 
15 *> [TXT]</a>
16 *> \endhtmlonly 
17 *
18 *  Definition:
19 *  ===========
20 *
21 *       SUBROUTINE ZLAQPS( M, N, OFFSET, NB, KB, A, LDA, JPVT, TAU, VN1,
22 *                          VN2, AUXV, F, LDF )
23
24 *       .. Scalar Arguments ..
25 *       INTEGER            KB, LDA, LDF, M, N, NB, OFFSET
26 *       ..
27 *       .. Array Arguments ..
28 *       INTEGER            JPVT( * )
29 *       DOUBLE PRECISION   VN1( * ), VN2( * )
30 *       COMPLEX*16         A( LDA, * ), AUXV( * ), F( LDF, * ), TAU( * )
31 *       ..
32 *  
33 *
34 *> \par Purpose:
35 *  =============
36 *>
37 *> \verbatim
38 *>
39 *> ZLAQPS computes a step of QR factorization with column pivoting
40 *> of a complex M-by-N matrix A by using Blas-3.  It tries to factorize
41 *> NB columns from A starting from the row OFFSET+1, and updates all
42 *> of the matrix with Blas-3 xGEMM.
43 *>
44 *> In some cases, due to catastrophic cancellations, it cannot
45 *> factorize NB columns.  Hence, the actual number of factorized
46 *> columns is returned in KB.
47 *>
48 *> Block A(1:OFFSET,1:N) is accordingly pivoted, but not factorized.
49 *> \endverbatim
50 *
51 *  Arguments:
52 *  ==========
53 *
54 *> \param[in] M
55 *> \verbatim
56 *>          M is INTEGER
57 *>          The number of rows of the matrix A. M >= 0.
58 *> \endverbatim
59 *>
60 *> \param[in] N
61 *> \verbatim
62 *>          N is INTEGER
63 *>          The number of columns of the matrix A. N >= 0
64 *> \endverbatim
65 *>
66 *> \param[in] OFFSET
67 *> \verbatim
68 *>          OFFSET is INTEGER
69 *>          The number of rows of A that have been factorized in
70 *>          previous steps.
71 *> \endverbatim
72 *>
73 *> \param[in] NB
74 *> \verbatim
75 *>          NB is INTEGER
76 *>          The number of columns to factorize.
77 *> \endverbatim
78 *>
79 *> \param[out] KB
80 *> \verbatim
81 *>          KB is INTEGER
82 *>          The number of columns actually factorized.
83 *> \endverbatim
84 *>
85 *> \param[in,out] A
86 *> \verbatim
87 *>          A is COMPLEX*16 array, dimension (LDA,N)
88 *>          On entry, the M-by-N matrix A.
89 *>          On exit, block A(OFFSET+1:M,1:KB) is the triangular
90 *>          factor obtained and block A(1:OFFSET,1:N) has been
91 *>          accordingly pivoted, but no factorized.
92 *>          The rest of the matrix, block A(OFFSET+1:M,KB+1:N) has
93 *>          been updated.
94 *> \endverbatim
95 *>
96 *> \param[in] LDA
97 *> \verbatim
98 *>          LDA is INTEGER
99 *>          The leading dimension of the array A. LDA >= max(1,M).
100 *> \endverbatim
101 *>
102 *> \param[in,out] JPVT
103 *> \verbatim
104 *>          JPVT is INTEGER array, dimension (N)
105 *>          JPVT(I) = K <==> Column K of the full matrix A has been
106 *>          permuted into position I in AP.
107 *> \endverbatim
108 *>
109 *> \param[out] TAU
110 *> \verbatim
111 *>          TAU is COMPLEX*16 array, dimension (KB)
112 *>          The scalar factors of the elementary reflectors.
113 *> \endverbatim
114 *>
115 *> \param[in,out] VN1
116 *> \verbatim
117 *>          VN1 is DOUBLE PRECISION array, dimension (N)
118 *>          The vector with the partial column norms.
119 *> \endverbatim
120 *>
121 *> \param[in,out] VN2
122 *> \verbatim
123 *>          VN2 is DOUBLE PRECISION array, dimension (N)
124 *>          The vector with the exact column norms.
125 *> \endverbatim
126 *>
127 *> \param[in,out] AUXV
128 *> \verbatim
129 *>          AUXV is COMPLEX*16 array, dimension (NB)
130 *>          Auxiliar vector.
131 *> \endverbatim
132 *>
133 *> \param[in,out] F
134 *> \verbatim
135 *>          F is COMPLEX*16 array, dimension (LDF,NB)
136 *>          Matrix F**H = L * Y**H * A.
137 *> \endverbatim
138 *>
139 *> \param[in] LDF
140 *> \verbatim
141 *>          LDF is INTEGER
142 *>          The leading dimension of the array F. LDF >= max(1,N).
143 *> \endverbatim
144 *
145 *  Authors:
146 *  ========
147 *
148 *> \author Univ. of Tennessee 
149 *> \author Univ. of California Berkeley 
150 *> \author Univ. of Colorado Denver 
151 *> \author NAG Ltd. 
152 *
153 *> \date September 2012
154 *
155 *> \ingroup complex16OTHERauxiliary
156 *
157 *> \par Contributors:
158 *  ==================
159 *>
160 *>    G. Quintana-Orti, Depto. de Informatica, Universidad Jaime I, Spain
161 *>    X. Sun, Computer Science Dept., Duke University, USA
162 *> \n
163 *>  Partial column norm updating strategy modified on April 2011
164 *>    Z. Drmac and Z. Bujanovic, Dept. of Mathematics,
165 *>    University of Zagreb, Croatia.
166 *
167 *> \par References:
168 *  ================
169 *>
170 *> LAPACK Working Note 176
171 *
172 *> \htmlonly
173 *> <a href="http://www.netlib.org/lapack/lawnspdf/lawn176.pdf">[PDF]</a> 
174 *> \endhtmlonly 
175 *
176 *  =====================================================================
177       SUBROUTINE ZLAQPS( M, N, OFFSET, NB, KB, A, LDA, JPVT, TAU, VN1,
178      $                   VN2, AUXV, F, LDF )
179 *
180 *  -- LAPACK auxiliary routine (version 3.4.2) --
181 *  -- LAPACK is a software package provided by Univ. of Tennessee,    --
182 *  -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
183 *     September 2012
184 *
185 *     .. Scalar Arguments ..
186       INTEGER            KB, LDA, LDF, M, N, NB, OFFSET
187 *     ..
188 *     .. Array Arguments ..
189       INTEGER            JPVT( * )
190       DOUBLE PRECISION   VN1( * ), VN2( * )
191       COMPLEX*16         A( LDA, * ), AUXV( * ), F( LDF, * ), TAU( * )
192 *     ..
193 *
194 *  =====================================================================
195 *
196 *     .. Parameters ..
197       DOUBLE PRECISION   ZERO, ONE
198       COMPLEX*16         CZERO, CONE
199       PARAMETER          ( ZERO = 0.0D+0, ONE = 1.0D+0,
200      $                   CZERO = ( 0.0D+0, 0.0D+0 ),
201      $                   CONE = ( 1.0D+0, 0.0D+0 ) )
202 *     ..
203 *     .. Local Scalars ..
204       INTEGER            ITEMP, J, K, LASTRK, LSTICC, PVT, RK
205       DOUBLE PRECISION   TEMP, TEMP2, TOL3Z
206       COMPLEX*16         AKK
207 *     ..
208 *     .. External Subroutines ..
209       EXTERNAL           ZGEMM, ZGEMV, ZLARFG, ZSWAP
210 *     ..
211 *     .. Intrinsic Functions ..
212       INTRINSIC          ABS, DBLE, DCONJG, MAX, MIN, NINT, SQRT
213 *     ..
214 *     .. External Functions ..
215       INTEGER            IDAMAX
216       DOUBLE PRECISION   DLAMCH, DZNRM2
217       EXTERNAL           IDAMAX, DLAMCH, DZNRM2
218 *     ..
219 *     .. Executable Statements ..
220 *
221       LASTRK = MIN( M, N+OFFSET )
222       LSTICC = 0
223       K = 0
224       TOL3Z = SQRT(DLAMCH('Epsilon'))
225 *
226 *     Beginning of while loop.
227 *
228    10 CONTINUE
229       IF( ( K.LT.NB ) .AND. ( LSTICC.EQ.0 ) ) THEN
230          K = K + 1
231          RK = OFFSET + K
232 *
233 *        Determine ith pivot column and swap if necessary
234 *
235          PVT = ( K-1 ) + IDAMAX( N-K+1, VN1( K ), 1 )
236          IF( PVT.NE.K ) THEN
237             CALL ZSWAP( M, A( 1, PVT ), 1, A( 1, K ), 1 )
238             CALL ZSWAP( K-1, F( PVT, 1 ), LDF, F( K, 1 ), LDF )
239             ITEMP = JPVT( PVT )
240             JPVT( PVT ) = JPVT( K )
241             JPVT( K ) = ITEMP
242             VN1( PVT ) = VN1( K )
243             VN2( PVT ) = VN2( K )
244          END IF
245 *
246 *        Apply previous Householder reflectors to column K:
247 *        A(RK:M,K) := A(RK:M,K) - A(RK:M,1:K-1)*F(K,1:K-1)**H.
248 *
249          IF( K.GT.1 ) THEN
250             DO 20 J = 1, K - 1
251                F( K, J ) = DCONJG( F( K, J ) )
252    20       CONTINUE
253             CALL ZGEMV( 'No transpose', M-RK+1, K-1, -CONE, A( RK, 1 ),
254      $                  LDA, F( K, 1 ), LDF, CONE, A( RK, K ), 1 )
255             DO 30 J = 1, K - 1
256                F( K, J ) = DCONJG( F( K, J ) )
257    30       CONTINUE
258          END IF
259 *
260 *        Generate elementary reflector H(k).
261 *
262          IF( RK.LT.M ) THEN
263             CALL ZLARFG( M-RK+1, A( RK, K ), A( RK+1, K ), 1, TAU( K ) )
264          ELSE
265             CALL ZLARFG( 1, A( RK, K ), A( RK, K ), 1, TAU( K ) )
266          END IF
267 *
268          AKK = A( RK, K )
269          A( RK, K ) = CONE
270 *
271 *        Compute Kth column of F:
272 *
273 *        Compute  F(K+1:N,K) := tau(K)*A(RK:M,K+1:N)**H*A(RK:M,K).
274 *
275          IF( K.LT.N ) THEN
276             CALL ZGEMV( 'Conjugate transpose', M-RK+1, N-K, TAU( K ),
277      $                  A( RK, K+1 ), LDA, A( RK, K ), 1, CZERO,
278      $                  F( K+1, K ), 1 )
279          END IF
280 *
281 *        Padding F(1:K,K) with zeros.
282 *
283          DO 40 J = 1, K
284             F( J, K ) = CZERO
285    40    CONTINUE
286 *
287 *        Incremental updating of F:
288 *        F(1:N,K) := F(1:N,K) - tau(K)*F(1:N,1:K-1)*A(RK:M,1:K-1)**H
289 *                    *A(RK:M,K).
290 *
291          IF( K.GT.1 ) THEN
292             CALL ZGEMV( 'Conjugate transpose', M-RK+1, K-1, -TAU( K ),
293      $                  A( RK, 1 ), LDA, A( RK, K ), 1, CZERO,
294      $                  AUXV( 1 ), 1 )
295 *
296             CALL ZGEMV( 'No transpose', N, K-1, CONE, F( 1, 1 ), LDF,
297      $                  AUXV( 1 ), 1, CONE, F( 1, K ), 1 )
298          END IF
299 *
300 *        Update the current row of A:
301 *        A(RK,K+1:N) := A(RK,K+1:N) - A(RK,1:K)*F(K+1:N,1:K)**H.
302 *
303          IF( K.LT.N ) THEN
304             CALL ZGEMM( 'No transpose', 'Conjugate transpose', 1, N-K,
305      $                  K, -CONE, A( RK, 1 ), LDA, F( K+1, 1 ), LDF,
306      $                  CONE, A( RK, K+1 ), LDA )
307          END IF
308 *
309 *        Update partial column norms.
310 *
311          IF( RK.LT.LASTRK ) THEN
312             DO 50 J = K + 1, N
313                IF( VN1( J ).NE.ZERO ) THEN
314 *
315 *                 NOTE: The following 4 lines follow from the analysis in
316 *                 Lapack Working Note 176.
317 *
318                   TEMP = ABS( A( RK, J ) ) / VN1( J )
319                   TEMP = MAX( ZERO, ( ONE+TEMP )*( ONE-TEMP ) )
320                   TEMP2 = TEMP*( VN1( J ) / VN2( J ) )**2
321                   IF( TEMP2 .LE. TOL3Z ) THEN
322                      VN2( J ) = DBLE( LSTICC )
323                      LSTICC = J
324                   ELSE
325                      VN1( J ) = VN1( J )*SQRT( TEMP )
326                   END IF
327                END IF
328    50       CONTINUE
329          END IF
330 *
331          A( RK, K ) = AKK
332 *
333 *        End of while loop.
334 *
335          GO TO 10
336       END IF
337       KB = K
338       RK = OFFSET + KB
339 *
340 *     Apply the block reflector to the rest of the matrix:
341 *     A(OFFSET+KB+1:M,KB+1:N) := A(OFFSET+KB+1:M,KB+1:N) -
342 *                         A(OFFSET+KB+1:M,1:KB)*F(KB+1:N,1:KB)**H.
343 *
344       IF( KB.LT.MIN( N, M-OFFSET ) ) THEN
345          CALL ZGEMM( 'No transpose', 'Conjugate transpose', M-RK, N-KB,
346      $               KB, -CONE, A( RK+1, 1 ), LDA, F( KB+1, 1 ), LDF,
347      $               CONE, A( RK+1, KB+1 ), LDA )
348       END IF
349 *
350 *     Recomputation of difficult columns.
351 *
352    60 CONTINUE
353       IF( LSTICC.GT.0 ) THEN
354          ITEMP = NINT( VN2( LSTICC ) )
355          VN1( LSTICC ) = DZNRM2( M-RK, A( RK+1, LSTICC ), 1 )
356 *
357 *        NOTE: The computation of VN1( LSTICC ) relies on the fact that 
358 *        SNRM2 does not fail on vectors with norm below the value of
359 *        SQRT(DLAMCH('S')) 
360 *
361          VN2( LSTICC ) = VN1( LSTICC )
362          LSTICC = ITEMP
363          GO TO 60
364       END IF
365 *
366       RETURN
367 *
368 *     End of ZLAQPS
369 *
370       END