LAPACK 3.12.1
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches

◆ dgehrd()

subroutine dgehrd ( integer n,
integer ilo,
integer ihi,
double precision, dimension( lda, * ) a,
integer lda,
double precision, dimension( * ) tau,
double precision, dimension( * ) work,
integer lwork,
integer info )

DGEHRD

Download DGEHRD + dependencies [TGZ] [ZIP] [TXT]

Purpose:
!>
!> DGEHRD reduces a real general matrix A to upper Hessenberg form H by
!> an orthogonal similarity transformation:  Q**T * A * Q = H .
!> 
Parameters
[in]N
!>          N is INTEGER
!>          The order of the matrix A.  N >= 0.
!> 
[in]ILO
!>          ILO is INTEGER
!> 
[in]IHI
!>          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.
!> 
[in,out]A
!>          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.
!> 
[in]LDA
!>          LDA is INTEGER
!>          The leading dimension of the array A.  LDA >= max(1,N).
!> 
[out]TAU
!>          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.
!> 
[out]WORK
!>          WORK is DOUBLE PRECISION array, dimension (MAX(1,LWORK))
!>          On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
!> 
[in]LWORK
!>          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.
!> 
[out]INFO
!>          INFO is INTEGER
!>          = 0:  successful exit
!>          < 0:  if INFO = -i, the i-th argument had an illegal value.
!> 
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Further Details:
!>
!>  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.)
!> 

Definition at line 164 of file dgehrd.f.

166*
167* -- LAPACK computational routine --
168* -- LAPACK is a software package provided by Univ. of Tennessee, --
169* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
170*
171* .. Scalar Arguments ..
172 INTEGER IHI, ILO, INFO, LDA, LWORK, N
173* ..
174* .. Array Arguments ..
175 DOUBLE PRECISION A( LDA, * ), TAU( * ), WORK( * )
176* ..
177*
178* =====================================================================
179*
180* .. Parameters ..
181 INTEGER NBMAX, LDT, TSIZE
182 parameter( nbmax = 64, ldt = nbmax+1,
183 $ tsize = ldt*nbmax )
184 DOUBLE PRECISION ZERO, ONE
185 parameter( zero = 0.0d+0,
186 $ one = 1.0d+0 )
187* ..
188* .. Local Scalars ..
189 LOGICAL LQUERY
190 INTEGER I, IB, IINFO, IWT, J, LDWORK, LWKOPT, NB,
191 $ NBMIN, NH, NX
192 DOUBLE PRECISION EI
193* ..
194* .. External Subroutines ..
195 EXTERNAL daxpy, dgehd2, dgemm, dlahr2, dlarfb,
196 $ dtrmm,
197 $ xerbla
198* ..
199* .. Intrinsic Functions ..
200 INTRINSIC max, min
201* ..
202* .. External Functions ..
203 INTEGER ILAENV
204 EXTERNAL ilaenv
205* ..
206* .. Executable Statements ..
207*
208* Test the input parameters
209*
210 info = 0
211 lquery = ( lwork.EQ.-1 )
212 IF( n.LT.0 ) THEN
213 info = -1
214 ELSE IF( ilo.LT.1 .OR. ilo.GT.max( 1, n ) ) THEN
215 info = -2
216 ELSE IF( ihi.LT.min( ilo, n ) .OR. ihi.GT.n ) THEN
217 info = -3
218 ELSE IF( lda.LT.max( 1, n ) ) THEN
219 info = -5
220 ELSE IF( lwork.LT.max( 1, n ) .AND. .NOT.lquery ) THEN
221 info = -8
222 END IF
223*
224 nh = ihi - ilo + 1
225 IF( info.EQ.0 ) THEN
226*
227* Compute the workspace requirements
228*
229 IF( nh.LE.1 ) THEN
230 lwkopt = 1
231 ELSE
232 nb = min( nbmax, ilaenv( 1, 'DGEHRD', ' ', n, ilo, ihi,
233 $ -1 ) )
234 lwkopt = n*nb + tsize
235 ENDIF
236 work( 1 ) = lwkopt
237 END IF
238*
239 IF( info.NE.0 ) THEN
240 CALL xerbla( 'DGEHRD', -info )
241 RETURN
242 ELSE IF( lquery ) THEN
243 RETURN
244 END IF
245*
246* Set elements 1:ILO-1 and IHI:N-1 of TAU to zero
247*
248 DO 10 i = 1, ilo - 1
249 tau( i ) = zero
250 10 CONTINUE
251 DO 20 i = max( 1, ihi ), n - 1
252 tau( i ) = zero
253 20 CONTINUE
254*
255* Quick return if possible
256*
257 IF( nh.LE.1 ) THEN
258 work( 1 ) = 1
259 RETURN
260 END IF
261*
262* Determine the block size
263*
264 nb = min( nbmax, ilaenv( 1, 'DGEHRD', ' ', n, ilo, ihi, -1 ) )
265 nbmin = 2
266 IF( nb.GT.1 .AND. nb.LT.nh ) THEN
267*
268* Determine when to cross over from blocked to unblocked code
269* (last block is always handled by unblocked code)
270*
271 nx = max( nb, ilaenv( 3, 'DGEHRD', ' ', n, ilo, ihi, -1 ) )
272 IF( nx.LT.nh ) THEN
273*
274* Determine if workspace is large enough for blocked code
275*
276 IF( lwork.LT.lwkopt ) THEN
277*
278* Not enough workspace to use optimal NB: determine the
279* minimum value of NB, and reduce NB or force use of
280* unblocked code
281*
282 nbmin = max( 2, ilaenv( 2, 'DGEHRD', ' ', n, ilo, ihi,
283 $ -1 ) )
284 IF( lwork.GE.(n*nbmin + tsize) ) THEN
285 nb = (lwork-tsize) / n
286 ELSE
287 nb = 1
288 END IF
289 END IF
290 END IF
291 END IF
292 ldwork = n
293*
294 IF( nb.LT.nbmin .OR. nb.GE.nh ) THEN
295*
296* Use unblocked code below
297*
298 i = ilo
299*
300 ELSE
301*
302* Use blocked code
303*
304 iwt = 1 + n*nb
305 DO 40 i = ilo, ihi - 1 - nx, nb
306 ib = min( nb, ihi-i )
307*
308* Reduce columns i:i+ib-1 to Hessenberg form, returning the
309* matrices V and T of the block reflector H = I - V*T*V**T
310* which performs the reduction, and also the matrix Y = A*V*T
311*
312 CALL dlahr2( ihi, i, ib, a( 1, i ), lda, tau( i ),
313 $ work( iwt ), ldt, work, ldwork )
314*
315* Apply the block reflector H to A(1:ihi,i+ib:ihi) from the
316* right, computing A := A - Y * V**T. V(i+ib,ib-1) must be set
317* to 1
318*
319 ei = a( i+ib, i+ib-1 )
320 a( i+ib, i+ib-1 ) = one
321 CALL dgemm( 'No transpose', 'Transpose',
322 $ ihi, ihi-i-ib+1,
323 $ ib, -one, work, ldwork, a( i+ib, i ), lda, one,
324 $ a( 1, i+ib ), lda )
325 a( i+ib, i+ib-1 ) = ei
326*
327* Apply the block reflector H to A(1:i,i+1:i+ib-1) from the
328* right
329*
330 CALL dtrmm( 'Right', 'Lower', 'Transpose',
331 $ 'Unit', i, ib-1,
332 $ one, a( i+1, i ), lda, work, ldwork )
333 DO 30 j = 0, ib-2
334 CALL daxpy( i, -one, work( ldwork*j+1 ), 1,
335 $ a( 1, i+j+1 ), 1 )
336 30 CONTINUE
337*
338* Apply the block reflector H to A(i+1:ihi,i+ib:n) from the
339* left
340*
341 CALL dlarfb( 'Left', 'Transpose', 'Forward',
342 $ 'Columnwise',
343 $ ihi-i, n-i-ib+1, ib, a( i+1, i ), lda,
344 $ work( iwt ), ldt, a( i+1, i+ib ), lda,
345 $ work, ldwork )
346 40 CONTINUE
347 END IF
348*
349* Use unblocked code to reduce the rest of the matrix
350*
351 CALL dgehd2( n, i, ihi, a, lda, tau, work, iinfo )
352*
353 work( 1 ) = lwkopt
354*
355 RETURN
356*
357* End of DGEHRD
358*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine daxpy(n, da, dx, incx, dy, incy)
DAXPY
Definition daxpy.f:89
subroutine dgehd2(n, ilo, ihi, a, lda, tau, work, info)
DGEHD2 reduces a general square matrix to upper Hessenberg form using an unblocked algorithm.
Definition dgehd2.f:147
subroutine dgemm(transa, transb, m, n, k, alpha, a, lda, b, ldb, beta, c, ldc)
DGEMM
Definition dgemm.f:188
integer function ilaenv(ispec, name, opts, n1, n2, n3, n4)
ILAENV
Definition ilaenv.f:160
subroutine dlahr2(n, k, nb, a, lda, tau, t, ldt, y, ldy)
DLAHR2 reduces the specified number of first columns of a general rectangular matrix A so that elemen...
Definition dlahr2.f:179
subroutine dlarfb(side, trans, direct, storev, m, n, k, v, ldv, t, ldt, c, ldc, work, ldwork)
DLARFB applies a block reflector or its transpose to a general rectangular matrix.
Definition dlarfb.f:195
subroutine dtrmm(side, uplo, transa, diag, m, n, alpha, a, lda, b, ldb)
DTRMM
Definition dtrmm.f:177
Here is the call graph for this function:
Here is the caller graph for this function: