-
Notifications
You must be signed in to change notification settings - Fork 453
Expand file tree
/
Copy pathmateqn.py
More file actions
659 lines (518 loc) · 21.4 KB
/
mateqn.py
File metadata and controls
659 lines (518 loc) · 21.4 KB
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
# mateqn.py - matrix equation solvers (Lyapunov, Riccati)
#
# Initial author: Bjorn Olofsson
# Creation date: 2011
"""Matrix equation solvers (Lyapunov, Riccati).
This module contains implementation of the functions lyap, dlyap, care
and dare for solution of Lyapunov and Riccati equations.
"""
import warnings
import numpy as np
import scipy as sp
from numpy import eye, finfo, inexact
from scipy.linalg import eigvals, solve
from .exception import ControlArgument, ControlDimension, ControlSlycot, \
slycot_check
# Make sure we have access to the right Slycot routines
try:
from slycot.exceptions import SlycotResultWarning
except ImportError:
SlycotResultWarning = UserWarning
try:
from slycot import sb03md57
# wrap without the deprecation warning
def sb03md(n, C, A, U, dico, job='X', fact='N', trana='N', ldwork=None):
ret = sb03md57(A, U, C, dico, job, fact, trana, ldwork)
return ret[2:]
except ImportError:
try:
from slycot import sb03md
except ImportError:
sb03md = None
try:
from slycot import sb04md
except ImportError:
sb04md = None
try:
from slycot import sb04qd
except ImportError:
sb0qmd = None
try:
from slycot import sg03ad
except ImportError:
sb04ad = None
__all__ = ['lyap', 'dlyap', 'dare', 'care']
#
# Lyapunov equation solvers lyap and dlyap
#
def lyap(A, Q, C=None, E=None, method=None):
"""Solves the continuous-time Lyapunov equation.
X = lyap(A, Q) solves
:math:`A X + X A^T + Q = 0`
where A and Q are square matrices of the same dimension. Q must be
symmetric.
X = lyap(A, Q, C) solves the Sylvester equation
:math:`A X + X Q + C = 0`
where A and Q are square matrices.
X = lyap(A, Q, None, E) solves the generalized continuous-time
Lyapunov equation
:math:`A X E^T + E X A^T + Q = 0`
where Q is a symmetric matrix and A, Q and E are square matrices of the
same dimension.
Parameters
----------
A, Q : 2D array_like
Input matrices for the Lyapunov or Sylvestor equation.
C : 2D array_like, optional
If present, solve the Sylvester equation.
E : 2D array_like, optional
If present, solve the generalized Lyapunov equation.
method : str, optional
Set the method used for computing the result. Current methods are
'slycot' and 'scipy'. If set to None (default), try 'slycot' first
and then 'scipy'.
Returns
-------
X : 2D array
Solution to the Lyapunov or Sylvester equation.
"""
# Decide what method to use
method = _slycot_or_scipy(method)
if method == 'slycot':
if sb03md is None:
raise ControlSlycot("Can't find slycot module 'sb03md'")
if sb04md is None:
raise ControlSlycot("Can't find slycot module 'sb04md'")
# Reshape input arrays
A = np.array(A, ndmin=2)
Q = np.array(Q, ndmin=2)
if C is not None:
C = np.array(C, ndmin=2)
if E is not None:
E = np.array(E, ndmin=2)
# Determine main dimensions
n = A.shape[0]
m = Q.shape[0]
# Check to make sure input matrices are the right shape and type
_check_shape(A, n, n, square=True, name="A")
# Solve standard Lyapunov equation
if C is None and E is None:
# Check to make sure input matrices are the right shape and type
_check_shape(Q, n, n, square=True, symmetric=True, name="Q")
if method == 'scipy':
# Solve the Lyapunov equation using SciPy
return sp.linalg.solve_continuous_lyapunov(A, -Q)
# Solve the Lyapunov equation by calling Slycot function sb03md
with warnings.catch_warnings():
warnings.simplefilter("error", category=SlycotResultWarning)
X, scale, sep, ferr, w = \
sb03md(n, -Q, A, eye(n, n), 'C', trana='T')
# Solve the Sylvester equation
elif C is not None and E is None:
# Check to make sure input matrices are the right shape and type
_check_shape(Q, m, m, square=True, name="Q")
_check_shape(C, n, m, name="C")
if method == 'scipy':
# Solve the Sylvester equation using SciPy
return sp.linalg.solve_sylvester(A, Q, -C)
# Solve the Sylvester equation by calling the Slycot function sb04md
X = sb04md(n, m, A, Q, -C)
# Solve the generalized Lyapunov equation
elif C is None and E is not None:
# Check to make sure input matrices are the right shape and type
_check_shape(Q, n, n, square=True, symmetric=True, name="Q")
_check_shape(E, n, n, square=True, name="E")
if method == 'scipy':
raise ControlArgument(
"method='scipy' not valid for generalized Lyapunov equation")
# Make sure we have access to the write Slycot routine
try:
from slycot import sg03ad
except ImportError:
raise ControlSlycot("Can't find slycot module 'sg03ad'")
# Solve the generalized Lyapunov equation by calling Slycot
# function sg03ad
with warnings.catch_warnings():
warnings.simplefilter("error", category=SlycotResultWarning)
A, E, Q, Z, X, scale, sep, ferr, alphar, alphai, beta = \
sg03ad('C', 'B', 'N', 'T', 'L', n,
A, E, eye(n, n), eye(n, n), -Q)
# Invalid set of input parameters (C and E specified)
else:
raise ControlArgument("Invalid set of input parameters")
return X
def dlyap(A, Q, C=None, E=None, method=None):
"""Solves the discrete-time Lyapunov equation.
X = dlyap(A, Q) solves
:math:`A X A^T - X + Q = 0`
where A and Q are square matrices of the same dimension. Further
Q must be symmetric.
dlyap(A, Q, C) solves the Sylvester equation
:math:`A X Q^T - X + C = 0`
where A and Q are square matrices.
dlyap(A, Q, None, E) solves the generalized discrete-time Lyapunov
equation
:math:`A X A^T - E X E^T + Q = 0`
where Q is a symmetric matrix and A, Q and E are square matrices of the
same dimension.
Parameters
----------
A, Q : 2D array_like
Input matrices for the Lyapunov or Sylvestor equation.
C : 2D array_like, optional
If present, solve the Sylvester equation.
E : 2D array_like, optional
If present, solve the generalized Lyapunov equation.
method : str, optional
Set the method used for computing the result. Current methods are
'slycot' and 'scipy'. If set to None (default), try 'slycot' first
and then 'scipy'.
Returns
-------
X : 2D array (or matrix)
Solution to the Lyapunov or Sylvester equation.
"""
# Decide what method to use
method = _slycot_or_scipy(method)
if method == 'slycot':
# Make sure we have access to the right slycot routines
if sb03md is None:
raise ControlSlycot("Can't find slycot module 'sb03md'")
if sb04qd is None:
raise ControlSlycot("Can't find slycot module 'sb04qd'")
if sg03ad is None:
raise ControlSlycot("Can't find slycot module 'sg03ad'")
# Reshape input arrays
A = np.array(A, ndmin=2)
Q = np.array(Q, ndmin=2)
if C is not None:
C = np.array(C, ndmin=2)
if E is not None:
E = np.array(E, ndmin=2)
# Determine main dimensions
n = A.shape[0]
m = Q.shape[0]
# Check to make sure input matrices are the right shape and type
_check_shape(A, n, n, square=True, name="A")
# Solve standard Lyapunov equation
if C is None and E is None:
# Check to make sure input matrices are the right shape and type
_check_shape(Q, n, n, square=True, symmetric=True, name="Q")
if method == 'scipy':
# Solve the Lyapunov equation using SciPy
return sp.linalg.solve_discrete_lyapunov(A, Q)
# Solve the Lyapunov equation by calling the Slycot function sb03md
with warnings.catch_warnings():
warnings.simplefilter("error", category=SlycotResultWarning)
X, scale, sep, ferr, w = \
sb03md(n, -Q, A, eye(n, n), 'D', trana='T')
# Solve the Sylvester equation
elif C is not None and E is None:
# Check to make sure input matrices are the right shape and type
_check_shape(Q, m, m, square=True, name="Q")
_check_shape(C, n, m, name="C")
if method == 'scipy':
raise ControlArgument(
"method='scipy' not valid for Sylvester equation")
# Solve the Sylvester equation by calling Slycot function sb04qd
X = sb04qd(n, m, -A, Q.T, C)
# Solve the generalized Lyapunov equation
elif C is None and E is not None:
# Check to make sure input matrices are the right shape and type
_check_shape(Q, n, n, square=True, symmetric=True, name="Q")
_check_shape(E, n, n, square=True, name="E")
if method == 'scipy':
raise ControlArgument(
"method='scipy' not valid for generalized Lyapunov equation")
# Solve the generalized Lyapunov equation by calling Slycot
# function sg03ad
with warnings.catch_warnings():
warnings.simplefilter("error", category=SlycotResultWarning)
A, E, Q, Z, X, scale, sep, ferr, alphar, alphai, beta = \
sg03ad('D', 'B', 'N', 'T', 'L', n,
A, E, eye(n, n), eye(n, n), -Q)
# Invalid set of input parameters (C and E specified)
else:
raise ControlArgument("Invalid set of input parameters")
return X
#
# Riccati equation solvers care and dare
#
def care(A, B, Q, R=None, S=None, E=None, stabilizing=True, method=None,
_As="A", _Bs="B", _Qs="Q", _Rs="R", _Ss="S", _Es="E"):
"""Solves the continuous-time algebraic Riccati equation.
X, L, G = care(A, B, Q, R=None) solves
:math:`A^T X + X A - X B R^{-1} B^T X + Q = 0`
where A and Q are square matrices of the same dimension. Further,
Q and R are a symmetric matrices. If R is None, it is set to the
identity matrix. The function returns the solution X, the gain
matrix G = B^T X and the closed loop eigenvalues L, i.e., the
eigenvalues of A - B G.
X, L, G = care(A, B, Q, R, S, E) solves the generalized
continuous-time algebraic Riccati equation
:math:`A^T X E + E^T X A - (E^T X B + S) R^{-1} (B^T X E + S^T) + Q = 0`
where A, Q and E are square matrices of the same dimension. Further, Q
and R are symmetric matrices. If R is None, it is set to the identity
matrix. The function returns the solution X, the gain matrix G = R^-1
(B^T X E + S^T) and the closed loop eigenvalues L, i.e., the eigenvalues
of A - B G , E.
Parameters
----------
A, B, Q : 2D array_like
Input matrices for the Riccati equation.
R, S, E : 2D array_like, optional
Input matrices for generalized Riccati equation.
method : str, optional
Set the method used for computing the result. Current methods are
'slycot' and 'scipy'. If set to None (default), try 'slycot' first
and then 'scipy'.
stabilizing : bool, optional
If `method` is 'slycot', unstabilized eigenvalues will be returned
in the initial elements of `L`. Not supported for 'scipy'.
Returns
-------
X : 2D array (or matrix)
Solution to the Riccati equation.
L : 1D array
Closed loop eigenvalues.
G : 2D array (or matrix)
Gain matrix.
"""
# Decide what method to use
method = _slycot_or_scipy(method)
# Reshape input arrays
A = np.array(A, ndmin=2)
B = np.array(B, ndmin=2)
Q = np.array(Q, ndmin=2)
R = np.eye(B.shape[1]) if R is None else np.array(R, ndmin=2)
if S is not None:
S = np.array(S, ndmin=2)
if E is not None:
E = np.array(E, ndmin=2)
# Determine main dimensions
n = A.shape[0]
m = B.shape[1]
# Check to make sure input matrices are the right shape and type
_check_shape(A, n, n, square=True, name=_As)
_check_shape(B, n, m, name=_Bs)
_check_shape(Q, n, n, square=True, symmetric=True, name=_Qs)
_check_shape(R, m, m, square=True, symmetric=True, name=_Rs)
# Solve the standard algebraic Riccati equation
if S is None and E is None:
# See if we should solve this using SciPy
if method == 'scipy':
if not stabilizing:
raise ControlArgument(
"method='scipy' not valid when stabilizing is not True")
X = sp.linalg.solve_continuous_are(A, B, Q, R)
K = np.linalg.solve(R, B.T @ X)
E, _ = np.linalg.eig(A - B @ K)
return X, E, K
# Make sure we can import required Slycot routines
try:
from slycot import sb02md
except ImportError:
raise ControlSlycot("Can't find slycot module 'sb02md'")
try:
from slycot import sb02mt
except ImportError:
raise ControlSlycot("Can't find slycot module 'sb02mt'")
# Solve the standard algebraic Riccati equation by calling Slycot
# functions sb02mt and sb02md
A_b, B_b, Q_b, R_b, L_b, ipiv, oufact, G = sb02mt(n, m, B, R)
sort = 'S' if stabilizing else 'U'
X, rcond, w, S_o, U, A_inv = sb02md(n, A, G, Q, 'C', sort=sort)
# Calculate the gain matrix G
G = solve(R, B.T) @ X
# Return the solution X, the closed-loop eigenvalues L and
# the gain matrix G
return X, w[:n], G
# Solve the generalized algebraic Riccati equation
else:
# Initialize optional matrices
S = np.zeros((n, m)) if S is None else np.array(S, ndmin=2)
E = np.eye(A.shape[0]) if E is None else np.array(E, ndmin=2)
# Check to make sure input matrices are the right shape and type
_check_shape(E, n, n, square=True, name=_Es)
_check_shape(S, n, m, name=_Ss)
# See if we should solve this using SciPy
if method == 'scipy':
if not stabilizing:
raise ControlArgument(
"method='scipy' not valid when stabilizing is not True")
X = sp.linalg.solve_continuous_are(A, B, Q, R, s=S, e=E)
K = np.linalg.solve(R, B.T @ X @ E + S.T)
eigs, _ = sp.linalg.eig(A - B @ K, E)
return X, eigs, K
# Make sure we can find the required Slycot routine
try:
from slycot import sg02ad
except ImportError:
raise ControlSlycot("Can't find slycot module sg02ad")
# Solve the generalized algebraic Riccati equation by calling the
# Slycot function sg02ad
with warnings.catch_warnings():
sort = 'S' if stabilizing else 'U'
warnings.simplefilter("error", category=SlycotResultWarning)
rcondu, X, alfar, alfai, beta, S_o, T, U, iwarn = \
sg02ad('C', 'B', 'N', 'U', 'N', 'N', sort,
'R', n, m, 0, A, E, B, Q, R, S)
# Calculate the closed-loop eigenvalues L
L = np.array([(alfar[i] + alfai[i]*1j) / beta[i] for i in range(n)])
# Calculate the gain matrix G
G = solve(R, B.T @ X @ E + S.T)
# Return the solution X, the closed-loop eigenvalues L and
# the gain matrix G
return X, L, G
def dare(A, B, Q, R, S=None, E=None, stabilizing=True, method=None,
_As="A", _Bs="B", _Qs="Q", _Rs="R", _Ss="S", _Es="E"):
"""Solves the discrete-time algebraic Riccati equation.
X, L, G = dare(A, B, Q, R) solves
:math:`A^T X A - X - A^T X B (B^T X B + R)^{-1} B^T X A + Q = 0`
where A and Q are square matrices of the same dimension. Further, Q
is a symmetric matrix. The function returns the solution X, the gain
matrix G = (B^T X B + R)^-1 B^T X A and the closed loop eigenvalues L,
i.e., the eigenvalues of A - B G.
X, L, G = dare(A, B, Q, R, S, E) solves the generalized discrete-time
algebraic Riccati equation
: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`
where A, Q and E are square matrices of the same dimension. Further, Q
and R are symmetric matrices. If R is None, it is set to the identity
matrix. The function returns the solution X, the gain matrix :math:`G =
(B^T X B + R)^{-1} (B^T X A + S^T)` and the closed loop eigenvalues L,
i.e., the (generalized) eigenvalues of A - B G (with respect to E, if
specified).
Parameters
----------
A, B, Q : 2D arrays
Input matrices for the Riccati equation.
R, S, E : 2D arrays, optional
Input matrices for generalized Riccati equation.
method : str, optional
Set the method used for computing the result. Current methods are
'slycot' and 'scipy'. If set to None (default), try 'slycot' first
and then 'scipy'.
stabilizing : bool, optional
If `method` is 'slycot', unstabilized eigenvalues will be returned
in the initial elements of `L`. Not supported for 'scipy'.
Returns
-------
X : 2D array (or matrix)
Solution to the Riccati equation.
L : 1D array
Closed loop eigenvalues.
G : 2D array (or matrix)
Gain matrix.
"""
# Decide what method to use
method = _slycot_or_scipy(method)
# Reshape input arrays
A = np.array(A, ndmin=2)
B = np.array(B, ndmin=2)
Q = np.array(Q, ndmin=2)
R = np.eye(B.shape[1]) if R is None else np.array(R, ndmin=2)
if S is not None:
S = np.array(S, ndmin=2)
if E is not None:
E = np.array(E, ndmin=2)
# Determine main dimensions
n = A.shape[0]
m = B.shape[1]
# Check to make sure input matrices are the right shape and type
_check_shape(A, n, n, square=True, name=_As)
_check_shape(B, n, m, name=_Bs)
_check_shape(Q, n, n, square=True, symmetric=True, name=_Qs)
_check_shape(R, m, m, square=True, symmetric=True, name=_Rs)
if E is not None:
_check_shape(E, n, n, square=True, name=_Es)
if S is not None:
_check_shape(S, n, m, name=_Ss)
# Figure out how to solve the problem
if method == 'scipy':
if not stabilizing:
raise ControlArgument(
"method='scipy' not valid when stabilizing is not True")
X = sp.linalg.solve_discrete_are(A, B, Q, R, e=E, s=S)
if S is None:
G = solve(B.T @ X @ B + R, B.T @ X @ A)
else:
G = solve(B.T @ X @ B + R, B.T @ X @ A + S.T)
if E is None:
L = eigvals(A - B @ G)
else:
L, _ = sp.linalg.eig(A - B @ G, E)
return X, L, G
# Make sure we can import required Slycot routine
try:
from slycot import sg02ad
except ImportError:
raise ControlSlycot("Can't find slycot module sg02ad")
# Initialize optional matrices
S = np.zeros((n, m)) if S is None else np.array(S, ndmin=2)
E = np.eye(A.shape[0]) if E is None else np.array(E, ndmin=2)
# Solve the generalized algebraic Riccati equation by calling the
# Slycot function sg02ad
sort = 'S' if stabilizing else 'U'
with warnings.catch_warnings():
warnings.simplefilter("error", category=SlycotResultWarning)
rcondu, X, alfar, alfai, beta, S_o, T, U, iwarn = \
sg02ad('D', 'B', 'N', 'U', 'N', 'N', sort,
'R', n, m, 0, A, E, B, Q, R, S)
# Calculate the closed-loop eigenvalues L
L = np.array([(alfar[i] + alfai[i]*1j) / beta[i] for i in range(n)])
# Calculate the gain matrix G
G = solve(B.T @ X @ B + R, B.T @ X @ A + S.T)
# Return the solution X, the closed-loop eigenvalues L and
# the gain matrix G
return X, L, G
# Utility function to decide on method to use
def _slycot_or_scipy(method):
if method == 'slycot' or (method is None and slycot_check()):
return 'slycot'
elif method == 'scipy' or (method is None and not slycot_check()):
return 'scipy'
else:
raise ControlArgument("Unknown method %s" % method)
# Utility function to check matrix dimensions
def _check_shape(M, n, m, square=False, symmetric=False, name="??"):
"""Check the shape and properties of a 2D array.
This function can be used to check to make sure a 2D array_like has the
right shape, along with other properties. If not, an appropriate error
message is generated.
Parameters
----------
M : array_like
Array to be checked.
n : int
Expected number of rows.
m : int
Expected number of columns.
square : bool, optional
If True, check to make sure the matrix is square.
symmetric : bool, optional
If True, check to make sure the matrix is symmetric.
name : str
Name of the matrix (for use in error messages).
Returns
-------
M : 2D array
Input array, converted to 2D if needed.
"""
M = np.atleast_2d(M)
if (square or symmetric) and M.shape[0] != M.shape[1]:
raise ControlDimension("%s must be a square matrix" % name)
if symmetric and not _is_symmetric(M):
raise ControlArgument("%s must be a symmetric matrix" % name)
if M.shape[0] != n or M.shape[1] != m:
raise ControlDimension(
f"Incompatible dimensions of {name} matrix; "
f"expected ({n}, {m}) but found {M.shape}")
return M
# Utility function to check if a matrix is symmetric
def _is_symmetric(M):
M = np.atleast_2d(M)
if isinstance(M[0, 0], inexact):
eps = finfo(M.dtype).eps
return ((M - M.T) < eps).all()
else:
return (M == M.T).all()