LAPACK  3.9.0
LAPACK: Linear Algebra PACKage

◆ clamtsqr()

subroutine clamtsqr ( character  SIDE,
character  TRANS,
integer  M,
integer  N,
integer  K,
integer  MB,
integer  NB,
complex, dimension( lda, * )  A,
integer  LDA,
complex, dimension( ldt, * )  T,
integer  LDT,
complex, dimension(ldc, * )  C,
integer  LDC,
complex, dimension( * )  WORK,
integer  LWORK,
integer  INFO 
)

CLAMTSQR

Purpose:
      CLAMTSQR overwrites the general complex M-by-N matrix C with


                 SIDE = 'L'     SIDE = 'R'
 TRANS = 'N':      Q * C          C * Q
 TRANS = 'C':      Q**H * C       C * Q**H
      where Q is a real orthogonal matrix defined as the product
      of blocked elementary reflectors computed by tall skinny
      QR factorization (CLATSQR)
Parameters
[in]SIDE
          SIDE is CHARACTER*1
          = 'L': apply Q or Q**H from the Left;
          = 'R': apply Q or Q**H from the Right.
[in]TRANS
          TRANS is CHARACTER*1
          = 'N':  No transpose, apply Q;
          = 'C':  Conjugate Transpose, apply Q**H.
[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 C. M >= N >= 0.
[in]K
          K is INTEGER
          The number of elementary reflectors whose product defines
          the matrix Q.
          N >= K >= 0;
[in]MB
          MB is INTEGER
          The block size to be used in the blocked QR.
          MB > N. (must be the same as DLATSQR)
[in]NB
          NB is INTEGER
          The column block size to be used in the blocked QR.
          N >= NB >= 1.
[in]A
          A is COMPLEX array, dimension (LDA,K)
          The i-th column must contain the vector which defines the
          blockedelementary reflector H(i), for i = 1,2,...,k, as
          returned by DLATSQR in the first k columns of
          its array argument A.
[in]LDA
          LDA is INTEGER
          The leading dimension of the array A.
          If SIDE = 'L', LDA >= max(1,M);
          if SIDE = 'R', LDA >= max(1,N).
[in]T
          T is COMPLEX array, dimension
          ( N * Number of blocks(CEIL(M-K/MB-K)),
          The blocked upper triangular block reflectors stored in compact form
          as a sequence of upper triangular blocks.  See below
          for further details.
[in]LDT
          LDT is INTEGER
          The leading dimension of the array T.  LDT >= NB.
[in,out]C
          C is COMPLEX array, dimension (LDC,N)
          On entry, the M-by-N matrix C.
          On exit, C is overwritten by Q*C or Q**H*C or C*Q**H or C*Q.
[in]LDC
          LDC is INTEGER
          The leading dimension of the array C. LDC >= max(1,M).
[out]WORK
         (workspace) COMPLEX array, dimension (MAX(1,LWORK))
[in]LWORK
          LWORK is INTEGER
          The dimension of the array WORK.

          If SIDE = 'L', LWORK >= max(1,N)*NB;
          if SIDE = 'R', LWORK >= max(1,MB)*NB.
          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.
Further Details:
 Tall-Skinny QR (TSQR) performs QR by a sequence of orthogonal transformations,
 representing Q as a product of other orthogonal matrices
   Q = Q(1) * Q(2) * . . . * Q(k)
 where each Q(i) zeros out subdiagonal entries of a block of MB rows of A:
   Q(1) zeros out the subdiagonal entries of rows 1:MB of A
   Q(2) zeros out the bottom MB-N rows of rows [1:N,MB+1:2*MB-N] of A
   Q(3) zeros out the bottom MB-N rows of rows [1:N,2*MB-N+1:3*MB-2*N] of A
   . . .

 Q(1) is computed by GEQRT, which represents Q(1) by Householder vectors
 stored under the diagonal of rows 1:MB of A, and by upper triangular
 block reflectors, stored in array T(1:LDT,1:N).
 For more information see Further Details in GEQRT.

 Q(i) for i>1 is computed by TPQRT, which represents Q(i) by Householder vectors
 stored in rows [(i-1)*(MB-N)+N+1:i*(MB-N)+N] of A, and by upper triangular
 block reflectors, stored in array T(1:LDT,(i-1)*N+1:i*N).
 The last Q(k) may use fewer rows.
 For more information see Further Details in TPQRT.

 For more details of the overall algorithm, see the description of
 Sequential TSQR in Section 2.2 of [1].

 [1] “Communication-Optimal Parallel and Sequential QR and LU Factorizations,”
     J. Demmel, L. Grigori, M. Hoemmen, J. Langou,
     SIAM J. Sci. Comput, vol. 34, no. 1, 2012

Definition at line 198 of file clamtsqr.f.

198 *
199 * -- LAPACK computational routine (version 3.7.1) --
200 * -- LAPACK is a software package provided by Univ. of Tennessee, --
201 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
202 * June 2017
203 *
204 * .. Scalar Arguments ..
205  CHARACTER SIDE, TRANS
206  INTEGER INFO, LDA, M, N, K, MB, NB, LDT, LWORK, LDC
207 * ..
208 * .. Array Arguments ..
209  COMPLEX A( LDA, * ), WORK( * ), C(LDC, * ),
210  $ T( LDT, * )
211 * ..
212 *
213 * =====================================================================
214 *
215 * ..
216 * .. Local Scalars ..
217  LOGICAL LEFT, RIGHT, TRAN, NOTRAN, LQUERY
218  INTEGER I, II, KK, LW, CTR
219 * ..
220 * .. External Functions ..
221  LOGICAL LSAME
222  EXTERNAL lsame
223 * .. External Subroutines ..
224  EXTERNAL cgemqrt, ctpmqrt, xerbla
225 * ..
226 * .. Executable Statements ..
227 *
228 * Test the input arguments
229 *
230  lquery = lwork.LT.0
231  notran = lsame( trans, 'N' )
232  tran = lsame( trans, 'C' )
233  left = lsame( side, 'L' )
234  right = lsame( side, 'R' )
235  IF (left) THEN
236  lw = n * nb
237  ELSE
238  lw = m * nb
239  END IF
240 *
241  info = 0
242  IF( .NOT.left .AND. .NOT.right ) THEN
243  info = -1
244  ELSE IF( .NOT.tran .AND. .NOT.notran ) THEN
245  info = -2
246  ELSE IF( m.LT.0 ) THEN
247  info = -3
248  ELSE IF( n.LT.0 ) THEN
249  info = -4
250  ELSE IF( k.LT.0 ) THEN
251  info = -5
252  ELSE IF( lda.LT.max( 1, k ) ) THEN
253  info = -9
254  ELSE IF( ldt.LT.max( 1, nb) ) THEN
255  info = -11
256  ELSE IF( ldc.LT.max( 1, m ) ) THEN
257  info = -13
258  ELSE IF(( lwork.LT.max(1,lw)).AND.(.NOT.lquery)) THEN
259  info = -15
260  END IF
261 *
262 * Determine the block size if it is tall skinny or short and wide
263 *
264  IF( info.EQ.0) THEN
265  work(1) = lw
266  END IF
267 *
268  IF( info.NE.0 ) THEN
269  CALL xerbla( 'CLAMTSQR', -info )
270  RETURN
271  ELSE IF (lquery) THEN
272  RETURN
273  END IF
274 *
275 * Quick return if possible
276 *
277  IF( min(m,n,k).EQ.0 ) THEN
278  RETURN
279  END IF
280 *
281  IF((mb.LE.k).OR.(mb.GE.max(m,n,k))) THEN
282  CALL cgemqrt( side, trans, m, n, k, nb, a, lda,
283  $ t, ldt, c, ldc, work, info)
284  RETURN
285  END IF
286 *
287  IF(left.AND.notran) THEN
288 *
289 * Multiply Q to the last block of C
290 *
291  kk = mod((m-k),(mb-k))
292  ctr = (m-k)/(mb-k)
293  IF (kk.GT.0) THEN
294  ii=m-kk+1
295  CALL ctpmqrt('L','N',kk , n, k, 0, nb, a(ii,1), lda,
296  $ t(1, ctr*k+1),ldt , c(1,1), ldc,
297  $ c(ii,1), ldc, work, info )
298  ELSE
299  ii=m+1
300  END IF
301 *
302  DO i=ii-(mb-k),mb+1,-(mb-k)
303 *
304 * Multiply Q to the current block of C (I:I+MB,1:N)
305 *
306  ctr = ctr - 1
307  CALL ctpmqrt('L','N',mb-k , n, k, 0,nb, a(i,1), lda,
308  $ t(1,ctr*k+1),ldt, c(1,1), ldc,
309  $ c(i,1), ldc, work, info )
310 
311  END DO
312 *
313 * Multiply Q to the first block of C (1:MB,1:N)
314 *
315  CALL cgemqrt('L','N',mb , n, k, nb, a(1,1), lda, t
316  $ ,ldt ,c(1,1), ldc, work, info )
317 *
318  ELSE IF (left.AND.tran) THEN
319 *
320 * Multiply Q to the first block of C
321 *
322  kk = mod((m-k),(mb-k))
323  ii=m-kk+1
324  ctr = 1
325  CALL cgemqrt('L','C',mb , n, k, nb, a(1,1), lda, t
326  $ ,ldt ,c(1,1), ldc, work, info )
327 *
328  DO i=mb+1,ii-mb+k,(mb-k)
329 *
330 * Multiply Q to the current block of C (I:I+MB,1:N)
331 *
332  CALL ctpmqrt('L','C',mb-k , n, k, 0,nb, a(i,1), lda,
333  $ t(1, ctr*k+1),ldt, c(1,1), ldc,
334  $ c(i,1), ldc, work, info )
335  ctr = ctr + 1
336 *
337  END DO
338  IF(ii.LE.m) THEN
339 *
340 * Multiply Q to the last block of C
341 *
342  CALL ctpmqrt('L','C',kk , n, k, 0,nb, a(ii,1), lda,
343  $ t(1,ctr*k+1), ldt, c(1,1), ldc,
344  $ c(ii,1), ldc, work, info )
345 *
346  END IF
347 *
348  ELSE IF(right.AND.tran) THEN
349 *
350 * Multiply Q to the last block of C
351 *
352  kk = mod((n-k),(mb-k))
353  ctr = (n-k)/(mb-k)
354  IF (kk.GT.0) THEN
355  ii=n-kk+1
356  CALL ctpmqrt('R','C',m , kk, k, 0, nb, a(ii,1), lda,
357  $ t(1, ctr*k+1), ldt, c(1,1), ldc,
358  $ c(1,ii), ldc, work, info )
359  ELSE
360  ii=n+1
361  END IF
362 *
363  DO i=ii-(mb-k),mb+1,-(mb-k)
364 *
365 * Multiply Q to the current block of C (1:M,I:I+MB)
366 *
367  ctr = ctr - 1
368  CALL ctpmqrt('R','C',m , mb-k, k, 0,nb, a(i,1), lda,
369  $ t(1,ctr*k+1), ldt, c(1,1), ldc,
370  $ c(1,i), ldc, work, info )
371  END DO
372 *
373 * Multiply Q to the first block of C (1:M,1:MB)
374 *
375  CALL cgemqrt('R','C',m , mb, k, nb, a(1,1), lda, t
376  $ ,ldt ,c(1,1), ldc, work, info )
377 *
378  ELSE IF (right.AND.notran) THEN
379 *
380 * Multiply Q to the first block of C
381 *
382  kk = mod((n-k),(mb-k))
383  ii=n-kk+1
384  ctr = 1
385  CALL cgemqrt('R','N', m, mb , k, nb, a(1,1), lda, t
386  $ ,ldt ,c(1,1), ldc, work, info )
387 *
388  DO i=mb+1,ii-mb+k,(mb-k)
389 *
390 * Multiply Q to the current block of C (1:M,I:I+MB)
391 *
392  CALL ctpmqrt('R','N', m, mb-k, k, 0,nb, a(i,1), lda,
393  $ t(1,ctr*k+1),ldt, c(1,1), ldc,
394  $ c(1,i), ldc, work, info )
395  ctr = ctr + 1
396 *
397  END DO
398  IF(ii.LE.n) THEN
399 *
400 * Multiply Q to the last block of C
401 *
402  CALL ctpmqrt('R','N', m, kk , k, 0,nb, a(ii,1), lda,
403  $ t(1,ctr*k+1),ldt, c(1,1), ldc,
404  $ c(1,ii), ldc, work, info )
405 *
406  END IF
407 *
408  END IF
409 *
410  work(1) = lw
411  RETURN
412 *
413 * End of CLAMTSQR
414 *
Here is the call graph for this function:
Here is the caller graph for this function:
ctpmqrt
subroutine ctpmqrt(SIDE, TRANS, M, N, K, L, NB, V, LDV, T, LDT, A, LDA, B, LDB, WORK, INFO)
CTPMQRT
Definition: ctpmqrt.f:218
xerbla
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
lsame
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:55
cgemqrt
subroutine cgemqrt(SIDE, TRANS, M, N, K, NB, V, LDV, T, LDT, C, LDC, WORK, INFO)
CGEMQRT
Definition: cgemqrt.f:170