control_systems_torbox 0.2.1

Control systems toolbox
Documentation
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
      SUBROUTINE SG03BU( TRANS, N, A, LDA, E, LDE, B, LDB, SCALE, DWORK,
     $                   INFO )
C
C     PURPOSE
C
C     To compute the Cholesky factor U of the matrix X, X = U**T * U or
C     X = U * U**T, which is the solution of the generalized d-stable
C     discrete-time Lyapunov equation
C
C         T            T                  2    T
C        A  * X * A - E  * X * E = - SCALE  * B  * B,                (1)
C
C     or the transposed equation
C
C                 T            T          2        T
C        A * X * A  - E * X * E  = - SCALE  * B * B ,                (2)
C
C     respectively, where A, E, B, and U are real N-by-N matrices. The
C     Cholesky factor U of the solution is computed without first
C     finding X. The pencil A - lambda * E must be in generalized Schur
C     form ( A upper quasitriangular, E upper triangular ). Moreover, it
C     must be d-stable, i.e., the moduli of its eigenvalues must be less
C     than one. B must be an upper triangular matrix with non-negative
C     entries on its main diagonal.
C
C     The resulting matrix U is upper triangular. The entries on its
C     main diagonal are non-negative. SCALE is an output scale factor
C     set to avoid overflow in U.
C
C     ARGUMENTS
C
C     Mode Parameters
C
C     TRANS   CHARACTER*1
C             Specifies whether equation (1) or equation (2) is to be
C             solved:
C             = 'N':  Solve equation (1);
C             = 'T':  Solve equation (2).
C
C     Input/Output Parameters
C
C     N       (input) INTEGER
C             The order of the matrices.  N >= 0.
C
C     A       (input) DOUBLE PRECISION array, dimension (LDA,N)
C             The leading N-by-N upper Hessenberg part of this array
C             must contain the quasitriangular matrix A. The elements
C             below the upper Hessenberg part are not referenced.
C
C     LDA     INTEGER
C             The leading dimension of the array A.  LDA >= MAX(1,N).
C
C     E       (input) DOUBLE PRECISION array, dimension (LDE,N)
C             The leading N-by-N upper triangular part of this array
C             must contain the triangular matrix E. The elements below
C             the main diagonal are not referenced.
C
C     LDE     INTEGER
C             The leading dimension of the array E.  LDE >= MAX(1,N).
C
C     B       (input/output) DOUBLE PRECISION array, dimension (LDB,N)
C             On entry, the leading N-by-N upper triangular part of this
C             array must contain the matrix B.
C             On exit, the leading N-by-N upper triangular part of this
C             array contains the solution matrix U. The elements below
C             the main diagonal are not referenced.
C
C     LDB     INTEGER
C             The leading dimension of the array B.  LDB >= MAX(1,N).
C
C     SCALE   (output) DOUBLE PRECISION
C             The scale factor set to avoid overflow in U.
C             0 < SCALE <= 1.
C
C     Workspace
C
C     DWORK   DOUBLE PRECISION array, dimension (6*N-6)
C
C     Error indicator
C
C     INFO    INTEGER
C             = 0:  successful exit;
C             < 0:  if INFO = -i, the i-th argument had an illegal
C                   value;
C             = 1:  the generalized Sylvester equation to be solved in
C                   step II (see METHOD) is (nearly) singular to working
C                   precision;  perturbed values were used to solve the
C                   equation (but the matrices A and E are unchanged);
C             = 2:  the generalized Schur form of the pencil
C                   A - lambda * E contains a 2-by-2 main diagonal block
C                   whose eigenvalues are not a pair of complex
C                   conjugate numbers;
C             = 3:  the pencil A - lambda * E is not d-stable, i.e.,
C                   there are eigenvalues outside the open unit circle;
C             = 4:  the LAPACK routine DSYEVX utilized to factorize M3
C                   failed to converge. This error is unlikely to occur.
C
C     METHOD
C
C     The method [2] used by the routine is an extension of Hammarling's
C     algorithm [1] to generalized Lyapunov equations.
C
C     We present the method for solving equation (1). Equation (2) can
C     be treated in a similar fashion. For simplicity, assume SCALE = 1.
C
C     The matrix A is an upper quasitriangular matrix, i.e., it is a
C     block triangular matrix with square blocks on the main diagonal
C     and the block order at most 2. We use the following partitioning
C     for the matrices A, E, B and the solution matrix U
C
C               ( A11   A12 )        ( E11   E12 )
C           A = (           ),   E = (           ),
C               (   0   A22 )        (   0   E22 )
C
C               ( B11   B12 )        ( U11   U12 )
C           B = (           ),   U = (           ).                  (3)
C               (   0   B22 )        (   0   U22 )
C
C     The size of the (1,1)-blocks is 1-by-1 (iff A(2,1) = 0.0) or
C     2-by-2.
C
C     We compute U11, U12**T, and U22 in three steps.
C
C     Step I:
C
C        From (1) and (3) we get the 1-by-1 or 2-by-2 equation
C
C                T      T                  T      T
C             A11  * U11  * U11 * A11 - E11  * U11  * U11 * E11
C
C                    T
C             = - B11  * B11.
C
C        For brevity, details are omitted here. The technique for
C        computing U11 is similar to those applied to standard Lyapunov
C        equations in Hammarling's algorithm ([1], section 6).
C
C        Furthermore, the auxiliary matrices M1 and M2 defined as
C        follows
C
C                               -1      -1
C           M1 = U11 * A11 * E11   * U11  ,
C
C                         -1      -1
C           M2 = B11 * E11   * U11  ,
C
C        are computed in a numerically reliable way.
C
C     Step II:
C
C        We solve for U12**T the generalized Sylvester equation
C
C              T      T           T      T
C           A22  * U12  * M1 - E22  * U12
C
C                  T           T      T      T      T
C           = - B12  * M2 + E12  * U11  - A12  * U11  * M1.
C
C     Step III:
C
C        One can show that
C
C              T      T                  T      T
C           A22  * U22  * U22 * A22 - E22  * U22  * U22 * E22  =
C
C                T              T
C           - B22  * B22 - y * y                                     (4)
C
C        holds, where y is defined as follows
C
C                  T      T      T      T
C           w = A12  * U11  + A22  * U12 ,
C
C                    T
C           y = ( B12   w ) * M3EV,
C
C        where M3EV is a matrix which fulfils
C
C                ( I-M2*M2**T   -M2*M1**T )              T
C           M3 = (                        ) = M3EV * M3EV .
C                (  -M1*M2**T  I-M1*M1**T )
C
C        M3 is positive semidefinite and its rank is equal to the size
C        of U11. Therefore, a matrix M3EV can be found by solving the
C        symmetric eigenvalue problem for M3 such that y consists of
C        either 1 or 2 rows.
C
C        If B22_tilde is the square triangular matrix arising from the
C        QR-factorization
C
C               ( B22_tilde )     ( B22  )
C           Q * (           )  =  (      ),
C               (     0     )     ( y**T )
C
C        then
C
C                T              T                T
C           - B22  * B22 - y * y   =  - B22_tilde  * B22_tilde.
C
C        Replacing the right hand side in (4) by the term
C        - B22_tilde**T * B22_tilde leads to a generalized Lyapunov
C        equation of lower dimension compared to (1).
C
C     The solution U of the equation (1) can be obtained by recursive
C     application of the steps I to III.
C
C     REFERENCES
C
C     [1] Hammarling, S.J.
C         Numerical solution of the stable, non-negative definite
C         Lyapunov equation.
C         IMA J. Num. Anal., 2, pp. 303-323, 1982.
C
C     [2] Penzl, T.
C         Numerical solution of generalized Lyapunov equations.
C         Advances in Comp. Math., vol. 8, pp. 33-48, 1998.
C
C     NUMERICAL ASPECTS
C
C     The routine requires 2*N**3 flops. Note that we count a single
C     floating point arithmetic operation as one flop.
C
C     FURTHER COMMENTS
C
C     The Lyapunov equation may be very ill-conditioned. In particular,
C     if the pencil A - lambda * E has a pair of almost reciprocal
C     eigenvalues, then the Lyapunov equation will be ill-conditioned.
C     Perturbed values were used to solve the equation.
C     A condition estimate can be obtained from the routine SG03AD.
C     When setting the error indicator INFO, the routine does not test
C     for near instability in the equation but only for exact
C     instability.
C
C     CONTRIBUTOR
C
C     T. Penzl, Technical University Chemnitz, Germany, Aug. 1998.
C
C     REVISIONS
C
C     Sep. 1998, Dec. 2021 (V. Sima).
C
C     KEYWORDS
C
C     Lyapunov equation
C
C     ******************************************************************
C
C     .. Parameters ..
      DOUBLE PRECISION  HALF, MONE, ONE, TWO, ZERO
      PARAMETER         ( HALF = 0.5D+0, MONE = -1.0D0, ONE = 1.0D+0,
     $                    TWO = 2.0D+0,  ZERO = 0.0D+0 )
C     .. Scalar Arguments ..
      CHARACTER         TRANS
      DOUBLE PRECISION  SCALE
      INTEGER           INFO, LDA, LDB, LDE, N
C     .. Array Arguments ..
      DOUBLE PRECISION  A(LDA,*), B(LDB,*), DWORK(*), E(LDE,*)
C     .. Local Scalars ..
      DOUBLE PRECISION  BIGNUM, C, DELTA1, EPS, R, S, SCALE1, SMLNUM, T,
     $                  UFLT, X
      INTEGER           I, INFO1, J, KB, KH, KL, KL1, L, LDWS, M, UIIPT,
     $                  WPT, YPT
      LOGICAL           NOTRNS
C     .. Local Arrays ..
      DOUBLE PRECISION  M1(2,2), M2(2,2), M3(4,4), M3C(4,4), M3EW(4),
     $                  RW(32), TM(2,2), UI(2,2)
      INTEGER           IW(24)
C     .. External Functions ..
      DOUBLE PRECISION  DLAMCH
      LOGICAL           LSAME
      EXTERNAL          DLAMCH, LSAME
C     .. External Subroutines ..
      EXTERNAL          DCOPY, DGEMM, DGEMV, DLABAD, DLACPY, DLARTG,
     $                  DLASCL, DLASET, DROT, DSCAL, DSYEVX, DSYRK,
     $                  SG03BW, SG03BX, XERBLA
C     .. Intrinsic Functions ..
      INTRINSIC         ABS, MAX, MIN, SQRT
C     .. Executable Statements ..
C
C     Decode input parameter.
C
      NOTRNS = LSAME( TRANS, 'N' )
C
C     Check the scalar input parameters.
C
      IF ( .NOT.( NOTRNS .OR. LSAME( TRANS, 'T' ) ) ) THEN
         INFO = -1
      ELSEIF ( N.LT.0 ) THEN
         INFO = -2
      ELSEIF ( LDA.LT.MAX( 1, N ) ) THEN
         INFO = -4
      ELSEIF ( LDE.LT.MAX( 1, N ) ) THEN
         INFO = -6
      ELSEIF ( LDB.LT.MAX( 1, N ) ) THEN
         INFO = -8
      ELSE
         INFO = 0
      END IF
      IF ( INFO.NE.0 ) THEN
         CALL XERBLA( 'SG03BU', -INFO )
         RETURN
      END IF
C
      SCALE = ONE
C
C     Quick return if possible.
C
      IF ( N.EQ.0 )
     $    RETURN
C
C     Set constants to control overflow.
C
      EPS    = DLAMCH( 'P' )
      UFLT   = DLAMCH( 'S' )
      SMLNUM = UFLT/EPS
      BIGNUM = ONE/SMLNUM
      CALL DLABAD( SMLNUM, BIGNUM )
C
C     Set workspace pointers and leading dimension of matrices in the
C     workspace.
C
      UIIPT = 1
      WPT   = 2*N-1
      YPT   = 4*N-3
      LDWS  = N-1
C
      IF ( NOTRNS ) THEN
C
C        Solve equation (1).
C
C        Main Loop. Compute block row U(KL:KH,KL:N). KB denotes the
C        number of rows in this block row.
C
         KH = 0
C        WHILE ( KH.LT.N ) DO
   20    CONTINUE
         IF ( KH.LT.N ) THEN
            KL = KH + 1
            IF ( KL.EQ.N ) THEN
               KH = N
               KB = 1
            ELSE
               IF ( A(KL+1,KL).EQ.ZERO ) THEN
                  KH = KL
                  KB = 1
               ELSE
                  KH = KL + 1
                  KB = 2
               END IF
            END IF
C
C           STEP I: Compute block U(KL:KH,KL:KH) and the auxiliary
C                   matrices M1 and M2. (For the moment the result
C                   U(KL:KH,KL:KH) is stored in UI).
C
            IF ( KB.EQ.1 ) THEN
               DELTA1 = E(KL,KL)
               T      = ABS( A(KL,KL) )
               X      = MAX( DELTA1, T )
               DELTA1 = DELTA1/X
               T      =      T/X
               IF ( DELTA1.LE.T ) THEN
                  INFO = 3
                  RETURN
               END IF
               DELTA1 = SQRT( ONE - T )*SQRT( ONE + T )*X
               T = B(KL,KL)*SMLNUM
               IF ( T.GT.DELTA1 ) THEN
                  SCALE1 = DELTA1/T
                  SCALE  = SCALE1*SCALE
                  CALL DLASCL( 'Upper', 0, 0, ONE, SCALE1, N, N, B, LDB,
     $                         INFO1 )
               END IF
C
               UI(1,1) = B(KL,KL)/DELTA1
               M1(1,1) = A(KL,KL)/E(KL,KL)
               M2(1,1) =   DELTA1/E(KL,KL)
C                
            ELSE
C
C              If a pair of complex conjugate eigenvalues occurs, apply
C              (complex) Hammarling algorithm for the 2-by-2 problem.
C
               CALL SG03BX( 'D', 'N', A(KL,KL), LDA, E(KL,KL), LDE,
     $                      B(KL,KL), LDB, UI, 2, SCALE1, M1, 2, M2, 2,
     $                      INFO )
               IF ( INFO.NE.0 )
     $            RETURN
               IF ( SCALE1.NE.ONE ) THEN
                  SCALE = SCALE1*SCALE
                  CALL DLASCL( 'Upper', 0, 0, ONE, SCALE1, N, N, B, LDB,
     $                         INFO1 )
               END IF
            END IF
C
            IF ( KH.LT.N ) THEN
C
C              STEP II: Compute U(KL:KH,KH+1:N) by solving a generalized
C                       Sylvester equation. (For the moment the result
C                       U(KL:KH,KH+1:N) is stored in the workspace.)
C
C              Form right hand side of the Sylvester equation.
C
               CALL DGEMM( 'T', 'N', N-KH, KB, KB, MONE, B(KL,KH+1),
     $                     LDB, M2, 2, ZERO, DWORK(UIIPT), LDWS )
               CALL DGEMM( 'T', 'T', N-KH, KB, KB, ONE, E(KL,KH+1),
     $                     LDE, UI, 2, ONE, DWORK(UIIPT), LDWS )
               CALL DGEMM( 'T', 'N', KB, KB, KB, ONE, UI, 2, M1, 2,
     $                     ZERO, TM, 2 )
               CALL DGEMM( 'T', 'N', N-KH, KB, KB, MONE, A(KL,KH+1),
     $                     LDA, TM, 2, ONE, DWORK(UIIPT), LDWS )
C
C              Solve generalized Sylvester equation.
C
               CALL DLASET( 'A', KB, KB, ZERO, MONE, TM, 2 )
               CALL SG03BW( 'N', N-KH, KB, A(KH+1,KH+1), LDA, M1, 2,
     $                      E(KH+1,KH+1), LDE, TM, 2, DWORK(UIIPT),
     $                      LDWS, SCALE1, INFO )
               IF ( SCALE1.NE.ONE ) THEN
                  SCALE = SCALE1*SCALE
                  CALL DLASCL( 'Upper', 0, 0, ONE, SCALE1, N, N, B, LDB,
     $                         INFO1 )
                  CALL DSCAL( 4, SCALE1, UI, 1 )
               END IF
C
C              STEP III: Form the right hand side matrix
C                        B(KH+1:N,KH+1:N) of the (smaller) Lyapunov
C                        equation to be solved during the next pass of
C                        the main loop.
C
C              Compute auxiliary matrices M3 and Y. The factorization
C              M3 = M3C * M3C**T is found by solving the symmetric
C              eigenvalue problem.
C
               CALL DLASET( 'U', 2*KB, 2*KB, ZERO, ONE, M3, 4 )
               CALL DSYRK(  'U', 'N', KB, KB, MONE, M2, 2, ONE, M3, 4 )
               CALL DGEMM(  'N', 'T', KB, KB, KB, MONE, M2, 2, M1, 2,
     $                      ZERO, M3(1,KB+1), 4 )
               CALL DSYRK(  'U', 'N', KB, KB, MONE, M1, 2, ONE,
     $                      M3(KB+1,KB+1), 4 )
               CALL DSYEVX( 'V', 'V', 'U', 2*KB, M3, 4, HALF, TWO, 1, 4,
     $                      TWO*UFLT, M, M3EW, M3C, 4, RW, 32, IW(5),
     $                      IW, INFO1 )
               IF ( INFO1.NE.0 ) THEN
                  INFO = 4
                  RETURN
               END IF
               CALL DGEMM( 'T', 'N', N-KH, KB, KB, ONE, B(KL,KH+1), LDB,
     $                     M3C, 4, ZERO, DWORK(YPT), LDWS )
               CALL DGEMM( 'T', 'T', N-KH, KB, KB, ONE, A(KL,KH+1), LDA,
     $                     UI, 2, ZERO, DWORK(WPT), LDWS )
               DO 40 I = 1, N-KH
                  CALL DGEMV( 'T', MIN( I+1, N-KH ), KB, ONE,
     $                        DWORK(UIIPT), LDWS, A(KH+1,KH+I), 1, ONE,
     $                        DWORK(WPT+I-1), LDWS )
   40          CONTINUE
               CALL DGEMM( 'N', 'N', N-KH, KB, KB, ONE, DWORK(WPT),
     $                     LDWS, M3C(KB+1,1), 4, ONE, DWORK(YPT), LDWS )
C
C              Overwrite B(KH+1:N,KH+1:N) with the triangular matrix
C              from the QR-factorization of the (N-KH+KB)-by-(N-KH)
C              matrix
C
C                          (  B(KH+1:N,KH+1:N)  )
C                          (                    )
C                          (       Y**T         ) .
C
               L = YPT - 1
               DO 80 J = 1, KB
                  DO 60 I = 1, N-KH
                     X = B(KH+I,KH+I)
                     T = DWORK(L+I)
                     CALL DLARTG( X, T, C, S, R )
                     B(KH+I,KH+I) = R
                     IF ( I.LT.N-KH )
     $                  CALL DROT( N-KH-I, B(KH+I,KH+I+1), LDB,
     $                             DWORK(L+I+1), 1, C, S )
   60             CONTINUE
                  L = L + LDWS
   80          CONTINUE
C
C              Make main diagonal elements of B(KH+1:N,KH+1:N) positive.
C
               DO 100 I = KH+1, N
                  IF ( B(I,I).LT.ZERO )
     $               CALL DSCAL( N-I+1, MONE, B(I,I), LDB )
  100          CONTINUE
C
C              Overwrite right hand side with the part of the solution
C              computed in step II.
C
               CALL DCOPY( N-KH, DWORK(UIIPT), 1, B(KL,KH+1), LDB )
               IF ( KH.GT.KL )
     $            CALL DCOPY( N-KH, DWORK(UIIPT+LDWS), 1, B(KH,KH+1),
     $                        LDB )
            END IF
C
C           Overwrite right hand side with the part of the solution
C           computed in step I.
C
            CALL DLACPY( 'U', KB, KB, UI, 2, B(KL,KL), LDB )
C
         GOTO 20
         END IF
C        END WHILE 20
C
      ELSE
C
C        Solve equation (2).
C
C        Main Loop. Compute block column U(1:KH,KL:KH). KB denotes the
C        number of columns in this block column.
C
         KL = N + 1
C        WHILE ( KL.GT.1 ) DO
  120    CONTINUE
         IF ( KL.GT.1 ) THEN
            KH = KL - 1
            IF ( KH.EQ.1 ) THEN
               KL = 1
               KB = 1
            ELSE
               IF ( A(KH,KH-1).EQ.ZERO ) THEN
                  KL = KH
                  KB = 1
               ELSE
                  KL = KH - 1
                  KB = 2
               END IF
            END IF
            KL1 = KL - 1
C
C           STEP I: Compute block U(KL:KH,KL:KH) and the auxiliary
C                   matrices M1 and M2. (For the moment the result
C                   U(KL:KH,KL:KH) is stored in UI).
C
            IF ( KB.EQ.1 ) THEN
               DELTA1 = E(KL,KL)
               T      = ABS( A(KL,KL) )
               X      = MAX( DELTA1, T )
               DELTA1 = DELTA1/X
               T      =      T/X
               IF ( DELTA1.LE.T ) THEN
                  INFO = 3
                  RETURN
               END IF
               DELTA1 = SQRT( ONE - T )*SQRT( ONE + T )*X
               T = B(KL,KL)*SMLNUM
               IF ( T.GT.DELTA1 ) THEN
                  SCALE1 = DELTA1/T
                  SCALE  = SCALE1*SCALE
                  CALL DLASCL( 'Upper', 0, 0, ONE, SCALE1, N, N, B, LDB,
     $                         INFO1 )
               END IF
C
               UI(1,1) = B(KL,KL)/DELTA1
               M1(1,1) = A(KL,KL)/E(KL,KL)
               M2(1,1) =   DELTA1/E(KL,KL)
C                
            ELSE
C
C              If a pair of complex conjugate eigenvalues occurs, apply
C              (complex) Hammarling algorithm for the 2-by-2 problem.
C
               CALL SG03BX( 'D', 'T', A(KL,KL), LDA, E(KL,KL), LDE,
     $                      B(KL,KL), LDB, UI, 2, SCALE1, M1, 2, M2, 2,
     $                      INFO )
               IF ( INFO.NE.0 )
     $            RETURN
               IF ( SCALE1.NE.ONE ) THEN
                  SCALE = SCALE1*SCALE
                  CALL DLASCL( 'Upper', 0, 0, ONE, SCALE1, N, N, B, LDB,
     $                         INFO1 )
               END IF
            END IF
C
            IF ( KL.GT.1 ) THEN
C
C              STEP II: Compute U(1:KL-1,KL:KH) by solving a generalized
C                       Sylvester equation. (For the moment the result
C                       U(1:KL-1,KL:KH) is stored in the workspace.)
C
C              Form right hand side of the Sylvester equation.
C
               CALL DGEMM( 'N', 'T', KL1, KB, KB, MONE, B(1,KL), LDB,
     $                     M2, 2, ZERO, DWORK(UIIPT), LDWS )
               CALL DGEMM( 'N', 'N', KL1, KB, KB, ONE, E(1,KL), LDE, UI,
     $                     2, ONE, DWORK(UIIPT), LDWS )
               CALL DGEMM( 'N', 'T', KB, KB, KB, ONE, UI, 2, M1, 2,
     $                     ZERO, TM, 2 )
               CALL DGEMM( 'N', 'N', KL1, KB, KB, MONE, A(1,KL), LDA,
     $                     TM, 2, ONE, DWORK(UIIPT), LDWS )
C
C              Solve generalized Sylvester equation.
C
               CALL DLASET( 'A', KB, KB, ZERO, MONE, TM, 2 )
               CALL SG03BW( 'T', KL1, KB, A, LDA, M1, 2, E, LDE, TM, 2,
     $                      DWORK(UIIPT), LDWS, SCALE1, INFO )
               IF ( SCALE1.NE.ONE ) THEN
                  SCALE = SCALE1*SCALE
                  CALL DLASCL( 'Upper', 0, 0, ONE, SCALE1, N, N, B, LDB,
     $                         INFO1 )
                  CALL DSCAL( 4, SCALE1, UI, 1 )
               END IF
C
C              STEP III: Form the right hand side matrix
C                        B(1:KL-1,1:KL-1) of the (smaller) Lyapunov
C                        equation to be solved during the next pass of
C                        the main loop.
C
C              Compute auxiliary matrices M3 and Y. The factorization
C              M3 = M3C * M3C**T is found by solving the symmetric
C              eigenvalue problem.
C
               CALL DLASET( 'U', 2*KB, 2*KB, ZERO, ONE, M3, 4 )
               CALL DSYRK(  'U', 'T', KB, KB, MONE, M2, 2, ONE, M3, 4 )
               CALL DGEMM(  'T', 'N', KB, KB, KB, MONE, M2, 2, M1, 2,
     $                      ZERO, M3(1,KB+1), 4 )
               CALL DSYRK(  'U', 'T', KB, KB, MONE, M1, 2, ONE,
     $                      M3(KB+1,KB+1), 4 )
               CALL DSYEVX( 'V', 'V', 'U', 2*KB, M3, 4, HALF, TWO, 1, 4,
     $                      TWO*UFLT, M, M3EW, M3C, 4, RW, 32, IW(5),
     $                      IW, INFO1 )
               IF ( INFO1.NE.0 ) THEN
                  INFO = 4
                  RETURN
               END IF
               CALL DGEMM( 'N', 'N', KL1, KB, KB, ONE, B(1,KL), LDB,
     $                     M3C, 4, ZERO, DWORK(YPT), LDWS )
               CALL DGEMM( 'N', 'N', KL1, KB, KB, ONE, A(1,KL), LDA, UI,
     $                     2, ZERO, DWORK(WPT), LDWS )
               DO 140 I = 1, KL1
                  CALL DGEMV( 'T', MIN( KL-I+1, KL1 ), KB, ONE,
     $                        DWORK(MAX( UIIPT, UIIPT+I-2 )), LDWS,
     $                        A(I,MAX( I-1, 1 )), LDA, ONE,
     $                        DWORK(WPT+I-1), LDWS )
  140          CONTINUE
               CALL DGEMM( 'N', 'N', KL1, KB, KB, ONE, DWORK(WPT), LDWS,
     $                     M3C(KB+1,1), 4, ONE, DWORK(YPT), LDWS )
C
C              Overwrite B(1:KL-1,1:KL-1) with the triangular matrix
C              from the RQ-factorization of the (KL-1)-by-KH matrix
C
C                          (                        )
C                          (  B(1:KL-1,1:KL-1)   Y  )
C                          (                        ).
C
               L = YPT - 1
               DO 180 J = 1, KB
                  DO 160 I = KL1, 1, -1
                     X = B(I,I)
                     T = DWORK(L+I)
                     CALL DLARTG( X, T, C, S, R )
                     B(I,I) = R
                     IF ( I.GT.1 )
     $                  CALL DROT( I-1, B(1,I), 1, DWORK(L+1), 1, C, S )
  160             CONTINUE
                  L = L + LDWS
  180          CONTINUE
C
C              Make main diagonal elements of B(1:KL-1,1:KL-1) positive.
C
               DO 200 I = 1, KL1
                  IF ( B(I,I).LT.ZERO )
     $               CALL DSCAL( I, MONE, B(1,I), 1 )
  200          CONTINUE
C
C              Overwrite right hand side with the part of the solution
C              computed in step II.
C
               CALL DLACPY( 'A', KL1, KB, DWORK(UIIPT), LDWS, B(1,KL),
     $                      LDB )
C
            END IF
C
C           Overwrite right hand side with the part of the solution
C           computed in step I.
C
            CALL DLACPY( 'U', KB, KB, UI, 2, B(KL,KL), LDB )
C
         GOTO 120
         END IF
C        END WHILE 120
C
      END IF
C
      RETURN
C *** Last line of SG03BU ***
      END