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