LAPACK 3.12.1
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 128 of file dgelqt3.f.

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