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

◆ cgeqp3()

subroutine cgeqp3 ( integer  m,
integer  n,
complex, dimension( lda, * )  a,
integer  lda,
integer, dimension( * )  jpvt,
complex, dimension( * )  tau,
complex, dimension( * )  work,
integer  lwork,
real, dimension( * )  rwork,
integer  info 
)

CGEQP3

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

Purpose:
 CGEQP3 computes a QR factorization with column pivoting of a
 matrix A:  A*P = Q*R  using Level 3 BLAS.
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 trapezoidal 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(J).ne.0, the J-th column of A is permuted
          to the front of A*P (a leading column); if JPVT(J)=0,
          the J-th column of A is a free column.
          On exit, if JPVT(J)=K, then the J-th column of A*P was the
          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 (MAX(1,LWORK))
          On exit, if INFO=0, WORK(1) returns the optimal LWORK.
[in]LWORK
          LWORK is INTEGER
          The dimension of the array WORK. LWORK >= N+1.
          For optimal performance LWORK >= ( N+1 )*NB, where NB
          is the optimal blocksize.

          If LWORK = -1, then a workspace query is assumed; the routine
          only calculates the optimal size of the WORK array, returns
          this value as the first entry of the WORK array, and no error
          message related to LWORK is issued by XERBLA.
[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(k), where k = min(m,n).

  Each H(i) has the form

     H(i) = I - tau * v * v**H

  where tau is a complex scalar, and v is a real/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), and tau in TAU(i).
Contributors:
G. Quintana-Orti, Depto. de Informatica, Universidad Jaime I, Spain X. Sun, Computer Science Dept., Duke University, USA

Definition at line 157 of file cgeqp3.f.

159*
160* -- LAPACK computational routine --
161* -- LAPACK is a software package provided by Univ. of Tennessee, --
162* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
163*
164* .. Scalar Arguments ..
165 INTEGER INFO, LDA, LWORK, M, N
166* ..
167* .. Array Arguments ..
168 INTEGER JPVT( * )
169 REAL RWORK( * )
170 COMPLEX A( LDA, * ), TAU( * ), WORK( * )
171* ..
172*
173* =====================================================================
174*
175* .. Parameters ..
176 INTEGER INB, INBMIN, IXOVER
177 parameter( inb = 1, inbmin = 2, ixover = 3 )
178* ..
179* .. Local Scalars ..
180 LOGICAL LQUERY
181 INTEGER FJB, IWS, J, JB, LWKOPT, MINMN, MINWS, NA, NB,
182 $ NBMIN, NFXD, NX, SM, SMINMN, SN, TOPBMN
183* ..
184* .. External Subroutines ..
185 EXTERNAL cgeqrf, claqp2, claqps, cswap, cunmqr, xerbla
186* ..
187* .. External Functions ..
188 INTEGER ILAENV
189 REAL SCNRM2
190 EXTERNAL ilaenv, scnrm2
191* ..
192* .. Intrinsic Functions ..
193 INTRINSIC int, max, min
194* ..
195* .. Executable Statements ..
196*
197* Test input arguments
198* ====================
199*
200 info = 0
201 lquery = ( lwork.EQ.-1 )
202 IF( m.LT.0 ) THEN
203 info = -1
204 ELSE IF( n.LT.0 ) THEN
205 info = -2
206 ELSE IF( lda.LT.max( 1, m ) ) THEN
207 info = -4
208 END IF
209*
210 IF( info.EQ.0 ) THEN
211 minmn = min( m, n )
212 IF( minmn.EQ.0 ) THEN
213 iws = 1
214 lwkopt = 1
215 ELSE
216 iws = n + 1
217 nb = ilaenv( inb, 'CGEQRF', ' ', m, n, -1, -1 )
218 lwkopt = ( n + 1 )*nb
219 END IF
220 work( 1 ) = cmplx( lwkopt )
221*
222 IF( ( lwork.LT.iws ) .AND. .NOT.lquery ) THEN
223 info = -8
224 END IF
225 END IF
226*
227 IF( info.NE.0 ) THEN
228 CALL xerbla( 'CGEQP3', -info )
229 RETURN
230 ELSE IF( lquery ) THEN
231 RETURN
232 END IF
233*
234* Move initial columns up front.
235*
236 nfxd = 1
237 DO 10 j = 1, n
238 IF( jpvt( j ).NE.0 ) THEN
239 IF( j.NE.nfxd ) THEN
240 CALL cswap( m, a( 1, j ), 1, a( 1, nfxd ), 1 )
241 jpvt( j ) = jpvt( nfxd )
242 jpvt( nfxd ) = j
243 ELSE
244 jpvt( j ) = j
245 END IF
246 nfxd = nfxd + 1
247 ELSE
248 jpvt( j ) = j
249 END IF
250 10 CONTINUE
251 nfxd = nfxd - 1
252*
253* Factorize fixed columns
254* =======================
255*
256* Compute the QR factorization of fixed columns and update
257* remaining columns.
258*
259 IF( nfxd.GT.0 ) THEN
260 na = min( m, nfxd )
261*CC CALL CGEQR2( M, NA, A, LDA, TAU, WORK, INFO )
262 CALL cgeqrf( m, na, a, lda, tau, work, lwork, info )
263 iws = max( iws, int( work( 1 ) ) )
264 IF( na.LT.n ) THEN
265*CC CALL CUNM2R( 'Left', 'Conjugate Transpose', M, N-NA,
266*CC $ NA, A, LDA, TAU, A( 1, NA+1 ), LDA, WORK,
267*CC $ INFO )
268 CALL cunmqr( 'Left', 'Conjugate Transpose', m, n-na, na, a,
269 $ lda, tau, a( 1, na+1 ), lda, work, lwork,
270 $ info )
271 iws = max( iws, int( work( 1 ) ) )
272 END IF
273 END IF
274*
275* Factorize free columns
276* ======================
277*
278 IF( nfxd.LT.minmn ) THEN
279*
280 sm = m - nfxd
281 sn = n - nfxd
282 sminmn = minmn - nfxd
283*
284* Determine the block size.
285*
286 nb = ilaenv( inb, 'CGEQRF', ' ', sm, sn, -1, -1 )
287 nbmin = 2
288 nx = 0
289*
290 IF( ( nb.GT.1 ) .AND. ( nb.LT.sminmn ) ) THEN
291*
292* Determine when to cross over from blocked to unblocked code.
293*
294 nx = max( 0, ilaenv( ixover, 'CGEQRF', ' ', sm, sn, -1,
295 $ -1 ) )
296*
297*
298 IF( nx.LT.sminmn ) THEN
299*
300* Determine if workspace is large enough for blocked code.
301*
302 minws = ( sn+1 )*nb
303 iws = max( iws, minws )
304 IF( lwork.LT.minws ) THEN
305*
306* Not enough workspace to use optimal NB: Reduce NB and
307* determine the minimum value of NB.
308*
309 nb = lwork / ( sn+1 )
310 nbmin = max( 2, ilaenv( inbmin, 'CGEQRF', ' ', sm, sn,
311 $ -1, -1 ) )
312*
313*
314 END IF
315 END IF
316 END IF
317*
318* Initialize partial column norms. The first N elements of work
319* store the exact column norms.
320*
321 DO 20 j = nfxd + 1, n
322 rwork( j ) = scnrm2( sm, a( nfxd+1, j ), 1 )
323 rwork( n+j ) = rwork( j )
324 20 CONTINUE
325*
326 IF( ( nb.GE.nbmin ) .AND. ( nb.LT.sminmn ) .AND.
327 $ ( nx.LT.sminmn ) ) THEN
328*
329* Use blocked code initially.
330*
331 j = nfxd + 1
332*
333* Compute factorization: while loop.
334*
335*
336 topbmn = minmn - nx
337 30 CONTINUE
338 IF( j.LE.topbmn ) THEN
339 jb = min( nb, topbmn-j+1 )
340*
341* Factorize JB columns among columns J:N.
342*
343 CALL claqps( m, n-j+1, j-1, jb, fjb, a( 1, j ), lda,
344 $ jpvt( j ), tau( j ), rwork( j ),
345 $ rwork( n+j ), work( 1 ), work( jb+1 ),
346 $ n-j+1 )
347*
348 j = j + fjb
349 GO TO 30
350 END IF
351 ELSE
352 j = nfxd + 1
353 END IF
354*
355* Use unblocked code to factor the last or only block.
356*
357*
358 IF( j.LE.minmn )
359 $ CALL claqp2( m, n-j+1, j-1, a( 1, j ), lda, jpvt( j ),
360 $ tau( j ), rwork( j ), rwork( n+j ), work( 1 ) )
361*
362 END IF
363*
364 work( 1 ) = cmplx( lwkopt )
365 RETURN
366*
367* End of CGEQP3
368*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine cgeqrf(m, n, a, lda, tau, work, lwork, info)
CGEQRF
Definition cgeqrf.f:146
integer function ilaenv(ispec, name, opts, n1, n2, n3, n4)
ILAENV
Definition ilaenv.f:162
subroutine claqp2(m, n, offset, a, lda, jpvt, tau, vn1, vn2, work)
CLAQP2 computes a QR factorization with column pivoting of the matrix block.
Definition claqp2.f:149
subroutine claqps(m, n, offset, nb, kb, a, lda, jpvt, tau, vn1, vn2, auxv, f, ldf)
CLAQPS computes a step of QR factorization with column pivoting of a real m-by-n matrix A by using BL...
Definition claqps.f:178
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 cunmqr(side, trans, m, n, k, a, lda, tau, c, ldc, work, lwork, info)
CUNMQR
Definition cunmqr.f:168
Here is the call graph for this function:
Here is the caller graph for this function: