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

◆ cgeqpf()

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.
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 145 of file cgeqpf.f.

146*
147* -- LAPACK computational routine --
148* -- LAPACK is a software package provided by Univ. of Tennessee, --
149* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
150*
151* .. Scalar Arguments ..
152 INTEGER INFO, LDA, M, N
153* ..
154* .. Array Arguments ..
155 INTEGER JPVT( * )
156 REAL RWORK( * )
157 COMPLEX A( LDA, * ), TAU( * ), WORK( * )
158* ..
159*
160* =====================================================================
161*
162* .. Parameters ..
163 REAL ZERO, ONE
164 parameter( zero = 0.0e+0, one = 1.0e+0 )
165* ..
166* .. Local Scalars ..
167 INTEGER I, ITEMP, J, MA, MN, PVT
168 REAL TEMP, TEMP2, TOL3Z
169 COMPLEX AII
170* ..
171* .. External Subroutines ..
172 EXTERNAL cgeqr2, clarf, clarfg, cswap, cunm2r, xerbla
173* ..
174* .. Intrinsic Functions ..
175 INTRINSIC abs, cmplx, conjg, max, min, sqrt
176* ..
177* .. External Functions ..
178 INTEGER ISAMAX
179 REAL SCNRM2, SLAMCH
180 EXTERNAL isamax, scnrm2, slamch
181* ..
182* .. Executable Statements ..
183*
184* Test the input arguments
185*
186 info = 0
187 IF( m.LT.0 ) THEN
188 info = -1
189 ELSE IF( n.LT.0 ) THEN
190 info = -2
191 ELSE IF( lda.LT.max( 1, m ) ) THEN
192 info = -4
193 END IF
194 IF( info.NE.0 ) THEN
195 CALL xerbla( 'CGEQPF', -info )
196 RETURN
197 END IF
198*
199 mn = min( m, n )
200 tol3z = sqrt(slamch('Epsilon'))
201*
202* Move initial columns up front
203*
204 itemp = 1
205 DO 10 i = 1, n
206 IF( jpvt( i ).NE.0 ) THEN
207 IF( i.NE.itemp ) THEN
208 CALL cswap( m, a( 1, i ), 1, a( 1, itemp ), 1 )
209 jpvt( i ) = jpvt( itemp )
210 jpvt( itemp ) = i
211 ELSE
212 jpvt( i ) = i
213 END IF
214 itemp = itemp + 1
215 ELSE
216 jpvt( i ) = i
217 END IF
218 10 CONTINUE
219 itemp = itemp - 1
220*
221* Compute the QR factorization and update remaining columns
222*
223 IF( itemp.GT.0 ) THEN
224 ma = min( itemp, m )
225 CALL cgeqr2( m, ma, a, lda, tau, work, info )
226 IF( ma.LT.n ) THEN
227 CALL cunm2r( 'Left', 'Conjugate transpose', m, n-ma, ma, a,
228 $ lda, tau, a( 1, ma+1 ), lda, work, info )
229 END IF
230 END IF
231*
232 IF( itemp.LT.mn ) THEN
233*
234* Initialize partial column norms. The first n elements of
235* work store the exact column norms.
236*
237 DO 20 i = itemp + 1, n
238 rwork( i ) = scnrm2( m-itemp, a( itemp+1, i ), 1 )
239 rwork( n+i ) = rwork( i )
240 20 CONTINUE
241*
242* Compute factorization
243*
244 DO 40 i = itemp + 1, mn
245*
246* Determine ith pivot column and swap if necessary
247*
248 pvt = ( i-1 ) + isamax( n-i+1, rwork( i ), 1 )
249*
250 IF( pvt.NE.i ) THEN
251 CALL cswap( m, a( 1, pvt ), 1, a( 1, i ), 1 )
252 itemp = jpvt( pvt )
253 jpvt( pvt ) = jpvt( i )
254 jpvt( i ) = itemp
255 rwork( pvt ) = rwork( i )
256 rwork( n+pvt ) = rwork( n+i )
257 END IF
258*
259* Generate elementary reflector H(i)
260*
261 aii = a( i, i )
262 CALL clarfg( m-i+1, aii, a( min( i+1, m ), i ), 1,
263 $ tau( i ) )
264 a( i, i ) = aii
265*
266 IF( i.LT.n ) THEN
267*
268* Apply H(i) to A(i:m,i+1:n) from the left
269*
270 aii = a( i, i )
271 a( i, i ) = cmplx( one )
272 CALL clarf( 'Left', m-i+1, n-i, a( i, i ), 1,
273 $ conjg( tau( i ) ), a( i, i+1 ), lda, work )
274 a( i, i ) = aii
275 END IF
276*
277* Update partial column norms
278*
279 DO 30 j = i + 1, n
280 IF( rwork( j ).NE.zero ) THEN
281*
282* NOTE: The following 4 lines follow from the analysis in
283* Lapack Working Note 176.
284*
285 temp = abs( a( i, j ) ) / rwork( j )
286 temp = max( zero, ( one+temp )*( one-temp ) )
287 temp2 = temp*( rwork( j ) / rwork( n+j ) )**2
288 IF( temp2 .LE. tol3z ) THEN
289 IF( m-i.GT.0 ) THEN
290 rwork( j ) = scnrm2( m-i, a( i+1, j ), 1 )
291 rwork( n+j ) = rwork( j )
292 ELSE
293 rwork( j ) = zero
294 rwork( n+j ) = zero
295 END IF
296 ELSE
297 rwork( j ) = rwork( j )*sqrt( temp )
298 END IF
299 END IF
300 30 CONTINUE
301*
302 40 CONTINUE
303 END IF
304 RETURN
305*
306* End of CGEQPF
307*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
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:128
integer function isamax(n, sx, incx)
ISAMAX
Definition isamax.f:71
real function slamch(cmach)
SLAMCH
Definition slamch.f:68
subroutine clarf(side, m, n, v, incv, tau, c, ldc, work)
CLARF applies an elementary reflector to a general rectangular matrix.
Definition clarf.f:126
subroutine clarfg(n, alpha, x, incx, tau)
CLARFG generates an elementary reflector (Householder matrix).
Definition clarfg.f:104
real(wp) function scnrm2(n, x, incx)
SCNRM2
Definition scnrm2.f90:90
subroutine cswap(n, cx, incx, cy, incy)
CSWAP
Definition cswap.f:81
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:157
Here is the call graph for this function:
Here is the caller graph for this function: