3
3
# Implementation of the functions lyap, dlyap, care and dare
4
4
# for solution of Lyapunov and Riccati equations.
5
5
#
6
- # Author : Bjorn Olofsson
6
+ # Original author : Bjorn Olofsson
7
7
8
8
# Copyright (c) 2011, All rights reserved.
9
9
@@ -162,6 +162,7 @@ def lyap(A, Q, C=None, E=None, method=None):
162
162
_check_shape ("Q" , Q , n , n , square = True , symmetric = True )
163
163
164
164
if method == 'scipy' :
165
+ # Solve the Lyapunov equation using SciPy
165
166
return sp .linalg .solve_continuous_lyapunov (A , - Q )
166
167
167
168
# Solve the Lyapunov equation by calling Slycot function sb03md
@@ -177,6 +178,7 @@ def lyap(A, Q, C=None, E=None, method=None):
177
178
_check_shape ("C" , C , n , m )
178
179
179
180
if method == 'scipy' :
181
+ # Solve the Sylvester equation using SciPy
180
182
return sp .linalg .solve_sylvester (A , Q , - C )
181
183
182
184
# Solve the Sylvester equation by calling the Slycot function sb04md
@@ -293,6 +295,7 @@ def dlyap(A, Q, C=None, E=None, method=None):
293
295
_check_shape ("Q" , Q , n , n , square = True , symmetric = True )
294
296
295
297
if method == 'scipy' :
298
+ # Solve the Lyapunov equation using SciPy
296
299
return sp .linalg .solve_discrete_lyapunov (A , Q )
297
300
298
301
# Solve the Lyapunov equation by calling the Slycot function sb03md
@@ -343,7 +346,8 @@ def dlyap(A, Q, C=None, E=None, method=None):
343
346
# Riccati equation solvers care and dare
344
347
#
345
348
346
- def care (A , B , Q , R = None , S = None , E = None , stabilizing = True , method = None ):
349
+ def care (A , B , Q , R = None , S = None , E = None , stabilizing = True , method = None ,
350
+ A_s = "A" , B_s = "B" , Q_s = "Q" , R_s = "R" , S_s = "S" , E_s = "E" ):
347
351
"""X, L, G = care(A, B, Q, R=None) solves the continuous-time
348
352
algebraic Riccati equation
349
353
@@ -395,24 +399,6 @@ def care(A, B, Q, R=None, S=None, E=None, stabilizing=True, method=None):
395
399
# Decide what method to use
396
400
method = _slycot_or_scipy (method )
397
401
398
- if method == 'slycot' :
399
- # Make sure we can import required slycot routines
400
- try :
401
- from slycot import sb02md
402
- except ImportError :
403
- raise ControlSlycot ("Can't find slycot module 'sb02md'" )
404
-
405
- try :
406
- from slycot import sb02mt
407
- except ImportError :
408
- raise ControlSlycot ("Can't find slycot module 'sb02mt'" )
409
-
410
- # Make sure we can find the required slycot routine
411
- try :
412
- from slycot import sg02ad
413
- except ImportError :
414
- raise ControlSlycot ("Can't find slycot module 'sg02ad'" )
415
-
416
402
# Reshape input arrays
417
403
A = np .array (A , ndmin = 2 )
418
404
B = np .array (B , ndmin = 2 )
@@ -428,10 +414,10 @@ def care(A, B, Q, R=None, S=None, E=None, stabilizing=True, method=None):
428
414
m = B .shape [1 ]
429
415
430
416
# Check to make sure input matrices are the right shape and type
431
- _check_shape ("A" , A , n , n , square = True )
432
- _check_shape ("B" , B , n , m )
433
- _check_shape ("Q" , Q , n , n , square = True , symmetric = True )
434
- _check_shape ("R" , R , m , m , square = True , symmetric = True )
417
+ _check_shape (A_s , A , n , n , square = True )
418
+ _check_shape (B_s , B , n , m )
419
+ _check_shape (Q_s , Q , n , n , square = True , symmetric = True )
420
+ _check_shape (R_s , R , m , m , square = True , symmetric = True )
435
421
436
422
# Solve the standard algebraic Riccati equation
437
423
if S is None and E is None :
@@ -446,9 +432,16 @@ def care(A, B, Q, R=None, S=None, E=None, stabilizing=True, method=None):
446
432
E , _ = np .linalg .eig (A - B @ K )
447
433
return _ssmatrix (X ), E , _ssmatrix (K )
448
434
449
- # Create back-up of arrays needed for later computations
450
- R_ba = copy (R )
451
- B_ba = copy (B )
435
+ # Make sure we can import required slycot routines
436
+ try :
437
+ from slycot import sb02md
438
+ except ImportError :
439
+ raise ControlSlycot ("Can't find slycot module 'sb02md'" )
440
+
441
+ try :
442
+ from slycot import sb02mt
443
+ except ImportError :
444
+ raise ControlSlycot ("Can't find slycot module 'sb02mt'" )
452
445
453
446
# Solve the standard algebraic Riccati equation by calling Slycot
454
447
# functions sb02mt and sb02md
@@ -458,7 +451,7 @@ def care(A, B, Q, R=None, S=None, E=None, stabilizing=True, method=None):
458
451
X , rcond , w , S_o , U , A_inv = sb02md (n , A , G , Q , 'C' , sort = sort )
459
452
460
453
# Calculate the gain matrix G
461
- G = solve (R_ba , B_ba .T ) @ X
454
+ G = solve (R , B .T ) @ X
462
455
463
456
# Return the solution X, the closed-loop eigenvalues L and
464
457
# the gain matrix G
@@ -471,8 +464,8 @@ def care(A, B, Q, R=None, S=None, E=None, stabilizing=True, method=None):
471
464
E = np .eye (A .shape [0 ]) if E is None else np .array (E , ndmin = 2 )
472
465
473
466
# Check to make sure input matrices are the right shape and type
474
- _check_shape ("E" , E , n , n , square = True )
475
- _check_shape ("S" , S , n , m )
467
+ _check_shape (E_s , E , n , n , square = True )
468
+ _check_shape (S_s , S , n , m )
476
469
477
470
# See if we should solve this using SciPy
478
471
if method == 'scipy' :
@@ -485,11 +478,11 @@ def care(A, B, Q, R=None, S=None, E=None, stabilizing=True, method=None):
485
478
eigs , _ = sp .linalg .eig (A - B @ K , E )
486
479
return _ssmatrix (X ), eigs , _ssmatrix (K )
487
480
488
- # Create back-up of arrays needed for later computations
489
- R_b = copy ( R )
490
- B_b = copy ( B )
491
- E_b = copy ( E )
492
- S_b = copy ( S )
481
+ # Make sure we can find the required slycot routine
482
+ try :
483
+ from slycot import sg02ad
484
+ except ImportError :
485
+ raise ControlSlycot ( "Can't find slycot module 'sg02ad'" )
493
486
494
487
# Solve the generalized algebraic Riccati equation by calling the
495
488
# Slycot function sg02ad
@@ -504,14 +497,15 @@ def care(A, B, Q, R=None, S=None, E=None, stabilizing=True, method=None):
504
497
L = np .array ([(alfar [i ] + alfai [i ]* 1j ) / beta [i ] for i in range (n )])
505
498
506
499
# Calculate the gain matrix G
507
- G = solve (R_b , B_b .T @ X @ E_b + S_b .T )
500
+ G = solve (R , B .T @ X @ E + S .T )
508
501
509
502
# Return the solution X, the closed-loop eigenvalues L and
510
503
# the gain matrix G
511
504
return _ssmatrix (X ), L , _ssmatrix (G )
512
505
513
- def dare (A , B , Q , R , S = None , E = None , stabilizing = True , method = None ):
514
- """(X, L, G) = dare(A, B, Q, R) solves the discrete-time algebraic Riccati
506
+ def dare (A , B , Q , R , S = None , E = None , stabilizing = True , method = None ,
507
+ A_s = "A" , B_s = "B" , Q_s = "Q" , R_s = "R" , S_s = "S" , E_s = "E" ):
508
+ """X, L, G = dare(A, B, Q, R) solves the discrete-time algebraic Riccati
515
509
equation
516
510
517
511
:math:`A^T X A - X - A^T X B (B^T X B + R)^{-1} B^T X A + Q = 0`
@@ -521,15 +515,17 @@ def dare(A, B, Q, R, S=None, E=None, stabilizing=True, method=None):
521
515
matrix G = (B^T X B + R)^-1 B^T X A and the closed loop eigenvalues L,
522
516
i.e., the eigenvalues of A - B G.
523
517
524
- ( X, L, G) = dare(A, B, Q, R, S, E) solves the generalized discrete-time
518
+ X, L, G = dare(A, B, Q, R, S, E) solves the generalized discrete-time
525
519
algebraic Riccati equation
526
520
527
521
:math:`A^T X A - E^T X E - (A^T X B + S) (B^T X B + R)^{-1} (B^T X A + S^T) + Q = 0`
528
522
529
- where A, Q and E are square matrices of the same dimension. Further, Q and
530
- R are symmetric matrices. The function returns the solution X, the gain
531
- matrix :math:`G = (B^T X B + R)^{-1} (B^T X A + S^T)` and the closed loop
532
- eigenvalues L, i.e., the eigenvalues of A - B G , E.
523
+ where A, Q and E are square matrices of the same dimension. Further, Q
524
+ and R are symmetric matrices. If R is None, it is set to the identity
525
+ matrix. The function returns the solution X, the gain matrix :math:`G =
526
+ (B^T X B + R)^{-1} (B^T X A + S^T)` and the closed loop eigenvalues L,
527
+ i.e., the (generalized) eigenvalues of A - B G (with respect to E, if
528
+ specified).
533
529
534
530
Parameters
535
531
----------
@@ -575,87 +571,43 @@ def dare(A, B, Q, R, S=None, E=None, stabilizing=True, method=None):
575
571
m = B .shape [1 ]
576
572
577
573
# Check to make sure input matrices are the right shape and type
578
- _check_shape ("A" , A , n , n , square = True )
574
+ _check_shape (A_s , A , n , n , square = True )
575
+ _check_shape (B_s , B , n , m )
576
+ _check_shape (Q_s , Q , n , n , square = True , symmetric = True )
577
+ _check_shape (R_s , R , m , m , square = True , symmetric = True )
578
+ if E is not None :
579
+ _check_shape (E_s , E , n , n , square = True )
580
+ if S is not None :
581
+ _check_shape (S_s , S , n , m )
579
582
580
583
# Figure out how to solve the problem
581
- if method == 'scipy' and not stabilizing :
582
- raise ControlArgument (
583
- "method='scipy' not valid when stabilizing is not True" )
584
-
585
- elif method == 'slycot' :
586
- return _dare_slycot (A , B , Q , R , S , E , stabilizing )
584
+ if method == 'scipy' :
585
+ if not stabilizing :
586
+ raise ControlArgument (
587
+ "method='scipy' not valid when stabilizing is not True" )
587
588
588
- else :
589
- _check_shape ("B" , B , n , m )
590
- _check_shape ("Q" , Q , n , n , square = True , symmetric = True )
591
- _check_shape ("R" , R , m , m , square = True , symmetric = True )
592
- if E is not None :
593
- _check_shape ("E" , E , n , n , square = True )
594
- if S is not None :
595
- _check_shape ("S" , S , n , m )
596
-
597
- Rmat = _ssmatrix (R )
598
- Qmat = _ssmatrix (Q )
599
- X = sp .linalg .solve_discrete_are (A , B , Qmat , Rmat , e = E , s = S )
589
+ X = sp .linalg .solve_discrete_are (A , B , Q , R , e = E , s = S )
600
590
if S is None :
601
- G = solve (B .T @ X @ B + Rmat , B .T @ X @ A )
591
+ G = solve (B .T @ X @ B + R , B .T @ X @ A )
602
592
else :
603
- G = solve (B .T @ X @ B + Rmat , B .T @ X @ A + S .T )
593
+ G = solve (B .T @ X @ B + R , B .T @ X @ A + S .T )
604
594
if E is None :
605
595
L = eigvals (A - B @ G )
606
596
else :
607
597
L , _ = sp .linalg .eig (A - B @ G , E )
608
598
609
599
return _ssmatrix (X ), L , _ssmatrix (G )
610
600
611
-
612
- def _dare_slycot (A , B , Q , R , S = None , E = None , stabilizing = True ):
613
601
# Make sure we can import required slycot routine
614
- try :
615
- from slycot import sb02md
616
- except ImportError :
617
- raise ControlSlycot ("Can't find slycot module 'sb02md'" )
618
-
619
- try :
620
- from slycot import sb02mt
621
- except ImportError :
622
- raise ControlSlycot ("Can't find slycot module 'sb02mt'" )
623
-
624
- # Make sure we can find the required slycot routine
625
602
try :
626
603
from slycot import sg02ad
627
604
except ImportError :
628
605
raise ControlSlycot ("Can't find slycot module 'sg02ad'" )
629
606
630
- # Reshape input arrays
631
- A = np .array (A , ndmin = 2 )
632
- B = np .array (B , ndmin = 2 )
633
- Q = np .array (Q , ndmin = 2 )
634
- R = np .eye (B .shape [1 ]) if R is None else np .array (R , ndmin = 2 )
635
-
636
- # Determine main dimensions
637
- n = A .shape [0 ]
638
- m = B .shape [1 ]
639
-
640
607
# Initialize optional matrices
641
608
S = np .zeros ((n , m )) if S is None else np .array (S , ndmin = 2 )
642
609
E = np .eye (A .shape [0 ]) if E is None else np .array (E , ndmin = 2 )
643
610
644
- # Check to make sure input matrices are the right shape and type
645
- _check_shape ("A" , A , n , n , square = True )
646
- _check_shape ("B" , B , n , m )
647
- _check_shape ("Q" , Q , n , n , square = True , symmetric = True )
648
- _check_shape ("R" , R , m , m , square = True , symmetric = True )
649
- _check_shape ("E" , E , n , n , square = True )
650
- _check_shape ("S" , S , n , m )
651
-
652
- # Create back-up of arrays needed for later computations
653
- A_b = copy (A )
654
- R_b = copy (R )
655
- B_b = copy (B )
656
- E_b = copy (E )
657
- S_b = copy (S )
658
-
659
611
# Solve the generalized algebraic Riccati equation by calling the
660
612
# Slycot function sg02ad
661
613
sort = 'S' if stabilizing else 'U'
@@ -669,7 +621,7 @@ def _dare_slycot(A, B, Q, R, S=None, E=None, stabilizing=True):
669
621
L = np .array ([(alfar [i ] + alfai [i ]* 1j ) / beta [i ] for i in range (n )])
670
622
671
623
# Calculate the gain matrix G
672
- G = solve (B_b .T @ X @ B_b + R_b , B_b .T @ X @ A_b + S_b .T )
624
+ G = solve (B .T @ X @ B + R , B .T @ X @ A + S .T )
673
625
674
626
# Return the solution X, the closed-loop eigenvalues L and
675
627
# the gain matrix G
0 commit comments