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

◆ 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 complex unitary 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. 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 block size to be used in the blocked QR.
          MB > N. (must be the same as CLATSQR)
[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 CLATSQR 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 unitary transformations,
 representing Q as a product of other unitary 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 195 of file clamtsqr.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, Q
217* ..
218* .. External Functions ..
219 LOGICAL LSAME
220 EXTERNAL lsame
221* .. External Subroutines ..
222 EXTERNAL cgemqrt, ctpmqrt, xerbla
223* ..
224* .. Executable Statements ..
225*
226* Test the input arguments
227*
228 lquery = lwork.LT.0
229 notran = lsame( trans, 'N' )
230 tran = lsame( trans, 'C' )
231 left = lsame( side, 'L' )
232 right = lsame( side, 'R' )
233 IF (left) THEN
234 lw = n * nb
235 q = m
236 ELSE
237 lw = m * nb
238 q = n
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.k ) 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( k.LT.nb .OR. nb.LT.1 ) THEN
253 info = -7
254 ELSE IF( lda.LT.max( 1, q ) ) THEN
255 info = -9
256 ELSE IF( ldt.LT.max( 1, nb) ) THEN
257 info = -11
258 ELSE IF( ldc.LT.max( 1, m ) ) THEN
259 info = -13
260 ELSE IF(( lwork.LT.max(1,lw)).AND.(.NOT.lquery)) THEN
261 info = -15
262 END IF
263*
264* Determine the block size if it is tall skinny or short and wide
265*
266 IF( info.EQ.0) THEN
267 work(1) = lw
268 END IF
269*
270 IF( info.NE.0 ) THEN
271 CALL xerbla( 'CLAMTSQR', -info )
272 RETURN
273 ELSE IF (lquery) THEN
274 RETURN
275 END IF
276*
277* Quick return if possible
278*
279 IF( min(m,n,k).EQ.0 ) THEN
280 RETURN
281 END IF
282*
283 IF((mb.LE.k).OR.(mb.GE.max(m,n,k))) THEN
284 CALL cgemqrt( side, trans, m, n, k, nb, a, lda,
285 $ t, ldt, c, ldc, work, info)
286 RETURN
287 END IF
288*
289 IF(left.AND.notran) THEN
290*
291* Multiply Q to the last block of C
292*
293 kk = mod((m-k),(mb-k))
294 ctr = (m-k)/(mb-k)
295 IF (kk.GT.0) THEN
296 ii=m-kk+1
297 CALL ctpmqrt('L','N',kk , n, k, 0, nb, a(ii,1), lda,
298 $ t(1, ctr*k+1),ldt , c(1,1), ldc,
299 $ c(ii,1), ldc, work, info )
300 ELSE
301 ii=m+1
302 END IF
303*
304 DO i=ii-(mb-k),mb+1,-(mb-k)
305*
306* Multiply Q to the current block of C (I:I+MB,1:N)
307*
308 ctr = ctr - 1
309 CALL ctpmqrt('L','N',mb-k , n, k, 0,nb, a(i,1), lda,
310 $ t(1,ctr*k+1),ldt, c(1,1), ldc,
311 $ c(i,1), ldc, work, info )
312
313 END DO
314*
315* Multiply Q to the first block of C (1:MB,1:N)
316*
317 CALL cgemqrt('L','N',mb , n, k, nb, a(1,1), lda, t
318 $ ,ldt ,c(1,1), ldc, work, info )
319*
320 ELSE IF (left.AND.tran) THEN
321*
322* Multiply Q to the first block of C
323*
324 kk = mod((m-k),(mb-k))
325 ii=m-kk+1
326 ctr = 1
327 CALL cgemqrt('L','C',mb , n, k, nb, a(1,1), lda, t
328 $ ,ldt ,c(1,1), ldc, work, info )
329*
330 DO i=mb+1,ii-mb+k,(mb-k)
331*
332* Multiply Q to the current block of C (I:I+MB,1:N)
333*
334 CALL ctpmqrt('L','C',mb-k , n, k, 0,nb, a(i,1), lda,
335 $ t(1, ctr*k+1),ldt, c(1,1), ldc,
336 $ c(i,1), ldc, work, info )
337 ctr = ctr + 1
338*
339 END DO
340 IF(ii.LE.m) THEN
341*
342* Multiply Q to the last block of C
343*
344 CALL ctpmqrt('L','C',kk , n, k, 0,nb, a(ii,1), lda,
345 $ t(1,ctr*k+1), ldt, c(1,1), ldc,
346 $ c(ii,1), ldc, work, info )
347*
348 END IF
349*
350 ELSE IF(right.AND.tran) THEN
351*
352* Multiply Q to the last block of C
353*
354 kk = mod((n-k),(mb-k))
355 ctr = (n-k)/(mb-k)
356 IF (kk.GT.0) THEN
357 ii=n-kk+1
358 CALL ctpmqrt('R','C',m , kk, k, 0, nb, a(ii,1), lda,
359 $ t(1, ctr*k+1), ldt, c(1,1), ldc,
360 $ c(1,ii), ldc, work, info )
361 ELSE
362 ii=n+1
363 END IF
364*
365 DO i=ii-(mb-k),mb+1,-(mb-k)
366*
367* Multiply Q to the current block of C (1:M,I:I+MB)
368*
369 ctr = ctr - 1
370 CALL ctpmqrt('R','C',m , mb-k, k, 0,nb, a(i,1), lda,
371 $ t(1,ctr*k+1), ldt, c(1,1), ldc,
372 $ c(1,i), ldc, work, info )
373 END DO
374*
375* Multiply Q to the first block of C (1:M,1:MB)
376*
377 CALL cgemqrt('R','C',m , mb, k, nb, a(1,1), lda, t
378 $ ,ldt ,c(1,1), ldc, work, info )
379*
380 ELSE IF (right.AND.notran) THEN
381*
382* Multiply Q to the first block of C
383*
384 kk = mod((n-k),(mb-k))
385 ii=n-kk+1
386 ctr = 1
387 CALL cgemqrt('R','N', m, mb , k, nb, a(1,1), lda, t
388 $ ,ldt ,c(1,1), ldc, work, info )
389*
390 DO i=mb+1,ii-mb+k,(mb-k)
391*
392* Multiply Q to the current block of C (1:M,I:I+MB)
393*
394 CALL ctpmqrt('R','N', m, mb-k, k, 0,nb, a(i,1), 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 ctpmqrt('R','N', m, kk , k, 0,nb, a(ii,1), 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 CLAMTSQR
416*
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:60
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:53
subroutine cgemqrt(SIDE, TRANS, M, N, K, NB, V, LDV, T, LDT, C, LDC, WORK, INFO)
CGEMQRT
Definition: cgemqrt.f:168
subroutine ctpmqrt(SIDE, TRANS, M, N, K, L, NB, V, LDV, T, LDT, A, LDA, B, LDB, WORK, INFO)
CTPMQRT
Definition: ctpmqrt.f:216
Here is the call graph for this function:
Here is the caller graph for this function: