LAPACK 3.12.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 197 of file clamtsqr.f.

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