kontrol.core
kontrol.core.controlutils module
Utility function related to python-control library.
- kontrol.core.controlutils.check_tf_equal(tf1, tf2, allclose_kwargs={}, minreal=True)
Check if two transfer functions are approximatedly equal by np.allclose.
- Parameters
tf1 (control.xferfcn.TransferFunction) – Transfer function 1.
tf2 (control.xferfcn.TrasnferFunction) – Transfer function 2.
allclose_kwargs (dict, optional) – Keyword arguments passed to np.allclose(), which is used to compare the list of poles and zeros and dcgain. Defaults {}.
minreal (boolean) – Use control.minreal to remove canceling zeros and poles before comparison.
- Return type
boolean
- kontrol.core.controlutils.clean_tf(tf, tol_order=5, small_number=1e-25)
Remove numerator/denominator coefficients that are small outliers
- Parameters
tf (TransferFunction) – The transfer function to be cleaned
tol_order (float, optional) – If the coefficient is
tol_order
order smaller than the rest of the coefficients, then this coefficient is an outlier. Defaults 5.small_number (float, optional) – A small number to be added to the log10 in case 0 is encountered.
- Returns
tf_cleaned – The cleaned transfer function.
- Return type
- kontrol.core.controlutils.clean_tf2(tf, tol_order=5, small_number=1e-25)
Remove zeros/poles that are outliers.
- Parameters
tf (TransferFunction) – The transfer function to be cleaned.
tol_order (float, optional) – If the frequency of the zero/pole is ``tol_order’’ order away from the mean order, then this zero/pole is an outlier. Defaults 5.
small_number (float, optional) – A small number to be added to the log10() in case 0 is encountered
- Returns
tf_cleaned – The cleaned transfer function.
- Return type
- kontrol.core.controlutils.clean_tf3(tf, tol_order=5, small_number=1e-25)
Remove first coefficient if it is much smaller than the second
- Parameters
tf (TransferFunction) – Transfer function to be cleaned
tol_order (float, optional) – First coefficient is removed if it is ``tol_order’’ order smaller than the second coefficient. Remove only if there exists such coefficient in both numerator and denominator.
small_number (float, optional) – A small number to be added to the log10() in case 0 is encountered.
- Returns
tf_cleaned – The cleaned TransferFunction
- Return type
- kontrol.core.controlutils.complex2wq(complex_frequency)
Convert a complex frequency to frequency and Q-factor.
- Parameters
complex_frequency (complex) – Complex frequency in rad/s.
- Returns
wn (float) – The frequency in rad/s
q (float) – The Q factor.
- kontrol.core.controlutils.convert_unstable_tf(control_tf)
Convert transfer function with unstable zeros and poles to tf without.
- Parameters
control_tf (control.xferfcn.TransferFunction) – The transfer function to be converted.
- Returns
tf_new – control_tf with unstable zeros and poles flipped.
- Return type
control.xferfcn.TransferFunction
Note
Warning can occur when the order gets high, e.g. 100. This happens due to numerical accuracy, i.e. coefficients get astronomically high when order gets high.
- kontrol.core.controlutils.generic_tf(zeros=None, poles=None, zeros_wn=None, zeros_q=None, poles_wn=None, poles_q=None, dcgain=1.0, unit='f')
Construct a generic transfer function object.
- Parameters
zeros (array, optional) – List of zeros in frequencies. Defaults [].
poles (array, optional) – List of poles. Defaults [].
zeros_wn (array, optional) – List of natural frequencies of numerator second-order sections.
zeros_q (array, optional.) – List of Q-value of numerator second-order sections.
poles_wn (array, optional) – List of natural frequencies of denominator second-order sections.
poles_q (array, optional.) – List of Q-value of denominator second-order sections.
dcgain (float, optional) – The DC gain of the transfer function. Defaults 1.
unit (str, optional) – The unit of the zeros, poles and the natural frequencies. Choose from [“f”, “s”, “Hz”, “omega”]. Defaults “f”.
- Returns
tf – The transfer function.
- Return type
control.xferfcn.TransferFunction
Note
No negative sign is needed for negative zeros and poles. For instance, generic_tf(poles=[1], unit=”s”) means \(\frac{1}{s+1}\).
- kontrol.core.controlutils.outlier_exists(tf, f, unit='f')
Checks for zeros and poles outside the frequency range.
- Parameters
tf (control.xferfcn.TransferFunction) – The transfer function.
f (array) – The frequency axis of interest.
unit (str, optional) – The unit of the zeros, poles and the natural frequencies. Choose from [“f”, “s”, “Hz”, “omega”]. Defaults “f”.
- Returns
If there’s any zero or pole outside the frequency range.
- Return type
boolean
- kontrol.core.controlutils.outliers(tf, f, unit='f')
Returns a list of zeros and poles outside the frequency range.
- Parameters
tf (control.xferfcn.TransferFunction) – The transfer function.
f (array) – The frequency axis of interest.
unit (str, optional) – The unit of the zeros, poles and the natural frequencies. Choose from [“f”, “s”, “Hz”, “omega”]. Defaults “f”.
- Returns
outlier_zeros (array) – Zeros outside the frequency range.
outlier_poles (array) – Poles outside the frequency range.
Note
We use the python-control package convention for the returned zeros and poles, so to preserve the type and format as much as possible for further processes.
- kontrol.core.controlutils.sos(natural_frequencies=None, quality_factors=None, gain=1.0, unit='f')
Second-order sections transfer fucntion.
- Parameters
natural_freuqnecies (array, optional) – List of natural frequencies. Defaults [].
quality_factors (list, optional) – List of quality factors. Defaults [].
gain (float, optional) – The gain of the transfer function. Defaults 1.
unit (str, optional) – The unit of the zeros, poles and the natural frequencies. Choose from [“f”, “s”, “Hz”, “omega”]. Defaults “f”.
- Returns
sos – The product of all second-order sections.
- Return type
control.xferfcn.TransferFunction
Note
\(K\prod\left(s^2+\omega_n/Q\,s+\omega_n^2\right)\).
- kontrol.core.controlutils.tf_order_split(tf, max_order=20)
Split TransferFunction objects into multiple ones with fewer order.
- Parameters
tf (TransferFunction) – The transfer function
max_order (int, optional) – The maxmium order of the split transfer functions. Defaults to 20 (the foton limit).
- Returns
The list of splitted transfer functions.
- Return type
list of TransferFunction.
- kontrol.core.controlutils.tfmatrix2tf(sys)
Convert a matrix of transfer functions to a MIMO transfer function.
- Parameters
sys (list of (list of control.xferfcn.TransferFunction)) – The transfer function matrix representing a MIMO system. sys[i][j] is the transfer function from the i+1 input port to the j+1 output port.
- Returns
The transfer function of the MIMO system.
- Return type
control.xferfcn.TransferFunction
- kontrol.core.controlutils.zpk(zeros, poles, gain, unit='f', negate=True)
Zero-pole-gain definition of transfer function.
- Parameters
zeros (list of floats) – A list of the location of the zeros
poles (list of floats) – A list of the location of the poles
gain (float) – The static gain of the transfer function
unit (str, optional) – The unit of the zeros, poles and the natural frequencies. Choose from [“f”, “s”, “Hz”, “omega”]. Defaults “f”.
negate (boolean, optional) – Negate zeros and poles in specification so negative sign is not needed for stable zeros and poles. Default to be True.
- Returns
zpk_tf – The zpk defined transfer function
- Return type
control.xferfcn.TransferFunction
Notes
Refrain from specifying imaginary zeros and poles. Use kontrol.utils.sos() for second-order sections instead. The zero and poles are negated by default.
kontrol.core.math module
Common math library
- kontrol.core.math.log_mse(x1, x2, weight=None, small_multiplier=1e-06)
Logarithmic mean square error/Mean square log error
- Parameters
x1 (array) – Array
x2 (array) – Array
weight (array or None, optional) – weighting function. If None, defaults to np.ones_like(x1)
small_multiplier (float, optional) – x1 will be modified by x1`+`small_multiplier`*`np.min(x2) if 0 exists in x1. Same goes to x2. If 0 both exists in x1 and x2, raise.
- Returns
Logarithmic mean square error between arrays x1 and x2.
- Return type
float
- kontrol.core.math.mse(x1, x2, weight=None)
Mean square error
- Parameters
x1 (array) – Array
x2 (array) – Array
weight (array or None, optional) – weighting function. If None, defaults to np.ones_like(x1)
- Returns
Mean square error between arrays x1 and x2.
- Return type
float
- kontrol.core.math.polyval(p, x)
- kontrol.core.math.quad_sum(*spectra)
Takes any number of same length spectrum and returns the quadrature sum.
- Parameters
*spectra – Variable length argument list of the spectra.
- Returns
qs – The quadrature sum of the spectra.
- Return type
np.ndarray
kontrol.core.spectral module
Spectral analysis related function library.
- kontrol.core.spectral.asd2rms(asd, f=None, df=1.0, return_series=True)
Calculate root-mean-squared value from amplitude spectral density
- Parameters
asd (array) – The amplitude spectral density
f (array, optional) – The frequency axis. Defaults
None
.df (float, optional) – The frequency spacing. Defaults 1.
return_series (bool, optional) – Returns the RMS as a frequency series, where each value is the integrated RMS from highest frequency. If False, returns a single RMS value for the whole spectrum. Defaults True.
- Returns
array – The integrated RMS series, if
return_series==True
.float – The integrated RMS value, if
return_series==False
.
Notes
The integrated RMS series is defined as
\[x_\mathrm{RMS}(f) = \int_\infty^{f}\,x(f')^2\,df'\,,\]where \(x(f)\) is the amplitude spectral density.
When
return_series
is False, only \(x_\mathrm{RMS}(0)\) is returned.
- kontrol.core.spectral.asd2ts(asd, f=None, fs=None, t=None, window=None, zero_mean=True)
Simulate time series from amplitude spectral density
- Parameters
asd (array) – The amplitude spectral density.
f (array, optional) – The frequency axis of the amplitude spectral density. This is used to calculate the sampling frequency. If
fs
is not specified,f
must be specified. DefaultsNone
.fs (float, optional) – The sampling frequency in Hz. If
f
is specified, the last element off
will be treated as the Nyquist frequency so the double of it would be the sampling frequency. DefaultsNone
.t (array, optional) – The desired time axis for the time series. If
t
is specified, then the time series will be interpolated usingnumpy.interp()
assuming that the time series is periodic. DefaultNone
.window (array, optional) – The FFT window used to compute the amplitude spectral density. Defaults
None
.zero_mean (boolean, optional) – Returns a time series with zero mean. Defaults
True
.
- Returns
t (array) – Time axis.
time_series (array) – The time series
Notes
Using a custom time axis is not recommended as interpolation leads to distortion of the signal. To extend the signal in time, it’s recommended to repeat the original time series instead.
- kontrol.core.spectral.pad_above_maxima(series, pad_index=-1, **kwargs)
Pad series above local maxima
- Parameters
series (array) – The series to be padded.
pad_index (int, optional) – The index of the local maxima. Defaults to -1.
**kwargs – Keyword arguments passed to scipy.signal.argrelmin()
- Returns
series_padded – The padded series
- Return type
array
- kontrol.core.spectral.pad_above_minima(series, pad_index=-1, **kwargs)
Pad series above local minima
- Parameters
series (array) – The series to be padded.
pad_index (int, optional) – The index of the local minima. Defaults to -1.
**kwargs – Keyword arguments passed to scipy.signal.argrelmin()
- Returns
series_padded – The padded series
- Return type
array
- kontrol.core.spectral.pad_below_maxima(series, pad_index=0, **kwargs)
Pad series below local maxima
- Parameters
series (array) – The series to be padded.
pad_index (int, optional) – The index of the local maxima. Defaults to 0.
**kwargs – Keyword arguments passed to scipy.signal.argrelmin()
- Returns
series_padded – The padded series
- Return type
array
- kontrol.core.spectral.pad_below_minima(series, pad_index=0, **kwargs)
Pad series below local minima
- Parameters
series (array) – The series to be padded.
pad_index (int, optional) – The index of the local minima. Defaults to 0.
**kwargs – Keyword arguments passed to scipy.signal.argrelmin()
- Returns
series_padded – The padded series
- Return type
array
- kontrol.core.spectral.three_channel_correlation(psd1, psd2=None, psd3=None, csd12=None, csd13=None, csd21=None, csd23=None, csd31=None, csd32=None, returnall=True)
Noise estimation from three sensors’ readouts.
- Parameters
psd1 (array) – The power spectral density of the readout of the first sensor.
psd2 (array, optional) – The power spectral density of the readout of the second sensor. Defaults
None
.psd3 (array, optional) – The power spectral density of the readout of the third sensor. Defaults
None
.csd12 (array, optional) – Cross power spectral density between readout 1 and 2. If not specified, this will be estimated as
psd1*psd2/csd21
. DefaultNone
.csd13 (array, optional) – Cross power spectral density between readout 1 and 3. If not specified, this will be estimated as
psd1*psd3/csd31
. DefaultNone
.csd21 (array, optional) – Cross power spectral density between readout 2 and 1. If not specified, this will be estimated as
psd2*psd1/csd12
. DefaultNone
.csd23 (array, optional) – Cross power spectral density between readout 2 and 3. If not specified, this will be estimated as
psd2*psd3/csd32
. DefaultNone
.csd31 (array, optional) – Cross power spectral density between readout 3 and 1. If not specified, this will be estimated as
psd1*psd3/csd13
. DefaultNone
.csd32 (array, optional) – Cross power spectral density between readout 3 and 1. If not specified, this will be estimated as
psd3*psd2/csd23
. DefaultNone
.returnall (boolean, optional) – If
True
, return all three noise estimations. IfFalse
, return noise estimation of first sensor only. Defaults True.
- Returns
noise1 (array) – Power spectral density of the estimated noise in psd1.
noise2 (array, optional) – Power spectral density of the estimated noise in psd2. Returns only if
returnall==True
noise3 (array, optional) – Power spectral density of the estimated noise in psd3. Returns only if
returnall==True
Notes
The PSD of the noise in psd1 is then computed as
\[P_{n_1n_1}(f) = \left\lvert P_{x_1x_1}(f) - \frac{P_{x_1x_3}(f)}{P_{x_2x_3}(f)}P_{x_2x_1}\right\vert\,,\]If
returnall
isTrue
, at leastpsd1
,psd2
,psd3
, (csd13
orcsd31
), (csd23
orcsd32
), and (csd12
andcsd21
) must be provided.References
- 2
R. Sleeman, A. Wettum, and J. Trampert. Three-channel correlation analysis: A new technique to measure instrumental noise of digitizers and seismic sensors. Bulletin of the Seismological Society of America, 96:258–271, 2006.
- kontrol.core.spectral.two_channel_correlation(psd, coh)
Noise estimation from two identical sensors’ readouts.
- Parameters
psd (array) – The power spectral density of the readout of the sensor
coh (array) – The coherence between readout of the sensor and another identical sensor.
- Returns
noise – Power spectral density of the estimated noise.
- Return type
array
Notes
The PSD of the noise is computed as
\[P_{nn}(f) = P_{x_1x_1}(f)\left(1-C_{x_1x_2}(f)^{\frac{1}{2}}\right)\,,\]where \(P_{x_1x_1}(f)\) is the power spectral density of the readout and \(C_{x_1x_2}(f)\) is the coherence between the two sensor readouts.
References
- 1
Aaron Barzilai, Tom VanZandt, and Tom Kenny. Technique for measurement of the noise of a sensor in the presence of large background signals. Review of Scientific Instruments, 69:2767–2772, 07 1998.
kontrol.core.foton module
KAGRA/LIGO Foton related utilities.
- kontrol.core.foton.foton2tf(foton_string)
Convert a Foton string to TransferFunction
- Parameters
foton_string (str) – The Foton string, e.g. zpk([0], [1; 1], 1, “n”).
- Returns
tf – The transfer function.
- Return type
- kontrol.core.foton.get_zpk2tf(get_zpk, plane='s')
Converts Foton’s get_zpk() output to TransferFunction
- Parameters
get_zpk (tuple(array, array, float)) – The output from Foton’s get_zpk() method of a filter instance.
plane (str, optional) – “s”, “f”, or “n”. Default “s”.
- Returns
tf – The converted transfer function.
- Return type
- kontrol.core.foton.notch(frequency, q, depth, significant_figures=6)
Returns the foton expression of a notch filter.
- Parameters
frequency (float) – The notch frequency (Hz).
q (float) – The quality factor.
depth (float) – The depth of the notch filter (magnitude).
significant_figures (int, optional) – Number of significant figures to print out. Defaults to 6.
- Returns
The foton representation of this notch filter.
- Return type
str
- kontrol.core.foton.rpoly2tf(foton_string)
Converts rpoly Foton strings to TransferFunction
- Parameters
foton_string (str) – The rpoly Foton string. E.g. rpoly([1; 2; 3], [2; 3; 4], 5)
- kontrol.core.foton.tf2foton(tf, expression='zpk', root_location='s', significant_figures=6, itol=1e-25, epsilon=1e-25)
Convert a single transfer function to foton expression.
- Parameters
tf (TransferFunction) – The transfer function object.
expression (str, optional) – Format of the foton expression. Choose from [“zpk”, “rpoly”]. Defaults to “zpk”.
root_location (str, optional) – Root location of the zeros and poles for expression==”zpk”. Choose from [“s”, “f”, “n”]. “s”: roots in s-plane, i.e. zpk([…], […], …, “s”). “f”: roots in frequency plane, i.e. zpk([…], [,,,], …, “f”). “n”: roots in frequency plane but negated and gains are normalized, i.e. real parts are positive zpk([…], […], …, “n”). Defaults to “s”.
significant_figures (int, optional) – Number of significant figures to print out. Defaults to 6.
itol (float, optional) – Treating complex roots as real roots if the ratio of the imaginary part and the real part is smaller than this tolerance Defaults to 1e-25.
epsilon (float, optional) – Small number to add to denominator to prevent division error. Defaults to 1e-25.
- Returns
foton_expression – The foton expression in selected format.
- Return type
str
- kontrol.core.foton.tf2rpoly(tf, significant_figures=6)
Convert a transfer function to foton rpoly expression.
- Parameters
tf (TransferFunction) – The transfer function object
significant_figures (int, optional) – Number of significant figures to print out. Defaults to 6.
- Returns
Foton express in foton rpoly expression.
- Return type
str
Notes
Only works for transfer functions with less than 20 orders.
- kontrol.core.foton.tf2zpk(tf, root_location='s', significant_figures=6, itol=1e-25, epsilon=1e-25)
Convert a single transfer function to foton zpk expression.
- Parameters
tf (TransferFunction) – The transfer function object.
root_location (str, optional) – Root location of the zeros and poles. Choose from [“s”, “f”, “n”]. “s”: roots in s-plane, i.e. zpk([…], […], …, “s”). “f”: roots in frequency plane, i.e. zpk([…], [,,,], …, “f”). “n”: roots in frequency plane but negated and gains are normalized, i.e. real parts are positive zpk([…], […], …, “n”). Defaults to “s”.
significant_figures (int, optional) – Number of significant figures to print out. Defaults to 6.
itol (float, optional) – Treating complex roots as real roots if the ratio of the imaginary part and the real part is smaller than this tolerance Defaults to 1e-25.
epsilon (float, optional) – Small number to add to denominator to prevent division error. Defaults to 1e-25.
- Returns
The foton zpk expression in selected format.
- Return type
str
Notes
Only works for transfer functions with less than 20 orders.