eqsig package

The eqsig package is based around the AccSignal and Signal objects to represent acceleration time series and other time series data. These objects hold commonly used signal parameters and use caching to avoid expensive recalculations.

Single time series objects

class eqsig.single.AccSignal(values, dt, label='m1', smooth_freq_range=(0.1, 30), smooth_fa_freqs=None, verbose=0, response_period_range=(0.1, 5), response_times=None, ccbox=0)[source]

Bases: Signal

A time series object

Parameters
  • values (array_like) – Acceleration time series - should be in m/s2

  • dt (float) – Time step

  • label (str (default='m1')) – Used for plotting

  • smooth_freq_range (tuple [floats] (default=(0.1, 30))) – lower and upper bound of frequency range to compute the smoothed Fourier amplitude spectrum

  • verbose (int (default=0)) – Level of output verbosity

  • response_times (tuple of floats (default=(0.1, 5))) –

  • ccbox (int) – colour index for plotting

clear_cache()[source]

Resets the dynamically calculated properties.

correct_me()[source]

This provides a correction to an acceleration time series

property displacement

Displacement time series

gen_response_spectrum(response_times=None, xi=-1, min_dt_ratio=4)[source]

Generate the response spectrum for the response_times for a given damping (xi).

Parameters
  • response_times (array_like) – Time periods of SDOFs that response should be computed for

  • xi (float (default=0.05 or previously set)) – Ratio of critical damping

  • min_dt_ratio (float) – Smallest ratio between response period and integration time step, acceleration time series will be interpolated to reach appropriate integration time step if not achieved using the acceleration time series time step.

generate_all_motion_stats()[source]

Deprecated: Use eqsig.im functions to calculate stats

generate_cumulative_stats()[source]

Deprecated: Use eqsig.im.calc_arias_intensity, eqsig.im.calc_cav

generate_displacement_and_velocity_series(trap=True)[source]

Calculates the displacement and velocity time series

generate_duration_stats()[source]

Deprecated: Use eqsig.im.calc_sig_dur or eqsig.im.calc_sig_dur, eqsig.im.calc_brac_dur or esig.im.calc_a_rms

generate_peak_values()[source]

DEPRECATED: all peak values are lazy loaded - this method does not need to be run.

generate_response_spectrum(response_times=None, xi=-1, min_dt_ratio=4)[source]
property pga

Absolute peak ground acceleration

property pgd

Absolute peak ground displacement

property pgv

Absolute peak ground velocity

rebase_displacement()[source]

Correction to make the displacement zero at the end of the record

remove_rolling_average(mtype='velocity', freq_window=5)[source]

Removes the rolling average

Parameters
  • mtype (str) – motion type to apply method to

  • freq_window (int) – window for applying the rolling average

reset_all_motion_stats()[source]
response_series(response_times=None, xi=-1)[source]

Generate the response time series for a set of elastic SDOFs for a given damping (xi).

property s_a

Pseudo maximum response accelerations of linear SDOFs

property s_d

Maximum response displacements of linear SDOFs

property s_v

Pseudo maximum response velocities of linear SDOFs

set_zero_residual_displacement(timezone=None)[source]
set_zero_residual_displacement_and_velocity(timezone=None)[source]
set_zero_residual_velocity(timezone=None)[source]
property velocity

Velocity time series

class eqsig.single.Signal(values, dt, label='m1', smooth_freq_range=(0.1, 30), smooth_fa_freqs=None, verbose=0, ccbox=0)[source]

Bases: object

A time series object

Parameters
  • values (array_like) – values of the time series

  • dt (float) – the time step

  • label (str, optional) – A name for the signal

  • smooth_freq_range (tuple, optional) – The frequency range for computing the smooth FAS

  • verbose (int, optional) – Level of console output

  • ccbox (int, optional) – colour id for plotting

add_constant(constant)[source]

Adds a single value to every value in the signal.

Parameters

constant (float) – Value to be added to time series

add_series(series)[source]

Adds a series of values to the values in the signal.

Parameters

series (array_like) – A series of values

add_signal(new_signal)[source]

Combines a signal.

Parameters

new_signal (Signal object) –

butter_pass(cut_off=(0.1, 15), **kwargs)[source]

Performs a Butterworth filter

Notes

Wraps the scipy ‘butter’ filter

Parameters
  • cut_off (tuple) – Lower and upper cut off frequencies for the filter, if None then no filter. e.g. (None, 15) applies a lowpass filter at 15Hz, whereas (0.1, 10) applies a bandpass filter at 0.1Hz to 10Hz.

  • filter_order (int) – Order of the Butterworth filter

  • remove_gibbs (str (default=None)) – Pads with zeros to remove the Gibbs filter effect, if = ‘start’ then pads at start, if = ‘end’ then pads at end, if = ‘mid’ then pads half at start and half at end

  • gibbs_extra (int) – each increment of the value doubles the record length using zero padding.

:param : :type : return:

clear_cache()[source]

Resets the dynamically calculated properties.

property dt

The time step

property fa_freqs
property fa_frequencies
property fa_spectrum

The one-sided Fourier Amplitude spectrum

property fa_spectrum_abs

The one-sided Fourier Amplitude spectrum

gen_fa_spectrum(p2_plus=0, n=None)[source]

Sets the one-sided Fourier Amplitude spectrum and frequencies

Parameters
  • p2_plus (int (default=0)) – Additional increments of n such that length of signal used is n and n=ceil(log2(self.npts) + p2_plus)).

  • n (int (optional)) – If specified the signal length is padded with zeros before computing FAS

gen_smooth_fa_spectrum(smooth_fa_freqs=None, band=40)[source]

Calculates the smoothed Fourier Amplitude Spectrum using the method by Konno and Ohmachi 1998

Parameters

band (int) – range to smooth over

generate_fa_spectrum()[source]
generate_smooth_fa_spectrum(band=40)[source]
get_section_average(start=0, end=-1, index=False)[source]

Gets the average value of a part of series.

Common use is so that it can be patched with another record.

Parameters
  • start (int or float, optional) – Section start point

  • end (int or float, optional) – Section end point

  • index (bool, optional) – if False then start and end are considered values in time.

Returns

The mean value of the section.

Return type

float

property npts

The number of points in the time series

remove_average(section=-1, verbose=-1)[source]

Calculates the average and removes it from the record

remove_poly(poly_fit=0)[source]

Calculates best fit polynomial and removes it from the record

reset_values(new_values)[source]
running_average(width=1)[source]

Averages the values over a width (chunk of numbers) replaces the value in the centre of the chunk.

Parameters

width (int) – The number of points over which values are averaged

set_smooth_fa_frequecies_by_range(limits, n_points)[source]
property smooth_fa_freqs
property smooth_fa_frequencies
property smooth_fa_spectrum
property smooth_freq_points
property smooth_freq_range
property time

An array of time of equal length to the time series

property values

Multiple time series objects

class eqsig.multiple.Cluster(values, dt, names=None, master_index=0, stypes='custom', **kwargs)[source]

Bases: object

This object represents a group of Signals or AccSignals that have the same time step

Parameters
  • values (2d_array_like) – An array containing arrays of values that represent time series

  • dt (float) – Time step of time series

  • names (list) – A list of names for each time series

  • master_index (int) – The index of the master time series - used in some signal matching methods

  • stypes (str or list of str) – The signal type for each time series, if ‘acc’ then AccSignal, else Signal.

  • kwargs

calculate_ratios()[source]
combine_motions(f_ch, low_index=0, high_index=1, **kwargs)[source]

This method combined the two ground motions by taking the high frequency of one values and the low frequency of the other. WARNING: records must have the same time step!!!

:param : :type : param f_ch: is the frequency change point :param : :type : param low_index: – :param : :type : param high_index: – :param : :type : param order: refers to the order of the BW filter :param combined values is returned as self.combo:

generate_response_spectrums()[source]
property n_signals
name_by_index(index)[source]
same_start(**kwargs)[source]

Starts both time series at the same y-value

signal_by_index(index)[source]
property time
time_match(**kwargs)[source]

This method determines the time lag between two records and removes the lag - the solution is based on a sum of the squared residuals. - The motions are both set to the base values, default = 0 - must have same time step and same magnitude

values(name)[source]
values_by_index(index)[source]
eqsig.multiple.combine_at_angle(acc_sig_ns, acc_sig_we, angle)[source]
eqsig.multiple.compute_rotated(acc_sig_ns, acc_sig_we, angle_off_ns=0.0, parameter=None, func=None, points=100)[source]

Computes the rotated value of a parameter.

Parameters
  • acc_sig_ns

  • acc_sig_we

  • angle_off_ns – Angle from North in degrees of the primary signal.

Returns

tuple, (angle to rotate, rotated)

Intensity measure functions

eqsig.im.calc_a_rms(asig, threshold)[source]

DEPRECATED

eqsig.im.calc_acc_rms(asig, threshold)[source]

Root mean squared acceleration

eqsig.im.calc_arias_intensity(acc_sig)[source]

Calculates the Arias Intensity [m/s]

Parameters

acc_sig (eqsig.AccSignal) –

Returns

A time series of the build up of Arias Intensity

Return type

array_like

eqsig.im.calc_bandwidth_f_max(asig, ratio=0.707)[source]

Upper frequency of smooth Fourier spectrum bandwidth

Parameters
  • asig (eqsig.AccSignal) – Acceleration time series object

  • ratio (float) – ratio of maximum value where bandwidth should be computed

Return type

float

eqsig.im.calc_bandwidth_f_min(asig, ratio=0.707)[source]

Lower frequency of smooth Fourier spectrum bandwidth

Parameters
  • asig (eqsig.AccSignal) – Acceleration time series object

  • ratio (float) – ratio of maximum value where bandwidth should be computed

Return type

float

eqsig.im.calc_bandwidth_freqs(asig, ratio=0.707)[source]

Lower and upper frequencies of smooth Fourier spectrum bandwidth

Parameters
  • asig (eqsig.AccSignal) – Acceleration time series object

  • ratio – ratio of maximum value where bandwidth should be computed

Returns

(lower, upper)

Return type

tuple

eqsig.im.calc_brac_dur(asig, threshold, se=False)[source]

Calculates the Bracketed duration based on some threshold

Parameters
  • asig (eqsig.AccSignal) – Acceleration time series object

  • threshold (float) – acceleration threshold to calculation duration start and end

  • se (bool, default=False) – If true then return the start and end times

Return type

float

eqsig.im.calc_bracketed_duration(asig, threshold)[source]

DEPRECATED: Use calc_brac_dur

eqsig.im.calc_cav(acc_sig)[source]

Calculates the Cumulative Absolute velocity

ref: Electrical Power Research Institute. Standardization of the Cumulative Absolute Velocity. 1991. EPRI TR-100082-1’2, Palo Alto, California.

eqsig.im.calc_cav_dp(asig)[source]

Calculates standardized cumulative absolute velocity

ref: Campbell KW, Bozorgnia Y. Predictive equations for the horizontal component of standardized cumulative absolute velocity as adapted for use in the shutdown of U.S. nuclear power plants. Nucl Eng Des 2011;241:2558-69.

Parameters

asig

Returns

eqsig.im.calc_cumulative_abs_displacement(asig)[source]

Cumulative absolute displacement - identical to integral of absolute velocity

eqsig.im.calc_cyc_amp_array_w_power_law(values, n_cyc, b)[source]

Calculate the equivalent uniform loading for a given number of cycles and power coefficient (b)

Parameters
  • values – array_like Time series of values

  • n_cyc

  • b

Returns

eqsig.im.calc_cyc_amp_combined_arrays_w_power_law(values0, values1, n_cyc, b)[source]

Computes the equivalent cyclic amplitude.

For a set number of cycles using a power law and assuming both components act equally

Parameters
  • values0 (array_like) – Time series of values

  • values1 (array_like) – Time series of values in orthogonal direction to values0

  • n_cyc (int or float) – Number of cycles

  • b (float) – Power law factor

Return type

array_like

eqsig.im.calc_cyc_amp_gm_arrays_w_power_law(values0, values1, n_cyc, b)[source]

Calculate the geometric mean equivalent uniform amplitude for a given number of cycles and power coefficient (b)

Parameters
  • values0 – array_like Time series of values

  • values1 – array_like Time series of values in orthogonal direction to values0

  • n_cycles

  • b

Returns

eqsig.im.calc_integral_of_abs_acceleration(asig)[source]

Integral of absolute acceleration

eqsig.im.calc_integral_of_abs_velocity(asig)[source]

Integral of absolute velocity - identical to cumulative absolute displacement

eqsig.im.calc_isv(acc_sig)[source]

Calculates the integral of the squared velocity

See Kokusho (2013) :return:

eqsig.im.calc_max_velocity_period(asig)[source]
eqsig.im.calc_n_cyc_array_w_power_law(values, a_ref, b, cut_off=0.01)[source]

Calculates the equivalent number of uniform amplitude cycles using a power law

Parameters
  • values (array_like) – Time series of values

  • a_ref (float) – Reference uniform amplitude

  • b (float or array_like) – Power law factor

  • cut_off (float) – Low amplitude cutoff value

Return type

array_like

eqsig.im.calc_peak(motion)[source]

Calculates the peak absolute response

eqsig.im.calc_sig_dur(asig, start=0.05, end=0.95, im=None, se=False)[source]

Computes the significant duration according to Trifunac and Brady (1975).

Parameters
  • asig (eqsig.AccSignal) – acceleration time series object

  • dt (float) – time step

  • start (float, default=0.05) – threshold to start the duration

  • end (float, default=0.95) – threshold to end the duration

  • im (function or None (default=None)) – A function that calculates a cumulative intensity measure, if =None, then use eqsig.im.calc_arias_intensity

  • se (bool, default=False) – If true then return the start and end times

Return type

tuple (start_time, end_time)

eqsig.im.calc_sig_dur_vals(motion, dt, start=0.05, end=0.95, se=False)[source]

Computes the significant duration using cumulative acceleration according to Trifunac and Brady (1975).

Parameters
  • motion (array-like) – acceleration time series

  • dt (float) – time step

  • start (float, default=0.05) – threshold to start the duration

  • end (float, default=0.95) – threshold to end the duration

  • se (bool, default=False) – If true then return the start and end times

Return type

tuple (start_time, end_time)

eqsig.im.calc_significant_duration(motion, dt, start=0.05, end=0.95)[source]

Deprecated. Use calc_sig_dur_vals

Parameters
  • motion

  • dt

  • start

  • end

eqsig.im.calc_sir(acc_sig)[source]

Calculates the shaking intensity rate

ref: Dashti, S., Bray, J. D., Pestana, J. S., Riemer, M., and Wilson, D., 2010. Centrifuge testing to evaluate and mitigate liquefaction-induced building settlement mechanisms, ASCE Journal of Geotechnical and Geoenv. Eng. 136, 918-929

Parameters

acc_sig (eqsig.AccSignal) – acceleration signal

Return type

float

eqsig.im.calc_unit_kinetic_energy(acc_signal)[source]

Calculates the cumulative absolute change in kinetic energy for a unit volume of soil with unit mass

Parameters

acc_signal (eqsig.AccSignal) –

eqsig.im.calculate_peak(motion)[source]

Calculates the peak absolute response

eqsig.im.cumulative_response_spectra(acc_signal, fun_name, periods=None, xi=None)[source]
eqsig.im.max_acceleration_period(asig)[source]
eqsig.im.max_fa_period(asig)[source]

Calculates the period corresponding to the maximum value in the Fourier amplitude spectrum

Signal processing and analysis functions

eqsig.functions.calc_fa_spectrum(sig, n=None, p2_plus=None)[source]

Produces the Fourier amplitude spectrum

Parameters

sig (eqsig.Signal) –

Returns

  • fa_spectrum (complex array_like) – Complex values of the spectrum

  • fa_frequencies (array_like) – Frequencies of the spectrum

eqsig.functions.calc_fourier_moment(asig, n)[source]

Original source unknown.

See :cite:`Rathje:2008va`

Parameters
  • asig

  • n

eqsig.functions.calc_roll_av_vals(values, steps, mode='forward')[source]

Calculates the rolling average of a series of values

Parameters
  • values (array_like) –

  • steps (int) – size of window to average over

  • mode (str (default='forward')) – if ‘forward’ current value at start of window if ‘backward’ current value at end of window if ‘centre’ or ‘center’ current value in centre of window

Return type

array_like (len same as input array)

eqsig.functions.calc_smooth_fa_spectrum(fa_frequencies, fa_spectrum, smooth_fa_frequencies=None, band=40)[source]

Calculates the smoothed Fourier Amplitude Spectrum using the method by Konno and Ohmachi (1998)

Note: different order of inputs than generate_smooth_fa_spectrum

Parameters
  • smooth_fa_frequencies (array_like) – Frequencies to compute the smoothed amplitude

  • fa_frequencies (array_like) – Frequencies of the Fourier amplitude spectrum

  • fa_spectrum (array_like) – Amplitudes of the Fourier amplitude spectrum

  • band – window parameter

Returns

smoothed_fa_spectrum – Amplitudes of smoothed Fourier spectrum at specified frequencies

Return type

array_like

eqsig.functions.calc_smooth_fa_spectrum_w_custom_matrix(asig, smooth_matrix)[source]

Calculates the smoothed Fourier Amplitude Spectrum using a custom filter

eqsig.functions.calc_smoothing_matrix_konno_1998(fa_frequencies, smooth_fa_frequencies=None, band=40)[source]
Calculates the smoothing matrix for computing the smoothed Fourier Amplitude Spectrum

using the method by Konno and Ohmachi 1998

Parameters
  • fa_frequencies (array_like) – Frequencies of FAS

  • smooth_fa_frequencies (array_like) – Frequencies that smooth FAS should be computed at

  • band (int) – Bandwidth of smoothing function

Return type

2d-array_like

eqsig.functions.calc_step_fn_steps_vals(values, ind=None)[source]
eqsig.functions.calc_step_fn_vals_error(values, pow=1, dir=None)[source]

Calculates the error function generated by fitting a step function to the values

Note: Assumes minimum error is at the minimum sum of the error, regardless of the pow. I.e. The best fit is taken as the mean of the values.

Parameters
  • values (array_like) –

  • pow (int) – The power that the error should be raised to

  • dir (str) – Desired direction of the step function if ‘down’, then all upward steps are set to 10x maximum error if ‘up’, then all downward steps are set to 10x maximum error else, no modification to error

Return type

array_like (len same as input array)

eqsig.functions.clean_out_non_changing(values)[source]

Takes an array removes all values that are the same as the previous value.

Parameters

values – array of floats

Returns

cleaned array, indices of clean values in original array

eqsig.functions.determine_indices_of_peaks_for_cleaned(values)[source]

DEPRECATED: Use determine_indices_of_peaks_for_cleaned_array()

eqsig.functions.determine_indices_of_peaks_for_cleaned_array(values)[source]

Determines the position of values that form a local peak in a signal.

Warning: data must be cleaned so that adjacent points have the same value

Parameters

values (array_like) – Array of values that peaks will be found in

Returns

peak_indices – Array of indices of peaks

Return type

array_like of int

eqsig.functions.determine_peak_only_delta_series_4_cleaned_data(values)[source]

Determines the

Note: array must not contain adjacent repeated values

Parameters

values

Returns

eqsig.functions.determine_peaks_only_delta_series(values)[source]

Creates an array with the changes between peak values and zeros for non-peak values.

:param : :type : param values: array_like, array of values :param : :type : return:

Examples

>>> values = np.array([0, 2, 1, 2, 0, 1, 0, -1, 0, 1, 0])
np.array([0, 2, 1, 2, 0.3, 1, 0.3, -1, 0.4, 1, 0])
>>> determine_peaks_only_delta_series(values)
array([0,  2, -1,  1,  0,  -1,  0,  -2,  0,  2,  0])
eqsig.functions.determine_pseudo_cyclic_peak_only_series(values)[source]

Creates an array with only peak values assuming an alternative sign and zeros for non-peak values.

:param : :type : param values: array_like, array of values :param : :type : return:

Examples

>>> values = np.array([0, 2, 1, 2, 0, 1, 0, -1, 0, 1, 0])
np.array([0, 2, 1, 2, 0.3, 1, 0.3, -1, 0.4, 1, 0])
>>> determine_pseudo_cyclic_peak_only_series(values)
array([0,  2, -1,  2,  0,  1,  0,  1,  0,  1,  0])
eqsig.functions.fas2signal(fas, dt, stype='signal')[source]

Convert a fourier spectrum to time series signal

Parameters
  • fas (array_like of img floats) – Positive part only

  • dt (float) – time step of time series

  • stype (str) – If ‘signal’ then return Signal, else return AccSignal

eqsig.functions.fas2values(fas, dt)[source]

Convert a fourier spectrum to time series signal

Parameters
  • fas (array_like of img floats) – Positive part only

  • dt (float) – time step of time series

  • stype (str) – If ‘signal’ then return Signal, else return AccSignal

eqsig.functions.gen_ricker_wavelet_asig(omega, t0, duration, dt)[source]

Generates an acceleration time series that is a Ricker wavelet

Parameters
  • omega

  • t0

  • duration (float) – Total duration of motion

  • dt (float) – Time step of motion

eqsig.functions.generate_fa_spectrum(sig, n_pad=True)[source]

Produces the Fourier amplitude spectrum

Parameters

sig (eqsig.Signal) –

Returns

  • fa_spectrum (complex array_like) – Complex values of the spectrum

  • fa_frequencies (array_like) – Frequencies of the spectrum

eqsig.functions.generate_smooth_fa_spectrum(smooth_fa_frequencies, fa_frequencies, fa_spectrum, band=40)[source]

Deprecated - use calc_smooth_fa_spectrum

eqsig.functions.get_bandwidth_boore_2003(asig)[source]
eqsig.functions.get_major_change_indices(y, rtol=1e-08, atol=1e-05, already_diff=False, dx=1)[source]

Get indices where a significant change in slope occurs

eqsig.functions.get_n_cyc_array(values, opt='all', start='origin')[source]

Given an array, create an array of the same length that numbers the peaks :param values: :param opt:

eqsig.functions.get_peak_array_indices(values, ptype='all')[source]

Find the indices for all of the local maxima and minima

:param : :type : param values: array_like, array of values :param : :type : return:

Examples

>>> values = np.array([0, 2, 1, 2, -1, 1, 1, 0.3, -1, 0.2, 1, 0.2])
np.array([0, 2, 1, 2, -1, 1, 1, 0.3, -1, 0.2, 1, 0.2])
>>> get_peak_array_indices(values)
np.array([0, 1, 2, 3, 4, 5, 8, 10, 11])
eqsig.functions.get_peak_indices(asig)[source]
eqsig.functions.get_section_average(series, start=0, end=-1, index=False)[source]

Gets the average value of a part of series.

Common use is so that it can be patched with another record.

Parameters
  • series – A TimeSeries object

  • start – int or float, optional, Section start point

  • end – int or float, optional, Section end point

  • index – bool, optional, if False then start and end are considered values in time.

:return float,

The mean value of the section.

eqsig.functions.get_sig_array_indexes_range(fas1_smooth, ratio=15)[source]
eqsig.functions.get_sig_freq_range(asig, ratio=15)[source]
eqsig.functions.get_switched_peak_array_indices(values)[source]

Find the indices for largest peak between each zero crossing

Parameters

values (array_like) –

Return type

array_like

eqsig.functions.get_switched_peak_indices(asig)[source]

Find the indices for largest peak between each zero crossing

Parameters

asig (eqsig.AccSignal) –

Return type

array_like

eqsig.functions.get_zero_and_peak_array_indices(pvals, zvals=None, min_step=0)[source]
Parameters
  • pvals

  • zvals

  • min_step – int number of extra steps between zero crossing and peak

Returns

eqsig.functions.get_zero_crossings_array_indices(values, keep_adj_zeros=False)[source]

Find the indices for values that are equal to zero or just passed through zero

Parameters
  • values (array_like) – array of values

  • keep_adj_zeros (bool,) – if false then if adjacent zeros are found, only first is included

:param : :type : return:

Examples

>>> values = np.array([0, 2, 1, 2, -1, 1, 0, 0, 1, 0.3, 0, -1, 0.2, 1, 0.2])
np.array([0, 2, 1, 2, -1, 1, 0, 0, 1, 0.3, 0, -1, 0.2, 1, 0.2])
>>> get_zero_crossings_array_indices(values, keep_adj_zeros=False)
np.array([0, 4, 5, 6, 10, 12])
eqsig.functions.get_zero_crossings_indices(asig)[source]
eqsig.functions.interp2d(x, xf, f)[source]

Can interpolate a table to get an array of values in 2D

Parameters
  • x (array_like) – 1d array of values to be interpolated

  • xf (1d array of values) –

  • f (array_like) – 2d array of function values size=(len(x), n)

Examples

>>> f = np.array([[0, 0, 0],
>>>              [0, 1, 4],
>>>              [2, 6, 2],
>>>              [10, 10, 10]
>>>              ])
>>> xf = np.array([0, 1, 2, 3])
>>> x = np.array([0.5, 1, 2.2, 2.5])
>>> f_interp = interp2d(x, xf, f)
>>> print(f_interp[0][0])
0.0
>>> print(f_interp[0][1])
0.5
>>> print(f_interp[0][2])
2.0
eqsig.functions.interp_array_to_approx_dt(values, dt, target_dt=0.01, even=True)[source]

Interpolate an array of values to a new time step

Similar to interp_to_approx_dt

Only a target time step is provided and the algorithm determines

what time step is best to minimise loss of data from aliasing

Parameters
  • values (array_like) – values of time series

  • dt (float) – Time step

  • target_dt (float) – Target time step

  • even (bool) – If true then forces the number of time steps to be an even number

Returns

  • new_values (array_like) – Interpolated value of time series

  • new_dt (float) – New time step of interpolate time series

eqsig.functions.interp_left(x0, x, y=None)[source]

Interpolation takes the lower value

Parameters
  • x0 (array_like) – Values to be interpolated on x-axis

  • x (array_like) – Existing values on x-axis

  • y (array_like) – Existing y-axis values

eqsig.functions.interp_to_approx_dt(asig, target_dt=0.01, even=True)[source]

Interpolate a signal to a new time step

Only a target time step is provided and the algorithm determines

what time step is best to minimise loss of data from aliasing

Parameters
  • asig (eqsig.AccSignal) – Acceleration time series object

  • target_dt (float) – Target time step

  • even (bool) – If true then forces the number of time steps to be an even number

Returns

new_asig – Acceleration time series object of interpolated time series

Return type

eqsig.AccSignal

eqsig.functions.join_sig_w_time_shift(sig, time_shifts, jtype='add')[source]

Zero pads values of a signal by an array of time shifts and joins it with the original

Parameters
  • sig (eqsig.Signal) – signal to be shifted

  • time_shifts (array_like) – Time shifts to be performed

  • jtype (str (default='add')) – if = ‘add’ then shifted and original signals are added, if =’sub’ then subtracted

Returns

shifted_values

Return type

array_like [2D shape(len(sig.values), len(shift))]

eqsig.functions.join_values_w_shifts(values, shifts, jtype='add')[source]

Zero pads values by an array of shifts and joins it with the original values

Parameters
  • values (array_like) – values to be shifted

  • shifts (array_like [int]) – Shifts to be performed

  • jtype (str (default='add')) – if = ‘add’ then shifted and original values are added, if =’sub’ then subtracted

Returns

shifted_values

Return type

array_like [2D shape(len(values), len(shift))]

eqsig.functions.put_array_in_2d_array(values, shifts, clip='none')[source]

Creates a 2D array where values appear on each line, shifted by a set of indices

Parameters
  • values (array_like (1d)) – Values to be shifted

  • shifts (array_like (int)) – Indices to shift values

  • clip (str or none) – if ‘end’ then returned 2D array trims values that overlap end of input values array

Return type

array_like (2D)

eqsig.functions.remove_poly(values, poly_fit=0)[source]

Calculates best fit polynomial and removes it from the record

eqsig.functions.resample_to_approx_dt(asig, target_dt=0.01, even=True)[source]

Resample a signal assuming periodic to a new time step

Only a target time step is provided and the algorithm determines

what time step is best to minimise loss of data from aliasing

Parameters
  • asig (eqsig.AccSignal) – Acceleration time series object

  • target_dt (float) – Target time step

  • even (bool) – If true then forces the number of time steps to be an even number

eqsig.functions.time_indices(npts, dt, start, end, index)[source]

Determine the new start and end indices of the time series.

Parameters
  • npts – Number of points in original time series

  • dt – Time step of original time series

  • start – int or float, optional, New start point

  • end – int or float, optional, New end point

  • index – bool, optional, if False then start and end are considered values in time.

Returns

tuple, start index, end index

eqsig.functions.time_series_from_motion(motion, dt)[source]

Loading and saving functions

The eqsig format is:
  • First line: name of signal

  • Second line: Number of points and time steps, separated by a space

  • Remaining lines: Each line contains one value of time series

eqsig.loader.load_3_comp_values_and_dt_from_v2a(ffp)[source]

Loads a ground motion file stored in the V2A format

Parameters

ffp

eqsig.loader.load_asig(ffp, load_label=False, m=1.0)[source]

Loads an AccSignal that was saved in eqsig input format.

Parameters
  • ffp (str) – Full file path to output file

  • load_label (bool) – if true then get label from file

  • m (float (default=1.0)) – Scale factor to apply to time series data when loading

Returns

asig

Return type

eqsig.AccSignal

eqsig.loader.load_sig(ffp, m=1.0)[source]

Loads a Signal that was saved in eqsig input format.

Parameters

ffp (str) – Full file path to output file

Returns

  • sig (eqsig.Signal)

  • m (float (default=1.0)) – Scale factor to apply to time series data when loading

eqsig.loader.load_signal(ffp, astype='sig')[source]
eqsig.loader.load_values_and_dt(ffp)[source]

Loads values and time step that were saved in eqsig input format.

Parameters

ffp (str) – Full file path to output file

Returns

  • values (array_like) – An array of values

  • dt (float) – Time step

eqsig.loader.save_signal(ffp, signal)[source]

Saves a Signal or AccSignal to the eqsig format

Parameters
  • ffp (str) – Full file path to output file

  • signal

eqsig.loader.save_values_and_dt(ffp, values, dt, label)[source]

Exports acceleration values to the eqsig format.

Parameters
  • ffp (str) – Full file path to output file

  • values (array_like) – An array of values

  • dt (float) – Time step

  • label (str) – A label of the data

eqsig.sdof module

eqsig.sdof.absmax(a, axis=None)[source]
eqsig.sdof.calc_input_energy_spectrum(acc_signal, periods=None, xi=None, series=False)[source]
eqsig.sdof.calc_resp_uke_spectrum(acc_signal, periods=None, xi=None)[source]

Calculates the sdof response (kinematic + stored) energy spectrum

Parameters
  • acc_signal

  • periods

  • xi

Returns

eqsig.sdof.compute_a_and_b(xi, w, dt)[source]

From the paper by Nigam and Jennings (1968), computes the two matrices.

Parameters
  • xi – critical damping ratio

  • w – angular frequencies

  • dt – time step

Returns

matrices A and B

eqsig.sdof.nigam_and_jennings_response(acc, dt, periods, xi)[source]

Implementation of the response spectrum calculation from Nigam and Jennings (1968).

Ref: Nigam, N. C., Jennings, P. C. (1968) Digital calculation of response spectra from strong-motion earthquake records. National Science Foundation.

Parameters
  • acc – acceleration in m/s2

  • periods – response periods of interest

  • dt – time step of the acceleration time series

  • xi – critical damping factor

Returns

response displacement, response velocity, response acceleration

eqsig.sdof.pseudo_response_spectra(motion, dt, periods, xi)[source]

Computes the maximum response displacement, pseudo velocity and pseudo acceleration.

Parameters
  • motion – array floats, acceleration in m/s2

  • dt – float, the time step

  • periods – array floats, The period of SDOF oscilator

  • xi – float, fraction of critical damping (e.g. 0.05)

Returns

tuple floats, (spectral displacement, pseudo spectral velocity, pseudo spectral acceleration)

eqsig.sdof.response_series(motion, dt, periods, xi)[source]

Computes the elastic response to the acceleration time series

Parameters
  • motion – array floats, acceleration in m/s2

  • dt – float, the time step

  • periods – array floats, The period of SDOF oscillator

  • xi – float, fraction of critical damping (e.g. 0.05)

Returns

tuple of float arrays, (response displacements, response velocities, response accelerations)

eqsig.sdof.single_elastic_response(motion, step, period, xi)[source]

Perform Duhamels integral to get the displacement. http://www.civil.utah.edu/~bartlett/CVEEN7330/Duhamel%27s_integral.pdf http://www1.aucegypt.edu/faculty/mharafa/MENG%20475/Forced%20Vibration.pdf

Parameters
  • motion – acceleration in m/s2

  • step – the time step

  • period – The period of SDOF oscillator

  • xi – fraction of critical damping (e.g. 0.05)

Returns

eqsig.sdof.slow_response_spectra(motion, step, periods, xis)[source]

Perform Duhamels integral to get the displacement. http://www.civil.utah.edu/~bartlett/CVEEN7330/Duhamel%27s_integral.pdf http://www1.aucegypt.edu/faculty/mharafa/MENG%20475/Forced%20Vibration.pdf

Parameters
  • motion – acceleration in m/s2

  • step – the time step

  • period – The period of SDOF oscilator

  • xi – fraction of critical damping (e.g. 0.05)

Returns

eqsig.sdof.time_the_generation_of_response_spectra()[source]
eqsig.sdof.true_response_spectra(motion, dt, periods, xi)[source]

Computes the actual maximum response values, not the pseudo values

Parameters
  • motion – array floats, acceleration in m/s2

  • dt – float, the time step

  • periods – array floats, The period of SDOF oscilator

  • xi – float, fraction of critical damping (e.g. 0.05)

Returns

tuple floats, (spectral displacement, spectral velocity, spectral acceleration)

eqsig.displacements module

eqsig.displacements.calc_velo_and_disp_from_accel_arr(acceleration, dt, trap=True)[source]

Computes the velocity and acceleration of an acceleration time series, using numerical integration.

Parameters
  • acceleration (array_like) – acceleration time series

  • dt (float) – time step

  • trap (bool (default=True)) – if True then uses trapezium integration

Returns

  • velocity (array_like (len=len(acceleration))) – velocity time series

  • displacement (array_like (len=len(acceleration))) – displacement time series

eqsig.displacements.velocity_and_displacement_from_acceleration(acceleration, dt, trap=True)[source]

Computes the velocity and acceleration of an acceleration time series, using numerical integration.

Parameters
  • acceleration (array_like) – acceleration time series

  • dt (float) – time step

  • trap (bool (default=True)) – if True then uses trapezium integration

Returns

  • velocity (array_like (len=len(acceleration))) – velocity time series

  • displacement (array_like (len=len(acceleration))) – displacement time series

eqsig.stockwell module

eqsig.stockwell.dep_itransform(stock)[source]

Performs an inverse Stockwell Transform

eqsig.stockwell.generate_gaussian(n_d2)[source]

create a gaussian distribution

eqsig.stockwell.get_max_stockwell_freq(asig)[source]
eqsig.stockwell.get_max_tifq_vals_freq(tifq_values, dt)[source]
eqsig.stockwell.get_stockwell_freqs(asig)[source]
eqsig.stockwell.get_stockwell_times(asig)[source]
eqsig.stockwell.itransform(stock)[source]

Performs an inverse Stockwell Transform

eqsig.stockwell.plot_fas_at_time(splot, asig, time)[source]

Plots the Fourier amplitude spectrum at a time based on a Stockwell transform

eqsig.stockwell.plot_max_freq_azimuth(splot, asig1, asig2, max_f=None, norm=False, r_steps=90)[source]
eqsig.stockwell.plot_stock(splot, asig, norm_x=False, norm_all=False, interp=False, cmap=None, vmin=None, vmax=None)[source]

Plots the Stockwell transform of an acceleration signal

Parameters
  • splot (matplotlib.pyplot.subplot) –

  • asig (eqsig.Signal) –

  • norm_x (bool, default=False) – If true then all values at a time step are normalised by the maximum value

  • norm_all (bool, default=False) – If true the all values are normalised by the maximum value

Return type

None

eqsig.stockwell.plot_tifq_vals(subplot, tifq_vals, dt, norm_all=False, norm_x=False, cmap=None)[source]

Plots a time frequency spectrum (e.g. Stockwell transform)

Same as plot_stock but takes in values (tifq_vals) and time step (dt) instead of AccSignal object.

Parameters
  • subplot (matplotlib.pyplot.subplot) –

  • tifq_vals (2-d array) – time-frequency values

  • dt (float) – time step

  • norm_x (bool, default=False) – If true then all values at a time step are normalised by the maximum value

  • norm_all (bool, default=False) – If true the all values are normalised by the maximum value

Return type

None

eqsig.stockwell.plot_windowed_fas_at_time(splot, asig, time, time_window=3)[source]

Plots the time averaged Fourier amplitude spectrum at a time based on a Stockwell transform

eqsig.stockwell.transform(acc, interp=False)[source]

Performs a Stockwell transform on an array

Assumes m = 1, p = 1

Parameters

acc – array_like

Returns

eqsig.stockwell.transform_slow(acc, interp=False, ith=0)[source]

Performs a Stockwell transform on an array

Assumes m = 1, p = 1

Parameters

acc – array_like

Returns

eqsig.stockwell.transform_w_scipy_fft(acc, interp=False)[source]

Performs a Stockwell transform on an array

Assumes m = 1, p = 1

Parameters

acc – array_like

Returns

eqsig.stockwell module

eqsig.stockwell.dep_itransform(stock)[source]

Performs an inverse Stockwell Transform

eqsig.stockwell.generate_gaussian(n_d2)[source]

create a gaussian distribution

eqsig.stockwell.get_max_stockwell_freq(asig)[source]
eqsig.stockwell.get_max_tifq_vals_freq(tifq_values, dt)[source]
eqsig.stockwell.get_stockwell_freqs(asig)[source]
eqsig.stockwell.get_stockwell_times(asig)[source]
eqsig.stockwell.itransform(stock)[source]

Performs an inverse Stockwell Transform

eqsig.stockwell.plot_fas_at_time(splot, asig, time)[source]

Plots the Fourier amplitude spectrum at a time based on a Stockwell transform

eqsig.stockwell.plot_max_freq_azimuth(splot, asig1, asig2, max_f=None, norm=False, r_steps=90)[source]
eqsig.stockwell.plot_stock(splot, asig, norm_x=False, norm_all=False, interp=False, cmap=None, vmin=None, vmax=None)[source]

Plots the Stockwell transform of an acceleration signal

Parameters
  • splot (matplotlib.pyplot.subplot) –

  • asig (eqsig.Signal) –

  • norm_x (bool, default=False) – If true then all values at a time step are normalised by the maximum value

  • norm_all (bool, default=False) – If true the all values are normalised by the maximum value

Return type

None

eqsig.stockwell.plot_tifq_vals(subplot, tifq_vals, dt, norm_all=False, norm_x=False, cmap=None)[source]

Plots a time frequency spectrum (e.g. Stockwell transform)

Same as plot_stock but takes in values (tifq_vals) and time step (dt) instead of AccSignal object.

Parameters
  • subplot (matplotlib.pyplot.subplot) –

  • tifq_vals (2-d array) – time-frequency values

  • dt (float) – time step

  • norm_x (bool, default=False) – If true then all values at a time step are normalised by the maximum value

  • norm_all (bool, default=False) – If true the all values are normalised by the maximum value

Return type

None

eqsig.stockwell.plot_windowed_fas_at_time(splot, asig, time, time_window=3)[source]

Plots the time averaged Fourier amplitude spectrum at a time based on a Stockwell transform

eqsig.stockwell.transform(acc, interp=False)[source]

Performs a Stockwell transform on an array

Assumes m = 1, p = 1

Parameters

acc – array_like

Returns

eqsig.stockwell.transform_slow(acc, interp=False, ith=0)[source]

Performs a Stockwell transform on an array

Assumes m = 1, p = 1

Parameters

acc – array_like

Returns

eqsig.stockwell.transform_w_scipy_fft(acc, interp=False)[source]

Performs a Stockwell transform on an array

Assumes m = 1, p = 1

Parameters

acc – array_like

Returns

eqsig.exceptions module

exception eqsig.exceptions.SignalProcessingError[source]

Bases: Exception

exception eqsig.exceptions.SignalProcessingWarning[source]

Bases: Warning

eqsig.exceptions.deprecation(message)[source]