LAPACK  3.9.0
LAPACK: Linear Algebra PACKage

◆ zungtsqr()

subroutine zungtsqr ( integer  M,
integer  N,
integer  MB,
integer  NB,
complex*16, dimension( lda, * )  A,
integer  LDA,
complex*16, dimension( ldt, * )  T,
integer  LDT,
complex*16, dimension( * )  WORK,
integer  LWORK,
integer  INFO 
)

ZUNGTSQR

Download ZUNGTSQR + dependencies [TGZ] [ZIP] [TXT] \par Purpose: @verbatim ZUNGTSQR 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 ZLATSQR Q_out = first_N_columns_of( Q(1)_in * Q(2)_in * ... * Q(k)_in ). See the documentation for ZLATSQR. \endverbatim @param [in] M @verbatim M is INTEGER The number of rows of the matrix A. M >= 0. \endverbatim @param [in] N @verbatim N is INTEGER The number of columns of the matrix A. M >= N >= 0. \endverbatim @param [in] MB @verbatim MB is INTEGER The row block size used by DLATSQR 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). \endverbatim @param [in] NB @verbatim NB is INTEGER The column block size used by ZLATSQR 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). \endverbatim @param [in,out] A @verbatim A is COMPLEX*16 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 ZLATSQR 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 ZLATSQR). On exit: The array A contains an M-by-N orthonormal matrix Q_out, i.e the columns of A are orthogonal unit vectors. \endverbatim @param [in] LDA @verbatim LDA is INTEGER The leading dimension of the array A. LDA >= max(1,M). \endverbatim @param [in] T @verbatim T is COMPLEX*16 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 ZLATSQR). \endverbatim @param [in] LDT @verbatim LDT is INTEGER The leading dimension of the array T. LDT >= max(1,min(NB1,N)). \endverbatim @param [out] WORK @verbatim (workspace) COMPLEX*16 array, dimension (MAX(2,LWORK)) On exit, if INFO = 0, WORK(1) returns the optimal LWORK. \endverbatim @param [in] LWORK @verbatim 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. \endverbatim @param [out] INFO @verbatim INFO is INTEGER = 0: successful exit < 0: if INFO = -i, the i-th argument had an illegal value \endverbatim \author Univ. of Tennessee \author Univ. of California Berkeley \author Univ. of Colorado Denver \author NAG Ltd. \date November 2019 \par Contributors: @verbatim November 2019, Igor Kozachenko, Computer Science Division, University of California, Berkeley \endverbatim

Definition at line 176 of file zungtsqr.f.

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