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.
3 * =========== DOCUMENTATION ===========
5 * Online html documentation available at
6 * http://www.netlib.org/lapack/explore-html/
9 *> Download DLAHR2 + dependencies
10 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dlahr2.f">
12 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dlahr2.f">
14 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dlahr2.f">
21 * SUBROUTINE DLAHR2( N, K, NB, A, LDA, TAU, T, LDT, Y, LDY )
23 * .. Scalar Arguments ..
24 * INTEGER K, LDA, LDT, LDY, N, NB
26 * .. Array Arguments ..
27 * DOUBLE PRECISION A( LDA, * ), T( LDT, NB ), TAU( NB ),
37 *> DLAHR2 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.
43 *> This is an auxiliary routine called by DGEHRD.
52 *> The order of the matrix A.
58 *> The offset for the reduction. Elements below the k-th
59 *> subdiagonal in the first NB columns are reduced to zero.
66 *> The number of columns to be reduced.
71 *> A is DOUBLE PRECISION array, dimension (LDA,N-K+1)
72 *> On entry, the n-by-(n-k+1) general matrix A.
73 *> On exit, the elements on and above the k-th subdiagonal in
74 *> the first NB columns are overwritten with the corresponding
75 *> elements of the reduced matrix; the elements below the k-th
76 *> subdiagonal, with the array TAU, represent the matrix Q as a
77 *> product of elementary reflectors. The other columns of A are
78 *> unchanged. See Further Details.
84 *> The leading dimension of the array A. LDA >= max(1,N).
89 *> TAU is DOUBLE PRECISION array, dimension (NB)
90 *> The scalar factors of the elementary reflectors. See Further
96 *> T is DOUBLE PRECISION array, dimension (LDT,NB)
97 *> The upper triangular matrix T.
103 *> The leading dimension of the array T. LDT >= NB.
108 *> Y is DOUBLE PRECISION array, dimension (LDY,NB)
109 *> The n-by-nb matrix Y.
115 *> The leading dimension of the array Y. LDY >= N.
121 *> \author Univ. of Tennessee
122 *> \author Univ. of California Berkeley
123 *> \author Univ. of Colorado Denver
126 *> \date September 2012
128 *> \ingroup doubleOTHERauxiliary
130 *> \par Further Details:
131 * =====================
135 *> The matrix Q is represented as a product of nb elementary reflectors
137 *> Q = H(1) H(2) . . . H(nb).
139 *> Each H(i) has the form
141 *> H(i) = I - tau * v * v**T
143 *> where tau is a real scalar, and v is a real vector with
144 *> v(1:i+k-1) = 0, v(i+k) = 1; v(i+k+1:n) is stored on exit in
145 *> A(i+k+1:n,i), and tau in TAU(i).
147 *> The elements of the vectors v together form the (n-k+1)-by-nb matrix
148 *> V which is needed, with T and Y, to apply the transformation to the
149 *> unreduced part of the matrix, using an update of the form:
150 *> A := (I - V*T*V**T) * (A - Y*V**T).
152 *> The contents of A on exit are illustrated by the following example
153 *> with n = 7, k = 3 and nb = 2:
163 *> where a denotes an element of the original matrix A, h denotes a
164 *> modified element of the upper Hessenberg matrix H, and vi denotes an
165 *> element of the vector defining H(i).
167 *> This subroutine is a slight modification of LAPACK-3.0's DLAHRD
168 *> incorporating improvements proposed by Quintana-Orti and Van de
169 *> Gejin. Note that the entries of A(1:K,2:NB) differ from those
170 *> returned by the original LAPACK-3.0's DLAHRD routine. (This
171 *> subroutine is not backward compatible with LAPACK-3.0's DLAHRD.)
177 *> Gregorio Quintana-Orti and Robert van de Geijn, "Improving the
178 *> performance of reduction to Hessenberg form," ACM Transactions on
179 *> Mathematical Software, 32(2):180-194, June 2006.
181 * =====================================================================
182 SUBROUTINE DLAHR2( N, K, NB, A, LDA, TAU, T, LDT, Y, LDY )
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..--
189 * .. Scalar Arguments ..
190 INTEGER K, LDA, LDT, LDY, N, NB
192 * .. Array Arguments ..
193 DOUBLE PRECISION A( LDA, * ), T( LDT, NB ), TAU( NB ),
197 * =====================================================================
200 DOUBLE PRECISION ZERO, ONE
201 PARAMETER ( ZERO = 0.0D+0,
204 * .. Local Scalars ..
208 * .. External Subroutines ..
209 EXTERNAL DAXPY, DCOPY, DGEMM, DGEMV, DLACPY,
210 $ DLARFG, DSCAL, DTRMM, DTRMV
212 * .. Intrinsic Functions ..
215 * .. Executable Statements ..
217 * Quick return if possible
227 * Update I-th column of A - Y * V**T
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 )
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
235 * Let V = ( V1 ) and b = ( b1 ) (first I-1 rows)
238 * where V1 is unit lower triangular
242 CALL DCOPY( I-1, A( K+1, I ), 1, T( 1, NB ), 1 )
243 CALL DTRMV( 'Lower', 'Transpose', 'UNIT',
245 $ LDA, T( 1, NB ), 1 )
247 * w := w + V2**T * b2
249 CALL DGEMV( 'Transpose', N-K-I+1, I-1,
251 $ LDA, A( K+I, I ), 1, ONE, T( 1, NB ), 1 )
255 CALL DTRMV( 'Upper', 'Transpose', 'NON-UNIT',
261 CALL DGEMV( 'NO TRANSPOSE', N-K-I+1, I-1, -ONE,
263 $ LDA, T( 1, NB ), 1, ONE, A( K+I, I ), 1 )
267 CALL DTRMV( 'Lower', 'NO TRANSPOSE',
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 )
275 * Generate the elementary reflector H(I) to annihilate
278 CALL DLARFG( N-K-I+1, A( K+I, I ), A( MIN( K+I+1, N ), I ), 1,
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,
293 $ T( 1, I ), 1, ONE, Y( K+1, I ), 1 )
294 CALL DSCAL( N-K, TAU( I ), Y( K+1, I ), 1 )
298 CALL DSCAL( I-1, -TAU( I ), T( 1, I ), 1 )
299 CALL DTRMV( 'Upper', 'No Transpose', 'NON-UNIT',
307 * Compute Y(1:K,1:NB)
309 CALL DLACPY( 'ALL', K, NB, A( 1, 2 ), LDA, Y, LDY )
310 CALL DTRMM( 'RIGHT', 'Lower', 'NO TRANSPOSE',
312 $ ONE, A( K+1, 1 ), LDA, Y, LDY )
314 $ CALL DGEMM( 'NO TRANSPOSE', 'NO TRANSPOSE', K,
316 $ A( 1, 2+NB ), LDA, A( K+1+NB, 1 ), LDA, ONE, Y,
318 CALL DTRMM( 'RIGHT', 'Upper', 'NO TRANSPOSE',
320 $ ONE, T, LDT, Y, LDY )