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