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

◆ dgelqt3()

recursive subroutine dgelqt3 ( integer  M,
integer  N,
double precision, dimension( lda, * )  A,
integer  LDA,
double precision, dimension( ldt, * )  T,
integer  LDT,
integer  INFO 
)

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

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

Purpose:
 DGELQT3 recursively computes a LQ factorization of a real 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 DOUBLE PRECISION array, dimension (LDA,N)
          On entry, the real M-by-N matrix A.  On exit, the elements on and
          below the diagonal contain the N-by-N lower triangular matrix L; the
          elements above the diagonal are the rows 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 DOUBLE PRECISION 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.
Further Details:
  The matrix V stores the elementary reflectors H(i) in the i-th row
  above the diagonal. For example, if M=5 and N=3, the matrix V is

               V = (  1  v1 v1 v1 v1 )
                   (     1  v2 v2 v2 )
                   (     1  v3 v3 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**T

  where V**T is the transpose of V.

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

Definition at line 130 of file dgelqt3.f.

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