LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
stzrzf.f
Go to the documentation of this file.
1*> \brief \b STZRZF
2*
3* =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6* http://www.netlib.org/lapack/explore-html/
7*
8*> \htmlonly
9*> Download STZRZF + dependencies
10*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/stzrzf.f">
11*> [TGZ]</a>
12*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/stzrzf.f">
13*> [ZIP]</a>
14*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/stzrzf.f">
15*> [TXT]</a>
16*> \endhtmlonly
17*
18* Definition:
19* ===========
20*
21* SUBROUTINE STZRZF( M, N, A, LDA, TAU, WORK, LWORK, INFO )
22*
23* .. Scalar Arguments ..
24* INTEGER INFO, LDA, LWORK, M, N
25* ..
26* .. Array Arguments ..
27* REAL A( LDA, * ), TAU( * ), WORK( * )
28* ..
29*
30*
31*> \par Purpose:
32* =============
33*>
34*> \verbatim
35*>
36*> STZRZF reduces the M-by-N ( M<=N ) real upper trapezoidal matrix A
37*> to upper triangular form by means of orthogonal transformations.
38*>
39*> The upper trapezoidal matrix A is factored as
40*>
41*> A = ( R 0 ) * Z,
42*>
43*> where Z is an N-by-N orthogonal matrix and R is an M-by-M upper
44*> triangular matrix.
45*> \endverbatim
46*
47* Arguments:
48* ==========
49*
50*> \param[in] M
51*> \verbatim
52*> M is INTEGER
53*> The number of rows of the matrix A. M >= 0.
54*> \endverbatim
55*>
56*> \param[in] N
57*> \verbatim
58*> N is INTEGER
59*> The number of columns of the matrix A. N >= M.
60*> \endverbatim
61*>
62*> \param[in,out] A
63*> \verbatim
64*> A is REAL array, dimension (LDA,N)
65*> On entry, the leading M-by-N upper trapezoidal part of the
66*> array A must contain the matrix to be factorized.
67*> On exit, the leading M-by-M upper triangular part of A
68*> contains the upper triangular matrix R, and elements M+1 to
69*> N of the first M rows of A, with the array TAU, represent the
70*> orthogonal matrix Z as a product of M elementary reflectors.
71*> \endverbatim
72*>
73*> \param[in] LDA
74*> \verbatim
75*> LDA is INTEGER
76*> The leading dimension of the array A. LDA >= max(1,M).
77*> \endverbatim
78*>
79*> \param[out] TAU
80*> \verbatim
81*> TAU is REAL array, dimension (M)
82*> The scalar factors of the elementary reflectors.
83*> \endverbatim
84*>
85*> \param[out] WORK
86*> \verbatim
87*> WORK is REAL array, dimension (MAX(1,LWORK))
88*> On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
89*> \endverbatim
90*>
91*> \param[in] LWORK
92*> \verbatim
93*> LWORK is INTEGER
94*> The dimension of the array WORK. LWORK >= max(1,M).
95*> For optimum performance LWORK >= M*NB, where NB is
96*> the optimal blocksize.
97*>
98*> If LWORK = -1, then a workspace query is assumed; the routine
99*> only calculates the optimal size of the WORK array, returns
100*> this value as the first entry of the WORK array, and no error
101*> message related to LWORK is issued by XERBLA.
102*> \endverbatim
103*>
104*> \param[out] INFO
105*> \verbatim
106*> INFO is INTEGER
107*> = 0: successful exit
108*> < 0: if INFO = -i, the i-th argument had an illegal value
109*> \endverbatim
110*
111* Authors:
112* ========
113*
114*> \author Univ. of Tennessee
115*> \author Univ. of California Berkeley
116*> \author Univ. of Colorado Denver
117*> \author NAG Ltd.
118*
119*> \ingroup tzrzf
120*
121*> \par Contributors:
122* ==================
123*>
124*> A. Petitet, Computer Science Dept., Univ. of Tenn., Knoxville, USA
125*
126*> \par Further Details:
127* =====================
128*>
129*> \verbatim
130*>
131*> The N-by-N matrix Z can be computed by
132*>
133*> Z = Z(1)*Z(2)* ... *Z(M)
134*>
135*> where each N-by-N Z(k) is given by
136*>
137*> Z(k) = I - tau(k)*v(k)*v(k)**T
138*>
139*> with v(k) is the kth row vector of the M-by-N matrix
140*>
141*> V = ( I A(:,M+1:N) )
142*>
143*> I is the M-by-M identity matrix, A(:,M+1:N)
144*> is the output stored in A on exit from STZRZF,
145*> and tau(k) is the kth element of the array TAU.
146*>
147*> \endverbatim
148*>
149* =====================================================================
150 SUBROUTINE stzrzf( M, N, A, LDA, TAU, WORK, LWORK, INFO )
151*
152* -- LAPACK computational routine --
153* -- LAPACK is a software package provided by Univ. of Tennessee, --
154* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
155*
156* .. Scalar Arguments ..
157 INTEGER INFO, LDA, LWORK, M, N
158* ..
159* .. Array Arguments ..
160 REAL A( LDA, * ), TAU( * ), WORK( * )
161* ..
162*
163* =====================================================================
164*
165* .. Parameters ..
166 REAL ZERO
167 parameter( zero = 0.0e+0 )
168* ..
169* .. Local Scalars ..
170 LOGICAL LQUERY
171 INTEGER I, IB, IWS, KI, KK, LDWORK, LWKMIN, LWKOPT,
172 $ M1, MU, NB, NBMIN, NX
173* ..
174* .. External Subroutines ..
175 EXTERNAL xerbla, slarzb, slarzt, slatrz
176* ..
177* .. Intrinsic Functions ..
178 INTRINSIC max, min
179* ..
180* .. External Functions ..
181 INTEGER ILAENV
182 REAL SROUNDUP_LWORK
183 EXTERNAL ilaenv, sroundup_lwork
184* ..
185* .. Executable Statements ..
186*
187* Test the input arguments
188*
189 info = 0
190 lquery = ( lwork.EQ.-1 )
191 IF( m.LT.0 ) THEN
192 info = -1
193 ELSE IF( n.LT.m ) THEN
194 info = -2
195 ELSE IF( lda.LT.max( 1, m ) ) THEN
196 info = -4
197 END IF
198*
199 IF( info.EQ.0 ) THEN
200 IF( m.EQ.0 .OR. m.EQ.n ) THEN
201 lwkopt = 1
202 lwkmin = 1
203 ELSE
204*
205* Determine the block size.
206*
207 nb = ilaenv( 1, 'SGERQF', ' ', m, n, -1, -1 )
208 lwkopt = m*nb
209 lwkmin = max( 1, m )
210 END IF
211 work( 1 ) = sroundup_lwork(lwkopt)
212*
213 IF( lwork.LT.lwkmin .AND. .NOT.lquery ) THEN
214 info = -7
215 END IF
216 END IF
217*
218 IF( info.NE.0 ) THEN
219 CALL xerbla( 'STZRZF', -info )
220 RETURN
221 ELSE IF( lquery ) THEN
222 RETURN
223 END IF
224*
225* Quick return if possible
226*
227 IF( m.EQ.0 ) THEN
228 RETURN
229 ELSE IF( m.EQ.n ) THEN
230 DO 10 i = 1, n
231 tau( i ) = zero
232 10 CONTINUE
233 RETURN
234 END IF
235*
236 nbmin = 2
237 nx = 1
238 iws = m
239 IF( nb.GT.1 .AND. nb.LT.m ) THEN
240*
241* Determine when to cross over from blocked to unblocked code.
242*
243 nx = max( 0, ilaenv( 3, 'SGERQF', ' ', m, n, -1, -1 ) )
244 IF( nx.LT.m ) THEN
245*
246* Determine if workspace is large enough for blocked code.
247*
248 ldwork = m
249 iws = ldwork*nb
250 IF( lwork.LT.iws ) THEN
251*
252* Not enough workspace to use optimal NB: reduce NB and
253* determine the minimum value of NB.
254*
255 nb = lwork / ldwork
256 nbmin = max( 2, ilaenv( 2, 'SGERQF', ' ', m, n, -1,
257 $ -1 ) )
258 END IF
259 END IF
260 END IF
261*
262 IF( nb.GE.nbmin .AND. nb.LT.m .AND. nx.LT.m ) THEN
263*
264* Use blocked code initially.
265* The last kk rows are handled by the block method.
266*
267 m1 = min( m+1, n )
268 ki = ( ( m-nx-1 ) / nb )*nb
269 kk = min( m, ki+nb )
270*
271 DO 20 i = m - kk + ki + 1, m - kk + 1, -nb
272 ib = min( m-i+1, nb )
273*
274* Compute the TZ factorization of the current block
275* A(i:i+ib-1,i:n)
276*
277 CALL slatrz( ib, n-i+1, n-m, a( i, i ), lda, tau( i ),
278 $ work )
279 IF( i.GT.1 ) THEN
280*
281* Form the triangular factor of the block reflector
282* H = H(i+ib-1) . . . H(i+1) H(i)
283*
284 CALL slarzt( 'Backward', 'Rowwise', n-m, ib, a( i, m1 ),
285 $ lda, tau( i ), work, ldwork )
286*
287* Apply H to A(1:i-1,i:n) from the right
288*
289 CALL slarzb( 'Right', 'No transpose', 'Backward',
290 $ 'Rowwise', i-1, n-i+1, ib, n-m, a( i, m1 ),
291 $ lda, work, ldwork, a( 1, i ), lda,
292 $ work( ib+1 ), ldwork )
293 END IF
294 20 CONTINUE
295 mu = i + nb - 1
296 ELSE
297 mu = m
298 END IF
299*
300* Use unblocked code to factor the last or only block
301*
302 IF( mu.GT.0 )
303 $ CALL slatrz( mu, n, n-m, a, lda, tau, work )
304*
305 work( 1 ) = sroundup_lwork(lwkopt)
306*
307 RETURN
308*
309* End of STZRZF
310*
311 END
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine slarzb(side, trans, direct, storev, m, n, k, l, v, ldv, t, ldt, c, ldc, work, ldwork)
SLARZB applies a block reflector or its transpose to a general matrix.
Definition slarzb.f:183
subroutine slarzt(direct, storev, n, k, v, ldv, tau, t, ldt)
SLARZT forms the triangular factor T of a block reflector H = I - vtvH.
Definition slarzt.f:185
subroutine slatrz(m, n, l, a, lda, tau, work)
SLATRZ factors an upper trapezoidal matrix by means of orthogonal transformations.
Definition slatrz.f:140
subroutine stzrzf(m, n, a, lda, tau, work, lwork, info)
STZRZF
Definition stzrzf.f:151