LAPACK 3.12.1
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
dgghrd.f
Go to the documentation of this file.
1*> \brief \b DGGHRD
2*
3* =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6* http://www.netlib.org/lapack/explore-html/
7*
8*> Download DGGHRD + dependencies
9*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dgghrd.f">
10*> [TGZ]</a>
11*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dgghrd.f">
12*> [ZIP]</a>
13*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dgghrd.f">
14*> [TXT]</a>
15*
16* Definition:
17* ===========
18*
19* SUBROUTINE DGGHRD( COMPQ, COMPZ, N, ILO, IHI, A, LDA, B, LDB, Q,
20* LDQ, Z, LDZ, INFO )
21*
22* .. Scalar Arguments ..
23* CHARACTER COMPQ, COMPZ
24* INTEGER IHI, ILO, INFO, LDA, LDB, LDQ, LDZ, N
25* ..
26* .. Array Arguments ..
27* DOUBLE PRECISION A( LDA, * ), B( LDB, * ), Q( LDQ, * ),
28* $ Z( LDZ, * )
29* ..
30*
31*
32*> \par Purpose:
33* =============
34*>
35*> \verbatim
36*>
37*> DGGHRD reduces a pair of real matrices (A,B) to generalized upper
38*> Hessenberg form using orthogonal transformations, where A is a
39*> general matrix and B is upper triangular. The form of the
40*> generalized eigenvalue problem is
41*> A*x = lambda*B*x,
42*> and B is typically made upper triangular by computing its QR
43*> factorization and moving the orthogonal matrix Q to the left side
44*> of the equation.
45*>
46*> This subroutine simultaneously reduces A to a Hessenberg matrix H:
47*> Q**T*A*Z = H
48*> and transforms B to another upper triangular matrix T:
49*> Q**T*B*Z = T
50*> in order to reduce the problem to its standard form
51*> H*y = lambda*T*y
52*> where y = Z**T*x.
53*>
54*> The orthogonal matrices Q and Z are determined as products of Givens
55*> rotations. They may either be formed explicitly, or they may be
56*> postmultiplied into input matrices Q1 and Z1, so that
57*>
58*> Q1 * A * Z1**T = (Q1*Q) * H * (Z1*Z)**T
59*>
60*> Q1 * B * Z1**T = (Q1*Q) * T * (Z1*Z)**T
61*>
62*> If Q1 is the orthogonal matrix from the QR factorization of B in the
63*> original equation A*x = lambda*B*x, then DGGHRD reduces the original
64*> problem to generalized Hessenberg form.
65*> \endverbatim
66*
67* Arguments:
68* ==========
69*
70*> \param[in] COMPQ
71*> \verbatim
72*> COMPQ is CHARACTER*1
73*> = 'N': do not compute Q;
74*> = 'I': Q is initialized to the unit matrix, and the
75*> orthogonal matrix Q is returned;
76*> = 'V': Q must contain an orthogonal matrix Q1 on entry,
77*> and the product Q1*Q is returned.
78*> \endverbatim
79*>
80*> \param[in] COMPZ
81*> \verbatim
82*> COMPZ is CHARACTER*1
83*> = 'N': do not compute Z;
84*> = 'I': Z is initialized to the unit matrix, and the
85*> orthogonal matrix Z is returned;
86*> = 'V': Z must contain an orthogonal matrix Z1 on entry,
87*> and the product Z1*Z is returned.
88*> \endverbatim
89*>
90*> \param[in] N
91*> \verbatim
92*> N is INTEGER
93*> The order of the matrices A and B. N >= 0.
94*> \endverbatim
95*>
96*> \param[in] ILO
97*> \verbatim
98*> ILO is INTEGER
99*> \endverbatim
100*>
101*> \param[in] IHI
102*> \verbatim
103*> IHI is INTEGER
104*>
105*> ILO and IHI mark the rows and columns of A which are to be
106*> reduced. It is assumed that A is already upper triangular
107*> in rows and columns 1:ILO-1 and IHI+1:N. ILO and IHI are
108*> normally set by a previous call to DGGBAL; otherwise they
109*> should be set to 1 and N respectively.
110*> 1 <= ILO <= IHI <= N, if N > 0; ILO=1 and IHI=0, if N=0.
111*> \endverbatim
112*>
113*> \param[in,out] A
114*> \verbatim
115*> A is DOUBLE PRECISION array, dimension (LDA, N)
116*> On entry, the N-by-N general matrix to be reduced.
117*> On exit, the upper triangle and the first subdiagonal of A
118*> are overwritten with the upper Hessenberg matrix H, and the
119*> rest is set to zero.
120*> \endverbatim
121*>
122*> \param[in] LDA
123*> \verbatim
124*> LDA is INTEGER
125*> The leading dimension of the array A. LDA >= max(1,N).
126*> \endverbatim
127*>
128*> \param[in,out] B
129*> \verbatim
130*> B is DOUBLE PRECISION array, dimension (LDB, N)
131*> On entry, the N-by-N upper triangular matrix B.
132*> On exit, the upper triangular matrix T = Q**T B Z. The
133*> elements below the diagonal are set to zero.
134*> \endverbatim
135*>
136*> \param[in] LDB
137*> \verbatim
138*> LDB is INTEGER
139*> The leading dimension of the array B. LDB >= max(1,N).
140*> \endverbatim
141*>
142*> \param[in,out] Q
143*> \verbatim
144*> Q is DOUBLE PRECISION array, dimension (LDQ, N)
145*> On entry, if COMPQ = 'V', the orthogonal matrix Q1,
146*> typically from the QR factorization of B.
147*> On exit, if COMPQ='I', the orthogonal matrix Q, and if
148*> COMPQ = 'V', the product Q1*Q.
149*> Not referenced if COMPQ='N'.
150*> \endverbatim
151*>
152*> \param[in] LDQ
153*> \verbatim
154*> LDQ is INTEGER
155*> The leading dimension of the array Q.
156*> LDQ >= N if COMPQ='V' or 'I'; LDQ >= 1 otherwise.
157*> \endverbatim
158*>
159*> \param[in,out] Z
160*> \verbatim
161*> Z is DOUBLE PRECISION array, dimension (LDZ, N)
162*> On entry, if COMPZ = 'V', the orthogonal matrix Z1.
163*> On exit, if COMPZ='I', the orthogonal matrix Z, and if
164*> COMPZ = 'V', the product Z1*Z.
165*> Not referenced if COMPZ='N'.
166*> \endverbatim
167*>
168*> \param[in] LDZ
169*> \verbatim
170*> LDZ is INTEGER
171*> The leading dimension of the array Z.
172*> LDZ >= N if COMPZ='V' or 'I'; LDZ >= 1 otherwise.
173*> \endverbatim
174*>
175*> \param[out] INFO
176*> \verbatim
177*> INFO is INTEGER
178*> = 0: successful exit.
179*> < 0: if INFO = -i, the i-th argument had an illegal value.
180*> \endverbatim
181*
182* Authors:
183* ========
184*
185*> \author Univ. of Tennessee
186*> \author Univ. of California Berkeley
187*> \author Univ. of Colorado Denver
188*> \author NAG Ltd.
189*
190*> \ingroup gghrd
191*
192*> \par Further Details:
193* =====================
194*>
195*> \verbatim
196*>
197*> This routine reduces A to Hessenberg and B to triangular form by
198*> an unblocked reduction, as described in _Matrix_Computations_,
199*> by Golub and Van Loan (Johns Hopkins Press.)
200*> \endverbatim
201*>
202* =====================================================================
203 SUBROUTINE dgghrd( COMPQ, COMPZ, N, ILO, IHI, A, LDA, B, LDB,
204 $ Q,
205 $ LDQ, Z, LDZ, INFO )
206*
207* -- LAPACK computational routine --
208* -- LAPACK is a software package provided by Univ. of Tennessee, --
209* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
210*
211* .. Scalar Arguments ..
212 CHARACTER COMPQ, COMPZ
213 INTEGER IHI, ILO, INFO, LDA, LDB, LDQ, LDZ, N
214* ..
215* .. Array Arguments ..
216 DOUBLE PRECISION A( LDA, * ), B( LDB, * ), Q( LDQ, * ),
217 $ Z( LDZ, * )
218* ..
219*
220* =====================================================================
221*
222* .. Parameters ..
223 DOUBLE PRECISION ONE, ZERO
224 PARAMETER ( ONE = 1.0d+0, zero = 0.0d+0 )
225* ..
226* .. Local Scalars ..
227 LOGICAL ILQ, ILZ
228 INTEGER ICOMPQ, ICOMPZ, JCOL, JROW
229 DOUBLE PRECISION C, S, TEMP
230* ..
231* .. External Functions ..
232 LOGICAL LSAME
233 EXTERNAL LSAME
234* ..
235* .. External Subroutines ..
236 EXTERNAL dlartg, dlaset, drot, xerbla
237* ..
238* .. Intrinsic Functions ..
239 INTRINSIC max
240* ..
241* .. Executable Statements ..
242*
243* Decode COMPQ
244*
245 IF( lsame( compq, 'N' ) ) THEN
246 ilq = .false.
247 icompq = 1
248 ELSE IF( lsame( compq, 'V' ) ) THEN
249 ilq = .true.
250 icompq = 2
251 ELSE IF( lsame( compq, 'I' ) ) THEN
252 ilq = .true.
253 icompq = 3
254 ELSE
255 icompq = 0
256 END IF
257*
258* Decode COMPZ
259*
260 IF( lsame( compz, 'N' ) ) THEN
261 ilz = .false.
262 icompz = 1
263 ELSE IF( lsame( compz, 'V' ) ) THEN
264 ilz = .true.
265 icompz = 2
266 ELSE IF( lsame( compz, 'I' ) ) THEN
267 ilz = .true.
268 icompz = 3
269 ELSE
270 icompz = 0
271 END IF
272*
273* Test the input parameters.
274*
275 info = 0
276 IF( icompq.LE.0 ) THEN
277 info = -1
278 ELSE IF( icompz.LE.0 ) THEN
279 info = -2
280 ELSE IF( n.LT.0 ) THEN
281 info = -3
282 ELSE IF( ilo.LT.1 ) THEN
283 info = -4
284 ELSE IF( ihi.GT.n .OR. ihi.LT.ilo-1 ) THEN
285 info = -5
286 ELSE IF( lda.LT.max( 1, n ) ) THEN
287 info = -7
288 ELSE IF( ldb.LT.max( 1, n ) ) THEN
289 info = -9
290 ELSE IF( ( ilq .AND. ldq.LT.n ) .OR. ldq.LT.1 ) THEN
291 info = -11
292 ELSE IF( ( ilz .AND. ldz.LT.n ) .OR. ldz.LT.1 ) THEN
293 info = -13
294 END IF
295 IF( info.NE.0 ) THEN
296 CALL xerbla( 'DGGHRD', -info )
297 RETURN
298 END IF
299*
300* Initialize Q and Z if desired.
301*
302 IF( icompq.EQ.3 )
303 $ CALL dlaset( 'Full', n, n, zero, one, q, ldq )
304 IF( icompz.EQ.3 )
305 $ CALL dlaset( 'Full', n, n, zero, one, z, ldz )
306*
307* Quick return if possible
308*
309 IF( n.LE.1 )
310 $ RETURN
311*
312* Zero out lower triangle of B
313*
314 DO 20 jcol = 1, n - 1
315 DO 10 jrow = jcol + 1, n
316 b( jrow, jcol ) = zero
317 10 CONTINUE
318 20 CONTINUE
319*
320* Reduce A and B
321*
322 DO 40 jcol = ilo, ihi - 2
323*
324 DO 30 jrow = ihi, jcol + 2, -1
325*
326* Step 1: rotate rows JROW-1, JROW to kill A(JROW,JCOL)
327*
328 temp = a( jrow-1, jcol )
329 CALL dlartg( temp, a( jrow, jcol ), c, s,
330 $ a( jrow-1, jcol ) )
331 a( jrow, jcol ) = zero
332 CALL drot( n-jcol, a( jrow-1, jcol+1 ), lda,
333 $ a( jrow, jcol+1 ), lda, c, s )
334 CALL drot( n+2-jrow, b( jrow-1, jrow-1 ), ldb,
335 $ b( jrow, jrow-1 ), ldb, c, s )
336 IF( ilq )
337 $ CALL drot( n, q( 1, jrow-1 ), 1, q( 1, jrow ), 1, c,
338 $ s )
339*
340* Step 2: rotate columns JROW, JROW-1 to kill B(JROW,JROW-1)
341*
342 temp = b( jrow, jrow )
343 CALL dlartg( temp, b( jrow, jrow-1 ), c, s,
344 $ b( jrow, jrow ) )
345 b( jrow, jrow-1 ) = zero
346 CALL drot( ihi, a( 1, jrow ), 1, a( 1, jrow-1 ), 1, c,
347 $ s )
348 CALL drot( jrow-1, b( 1, jrow ), 1, b( 1, jrow-1 ), 1, c,
349 $ s )
350 IF( ilz )
351 $ CALL drot( n, z( 1, jrow ), 1, z( 1, jrow-1 ), 1, c,
352 $ s )
353 30 CONTINUE
354 40 CONTINUE
355*
356 RETURN
357*
358* End of DGGHRD
359*
360 END
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine dgghrd(compq, compz, n, ilo, ihi, a, lda, b, ldb, q, ldq, z, ldz, info)
DGGHRD
Definition dgghrd.f:206
subroutine dlartg(f, g, c, s, r)
DLARTG generates a plane rotation with real cosine and real sine.
Definition dlartg.f90:111
subroutine dlaset(uplo, m, n, alpha, beta, a, lda)
DLASET initializes the off-diagonal elements and the diagonal elements of a matrix to given values.
Definition dlaset.f:108
subroutine drot(n, dx, incx, dy, incy, c, s)
DROT
Definition drot.f:92