Skip to content

Commit 66c0088

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

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
@@ -385,23 +385,100 @@ def horner(self, s):
385385
self.B) + self.D
386386
return array(resp)
387387

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

406483
# Compute poles and zeros
407484
def pole(self):

0 commit comments

Comments
 (0)