Skip to content

Commit fbcdc21

Browse files
committed
Use slycot's tb05ad for faster and more accurate frequency response evaluation of state-space systems. If slycot is unavailible, use the built in horner function (instead of converting to a transfer function, as was done before).
1 parent da0ddc1 commit fbcdc21

File tree

1 file changed

+88
-11
lines changed

1 file changed

+88
-11
lines changed

control/statesp.py

Lines changed: 88 additions & 11 deletions
Original file line numberDiff line numberDiff line change
@@ -378,23 +378,100 @@ def horner(self, s):
378378
self.B) + self.D
379379
return array(resp)
380380

381+
381382
# Method for generating the frequency response of the system
382-
def freqresp(self, omega):
383-
"""Evaluate the system's transfer func. at a list of ang. frequencies.
383+
def freqresp(self, omega_s):
384+
"""Evaluate the system's transfer func. at a list of freqs, omega.
384385
385386
mag, phase, omega = self.freqresp(omega)
386387
387-
reports the value of the magnitude, phase, and angular frequency of the
388-
system's transfer function matrix evaluated at s = i * omega, where
389-
omega is a list of angular frequencies, and is a sorted version of the
390-
input omega.
388+
Reports the frequency response of the system,
389+
390+
G(j*omega) = mag*exp(j*phase)
391+
392+
for continuous time. For discrete time systems, the response is
393+
evaluated around the unit circle such that
394+
395+
G(exp(j*omega*dt)) = mag*exp(j*phase).
396+
397+
Inputs:
398+
------
399+
omega_s: A list of frequencies in radians/sec at which the system
400+
should be evaluated. The list can be either a python list
401+
or a numpy array and will be sorted before evaluation.
402+
403+
Returns:
404+
-------
405+
mag: The magnitude (absolute value, not dB or log10) of the system
406+
frequency response.
407+
408+
phase: The wrapped phase in radians of the system frequency
409+
response.
410+
411+
omega_s: The list of sorted frequencies at which the response
412+
was evaluated.
391413
392414
"""
393-
# when evaluating at many frequencies, much faster to convert to
394-
# transfer function first and then evaluate, than to solve an
395-
# n-dimensional linear system at each frequency
396-
tf = _convertToTransferFunction(self)
397-
return tf.freqresp(omega)
415+
416+
# In case omega is passed in as a list, rather than a proper array.
417+
omega_s = np.asarray(omega_s)
418+
419+
numFreqs = len(omega_s)
420+
Gfrf = np.empty((self.outputs, self.inputs, numFreqs),
421+
dtype=np.complex128)
422+
423+
# Sort frequency and calculate complex frequencies on either imaginary
424+
# axis (continuous time) or unit circle (discrete time).
425+
omega_s.sort()
426+
if isdtime(self, strict=True):
427+
dt = timebase(self)
428+
cmplx_freqs = exp(1.j * omega_s * dt)
429+
if ((omega_s * dt).any() > pi):
430+
warn_message = ("evalfr: frequency evaluation"
431+
" above Nyquist frequency")
432+
warnings.warn(warn_message)
433+
else:
434+
cmplx_freqs = omega_s * 1.j
435+
436+
# Do the frequency response evaluation. Use TB05AD from Slycot
437+
# if it's available, otherwise use the built-in horners function.
438+
try:
439+
from slycot import tb05ad
440+
441+
n = np.shape(self.A)[0]
442+
m = self.inputs
443+
p = self.outputs
444+
# The first call both evalates C(sI-A)^-1 B and also returns
445+
# hessenberg transformed matrices at, bt, ct.
446+
result = tb05ad(n, m, p, cmplx_freqs[0], self.A,
447+
self.B, self.C, job='NG')
448+
# When job='NG', result = (at, bt, ct, g_i, hinvb, info)
449+
at = result[0]
450+
bt = result[1]
451+
ct = result[2]
452+
453+
# TB05AD freqency evaluation does not include direct feedthrough.
454+
Gfrf[:, :, 0] = result[3] + self.D
455+
456+
# Now, iterate through the remaining frequencies using the
457+
# transformed state matrices, at, bt, ct.
458+
459+
# Start at the second frequency, already have the first.
460+
for kk, cmplx_freqs_kk in enumerate(cmplx_freqs[1:numFreqs]):
461+
result = tb05ad(n, m, p, cmplx_freqs_kk, at,
462+
bt, ct, job='NH')
463+
# When job='NH', result = (g_i, hinvb, info)
464+
465+
# kk+1 because enumerate starts at kk = 0.
466+
# but zero-th spot is already filled.
467+
Gfrf[:, :, kk+1] = result[0] + self.D
468+
469+
except ImportError: # Slycot unavailable. Fall back to horner.
470+
for kk, cmplx_freqs_kk in enumerate(cmplx_freqs):
471+
Gfrf[:, :, kk] = self.horner(cmplx_freqs_kk)
472+
473+
# mag phase omega_s
474+
return np.abs(Gfrf), np.angle(Gfrf), omega_s
398475

399476
# Compute poles and zeros
400477
def pole(self):

0 commit comments

Comments
 (0)