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

◆ cungtsqr()

subroutine cungtsqr ( integer m,
integer n,
integer mb,
integer nb,
complex, dimension( lda, * ) a,
integer lda,
complex, dimension( ldt, * ) t,
integer ldt,
complex, dimension( * ) work,
integer lwork,
integer info )

CUNGTSQR

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

Purpose:
!>
!> CUNGTSQR generates an M-by-N complex matrix Q_out with orthonormal
!> columns, which are the first N columns of a product of comlpex unitary
!> matrices of order M which are returned by CLATSQR
!>
!>      Q_out = first_N_columns_of( Q(1)_in * Q(2)_in * ... * Q(k)_in ).
!>
!> See the documentation for CLATSQR.
!> 
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. M >= N >= 0.
!> 
[in]MB
!>          MB is INTEGER
!>          The row block size used by CLATSQR to return
!>          arrays A and T. MB > N.
!>          (Note that if MB > M, then M is used instead of MB
!>          as the row block size).
!> 
[in]NB
!>          NB is INTEGER
!>          The column block size used by CLATSQR to return
!>          arrays A and T. NB >= 1.
!>          (Note that if NB > N, then N is used instead of NB
!>          as the column block size).
!> 
[in,out]A
!>          A is COMPLEX array, dimension (LDA,N)
!>
!>          On entry:
!>
!>             The elements on and above the diagonal are not accessed.
!>             The elements below the diagonal represent the unit
!>             lower-trapezoidal blocked matrix V computed by CLATSQR
!>             that defines the input matrices Q_in(k) (ones on the
!>             diagonal are not stored) (same format as the output A
!>             below the diagonal in CLATSQR).
!>
!>          On exit:
!>
!>             The array A contains an M-by-N orthonormal matrix Q_out,
!>             i.e the columns of A are orthogonal unit vectors.
!> 
[in]LDA
!>          LDA is INTEGER
!>          The leading dimension of the array A.  LDA >= max(1,M).
!> 
[in]T
!>          T is COMPLEX array,
!>          dimension (LDT, N * NIRB)
!>          where NIRB = Number_of_input_row_blocks
!>                     = MAX( 1, CEIL((M-N)/(MB-N)) )
!>          Let NICB = Number_of_input_col_blocks
!>                   = CEIL(N/NB)
!>
!>          The upper-triangular block reflectors used to define the
!>          input matrices Q_in(k), k=(1:NIRB*NICB). The block
!>          reflectors are stored in compact form in NIRB block
!>          reflector sequences. Each of NIRB block reflector sequences
!>          is stored in a larger NB-by-N column block of T and consists
!>          of NICB smaller NB-by-NB upper-triangular column blocks.
!>          (same format as the output T in CLATSQR).
!> 
[in]LDT
!>          LDT is INTEGER
!>          The leading dimension of the array T.
!>          LDT >= max(1,min(NB1,N)).
!> 
[out]WORK
!>          (workspace) COMPLEX array, dimension (MAX(2,LWORK))
!>          On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
!> 
[in]LWORK
!>          LWORK is INTEGER
!>          The dimension of the array WORK.  LWORK >= (M+NB)*N.
!>          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.
Contributors:
!>
!> November 2019, Igor Kozachenko,
!>                Computer Science Division,
!>                University of California, Berkeley
!>
!> 

Definition at line 172 of file cungtsqr.f.

174 IMPLICIT NONE
175*
176* -- LAPACK computational routine --
177* -- LAPACK is a software package provided by Univ. of Tennessee, --
178* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
179*
180* .. Scalar Arguments ..
181 INTEGER INFO, LDA, LDT, LWORK, M, N, MB, NB
182* ..
183* .. Array Arguments ..
184 COMPLEX A( LDA, * ), T( LDT, * ), WORK( * )
185* ..
186*
187* =====================================================================
188*
189* .. Parameters ..
190 COMPLEX CONE, CZERO
191 parameter( cone = ( 1.0e+0, 0.0e+0 ),
192 $ czero = ( 0.0e+0, 0.0e+0 ) )
193* ..
194* .. Local Scalars ..
195 LOGICAL LQUERY
196 INTEGER IINFO, LDC, LWORKOPT, LC, LW, NBLOCAL, J
197* ..
198* .. External Subroutines ..
199 EXTERNAL ccopy, clamtsqr, claset, xerbla
200* ..
201* .. Intrinsic Functions ..
202 INTRINSIC cmplx, max, min
203* ..
204* .. Executable Statements ..
205*
206* Test the input parameters
207*
208 lquery = lwork.EQ.-1
209 info = 0
210 IF( m.LT.0 ) THEN
211 info = -1
212 ELSE IF( n.LT.0 .OR. m.LT.n ) THEN
213 info = -2
214 ELSE IF( mb.LE.n ) THEN
215 info = -3
216 ELSE IF( nb.LT.1 ) THEN
217 info = -4
218 ELSE IF( lda.LT.max( 1, m ) ) THEN
219 info = -6
220 ELSE IF( ldt.LT.max( 1, min( nb, n ) ) ) THEN
221 info = -8
222 ELSE
223*
224* Test the input LWORK for the dimension of the array WORK.
225* This workspace is used to store array C(LDC, N) and WORK(LWORK)
226* in the call to CLAMTSQR. See the documentation for CLAMTSQR.
227*
228 IF( lwork.LT.2 .AND. (.NOT.lquery) ) THEN
229 info = -10
230 ELSE
231*
232* Set block size for column blocks
233*
234 nblocal = min( nb, n )
235*
236* LWORK = -1, then set the size for the array C(LDC,N)
237* in CLAMTSQR call and set the optimal size of the work array
238* WORK(LWORK) in CLAMTSQR call.
239*
240 ldc = m
241 lc = ldc*n
242 lw = n * nblocal
243*
244 lworkopt = lc+lw
245*
246 IF( ( lwork.LT.max( 1, lworkopt ) ).AND.(.NOT.lquery) ) THEN
247 info = -10
248 END IF
249 END IF
250*
251 END IF
252*
253* Handle error in the input parameters and return workspace query.
254*
255 IF( info.NE.0 ) THEN
256 CALL xerbla( 'CUNGTSQR', -info )
257 RETURN
258 ELSE IF ( lquery ) THEN
259 work( 1 ) = cmplx( lworkopt )
260 RETURN
261 END IF
262*
263* Quick return if possible
264*
265 IF( min( m, n ).EQ.0 ) THEN
266 work( 1 ) = cmplx( lworkopt )
267 RETURN
268 END IF
269*
270* (1) Form explicitly the tall-skinny M-by-N left submatrix Q1_in
271* of M-by-M orthogonal matrix Q_in, which is implicitly stored in
272* the subdiagonal part of input array A and in the input array T.
273* Perform by the following operation using the routine CLAMTSQR.
274*
275* Q1_in = Q_in * ( I ), where I is a N-by-N identity matrix,
276* ( 0 ) 0 is a (M-N)-by-N zero matrix.
277*
278* (1a) Form M-by-N matrix in the array WORK(1:LDC*N) with ones
279* on the diagonal and zeros elsewhere.
280*
281 CALL claset( 'F', m, n, czero, cone, work, ldc )
282*
283* (1b) On input, WORK(1:LDC*N) stores ( I );
284* ( 0 )
285*
286* On output, WORK(1:LDC*N) stores Q1_in.
287*
288 CALL clamtsqr( 'L', 'N', m, n, n, mb, nblocal, a, lda, t, ldt,
289 $ work, ldc, work( lc+1 ), lw, iinfo )
290*
291* (2) Copy the result from the part of the work array (1:M,1:N)
292* with the leading dimension LDC that starts at WORK(1) into
293* the output array A(1:M,1:N) column-by-column.
294*
295 DO j = 1, n
296 CALL ccopy( m, work( (j-1)*ldc + 1 ), 1, a( 1, j ), 1 )
297 END DO
298*
299 work( 1 ) = cmplx( lworkopt )
300 RETURN
301*
302* End of CUNGTSQR
303*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine ccopy(n, cx, incx, cy, incy)
CCOPY
Definition ccopy.f:81
subroutine clamtsqr(side, trans, m, n, k, mb, nb, a, lda, t, ldt, c, ldc, work, lwork, info)
CLAMTSQR
Definition clamtsqr.f:201
subroutine claset(uplo, m, n, alpha, beta, a, lda)
CLASET initializes the off-diagonal elements and the diagonal elements of a matrix to given values.
Definition claset.f:104
Here is the call graph for this function:
Here is the caller graph for this function: