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

A collection of functions to deal with second order dynamic systems

This module is used throughout PyDynamic and is specialized for second order dynamic systems, such as the ones used for high-class accelerometers.

This module contains the following functions:

  • sos_absphase(): Propagation of uncertainty from physical parameters to real and imaginary part of system’s transfer function using GUM S2 Monte Carlo

  • sos_FreqResp(): Calculation of the system frequency response

  • sos_phys2filter(): Calculation of continuous filter coefficients from physical parameters

  • 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_absphase(S, d, f0, uS, ud, uf0, f, runs=10000)[source]

Propagation of uncertainty from physical parameters to amplitude and phase

Propagation of uncertainties from physical parameters to amplitude and phase of system’s transfer function is performed 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 amplitude and phase

  • runs (int, optional) – number of Monte Carlo runs

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_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_realimag(S, d, f0, uS, ud, uf0, f, runs=10000)[source]

Propagation of uncertainty from physical parameters to real and imaginary part

Propagation of uncertainties from physical parameters to real and imaginary part of system’s transfer function is performed 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

  • runs (int, optional) – number of Monte Carlo runs

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

This module is a collection of functions 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 an IIR filter with certain coefficients 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.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.isstable(b, a, ftype='digital')[source]

Determine whether an IIR filter with certain coefficients is stable

Determine whether IIR filter with coefficients b and a is stable by checking the roots of the polynomial a.

Parameters
  • b (ndarray) – Filter numerator coefficients. These are only part of the input parameters for compatibility reasons (especially with MATLAB code). During the computation they are actually not used.

  • a (ndarray) – Filter denominator coefficients.

  • ftype (string, optional) – Filter type. ‘digital’ if in discrete-time (default) and ‘analog’ if in continuous-time.

Returns

stable – Whether filter is stable or not.

Return type

bool

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.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.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

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

Shortcut for calculation of unwrapped angle of complex values

Test signals

A collection of test signals which can be used to simulate dynamic measurements

This module contains the following functions:

  • GaussianPulse(): Generate a Gaussian pulse at t0 with height m0 and std sigma

  • multi_sine(): Generate a multi-sine signal as summation of single sine signals

  • rect(): Rectangular signal of given height and width \(t_1 - t_0\)

  • shocklikeGaussian(): Generate a signal that resembles a shock excitation as a Gaussian

  • sine(): Generate a sine signal

  • squarepulse(): Generates a series of rect functions to represent a square pulse signal

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

Generate 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,)

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

Base class for generation of a correlated noise process

PyDynamic.misc.testsignals.multi_sine(time, amps, freqs, noise=0.0)[source]

Generate a batch of a summation of sine signals with normally distributed noise

Parameters
  • time (np.ndarray of shape (N,)) – time instants

  • amps (list or np.ndarray of shape (M,) of floating point values) – amplitudes of the sine signals

  • freqs (list or np.ndarray of shape (M,) of floating point values) – frequencies of the sine signals in Hz

  • 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,)

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.shocklikeGaussian(time, t0, m0, sigma, noise=0.0)[source]

Generate a signal that resembles a shock excitation as a Gaussian

The main shock is 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.sine(time, amp=1.0, freq=1.0, noise=0.0)[source]

Generate a batch of a sine signal with normally distributed noise

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 = 1.0)

  • 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,)

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,)

Noise related functions

Collection of noise-signals

This module contains the following functions:

  • ARMA(): autoregressive moving average noise process

  • get_alpha(): normal distributed signal amplitudes with equal power spectral density

  • power_law_acf(): The theoretic right-sided autocorrelation (Rww) of different colors of noise

  • power_law_noise(): normal distributed signal amplitudes with power spectrum \(f^\alpha\)

  • power_law_acf(): (theoretical) autocorrelation function of power law noise

  • white_gaussian(): Draw random samples from a normal (Gaussian) distribution

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

Generate time-series of a predefined ARMA-process

The process is generated 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
Returns

e – time-series of the predefined ARMA-process

Return type

np.ndarray of shape (length, )

References

PyDynamic.misc.noise.get_alpha(color_value=0)[source]

Return the matching alpha for the provided color value

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 color names -> return alpha

  • if numeric -> directly return value

Returns

alpha – exponent alpha or directly numeric value

Return type

float

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

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 math:f^alpha. Sww and Rww form a Fourier-pair. Therefore Rww = ifft(Sww).

PyDynamic.misc.noise.power_law_noise(N=None, w=None, color_value='white', mean=0.0, std=1.0)[source]

Generate colored noise

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] .

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 color names 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 – filtered noise signal

Return type

np.ndarray

PyDynamic.misc.noise.white_gaussian(N, mean=0, std=1)[source]

Draw random samples from a normal (Gaussian) distribution

Parameters
  • N (int or tuple of ints) – Output shape. If the given shape is, e.g., (m, n, k), then m * n * k samples are drawn.

  • mean (float or array_like of floats) – Mean (“centre”) of the distribution.

  • std (float or array_like of floats) – Standard deviation (spread or “width”) of the distribution. Must be non-negative.

Returns

Drawn samples from the parameterized normal distribution.

Return type

np.ndarray

Miscellaneous useful helper functions

A collection of miscellaneous helper functions

This module contains the following functions:

PyDynamic.misc.tools.FreqResp2RealImag(Abs: numpy.ndarray, Phase: numpy.ndarray, Unc: numpy.ndarray, MCruns: Optional[int] = 1000)[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 either as full covariance matrix or as its main diagonal

  • MCruns (int, optional) – number of iterations for Monte Carlo simulation, defaults to 1000

Returns

  • Re, Im ((N,) array_like) – best estimate of real and imaginary parts

  • URI ((2N, 2N) array_like) – uncertainties assoc. with Re and Im

PyDynamic.misc.tools.complex_2_real_imag(array: numpy.ndarray) numpy.ndarray[source]

Take an array of any non-flexible scalar dtype to return real and imaginary part

The input array \(x \in \mathbb R^m\) is reassembled to the form of the expected input of some of the functions in the modules propagate_DFT and fit_filter: \(y = \left( \operatorname{Re}(x), \operatorname{Im}(x) \right)\).

Parameters

array (np.ndarray of shape (M,)) – the array to assemble the version with real and imaginary parts from

Returns

the array of real and imaginary parts

Return type

np.ndarray of shape (2M,)

PyDynamic.misc.tools.is_2d_matrix(ndarray: numpy.ndarray) bool[source]

Check if a np.ndarray is a matrix, i.e. is of shape (n,m)

Parameters

ndarray (np.ndarray) – the array to check

Returns

True, if the array expands over exactly two dimensions, False otherwise

Return type

bool

PyDynamic.misc.tools.is_2d_square_matrix(ndarray: numpy.ndarray) bool[source]

Check if a np.ndarray is a two-dimensional square matrix, i.e. is of shape (n,n)

Parameters

ndarray (np.ndarray) – the array to check

Returns

True, if the array expands over exactly two dimensions of similar size, False otherwise

Return type

bool

PyDynamic.misc.tools.is_vector(ndarray: numpy.ndarray) bool[source]

Check if a np.ndarray is a vector, i.e. is of shape (n,)

Parameters

ndarray (np.ndarray) – the array to check

Returns

True, if the array expands over one dimension only, False otherwise

Return type

bool

PyDynamic.misc.tools.make_equidistant(*args, **kwargs)[source]

Deprecated since version 2.0.0: Please use PyDynamic.uncertainty.interpolate.interp1d_unc()

PyDynamic.misc.tools.make_semiposdef(matrix: numpy.ndarray, maxiter: Optional[int] = 10, tol: Optional[float] = 1e-12, verbose: Optional[bool] = False) numpy.ndarray[source]

Make quadratic matrix positive semi-definite by increasing its eigenvalues

Parameters
  • matrix (array_like of shape (N,N)) – the matrix to process

  • maxiter (int, optional) – the maximum number of iterations for increasing the eigenvalues, defaults to 10

  • tol (float, optional) – tolerance for deciding if pos. semi-def., defaults to 1e-12

  • verbose (bool, optional) – If True print smallest eigenvalue of the resulting matrix, if False (default) be quiet

Returns

quadratic positive semi-definite matrix

Return type

(N,N) array_like

Raises

ValueError – If matrix is not square.

PyDynamic.misc.tools.normalize_vector_or_matrix(numbers: numpy.ndarray) numpy.ndarray[source]

Scale an array of numbers to the interval between zero and one

If all values in the array are the same, the output array will be constant zero.

Parameters

numbers (np.ndarray) – the numpy.ndarray to normalize

Returns

the normalized array

Return type

np.ndarray

PyDynamic.misc.tools.number_of_rows_equals_vector_dim(matrix: numpy.ndarray, vector: numpy.ndarray) bool[source]

Check if a matrix has the same number of rows as a vector

Parameters
  • matrix (np.ndarray) – the matrix, that is supposed to have the same number of rows

  • vector (np.ndarray) – the vector, that is supposed to have the same number of elements

Returns

True, if the number of rows coincide, False otherwise

Return type

bool

PyDynamic.misc.tools.plot_vectors_and_covariances_comparison(vector_1: numpy.ndarray, vector_2: numpy.ndarray, covariance_1: numpy.ndarray, covariance_2: numpy.ndarray, title: Optional[str] = 'Comparison between two vectors and corresponding uncertainties', label_1: Optional[str] = 'vector_1', label_2: Optional[str] = 'vector_2')[source]

Plot two vectors and their covariances side-by-side for visual comparison

Parameters
  • vector_1 (np.ndarray) – the first vector to compare

  • vector_2 (np.ndarray) – the second vector to compare

  • covariance_1 (np.ndarray) – the first covariance matrix to compare

  • covariance_2 (np.ndarray) – the second covariance matrix to compare

  • title (str, optional) – the title for the comparison plot, defaults to “Comparison between two vectors and corresponding uncertainties”

  • label_1 (str, optional) – the label for the first vector in the legend and title for the first covariance plot, defaults to “vector_1”

  • label_2 (str, optional) – the label for the second vector in the legend and title for the second covariance plot, defaults to “vector_2”

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.progress_bar(count, count_max, width: Optional[int] = 30, prefix: Optional[str] = '', done_indicator: Optional[str] = '#', todo_indicator: Optional[str] = '.', fout: Optional = None)[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), default to 30

  • prefix (str, optional) – some text that will be printed in front of the bar (i.e. “Progress of ABC:”), if not given only progressbar itself will be printed

  • done_indicator (str, optional) – what character is used as “already-done”-indicator, defaults to “#”

  • todo_indicator (str, optional) – what character is used as “not-done-yet”-indicator, defaults to “.”

  • fout (file-object, optional) – where the progress-bar should be written/printed to, defaults to direct print to stdout

PyDynamic.misc.tools.real_imag_2_complex(array: numpy.ndarray) numpy.ndarray[source]

Take a np.ndarray with real and imaginary parts and return dtype complex ndarray

The input array \(x \in \mathbb R^{2m}\) representing a complex vector \(y \in \mathbb C^m\) has the form of the expected input of some of the functions in the modules propagate_DFT and fit_filter: \(x = \left( \operatorname{Re}(y), \operatorname{Im}(y) \right)\) or a np.ndarray containing several of these.

Parameters

array (np.ndarray of shape (N,2M) or of shape (2M,)) – the array of any integer or floating dtype to assemble the complex version of

Returns

the complex array

Return type

np.ndarray of shape (N,M) or of shape (M,)

PyDynamic.misc.tools.separate_real_imag_of_mc_samples(array: numpy.ndarray) List[numpy.ndarray][source]

Split a np.ndarray containing MonteCarlo samples real and imaginary parts

The input array \(x \in \mathbb R^{n \times 2m}\) representing an n-elemental array of complex vectors \(y_i \in \mathbb C^m\) has the form of the expected input of some of the functions in the modules propagate_DFT and fit_filter: \(x = \left( \operatorname{Re}(y_i), \operatorname{Im}(y_i) \right)_{i=1,\ldots,n}\).

Parameters

array (np.ndarray of shape (N,2M)) – the array of any integer or floating dtype to assemble the complex version of

Returns

two-element list of the two arrays containing the real and imaginary parts

Return type

list of two np.ndarrays of shape (N,M)

PyDynamic.misc.tools.separate_real_imag_of_vector(vector: numpy.ndarray) List[numpy.ndarray][source]

Split a np.ndarray containing real and imaginary parts into half

The input array \(x \in \mathbb R^{2m}\) representing a complex vector \(y \in \mathbb C^m\) has the form of the expected input of some of the functions in the modules propagate_DFT and fit_filter: \(x = \left( \operatorname{Re}(y), \operatorname{Im}(y) \right)\).

Parameters

vector (np.ndarray of shape (2M,)) – the array of any integer or floating dtype to assemble the complex version of

Returns

two-element list of the two arrays containing the real and imaginary parts

Return type

list of two np.ndarrays of shape (M,)

PyDynamic.misc.tools.shift_uncertainty(x: numpy.ndarray, ux: numpy.ndarray, shift: int)[source]

Shift the elements in the vector x and associated uncertainties ux

This function uses numpy.roll() to shift the elements in x and ux. See the linked official documentation for details.

Parameters
  • x (np.ndarray of shape (N,)) – vector of estimates

  • ux (float, np.ndarray of shape (N,) or of shape (N,N)) – uncertainty associated with the vector of estimates

  • shift (int) – amount of shift

Returns

  • shifted_x ((N,) np.ndarray) – shifted vector of estimates

  • shifted_ux (float, np.ndarray of shape (N,) or of shape (N,N)) – uncertainty associated with the shifted vector of estimates

Raises

ValueError – If shift, x or ux are of unexpected type, dimensions of x and ux do not fit or ux is of unexpected shape

PyDynamic.misc.tools.trimOrPad(array: Union[List, numpy.ndarray], length: int, mode: Optional[str] = 'constant')[source]

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

Either trim or zero-pad an array to achieve the required length. Both actions are applied to the end of the array.

Parameters
  • array (list, 1D np.ndarray) – original data

  • length (int) – length of output

  • mode (str, optional) – handed over to np.pad, default “constant”

Returns

array_modified – An either trimmed or zero-padded array

Return type

np.ndarray of shape (length,)