00001 SUBROUTINE DLAHR2( N, K, NB, A, LDA, TAU, T, LDT, Y, LDY ) 00002 * 00003 * -- LAPACK auxiliary routine (version 3.2.1) -- 00004 * -- LAPACK is a software package provided by Univ. of Tennessee, -- 00005 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..-- 00006 * -- April 2009 -- 00007 * 00008 * .. Scalar Arguments .. 00009 INTEGER K, LDA, LDT, LDY, N, NB 00010 * .. 00011 * .. Array Arguments .. 00012 DOUBLE PRECISION A( LDA, * ), T( LDT, NB ), TAU( NB ), 00013 $ Y( LDY, NB ) 00014 * .. 00015 * 00016 * Purpose 00017 * ======= 00018 * 00019 * DLAHR2 reduces the first NB columns of A real general n-BY-(n-k+1) 00020 * matrix A so that elements below the k-th subdiagonal are zero. The 00021 * reduction is performed by an orthogonal similarity transformation 00022 * Q' * A * Q. The routine returns the matrices V and T which determine 00023 * Q as a block reflector I - V*T*V', and also the matrix Y = A * V * T. 00024 * 00025 * This is an auxiliary routine called by DGEHRD. 00026 * 00027 * Arguments 00028 * ========= 00029 * 00030 * N (input) INTEGER 00031 * The order of the matrix A. 00032 * 00033 * K (input) INTEGER 00034 * The offset for the reduction. Elements below the k-th 00035 * subdiagonal in the first NB columns are reduced to zero. 00036 * K < N. 00037 * 00038 * NB (input) INTEGER 00039 * The number of columns to be reduced. 00040 * 00041 * A (input/output) DOUBLE PRECISION array, dimension (LDA,N-K+1) 00042 * On entry, the n-by-(n-k+1) general matrix A. 00043 * On exit, the elements on and above the k-th subdiagonal in 00044 * the first NB columns are overwritten with the corresponding 00045 * elements of the reduced matrix; the elements below the k-th 00046 * subdiagonal, with the array TAU, represent the matrix Q as a 00047 * product of elementary reflectors. The other columns of A are 00048 * unchanged. See Further Details. 00049 * 00050 * LDA (input) INTEGER 00051 * The leading dimension of the array A. LDA >= max(1,N). 00052 * 00053 * TAU (output) DOUBLE PRECISION array, dimension (NB) 00054 * The scalar factors of the elementary reflectors. See Further 00055 * Details. 00056 * 00057 * T (output) DOUBLE PRECISION array, dimension (LDT,NB) 00058 * The upper triangular matrix T. 00059 * 00060 * LDT (input) INTEGER 00061 * The leading dimension of the array T. LDT >= NB. 00062 * 00063 * Y (output) DOUBLE PRECISION array, dimension (LDY,NB) 00064 * The n-by-nb matrix Y. 00065 * 00066 * LDY (input) INTEGER 00067 * The leading dimension of the array Y. LDY >= N. 00068 * 00069 * Further Details 00070 * =============== 00071 * 00072 * The matrix Q is represented as a product of nb elementary reflectors 00073 * 00074 * Q = H(1) H(2) . . . H(nb). 00075 * 00076 * Each H(i) has the form 00077 * 00078 * H(i) = I - tau * v * v' 00079 * 00080 * where tau is a real scalar, and v is a real vector with 00081 * v(1:i+k-1) = 0, v(i+k) = 1; v(i+k+1:n) is stored on exit in 00082 * A(i+k+1:n,i), and tau in TAU(i). 00083 * 00084 * The elements of the vectors v together form the (n-k+1)-by-nb matrix 00085 * V which is needed, with T and Y, to apply the transformation to the 00086 * unreduced part of the matrix, using an update of the form: 00087 * A := (I - V*T*V') * (A - Y*V'). 00088 * 00089 * The contents of A on exit are illustrated by the following example 00090 * with n = 7, k = 3 and nb = 2: 00091 * 00092 * ( a a a a a ) 00093 * ( a a a a a ) 00094 * ( a a a a a ) 00095 * ( h h a a a ) 00096 * ( v1 h a a a ) 00097 * ( v1 v2 a a a ) 00098 * ( v1 v2 a a a ) 00099 * 00100 * where a denotes an element of the original matrix A, h denotes a 00101 * modified element of the upper Hessenberg matrix H, and vi denotes an 00102 * element of the vector defining H(i). 00103 * 00104 * This subroutine is a slight modification of LAPACK-3.0's DLAHRD 00105 * incorporating improvements proposed by Quintana-Orti and Van de 00106 * Gejin. Note that the entries of A(1:K,2:NB) differ from those 00107 * returned by the original LAPACK-3.0's DLAHRD routine. (This 00108 * subroutine is not backward compatible with LAPACK-3.0's DLAHRD.) 00109 * 00110 * References 00111 * ========== 00112 * 00113 * Gregorio Quintana-Orti and Robert van de Geijn, "Improving the 00114 * performance of reduction to Hessenberg form," ACM Transactions on 00115 * Mathematical Software, 32(2):180-194, June 2006. 00116 * 00117 * ===================================================================== 00118 * 00119 * .. Parameters .. 00120 DOUBLE PRECISION ZERO, ONE 00121 PARAMETER ( ZERO = 0.0D+0, 00122 $ ONE = 1.0D+0 ) 00123 * .. 00124 * .. Local Scalars .. 00125 INTEGER I 00126 DOUBLE PRECISION EI 00127 * .. 00128 * .. External Subroutines .. 00129 EXTERNAL DAXPY, DCOPY, DGEMM, DGEMV, DLACPY, 00130 $ DLARFG, DSCAL, DTRMM, DTRMV 00131 * .. 00132 * .. Intrinsic Functions .. 00133 INTRINSIC MIN 00134 * .. 00135 * .. Executable Statements .. 00136 * 00137 * Quick return if possible 00138 * 00139 IF( N.LE.1 ) 00140 $ RETURN 00141 * 00142 DO 10 I = 1, NB 00143 IF( I.GT.1 ) THEN 00144 * 00145 * Update A(K+1:N,I) 00146 * 00147 * Update I-th column of A - Y * V' 00148 * 00149 CALL DGEMV( 'NO TRANSPOSE', N-K, I-1, -ONE, Y(K+1,1), LDY, 00150 $ A( K+I-1, 1 ), LDA, ONE, A( K+1, I ), 1 ) 00151 * 00152 * Apply I - V * T' * V' to this column (call it b) from the 00153 * left, using the last column of T as workspace 00154 * 00155 * Let V = ( V1 ) and b = ( b1 ) (first I-1 rows) 00156 * ( V2 ) ( b2 ) 00157 * 00158 * where V1 is unit lower triangular 00159 * 00160 * w := V1' * b1 00161 * 00162 CALL DCOPY( I-1, A( K+1, I ), 1, T( 1, NB ), 1 ) 00163 CALL DTRMV( 'Lower', 'Transpose', 'UNIT', 00164 $ I-1, A( K+1, 1 ), 00165 $ LDA, T( 1, NB ), 1 ) 00166 * 00167 * w := w + V2'*b2 00168 * 00169 CALL DGEMV( 'Transpose', N-K-I+1, I-1, 00170 $ ONE, A( K+I, 1 ), 00171 $ LDA, A( K+I, I ), 1, ONE, T( 1, NB ), 1 ) 00172 * 00173 * w := T'*w 00174 * 00175 CALL DTRMV( 'Upper', 'Transpose', 'NON-UNIT', 00176 $ I-1, T, LDT, 00177 $ T( 1, NB ), 1 ) 00178 * 00179 * b2 := b2 - V2*w 00180 * 00181 CALL DGEMV( 'NO TRANSPOSE', N-K-I+1, I-1, -ONE, 00182 $ A( K+I, 1 ), 00183 $ LDA, T( 1, NB ), 1, ONE, A( K+I, I ), 1 ) 00184 * 00185 * b1 := b1 - V1*w 00186 * 00187 CALL DTRMV( 'Lower', 'NO TRANSPOSE', 00188 $ 'UNIT', I-1, 00189 $ A( K+1, 1 ), LDA, T( 1, NB ), 1 ) 00190 CALL DAXPY( I-1, -ONE, T( 1, NB ), 1, A( K+1, I ), 1 ) 00191 * 00192 A( K+I-1, I-1 ) = EI 00193 END IF 00194 * 00195 * Generate the elementary reflector H(I) to annihilate 00196 * A(K+I+1:N,I) 00197 * 00198 CALL DLARFG( N-K-I+1, A( K+I, I ), A( MIN( K+I+1, N ), I ), 1, 00199 $ TAU( I ) ) 00200 EI = A( K+I, I ) 00201 A( K+I, I ) = ONE 00202 * 00203 * Compute Y(K+1:N,I) 00204 * 00205 CALL DGEMV( 'NO TRANSPOSE', N-K, N-K-I+1, 00206 $ ONE, A( K+1, I+1 ), 00207 $ LDA, A( K+I, I ), 1, ZERO, Y( K+1, I ), 1 ) 00208 CALL DGEMV( 'Transpose', N-K-I+1, I-1, 00209 $ ONE, A( K+I, 1 ), LDA, 00210 $ A( K+I, I ), 1, ZERO, T( 1, I ), 1 ) 00211 CALL DGEMV( 'NO TRANSPOSE', N-K, I-1, -ONE, 00212 $ Y( K+1, 1 ), LDY, 00213 $ T( 1, I ), 1, ONE, Y( K+1, I ), 1 ) 00214 CALL DSCAL( N-K, TAU( I ), Y( K+1, I ), 1 ) 00215 * 00216 * Compute T(1:I,I) 00217 * 00218 CALL DSCAL( I-1, -TAU( I ), T( 1, I ), 1 ) 00219 CALL DTRMV( 'Upper', 'No Transpose', 'NON-UNIT', 00220 $ I-1, T, LDT, 00221 $ T( 1, I ), 1 ) 00222 T( I, I ) = TAU( I ) 00223 * 00224 10 CONTINUE 00225 A( K+NB, NB ) = EI 00226 * 00227 * Compute Y(1:K,1:NB) 00228 * 00229 CALL DLACPY( 'ALL', K, NB, A( 1, 2 ), LDA, Y, LDY ) 00230 CALL DTRMM( 'RIGHT', 'Lower', 'NO TRANSPOSE', 00231 $ 'UNIT', K, NB, 00232 $ ONE, A( K+1, 1 ), LDA, Y, LDY ) 00233 IF( N.GT.K+NB ) 00234 $ CALL DGEMM( 'NO TRANSPOSE', 'NO TRANSPOSE', K, 00235 $ NB, N-K-NB, ONE, 00236 $ A( 1, 2+NB ), LDA, A( K+1+NB, 1 ), LDA, ONE, Y, 00237 $ LDY ) 00238 CALL DTRMM( 'RIGHT', 'Upper', 'NO TRANSPOSE', 00239 $ 'NON-UNIT', K, NB, 00240 $ ONE, T, LDT, Y, LDY ) 00241 * 00242 RETURN 00243 * 00244 * End of DLAHR2 00245 * 00246 END