1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361
|
*> \brief \b DGEHRD
*
* =========== DOCUMENTATION ===========
*
* Online html documentation available at
* http://www.netlib.org/lapack/explore-html/
*
*> \htmlonly
*> Download DGEHRD + dependencies
*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dgehrd.f">
*> [TGZ]</a>
*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dgehrd.f">
*> [ZIP]</a>
*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dgehrd.f">
*> [TXT]</a>
*> \endhtmlonly
*
* Definition:
* ===========
*
* SUBROUTINE DGEHRD( N, ILO, IHI, A, LDA, TAU, WORK, LWORK, INFO )
*
* .. Scalar Arguments ..
* INTEGER IHI, ILO, INFO, LDA, LWORK, N
* ..
* .. Array Arguments ..
* DOUBLE PRECISION A( LDA, * ), TAU( * ), WORK( * )
* ..
*
*
*> \par Purpose:
* =============
*>
*> \verbatim
*>
*> DGEHRD reduces a real general matrix A to upper Hessenberg form H by
*> an orthogonal similarity transformation: Q**T * A * Q = H .
*> \endverbatim
*
* Arguments:
* ==========
*
*> \param[in] N
*> \verbatim
*> N is INTEGER
*> The order of the matrix A. N >= 0.
*> \endverbatim
*>
*> \param[in] ILO
*> \verbatim
*> ILO is INTEGER
*> \endverbatim
*>
*> \param[in] IHI
*> \verbatim
*> IHI is INTEGER
*>
*> It is assumed that A is already upper triangular in rows
*> and columns 1:ILO-1 and IHI+1:N. ILO and IHI are normally
*> set by a previous call to DGEBAL; otherwise they should be
*> set to 1 and N respectively. See Further Details.
*> 1 <= ILO <= IHI <= N, if N > 0; ILO=1 and IHI=0, if N=0.
*> \endverbatim
*>
*> \param[in,out] A
*> \verbatim
*> A is DOUBLE PRECISION array, dimension (LDA,N)
*> On entry, the N-by-N general matrix to be reduced.
*> On exit, the upper triangle and the first subdiagonal of A
*> are overwritten with the upper Hessenberg matrix H, and the
*> elements below the first subdiagonal, with the array TAU,
*> represent the orthogonal matrix Q as a product of elementary
*> reflectors. See Further Details.
*> \endverbatim
*>
*> \param[in] LDA
*> \verbatim
*> LDA is INTEGER
*> The leading dimension of the array A. LDA >= max(1,N).
*> \endverbatim
*>
*> \param[out] TAU
*> \verbatim
*> TAU is DOUBLE PRECISION array, dimension (N-1)
*> The scalar factors of the elementary reflectors (see Further
*> Details). Elements 1:ILO-1 and IHI:N-1 of TAU are set to
*> zero.
*> \endverbatim
*>
*> \param[out] WORK
*> \verbatim
*> WORK is DOUBLE PRECISION array, dimension (MAX(1,LWORK))
*> On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
*> \endverbatim
*>
*> \param[in] LWORK
*> \verbatim
*> LWORK is INTEGER
*> The length of the array WORK. LWORK >= max(1,N).
*> For good performance, LWORK should generally be larger.
*>
*> If LWORK = -1, then a workspace query is assumed; the routine
*> only calculates the optimal size of the WORK array, returns
*> this value as the first entry of the WORK array, and no error
*> message related to LWORK is issued by XERBLA.
*> \endverbatim
*>
*> \param[out] INFO
*> \verbatim
*> INFO is INTEGER
*> = 0: successful exit
*> < 0: if INFO = -i, the i-th argument had an illegal value.
*> \endverbatim
*
* Authors:
* ========
*
*> \author Univ. of Tennessee
*> \author Univ. of California Berkeley
*> \author Univ. of Colorado Denver
*> \author NAG Ltd.
*
*> \ingroup gehrd
*
*> \par Further Details:
* =====================
*>
*> \verbatim
*>
*> The matrix Q is represented as a product of (ihi-ilo) elementary
*> reflectors
*>
*> Q = H(ilo) H(ilo+1) . . . H(ihi-1).
*>
*> Each H(i) has the form
*>
*> H(i) = I - tau * v * v**T
*>
*> where tau is a real scalar, and v is a real vector with
*> v(1:i) = 0, v(i+1) = 1 and v(ihi+1:n) = 0; v(i+2:ihi) is stored on
*> exit in A(i+2:ihi,i), and tau in TAU(i).
*>
*> The contents of A are illustrated by the following example, with
*> n = 7, ilo = 2 and ihi = 6:
*>
*> on entry, on exit,
*>
*> ( a a a a a a a ) ( a a h h h h a )
*> ( a a a a a a ) ( a h h h h a )
*> ( a a a a a a ) ( h h h h h h )
*> ( a a a a a a ) ( v2 h h h h h )
*> ( a a a a a a ) ( v2 v3 h h h h )
*> ( a a a a a a ) ( v2 v3 v4 h h h )
*> ( a ) ( a )
*>
*> where a denotes an element of the original matrix A, h denotes a
*> modified element of the upper Hessenberg matrix H, and vi denotes an
*> element of the vector defining H(i).
*>
*> This file is a slight modification of LAPACK-3.0's DGEHRD
*> subroutine incorporating improvements proposed by Quintana-Orti and
*> Van de Geijn (2006). (See DLAHR2.)
*> \endverbatim
*>
* =====================================================================
SUBROUTINE DGEHRD( N, ILO, IHI, A, LDA, TAU, WORK, LWORK,
$ INFO )
*
* -- LAPACK computational routine --
* -- LAPACK is a software package provided by Univ. of Tennessee, --
* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
*
* .. Scalar Arguments ..
INTEGER IHI, ILO, INFO, LDA, LWORK, N
* ..
* .. Array Arguments ..
DOUBLE PRECISION A( LDA, * ), TAU( * ), WORK( * )
* ..
*
* =====================================================================
*
* .. Parameters ..
INTEGER NBMAX, LDT, TSIZE
PARAMETER ( NBMAX = 64, LDT = NBMAX+1,
$ TSIZE = LDT*NBMAX )
DOUBLE PRECISION ZERO, ONE
PARAMETER ( ZERO = 0.0D+0,
$ ONE = 1.0D+0 )
* ..
* .. Local Scalars ..
LOGICAL LQUERY
INTEGER I, IB, IINFO, IWT, J, LDWORK, LWKOPT, NB,
$ NBMIN, NH, NX
DOUBLE PRECISION EI
* ..
* .. External Subroutines ..
EXTERNAL DAXPY, DGEHD2, DGEMM, DLAHR2, DLARFB,
$ DTRMM,
$ XERBLA
* ..
* .. Intrinsic Functions ..
INTRINSIC MAX, MIN
* ..
* .. External Functions ..
INTEGER ILAENV
EXTERNAL ILAENV
* ..
* .. Executable Statements ..
*
* Test the input parameters
*
INFO = 0
LQUERY = ( LWORK.EQ.-1 )
IF( N.LT.0 ) THEN
INFO = -1
ELSE IF( ILO.LT.1 .OR. ILO.GT.MAX( 1, N ) ) THEN
INFO = -2
ELSE IF( IHI.LT.MIN( ILO, N ) .OR. IHI.GT.N ) THEN
INFO = -3
ELSE IF( LDA.LT.MAX( 1, N ) ) THEN
INFO = -5
ELSE IF( LWORK.LT.MAX( 1, N ) .AND. .NOT.LQUERY ) THEN
INFO = -8
END IF
*
NH = IHI - ILO + 1
IF( INFO.EQ.0 ) THEN
*
* Compute the workspace requirements
*
IF( NH.LE.1 ) THEN
LWKOPT = 1
ELSE
NB = MIN( NBMAX, ILAENV( 1, 'DGEHRD', ' ', N, ILO, IHI,
$ -1 ) )
LWKOPT = N*NB + TSIZE
ENDIF
WORK( 1 ) = LWKOPT
END IF
*
IF( INFO.NE.0 ) THEN
CALL XERBLA( 'DGEHRD', -INFO )
RETURN
ELSE IF( LQUERY ) THEN
RETURN
END IF
*
* Set elements 1:ILO-1 and IHI:N-1 of TAU to zero
*
DO 10 I = 1, ILO - 1
TAU( I ) = ZERO
10 CONTINUE
DO 20 I = MAX( 1, IHI ), N - 1
TAU( I ) = ZERO
20 CONTINUE
*
* Quick return if possible
*
IF( NH.LE.1 ) THEN
WORK( 1 ) = 1
RETURN
END IF
*
* Determine the block size
*
NB = MIN( NBMAX, ILAENV( 1, 'DGEHRD', ' ', N, ILO, IHI, -1 ) )
NBMIN = 2
IF( NB.GT.1 .AND. NB.LT.NH ) THEN
*
* Determine when to cross over from blocked to unblocked code
* (last block is always handled by unblocked code)
*
NX = MAX( NB, ILAENV( 3, 'DGEHRD', ' ', N, ILO, IHI, -1 ) )
IF( NX.LT.NH ) THEN
*
* Determine if workspace is large enough for blocked code
*
IF( LWORK.LT.LWKOPT ) THEN
*
* Not enough workspace to use optimal NB: determine the
* minimum value of NB, and reduce NB or force use of
* unblocked code
*
NBMIN = MAX( 2, ILAENV( 2, 'DGEHRD', ' ', N, ILO, IHI,
$ -1 ) )
IF( LWORK.GE.(N*NBMIN + TSIZE) ) THEN
NB = (LWORK-TSIZE) / N
ELSE
NB = 1
END IF
END IF
END IF
END IF
LDWORK = N
*
IF( NB.LT.NBMIN .OR. NB.GE.NH ) THEN
*
* Use unblocked code below
*
I = ILO
*
ELSE
*
* Use blocked code
*
IWT = 1 + N*NB
DO 40 I = ILO, IHI - 1 - NX, NB
IB = MIN( NB, IHI-I )
*
* Reduce columns i:i+ib-1 to Hessenberg form, returning the
* matrices V and T of the block reflector H = I - V*T*V**T
* which performs the reduction, and also the matrix Y = A*V*T
*
CALL DLAHR2( IHI, I, IB, A( 1, I ), LDA, TAU( I ),
$ WORK( IWT ), LDT, WORK, LDWORK )
*
* Apply the block reflector H to A(1:ihi,i+ib:ihi) from the
* right, computing A := A - Y * V**T. V(i+ib,ib-1) must be set
* to 1
*
EI = A( I+IB, I+IB-1 )
A( I+IB, I+IB-1 ) = ONE
CALL DGEMM( 'No transpose', 'Transpose',
$ IHI, IHI-I-IB+1,
$ IB, -ONE, WORK, LDWORK, A( I+IB, I ), LDA, ONE,
$ A( 1, I+IB ), LDA )
A( I+IB, I+IB-1 ) = EI
*
* Apply the block reflector H to A(1:i,i+1:i+ib-1) from the
* right
*
CALL DTRMM( 'Right', 'Lower', 'Transpose',
$ 'Unit', I, IB-1,
$ ONE, A( I+1, I ), LDA, WORK, LDWORK )
DO 30 J = 0, IB-2
CALL DAXPY( I, -ONE, WORK( LDWORK*J+1 ), 1,
$ A( 1, I+J+1 ), 1 )
30 CONTINUE
*
* Apply the block reflector H to A(i+1:ihi,i+ib:n) from the
* left
*
CALL DLARFB( 'Left', 'Transpose', 'Forward',
$ 'Columnwise',
$ IHI-I, N-I-IB+1, IB, A( I+1, I ), LDA,
$ WORK( IWT ), LDT, A( I+1, I+IB ), LDA,
$ WORK, LDWORK )
40 CONTINUE
END IF
*
* Use unblocked code to reduce the rest of the matrix
*
CALL DGEHD2( N, I, IHI, A, LDA, TAU, WORK, IINFO )
*
WORK( 1 ) = LWKOPT
*
RETURN
*
* End of DGEHRD
*
END
|