LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
recursive subroutine cgeqrt3 ( integer  M,
integer  N,
complex, dimension( lda, * )  A,
integer  LDA,
complex, dimension( ldt, * )  T,
integer  LDT,
integer  INFO 
)

CGEQRT3 recursively computes a QR factorization of a general real or complex matrix using the compact WY representation of Q.

Download CGEQRT3 + dependencies [TGZ] [ZIP] [TXT]

Purpose:
 CGEQRT3 recursively computes a QR factorization of a complex M-by-N matrix A, 
 using the compact WY representation of Q. 

 Based on the algorithm of Elmroth and Gustavson, 
 IBM J. Res. Develop. Vol 44 No. 4 July 2000.
Parameters
[in]M
          M is INTEGER
          The number of rows of the matrix A.  M >= N.
[in]N
          N is INTEGER
          The number of columns of the matrix A.  N >= 0.
[in,out]A
          A is COMPLEX array, dimension (LDA,N)
          On entry, the complex M-by-N matrix A.  On exit, the elements on and
          above the diagonal contain the N-by-N upper triangular matrix R; the
          elements below the diagonal are the columns of V.  See below for
          further details.
[in]LDA
          LDA is INTEGER
          The leading dimension of the array A.  LDA >= max(1,M).
[out]T
          T is COMPLEX array, dimension (LDT,N)
          The N-by-N upper triangular factor of the block reflector.
          The elements on and above the diagonal contain the block
          reflector T; the elements below the diagonal are not used.
          See below for further details.
[in]LDT
          LDT is INTEGER
          The leading dimension of the array T.  LDT >= max(1,N).
[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.
Date
June 2016
Further Details:
  The matrix V stores the elementary reflectors H(i) in the i-th column
  below the diagonal. For example, if M=5 and N=3, the matrix V is

               V = (  1       )
                   ( v1  1    )
                   ( v1 v2  1 )
                   ( v1 v2 v3 )
                   ( v1 v2 v3 )

  where the vi's represent the vectors which define H(i), which are returned
  in the matrix A.  The 1's along the diagonal of V are not stored in A.  The
  block reflector H is then given by

               H = I - V * T * V**H

  where V**H is the conjugate transpose of V.

  For details of the algorithm, see Elmroth and Gustavson (cited above).

Definition at line 134 of file cgeqrt3.f.

134 *
135 * -- LAPACK computational routine (version 3.6.1) --
136 * -- LAPACK is a software package provided by Univ. of Tennessee, --
137 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
138 * June 2016
139 *
140 * .. Scalar Arguments ..
141  INTEGER info, lda, m, n, ldt
142 * ..
143 * .. Array Arguments ..
144  COMPLEX a( lda, * ), t( ldt, * )
145 * ..
146 *
147 * =====================================================================
148 *
149 * .. Parameters ..
150  COMPLEX one
151  parameter( one = (1.0,0.0) )
152 * ..
153 * .. Local Scalars ..
154  INTEGER i, i1, j, j1, n1, n2, iinfo
155 * ..
156 * .. External Subroutines ..
157  EXTERNAL clarfg, ctrmm, cgemm, xerbla
158 * ..
159 * .. Executable Statements ..
160 *
161  info = 0
162  IF( n .LT. 0 ) THEN
163  info = -2
164  ELSE IF( m .LT. n ) THEN
165  info = -1
166  ELSE IF( lda .LT. max( 1, m ) ) THEN
167  info = -4
168  ELSE IF( ldt .LT. max( 1, n ) ) THEN
169  info = -6
170  END IF
171  IF( info.NE.0 ) THEN
172  CALL xerbla( 'CGEQRT3', -info )
173  RETURN
174  END IF
175 *
176  IF( n.EQ.1 ) THEN
177 *
178 * Compute Householder transform when N=1
179 *
180  CALL clarfg( m, a(1,1), a( min( 2, m ), 1 ), 1, t(1,1) )
181 *
182  ELSE
183 *
184 * Otherwise, split A into blocks...
185 *
186  n1 = n/2
187  n2 = n-n1
188  j1 = min( n1+1, n )
189  i1 = min( n+1, m )
190 *
191 * Compute A(1:M,1:N1) <- (Y1,R1,T1), where Q1 = I - Y1 T1 Y1**H
192 *
193  CALL cgeqrt3( m, n1, a, lda, t, ldt, iinfo )
194 *
195 * Compute A(1:M,J1:N) = Q1**H A(1:M,J1:N) [workspace: T(1:N1,J1:N)]
196 *
197  DO j=1,n2
198  DO i=1,n1
199  t( i, j+n1 ) = a( i, j+n1 )
200  END DO
201  END DO
202  CALL ctrmm( 'L', 'L', 'C', 'U', n1, n2, one,
203  & a, lda, t( 1, j1 ), ldt )
204 *
205  CALL cgemm( 'C', 'N', n1, n2, m-n1, one, a( j1, 1 ), lda,
206  & a( j1, j1 ), lda, one, t( 1, j1 ), ldt)
207 *
208  CALL ctrmm( 'L', 'U', 'C', 'N', n1, n2, one,
209  & t, ldt, t( 1, j1 ), ldt )
210 *
211  CALL cgemm( 'N', 'N', m-n1, n2, n1, -one, a( j1, 1 ), lda,
212  & t( 1, j1 ), ldt, one, a( j1, j1 ), lda )
213 *
214  CALL ctrmm( 'L', 'L', 'N', 'U', n1, n2, one,
215  & a, lda, t( 1, j1 ), ldt )
216 *
217  DO j=1,n2
218  DO i=1,n1
219  a( i, j+n1 ) = a( i, j+n1 ) - t( i, j+n1 )
220  END DO
221  END DO
222 *
223 * Compute A(J1:M,J1:N) <- (Y2,R2,T2) where Q2 = I - Y2 T2 Y2**H
224 *
225  CALL cgeqrt3( m-n1, n2, a( j1, j1 ), lda,
226  & t( j1, j1 ), ldt, iinfo )
227 *
228 * Compute T3 = T(1:N1,J1:N) = -T1 Y1**H Y2 T2
229 *
230  DO i=1,n1
231  DO j=1,n2
232  t( i, j+n1 ) = conjg(a( j+n1, i ))
233  END DO
234  END DO
235 *
236  CALL ctrmm( 'R', 'L', 'N', 'U', n1, n2, one,
237  & a( j1, j1 ), lda, t( 1, j1 ), ldt )
238 *
239  CALL cgemm( 'C', 'N', n1, n2, m-n, one, a( i1, 1 ), lda,
240  & a( i1, j1 ), lda, one, t( 1, j1 ), ldt )
241 *
242  CALL ctrmm( 'L', 'U', 'N', 'N', n1, n2, -one, t, ldt,
243  & t( 1, j1 ), ldt )
244 *
245  CALL ctrmm( 'R', 'U', 'N', 'N', n1, n2, one,
246  & t( j1, j1 ), ldt, t( 1, j1 ), ldt )
247 *
248 * Y = (Y1,Y2); R = [ R1 A(1:N1,J1:N) ]; T = [T1 T3]
249 * [ 0 R2 ] [ 0 T2]
250 *
251  END IF
252 *
253  RETURN
254 *
255 * End of CGEQRT3
256 *
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine ctrmm(SIDE, UPLO, TRANSA, DIAG, M, N, ALPHA, A, LDA, B, LDB)
CTRMM
Definition: ctrmm.f:179
recursive subroutine cgeqrt3(M, N, A, LDA, T, LDT, INFO)
CGEQRT3 recursively computes a QR factorization of a general real or complex matrix using the compact...
Definition: cgeqrt3.f:134
subroutine cgemm(TRANSA, TRANSB, M, N, K, ALPHA, A, LDA, B, LDB, BETA, C, LDC)
CGEMM
Definition: cgemm.f:189
subroutine clarfg(N, ALPHA, X, INCX, TAU)
CLARFG generates an elementary reflector (Householder matrix).
Definition: clarfg.f:108

Here is the call graph for this function:

Here is the caller graph for this function: