LAPACK  3.7.0
LAPACK: Linear Algebra PACKage
subroutine zlamswlq ( character  SIDE,
character  TRANS,
integer  M,
integer  N,
integer  K,
integer  MB,
integer  NB,
complex*16, dimension( lda, * )  A,
integer  LDA,
complex*16, dimension( ldt, * )  T,
integer  LDT,
complex*16, dimension(ldc, * )  C,
integer  LDC,
complex*16, dimension( * )  WORK,
integer  LWORK,
integer  INFO 
)
Purpose:

ZLAMQRTS overwrites the general real M-by-N matrix C with

SIDE = 'L' SIDE = 'R' TRANS = 'N': Q * C C * Q TRANS = 'T': Q**T * C C * Q**T where Q is a real orthogonal matrix defined as the product of blocked elementary reflectors computed by short wide LQ factorization (ZLASWLQ)

Parameters
[in]SIDE
          SIDE is CHARACTER*1
          = 'L': apply Q or Q**T from the Left;
          = 'R': apply Q or Q**T from the Right.
[in]TRANS
          TRANS is CHARACTER*1
          = 'N':  No transpose, apply Q;
          = 'T':  Transpose, apply Q**T.
[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. N >= M.
[in]K
          K is INTEGER
          The number of elementary reflectors whose product defines
          the matrix Q.
          M >= K >= 0;
[in]MB
          MB is INTEGER
          The row block size to be used in the blocked QR.
          M >= MB >= 1
[in]NB
          NB is INTEGER
          The column block size to be used in the blocked QR.
          NB > M.
[in]NB
          NB is INTEGER
          The block size to be used in the blocked QR.
                MB > M.
[in,out]A
          A is COMPLEX*16 array, dimension (LDA,K)
          The i-th row must contain the vector which defines the blocked
          elementary reflector H(i), for i = 1,2,...,k, as returned by
          DLASWLQ in the first k rows 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*16 array, dimension
          ( M * Number of blocks(CEIL(N-K/NB-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 >= MB.
[in,out]C
          C is COMPLEX*16 array, dimension (LDC,N)
          On entry, the M-by-N matrix C.
          On exit, C is overwritten by Q*C or Q**T*C or C*Q**T or C*Q.
[in]LDC
          LDC is INTEGER
          The leading dimension of the array C. LDC >= max(1,M).
[out]WORK
         (workspace) COMPLEX*16 array, dimension (MAX(1,LWORK))
[in]LWORK
          LWORK is INTEGER
          The dimension of the array WORK.
          If SIDE = 'L', LWORK >= max(1,NB) * MB;
          if SIDE = 'R', LWORK >= max(1,M) * MB.
          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:
Short-Wide LQ (SWLQ) performs LQ 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 upper diagonal entries of a block of NB rows of A: Q(1) zeros out the upper diagonal entries of rows 1:NB of A Q(2) zeros out the bottom MB-N rows of rows [1:M,NB+1:2*NB-M] of A Q(3) zeros out the bottom MB-N rows of rows [1:M,2*NB-M+1:3*NB-2*M] of A . . .

Q(1) is computed by GELQT, 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 GELQT.

Q(i) for i>1 is computed by TPLQT, which represents Q(i) by Householder vectors stored in columns [(i-1)*(NB-M)+M+1:i*(NB-M)+M] of A, and by upper triangular block reflectors, stored in array T(1:LDT,(i-1)*M+1:i*M). 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 202 of file zlamswlq.f.

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

Here is the call graph for this function:

Here is the caller graph for this function: