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

◆ cgghrd()

subroutine cgghrd ( character  compq,
character  compz,
integer  n,
integer  ilo,
integer  ihi,
complex, dimension( lda, * )  a,
integer  lda,
complex, dimension( ldb, * )  b,
integer  ldb,
complex, dimension( ldq, * )  q,
integer  ldq,
complex, dimension( ldz, * )  z,
integer  ldz,
integer  info 
)

CGGHRD

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

Purpose:
 CGGHRD reduces a pair of complex matrices (A,B) to generalized upper
 Hessenberg form using unitary transformations, where A is a
 general matrix and B is upper triangular.  The form of the generalized
 eigenvalue problem is
    A*x = lambda*B*x,
 and B is typically made upper triangular by computing its QR
 factorization and moving the unitary matrix Q to the left side
 of the equation.

 This subroutine simultaneously reduces A to a Hessenberg matrix H:
    Q**H*A*Z = H
 and transforms B to another upper triangular matrix T:
    Q**H*B*Z = T
 in order to reduce the problem to its standard form
    H*y = lambda*T*y
 where y = Z**H*x.

 The unitary matrices Q and Z are determined as products of Givens
 rotations.  They may either be formed explicitly, or they may be
 postmultiplied into input matrices Q1 and Z1, so that
      Q1 * A * Z1**H = (Q1*Q) * H * (Z1*Z)**H
      Q1 * B * Z1**H = (Q1*Q) * T * (Z1*Z)**H
 If Q1 is the unitary matrix from the QR factorization of B in the
 original equation A*x = lambda*B*x, then CGGHRD reduces the original
 problem to generalized Hessenberg form.
Parameters
[in]COMPQ
          COMPQ is CHARACTER*1
          = 'N': do not compute Q;
          = 'I': Q is initialized to the unit matrix, and the
                 unitary matrix Q is returned;
          = 'V': Q must contain a unitary matrix Q1 on entry,
                 and the product Q1*Q is returned.
[in]COMPZ
          COMPZ is CHARACTER*1
          = 'N': do not compute Z;
          = 'I': Z is initialized to the unit matrix, and the
                 unitary matrix Z is returned;
          = 'V': Z must contain a unitary matrix Z1 on entry,
                 and the product Z1*Z is returned.
[in]N
          N is INTEGER
          The order of the matrices A and B.  N >= 0.
[in]ILO
          ILO is INTEGER
[in]IHI
          IHI is INTEGER

          ILO and IHI mark the rows and columns of A which are to be
          reduced.  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 CGGBAL; otherwise they
          should be set to 1 and N respectively.
          1 <= ILO <= IHI <= N, if N > 0; ILO=1 and IHI=0, if N=0.
[in,out]A
          A is COMPLEX 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
          rest is set to zero.
[in]LDA
          LDA is INTEGER
          The leading dimension of the array A.  LDA >= max(1,N).
[in,out]B
          B is COMPLEX array, dimension (LDB, N)
          On entry, the N-by-N upper triangular matrix B.
          On exit, the upper triangular matrix T = Q**H B Z.  The
          elements below the diagonal are set to zero.
[in]LDB
          LDB is INTEGER
          The leading dimension of the array B.  LDB >= max(1,N).
[in,out]Q
          Q is COMPLEX array, dimension (LDQ, N)
          On entry, if COMPQ = 'V', the unitary matrix Q1, typically
          from the QR factorization of B.
          On exit, if COMPQ='I', the unitary matrix Q, and if
          COMPQ = 'V', the product Q1*Q.
          Not referenced if COMPQ='N'.
[in]LDQ
          LDQ is INTEGER
          The leading dimension of the array Q.
          LDQ >= N if COMPQ='V' or 'I'; LDQ >= 1 otherwise.
[in,out]Z
          Z is COMPLEX array, dimension (LDZ, N)
          On entry, if COMPZ = 'V', the unitary matrix Z1.
          On exit, if COMPZ='I', the unitary matrix Z, and if
          COMPZ = 'V', the product Z1*Z.
          Not referenced if COMPZ='N'.
[in]LDZ
          LDZ is INTEGER
          The leading dimension of the array Z.
          LDZ >= N if COMPZ='V' or 'I'; LDZ >= 1 otherwise.
[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:
  This routine reduces A to Hessenberg and B to triangular form by
  an unblocked reduction, as described in _Matrix_Computations_,
  by Golub and van Loan (Johns Hopkins Press).

Definition at line 202 of file cgghrd.f.

204*
205* -- LAPACK computational routine --
206* -- LAPACK is a software package provided by Univ. of Tennessee, --
207* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
208*
209* .. Scalar Arguments ..
210 CHARACTER COMPQ, COMPZ
211 INTEGER IHI, ILO, INFO, LDA, LDB, LDQ, LDZ, N
212* ..
213* .. Array Arguments ..
214 COMPLEX A( LDA, * ), B( LDB, * ), Q( LDQ, * ),
215 $ Z( LDZ, * )
216* ..
217*
218* =====================================================================
219*
220* .. Parameters ..
221 COMPLEX CONE, CZERO
222 parameter( cone = ( 1.0e+0, 0.0e+0 ),
223 $ czero = ( 0.0e+0, 0.0e+0 ) )
224* ..
225* .. Local Scalars ..
226 LOGICAL ILQ, ILZ
227 INTEGER ICOMPQ, ICOMPZ, JCOL, JROW
228 REAL C
229 COMPLEX CTEMP, S
230* ..
231* .. External Functions ..
232 LOGICAL LSAME
233 EXTERNAL lsame
234* ..
235* .. External Subroutines ..
236 EXTERNAL clartg, claset, crot, xerbla
237* ..
238* .. Intrinsic Functions ..
239 INTRINSIC conjg, 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( 'CGGHRD', -info )
297 RETURN
298 END IF
299*
300* Initialize Q and Z if desired.
301*
302 IF( icompq.EQ.3 )
303 $ CALL claset( 'Full', n, n, czero, cone, q, ldq )
304 IF( icompz.EQ.3 )
305 $ CALL claset( 'Full', n, n, czero, cone, 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 ) = czero
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 ctemp = a( jrow-1, jcol )
329 CALL clartg( ctemp, a( jrow, jcol ), c, s,
330 $ a( jrow-1, jcol ) )
331 a( jrow, jcol ) = czero
332 CALL crot( n-jcol, a( jrow-1, jcol+1 ), lda,
333 $ a( jrow, jcol+1 ), lda, c, s )
334 CALL crot( n+2-jrow, b( jrow-1, jrow-1 ), ldb,
335 $ b( jrow, jrow-1 ), ldb, c, s )
336 IF( ilq )
337 $ CALL crot( n, q( 1, jrow-1 ), 1, q( 1, jrow ), 1, c,
338 $ conjg( s ) )
339*
340* Step 2: rotate columns JROW, JROW-1 to kill B(JROW,JROW-1)
341*
342 ctemp = b( jrow, jrow )
343 CALL clartg( ctemp, b( jrow, jrow-1 ), c, s,
344 $ b( jrow, jrow ) )
345 b( jrow, jrow-1 ) = czero
346 CALL crot( ihi, a( 1, jrow ), 1, a( 1, jrow-1 ), 1, c, 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, s )
351 30 CONTINUE
352 40 CONTINUE
353*
354 RETURN
355*
356* End of CGGHRD
357*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
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:106
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
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:103
Here is the call graph for this function:
Here is the caller graph for this function: