LAPACK 3.12.1
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches

◆ dlaqps()

subroutine dlaqps ( integer m,
integer n,
integer offset,
integer nb,
integer kb,
double precision, dimension( lda, * ) a,
integer lda,
integer, dimension( * ) jpvt,
double precision, dimension( * ) tau,
double precision, dimension( * ) vn1,
double precision, dimension( * ) vn2,
double precision, dimension( * ) auxv,
double precision, dimension( ldf, * ) f,
integer ldf )

DLAQPS computes a step of QR factorization with column pivoting of a real m-by-n matrix A by using BLAS level 3.

Download DLAQPS + dependencies [TGZ] [ZIP] [TXT]

Purpose:
!>
!> DLAQPS computes a step of QR factorization with column pivoting
!> of a real M-by-N matrix A by using Blas-3.  It tries to factorize
!> NB columns from A starting from the row OFFSET+1, and updates all
!> of the matrix with Blas-3 xGEMM.
!>
!> In some cases, due to catastrophic cancellations, it cannot
!> factorize NB columns.  Hence, the actual number of factorized
!> columns is returned in KB.
!>
!> Block A(1:OFFSET,1:N) is accordingly pivoted, but not factorized.
!> 
Parameters
[in]M
!>          M is INTEGER
!>          The number of rows of the matrix A. M >= 0.
!> 
[in]N
!>          N is INTEGER
!>          The number of columns of the matrix A. N >= 0
!> 
[in]OFFSET
!>          OFFSET is INTEGER
!>          The number of rows of A that have been factorized in
!>          previous steps.
!> 
[in]NB
!>          NB is INTEGER
!>          The number of columns to factorize.
!> 
[out]KB
!>          KB is INTEGER
!>          The number of columns actually factorized.
!> 
[in,out]A
!>          A is DOUBLE PRECISION array, dimension (LDA,N)
!>          On entry, the M-by-N matrix A.
!>          On exit, block A(OFFSET+1:M,1:KB) is the triangular
!>          factor obtained and block A(1:OFFSET,1:N) has been
!>          accordingly pivoted, but no factorized.
!>          The rest of the matrix, block A(OFFSET+1:M,KB+1:N) has
!>          been updated.
!> 
[in]LDA
!>          LDA is INTEGER
!>          The leading dimension of the array A. LDA >= max(1,M).
!> 
[in,out]JPVT
!>          JPVT is INTEGER array, dimension (N)
!>          JPVT(I) = K <==> Column K of the full matrix A has been
!>          permuted into position I in AP.
!> 
[out]TAU
!>          TAU is DOUBLE PRECISION array, dimension (KB)
!>          The scalar factors of the elementary reflectors.
!> 
[in,out]VN1
!>          VN1 is DOUBLE PRECISION array, dimension (N)
!>          The vector with the partial column norms.
!> 
[in,out]VN2
!>          VN2 is DOUBLE PRECISION array, dimension (N)
!>          The vector with the exact column norms.
!> 
[in,out]AUXV
!>          AUXV is DOUBLE PRECISION array, dimension (NB)
!>          Auxiliary vector.
!> 
[in,out]F
!>          F is DOUBLE PRECISION array, dimension (LDF,NB)
!>          Matrix F**T = L*Y**T*A.
!> 
[in]LDF
!>          LDF is INTEGER
!>          The leading dimension of the array F. LDF >= max(1,N).
!> 
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Contributors:
G. Quintana-Orti, Depto. de Informatica, Universidad Jaime I, Spain X. Sun, Computer Science Dept., Duke University, USA
Partial column norm updating strategy modified on April 2011 Z. Drmac and Z. Bujanovic, Dept. of Mathematics, University of Zagreb, Croatia.
References:
LAPACK Working Note 176 [PDF]

Definition at line 171 of file dlaqps.f.

174*
175* -- LAPACK auxiliary routine --
176* -- LAPACK is a software package provided by Univ. of Tennessee, --
177* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
178*
179* .. Scalar Arguments ..
180 INTEGER KB, LDA, LDF, M, N, NB, OFFSET
181* ..
182* .. Array Arguments ..
183 INTEGER JPVT( * )
184 DOUBLE PRECISION A( LDA, * ), AUXV( * ), F( LDF, * ), TAU( * ),
185 $ VN1( * ), VN2( * )
186* ..
187*
188* =====================================================================
189*
190* .. Parameters ..
191 DOUBLE PRECISION ZERO, ONE
192 parameter( zero = 0.0d+0, one = 1.0d+0 )
193* ..
194* .. Local Scalars ..
195 INTEGER ITEMP, J, K, LASTRK, LSTICC, PVT, RK
196 DOUBLE PRECISION AKK, TEMP, TEMP2, TOL3Z
197* ..
198* .. External Subroutines ..
199 EXTERNAL dgemm, dgemv, dlarfg, dswap
200* ..
201* .. Intrinsic Functions ..
202 INTRINSIC abs, dble, max, min, nint, sqrt
203* ..
204* .. External Functions ..
205 INTEGER IDAMAX
206 DOUBLE PRECISION DLAMCH, DNRM2
207 EXTERNAL idamax, dlamch, dnrm2
208* ..
209* .. Executable Statements ..
210*
211 lastrk = min( m, n+offset )
212 lsticc = 0
213 k = 0
214 tol3z = sqrt(dlamch('Epsilon'))
215*
216* Beginning of while loop.
217*
218 10 CONTINUE
219 IF( ( k.LT.nb ) .AND. ( lsticc.EQ.0 ) ) THEN
220 k = k + 1
221 rk = offset + k
222*
223* Determine ith pivot column and swap if necessary
224*
225 pvt = ( k-1 ) + idamax( n-k+1, vn1( k ), 1 )
226 IF( pvt.NE.k ) THEN
227 CALL dswap( m, a( 1, pvt ), 1, a( 1, k ), 1 )
228 CALL dswap( k-1, f( pvt, 1 ), ldf, f( k, 1 ), ldf )
229 itemp = jpvt( pvt )
230 jpvt( pvt ) = jpvt( k )
231 jpvt( k ) = itemp
232 vn1( pvt ) = vn1( k )
233 vn2( pvt ) = vn2( k )
234 END IF
235*
236* Apply previous Householder reflectors to column K:
237* A(RK:M,K) := A(RK:M,K) - A(RK:M,1:K-1)*F(K,1:K-1)**T.
238*
239 IF( k.GT.1 ) THEN
240 CALL dgemv( 'No transpose', m-rk+1, k-1, -one, a( rk,
241 $ 1 ),
242 $ lda, f( k, 1 ), ldf, one, a( rk, k ), 1 )
243 END IF
244*
245* Generate elementary reflector H(k).
246*
247 IF( rk.LT.m ) THEN
248 CALL dlarfg( m-rk+1, a( rk, k ), a( rk+1, k ), 1,
249 $ tau( k ) )
250 ELSE
251 CALL dlarfg( 1, a( rk, k ), a( rk, k ), 1, tau( k ) )
252 END IF
253*
254 akk = a( rk, k )
255 a( rk, k ) = one
256*
257* Compute Kth column of F:
258*
259* Compute F(K+1:N,K) := tau(K)*A(RK:M,K+1:N)**T*A(RK:M,K).
260*
261 IF( k.LT.n ) THEN
262 CALL dgemv( 'Transpose', m-rk+1, n-k, tau( k ),
263 $ a( rk, k+1 ), lda, a( rk, k ), 1, zero,
264 $ f( k+1, k ), 1 )
265 END IF
266*
267* Padding F(1:K,K) with zeros.
268*
269 DO 20 j = 1, k
270 f( j, k ) = zero
271 20 CONTINUE
272*
273* Incremental updating of F:
274* F(1:N,K) := F(1:N,K) - tau(K)*F(1:N,1:K-1)*A(RK:M,1:K-1)**T
275* *A(RK:M,K).
276*
277 IF( k.GT.1 ) THEN
278 CALL dgemv( 'Transpose', m-rk+1, k-1, -tau( k ), a( rk,
279 $ 1 ),
280 $ lda, a( rk, k ), 1, zero, auxv( 1 ), 1 )
281*
282 CALL dgemv( 'No transpose', n, k-1, one, f( 1, 1 ), ldf,
283 $ auxv( 1 ), 1, one, f( 1, k ), 1 )
284 END IF
285*
286* Update the current row of A:
287* A(RK,K+1:N) := A(RK,K+1:N) - A(RK,1:K)*F(K+1:N,1:K)**T.
288*
289 IF( k.LT.n ) THEN
290 CALL dgemv( 'No transpose', n-k, k, -one, f( k+1, 1 ),
291 $ ldf,
292 $ a( rk, 1 ), lda, one, a( rk, k+1 ), lda )
293 END IF
294*
295* Update partial column norms.
296*
297 IF( rk.LT.lastrk ) THEN
298 DO 30 j = k + 1, n
299 IF( vn1( j ).NE.zero ) THEN
300*
301* NOTE: The following 4 lines follow from the analysis in
302* Lapack Working Note 176.
303*
304 temp = abs( a( rk, j ) ) / vn1( j )
305 temp = max( zero, ( one+temp )*( one-temp ) )
306 temp2 = temp*( vn1( j ) / vn2( j ) )**2
307 IF( temp2 .LE. tol3z ) THEN
308 vn2( j ) = dble( lsticc )
309 lsticc = j
310 ELSE
311 vn1( j ) = vn1( j )*sqrt( temp )
312 END IF
313 END IF
314 30 CONTINUE
315 END IF
316*
317 a( rk, k ) = akk
318*
319* End of while loop.
320*
321 GO TO 10
322 END IF
323 kb = k
324 rk = offset + kb
325*
326* Apply the block reflector to the rest of the matrix:
327* A(OFFSET+KB+1:M,KB+1:N) := A(OFFSET+KB+1:M,KB+1:N) -
328* A(OFFSET+KB+1:M,1:KB)*F(KB+1:N,1:KB)**T.
329*
330 IF( kb.LT.min( n, m-offset ) ) THEN
331 CALL dgemm( 'No transpose', 'Transpose', m-rk, n-kb, kb,
332 $ -one,
333 $ a( rk+1, 1 ), lda, f( kb+1, 1 ), ldf, one,
334 $ a( rk+1, kb+1 ), lda )
335 END IF
336*
337* Recomputation of difficult columns.
338*
339 40 CONTINUE
340 IF( lsticc.GT.0 ) THEN
341 itemp = nint( vn2( lsticc ) )
342 vn1( lsticc ) = dnrm2( m-rk, a( rk+1, lsticc ), 1 )
343*
344* NOTE: The computation of VN1( LSTICC ) relies on the fact that
345* SNRM2 does not fail on vectors with norm below the value of
346* SQRT(DLAMCH('S'))
347*
348 vn2( lsticc ) = vn1( lsticc )
349 lsticc = itemp
350 GO TO 40
351 END IF
352*
353 RETURN
354*
355* End of DLAQPS
356*
subroutine dgemm(transa, transb, m, n, k, alpha, a, lda, b, ldb, beta, c, ldc)
DGEMM
Definition dgemm.f:188
subroutine dgemv(trans, m, n, alpha, a, lda, x, incx, beta, y, incy)
DGEMV
Definition dgemv.f:158
integer function idamax(n, dx, incx)
IDAMAX
Definition idamax.f:71
double precision function dlamch(cmach)
DLAMCH
Definition dlamch.f:69
subroutine dlarfg(n, alpha, x, incx, tau)
DLARFG generates an elementary reflector (Householder matrix).
Definition dlarfg.f:104
real(wp) function dnrm2(n, x, incx)
DNRM2
Definition dnrm2.f90:89
subroutine dswap(n, dx, incx, dy, incy)
DSWAP
Definition dswap.f:82
Here is the call graph for this function:
Here is the caller graph for this function: