LAPACK  3.7.0
LAPACK: Linear Algebra PACKage
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 DGEQRT3 + 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.
Date
December 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 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 133 of file dgelqt3.f.

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

Here is the call graph for this function:

Here is the caller graph for this function: