LAPACK 3.12.1
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
sgeqp3.f
Go to the documentation of this file.
1*> \brief \b SGEQP3
2*
3* =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6* http://www.netlib.org/lapack/explore-html/
7*
8*> Download SGEQP3 + dependencies
9*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/sgeqp3.f">
10*> [TGZ]</a>
11*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/sgeqp3.f">
12*> [ZIP]</a>
13*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/sgeqp3.f">
14*> [TXT]</a>
15*
16* Definition:
17* ===========
18*
19* SUBROUTINE SGEQP3( M, N, A, LDA, JPVT, TAU, WORK, LWORK, INFO )
20*
21* .. Scalar Arguments ..
22* INTEGER INFO, LDA, LWORK, M, N
23* ..
24* .. Array Arguments ..
25* INTEGER JPVT( * )
26* REAL A( LDA, * ), TAU( * ), WORK( * )
27* ..
28*
29*
30*> \par Purpose:
31* =============
32*>
33*> \verbatim
34*>
35*> SGEQP3 computes a QR factorization with column pivoting of a
36*> matrix A: A*P = Q*R using Level 3 BLAS.
37*> \endverbatim
38*
39* Arguments:
40* ==========
41*
42*> \param[in] M
43*> \verbatim
44*> M is INTEGER
45*> The number of rows of the matrix A. M >= 0.
46*> \endverbatim
47*>
48*> \param[in] N
49*> \verbatim
50*> N is INTEGER
51*> The number of columns of the matrix A. N >= 0.
52*> \endverbatim
53*>
54*> \param[in,out] A
55*> \verbatim
56*> A is REAL array, dimension (LDA,N)
57*> On entry, the M-by-N matrix A.
58*> On exit, the upper triangle of the array contains the
59*> min(M,N)-by-N upper trapezoidal matrix R; the elements below
60*> the diagonal, together with the array TAU, represent the
61*> orthogonal matrix Q as a product of min(M,N) elementary
62*> reflectors.
63*> \endverbatim
64*>
65*> \param[in] LDA
66*> \verbatim
67*> LDA is INTEGER
68*> The leading dimension of the array A. LDA >= max(1,M).
69*> \endverbatim
70*>
71*> \param[in,out] JPVT
72*> \verbatim
73*> JPVT is INTEGER array, dimension (N)
74*> On entry, if JPVT(J).ne.0, the J-th column of A is permuted
75*> to the front of A*P (a leading column); if JPVT(J)=0,
76*> the J-th column of A is a free column.
77*> On exit, if JPVT(J)=K, then the J-th column of A*P was the
78*> the K-th column of A.
79*> \endverbatim
80*>
81*> \param[out] TAU
82*> \verbatim
83*> TAU is REAL array, dimension (min(M,N))
84*> The scalar factors of the elementary reflectors.
85*> \endverbatim
86*>
87*> \param[out] WORK
88*> \verbatim
89*> WORK is REAL array, dimension (MAX(1,LWORK))
90*> On exit, if INFO=0, WORK(1) returns the optimal LWORK.
91*> \endverbatim
92*>
93*> \param[in] LWORK
94*> \verbatim
95*> LWORK is INTEGER
96*> The dimension of the array WORK. LWORK >= 3*N+1.
97*> For optimal performance LWORK >= 2*N+( N+1 )*NB, where NB
98*> is the optimal blocksize.
99*>
100*> If LWORK = -1, then a workspace query is assumed; the routine
101*> only calculates the optimal size of the WORK array, returns
102*> this value as the first entry of the WORK array, and no error
103*> message related to LWORK is issued by XERBLA.
104*> \endverbatim
105*>
106*> \param[out] INFO
107*> \verbatim
108*> INFO is INTEGER
109*> = 0: successful exit.
110*> < 0: if INFO = -i, the i-th argument had an illegal value.
111*> \endverbatim
112*
113* Authors:
114* ========
115*
116*> \author Univ. of Tennessee
117*> \author Univ. of California Berkeley
118*> \author Univ. of Colorado Denver
119*> \author NAG Ltd.
120*
121*> \ingroup geqp3
122*
123*> \par Further Details:
124* =====================
125*>
126*> \verbatim
127*>
128*> The matrix Q is represented as a product of elementary reflectors
129*>
130*> Q = H(1) H(2) . . . H(k), where k = min(m,n).
131*>
132*> Each H(i) has the form
133*>
134*> H(i) = I - tau * v * v**T
135*>
136*> where tau is a real scalar, and v is a real/complex vector
137*> with v(1:i-1) = 0 and v(i) = 1; v(i+1:m) is stored on exit in
138*> A(i+1:m,i), and tau in TAU(i).
139*> \endverbatim
140*
141*> \par Contributors:
142* ==================
143*>
144*> G. Quintana-Orti, Depto. de Informatica, Universidad Jaime I, Spain
145*> X. Sun, Computer Science Dept., Duke University, USA
146*>
147* =====================================================================
148 SUBROUTINE sgeqp3( M, N, A, LDA, JPVT, TAU, WORK, LWORK, INFO )
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 REAL 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 sgeqrf, slaqp2, slaqps, sormqr, sswap,
175 $ xerbla
176* ..
177* .. External Functions ..
178 INTEGER ILAENV
179 REAL SNRM2, SROUNDUP_LWORK
180 EXTERNAL ilaenv, snrm2, sroundup_lwork
181* ..
182* .. Intrinsic Functions ..
183 INTRINSIC int, max, min
184* Test input arguments
185* ====================
186*
187 info = 0
188 lquery = ( lwork.EQ.-1 )
189 IF( m.LT.0 ) THEN
190 info = -1
191 ELSE IF( n.LT.0 ) THEN
192 info = -2
193 ELSE IF( lda.LT.max( 1, m ) ) THEN
194 info = -4
195 END IF
196*
197 IF( info.EQ.0 ) THEN
198 minmn = min( m, n )
199 IF( minmn.EQ.0 ) THEN
200 iws = 1
201 lwkopt = 1
202 ELSE
203 iws = 3*n + 1
204 nb = ilaenv( inb, 'SGEQRF', ' ', m, n, -1, -1 )
205 lwkopt = 2*n + ( n + 1 )*nb
206 END IF
207 work( 1 ) = sroundup_lwork(lwkopt)
208*
209 IF( ( lwork.LT.iws ) .AND. .NOT.lquery ) THEN
210 info = -8
211 END IF
212 END IF
213*
214 IF( info.NE.0 ) THEN
215 CALL xerbla( 'SGEQP3', -info )
216 RETURN
217 ELSE IF( lquery ) THEN
218 RETURN
219 END IF
220*
221* Move initial columns up front.
222*
223 nfxd = 1
224 DO 10 j = 1, n
225 IF( jpvt( j ).NE.0 ) THEN
226 IF( j.NE.nfxd ) THEN
227 CALL sswap( m, a( 1, j ), 1, a( 1, nfxd ), 1 )
228 jpvt( j ) = jpvt( nfxd )
229 jpvt( nfxd ) = j
230 ELSE
231 jpvt( j ) = j
232 END IF
233 nfxd = nfxd + 1
234 ELSE
235 jpvt( j ) = j
236 END IF
237 10 CONTINUE
238 nfxd = nfxd - 1
239*
240* Factorize fixed columns
241* =======================
242*
243* Compute the QR factorization of fixed columns and update
244* remaining columns.
245*
246 IF( nfxd.GT.0 ) THEN
247 na = min( m, nfxd )
248*CC CALL SGEQR2( M, NA, A, LDA, TAU, WORK, INFO )
249 CALL sgeqrf( m, na, a, lda, tau, work, lwork, info )
250 iws = max( iws, int( work( 1 ) ) )
251 IF( na.LT.n ) THEN
252*CC CALL SORM2R( 'Left', 'Transpose', M, N-NA, NA, A, LDA,
253*CC $ TAU, A( 1, NA+1 ), LDA, WORK, INFO )
254 CALL sormqr( 'Left', 'Transpose', m, n-na, na, a, lda,
255 $ tau,
256 $ a( 1, na+1 ), lda, work, lwork, info )
257 iws = max( iws, int( work( 1 ) ) )
258 END IF
259 END IF
260*
261* Factorize free columns
262* ======================
263*
264 IF( nfxd.LT.minmn ) THEN
265*
266 sm = m - nfxd
267 sn = n - nfxd
268 sminmn = minmn - nfxd
269*
270* Determine the block size.
271*
272 nb = ilaenv( inb, 'SGEQRF', ' ', sm, sn, -1, -1 )
273 nbmin = 2
274 nx = 0
275*
276 IF( ( nb.GT.1 ) .AND. ( nb.LT.sminmn ) ) THEN
277*
278* Determine when to cross over from blocked to unblocked code.
279*
280 nx = max( 0, ilaenv( ixover, 'SGEQRF', ' ', sm, sn, -1,
281 $ -1 ) )
282*
283*
284 IF( nx.LT.sminmn ) THEN
285*
286* Determine if workspace is large enough for blocked code.
287*
288 minws = 2*sn + ( sn+1 )*nb
289 iws = max( iws, minws )
290 IF( lwork.LT.minws ) THEN
291*
292* Not enough workspace to use optimal NB: Reduce NB and
293* determine the minimum value of NB.
294*
295 nb = ( lwork-2*sn ) / ( sn+1 )
296 nbmin = max( 2, ilaenv( inbmin, 'SGEQRF', ' ', sm,
297 $ sn,
298 $ -1, -1 ) )
299*
300*
301 END IF
302 END IF
303 END IF
304*
305* Initialize partial column norms. The first N elements of work
306* store the exact column norms.
307*
308 DO 20 j = nfxd + 1, n
309 work( j ) = snrm2( sm, a( nfxd+1, j ), 1 )
310 work( n+j ) = work( j )
311 20 CONTINUE
312*
313 IF( ( nb.GE.nbmin ) .AND. ( nb.LT.sminmn ) .AND.
314 $ ( nx.LT.sminmn ) ) THEN
315*
316* Use blocked code initially.
317*
318 j = nfxd + 1
319*
320* Compute factorization: while loop.
321*
322*
323 topbmn = minmn - nx
324 30 CONTINUE
325 IF( j.LE.topbmn ) THEN
326 jb = min( nb, topbmn-j+1 )
327*
328* Factorize JB columns among columns J:N.
329*
330 CALL slaqps( m, n-j+1, j-1, jb, fjb, a( 1, j ), lda,
331 $ jpvt( j ), tau( j ), work( j ), work( n+j ),
332 $ work( 2*n+1 ), work( 2*n+jb+1 ), n-j+1 )
333*
334 j = j + fjb
335 GO TO 30
336 END IF
337 ELSE
338 j = nfxd + 1
339 END IF
340*
341* Use unblocked code to factor the last or only block.
342*
343*
344 IF( j.LE.minmn )
345 $ CALL slaqp2( m, n-j+1, j-1, a( 1, j ), lda, jpvt( j ),
346 $ tau( j ), work( j ), work( n+j ),
347 $ work( 2*n+1 ) )
348*
349 END IF
350*
351 work( 1 ) = sroundup_lwork(iws)
352 RETURN
353*
354* End of SGEQP3
355*
356 END
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine sgeqp3(m, n, a, lda, jpvt, tau, work, lwork, info)
SGEQP3
Definition sgeqp3.f:149
subroutine sgeqrf(m, n, a, lda, tau, work, lwork, info)
SGEQRF
Definition sgeqrf.f:144
subroutine slaqp2(m, n, offset, a, lda, jpvt, tau, vn1, vn2, work)
SLAQP2 computes a QR factorization with column pivoting of the matrix block.
Definition slaqp2.f:145
subroutine slaqps(m, n, offset, nb, kb, a, lda, jpvt, tau, vn1, vn2, auxv, f, ldf)
SLAQPS computes a step of QR factorization with column pivoting of a real m-by-n matrix A by using BL...
Definition slaqps.f:175
subroutine sswap(n, sx, incx, sy, incy)
SSWAP
Definition sswap.f:82
subroutine sormqr(side, trans, m, n, k, a, lda, tau, c, ldc, work, lwork, info)
SORMQR
Definition sormqr.f:166