LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine dlahr2 ( integer  N,
integer  K,
integer  NB,
double precision, dimension( lda, * )  A,
integer  LDA,
double precision, dimension( nb )  TAU,
double precision, dimension( ldt, nb )  T,
integer  LDT,
double precision, dimension( ldy, nb )  Y,
integer  LDY 
)

DLAHR2 reduces the specified number of first columns of a general rectangular matrix A so that elements below the specified subdiagonal are zero, and returns auxiliary matrices which are needed to apply the transformation to the unreduced part of A.

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

Purpose:
 DLAHR2 reduces the first NB columns of A real general n-BY-(n-k+1)
 matrix A so that elements below the k-th subdiagonal are zero. The
 reduction is performed by an orthogonal similarity transformation
 Q**T * A * Q. The routine returns the matrices V and T which determine
 Q as a block reflector I - V*T*V**T, and also the matrix Y = A * V * T.

 This is an auxiliary routine called by DGEHRD.
Parameters
[in]N
          N is INTEGER
          The order of the matrix A.
[in]K
          K is INTEGER
          The offset for the reduction. Elements below the k-th
          subdiagonal in the first NB columns are reduced to zero.
          K < N.
[in]NB
          NB is INTEGER
          The number of columns to be reduced.
[in,out]A
          A is DOUBLE PRECISION array, dimension (LDA,N-K+1)
          On entry, the n-by-(n-k+1) general matrix A.
          On exit, the elements on and above the k-th subdiagonal in
          the first NB columns are overwritten with the corresponding
          elements of the reduced matrix; the elements below the k-th
          subdiagonal, with the array TAU, represent the matrix Q as a
          product of elementary reflectors. The other columns of A are
          unchanged. See Further Details.
[in]LDA
          LDA is INTEGER
          The leading dimension of the array A.  LDA >= max(1,N).
[out]TAU
          TAU is DOUBLE PRECISION array, dimension (NB)
          The scalar factors of the elementary reflectors. See Further
          Details.
[out]T
          T is DOUBLE PRECISION array, dimension (LDT,NB)
          The upper triangular matrix T.
[in]LDT
          LDT is INTEGER
          The leading dimension of the array T.  LDT >= NB.
[out]Y
          Y is DOUBLE PRECISION array, dimension (LDY,NB)
          The n-by-nb matrix Y.
[in]LDY
          LDY is INTEGER
          The leading dimension of the array Y. LDY >= N.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
September 2012
Further Details:
  The matrix Q is represented as a product of nb elementary reflectors

     Q = H(1) H(2) . . . H(nb).

  Each H(i) has the form

     H(i) = I - tau * v * v**T

  where tau is a real scalar, and v is a real vector with
  v(1:i+k-1) = 0, v(i+k) = 1; v(i+k+1:n) is stored on exit in
  A(i+k+1:n,i), and tau in TAU(i).

  The elements of the vectors v together form the (n-k+1)-by-nb matrix
  V which is needed, with T and Y, to apply the transformation to the
  unreduced part of the matrix, using an update of the form:
  A := (I - V*T*V**T) * (A - Y*V**T).

  The contents of A on exit are illustrated by the following example
  with n = 7, k = 3 and nb = 2:

     ( a   a   a   a   a )
     ( a   a   a   a   a )
     ( a   a   a   a   a )
     ( h   h   a   a   a )
     ( v1  h   a   a   a )
     ( v1  v2  a   a   a )
     ( v1  v2  a   a   a )

  where a denotes an element of the original matrix A, h denotes a
  modified element of the upper Hessenberg matrix H, and vi denotes an
  element of the vector defining H(i).

  This subroutine is a slight modification of LAPACK-3.0's DLAHRD
  incorporating improvements proposed by Quintana-Orti and Van de
  Gejin. Note that the entries of A(1:K,2:NB) differ from those
  returned by the original LAPACK-3.0's DLAHRD routine. (This
  subroutine is not backward compatible with LAPACK-3.0's DLAHRD.)
References:
Gregorio Quintana-Orti and Robert van de Geijn, "Improving the performance of reduction to Hessenberg form," ACM Transactions on Mathematical Software, 32(2):180-194, June 2006.

Definition at line 183 of file dlahr2.f.

183 *
184 * -- LAPACK auxiliary routine (version 3.4.2) --
185 * -- LAPACK is a software package provided by Univ. of Tennessee, --
186 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
187 * September 2012
188 *
189 * .. Scalar Arguments ..
190  INTEGER k, lda, ldt, ldy, n, nb
191 * ..
192 * .. Array Arguments ..
193  DOUBLE PRECISION a( lda, * ), t( ldt, nb ), tau( nb ),
194  $ y( ldy, nb )
195 * ..
196 *
197 * =====================================================================
198 *
199 * .. Parameters ..
200  DOUBLE PRECISION zero, one
201  parameter ( zero = 0.0d+0,
202  $ one = 1.0d+0 )
203 * ..
204 * .. Local Scalars ..
205  INTEGER i
206  DOUBLE PRECISION ei
207 * ..
208 * .. External Subroutines ..
209  EXTERNAL daxpy, dcopy, dgemm, dgemv, dlacpy,
210  $ dlarfg, dscal, dtrmm, dtrmv
211 * ..
212 * .. Intrinsic Functions ..
213  INTRINSIC min
214 * ..
215 * .. Executable Statements ..
216 *
217 * Quick return if possible
218 *
219  IF( n.LE.1 )
220  $ RETURN
221 *
222  DO 10 i = 1, nb
223  IF( i.GT.1 ) THEN
224 *
225 * Update A(K+1:N,I)
226 *
227 * Update I-th column of A - Y * V**T
228 *
229  CALL dgemv( 'NO TRANSPOSE', n-k, i-1, -one, y(k+1,1), ldy,
230  $ a( k+i-1, 1 ), lda, one, a( k+1, i ), 1 )
231 *
232 * Apply I - V * T**T * V**T to this column (call it b) from the
233 * left, using the last column of T as workspace
234 *
235 * Let V = ( V1 ) and b = ( b1 ) (first I-1 rows)
236 * ( V2 ) ( b2 )
237 *
238 * where V1 is unit lower triangular
239 *
240 * w := V1**T * b1
241 *
242  CALL dcopy( i-1, a( k+1, i ), 1, t( 1, nb ), 1 )
243  CALL dtrmv( 'Lower', 'Transpose', 'UNIT',
244  $ i-1, a( k+1, 1 ),
245  $ lda, t( 1, nb ), 1 )
246 *
247 * w := w + V2**T * b2
248 *
249  CALL dgemv( 'Transpose', n-k-i+1, i-1,
250  $ one, a( k+i, 1 ),
251  $ lda, a( k+i, i ), 1, one, t( 1, nb ), 1 )
252 *
253 * w := T**T * w
254 *
255  CALL dtrmv( 'Upper', 'Transpose', 'NON-UNIT',
256  $ i-1, t, ldt,
257  $ t( 1, nb ), 1 )
258 *
259 * b2 := b2 - V2*w
260 *
261  CALL dgemv( 'NO TRANSPOSE', n-k-i+1, i-1, -one,
262  $ a( k+i, 1 ),
263  $ lda, t( 1, nb ), 1, one, a( k+i, i ), 1 )
264 *
265 * b1 := b1 - V1*w
266 *
267  CALL dtrmv( 'Lower', 'NO TRANSPOSE',
268  $ 'UNIT', i-1,
269  $ a( k+1, 1 ), lda, t( 1, nb ), 1 )
270  CALL daxpy( i-1, -one, t( 1, nb ), 1, a( k+1, i ), 1 )
271 *
272  a( k+i-1, i-1 ) = ei
273  END IF
274 *
275 * Generate the elementary reflector H(I) to annihilate
276 * A(K+I+1:N,I)
277 *
278  CALL dlarfg( n-k-i+1, a( k+i, i ), a( min( k+i+1, n ), i ), 1,
279  $ tau( i ) )
280  ei = a( k+i, i )
281  a( k+i, i ) = one
282 *
283 * Compute Y(K+1:N,I)
284 *
285  CALL dgemv( 'NO TRANSPOSE', n-k, n-k-i+1,
286  $ one, a( k+1, i+1 ),
287  $ lda, a( k+i, i ), 1, zero, y( k+1, i ), 1 )
288  CALL dgemv( 'Transpose', n-k-i+1, i-1,
289  $ one, a( k+i, 1 ), lda,
290  $ a( k+i, i ), 1, zero, t( 1, i ), 1 )
291  CALL dgemv( 'NO TRANSPOSE', n-k, i-1, -one,
292  $ y( k+1, 1 ), ldy,
293  $ t( 1, i ), 1, one, y( k+1, i ), 1 )
294  CALL dscal( n-k, tau( i ), y( k+1, i ), 1 )
295 *
296 * Compute T(1:I,I)
297 *
298  CALL dscal( i-1, -tau( i ), t( 1, i ), 1 )
299  CALL dtrmv( 'Upper', 'No Transpose', 'NON-UNIT',
300  $ i-1, t, ldt,
301  $ t( 1, i ), 1 )
302  t( i, i ) = tau( i )
303 *
304  10 CONTINUE
305  a( k+nb, nb ) = ei
306 *
307 * Compute Y(1:K,1:NB)
308 *
309  CALL dlacpy( 'ALL', k, nb, a( 1, 2 ), lda, y, ldy )
310  CALL dtrmm( 'RIGHT', 'Lower', 'NO TRANSPOSE',
311  $ 'UNIT', k, nb,
312  $ one, a( k+1, 1 ), lda, y, ldy )
313  IF( n.GT.k+nb )
314  $ CALL dgemm( 'NO TRANSPOSE', 'NO TRANSPOSE', k,
315  $ nb, n-k-nb, one,
316  $ a( 1, 2+nb ), lda, a( k+1+nb, 1 ), lda, one, y,
317  $ ldy )
318  CALL dtrmm( 'RIGHT', 'Upper', 'NO TRANSPOSE',
319  $ 'NON-UNIT', k, nb,
320  $ one, t, ldt, y, ldy )
321 *
322  RETURN
323 *
324 * End of DLAHR2
325 *
subroutine dcopy(N, DX, INCX, DY, INCY)
DCOPY
Definition: dcopy.f:53
subroutine dgemv(TRANS, M, N, ALPHA, A, LDA, X, INCX, BETA, Y, INCY)
DGEMV
Definition: dgemv.f:158
subroutine dtrmm(SIDE, UPLO, TRANSA, DIAG, M, N, ALPHA, A, LDA, B, LDB)
DTRMM
Definition: dtrmm.f:179
subroutine dlacpy(UPLO, M, N, A, LDA, B, LDB)
DLACPY copies all or part of one two-dimensional array to another.
Definition: dlacpy.f:105
subroutine daxpy(N, DA, DX, INCX, DY, INCY)
DAXPY
Definition: daxpy.f:54
subroutine dgemm(TRANSA, TRANSB, M, N, K, ALPHA, A, LDA, B, LDB, BETA, C, LDC)
DGEMM
Definition: dgemm.f:189
subroutine dscal(N, DA, DX, INCX)
DSCAL
Definition: dscal.f:55
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: