LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
zungtsqr.f
Go to the documentation of this file.
1*> \brief \b ZUNGTSQR
2*
3* =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6* http://www.netlib.org/lapack/explore-html/
7*
8*> \htmlonly
9*> Download ZUNGTSQR + dependencies
10*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/zuntsqr.f">
11*> [TGZ]</a>
12*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/zungtsqr.f">
13*> [ZIP]</a>
14*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/zungtsqr.f">
15*> [TXT]</a>
16*> \endhtmlonly
17*
18* Definition:
19* ===========
20*
21* SUBROUTINE ZUNGTSQR( M, N, MB, NB, A, LDA, T, LDT, WORK, LWORK,
22* $ INFO )
23*
24* .. Scalar Arguments ..
25* INTEGER INFO, LDA, LDT, LWORK, M, N, MB, NB
26* ..
27* .. Array Arguments ..
28* COMPLEX*16 A( LDA, * ), T( LDT, * ), WORK( * )
29* ..
30*
31*> \par Purpose:
32* =============
33*>
34*> \verbatim
35*>
36*> ZUNGTSQR generates an M-by-N complex matrix Q_out with orthonormal
37*> columns, which are the first N columns of a product of comlpex unitary
38*> matrices of order M which are returned by ZLATSQR
39*>
40*> Q_out = first_N_columns_of( Q(1)_in * Q(2)_in * ... * Q(k)_in ).
41*>
42*> See the documentation for ZLATSQR.
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. M >= N >= 0.
58*> \endverbatim
59*>
60*> \param[in] MB
61*> \verbatim
62*> MB is INTEGER
63*> The row block size used by ZLATSQR to return
64*> arrays A and T. MB > N.
65*> (Note that if MB > M, then M is used instead of MB
66*> as the row block size).
67*> \endverbatim
68*>
69*> \param[in] NB
70*> \verbatim
71*> NB is INTEGER
72*> The column block size used by ZLATSQR to return
73*> arrays A and T. NB >= 1.
74*> (Note that if NB > N, then N is used instead of NB
75*> as the column block size).
76*> \endverbatim
77*>
78*> \param[in,out] A
79*> \verbatim
80*> A is COMPLEX*16 array, dimension (LDA,N)
81*>
82*> On entry:
83*>
84*> The elements on and above the diagonal are not accessed.
85*> The elements below the diagonal represent the unit
86*> lower-trapezoidal blocked matrix V computed by ZLATSQR
87*> that defines the input matrices Q_in(k) (ones on the
88*> diagonal are not stored) (same format as the output A
89*> below the diagonal in ZLATSQR).
90*>
91*> On exit:
92*>
93*> The array A contains an M-by-N orthonormal matrix Q_out,
94*> i.e the columns of A are orthogonal unit vectors.
95*> \endverbatim
96*>
97*> \param[in] LDA
98*> \verbatim
99*> LDA is INTEGER
100*> The leading dimension of the array A. LDA >= max(1,M).
101*> \endverbatim
102*>
103*> \param[in] T
104*> \verbatim
105*> T is COMPLEX*16 array,
106*> dimension (LDT, N * NIRB)
107*> where NIRB = Number_of_input_row_blocks
108*> = MAX( 1, CEIL((M-N)/(MB-N)) )
109*> Let NICB = Number_of_input_col_blocks
110*> = CEIL(N/NB)
111*>
112*> The upper-triangular block reflectors used to define the
113*> input matrices Q_in(k), k=(1:NIRB*NICB). The block
114*> reflectors are stored in compact form in NIRB block
115*> reflector sequences. Each of NIRB block reflector sequences
116*> is stored in a larger NB-by-N column block of T and consists
117*> of NICB smaller NB-by-NB upper-triangular column blocks.
118*> (same format as the output T in ZLATSQR).
119*> \endverbatim
120*>
121*> \param[in] LDT
122*> \verbatim
123*> LDT is INTEGER
124*> The leading dimension of the array T.
125*> LDT >= max(1,min(NB1,N)).
126*> \endverbatim
127*>
128*> \param[out] WORK
129*> \verbatim
130*> (workspace) COMPLEX*16 array, dimension (MAX(2,LWORK))
131*> On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
132*> \endverbatim
133*>
134*> \param[in] LWORK
135*> \verbatim
136*> LWORK is INTEGER
137*> The dimension of the array WORK. LWORK >= (M+NB)*N.
138*> If LWORK = -1, then a workspace query is assumed.
139*> The routine only calculates the optimal size of the WORK
140*> array, returns this value as the first entry of the WORK
141*> array, and no error message related to LWORK is issued
142*> by XERBLA.
143*> \endverbatim
144*>
145*> \param[out] INFO
146*> \verbatim
147*> INFO is INTEGER
148*> = 0: successful exit
149*> < 0: if INFO = -i, the i-th argument had an illegal value
150*> \endverbatim
151*>
152* Authors:
153* ========
154*
155*> \author Univ. of Tennessee
156*> \author Univ. of California Berkeley
157*> \author Univ. of Colorado Denver
158*> \author NAG Ltd.
159*
160*> \ingroup ungtsqr
161*
162*> \par Contributors:
163* ==================
164*>
165*> \verbatim
166*>
167*> November 2019, Igor Kozachenko,
168*> Computer Science Division,
169*> University of California, Berkeley
170*>
171*> \endverbatim
172*
173* =====================================================================
174 SUBROUTINE zungtsqr( M, N, MB, NB, A, LDA, T, LDT, WORK, LWORK,
175 $ INFO )
176 IMPLICIT NONE
177*
178* -- LAPACK computational routine --
179* -- LAPACK is a software package provided by Univ. of Tennessee, --
180* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
181*
182* .. Scalar Arguments ..
183 INTEGER INFO, LDA, LDT, LWORK, M, N, MB, NB
184* ..
185* .. Array Arguments ..
186 COMPLEX*16 A( LDA, * ), T( LDT, * ), WORK( * )
187* ..
188*
189* =====================================================================
190*
191* .. Parameters ..
192 COMPLEX*16 CONE, CZERO
193 parameter( cone = ( 1.0d+0, 0.0d+0 ),
194 $ czero = ( 0.0d+0, 0.0d+0 ) )
195* ..
196* .. Local Scalars ..
197 LOGICAL LQUERY
198 INTEGER IINFO, LDC, LWORKOPT, LC, LW, NBLOCAL, J
199* ..
200* .. External Subroutines ..
201 EXTERNAL zcopy, zlamtsqr, zlaset, xerbla
202* ..
203* .. Intrinsic Functions ..
204 INTRINSIC dcmplx, max, min
205* ..
206* .. Executable Statements ..
207*
208* Test the input parameters
209*
210 lquery = lwork.EQ.-1
211 info = 0
212 IF( m.LT.0 ) THEN
213 info = -1
214 ELSE IF( n.LT.0 .OR. m.LT.n ) THEN
215 info = -2
216 ELSE IF( mb.LE.n ) THEN
217 info = -3
218 ELSE IF( nb.LT.1 ) THEN
219 info = -4
220 ELSE IF( lda.LT.max( 1, m ) ) THEN
221 info = -6
222 ELSE IF( ldt.LT.max( 1, min( nb, n ) ) ) THEN
223 info = -8
224 ELSE
225*
226* Test the input LWORK for the dimension of the array WORK.
227* This workspace is used to store array C(LDC, N) and WORK(LWORK)
228* in the call to ZLAMTSQR. See the documentation for ZLAMTSQR.
229*
230 IF( lwork.LT.2 .AND. (.NOT.lquery) ) THEN
231 info = -10
232 ELSE
233*
234* Set block size for column blocks
235*
236 nblocal = min( nb, n )
237*
238* LWORK = -1, then set the size for the array C(LDC,N)
239* in ZLAMTSQR call and set the optimal size of the work array
240* WORK(LWORK) in ZLAMTSQR call.
241*
242 ldc = m
243 lc = ldc*n
244 lw = n * nblocal
245*
246 lworkopt = lc+lw
247*
248 IF( ( lwork.LT.max( 1, lworkopt ) ).AND.(.NOT.lquery) ) THEN
249 info = -10
250 END IF
251 END IF
252*
253 END IF
254*
255* Handle error in the input parameters and return workspace query.
256*
257 IF( info.NE.0 ) THEN
258 CALL xerbla( 'ZUNGTSQR', -info )
259 RETURN
260 ELSE IF ( lquery ) THEN
261 work( 1 ) = dcmplx( lworkopt )
262 RETURN
263 END IF
264*
265* Quick return if possible
266*
267 IF( min( m, n ).EQ.0 ) THEN
268 work( 1 ) = dcmplx( lworkopt )
269 RETURN
270 END IF
271*
272* (1) Form explicitly the tall-skinny M-by-N left submatrix Q1_in
273* of M-by-M orthogonal matrix Q_in, which is implicitly stored in
274* the subdiagonal part of input array A and in the input array T.
275* Perform by the following operation using the routine ZLAMTSQR.
276*
277* Q1_in = Q_in * ( I ), where I is a N-by-N identity matrix,
278* ( 0 ) 0 is a (M-N)-by-N zero matrix.
279*
280* (1a) Form M-by-N matrix in the array WORK(1:LDC*N) with ones
281* on the diagonal and zeros elsewhere.
282*
283 CALL zlaset( 'F', m, n, czero, cone, work, ldc )
284*
285* (1b) On input, WORK(1:LDC*N) stores ( I );
286* ( 0 )
287*
288* On output, WORK(1:LDC*N) stores Q1_in.
289*
290 CALL zlamtsqr( 'L', 'N', m, n, n, mb, nblocal, a, lda, t, ldt,
291 $ work, ldc, work( lc+1 ), lw, iinfo )
292*
293* (2) Copy the result from the part of the work array (1:M,1:N)
294* with the leading dimension LDC that starts at WORK(1) into
295* the output array A(1:M,1:N) column-by-column.
296*
297 DO j = 1, n
298 CALL zcopy( m, work( (j-1)*ldc + 1 ), 1, a( 1, j ), 1 )
299 END DO
300*
301 work( 1 ) = dcmplx( lworkopt )
302 RETURN
303*
304* End of ZUNGTSQR
305*
306 END
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine zcopy(n, zx, incx, zy, incy)
ZCOPY
Definition zcopy.f:81
subroutine zlamtsqr(side, trans, m, n, k, mb, nb, a, lda, t, ldt, c, ldc, work, lwork, info)
ZLAMTSQR
Definition zlamtsqr.f:199
subroutine zlaset(uplo, m, n, alpha, beta, a, lda)
ZLASET initializes the off-diagonal elements and the diagonal elements of a matrix to given values.
Definition zlaset.f:106
subroutine zungtsqr(m, n, mb, nb, a, lda, t, ldt, work, lwork, info)
ZUNGTSQR
Definition zungtsqr.f:176