Skip to content

Commit f630b4a

Browse files
roryyorkemurrayrm
authored andcommitted
Change docstrings to numpydoc conventions; update doc func index
1 parent 14079fb commit f630b4a

2 files changed

Lines changed: 42 additions & 36 deletions

File tree

control/canonical.py

Lines changed: 40 additions & 36 deletions
Original file line numberDiff line numberDiff line change
@@ -24,7 +24,7 @@ def canonical_form(xsys, form='reachable'):
2424
----------
2525
xsys : StateSpace object
2626
System to be transformed, with state 'x'
27-
form : String
27+
form : str
2828
Canonical form for transformation. Chosen from:
2929
* 'reachable' - reachable canonical form
3030
* 'observable' - observable canonical form
@@ -34,7 +34,7 @@ def canonical_form(xsys, form='reachable'):
3434
-------
3535
zsys : StateSpace object
3636
System in desired canonical form, with state 'z'
37-
T : matrix
37+
T : (M, M) real ndarray
3838
Coordinate transformation matrix, z = T * x
3939
"""
4040

@@ -63,7 +63,7 @@ def reachable_form(xsys):
6363
-------
6464
zsys : StateSpace object
6565
System in reachable canonical form, with state `z`
66-
T : matrix
66+
T : (M, M) real ndarray
6767
Coordinate transformation: z = T * x
6868
"""
6969
# Check to make sure we have a SISO system
@@ -117,7 +117,7 @@ def observable_form(xsys):
117117
-------
118118
zsys : StateSpace object
119119
System in observable canonical form, with state `z`
120-
T : matrix
120+
T : (M, M) real ndarray
121121
Coordinate transformation: z = T * x
122122
"""
123123
# Check to make sure we have a SISO system
@@ -164,11 +164,11 @@ def similarity_transform(xsys, T, timescale=1, inverse=False):
164164
----------
165165
xsys : StateSpace object
166166
System to transform
167-
T : 2D invertible array
167+
T : (M, M) array_like
168168
The matrix `T` defines the new set of coordinates z = T x.
169-
timescale : float
169+
timescale : float, optional
170170
If present, also rescale the time unit to tau = timescale * t
171-
inverse: boolean
171+
inverse: boolean, optional
172172
If True (default), transform so z = T x. If False, transform
173173
so x = T z.
174174
@@ -181,6 +181,8 @@ def similarity_transform(xsys, T, timescale=1, inverse=False):
181181
# Create a new system, starting with a copy of the old one
182182
zsys = StateSpace(xsys)
183183

184+
T = np.atleast_2d(T)
185+
184186
# Define a function to compute the right inverse (solve x M = y)
185187
def rsolve(M, y):
186188
return transpose(solve(transpose(M), transpose(y)))
@@ -207,14 +209,16 @@ def _bdschur_defective(blksizes, eigvals):
207209
208210
Parameters
209211
----------
210-
blksizes: size of Schur blocks
211-
eigvals: eigenvalues
212+
blksizes: (N,) int ndarray
213+
size of Schur blocks
214+
eigvals: (M,) real or complex ndarray
215+
Eigenvalues
212216
213217
Returns
214218
-------
215-
True iff Schur blocks are defective
219+
True iff Schur blocks are defective.
216220
217-
blksizes, eigvals are 3rd and 4th results returned by mb03rd.
221+
blksizes, eigvals are the 3rd and 4th results returned by mb03rd.
218222
"""
219223
if any(blksizes > 2):
220224
return True
@@ -245,22 +249,22 @@ def _bdschur_condmax_search(aschur, tschur, condmax):
245249
246250
Parameters
247251
----------
248-
aschur: (n, n) array
249-
real Schur-form matrix
250-
tschur: (n, n) array
251-
orthogonal transformation giving aschur from some initial matrix a
252-
condmax: positive scalar >= 1
253-
maximum condition number of final transformation
252+
aschur: (N, N) real ndarray
253+
Real Schur-form matrix
254+
tschur: (N, N) real ndarray
255+
Orthogonal transformation giving aschur from some initial matrix a
256+
condmax: float
257+
Maximum condition number of final transformation. Must be >= 1.
254258
255259
Returns
256260
-------
257-
amodal: n, n array
261+
amodal: (N, N) real ndarray
258262
block diagonal Schur form
259-
tmodal:
263+
tmodal: (N, N) real ndarray
260264
similarity transformation give amodal from aschur
261-
blksizes:
265+
blksizes: (M,) int ndarray
262266
Array of Schur block sizes
263-
eigvals:
267+
eigvals: (N,) real or complex ndarray
264268
Eigenvalues of amodal (and a, etc.)
265269
266270
Notes
@@ -338,20 +342,20 @@ def bdschur(a, condmax=None, sort=None):
338342
339343
Parameters
340344
----------
341-
a: real (n, n) array
342-
Matrix to decompose
343-
condmax: real scalar >= 1
344-
If None (default), use `1/sqrt(eps)`, which is approximately 2e-8
345-
sort: None, 'continuous', or 'discrete'
346-
See below
345+
a : (M, M) array_like
346+
Real matrix to decompose
347+
condmax : None or float, optional
348+
If None (default), use 1/sqrt(eps), which is approximately 1e8
349+
sort : {None, 'continuous', 'discrete'}
350+
Block sorting; see below.
347351
348352
Returns
349353
-------
350-
amodal: (n, n) array, dtype `np.double`
354+
amodal : (M, M) real ndarray
351355
Block-diagonal Schur decomposition of `a`
352-
tmodal: (n, n) array
353-
similarity transform relating `a` and `amodal`
354-
blksizes:
356+
tmodal : (M, M) real ndarray
357+
Similarity transform relating `a` and `amodal`
358+
blksizes : (N,) int ndarray
355359
Array of Schur block sizes
356360
357361
Notes
@@ -365,7 +369,7 @@ def bdschur(a, condmax=None, sort=None):
365369
366370
If `sort` is 'discrete', the blocks are sorted as for
367371
'continuous', but applied to log of eigenvalues
368-
(continuous-equivalent).
372+
(i.e., continuous-equivalent eigenvalues).
369373
"""
370374
if condmax is None:
371375
condmax = np.finfo(np.float64).eps ** -0.5
@@ -421,16 +425,16 @@ def modal_form(xsys, condmax=None, sort=False):
421425
----------
422426
xsys : StateSpace object
423427
System to be transformed, with state `x`
424-
condmax: real scalar >= 1
428+
condmax : None or float, optional
425429
An upper bound on individual transformations. If None, use `bdschur` default.
426-
sort: False (default)
427-
If true, Schur blocks will be sorted. See `bdschur` for sort order.
430+
sort : bool, optional
431+
If False (default), Schur blocks will not be sorted. See `bdschur` for sort order.
428432
429433
Returns
430434
-------
431435
zsys : StateSpace object
432436
System in modal canonical form, with state `z`
433-
T : matrix
437+
T : (M, M) ndarray
434438
Coordinate transformation: z = T * x
435439
"""
436440

doc/control.rst

Lines changed: 2 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -155,6 +155,7 @@ Utility functions and conversions
155155
:toctree: generated/
156156

157157
augw
158+
bdschur
158159
canonical_form
159160
damp
160161
db2mag
@@ -163,6 +164,7 @@ Utility functions and conversions
163164
issiso
164165
issys
165166
mag2db
167+
modal_form
166168
observable_form
167169
pade
168170
reachable_form

0 commit comments

Comments
 (0)