@@ -197,10 +197,48 @@ def modred(sys, ELIM, method='matchdc'):
197197 rsys = StateSpace (Ar ,Br ,Cr ,Dr )
198198 return rsys
199199
200+ def stabsep (T_schur , Z_schur , sys , ldim , no_in , no_out ):
201+ """
202+ Performs stable/unstabe decomposition of sys after Schur forms have been computed for system matrix.
203+
204+ Reference: Hsu,C.S., and Hou,D., 1991,
205+ Reducing unstable linear control systems via real Schur transformation.
206+ Electronics Letters, 27, 984-986.
207+
208+ """
209+ #Author: M. Clement (mdclemen@eng.ucsd.edu) 2016
210+ As = np .asmatrix (T_schur )
211+ Bs = Z_schur .T * sys .B
212+ Cs = sys .C * Z_schur
213+ #from ref 1 eq(1) As = [A_ Ac], Bs = [B_], and Cs = [C_ C+]; _ denotes stable subsystem
214+ # [0 A+] [B+]
215+ A_ = As [0 :ldim ,0 :ldim ]
216+ Ac = As [0 :ldim ,ldim ::]
217+ Ap = As [ldim ::,ldim ::]
218+
219+ B_ = Bs [0 :ldim ,:]
220+ Bp = Bs [ldim ::,:]
221+
222+ C_ = Cs [:,0 :ldim ]
223+ Cp = Cs [:,ldim ::]
224+ #do some more tricky math IAW ref 1 eq(3)
225+ B_tilde = np .bmat ([[B_ , Ac ]])
226+ D_tilde = np .bmat ([[np .zeros ((no_out , no_in )), Cp ]])
227+
228+ return A_ , B_tilde , C_ , D_tilde , Ap , Bp , Cp
229+
230+
200231def balred (sys , orders , method = 'truncate' ):
201232 """
202233 Balanced reduced order model of sys of a given order.
203234 States are eliminated based on Hankel singular value.
235+ If sys has unstable modes, they are removed, the
236+ balanced realization is done on the stable part, then
237+ reinserted IAW reference below.
238+
239+ Reference: Hsu,C.S., and Hou,D., 1991,
240+ Reducing unstable linear control systems via real Schur transformation.
241+ Electronics Letters, 27, 984-986.
204242
205243 Parameters
206244 ----------
@@ -215,23 +253,39 @@ def balred(sys, orders, method='truncate'):
215253 Returns
216254 -------
217255 rsys: StateSpace
218- A reduced order model
256+ A reduced order model or a list of reduced order models if orders is a list
219257
220258 Raises
221259 ------
222260 ValueError
223- * if `method` is not ``'truncate'``
224- * if eigenvalues of `sys.A` are not all in left half plane
225- (`sys` must be stable)
261+ * if `method` is not ``'truncate'`` or ``'matchdc'``
226262 ImportError
227- if slycot routine ab09ad is not found
263+ if slycot routine ab09ad or ab09bd is not found
264+
265+ ValueError
266+ if there are more unstable modes than any value in orders
228267
229268 Examples
230269 --------
231- >>> rsys = balred(sys, order , method='truncate')
270+ >>> rsys = balred(sys, orders , method='truncate')
232271
233272 """
273+ if method != 'truncate' and method != 'matchdc' :
274+ raise ValueError ("supported methods are 'truncate' or 'matchdc'" )
275+ elif method == 'truncate' :
276+ try :
277+ from slycot import ab09ad
278+ except ImportError :
279+ raise ControlSlycot ("can't find slycot subroutine ab09ad" )
280+ elif method == 'matchdc' :
281+ try :
282+ from slycot import ab09bd
283+ except ImportError :
284+ raise ControlSlycot ("can't find slycot subroutine ab09bd" )
285+
234286
287+ from scipy .linalg import schur #, cholesky, svd
288+ from numpy .linalg import cholesky , svd
235289 #Check for ss system object, need a utility for this?
236290
237291 #TODO: Check for continous or discrete, only continuous supported right now
@@ -241,34 +295,81 @@ def balred(sys, orders, method='truncate'):
241295 # dico = 'D'
242296 # else:
243297 dico = 'C'
244-
245- #Check system is stable
246- D ,V = np .linalg .eig (sys .A )
247- # print D.shape
248- # print D
249- for e in D :
250- if e .real >= 0 :
251- raise ValueError ("Oops, the system is unstable!" )
252-
253- if method == 'matchdc' :
254- raise ValueError ("MatchDC not yet supported!" )
255- elif method == 'truncate' :
256- try :
257- from slycot import ab09ad
258- except ImportError :
259- raise ControlSlycot ("can't find slycot subroutine ab09ad" )
260- job = 'B' # balanced (B) or not (N)
261- equil = 'N' # scale (S) or not (N)
262- n = np .size (sys .A ,0 )
263- m = np .size (sys .B ,1 )
264- p = np .size (sys .C ,0 )
265- Nr , Ar , Br , Cr , hsv = ab09ad (dico ,job ,equil ,n ,m ,p ,sys .A ,sys .B ,sys .C ,nr = orders ,tol = 0.0 )
266-
267- rsys = StateSpace (Ar , Br , Cr , sys .D )
298+ job = 'B' # balanced (B) or not (N)
299+ equil = 'N' # scale (S) or not (N)
300+
301+ rsys = [] #empty list for reduced systems
302+
303+ #check if orders is a list or a scalar
304+ try :
305+ order = iter (orders )
306+ except TypeError : #if orders is a scalar
307+ orders = [orders ]
308+
309+ #first get original system order
310+ nn = sys .A .shape [0 ] #no. of states
311+ mm = sys .B .shape [1 ] #no. of inputs
312+ rr = sys .C .shape [0 ] #no. of outputs
313+ #first do the schur decomposition
314+ T , V , l = schur (sys .A , sort = 'lhp' ) #l will contain the number of eigenvalues in the open left half plane, i.e. no. of stable eigenvalues
315+
316+ for i in orders :
317+ rorder = i - (nn - l )
318+ if rorder <= 0 :
319+ raise ValueError ("System has %i unstable states which is more than ORDER(%i)" % (nn - l , i ))
320+
321+ for i in orders :
322+ if (nn - l ) > 0 : #handles the stable/unstable decomposition if unstable eigenvalues are found, nn - l is the number of ustable eigenvalues
323+ #Author: M. Clement (mdclemen@eng.ucsd.edu) 2016
324+ print ("Unstable eigenvalues found, performing stable/unstable decomposition" )
325+
326+ rorder = i - (nn - l )
327+ A_ , B_tilde , C_ , D_tilde , Ap , Bp , Cp = stabsep (T , V , sys , l , mm , rr )
328+
329+ subSys = StateSpace (A_ , B_tilde , C_ , D_tilde )
330+ n = np .size (subSys .A ,0 )
331+ m = np .size (subSys .B ,1 )
332+ p = np .size (subSys .C ,0 )
333+
334+ if method == 'truncate' :
335+ Nr , Ar , Br , Cr , hsv = ab09ad (dico ,job ,equil ,n ,m ,p ,subSys .A ,subSys .B ,subSys .C ,nr = rorder ,tol = 0.0 )
336+ rsubSys = StateSpace (Ar , Br , Cr , np .zeros ((p ,m )))
337+
338+ elif method == 'matchdc' :
339+ Nr , Ar , Br , Cr , Dr , hsv = ab09bd (dico ,job ,equil ,n ,m ,p ,subSys .A ,subSys .B ,subSys .C ,subSys .D ,nr = rorder ,tol1 = 0.0 ,tol2 = 0.0 )
340+ rsubSys = StateSpace (Ar , Br , Cr , Dr )
341+
342+ A_r = rsubSys .A
343+ #IAW ref 1 eq(4) B^{tilde}_r = [B_r, Acr]
344+ B_r = rsubSys .B [:,0 :mm ]
345+ Acr = rsubSys .B [:,mm :mm + (nn - l )]
346+ C_r = rsubSys .C
347+
348+ #now put the unstable subsystem back in
349+ Ar = np .bmat ([[A_r , Acr ], [np .zeros ((nn - l ,rorder )), Ap ]])
350+ Br = np .bmat ([[B_r ], [Bp ]])
351+ Cr = np .bmat ([[C_r , Cp ]])
352+
353+ rsys .append (StateSpace (Ar , Br , Cr , sys .D ))
354+
355+ else : #stable system branch
356+ n = np .size (sys .A ,0 )
357+ m = np .size (sys .B ,1 )
358+ p = np .size (sys .C ,0 )
359+ if method == 'truncate' :
360+ Nr , Ar , Br , Cr , hsv = ab09ad (dico ,job ,equil ,n ,m ,p ,sys .A ,sys .B ,sys .C ,nr = i ,tol = 0.0 )
361+ rsys .append (StateSpace (Ar , Br , Cr , sys .D ))
362+
363+ elif method == 'matchdc' :
364+ Nr , Ar , Br , Cr , Dr , hsv = ab09bd (dico ,job ,equil ,n ,m ,p ,sys .A ,sys .B ,sys .C ,sys .D ,nr = rorder ,tol1 = 0.0 ,tol2 = 0.0 )
365+ rsys .append (StateSpace (Ar , Br , Cr , Dr ))
366+
367+ #if orders was a scalar, just return the single reduced model, not a list
368+ if len (orders ) == 1 :
369+ return rsys [0 ]
370+ #if orders was a list/vector, return a list/vector of systems
268371 else :
269- raise ValueError ("Oops, method is not supported!" )
270-
271- return rsys
372+ return rsys
272373
273374def minreal (sys , tol = None , verbose = True ):
274375 '''
0 commit comments