LAPACK 3.12.1
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
dlahr2.f
Go to the documentation of this file.
1*> \brief \b 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.
2*
3* =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6* http://www.netlib.org/lapack/explore-html/
7*
8*> Download DLAHR2 + dependencies
9*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dlahr2.f">
10*> [TGZ]</a>
11*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dlahr2.f">
12*> [ZIP]</a>
13*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dlahr2.f">
14*> [TXT]</a>
15*
16* Definition:
17* ===========
18*
19* SUBROUTINE DLAHR2( N, K, NB, A, LDA, TAU, T, LDT, Y, LDY )
20*
21* .. Scalar Arguments ..
22* INTEGER K, LDA, LDT, LDY, N, NB
23* ..
24* .. Array Arguments ..
25* DOUBLE PRECISION A( LDA, * ), T( LDT, NB ), TAU( NB ),
26* $ Y( LDY, NB )
27* ..
28*
29*
30*> \par Purpose:
31* =============
32*>
33*> \verbatim
34*>
35*> DLAHR2 reduces the first NB columns of A real general n-BY-(n-k+1)
36*> matrix A so that elements below the k-th subdiagonal are zero. The
37*> reduction is performed by an orthogonal similarity transformation
38*> Q**T * A * Q. The routine returns the matrices V and T which determine
39*> Q as a block reflector I - V*T*V**T, and also the matrix Y = A * V * T.
40*>
41*> This is an auxiliary routine called by DGEHRD.
42*> \endverbatim
43*
44* Arguments:
45* ==========
46*
47*> \param[in] N
48*> \verbatim
49*> N is INTEGER
50*> The order of the matrix A.
51*> \endverbatim
52*>
53*> \param[in] K
54*> \verbatim
55*> K is INTEGER
56*> The offset for the reduction. Elements below the k-th
57*> subdiagonal in the first NB columns are reduced to zero.
58*> K < N.
59*> \endverbatim
60*>
61*> \param[in] NB
62*> \verbatim
63*> NB is INTEGER
64*> The number of columns to be reduced.
65*> \endverbatim
66*>
67*> \param[in,out] A
68*> \verbatim
69*> A is DOUBLE PRECISION array, dimension (LDA,N-K+1)
70*> On entry, the n-by-(n-k+1) general matrix A.
71*> On exit, the elements on and above the k-th subdiagonal in
72*> the first NB columns are overwritten with the corresponding
73*> elements of the reduced matrix; the elements below the k-th
74*> subdiagonal, with the array TAU, represent the matrix Q as a
75*> product of elementary reflectors. The other columns of A are
76*> unchanged. See Further Details.
77*> \endverbatim
78*>
79*> \param[in] LDA
80*> \verbatim
81*> LDA is INTEGER
82*> The leading dimension of the array A. LDA >= max(1,N).
83*> \endverbatim
84*>
85*> \param[out] TAU
86*> \verbatim
87*> TAU is DOUBLE PRECISION array, dimension (NB)
88*> The scalar factors of the elementary reflectors. See Further
89*> Details.
90*> \endverbatim
91*>
92*> \param[out] T
93*> \verbatim
94*> T is DOUBLE PRECISION array, dimension (LDT,NB)
95*> The upper triangular matrix T.
96*> \endverbatim
97*>
98*> \param[in] LDT
99*> \verbatim
100*> LDT is INTEGER
101*> The leading dimension of the array T. LDT >= NB.
102*> \endverbatim
103*>
104*> \param[out] Y
105*> \verbatim
106*> Y is DOUBLE PRECISION array, dimension (LDY,NB)
107*> The n-by-nb matrix Y.
108*> \endverbatim
109*>
110*> \param[in] LDY
111*> \verbatim
112*> LDY is INTEGER
113*> The leading dimension of the array Y. LDY >= N.
114*> \endverbatim
115*
116* Authors:
117* ========
118*
119*> \author Univ. of Tennessee
120*> \author Univ. of California Berkeley
121*> \author Univ. of Colorado Denver
122*> \author NAG Ltd.
123*
124*> \ingroup lahr2
125*
126*> \par Further Details:
127* =====================
128*>
129*> \verbatim
130*>
131*> The matrix Q is represented as a product of nb elementary reflectors
132*>
133*> Q = H(1) H(2) . . . H(nb).
134*>
135*> Each H(i) has the form
136*>
137*> H(i) = I - tau * v * v**T
138*>
139*> where tau is a real scalar, and v is a real vector with
140*> v(1:i+k-1) = 0, v(i+k) = 1; v(i+k+1:n) is stored on exit in
141*> A(i+k+1:n,i), and tau in TAU(i).
142*>
143*> The elements of the vectors v together form the (n-k+1)-by-nb matrix
144*> V which is needed, with T and Y, to apply the transformation to the
145*> unreduced part of the matrix, using an update of the form:
146*> A := (I - V*T*V**T) * (A - Y*V**T).
147*>
148*> The contents of A on exit are illustrated by the following example
149*> with n = 7, k = 3 and nb = 2:
150*>
151*> ( a a a a a )
152*> ( a a a a a )
153*> ( a a a a a )
154*> ( h h a a a )
155*> ( v1 h a a a )
156*> ( v1 v2 a a a )
157*> ( v1 v2 a a a )
158*>
159*> where a denotes an element of the original matrix A, h denotes a
160*> modified element of the upper Hessenberg matrix H, and vi denotes an
161*> element of the vector defining H(i).
162*>
163*> This subroutine is a slight modification of LAPACK-3.0's DLAHRD
164*> incorporating improvements proposed by Quintana-Orti and Van de
165*> Gejin. Note that the entries of A(1:K,2:NB) differ from those
166*> returned by the original LAPACK-3.0's DLAHRD routine. (This
167*> subroutine is not backward compatible with LAPACK-3.0's DLAHRD.)
168*> \endverbatim
169*
170*> \par References:
171* ================
172*>
173*> Gregorio Quintana-Orti and Robert van de Geijn, "Improving the
174*> performance of reduction to Hessenberg form," ACM Transactions on
175*> Mathematical Software, 32(2):180-194, June 2006.
176*>
177* =====================================================================
178 SUBROUTINE dlahr2( N, K, NB, A, LDA, TAU, T, LDT, Y, LDY )
179*
180* -- LAPACK auxiliary routine --
181* -- LAPACK is a software package provided by Univ. of Tennessee, --
182* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
183*
184* .. Scalar Arguments ..
185 INTEGER K, LDA, LDT, LDY, N, NB
186* ..
187* .. Array Arguments ..
188 DOUBLE PRECISION A( LDA, * ), T( LDT, NB ), TAU( NB ),
189 $ Y( LDY, NB )
190* ..
191*
192* =====================================================================
193*
194* .. Parameters ..
195 DOUBLE PRECISION ZERO, ONE
196 parameter( zero = 0.0d+0,
197 $ one = 1.0d+0 )
198* ..
199* .. Local Scalars ..
200 INTEGER I
201 DOUBLE PRECISION EI
202* ..
203* .. External Subroutines ..
204 EXTERNAL daxpy, dcopy, dgemm, dgemv, dlacpy,
206* ..
207* .. Intrinsic Functions ..
208 INTRINSIC min
209* ..
210* .. Executable Statements ..
211*
212* Quick return if possible
213*
214 IF( n.LE.1 )
215 $ RETURN
216*
217 DO 10 i = 1, nb
218 IF( i.GT.1 ) THEN
219*
220* Update A(K+1:N,I)
221*
222* Update I-th column of A - Y * V**T
223*
224 CALL dgemv( 'NO TRANSPOSE', n-k, i-1, -one, y(k+1,1),
225 $ ldy,
226 $ a( k+i-1, 1 ), lda, one, a( k+1, i ), 1 )
227*
228* Apply I - V * T**T * V**T to this column (call it b) from the
229* left, using the last column of T as workspace
230*
231* Let V = ( V1 ) and b = ( b1 ) (first I-1 rows)
232* ( V2 ) ( b2 )
233*
234* where V1 is unit lower triangular
235*
236* w := V1**T * b1
237*
238 CALL dcopy( i-1, a( k+1, i ), 1, t( 1, nb ), 1 )
239 CALL dtrmv( 'Lower', 'Transpose', 'UNIT',
240 $ i-1, a( k+1, 1 ),
241 $ lda, t( 1, nb ), 1 )
242*
243* w := w + V2**T * b2
244*
245 CALL dgemv( 'Transpose', n-k-i+1, i-1,
246 $ one, a( k+i, 1 ),
247 $ lda, a( k+i, i ), 1, one, t( 1, nb ), 1 )
248*
249* w := T**T * w
250*
251 CALL dtrmv( 'Upper', 'Transpose', 'NON-UNIT',
252 $ i-1, t, ldt,
253 $ t( 1, nb ), 1 )
254*
255* b2 := b2 - V2*w
256*
257 CALL dgemv( 'NO TRANSPOSE', n-k-i+1, i-1, -one,
258 $ a( k+i, 1 ),
259 $ lda, t( 1, nb ), 1, one, a( k+i, i ), 1 )
260*
261* b1 := b1 - V1*w
262*
263 CALL dtrmv( 'Lower', 'NO TRANSPOSE',
264 $ 'UNIT', i-1,
265 $ a( k+1, 1 ), lda, t( 1, nb ), 1 )
266 CALL daxpy( i-1, -one, t( 1, nb ), 1, a( k+1, i ), 1 )
267*
268 a( k+i-1, i-1 ) = ei
269 END IF
270*
271* Generate the elementary reflector H(I) to annihilate
272* A(K+I+1:N,I)
273*
274 CALL dlarfg( n-k-i+1, a( k+i, i ), a( min( k+i+1, n ), i ),
275 $ 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*
323 END
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:101
subroutine dlahr2(n, k, nb, a, lda, tau, t, ldt, y, ldy)
DLAHR2 reduces the specified number of first columns of a general rectangular matrix A so that elemen...
Definition dlahr2.f:179
subroutine dlarfg(n, alpha, x, incx, tau)
DLARFG generates an elementary reflector (Householder matrix).
Definition dlarfg.f:104
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