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

◆ dlahr2()

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.
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 180 of file dlahr2.f.

181*
182* -- LAPACK auxiliary routine --
183* -- LAPACK is a software package provided by Univ. of Tennessee, --
184* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
185*
186* .. Scalar Arguments ..
187 INTEGER K, LDA, LDT, LDY, N, NB
188* ..
189* .. Array Arguments ..
190 DOUBLE PRECISION A( LDA, * ), T( LDT, NB ), TAU( NB ),
191 $ Y( LDY, NB )
192* ..
193*
194* =====================================================================
195*
196* .. Parameters ..
197 DOUBLE PRECISION ZERO, ONE
198 parameter( zero = 0.0d+0,
199 $ one = 1.0d+0 )
200* ..
201* .. Local Scalars ..
202 INTEGER I
203 DOUBLE PRECISION EI
204* ..
205* .. External Subroutines ..
206 EXTERNAL daxpy, dcopy, dgemm, dgemv, dlacpy,
208* ..
209* .. Intrinsic Functions ..
210 INTRINSIC min
211* ..
212* .. Executable Statements ..
213*
214* Quick return if possible
215*
216 IF( n.LE.1 )
217 $ RETURN
218*
219 DO 10 i = 1, nb
220 IF( i.GT.1 ) THEN
221*
222* Update A(K+1:N,I)
223*
224* Update I-th column of A - Y * V**T
225*
226 CALL dgemv( 'NO TRANSPOSE', n-k, i-1, -one, y(k+1,1), ldy,
227 $ a( k+i-1, 1 ), lda, one, a( k+1, i ), 1 )
228*
229* Apply I - V * T**T * V**T to this column (call it b) from the
230* left, using the last column of T as workspace
231*
232* Let V = ( V1 ) and b = ( b1 ) (first I-1 rows)
233* ( V2 ) ( b2 )
234*
235* where V1 is unit lower triangular
236*
237* w := V1**T * b1
238*
239 CALL dcopy( i-1, a( k+1, i ), 1, t( 1, nb ), 1 )
240 CALL dtrmv( 'Lower', 'Transpose', 'UNIT',
241 $ i-1, a( k+1, 1 ),
242 $ lda, t( 1, nb ), 1 )
243*
244* w := w + V2**T * b2
245*
246 CALL dgemv( 'Transpose', n-k-i+1, i-1,
247 $ one, a( k+i, 1 ),
248 $ lda, a( k+i, i ), 1, one, t( 1, nb ), 1 )
249*
250* w := T**T * w
251*
252 CALL dtrmv( 'Upper', 'Transpose', 'NON-UNIT',
253 $ i-1, t, ldt,
254 $ t( 1, nb ), 1 )
255*
256* b2 := b2 - V2*w
257*
258 CALL dgemv( 'NO TRANSPOSE', n-k-i+1, i-1, -one,
259 $ a( k+i, 1 ),
260 $ lda, t( 1, nb ), 1, one, a( k+i, i ), 1 )
261*
262* b1 := b1 - V1*w
263*
264 CALL dtrmv( 'Lower', 'NO TRANSPOSE',
265 $ 'UNIT', i-1,
266 $ a( k+1, 1 ), lda, t( 1, nb ), 1 )
267 CALL daxpy( i-1, -one, t( 1, nb ), 1, a( k+1, i ), 1 )
268*
269 a( k+i-1, i-1 ) = ei
270 END IF
271*
272* Generate the elementary reflector H(I) to annihilate
273* A(K+I+1:N,I)
274*
275 CALL dlarfg( n-k-i+1, a( k+i, i ), a( min( k+i+1, n ), i ), 1,
276 $ tau( i ) )
277 ei = a( k+i, i )
278 a( k+i, i ) = one
279*
280* Compute Y(K+1:N,I)
281*
282 CALL dgemv( 'NO TRANSPOSE', n-k, n-k-i+1,
283 $ one, a( k+1, i+1 ),
284 $ lda, a( k+i, i ), 1, zero, y( k+1, i ), 1 )
285 CALL dgemv( 'Transpose', n-k-i+1, i-1,
286 $ one, a( k+i, 1 ), lda,
287 $ a( k+i, i ), 1, zero, t( 1, i ), 1 )
288 CALL dgemv( 'NO TRANSPOSE', n-k, i-1, -one,
289 $ y( k+1, 1 ), ldy,
290 $ t( 1, i ), 1, one, y( k+1, i ), 1 )
291 CALL dscal( n-k, tau( i ), y( k+1, i ), 1 )
292*
293* Compute T(1:I,I)
294*
295 CALL dscal( i-1, -tau( i ), t( 1, i ), 1 )
296 CALL dtrmv( 'Upper', 'No Transpose', 'NON-UNIT',
297 $ i-1, t, ldt,
298 $ t( 1, i ), 1 )
299 t( i, i ) = tau( i )
300*
301 10 CONTINUE
302 a( k+nb, nb ) = ei
303*
304* Compute Y(1:K,1:NB)
305*
306 CALL dlacpy( 'ALL', k, nb, a( 1, 2 ), lda, y, ldy )
307 CALL dtrmm( 'RIGHT', 'Lower', 'NO TRANSPOSE',
308 $ 'UNIT', k, nb,
309 $ one, a( k+1, 1 ), lda, y, ldy )
310 IF( n.GT.k+nb )
311 $ CALL dgemm( 'NO TRANSPOSE', 'NO TRANSPOSE', k,
312 $ nb, n-k-nb, one,
313 $ a( 1, 2+nb ), lda, a( k+1+nb, 1 ), lda, one, y,
314 $ ldy )
315 CALL dtrmm( 'RIGHT', 'Upper', 'NO TRANSPOSE',
316 $ 'NON-UNIT', k, nb,
317 $ one, t, ldt, y, ldy )
318*
319 RETURN
320*
321* End of DLAHR2
322*
subroutine daxpy(n, da, dx, incx, dy, incy)
DAXPY
Definition daxpy.f:89
subroutine dcopy(n, dx, incx, dy, incy)
DCOPY
Definition dcopy.f:82
subroutine dgemm(transa, transb, m, n, k, alpha, a, lda, b, ldb, beta, c, ldc)
DGEMM
Definition dgemm.f:188
subroutine dgemv(trans, m, n, alpha, a, lda, x, incx, beta, y, incy)
DGEMV
Definition dgemv.f:158
subroutine dlacpy(uplo, m, n, a, lda, b, ldb)
DLACPY copies all or part of one two-dimensional array to another.
Definition dlacpy.f:103
subroutine dlarfg(n, alpha, x, incx, tau)
DLARFG generates an elementary reflector (Householder matrix).
Definition dlarfg.f:106
subroutine dscal(n, da, dx, incx)
DSCAL
Definition dscal.f:79
subroutine dtrmm(side, uplo, transa, diag, m, n, alpha, a, lda, b, ldb)
DTRMM
Definition dtrmm.f:177
subroutine dtrmv(uplo, trans, diag, n, a, lda, x, incx)
DTRMV
Definition dtrmv.f:147
Here is the call graph for this function:
Here is the caller graph for this function: