pyOMA.core.PreProcessingTools.PreProcessSignals#

class pyOMA.core.PreProcessingTools.PreProcessSignals(signals, sampling_rate, ref_channels=None, accel_channels=None, velo_channels=None, disp_channels=None, setup_name=None, channel_headers=None, start_time=None, F=None, **kwargs)[source]#

Bases: object

A simple pre-processor for signals * load ascii datafiles * specify sampling rate, reference channels * specify geometry, channel-dof-assignments * specify channel quantities such as acceleration, velocity, etc * remove (constant) offsets from signals * decimate signals * compute correlation functions and power spectral densities * filter signals

Subsequent modules of pyOMA (SysId, Modal Analysis, Stabilization, Mode shape visualization) rely on the variables and methods provided by this class.

  • time-step integration of signals

  • Multi-block Blackman-Tukey PSD

__init__(signals, sampling_rate, ref_channels=None, accel_channels=None, velo_channels=None, disp_channels=None, setup_name=None, channel_headers=None, start_time=None, F=None, **kwargs)[source]#

Methods

__init__(signals, sampling_rate[, ...])

add_chan_dofs(chan_dofs)

chan_dofs = [ (chan_num, node_name, az, elev, chan_name) , .

add_noise([amplitude, snr])

corr_blackman_tukey(m_lags[, n_segments, ...])

Estimate the (cross- and auto-) correlation functions (C/ACF), by direct computation of the standard un-biased estimator:

corr_welch([m_lags, n_segments, refs_only])

Estimate the (cross- and auto-) correlation functions (C/ACF), by the inverse Fourier Transform of Power Spectral Densities, estimated according to Welch's method.

correct_offset()

corrects a constant offset from measured signals

correlation([m_lags, method])

A convenience method for obtaining the correlation sequence by the default or any specified estimation method.

decimate_signals(decimate_factor[, nyq_rat, ...])

decimates signals data filter type and order are choosable (order 8 and type cheby1 are standard for scipy signal.decimate function) maximum ripple in the passband (rp) and minimum attenuation in the stop band (rs) are modifiable

filter_signals([lowpass, highpass, ...])

init_from_config(conf_file, meas_file[, ...])

initializes the PreProcessor object with a configuration file

load_chan_dofs(fname)

chan_dofs[i] = (chan_num, node_name, az, elev, chan_name)

load_measurement_file(fname, **kwargs)

A method for loading a signals file

load_state(fname)

precondition_signals([method])

psd([n_lines, method])

A convenience method for obtaining the PSD by the default or any specified estimation method.

psd_blackman_tukey([n_lines, refs_only, window])

Estimate the (cross- and auto-) power spectral densities (PSD), by Fourier Transform of correlation functions estimated according to Blackman-Tukey's method.

psd_welch([n_lines, n_segments, refs_only, ...])

Estimate the (cross- and auto-) power spectral densities (PSD), according to Welch's method.

save_state(fname)

sv_psd([n_lines])

Compute the singular values of the power spectral density matrices, for which the complete (all cross spectral densities) matrices are used.

take_chan_dof(chan, node, dof)

validate_channels(channels[, quant_check])

welch(n_lines, **kwargs)

Attributes

accel_channels

corr_matrices

corr_matrix

disp_channels

dt

duration

freqs

returns: freqs -- Array with the frequency lines corresponding to the spectral values :rtype: np.ndarray (n_lines, )

freqs_bt

returns: freqs -- Array with the frequency lines corresponding to the spectral values :rtype: np.ndarray (n_lines, )

freqs_wl

returns: freqs -- Array with the frequency lines corresponding to the spectral values :rtype: np.ndarray (n_lines, )

lags

lags_bt

lags_wl

m_lags

n_lines

n_segments

num_analised_channels

num_ref_channels

psd_matrix

ref_channels

signal_power

signal_rms

t

total_time_steps

velo_channels

add_chan_dofs(chan_dofs)[source]#

chan_dofs = [ (chan_num, node_name, az, elev, chan_name) , … ] This function is not checking if channels or nodes actually exist the former should be added the latter might only be possible, if the geometry object is known to the class

corr_blackman_tukey(m_lags, n_segments=None, refs_only=True, **kwargs)[source]#

Estimate the (cross- and auto-) correlation functions (C/ACF), by direct computation of the standard un-biased estimator:

\[\hat{R}_{fg}[m] = \frac{1}{N - m}\sum_{n=0}^{N - m - 1} f[n] g[n + m]\]

Computes correlation functions of all channels with selected reference channels up to, but excluding, a time lag of m_lags. Normalization is done according to the unbiased estimator, i.e. 0-lag correlation value must be multiplied by n_lines to get the signals cross-power.

Variance estimation for each time lag is performed by dividing the signals into n_segments non-overlapping blocks for individual estimation of correlation functions. With increasing numbers of non-overlapping blocks the confidence intervals of the correlation functions increase (“worsen”), especially at higher lags and short block lengths, due to a larger number of time steps being discarded.

Note that:

m_lags = n_lines // 2 + 1

n_lines = (m_lags - 1) * 2

Parameters:
  • m_lags (integer, optional) – Number of lags (positive). Note: this includes the 0-lag, therefore excludes the “m_lags”-lag.

  • n_segments (integer, optional) – Number of blocks to perform averaging over. If blocks are shorter than m_lags it raises a ValueError.

  • refs_only (bool, optional) – Compute cross-ACFss only with reference channels

  • kwargs – Additional kwargs are currently not used

Returns:

corr_matrix – Array of shape (num_channels, num_ref_channels, m_lags) containing the correlation values of the respective channels and lags

Return type:

np.ndarray

See also

corr_welch

Correlation function estimation by Welch’s method, possibly faster, but distorted for short segments and biased through windowing.

corr_welch(m_lags=None, n_segments=None, refs_only=True, **kwargs)[source]#

Estimate the (cross- and auto-) correlation functions (C/ACF), by the inverse Fourier Transform of Power Spectral Densities, estimated according to Welch’s method. Bias due to windowing of the underlying PSD persists. Normalization is done according to the unbiased estimator, i.e. 0-lag correlation value must be multiplied by n_lines to get the signals cross-power.

Note that:

m_lags = n_lines // 2 + 1

n_lines = (m_lags - 1) * 2

N_segment = N // n_segments

Parameters:
  • m_lags (integer, optional) – Total number of lags (positive). Note: this includes the 0-lag, therefore excludes the m_lags-lag.

  • n_segments (integer, optional) – Number of segments to perform averaging over resulting segment length must be smaller or equal n_lines

  • refs_only (bool, optional) – Compute cross-ACFss only with reference channels

  • kwargs – Additional kwargs are passed to self.psd_welch and further

Returns:

corr_matrix – Array of shape (num_channels, num_ref_channels, m_lags) containing the correlation values of the respective channels and lags

Return type:

np.ndarray

See also

psd_welch

PSD estimation algorithm used by this method.

Todo

  • deconvolve window (if possible)

correct_offset()[source]#

corrects a constant offset from measured signals

..TODO::
  • remove linear, … ofsets as well

correlation(m_lags=None, method=None, **kwargs)[source]#

A convenience method for obtaining the correlation sequence by the default or any specified estimation method.

Parameters:
  • m_lags (integer, optional) – Number of lags (positive). Note: this includes the 0-lag, therefore excludes the m_lags-lag.

  • method (str, optional) – The method to use for spectral estimation

  • kwargs – Additional parameters are passed to the spectral estimation method

Returns:

corr_matrix – Array of shape (num_channels, num_ref_channels, m_lags) containing the correlation values of the respective channels and lags

Return type:

np.ndarray

decimate_signals(decimate_factor, nyq_rat=2.5, highpass=None, order=None, filter_type='cheby1')[source]#

decimates signals data filter type and order are choosable (order 8 and type cheby1 are standard for scipy signal.decimate function) maximum ripple in the passband (rp) and minimum attenuation in the stop band (rs) are modifiable

property freqs#

returns: freqs – Array with the frequency lines corresponding to the spectral values :rtype: np.ndarray (n_lines, )

property freqs_bt#

returns: freqs – Array with the frequency lines corresponding to the spectral values :rtype: np.ndarray (n_lines, )

property freqs_wl#

returns: freqs – Array with the frequency lines corresponding to the spectral values :rtype: np.ndarray (n_lines, )

classmethod init_from_config(conf_file, meas_file, chan_dofs_file=None, **kwargs)[source]#

initializes the PreProcessor object with a configuration file

to remove channels at loading time use ‘usecols’ keyword argument if delete_channels are specified, these will be checked against all other channel definitions, which will be adjusted accordingly

static load_chan_dofs(fname)[source]#
chan_dofs[i] = (chan_num, node_name, az, elev, chan_name)

= (int, str, float,float, str)

azimuth angle starting from x axis towards y axis elevation defined from x-y plane up x: 0.0, 0.0 y: 90.0, 0.0 z: 0.0, 90.0 channels not present in the file will be removed later nodes do not have to, but should exist, as this information is also used for merging multiple setups, which does not rely on any “real” geometry

static load_measurement_file(fname, **kwargs)[source]#

A method for loading a signals file

Parameters:

fname (str) – The full path of the signals file

Returns:

  • headers (list of str) – The names of all channels

  • units (list of str) – The units of all channels

  • start_time (datetime.datetime) – The starting time of the measured signal

  • sample_rate (float) – The sample rate, at wich the signal was acquired

  • signals (ndarray) – Array of shape (num_timesteps, num_channels) which contains the acquired signal

psd(n_lines=None, method=None, **kwargs)[source]#

A convenience method for obtaining the PSD by the default or any specified estimation method.

Parameters:
  • n_lines (integer, optional) – Number of frequency lines (positive + negative)

  • method – The method to use for spectral estimation

  • **kwargs – Additional parameters are passed to the spectral estimation method

Returns:

psd_matrix – Array of shape (num_channels, num_ref_channels, n_lines // 2 + 1) containing the power density values of the respective channels and frequencies

Return type:

np.ndarray

psd_blackman_tukey(n_lines=None, refs_only=True, window='hamming', **kwargs)[source]#

Estimate the (cross- and auto-) power spectral densities (PSD), by Fourier Transform of correlation functions estimated according to Blackman-Tukey’s method. Non-negativeness of the PSD is ensured by using a lag window, i.e. convolving the temporal window with itself. Normalization is applied w.r.t. conservation of energy, i.e. magnitudes will change with n_lines but power stays constant.

Note that:

m_lags = n_lines // 2 + 1 n_lines = (m_lags - 1) * 2

Parameters:
  • n_lines (integer, optional) – Number of frequency lines (positive + negative)

  • refs_only (bool, optional) – Compute cross-PDSs only with reference channels

  • window (str or tuple or array_like, optional) – Desired temporal window to be applied to the correlation sequence after conversion to a lag window by “self-convolution” See scipy.signal.get_window() for more information

  • kwargs – Additional kwargs are passed to self.corr_blackman_tukey

Returns:

psd_matrix – Array of shape (num_channels, num_ref_channels, n_lines // 2 + 1) containing the power density values of the respective channels and frequencies

Return type:

np.ndarray

psd_welch(n_lines=None, n_segments=None, refs_only=True, window='hamming', **kwargs)[source]#

Estimate the (cross- and auto-) power spectral densities (PSD), according to Welch’s method. No overlapping is allowed (deliberately to ensure statistical independence of blocks for variance estimation). Segments are n_lines // 2 long and zero padded to n_lines to allow estimation of the full correlation sequence, which is twice as long as the input signal. Normalization is applied w.r.t. conservation of energy, i.e. magnitudes will change with n_lines but power stays constant.

Parameters:
  • n_lines (integer, optional) – Number of frequency lines (positive + negative)

  • n_segments (integer, optional) – Number of segments to perform averaging over resulting segment length must be smaller or equal n_lines

  • refs_only (bool, optional) – Compute cross-PDSs only with reference channels

  • window (str or tuple or array_like, optional) – Desired window to use. See scipy.signal.get_window() for more information

  • kwargs – Additional kwargs are passed to scipy.signals.csd

Returns:

psd_matrix – Array of shape (num_channels, num_ref_channels, n_lines // 2 + 1) containing the power density values of the respective channels and frequencies

Return type:

np.ndarray

sv_psd(n_lines=None, **kwargs)[source]#

Compute the singular values of the power spectral density matrices, for which the complete (all cross spectral densities) matrices are used.

Parameters:
  • n_lines (integer, optional) – Number of frequency lines (positive + negative)

  • kwargs – Additional parameters are passed to the spectral estimation method