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

◆ zgghrd()

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

ZGGHRD

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

Purpose:
!>
!> ZGGHRD 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 ZGGHRD 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 ZGGBAL; 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*16 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*16 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*16 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*16 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 200 of file zgghrd.f.

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*16 A( LDA, * ), B( LDB, * ), Q( LDQ, * ),
214 $ Z( LDZ, * )
215* ..
216*
217* =====================================================================
218*
219* .. Parameters ..
220 COMPLEX*16 CONE, CZERO
221 parameter( cone = ( 1.0d+0, 0.0d+0 ),
222 $ czero = ( 0.0d+0, 0.0d+0 ) )
223* ..
224* .. Local Scalars ..
225 LOGICAL ILQ, ILZ
226 INTEGER ICOMPQ, ICOMPZ, JCOL, JROW
227 DOUBLE PRECISION C
228 COMPLEX*16 CTEMP, S
229* ..
230* .. External Functions ..
231 LOGICAL LSAME
232 EXTERNAL lsame
233* ..
234* .. External Subroutines ..
235 EXTERNAL xerbla, zlartg, zlaset, zrot
236* ..
237* .. Intrinsic Functions ..
238 INTRINSIC dconjg, 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( 'ZGGHRD', -info )
296 RETURN
297 END IF
298*
299* Initialize Q and Z if desired.
300*
301 IF( icompq.EQ.3 )
302 $ CALL zlaset( 'Full', n, n, czero, cone, q, ldq )
303 IF( icompz.EQ.3 )
304 $ CALL zlaset( '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 zlartg( ctemp, a( jrow, jcol ), c, s,
329 $ a( jrow-1, jcol ) )
330 a( jrow, jcol ) = czero
331 CALL zrot( n-jcol, a( jrow-1, jcol+1 ), lda,
332 $ a( jrow, jcol+1 ), lda, c, s )
333 CALL zrot( n+2-jrow, b( jrow-1, jrow-1 ), ldb,
334 $ b( jrow, jrow-1 ), ldb, c, s )
335 IF( ilq )
336 $ CALL zrot( n, q( 1, jrow-1 ), 1, q( 1, jrow ), 1, c,
337 $ dconjg( s ) )
338*
339* Step 2: rotate columns JROW, JROW-1 to kill B(JROW,JROW-1)
340*
341 ctemp = b( jrow, jrow )
342 CALL zlartg( ctemp, b( jrow, jrow-1 ), c, s,
343 $ b( jrow, jrow ) )
344 b( jrow, jrow-1 ) = czero
345 CALL zrot( ihi, a( 1, jrow ), 1, a( 1, jrow-1 ), 1, c,
346 $ s )
347 CALL zrot( jrow-1, b( 1, jrow ), 1, b( 1, jrow-1 ), 1, c,
348 $ s )
349 IF( ilz )
350 $ CALL zrot( 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 ZGGHRD
358*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine zlartg(f, g, c, s, r)
ZLARTG generates a plane rotation with real cosine and complex sine.
Definition zlartg.f90:116
subroutine zlaset(uplo, m, n, alpha, beta, a, lda)
ZLASET initializes the off-diagonal elements and the diagonal elements of a matrix to given values.
Definition zlaset.f:104
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
subroutine zrot(n, cx, incx, cy, incy, c, s)
ZROT applies a plane rotation with real cosine and complex sine to a pair of complex vectors.
Definition zrot.f:101
Here is the call graph for this function:
Here is the caller graph for this function: