LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
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.
Date
September 2012
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 151 of file zlaqp2.f.

151 *
152 * -- LAPACK auxiliary routine (version 3.4.2) --
153 * -- LAPACK is a software package provided by Univ. of Tennessee, --
154 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
155 * September 2012
156 *
157 * .. Scalar Arguments ..
158  INTEGER lda, m, n, offset
159 * ..
160 * .. Array Arguments ..
161  INTEGER jpvt( * )
162  DOUBLE PRECISION vn1( * ), vn2( * )
163  COMPLEX*16 a( lda, * ), tau( * ), work( * )
164 * ..
165 *
166 * =====================================================================
167 *
168 * .. Parameters ..
169  DOUBLE PRECISION zero, one
170  COMPLEX*16 cone
171  parameter ( zero = 0.0d+0, one = 1.0d+0,
172  $ cone = ( 1.0d+0, 0.0d+0 ) )
173 * ..
174 * .. Local Scalars ..
175  INTEGER i, itemp, j, mn, offpi, pvt
176  DOUBLE PRECISION temp, temp2, tol3z
177  COMPLEX*16 aii
178 * ..
179 * .. External Subroutines ..
180  EXTERNAL zlarf, zlarfg, zswap
181 * ..
182 * .. Intrinsic Functions ..
183  INTRINSIC abs, dconjg, max, min, sqrt
184 * ..
185 * .. External Functions ..
186  INTEGER idamax
187  DOUBLE PRECISION dlamch, dznrm2
188  EXTERNAL idamax, dlamch, dznrm2
189 * ..
190 * .. Executable Statements ..
191 *
192  mn = min( m-offset, n )
193  tol3z = sqrt(dlamch('Epsilon'))
194 *
195 * Compute factorization.
196 *
197  DO 20 i = 1, mn
198 *
199  offpi = offset + i
200 *
201 * Determine ith pivot column and swap if necessary.
202 *
203  pvt = ( i-1 ) + idamax( n-i+1, vn1( i ), 1 )
204 *
205  IF( pvt.NE.i ) THEN
206  CALL zswap( m, a( 1, pvt ), 1, a( 1, i ), 1 )
207  itemp = jpvt( pvt )
208  jpvt( pvt ) = jpvt( i )
209  jpvt( i ) = itemp
210  vn1( pvt ) = vn1( i )
211  vn2( pvt ) = vn2( i )
212  END IF
213 *
214 * Generate elementary reflector H(i).
215 *
216  IF( offpi.LT.m ) THEN
217  CALL zlarfg( m-offpi+1, a( offpi, i ), a( offpi+1, i ), 1,
218  $ tau( i ) )
219  ELSE
220  CALL zlarfg( 1, a( m, i ), a( m, i ), 1, tau( i ) )
221  END IF
222 *
223  IF( i.LT.n ) THEN
224 *
225 * Apply H(i)**H to A(offset+i:m,i+1:n) from the left.
226 *
227  aii = a( offpi, i )
228  a( offpi, i ) = cone
229  CALL zlarf( 'Left', m-offpi+1, n-i, a( offpi, i ), 1,
230  $ dconjg( tau( i ) ), a( offpi, i+1 ), lda,
231  $ work( 1 ) )
232  a( offpi, i ) = aii
233  END IF
234 *
235 * Update partial column norms.
236 *
237  DO 10 j = i + 1, n
238  IF( vn1( j ).NE.zero ) THEN
239 *
240 * NOTE: The following 4 lines follow from the analysis in
241 * Lapack Working Note 176.
242 *
243  temp = one - ( abs( a( offpi, j ) ) / vn1( j ) )**2
244  temp = max( temp, zero )
245  temp2 = temp*( vn1( j ) / vn2( j ) )**2
246  IF( temp2 .LE. tol3z ) THEN
247  IF( offpi.LT.m ) THEN
248  vn1( j ) = dznrm2( m-offpi, a( offpi+1, j ), 1 )
249  vn2( j ) = vn1( j )
250  ELSE
251  vn1( j ) = zero
252  vn2( j ) = zero
253  END IF
254  ELSE
255  vn1( j ) = vn1( j )*sqrt( temp )
256  END IF
257  END IF
258  10 CONTINUE
259 *
260  20 CONTINUE
261 *
262  RETURN
263 *
264 * End of ZLAQP2
265 *
integer function idamax(N, DX, INCX)
IDAMAX
Definition: idamax.f:53
double precision function dlamch(CMACH)
DLAMCH
Definition: dlamch.f:65
subroutine zlarfg(N, ALPHA, X, INCX, TAU)
ZLARFG generates an elementary reflector (Householder matrix).
Definition: zlarfg.f:108
subroutine zswap(N, ZX, INCX, ZY, INCY)
ZSWAP
Definition: zswap.f:52
double precision function dznrm2(N, X, INCX)
DZNRM2
Definition: dznrm2.f:56
subroutine zlarf(SIDE, M, N, V, INCV, TAU, C, LDC, WORK)
ZLARF applies an elementary reflector to a general rectangular matrix.
Definition: zlarf.f:130

Here is the call graph for this function:

Here is the caller graph for this function: