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

◆ zlamswlq()

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 
)

ZLAMSWLQ

Purpose:
    ZLAMSWLQ 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 complex unitary 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**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 C.  M >=0.
[in]N
          N is INTEGER
          The number of columns of the matrix C. N >= 0.
[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 LQ.
          M >= MB >= 1
[in]NB
          NB is INTEGER
          The column block size to be used in the blocked LQ.
          NB > M.
[in]A
          A is COMPLEX*16 array, dimension
                               (LDA,M) if SIDE = 'L',
                               (LDA,N) if SIDE = 'R'
          The i-th row must contain the vector which defines the blocked
          elementary reflector H(i), for i = 1,2,...,k, as returned by
          ZLASWLQ in the first k rows of its array argument A.
[in]LDA
          LDA is INTEGER
          The leading dimension of the array A. LDA >= MAX(1,K).
[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**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*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 unitary transformations,
 representing Q as a product of other unitary 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 TPLQT.

 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 193 of file zlamswlq.f.

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