KTH framework for Nek5000 toolboxes; testing version  0.0.1
dgehrd.f
Go to the documentation of this file.
1  SUBROUTINE dgehrd( N, ILO, IHI, A, LDA, TAU, WORK, LWORK, INFO )
2 *
3 * -- LAPACK routine (version 3.0) --
4 * Univ. of Tennessee, Univ. of California Berkeley, NAG Ltd.,
5 * Courant Institute, Argonne National Lab, and Rice University
6 * June 30, 1999
7 *
8 * .. Scalar Arguments ..
9  INTEGER IHI, ILO, INFO, LDA, LWORK, N
10 * ..
11 * .. Array Arguments ..
12  DOUBLE PRECISION A( LDA, * ), TAU( * ), WORK( * )
13 * ..
14 *
15 * Purpose
16 * =======
17 *
18 * DGEHRD reduces a real general matrix A to upper Hessenberg form H by
19 * an orthogonal similarity transformation: Q' * A * Q = H .
20 *
21 * Arguments
22 * =========
23 *
24 * N (input) INTEGER
25 * The order of the matrix A. N >= 0.
26 *
27 * ILO (input) INTEGER
28 * IHI (input) INTEGER
29 * It is assumed that A is already upper triangular in rows
30 * and columns 1:ILO-1 and IHI+1:N. ILO and IHI are normally
31 * set by a previous call to DGEBAL; otherwise they should be
32 * set to 1 and N respectively. See Further Details.
33 * 1 <= ILO <= IHI <= N, if N > 0; ILO=1 and IHI=0, if N=0.
34 *
35 * A (input/output) DOUBLE PRECISION array, dimension (LDA,N)
36 * On entry, the N-by-N general matrix to be reduced.
37 * On exit, the upper triangle and the first subdiagonal of A
38 * are overwritten with the upper Hessenberg matrix H, and the
39 * elements below the first subdiagonal, with the array TAU,
40 * represent the orthogonal matrix Q as a product of elementary
41 * reflectors. See Further Details.
42 *
43 * LDA (input) INTEGER
44 * The leading dimension of the array A. LDA >= max(1,N).
45 *
46 * TAU (output) DOUBLE PRECISION array, dimension (N-1)
47 * The scalar factors of the elementary reflectors (see Further
48 * Details). Elements 1:ILO-1 and IHI:N-1 of TAU are set to
49 * zero.
50 *
51 * WORK (workspace/output) DOUBLE PRECISION array, dimension (LWORK)
52 * On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
53 *
54 * LWORK (input) INTEGER
55 * The length of the array WORK. LWORK >= max(1,N).
56 * For optimum performance LWORK >= N*NB, where NB is the
57 * optimal blocksize.
58 *
59 * If LWORK = -1, then a workspace query is assumed; the routine
60 * only calculates the optimal size of the WORK array, returns
61 * this value as the first entry of the WORK array, and no error
62 * message related to LWORK is issued by XERBLA.
63 *
64 * INFO (output) INTEGER
65 * = 0: successful exit
66 * < 0: if INFO = -i, the i-th argument had an illegal value.
67 *
68 * Further Details
69 * ===============
70 *
71 * The matrix Q is represented as a product of (ihi-ilo) elementary
72 * reflectors
73 *
74 * Q = H(ilo) H(ilo+1) . . . H(ihi-1).
75 *
76 * Each H(i) has the form
77 *
78 * H(i) = I - tau * v * v'
79 *
80 * where tau is a real scalar, and v is a real vector with
81 * v(1:i) = 0, v(i+1) = 1 and v(ihi+1:n) = 0; v(i+2:ihi) is stored on
82 * exit in A(i+2:ihi,i), and tau in TAU(i).
83 *
84 * The contents of A are illustrated by the following example, with
85 * n = 7, ilo = 2 and ihi = 6:
86 *
87 * on entry, on exit,
88 *
89 * ( a a a a a a a ) ( a a h h h h a )
90 * ( a a a a a a ) ( a h h h h a )
91 * ( a a a a a a ) ( h h h h h h )
92 * ( a a a a a a ) ( v2 h h h h h )
93 * ( a a a a a a ) ( v2 v3 h h h h )
94 * ( a a a a a a ) ( v2 v3 v4 h h h )
95 * ( a ) ( a )
96 *
97 * where a denotes an element of the original matrix A, h denotes a
98 * modified element of the upper Hessenberg matrix H, and vi denotes an
99 * element of the vector defining H(i).
100 *
101 * =====================================================================
102 *
103 * .. Parameters ..
104  INTEGER NBMAX, LDT
105  parameter( nbmax = 64, ldt = nbmax+1 )
106  DOUBLE PRECISION ZERO, ONE
107  parameter( zero = 0.0d+0, one = 1.0d+0 )
108 * ..
109 * .. Local Scalars ..
110  LOGICAL LQUERY
111  INTEGER I, IB, IINFO, IWS, LDWORK, LWKOPT, NB, NBMIN,
112  $ NH, NX
113  DOUBLE PRECISION EI
114 * ..
115 * .. Local Arrays ..
116  DOUBLE PRECISION T( LDT, NBMAX )
117 * ..
118 * .. External Subroutines ..
119  EXTERNAL dgehd2, dgemm, dlahrd, dlarfb, xerbla
120 * ..
121 * .. Intrinsic Functions ..
122  INTRINSIC max, min
123 * ..
124 * .. External Functions ..
125  INTEGER ILAENV
126  EXTERNAL ilaenv
127 * ..
128 * .. Executable Statements ..
129 *
130 * Test the input parameters
131 *
132  info = 0
133  nb = min( nbmax, ilaenv( 1, 'DGEHRD', ' ', n, ilo, ihi, -1 ) )
134  lwkopt = n*nb
135  work( 1 ) = lwkopt
136  lquery = ( lwork.EQ.-1 )
137  IF( n.LT.0 ) THEN
138  info = -1
139  ELSE IF( ilo.LT.1 .OR. ilo.GT.max( 1, n ) ) THEN
140  info = -2
141  ELSE IF( ihi.LT.min( ilo, n ) .OR. ihi.GT.n ) THEN
142  info = -3
143  ELSE IF( lda.LT.max( 1, n ) ) THEN
144  info = -5
145  ELSE IF( lwork.LT.max( 1, n ) .AND. .NOT.lquery ) THEN
146  info = -8
147  END IF
148  IF( info.NE.0 ) THEN
149  CALL xerbla( 'DGEHRD', -info )
150  RETURN
151  ELSE IF( lquery ) THEN
152  RETURN
153  END IF
154 *
155 * Set elements 1:ILO-1 and IHI:N-1 of TAU to zero
156 *
157  DO 10 i = 1, ilo - 1
158  tau( i ) = zero
159  10 CONTINUE
160  DO 20 i = max( 1, ihi ), n - 1
161  tau( i ) = zero
162  20 CONTINUE
163 *
164 * Quick return if possible
165 *
166  nh = ihi - ilo + 1
167  IF( nh.LE.1 ) THEN
168  work( 1 ) = 1
169  RETURN
170  END IF
171 *
172 * Determine the block size.
173 *
174  nb = min( nbmax, ilaenv( 1, 'DGEHRD', ' ', n, ilo, ihi, -1 ) )
175  nbmin = 2
176  iws = 1
177  IF( nb.GT.1 .AND. nb.LT.nh ) THEN
178 *
179 * Determine when to cross over from blocked to unblocked code
180 * (last block is always handled by unblocked code).
181 *
182  nx = max( nb, ilaenv( 3, 'DGEHRD', ' ', n, ilo, ihi, -1 ) )
183  IF( nx.LT.nh ) THEN
184 *
185 * Determine if workspace is large enough for blocked code.
186 *
187  iws = n*nb
188  IF( lwork.LT.iws ) THEN
189 *
190 * Not enough workspace to use optimal NB: determine the
191 * minimum value of NB, and reduce NB or force use of
192 * unblocked code.
193 *
194  nbmin = max( 2, ilaenv( 2, 'DGEHRD', ' ', n, ilo, ihi,
195  $ -1 ) )
196  IF( lwork.GE.n*nbmin ) THEN
197  nb = lwork / n
198  ELSE
199  nb = 1
200  END IF
201  END IF
202  END IF
203  END IF
204  ldwork = n
205 *
206  IF( nb.LT.nbmin .OR. nb.GE.nh ) THEN
207 *
208 * Use unblocked code below
209 *
210  i = ilo
211 *
212  ELSE
213 *
214 * Use blocked code
215 *
216  DO 30 i = ilo, ihi - 1 - nx, nb
217  ib = min( nb, ihi-i )
218 *
219 * Reduce columns i:i+ib-1 to Hessenberg form, returning the
220 * matrices V and T of the block reflector H = I - V*T*V'
221 * which performs the reduction, and also the matrix Y = A*V*T
222 *
223  CALL dlahrd( ihi, i, ib, a( 1, i ), lda, tau( i ), t, ldt,
224  $ work, ldwork )
225 *
226 * Apply the block reflector H to A(1:ihi,i+ib:ihi) from the
227 * right, computing A := A - Y * V'. V(i+ib,ib-1) must be set
228 * to 1.
229 *
230  ei = a( i+ib, i+ib-1 )
231  a( i+ib, i+ib-1 ) = one
232  CALL dgemm( 'No transpose', 'Transpose', ihi, ihi-i-ib+1,
233  $ ib, -one, work, ldwork, a( i+ib, i ), lda, one,
234  $ a( 1, i+ib ), lda )
235  a( i+ib, i+ib-1 ) = ei
236 *
237 * Apply the block reflector H to A(i+1:ihi,i+ib:n) from the
238 * left
239 *
240  CALL dlarfb( 'Left', 'Transpose', 'Forward', 'Columnwise',
241  $ ihi-i, n-i-ib+1, ib, a( i+1, i ), lda, t, ldt,
242  $ a( i+1, i+ib ), lda, work, ldwork )
243  30 CONTINUE
244  END IF
245 *
246 * Use unblocked code to reduce the rest of the matrix
247 *
248  CALL dgehd2( n, i, ihi, a, lda, tau, work, iinfo )
249  work( 1 ) = iws
250 *
251  RETURN
252 *
253 * End of DGEHRD
254 *
255  END
subroutine dgehd2(N, ILO, IHI, A, LDA, TAU, WORK, INFO)
Definition: dgehd2.f:2
subroutine dgehrd(N, ILO, IHI, A, LDA, TAU, WORK, LWORK, INFO)
Definition: dgehrd.f:2
subroutine dgemm(TRANSA, TRANSB, M, N, K, ALPHA, A, LDA, B, LDB, BETA, C, LDC)
Definition: dgemm.f:3
subroutine dlahrd(N, K, NB, A, LDA, TAU, T, LDT, Y, LDY)
Definition: dlahrd.f:2
subroutine dlarfb(SIDE, TRANS, DIRECT, STOREV, M, N, K, V, LDV, T, LDT, C, LDC, WORK, LDWORK)
Definition: dlarfb.f:3
subroutine xerbla(SRNAME, INFO)
Definition: xerbla.f:2