Miscellaneous

The Miscellaneous package provides various functions and methods which are used in the examples and in some of the other implemented routines.

The package contains the following modules:

Tools for 2nd order systems

The PyDynamic.misc.SecondOrderSystem module is a collection of methods that are used throughout the whole package, specialized for second order dynamic systems, such as the ones used for high-class accelerometers.

This module contains the following functions:

  • sos_FreqResp(): Calculation of the system frequency response
  • sos_phys2filter(): Calculation of continuous filter coefficients from physical parameters
  • sos_absphase(): Propagation of uncertainty from physical parameters to real and imaginary part of system’s transfer function using GUM S2 Monte Carlo
  • sos_realimag(): Propagation of uncertainty from physical parameters to real and imaginary part of system’s transfer function using GUM S2 Monte Carlo
PyDynamic.misc.SecondOrderSystem.sos_FreqResp(S, d, f0, freqs)[source]

Calculation of the system frequency response

The frequency response is calculated from the continuous physical model of a second order system given by

\(H(f) = \frac{4S\pi^2f_0^2}{(2\pi f_0)^2 + 2jd(2\pi f_0)f - f^2}\)

If the provided system parameters are vectors then \(H(f)\) is calculated for each set of parameters. This is helpful for Monte Carlo simulations by using draws from the model parameters

Parameters:
  • S (float or ndarray shape (K,)) – static gain
  • d (float or ndarray shape (K,)) – damping parameter
  • f0 (float or ndarray shape (K,)) – resonance frequency
  • freqs (ndarray shape (N,)) – frequencies at which to calculate the freq response
Returns:

H – complex frequency response values

Return type:

ndarray shape (N,) or ndarray shape (N,K)

PyDynamic.misc.SecondOrderSystem.sos_phys2filter(S, d, f0)[source]

Calculation of continuous filter coefficients from physical parameters.

If the provided system parameters are vectors then the filter coefficients are calculated for each set of parameters. This is helpful for Monte Carlo simulations by using draws from the model parameters

Parameters:
  • S (float) – static gain
  • d (float) – damping parameter
  • f0 (float) – resonance frequency
Returns:

b,a – analogue filter coefficients

Return type:

ndarray

PyDynamic.misc.SecondOrderSystem.sos_absphase(S, d, f0, uS, ud, uf0, f, runs=10000)[source]

Propagation of uncertainty from physical parameters to real and imaginary part of system’s transfer function using GUM S2 Monte Carlo.

Parameters:
  • S (float) – static gain
  • d (float) – damping
  • f0 (float) – resonance frequency
  • uS (float) – uncertainty associated with static gain
  • ud (float) – uncertainty associated with damping
  • uf0 (float) – uncertainty associated with resonance frequency
  • f (ndarray, shape (N,)) – frequency values at which to calculate amplitue and phase
Returns:

  • Hmean (ndarray, shape (N,)) – best estimate of complex frequency response values
  • Hcov (ndarray, shape (2N,2N)) – covariance matrix [ [u(abs,abs), u(abs,phase)], [u(phase,abs), u(phase,phase)] ]

PyDynamic.misc.SecondOrderSystem.sos_realimag(S, d, f0, uS, ud, uf0, f, runs=10000)[source]

Propagation of uncertainty from physical parameters to real and imaginary part of system’s transfer function using GUM S2 Monte Carlo.

Parameters:
  • S (float) – static gain
  • d (float) – damping
  • f0 (float) – resonance frequency
  • uS (float) – uncertainty associated with static gain
  • ud (float) – uncertainty associated with damping
  • uf0 (float) – uncertainty associated with resonance frequency
  • f (ndarray, shape (N,)) – frequency values at which to calculate real and imaginary part
Returns:

  • Hmean (ndarray, shape (N,)) – best estimate of complex frequency response values
  • Hcov (ndarray, shape (2N,2N)) – covariance matrix [ [u(real,real), u(real,imag)], [u(imag,real), u(imag,imag)] ]

Tools for digital filters

The PyDynamic.misc.filterstuff module is a collection of methods which are related to filter design.

This module contains the following functions:

  • db(): Calculation of decibel values \(20\log_{10}(x)\) for a vector of values
  • ua(): Shortcut for calculation of unwrapped angle of complex values
  • grpdelay(): Calculation of the group delay of a digital filter
  • mapinside(): Maps the roots of polynomial with coefficients \(a\) to the unit circle
  • kaiser_lowpass(): Design of a FIR lowpass filter using the window technique with a Kaiser window.
  • isstable(): Determine whether a given IIR filter is stable
  • savitzky_golay(): Smooth (and optionally differentiate) data with a Savitzky-Golay filter
PyDynamic.misc.filterstuff.db(vals)[source]

Calculation of decibel values \(20\log_{10}(x)\) for a vector of values

PyDynamic.misc.filterstuff.ua(vals)[source]

Shortcut for calculation of unwrapped angle of complex values

PyDynamic.misc.filterstuff.grpdelay(b, a, Fs, nfft=512)[source]

Calculation of the group delay of a digital filter

Parameters:
  • b (ndarray) – IIR filter numerator coefficients
  • a (ndarray) – IIR filter denominator coefficients
  • Fs (float) – sampling frequency of the filter
  • nfft (int) – number of FFT bins
Returns:

  • group_delay (np.ndarray) – group delay values
  • frequencies (ndarray) – frequencies at which the group delay is calculated

References

PyDynamic.misc.filterstuff.mapinside(a)[source]

Maps the roots of polynomial to the unit circle.

Maps the roots of polynomial with coefficients \(a\) to the unit circle.

Parameters:a (ndarray) – polynomial coefficients
Returns:a – polynomial coefficients with all roots inside or on the unit circle
Return type:ndarray
PyDynamic.misc.filterstuff.kaiser_lowpass(L, fcut, Fs, beta=8.0)[source]

Design of a FIR lowpass filter using the window technique with a Kaiser window.

This method uses a Kaiser window. Filters of that type are often used as FIR low-pass filters due to their linear phase.

Parameters:
  • L (int) – filter order (window length)
  • fcut (float) – desired cut-off frequency
  • Fs (float) – sampling frequency
  • beta (float) – scaling parameter for the Kaiser window
Returns:

  • blow (ndarray) – FIR filter coefficients
  • shift (int) – delay of the filter (in samples)

PyDynamic.misc.filterstuff.isstable(b, a, ftype='digital')[source]

Determine whether IIR filter (b,a) is stable

Determine whether IIR filter (b,a) is stable by checking roots of the polynomial ´a´.

Parameters:
  • b (ndarray) – filter numerator coefficients
  • a (ndarray) – filter denominator coefficients
  • ftype (string) – type of filter (digital or analog)
Returns:

stable – whether filter is stable or not

Return type:

bool

PyDynamic.misc.filterstuff.savitzky_golay(y, window_size, order, deriv=0, delta=1.0)[source]

Smooth (and optionally differentiate) data with a Savitzky-Golay filter

The Savitzky-Golay filter removes high frequency noise from data. It has the advantage of preserving the original shape and features of the signal better than other types of filtering approaches, such as moving averages techniques.

Source obtained from scipy cookbook (online), downloaded 2013-09-13

Parameters:
  • y (ndarray, shape (N,)) – the values of the time history of the signal
  • window_size (int) – the length of the window. Must be an odd integer number
  • order (int) – the order of the polynomial used in the filtering. Must be less then window_size - 1.
  • deriv (int, optional) – The order of the derivative to compute. This must be a nonnegative integer. The default is 0, which means to filter the data without differentiating.
  • delta (float, optional) – The spacing of the samples to which the filter will be applied. This is only used if deriv > 0. This includes a factor \(n! / h^n\), where \(n\) is represented by deriv and \(1/h\) by delta.
Returns:

ys – the smoothed signal (or it’s n-th derivative).

Return type:

ndarray, shape (N,)

Notes

The Savitzky-Golay is a type of low-pass filter, particularly suited for smoothing noisy data. The main idea behind this approach is to make for each point a least-square fit with a polynomial of high order over a odd-sized window centered at the point.

References

Test signals

The PyDynamic.misc.testsignals module is a collection of test signals which can be used to simulate dynamic measurements and test methods.

This module contains the following functions:

  • shocklikeGaussian(): signal that resembles a shock excitation as a Gaussian followed by a smaller Gaussian of opposite sign
  • GaussianPulse(): Generates a Gaussian pulse at \(t_0\) with height \(m_0\) and std \(sigma\)
  • rect(): Rectangular signal of given height and width \(t_1 - t_0\)
  • squarepulse(): Generates a series of rect functions to represent a square pulse signal
  • sine(): Generate a sine signal
PyDynamic.misc.testsignals.shocklikeGaussian(time, t0, m0, sigma, noise=0.0)[source]

Generates a signal that resembles a shock excitation as a Gaussian followed by a smaller Gaussian of opposite sign.

Parameters:
  • time (np.ndarray of shape (N,)) – time instants (equidistant)
  • t0 (float) – time instant of signal maximum
  • m0 (float) – signal maximum
  • sigma (float) – std of main pulse
  • noise (float, optional) – std of simulated signal noise
Returns:

x – signal amplitudes at time instants

Return type:

np.ndarray of shape (N,)

PyDynamic.misc.testsignals.GaussianPulse(time, t0, m0, sigma, noise=0.0)[source]

Generates a Gaussian pulse at t0 with height m0 and std sigma

Parameters:
  • time (np.ndarray of shape (N,)) – time instants (equidistant)
  • t0 (float) – time instant of signal maximum
  • m0 (float) – signal maximum
  • sigma (float) – std of pulse
  • noise (float, optional) – std of simulated signal noise
Returns:

x – signal amplitudes at time instants

Return type:

np.ndarray of shape (N,)

PyDynamic.misc.testsignals.rect(time, t0, t1, height=1, noise=0.0)[source]

Rectangular signal of given height and width t1-t0

Parameters:
  • time (np.ndarray of shape (N,)) – time instants (equidistant)
  • t0 (float) – time instant of rect lhs
  • t1 (float) – time instant of rect rhs
  • height (float) – signal maximum
  • noise (float or numpy.ndarray of shape (N,), optional) – float: standard deviation of additive white gaussian noise ndarray: user-defined additive noise
Returns:

x – signal amplitudes at time instants

Return type:

np.ndarray of shape (N,)

PyDynamic.misc.testsignals.squarepulse(time, height, numpulse=4, noise=0.0)[source]

Generates a series of rect functions to represent a square pulse signal

Parameters:
  • time (np.ndarray of shape (N,)) – time instants
  • height (float) – height of the rectangular pulses
  • numpulse (int) – number of pulses
  • noise (float, optional) – std of simulated signal noise
Returns:

x – signal amplitude at time instants

Return type:

np.ndarray of shape (N,)

class PyDynamic.misc.testsignals.corr_noise(w, sigma, seed=None)[source]

Base class for generation of a correlated noise process.

PyDynamic.misc.testsignals.sine(time, amp=1.0, freq=6.283185307179586, noise=0.0)[source]

Generate a sine signal

Parameters:
  • time (np.ndarray of shape (N,)) – time instants
  • amp (float, optional) – amplitude of the sine (default = 1.0)
  • freq (float, optional) – frequency of the sine in Hz (default = \(2 * \pi\))
  • noise (float, optional) – std of simulated signal noise (default = 0.0)
Returns:

x – signal amplitude at time instants

Return type:

np.ndarray of shape (N,)

Noise related functions

Collection of noise-signals

This module contains the following functions:

  • get_alpha(): normal distributed signal amplitudes with equal power spectral density
  • power_law_noise(): normal distributed signal amplitudes with power spectrum :math:f^lpha
  • power_law_acf(): (theoretical) autocorrelation function of power law noise
  • ARMA(): autoregressive moving average noise process
PyDynamic.misc.noise.get_alpha(color_value=0)[source]

Translate a color (given as string) into an exponent alpha or directly hand through a given numeric value of alpha.

Parameters:color_value (str, int or float) – if string -> check against known colornames -> return alpha if numeric -> directly return value
Returns:alpha
Return type:float
PyDynamic.misc.noise.power_law_noise(N=None, w=None, color_value='white', mean=0.0, std=1.0)[source]

Generate colored noise by * generate white gaussian noise * multiplying its Fourier-transform with f^(alpha/2) * inverse Fourier-transform to yield the colored/correlated noise * further adjustments to fit to specified mean/std

based on [Zhivomirov2018](A Method for Colored Noise Generation)

Parameters:
  • N (int) – length of noise to be generated
  • w (numpy.ndarray) – user-defined white noise if provided, N is ignored!
  • color_value (str, int or float) – if string -> check against known colornames if numeric -> used as alpha to shape PSD
  • mean (float) – mean of the output signal
  • std (float) – standard deviation of the output signal
Returns:

w_filt

Return type:

filtered noise signal

PyDynamic.misc.noise.power_law_acf(N, color_value='white', std=1.0)[source]

Return the theoretic right-sided autocorrelation (Rww) of different colors of noise.

Colors of noise are defined to have a power spectral density (Sww) proportional to f^lpha. Sww and Rww form a Fourier-pair. Therefore Rww = ifft(Sww).

PyDynamic.misc.noise.ARMA(length, phi=0.0, theta=0.0, std=1.0)[source]

Generate time-series of a predefined ARMA-process based on this equation: \(\sum_{j=1}^{\min(p,n-1)} \phi_j \epsilon[n-j] + \sum_{j=1}^{\min(q,n-1)} \theta_j w[n-j]\) where w is white gaussian noise. Equation and algorithm taken from [Eichst2012] .

Parameters:
  • length (int) – how long the drawn sample will be
  • phi (float, list or numpy.ndarray, shape (p, )) – AR-coefficients
  • theta (float, list or numpy.ndarray) – MA-coefficients
  • std (float) – std of the gaussian white noise that is feeded into the ARMA-model
Returns:

e – time-series of the predefined ARMA-process

Return type:

numpy.ndarray, shape (length, )

References

Miscellaneous useful helper functions

The PyDynamic.misc.tools module is a collection of miscellaneous helper functions.

This module contains the following functions:

PyDynamic.misc.tools.print_mat(matrix, prec=5, vertical=False, retS=False)[source]

Print matrix (2D array) to the console or return as formatted string

Parameters:
  • matrix ((M,N) array_like) –
  • prec (int) – the precision of the output
  • vertical (bool) – print out vertical or not
  • retS (bool) – print or return string
Returns:

s – if retS is True

Return type:

str

PyDynamic.misc.tools.print_vec(vector, prec=5, retS=False, vertical=False)[source]

Print vector (1D array) to the console or return as formatted string

Parameters:
  • vector ((M,) array_like) –
  • prec (int) – the precision of the output
  • vertical (bool) – print out vertical or not
  • retS (bool) – print or return string
Returns:

s – if retS is True

Return type:

str

PyDynamic.misc.tools.make_semiposdef(matrix, maxiter=10, tol=1e-12, verbose=False)[source]

Make quadratic matrix positive semi-definite by increasing its eigenvalues

Parameters:
  • matrix ((N,N) array_like) – the matrix to process
  • maxiter (int) – the maximum number of iterations for increasing the eigenvalues
  • tol (float) – tolerance for deciding if pos. semi-def.
  • verbose (bool) – If True print smallest eigenvalue of the resulting matrix
Returns:

quadratic positive semi-definite matrix

Return type:

(N,N) array_like

PyDynamic.misc.tools.FreqResp2RealImag(Abs, Phase, Unc, MCruns=10000.0)[source]

Calculate real and imaginary parts from frequency response

Calculate real and imaginary parts from amplitude and phase with associated uncertainties.

Parameters:
  • Abs ((N,) array_like) – amplitude values
  • Phase ((N,) array_like) – phase values in rad
  • Unc ((2N, 2N) or (2N,) array_like) – uncertainties
  • MCruns (bool) – Iterations for Monte Carlo simulation
Returns:

  • Re, Im ((N,) array_like) – real and imaginary parts (best estimate)
  • URI ((2N, 2N) array_like) – uncertainties assoc. with Re and Im

PyDynamic.misc.tools.make_equidistant(t, y, uy, dt=0.05, kind='linear')[source]

Interpolate non-equidistant time series to equidistant

Interpolate measurement values and propagate uncertainties accordingly.

Parameters:
  • t ((N,) array_like) – timestamps (or frequencies)
  • y ((N,) array_like) – corresponding measurement values
  • uy ((N,) array_like) – corresponding measurement values’ standard uncertainties
  • dt (float, optional) – desired interval length
  • kind (str, optional) – Specifies the kind of interpolation for the measurement values as a string (‘previous’, ‘next’, ‘nearest’ or ‘linear’).
Returns:

  • t_new ((M,) array_like) – interpolation timestamps (or frequencies)
  • y_new ((M,) array_like) – interpolated measurement values
  • uy_new ((M,) array_like) – interpolated measurement values’ standard uncertainties

References

PyDynamic.misc.tools.trimOrPad(array, length, mode='constant')[source]

Trim or pad (with zeros) a vector to desired length

PyDynamic.misc.tools.progress_bar(count, count_max, width=30, prefix='', done_indicator='#', todo_indicator='.', fout=<_io.TextIOWrapper name='<stdout>' mode='w' encoding='UTF-8'>)[source]

A simple and reusable progress-bar

Parameters:
  • count (int) – current status of iterations, assumed to be zero-based
  • count_max (int) – total number of iterations
  • width (int, optional) – width of the actual progressbar (actual printed line will be wider)
  • prefix (str, optional) – some text that will be printed in front of the bar (i.e. “Progress of ABC:”)
  • done_indicator (str, optional) – what character is used as “already-done”-indicator
  • todo_indicator (str, optional) – what character is used as “not-done-yet”-indicator
  • fout (file-object, optional) – where the progress-bar should be written/printed to