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

◆ dgeqp3()

subroutine dgeqp3 ( integer m,
integer n,
double precision, dimension( lda, * ) a,
integer lda,
integer, dimension( * ) jpvt,
double precision, dimension( * ) tau,
double precision, dimension( * ) work,
integer lwork,
integer info )

DGEQP3

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

Purpose:
!>
!> DGEQP3 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 DOUBLE PRECISION 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
!>          orthogonal 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 DOUBLE PRECISION array, dimension (min(M,N))
!>          The scalar factors of the elementary reflectors.
!> 
[out]WORK
!>          WORK is DOUBLE PRECISION 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 >= 3*N+1.
!>          For optimal performance LWORK >= 2*N+( 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]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**T
!>
!>  where tau is a real 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 148 of file dgeqp3.f.

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