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