LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine dgeqrt2 ( integer  M,
integer  N,
double precision, dimension( lda, * )  A,
integer  LDA,
double precision, dimension( ldt, * )  T,
integer  LDT,
integer  INFO 
)

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

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

Purpose:
 DGEQRT2 computes a QR factorization of a real M-by-N matrix A, 
 using the compact WY representation of Q. 
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
          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 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
September 2012
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**T

  where V**T is the transpose of V.

Definition at line 129 of file dgeqrt2.f.

129 *
130 * -- LAPACK computational routine (version 3.4.2) --
131 * -- LAPACK is a software package provided by Univ. of Tennessee, --
132 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
133 * September 2012
134 *
135 * .. Scalar Arguments ..
136  INTEGER info, lda, ldt, m, n
137 * ..
138 * .. Array Arguments ..
139  DOUBLE PRECISION a( lda, * ), t( ldt, * )
140 * ..
141 *
142 * =====================================================================
143 *
144 * .. Parameters ..
145  DOUBLE PRECISION one, zero
146  parameter( one = 1.0d+00, zero = 0.0d+00 )
147 * ..
148 * .. Local Scalars ..
149  INTEGER i, k
150  DOUBLE PRECISION aii, alpha
151 * ..
152 * .. External Subroutines ..
153  EXTERNAL dlarfg, dgemv, dger, dtrmv, xerbla
154 * ..
155 * .. Executable Statements ..
156 *
157 * Test the input arguments
158 *
159  info = 0
160  IF( m.LT.0 ) THEN
161  info = -1
162  ELSE IF( n.LT.0 ) THEN
163  info = -2
164  ELSE IF( lda.LT.max( 1, m ) ) THEN
165  info = -4
166  ELSE IF( ldt.LT.max( 1, n ) ) THEN
167  info = -6
168  END IF
169  IF( info.NE.0 ) THEN
170  CALL xerbla( 'DGEQRT2', -info )
171  RETURN
172  END IF
173 *
174  k = min( m, n )
175 *
176  DO i = 1, k
177 *
178 * Generate elem. refl. H(i) to annihilate A(i+1:m,i), tau(I) -> T(I,1)
179 *
180  CALL dlarfg( m-i+1, a( i, i ), a( min( i+1, m ), i ), 1,
181  $ t( i, 1 ) )
182  IF( i.LT.n ) THEN
183 *
184 * Apply H(i) to A(I:M,I+1:N) from the left
185 *
186  aii = a( i, i )
187  a( i, i ) = one
188 *
189 * W(1:N-I) := A(I:M,I+1:N)^H * A(I:M,I) [W = T(:,N)]
190 *
191  CALL dgemv( 'T',m-i+1, n-i, one, a( i, i+1 ), lda,
192  $ a( i, i ), 1, zero, t( 1, n ), 1 )
193 *
194 * A(I:M,I+1:N) = A(I:m,I+1:N) + alpha*A(I:M,I)*W(1:N-1)^H
195 *
196  alpha = -(t( i, 1 ))
197  CALL dger( m-i+1, n-i, alpha, a( i, i ), 1,
198  $ t( 1, n ), 1, a( i, i+1 ), lda )
199  a( i, i ) = aii
200  END IF
201  END DO
202 *
203  DO i = 2, n
204  aii = a( i, i )
205  a( i, i ) = one
206 *
207 * T(1:I-1,I) := alpha * A(I:M,1:I-1)**T * A(I:M,I)
208 *
209  alpha = -t( i, 1 )
210  CALL dgemv( 'T', m-i+1, i-1, alpha, a( i, 1 ), lda,
211  $ a( i, i ), 1, zero, t( 1, i ), 1 )
212  a( i, i ) = aii
213 *
214 * T(1:I-1,I) := T(1:I-1,1:I-1) * T(1:I-1,I)
215 *
216  CALL dtrmv( 'U', 'N', 'N', i-1, t, ldt, t( 1, i ), 1 )
217 *
218 * T(I,I) = tau(I)
219 *
220  t( i, i ) = t( i, 1 )
221  t( i, 1) = zero
222  END DO
223 
224 *
225 * End of DGEQRT2
226 *
subroutine dgemv(TRANS, M, N, ALPHA, A, LDA, X, INCX, BETA, Y, INCY)
DGEMV
Definition: dgemv.f:158
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine dger(M, N, ALPHA, X, INCX, Y, INCY, A, LDA)
DGER
Definition: dger.f:132
subroutine dlarfg(N, ALPHA, X, INCX, TAU)
DLARFG generates an elementary reflector (Householder matrix).
Definition: dlarfg.f:108
subroutine dtrmv(UPLO, TRANS, DIAG, N, A, LDA, X, INCX)
DTRMV
Definition: dtrmv.f:149

Here is the call graph for this function:

Here is the caller graph for this function: