LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
zlaqz0.f
Go to the documentation of this file.
1*> \brief \b ZLAQZ0
2*
3* =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6* http://www.netlib.org/lapack/explore-html/
7*
8*> \htmlonly
9*> Download ZLAQZ0 + dependencies
10*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/zlaqz0.f">
11*> [TGZ]</a>
12*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/zlaqz0.f">
13*> [ZIP]</a>
14*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/zlaqz0.f">
15*> [TXT]</a>
16*> \endhtmlonly
17*
18* Definition:
19* ===========
20*
21* SUBROUTINE ZLAQZ0( WANTS, WANTQ, WANTZ, N, ILO, IHI, A, LDA, B,
22* $ LDB, ALPHA, BETA, Q, LDQ, Z, LDZ, WORK, LWORK, RWORK, REC,
23* $ INFO )
24* IMPLICIT NONE
25*
26* Arguments
27* CHARACTER, INTENT( IN ) :: WANTS, WANTQ, WANTZ
28* INTEGER, INTENT( IN ) :: N, ILO, IHI, LDA, LDB, LDQ, LDZ, LWORK,
29* $ REC
30* INTEGER, INTENT( OUT ) :: INFO
31* COMPLEX*16, INTENT( INOUT ) :: A( LDA, * ), B( LDB, * ), Q( LDQ,
32* $ * ), Z( LDZ, * ), ALPHA( * ), BETA( * ), WORK( * )
33* DOUBLE PRECISION, INTENT( OUT ) :: RWORK( * )
34* ..
35*
36*
37*> \par Purpose:
38* =============
39*>
40*> \verbatim
41*>
42*> ZLAQZ0 computes the eigenvalues of a real matrix pair (H,T),
43*> where H is an upper Hessenberg matrix and T is upper triangular,
44*> using the double-shift QZ method.
45*> Matrix pairs of this type are produced by the reduction to
46*> generalized upper Hessenberg form of a real matrix pair (A,B):
47*>
48*> A = Q1*H*Z1**H, B = Q1*T*Z1**H,
49*>
50*> as computed by ZGGHRD.
51*>
52*> If JOB='S', then the Hessenberg-triangular pair (H,T) is
53*> also reduced to generalized Schur form,
54*>
55*> H = Q*S*Z**H, T = Q*P*Z**H,
56*>
57*> where Q and Z are unitary matrices, P and S are an upper triangular
58*> matrices.
59*>
60*> Optionally, the unitary matrix Q from the generalized Schur
61*> factorization may be postmultiplied into an input matrix Q1, and the
62*> unitary matrix Z may be postmultiplied into an input matrix Z1.
63*> If Q1 and Z1 are the unitary matrices from ZGGHRD that reduced
64*> the matrix pair (A,B) to generalized upper Hessenberg form, then the
65*> output matrices Q1*Q and Z1*Z are the unitary factors from the
66*> generalized Schur factorization of (A,B):
67*>
68*> A = (Q1*Q)*S*(Z1*Z)**H, B = (Q1*Q)*P*(Z1*Z)**H.
69*>
70*> To avoid overflow, eigenvalues of the matrix pair (H,T) (equivalently,
71*> of (A,B)) are computed as a pair of values (alpha,beta), where alpha is
72*> complex and beta real.
73*> If beta is nonzero, lambda = alpha / beta is an eigenvalue of the
74*> generalized nonsymmetric eigenvalue problem (GNEP)
75*> A*x = lambda*B*x
76*> and if alpha is nonzero, mu = beta / alpha is an eigenvalue of the
77*> alternate form of the GNEP
78*> mu*A*y = B*y.
79*> Eigenvalues can be read directly from the generalized Schur
80*> form:
81*> alpha = S(i,i), beta = P(i,i).
82*>
83*> Ref: C.B. Moler & G.W. Stewart, "An Algorithm for Generalized Matrix
84*> Eigenvalue Problems", SIAM J. Numer. Anal., 10(1973),
85*> pp. 241--256.
86*>
87*> Ref: B. Kagstrom, D. Kressner, "Multishift Variants of the QZ
88*> Algorithm with Aggressive Early Deflation", SIAM J. Numer.
89*> Anal., 29(2006), pp. 199--227.
90*>
91*> Ref: T. Steel, D. Camps, K. Meerbergen, R. Vandebril "A multishift,
92*> multipole rational QZ method with aggressive early deflation"
93*> \endverbatim
94*
95* Arguments:
96* ==========
97*
98*> \param[in] WANTS
99*> \verbatim
100*> WANTS is CHARACTER*1
101*> = 'E': Compute eigenvalues only;
102*> = 'S': Compute eigenvalues and the Schur form.
103*> \endverbatim
104*>
105*> \param[in] WANTQ
106*> \verbatim
107*> WANTQ is CHARACTER*1
108*> = 'N': Left Schur vectors (Q) are not computed;
109*> = 'I': Q is initialized to the unit matrix and the matrix Q
110*> of left Schur vectors of (A,B) is returned;
111*> = 'V': Q must contain an unitary matrix Q1 on entry and
112*> the product Q1*Q is returned.
113*> \endverbatim
114*>
115*> \param[in] WANTZ
116*> \verbatim
117*> WANTZ is CHARACTER*1
118*> = 'N': Right Schur vectors (Z) are not computed;
119*> = 'I': Z is initialized to the unit matrix and the matrix Z
120*> of right Schur vectors of (A,B) is returned;
121*> = 'V': Z must contain an unitary matrix Z1 on entry and
122*> the product Z1*Z is returned.
123*> \endverbatim
124*>
125*> \param[in] N
126*> \verbatim
127*> N is INTEGER
128*> The order of the matrices A, B, Q, and Z. N >= 0.
129*> \endverbatim
130*>
131*> \param[in] ILO
132*> \verbatim
133*> ILO is INTEGER
134*> \endverbatim
135*>
136*> \param[in] IHI
137*> \verbatim
138*> IHI is INTEGER
139*> ILO and IHI mark the rows and columns of A which are in
140*> Hessenberg form. It is assumed that A is already upper
141*> triangular in rows and columns 1:ILO-1 and IHI+1:N.
142*> If N > 0, 1 <= ILO <= IHI <= N; if N = 0, ILO=1 and IHI=0.
143*> \endverbatim
144*>
145*> \param[in,out] A
146*> \verbatim
147*> A is COMPLEX*16 array, dimension (LDA, N)
148*> On entry, the N-by-N upper Hessenberg matrix A.
149*> On exit, if JOB = 'S', A contains the upper triangular
150*> matrix S from the generalized Schur factorization.
151*> If JOB = 'E', the diagonal blocks of A match those of S, but
152*> the rest of A is unspecified.
153*> \endverbatim
154*>
155*> \param[in] LDA
156*> \verbatim
157*> LDA is INTEGER
158*> The leading dimension of the array A. LDA >= max( 1, N ).
159*> \endverbatim
160*>
161*> \param[in,out] B
162*> \verbatim
163*> B is COMPLEX*16 array, dimension (LDB, N)
164*> On entry, the N-by-N upper triangular matrix B.
165*> On exit, if JOB = 'S', B contains the upper triangular
166*> matrix P from the generalized Schur factorization;
167*> If JOB = 'E', the diagonal blocks of B match those of P, but
168*> the rest of B is unspecified.
169*> \endverbatim
170*>
171*> \param[in] LDB
172*> \verbatim
173*> LDB is INTEGER
174*> The leading dimension of the array B. LDB >= max( 1, N ).
175*> \endverbatim
176*>
177*> \param[out] ALPHA
178*> \verbatim
179*> ALPHA is COMPLEX*16 array, dimension (N)
180*> Each scalar alpha defining an eigenvalue
181*> of GNEP.
182*> \endverbatim
183*>
184*> \param[out] BETA
185*> \verbatim
186*> BETA is COMPLEX*16 array, dimension (N)
187*> The scalars beta that define the eigenvalues of GNEP.
188*> Together, the quantities alpha = ALPHA(j) and
189*> beta = BETA(j) represent the j-th eigenvalue of the matrix
190*> pair (A,B), in one of the forms lambda = alpha/beta or
191*> mu = beta/alpha. Since either lambda or mu may overflow,
192*> they should not, in general, be computed.
193*> \endverbatim
194*>
195*> \param[in,out] Q
196*> \verbatim
197*> Q is COMPLEX*16 array, dimension (LDQ, N)
198*> On entry, if COMPQ = 'V', the unitary matrix Q1 used in
199*> the reduction of (A,B) to generalized Hessenberg form.
200*> On exit, if COMPQ = 'I', the unitary matrix of left Schur
201*> vectors of (A,B), and if COMPQ = 'V', the unitary matrix
202*> of left Schur vectors of (A,B).
203*> Not referenced if COMPQ = 'N'.
204*> \endverbatim
205*>
206*> \param[in] LDQ
207*> \verbatim
208*> LDQ is INTEGER
209*> The leading dimension of the array Q. LDQ >= 1.
210*> If COMPQ='V' or 'I', then LDQ >= N.
211*> \endverbatim
212*>
213*> \param[in,out] Z
214*> \verbatim
215*> Z is COMPLEX*16 array, dimension (LDZ, N)
216*> On entry, if COMPZ = 'V', the unitary matrix Z1 used in
217*> the reduction of (A,B) to generalized Hessenberg form.
218*> On exit, if COMPZ = 'I', the unitary matrix of
219*> right Schur vectors of (H,T), and if COMPZ = 'V', the
220*> unitary matrix of right Schur vectors of (A,B).
221*> Not referenced if COMPZ = 'N'.
222*> \endverbatim
223*>
224*> \param[in] LDZ
225*> \verbatim
226*> LDZ is INTEGER
227*> The leading dimension of the array Z. LDZ >= 1.
228*> If COMPZ='V' or 'I', then LDZ >= N.
229*> \endverbatim
230*>
231*> \param[out] WORK
232*> \verbatim
233*> WORK is COMPLEX*16 array, dimension (MAX(1,LWORK))
234*> On exit, if INFO >= 0, WORK(1) returns the optimal LWORK.
235*> \endverbatim
236*>
237*> \param[out] RWORK
238*> \verbatim
239*> RWORK is DOUBLE PRECISION array, dimension (N)
240*> \endverbatim
241*>
242*> \param[in] LWORK
243*> \verbatim
244*> LWORK is INTEGER
245*> The dimension of the array WORK. LWORK >= max(1,N).
246*>
247*> If LWORK = -1, then a workspace query is assumed; the routine
248*> only calculates the optimal size of the WORK array, returns
249*> this value as the first entry of the WORK array, and no error
250*> message related to LWORK is issued by XERBLA.
251*> \endverbatim
252*>
253*> \param[in] REC
254*> \verbatim
255*> REC is INTEGER
256*> REC indicates the current recursion level. Should be set
257*> to 0 on first call.
258*> \endverbatim
259*>
260*> \param[out] INFO
261*> \verbatim
262*> INFO is INTEGER
263*> = 0: successful exit
264*> < 0: if INFO = -i, the i-th argument had an illegal value
265*> = 1,...,N: the QZ iteration did not converge. (A,B) is not
266*> in Schur form, but ALPHA(i) and
267*> BETA(i), i=INFO+1,...,N should be correct.
268*> \endverbatim
269*
270* Authors:
271* ========
272*
273*> \author Thijs Steel, KU Leuven
274*
275*> \date May 2020
276*
277*> \ingroup laqz0
278*>
279* =====================================================================
280 RECURSIVE SUBROUTINE zlaqz0( WANTS, WANTQ, WANTZ, N, ILO, IHI, A,
281 $ LDA, B, LDB, ALPHA, BETA, Q, LDQ, Z,
282 $ LDZ, WORK, LWORK, RWORK, REC,
283 $ INFO )
284 IMPLICIT NONE
285
286* Arguments
287 CHARACTER, INTENT( IN ) :: wants, wantq, wantz
288 INTEGER, INTENT( IN ) :: n, ilo, ihi, lda, ldb, ldq, ldz, lwork,
289 $ rec
290 INTEGER, INTENT( OUT ) :: info
291 COMPLEX*16, INTENT( INOUT ) :: a( lda, * ), b( ldb, * ), q( ldq,
292 $ * ), z( ldz, * ), alpha( * ), beta( * ), work( * )
293 DOUBLE PRECISION, INTENT( OUT ) :: rwork( * )
294
295* Parameters
296 COMPLEX*16 czero, cone
297 PARAMETER ( czero = ( 0.0d+0, 0.0d+0 ), cone = ( 1.0d+0,
298 $ 0.0d+0 ) )
299 DOUBLE PRECISION :: zero, one, half
300 parameter( zero = 0.0d0, one = 1.0d0, half = 0.5d0 )
301
302* Local scalars
303 DOUBLE PRECISION :: smlnum, ulp, safmin, safmax, c1, tempr,
304 $ bnorm, btol
305 COMPLEX*16 :: eshift, s1, temp
306 INTEGER :: istart, istop, iiter, maxit, istart2, k, ld, nshifts,
307 $ nblock, nw, nmin, nibble, n_undeflated, n_deflated,
308 $ ns, sweep_info, shiftpos, lworkreq, k2, istartm,
309 $ istopm, iwants, iwantq, iwantz, norm_info, aed_info,
310 $ nwr, nbr, nsr, itemp1, itemp2, rcost
311 LOGICAL :: ilschur, ilq, ilz
312 CHARACTER :: jbcmpz*3
313
314* External Functions
315 EXTERNAL :: xerbla, zhgeqz, zlaqz2, zlaqz3, zlaset,
316 $ zlartg, zrot
317 DOUBLE PRECISION, EXTERNAL :: dlamch, zlanhs
318 LOGICAL, EXTERNAL :: lsame
319 INTEGER, EXTERNAL :: ilaenv
320
321*
322* Decode wantS,wantQ,wantZ
323*
324 IF( lsame( wants, 'E' ) ) THEN
325 ilschur = .false.
326 iwants = 1
327 ELSE IF( lsame( wants, 'S' ) ) THEN
328 ilschur = .true.
329 iwants = 2
330 ELSE
331 iwants = 0
332 END IF
333
334 IF( lsame( wantq, 'N' ) ) THEN
335 ilq = .false.
336 iwantq = 1
337 ELSE IF( lsame( wantq, 'V' ) ) THEN
338 ilq = .true.
339 iwantq = 2
340 ELSE IF( lsame( wantq, 'I' ) ) THEN
341 ilq = .true.
342 iwantq = 3
343 ELSE
344 iwantq = 0
345 END IF
346
347 IF( lsame( wantz, 'N' ) ) THEN
348 ilz = .false.
349 iwantz = 1
350 ELSE IF( lsame( wantz, 'V' ) ) THEN
351 ilz = .true.
352 iwantz = 2
353 ELSE IF( lsame( wantz, 'I' ) ) THEN
354 ilz = .true.
355 iwantz = 3
356 ELSE
357 iwantz = 0
358 END IF
359*
360* Check Argument Values
361*
362 info = 0
363 IF( iwants.EQ.0 ) THEN
364 info = -1
365 ELSE IF( iwantq.EQ.0 ) THEN
366 info = -2
367 ELSE IF( iwantz.EQ.0 ) THEN
368 info = -3
369 ELSE IF( n.LT.0 ) THEN
370 info = -4
371 ELSE IF( ilo.LT.1 ) THEN
372 info = -5
373 ELSE IF( ihi.GT.n .OR. ihi.LT.ilo-1 ) THEN
374 info = -6
375 ELSE IF( lda.LT.n ) THEN
376 info = -8
377 ELSE IF( ldb.LT.n ) THEN
378 info = -10
379 ELSE IF( ldq.LT.1 .OR. ( ilq .AND. ldq.LT.n ) ) THEN
380 info = -15
381 ELSE IF( ldz.LT.1 .OR. ( ilz .AND. ldz.LT.n ) ) THEN
382 info = -17
383 END IF
384 IF( info.NE.0 ) THEN
385 CALL xerbla( 'ZLAQZ0', -info )
386 RETURN
387 END IF
388
389*
390* Quick return if possible
391*
392 IF( n.LE.0 ) THEN
393 work( 1 ) = dble( 1 )
394 RETURN
395 END IF
396
397*
398* Get the parameters
399*
400 jbcmpz( 1:1 ) = wants
401 jbcmpz( 2:2 ) = wantq
402 jbcmpz( 3:3 ) = wantz
403
404 nmin = ilaenv( 12, 'ZLAQZ0', jbcmpz, n, ilo, ihi, lwork )
405
406 nwr = ilaenv( 13, 'ZLAQZ0', jbcmpz, n, ilo, ihi, lwork )
407 nwr = max( 2, nwr )
408 nwr = min( ihi-ilo+1, ( n-1 ) / 3, nwr )
409
410 nibble = ilaenv( 14, 'ZLAQZ0', jbcmpz, n, ilo, ihi, lwork )
411
412 nsr = ilaenv( 15, 'ZLAQZ0', jbcmpz, n, ilo, ihi, lwork )
413 nsr = min( nsr, ( n+6 ) / 9, ihi-ilo )
414 nsr = max( 2, nsr-mod( nsr, 2 ) )
415
416 rcost = ilaenv( 17, 'ZLAQZ0', jbcmpz, n, ilo, ihi, lwork )
417 itemp1 = int( nsr/sqrt( 1+2*nsr/( dble( rcost )/100*n ) ) )
418 itemp1 = ( ( itemp1-1 )/4 )*4+4
419 nbr = nsr+itemp1
420
421 IF( n .LT. nmin .OR. rec .GE. 2 ) THEN
422 CALL zhgeqz( wants, wantq, wantz, n, ilo, ihi, a, lda, b, ldb,
423 $ alpha, beta, q, ldq, z, ldz, work, lwork, rwork,
424 $ info )
425 RETURN
426 END IF
427
428*
429* Find out required workspace
430*
431
432* Workspace query to ZLAQZ2
433 nw = max( nwr, nmin )
434 CALL zlaqz2( ilschur, ilq, ilz, n, ilo, ihi, nw, a, lda, b, ldb,
435 $ q, ldq, z, ldz, n_undeflated, n_deflated, alpha,
436 $ beta, work, nw, work, nw, work, -1, rwork, rec,
437 $ aed_info )
438 itemp1 = int( work( 1 ) )
439* Workspace query to ZLAQZ3
440 CALL zlaqz3( ilschur, ilq, ilz, n, ilo, ihi, nsr, nbr, alpha,
441 $ beta, a, lda, b, ldb, q, ldq, z, ldz, work, nbr,
442 $ work, nbr, work, -1, sweep_info )
443 itemp2 = int( work( 1 ) )
444
445 lworkreq = max( itemp1+2*nw**2, itemp2+2*nbr**2 )
446 IF ( lwork .EQ.-1 ) THEN
447 work( 1 ) = dble( lworkreq )
448 RETURN
449 ELSE IF ( lwork .LT. lworkreq ) THEN
450 info = -19
451 END IF
452 IF( info.NE.0 ) THEN
453 CALL xerbla( 'ZLAQZ0', info )
454 RETURN
455 END IF
456*
457* Initialize Q and Z
458*
459 IF( iwantq.EQ.3 ) CALL zlaset( 'FULL', n, n, czero, cone, q,
460 $ ldq )
461 IF( iwantz.EQ.3 ) CALL zlaset( 'FULL', n, n, czero, cone, z,
462 $ ldz )
463
464* Get machine constants
465 safmin = dlamch( 'SAFE MINIMUM' )
466 safmax = one/safmin
467 ulp = dlamch( 'PRECISION' )
468 smlnum = safmin*( dble( n )/ulp )
469
470 bnorm = zlanhs( 'F', ihi-ilo+1, b( ilo, ilo ), ldb, rwork )
471 btol = max( safmin, ulp*bnorm )
472
473 istart = ilo
474 istop = ihi
475 maxit = 30*( ihi-ilo+1 )
476 ld = 0
477
478 DO iiter = 1, maxit
479 IF( iiter .GE. maxit ) THEN
480 info = istop+1
481 GOTO 80
482 END IF
483 IF ( istart+1 .GE. istop ) THEN
484 istop = istart
485 EXIT
486 END IF
487
488* Check deflations at the end
489 IF ( abs( a( istop, istop-1 ) ) .LE. max( smlnum,
490 $ ulp*( abs( a( istop, istop ) )+abs( a( istop-1,
491 $ istop-1 ) ) ) ) ) THEN
492 a( istop, istop-1 ) = czero
493 istop = istop-1
494 ld = 0
495 eshift = czero
496 END IF
497* Check deflations at the start
498 IF ( abs( a( istart+1, istart ) ) .LE. max( smlnum,
499 $ ulp*( abs( a( istart, istart ) )+abs( a( istart+1,
500 $ istart+1 ) ) ) ) ) THEN
501 a( istart+1, istart ) = czero
502 istart = istart+1
503 ld = 0
504 eshift = czero
505 END IF
506
507 IF ( istart+1 .GE. istop ) THEN
508 EXIT
509 END IF
510
511* Check interior deflations
512 istart2 = istart
513 DO k = istop, istart+1, -1
514 IF ( abs( a( k, k-1 ) ) .LE. max( smlnum, ulp*( abs( a( k,
515 $ k ) )+abs( a( k-1, k-1 ) ) ) ) ) THEN
516 a( k, k-1 ) = czero
517 istart2 = k
518 EXIT
519 END IF
520 END DO
521
522* Get range to apply rotations to
523 IF ( ilschur ) THEN
524 istartm = 1
525 istopm = n
526 ELSE
527 istartm = istart2
528 istopm = istop
529 END IF
530
531* Check infinite eigenvalues, this is done without blocking so might
532* slow down the method when many infinite eigenvalues are present
533 k = istop
534 DO WHILE ( k.GE.istart2 )
535
536 IF( abs( b( k, k ) ) .LT. btol ) THEN
537* A diagonal element of B is negligible, move it
538* to the top and deflate it
539
540 DO k2 = k, istart2+1, -1
541 CALL zlartg( b( k2-1, k2 ), b( k2-1, k2-1 ), c1, s1,
542 $ temp )
543 b( k2-1, k2 ) = temp
544 b( k2-1, k2-1 ) = czero
545
546 CALL zrot( k2-2-istartm+1, b( istartm, k2 ), 1,
547 $ b( istartm, k2-1 ), 1, c1, s1 )
548 CALL zrot( min( k2+1, istop )-istartm+1, a( istartm,
549 $ k2 ), 1, a( istartm, k2-1 ), 1, c1, s1 )
550 IF ( ilz ) THEN
551 CALL zrot( n, z( 1, k2 ), 1, z( 1, k2-1 ), 1, c1,
552 $ s1 )
553 END IF
554
555 IF( k2.LT.istop ) THEN
556 CALL zlartg( a( k2, k2-1 ), a( k2+1, k2-1 ), c1,
557 $ s1, temp )
558 a( k2, k2-1 ) = temp
559 a( k2+1, k2-1 ) = czero
560
561 CALL zrot( istopm-k2+1, a( k2, k2 ), lda, a( k2+1,
562 $ k2 ), lda, c1, s1 )
563 CALL zrot( istopm-k2+1, b( k2, k2 ), ldb, b( k2+1,
564 $ k2 ), ldb, c1, s1 )
565 IF( ilq ) THEN
566 CALL zrot( n, q( 1, k2 ), 1, q( 1, k2+1 ), 1,
567 $ c1, dconjg( s1 ) )
568 END IF
569 END IF
570
571 END DO
572
573 IF( istart2.LT.istop )THEN
574 CALL zlartg( a( istart2, istart2 ), a( istart2+1,
575 $ istart2 ), c1, s1, temp )
576 a( istart2, istart2 ) = temp
577 a( istart2+1, istart2 ) = czero
578
579 CALL zrot( istopm-( istart2+1 )+1, a( istart2,
580 $ istart2+1 ), lda, a( istart2+1,
581 $ istart2+1 ), lda, c1, s1 )
582 CALL zrot( istopm-( istart2+1 )+1, b( istart2,
583 $ istart2+1 ), ldb, b( istart2+1,
584 $ istart2+1 ), ldb, c1, s1 )
585 IF( ilq ) THEN
586 CALL zrot( n, q( 1, istart2 ), 1, q( 1,
587 $ istart2+1 ), 1, c1, dconjg( s1 ) )
588 END IF
589 END IF
590
591 istart2 = istart2+1
592
593 END IF
594 k = k-1
595 END DO
596
597* istart2 now points to the top of the bottom right
598* unreduced Hessenberg block
599 IF ( istart2 .GE. istop ) THEN
600 istop = istart2-1
601 ld = 0
602 eshift = czero
603 cycle
604 END IF
605
606 nw = nwr
607 nshifts = nsr
608 nblock = nbr
609
610 IF ( istop-istart2+1 .LT. nmin ) THEN
611* Setting nw to the size of the subblock will make AED deflate
612* all the eigenvalues. This is slightly more efficient than just
613* using qz_small because the off diagonal part gets updated via BLAS.
614 IF ( istop-istart+1 .LT. nmin ) THEN
615 nw = istop-istart+1
616 istart2 = istart
617 ELSE
618 nw = istop-istart2+1
619 END IF
620 END IF
621
622*
623* Time for AED
624*
625 CALL zlaqz2( ilschur, ilq, ilz, n, istart2, istop, nw, a, lda,
626 $ b, ldb, q, ldq, z, ldz, n_undeflated, n_deflated,
627 $ alpha, beta, work, nw, work( nw**2+1 ), nw,
628 $ work( 2*nw**2+1 ), lwork-2*nw**2, rwork, rec,
629 $ aed_info )
630
631 IF ( n_deflated > 0 ) THEN
632 istop = istop-n_deflated
633 ld = 0
634 eshift = czero
635 END IF
636
637 IF ( 100*n_deflated > nibble*( n_deflated+n_undeflated ) .OR.
638 $ istop-istart2+1 .LT. nmin ) THEN
639* AED has uncovered many eigenvalues. Skip a QZ sweep and run
640* AED again.
641 cycle
642 END IF
643
644 ld = ld+1
645
646 ns = min( nshifts, istop-istart2 )
647 ns = min( ns, n_undeflated )
648 shiftpos = istop-n_undeflated+1
649
650 IF ( mod( ld, 6 ) .EQ. 0 ) THEN
651*
652* Exceptional shift. Chosen for no particularly good reason.
653*
654 IF( ( dble( maxit )*safmin )*abs( a( istop,
655 $ istop-1 ) ).LT.abs( a( istop-1, istop-1 ) ) ) THEN
656 eshift = a( istop, istop-1 )/b( istop-1, istop-1 )
657 ELSE
658 eshift = eshift+cone/( safmin*dble( maxit ) )
659 END IF
660 alpha( shiftpos ) = cone
661 beta( shiftpos ) = eshift
662 ns = 1
663 END IF
664
665*
666* Time for a QZ sweep
667*
668 CALL zlaqz3( ilschur, ilq, ilz, n, istart2, istop, ns, nblock,
669 $ alpha( shiftpos ), beta( shiftpos ), a, lda, b,
670 $ ldb, q, ldq, z, ldz, work, nblock, work( nblock**
671 $ 2+1 ), nblock, work( 2*nblock**2+1 ),
672 $ lwork-2*nblock**2, sweep_info )
673
674 END DO
675
676*
677* Call ZHGEQZ to normalize the eigenvalue blocks and set the eigenvalues
678* If all the eigenvalues have been found, ZHGEQZ will not do any iterations
679* and only normalize the blocks. In case of a rare convergence failure,
680* the single shift might perform better.
681*
682 80 CALL zhgeqz( wants, wantq, wantz, n, ilo, ihi, a, lda, b, ldb,
683 $ alpha, beta, q, ldq, z, ldz, work, lwork, rwork,
684 $ norm_info )
685
686 info = norm_info
687
688 END SUBROUTINE
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine zhgeqz(job, compq, compz, n, ilo, ihi, h, ldh, t, ldt, alpha, beta, q, ldq, z, ldz, work, lwork, rwork, info)
ZHGEQZ
Definition zhgeqz.f:284
integer function ilaenv(ispec, name, opts, n1, n2, n3, n4)
ILAENV
Definition ilaenv.f:162
double precision function dlamch(cmach)
DLAMCH
Definition dlamch.f:69
double precision function zlanhs(norm, n, a, lda, work)
ZLANHS returns the value of the 1-norm, Frobenius norm, infinity-norm, or the largest absolute value ...
Definition zlanhs.f:109
recursive subroutine zlaqz0(wants, wantq, wantz, n, ilo, ihi, a, lda, b, ldb, alpha, beta, q, ldq, z, ldz, work, lwork, rwork, rec, info)
ZLAQZ0
Definition zlaqz0.f:284
recursive subroutine zlaqz2(ilschur, ilq, ilz, n, ilo, ihi, nw, a, lda, b, ldb, q, ldq, z, ldz, ns, nd, alpha, beta, qc, ldqc, zc, ldzc, work, lwork, rwork, rec, info)
ZLAQZ2
Definition zlaqz2.f:234
subroutine zlaqz3(ilschur, ilq, ilz, n, ilo, ihi, nshifts, nblock_desired, alpha, beta, a, lda, b, ldb, q, ldq, z, ldz, qc, ldqc, zc, ldzc, work, lwork, info)
ZLAQZ3
Definition zlaqz3.f:208
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:106
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:103