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

◆ zlaqp2()

subroutine zlaqp2 ( integer  m,
integer  n,
integer  offset,
complex*16, dimension( lda, * )  a,
integer  lda,
integer, dimension( * )  jpvt,
complex*16, dimension( * )  tau,
double precision, dimension( * )  vn1,
double precision, dimension( * )  vn2,
complex*16, dimension( * )  work 
)

ZLAQP2 computes a QR factorization with column pivoting of the matrix block.

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

Purpose:
 ZLAQP2 computes a QR factorization with column pivoting of
 the block A(OFFSET+1:M,1:N).
 The 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 the matrix A that must be pivoted
          but no factorized. OFFSET >= 0.
[in,out]A
          A is COMPLEX*16 array, dimension (LDA,N)
          On entry, the M-by-N matrix A.
          On exit, the upper triangle of block A(OFFSET+1:M,1:N) is
          the triangular factor obtained; the elements in block
          A(OFFSET+1:M,1:N) below the diagonal, together with the
          array TAU, represent the orthogonal matrix Q as a product of
          elementary reflectors. Block A(1:OFFSET,1:N) has been
          accordingly pivoted, but no factorized.
[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)
          On entry, if JPVT(i) .ne. 0, the i-th column of A is permuted
          to the front of A*P (a leading column); if JPVT(i) = 0,
          the i-th column of A is a free column.
          On exit, if JPVT(i) = k, then the i-th column of A*P
          was the k-th column of A.
[out]TAU
          TAU is COMPLEX*16 array, dimension (min(M,N))
          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.
[out]WORK
          WORK is COMPLEX*16 array, dimension (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 147 of file zlaqp2.f.

149*
150* -- LAPACK auxiliary routine --
151* -- LAPACK is a software package provided by Univ. of Tennessee, --
152* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
153*
154* .. Scalar Arguments ..
155 INTEGER LDA, M, N, OFFSET
156* ..
157* .. Array Arguments ..
158 INTEGER JPVT( * )
159 DOUBLE PRECISION VN1( * ), VN2( * )
160 COMPLEX*16 A( LDA, * ), TAU( * ), WORK( * )
161* ..
162*
163* =====================================================================
164*
165* .. Parameters ..
166 DOUBLE PRECISION ZERO, ONE
167 COMPLEX*16 CONE
168 parameter( zero = 0.0d+0, one = 1.0d+0,
169 $ cone = ( 1.0d+0, 0.0d+0 ) )
170* ..
171* .. Local Scalars ..
172 INTEGER I, ITEMP, J, MN, OFFPI, PVT
173 DOUBLE PRECISION TEMP, TEMP2, TOL3Z
174 COMPLEX*16 AII
175* ..
176* .. External Subroutines ..
177 EXTERNAL zlarf, zlarfg, zswap
178* ..
179* .. Intrinsic Functions ..
180 INTRINSIC abs, dconjg, max, min, sqrt
181* ..
182* .. External Functions ..
183 INTEGER IDAMAX
184 DOUBLE PRECISION DLAMCH, DZNRM2
185 EXTERNAL idamax, dlamch, dznrm2
186* ..
187* .. Executable Statements ..
188*
189 mn = min( m-offset, n )
190 tol3z = sqrt(dlamch('Epsilon'))
191*
192* Compute factorization.
193*
194 DO 20 i = 1, mn
195*
196 offpi = offset + i
197*
198* Determine ith pivot column and swap if necessary.
199*
200 pvt = ( i-1 ) + idamax( n-i+1, vn1( i ), 1 )
201*
202 IF( pvt.NE.i ) THEN
203 CALL zswap( m, a( 1, pvt ), 1, a( 1, i ), 1 )
204 itemp = jpvt( pvt )
205 jpvt( pvt ) = jpvt( i )
206 jpvt( i ) = itemp
207 vn1( pvt ) = vn1( i )
208 vn2( pvt ) = vn2( i )
209 END IF
210*
211* Generate elementary reflector H(i).
212*
213 IF( offpi.LT.m ) THEN
214 CALL zlarfg( m-offpi+1, a( offpi, i ), a( offpi+1, i ), 1,
215 $ tau( i ) )
216 ELSE
217 CALL zlarfg( 1, a( m, i ), a( m, i ), 1, tau( i ) )
218 END IF
219*
220 IF( i.LT.n ) THEN
221*
222* Apply H(i)**H to A(offset+i:m,i+1:n) from the left.
223*
224 aii = a( offpi, i )
225 a( offpi, i ) = cone
226 CALL zlarf( 'Left', m-offpi+1, n-i, a( offpi, i ), 1,
227 $ dconjg( tau( i ) ), a( offpi, i+1 ), lda,
228 $ work( 1 ) )
229 a( offpi, i ) = aii
230 END IF
231*
232* Update partial column norms.
233*
234 DO 10 j = i + 1, n
235 IF( vn1( j ).NE.zero ) THEN
236*
237* NOTE: The following 4 lines follow from the analysis in
238* Lapack Working Note 176.
239*
240 temp = one - ( abs( a( offpi, j ) ) / vn1( j ) )**2
241 temp = max( temp, zero )
242 temp2 = temp*( vn1( j ) / vn2( j ) )**2
243 IF( temp2 .LE. tol3z ) THEN
244 IF( offpi.LT.m ) THEN
245 vn1( j ) = dznrm2( m-offpi, a( offpi+1, j ), 1 )
246 vn2( j ) = vn1( j )
247 ELSE
248 vn1( j ) = zero
249 vn2( j ) = zero
250 END IF
251 ELSE
252 vn1( j ) = vn1( j )*sqrt( temp )
253 END IF
254 END IF
255 10 CONTINUE
256*
257 20 CONTINUE
258*
259 RETURN
260*
261* End of ZLAQP2
262*
integer function idamax(n, dx, incx)
IDAMAX
Definition idamax.f:71
double precision function dlamch(cmach)
DLAMCH
Definition dlamch.f:69
subroutine zlarf(side, m, n, v, incv, tau, c, ldc, work)
ZLARF applies an elementary reflector to a general rectangular matrix.
Definition zlarf.f:128
subroutine zlarfg(n, alpha, x, incx, tau)
ZLARFG generates an elementary reflector (Householder matrix).
Definition zlarfg.f:106
real(wp) function dznrm2(n, x, incx)
DZNRM2
Definition dznrm2.f90:90
subroutine zswap(n, zx, incx, zy, incy)
ZSWAP
Definition zswap.f:81
Here is the call graph for this function:
Here is the caller graph for this function: