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

◆ clamswlq()

subroutine clamswlq ( 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 
)

CLAMSWLQ

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


                    SIDE = 'L'     SIDE = 'R'
    TRANS = 'N':      Q * C          C * Q
    TRANS = 'T':      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 (CLASWLQ)
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 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
          CLASWLQ 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 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 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,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 195 of file clamswlq.f.

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