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

DLAMQRTS 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 (DLASWLQ)

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 DOUBLE PRECISION 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 DOUBLE PRECISION 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 DOUBLE PRECISION 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) DOUBLE PRECISION 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 dlamswlq.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
211 * ..
212 * .. Array Arguments ..
213  DOUBLE PRECISION 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, lw
223 * ..
224 * .. External Functions ..
225  LOGICAL lsame
226  EXTERNAL lsame
227 * .. External Subroutines ..
228  EXTERNAL dtpmlqt, dgemlqt, 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, 'T' )
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( 'DLAMSWLQ', -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 dgemlqt( 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  IF (kk.GT.0) THEN
294  ii=m-kk+1
295  CALL dtpmlqt('L','T',kk , n, k, 0, mb, a(1,ii), 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-(nb-k),nb+1,-(nb-k)
303 *
304 * Multiply Q to the current block of C (1:M,I:I+NB)
305 *
306  ctr = ctr - 1
307  CALL dtpmlqt('L','T',nb-k , n, k, 0,mb, a(1,i), 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:M,1:NB)
314 *
315  CALL dgemlqt('L','T',nb , n, k, mb, a(1,1), lda, t
316  $ ,ldt ,c(1,1), ldc, work, info )
317 *
318  ELSE IF (left.AND.notran) THEN
319 *
320 * Multiply Q to the first block of C
321 *
322  kk = mod((m-k),(nb-k))
323  ii=m-kk+1
324  ctr = 1
325  CALL dgemlqt('L','N',nb , n, k, mb, a(1,1), lda, t
326  $ ,ldt ,c(1,1), ldc, work, info )
327 *
328  DO i=nb+1,ii-nb+k,(nb-k)
329 *
330 * Multiply Q to the current block of C (I:I+NB,1:N)
331 *
332  CALL dtpmlqt('L','N',nb-k , n, k, 0,mb, a(1,i), 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 dtpmlqt('L','N',kk , n, k, 0, mb, a(1,ii), 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.notran) THEN
349 *
350 * Multiply Q to the last block of C
351 *
352  kk = mod((n-k),(nb-k))
353  ctr = (n-k)/(nb-k)
354  IF (kk.GT.0) THEN
355  ii=n-kk+1
356  CALL dtpmlqt('R','N',m , kk, k, 0, mb, a(1, ii), 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-(nb-k),nb+1,-(nb-k)
364 *
365 * Multiply Q to the current block of C (1:M,I:I+MB)
366 *
367  ctr = ctr - 1
368  CALL dtpmlqt('R','N', m, nb-k, k, 0, mb, a(1, i), lda,
369  $ t(1,ctr*k+1), ldt, c(1,1), ldc,
370  $ c(1,i), ldc, work, info )
371 *
372  END DO
373 *
374 * Multiply Q to the first block of C (1:M,1:MB)
375 *
376  CALL dgemlqt('R','N',m , nb, k, mb, a(1,1), lda, t
377  $ ,ldt ,c(1,1), ldc, work, info )
378 *
379  ELSE IF (right.AND.tran) THEN
380 *
381 * Multiply Q to the first block of C
382 *
383  kk = mod((n-k),(nb-k))
384  ctr = 1
385  ii=n-kk+1
386  CALL dgemlqt('R','T',m , nb, k, mb, a(1,1), lda, t
387  $ ,ldt ,c(1,1), ldc, work, info )
388 *
389  DO i=nb+1,ii-nb+k,(nb-k)
390 *
391 * Multiply Q to the current block of C (1:M,I:I+MB)
392 *
393  CALL dtpmlqt('R','T',m , nb-k, k, 0,mb, a(1,i), lda,
394  $ t(1,ctr*k+1), ldt, c(1,1), ldc,
395  $ c(1,i), ldc, work, info )
396  ctr = ctr + 1
397 *
398  END DO
399  IF(ii.LE.n) THEN
400 *
401 * Multiply Q to the last block of C
402 *
403  CALL dtpmlqt('R','T',m , kk, k, 0,mb, a(1,ii), lda,
404  $ t(1,ctr*k+1),ldt, c(1,1), ldc,
405  $ c(1,ii), ldc, work, info )
406 *
407  END IF
408 *
409  END IF
410 *
411  work(1) = lw
412  RETURN
413 *
414 * End of DLAMSWLQ
415 *
subroutine dgemlqt(SIDE, TRANS, M, N, K, MB, V, LDV, T, LDT, C, LDC, WORK, INFO)
DGEMLQT
Definition: dgemlqt.f:170
subroutine dtpmlqt(SIDE, TRANS, M, N, K, L, MB, V, LDV, T, LDT, A, LDA, B, LDB, WORK, INFO)
DTPMLQT
Definition: dtpmlqt.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: