LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine cgeqpf ( integer  M,
integer  N,
complex, dimension( lda, * )  A,
integer  LDA,
integer, dimension( * )  JPVT,
complex, dimension( * )  TAU,
complex, dimension( * )  WORK,
real, dimension( * )  RWORK,
integer  INFO 
)

CGEQPF

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

Purpose:
 This routine is deprecated and has been replaced by routine CGEQP3.

 CGEQPF computes a QR factorization with column pivoting of a
 complex M-by-N matrix A: A*P = Q*R.
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,out]A
          A is COMPLEX array, dimension (LDA,N)
          On entry, the M-by-N matrix A.
          On exit, the upper triangle of the array contains the
          min(M,N)-by-N upper triangular matrix R; the elements
          below the diagonal, together with the array TAU,
          represent the unitary matrix Q as a product of
          min(m,n) elementary reflectors.
[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 array, dimension (min(M,N))
          The scalar factors of the elementary reflectors.
[out]WORK
          WORK is COMPLEX array, dimension (N)
[out]RWORK
          RWORK is REAL array, dimension (2*N)
[out]INFO
          INFO is INTEGER
          = 0:  successful exit
          < 0:  if INFO = -i, the i-th argument had an illegal value
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
November 2011
Further Details:
  The matrix Q is represented as a product of elementary reflectors

     Q = H(1) H(2) . . . H(n)

  Each H(i) has the form

     H = I - tau * v * v**H

  where tau is a complex scalar, and v is a complex vector with
  v(1:i-1) = 0 and v(i) = 1; v(i+1:m) is stored on exit in A(i+1:m,i).

  The matrix P is represented in jpvt as follows: If
     jpvt(j) = i
  then the jth column of P is the ith canonical unit vector.

  Partial column norm updating strategy modified by
    Z. Drmac and Z. Bujanovic, Dept. of Mathematics,
    University of Zagreb, Croatia.
  -- April 2011                                                      --
  For more details see LAPACK Working Note 176.

Definition at line 150 of file cgeqpf.f.

150 *
151 * -- LAPACK computational routine (version 3.4.0) --
152 * -- LAPACK is a software package provided by Univ. of Tennessee, --
153 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
154 * November 2011
155 *
156 * .. Scalar Arguments ..
157  INTEGER info, lda, m, n
158 * ..
159 * .. Array Arguments ..
160  INTEGER jpvt( * )
161  REAL rwork( * )
162  COMPLEX a( lda, * ), tau( * ), work( * )
163 * ..
164 *
165 * =====================================================================
166 *
167 * .. Parameters ..
168  REAL zero, one
169  parameter ( zero = 0.0e+0, one = 1.0e+0 )
170 * ..
171 * .. Local Scalars ..
172  INTEGER i, itemp, j, ma, mn, pvt
173  REAL temp, temp2, tol3z
174  COMPLEX aii
175 * ..
176 * .. External Subroutines ..
177  EXTERNAL cgeqr2, clarf, clarfg, cswap, cunm2r, xerbla
178 * ..
179 * .. Intrinsic Functions ..
180  INTRINSIC abs, cmplx, conjg, max, min, sqrt
181 * ..
182 * .. External Functions ..
183  INTEGER isamax
184  REAL scnrm2, slamch
185  EXTERNAL isamax, scnrm2, slamch
186 * ..
187 * .. Executable Statements ..
188 *
189 * Test the input arguments
190 *
191  info = 0
192  IF( m.LT.0 ) THEN
193  info = -1
194  ELSE IF( n.LT.0 ) THEN
195  info = -2
196  ELSE IF( lda.LT.max( 1, m ) ) THEN
197  info = -4
198  END IF
199  IF( info.NE.0 ) THEN
200  CALL xerbla( 'CGEQPF', -info )
201  RETURN
202  END IF
203 *
204  mn = min( m, n )
205  tol3z = sqrt(slamch('Epsilon'))
206 *
207 * Move initial columns up front
208 *
209  itemp = 1
210  DO 10 i = 1, n
211  IF( jpvt( i ).NE.0 ) THEN
212  IF( i.NE.itemp ) THEN
213  CALL cswap( m, a( 1, i ), 1, a( 1, itemp ), 1 )
214  jpvt( i ) = jpvt( itemp )
215  jpvt( itemp ) = i
216  ELSE
217  jpvt( i ) = i
218  END IF
219  itemp = itemp + 1
220  ELSE
221  jpvt( i ) = i
222  END IF
223  10 CONTINUE
224  itemp = itemp - 1
225 *
226 * Compute the QR factorization and update remaining columns
227 *
228  IF( itemp.GT.0 ) THEN
229  ma = min( itemp, m )
230  CALL cgeqr2( m, ma, a, lda, tau, work, info )
231  IF( ma.LT.n ) THEN
232  CALL cunm2r( 'Left', 'Conjugate transpose', m, n-ma, ma, a,
233  $ lda, tau, a( 1, ma+1 ), lda, work, info )
234  END IF
235  END IF
236 *
237  IF( itemp.LT.mn ) THEN
238 *
239 * Initialize partial column norms. The first n elements of
240 * work store the exact column norms.
241 *
242  DO 20 i = itemp + 1, n
243  rwork( i ) = scnrm2( m-itemp, a( itemp+1, i ), 1 )
244  rwork( n+i ) = rwork( i )
245  20 CONTINUE
246 *
247 * Compute factorization
248 *
249  DO 40 i = itemp + 1, mn
250 *
251 * Determine ith pivot column and swap if necessary
252 *
253  pvt = ( i-1 ) + isamax( n-i+1, rwork( i ), 1 )
254 *
255  IF( pvt.NE.i ) THEN
256  CALL cswap( m, a( 1, pvt ), 1, a( 1, i ), 1 )
257  itemp = jpvt( pvt )
258  jpvt( pvt ) = jpvt( i )
259  jpvt( i ) = itemp
260  rwork( pvt ) = rwork( i )
261  rwork( n+pvt ) = rwork( n+i )
262  END IF
263 *
264 * Generate elementary reflector H(i)
265 *
266  aii = a( i, i )
267  CALL clarfg( m-i+1, aii, a( min( i+1, m ), i ), 1,
268  $ tau( i ) )
269  a( i, i ) = aii
270 *
271  IF( i.LT.n ) THEN
272 *
273 * Apply H(i) to A(i:m,i+1:n) from the left
274 *
275  aii = a( i, i )
276  a( i, i ) = cmplx( one )
277  CALL clarf( 'Left', m-i+1, n-i, a( i, i ), 1,
278  $ conjg( tau( i ) ), a( i, i+1 ), lda, work )
279  a( i, i ) = aii
280  END IF
281 *
282 * Update partial column norms
283 *
284  DO 30 j = i + 1, n
285  IF( rwork( j ).NE.zero ) THEN
286 *
287 * NOTE: The following 4 lines follow from the analysis in
288 * Lapack Working Note 176.
289 *
290  temp = abs( a( i, j ) ) / rwork( j )
291  temp = max( zero, ( one+temp )*( one-temp ) )
292  temp2 = temp*( rwork( j ) / rwork( n+j ) )**2
293  IF( temp2 .LE. tol3z ) THEN
294  IF( m-i.GT.0 ) THEN
295  rwork( j ) = scnrm2( m-i, a( i+1, j ), 1 )
296  rwork( n+j ) = rwork( j )
297  ELSE
298  rwork( j ) = zero
299  rwork( n+j ) = zero
300  END IF
301  ELSE
302  rwork( j ) = rwork( j )*sqrt( temp )
303  END IF
304  END IF
305  30 CONTINUE
306 *
307  40 CONTINUE
308  END IF
309  RETURN
310 *
311 * End of CGEQPF
312 *
real function scnrm2(N, X, INCX)
SCNRM2
Definition: scnrm2.f:56
subroutine cgeqr2(M, N, A, LDA, TAU, WORK, INFO)
CGEQR2 computes the QR factorization of a general rectangular matrix using an unblocked algorithm...
Definition: cgeqr2.f:123
integer function isamax(N, SX, INCX)
ISAMAX
Definition: isamax.f:53
subroutine cunm2r(SIDE, TRANS, M, N, K, A, LDA, TAU, C, LDC, WORK, INFO)
CUNM2R multiplies a general matrix by the unitary matrix from a QR factorization determined by cgeqrf...
Definition: cunm2r.f:161
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine clarf(SIDE, M, N, V, INCV, TAU, C, LDC, WORK)
CLARF applies an elementary reflector to a general rectangular matrix.
Definition: clarf.f:130
subroutine cswap(N, CX, INCX, CY, INCY)
CSWAP
Definition: cswap.f:52
real function slamch(CMACH)
SLAMCH
Definition: slamch.f:69
subroutine clarfg(N, ALPHA, X, INCX, TAU)
CLARFG generates an elementary reflector (Householder matrix).
Definition: clarfg.f:108

Here is the call graph for this function:

Here is the caller graph for this function: