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