Code documentation

Definitions of signals

pyExSi.signals.burst_random(N, A=1.0, ratio=0.5, distribution='uniform', n_bursts=1, periodic_bursts=True, rg=None)[source]

Generate a zero-mean burst random excitation signal time series.

Parameters:
  • N – Number of time points.

  • A – Amplitude of the random signal. For ‘uniform’ distribution, this is the peak-to-peak amplitude, for ‘normal’ distribution this is the RMS.

  • ratio – The ratio of burst legth ot the total legth of the time series.

  • distribution – ‘uniform’, ‘normal’ or ‘pseudorandom’. Defaults to ‘uniform’.

  • n_bursts – Number of burst repetition. The output time series will have N*n_bursts points. Defaults to 1.

  • periodic_bursts – If True, bursts are periodically repeated n_bursts times, otherwise a uniquely random burst is generated for each repetition. Defaults to True.

  • rg (numpy.random._generator.Generator) – Initialized Generator object

Returns:

Burst random signal time series.

Example

>>> import numpy as np
>>> import matplotlib.pyplot as plt
>>> import pyExSi as es
>>> N = 1000
>>> amplitude = 5
>>> x = es.burst_random(N, A=amplitude, ratio=0.1, distribution='normal', n_bursts=3)
>>> plt.plot(x)
>>> plt.show()
pyExSi.signals.get_kurtosis(signal)[source]

Kurtosis of signal.

Parameters:

signal (array) – input signal.

Returns:

kurtosis

pyExSi.signals.get_psd(freq, freq_lower, freq_upper, variance=1)[source]

One-sided flat-shaped power spectral density (PSD).

Parameters:
  • freq (array) – Frequency vector [Hz]

  • freq_lower (float) – Lower frequency of PSD [Hz]

  • freq_upper (float) – Upper frequency of PSD [Hz]

  • variance (float) – Variance of random process, described by PSD [unit^2]

Returns:

one-sided flat-shaped PSD [unit^2/Hz]

Example

>>> import numpy as np
>>> import matplotlib.pyplot as plt
>>> import pyExSi as es
>>> N = 1000 # number of data points of time signal
>>> fs = 100 # sampling frequency [Hz]
>>> t = np.arange(0,N)/fs # time vector
>>> M = N // 2 + 1 # number of data points of frequency vector
>>> freq = np.arange(0, M, 1) * fs / N # frequency vector
>>> freq_lower = 10 # PSD lower frequency limit  [Hz]
>>> freq_upper = 20 # PSD upper frequency limit [Hz]
>>> PSD = es.get_psd(freq, freq_lower, freq_upper) # one-sided flat-shaped PSD
>>> plt.plot(freq,PSD)
>>> plt.xlabel(f [Hz])
>>> plt.ylabel(PSD [unit^2/Hz])
>>> plt.show()
pyExSi.signals.impulse(N, n_start=0, width=None, amplitude=1.0, window='sine')[source]

Impact impulse of the shape defined with the parameter window.

Parameters:
  • N (int) – Number of points in time signal.

  • width (int) – Number of points for pulse width, None results in width=N

  • amplitude (float) – Amplitude of pulse.

  • window (string, float, or tuple) – The type of window to create. See scipy.signal.windows for more details.

Returns:

impact pulse.

Example

>>> import numpy as np
>>> import matplotlib.pyplot as plt
>>> import pyExSi as es
>>> N = 1000
>>> n_start = 100
>>> width = 200
>>> amplitude = 3
>>> x_1 = es.impulse(N=N, n_start=n_start, width=width, amplitude=amplitude, window='triang')
>>> x_2 = es.impulse(N=N, n_start=n_start, width=width, amplitude=amplitude, window=('exponential',0,10))
>>> t = np.linspace(0,10,N)
>>> plt.plot(t,x_1, label='tringular')
>>> plt.plot(t,x_2, label='exponential')
>>> plt.legend()
>>> plt.show()
pyExSi.signals.nonstationary_signal(N, PSD, fs, k_u=3, modulating_signal=('PSD', None), param1_list=None, param2_list=None, seed=None, SQ=False)[source]

Non-stationary non-Gaussian realization of random process.

Non-stationarity random process is obtained by amplitude modulation of Gaussian random process[1]. Gaussian random process is obtained with IFFT of amplitude spectra with random phase [2]. Tuple modulating_signal selects the type of modulating signal: ‘PSD’ for random process realization [3], where PSD_modulating is power spectrum density of modulating signal, and ‘CSI’ for cubic spline interpolation [4,5], with sample step delta_n. The desired kurtosis k_u is obtained by iteration over lists param1_list and param2_list (for ‘PSD’ p and delta_m are needed, for ‘CSI’ alpha and beta are needed).

Parameters:
  • N ({int, float}) – Number of data points in returned signal

  • PSD (array) – One-sided power spectral density of carrier signal

  • fs ({int, float}) – sampling frequency

  • k_u (float) – Desired kurtosis value of returned signal. Defaults to 3 (Gaussian random process).

  • modulating_signal (tuple with name and parameter.) – Delects type of modulating signal and provides needed parameter.

  • param1_list (list of floats) – List of first parameter for modulating signal generation. Contains parameters p or alpha

  • param2_list (list of floats) – List of second parameter for modulating signal generation. Contains parameters delta_m or beta

  • seed ({None, int, array_like[ints], SeedSequence, BitGenerator, Generator}, optional) – A seed to initialize the BitGenerator. For details, see numpy.random.default_rng()

  • SQ (boolean) – If squeezing of signal [4] is required, set ‘True’. Defaults to ‘False’

Returns:

nonstationary signal. Optionally, stationary and modulating_signal are returned as well.

References

[1] Frederic Kihm, Stephen A. Rizzi, N. S. Ferguson, and Andrew Halfpenny. Understanding how kurtosis is transferred from input acceleration to stress response and it’s influence on fatigue life. In Proceedings of the XI International Conference on Recent Advances in Structural Dynamics, Pisa, Italy, 07 2013.

[2] D. E. Newland. An Introduction to Random Vibrations, Spectral & Wavelet Analysis. Dover Publications, 2005

[3] Arvid Trapp, Mafake James Makua, and Peter Wolfsteiner. Fatigue assessment of amplitude-modulated nonstationary random vibration loading. Procedia Structural Integrity, 17:379—-386, 2019.

[4] Lorenzo Capponi, Martin Česnik, Janko Slavič, Filippo Cianetti, and Miha Boltežar. Non-stationarity index in vibration fatigue: Theoretical and ex-perimental research.International Journal of Fatigue, 104:221–230, 2017.

[5] Janko Slavič, Matjaž Mršnik, Martin Česnik, Jaka Javh, Miha Boltežar. Vibration Fatigue by Spectral Methods, From Structural Dynamics to Fatigue Damage – Theory and Experiments, ISBN: 9780128221907, Elsevier, 1st September 2020

Example

>>> import numpy as np
>>> import matplotlib.pyplot as plt
>>> import pyExSi as es
>>> N = 1000 # number of data points of time signal
>>> fs = 100 # sampling frequency [Hz]
>>> t = np.arange(0,N)/fs # time vector
>>> M = N // 2 + 1 # number of data points of frequency vector
>>> freq = np.arange(0, M, 1) * fs / N # frequency vector
>>> freq_lower = 10 # PSD lower frequency limit  [Hz]
>>> freq_upper = 20 # PSD upper frequency limit [Hz]
>>> freq_lower_mod = 1 # modulating signals's PSD lower frequency limit  [Hz]
>>> freq_upper_mod = 2 # modulating signals's PSD upper frequency limit [Hz]

PSD of stationary and modulating signal

>>> PSD = es.get_psd(freq, freq_lower, freq_upper) # one-sided flat-shaped PSD
>>> PSD_modulating = es.get_psd(freq, freq_lower_mod, freq_upper_mod) # one-sided flat-shaped PSD

Specify kurtosis and return non-stationary signal

>>> k_u = 5
>>> x_nonstationary_1 = es.nonstationary_signal(N,PSD,fs,k_u=k_u,modulating_signal=('PSD',PSD_modulating))

Calculate kurtosis

>>> k_u_1 = es.get_kurtosis(x_nonstationary_1)
>>> print(f'desired kurtosis :{k_u:.3f}', actual kurtosis :{k_u_1:.3f}')

Refined array with amplitude modulation parameters

>>> delta_m_list = np.arange(.1,2.1,.1)
>>> p_list = np.arange(.1,2.1,.1)
>>> x_nonstationary_2 = es.nonstationary_signal(N,PSD,fs,k_u=k_u,modulating_signal=('PSD',PSD_modulating),
                                                param1_list=delta_m_list,param2_list=p_list)
>>> k_u_2 = es.get_kurtosis(x_nonstationary_2)
>>> print(f'desired kurtosis :{k_u:.3f}', actual kurtosis :{k_u_2:.3f}')

Define array of parameters alpha and beta

>>> alpha_list = np.arange(1,4,.5)
>>> beta_list = np.arange(1,4,.5)
>>> x_nonstationary_3 = es.nonstationary_signal(N,PSD,fs,k_u=10,modulating_signal=('CSI',delta_n),
                                                    param1_list=alpha_list,param2_list=beta_list)
>>> k_u_3 = es.get_kurtosis(x_nonstationary_3)
>>> print(f'desired kurtosis :{k_u:.3f}', actual kurtosis :{k_u_3:.3f}')
>>> plt.plot(t, x_nonstationary_2, label='PSD')
>>> plt.plot(t, x_nonstationary_3, label='CSI)
>>> plt.xlabel(t [s])
>>> plt.ylabel(x [unit])
>>> plt.legend()
>>> plt.show()
pyExSi.signals.normal_random(N, rg=None)[source]

Normal random distribution.

Parameters:
  • N (int) – Number of points.

  • rg (numpy.random._generator.Generator) – Initialized Generator object

Returns:

Random samples from the “standard normal” distribution

Example

>>> import numpy as np
>>> import matplotlib.pyplot as plt
>>> import pyExSi as es
>>> N = 100
>>> x = es.uniform_random(N=N)
>>> plt.plot(x)
>>> plt.show()
pyExSi.signals.pseudo_random(N, rg=None)[source]

Pseudorandom distribution.

Magnitudes are 1, phase is random.

Parameters:
  • N (int) – Number of points.

  • rg (numpy.random._generator.Generator) – Initialized Generator object

Returns:

Random samples from the “standard normal” distribution

Example

>>> import numpy as np
>>> import matplotlib.pyplot as plt
>>> import pyExSi as es
>>> N = 100
>>> x = es.pseudo_random(N=N)
>>> plt.plot(x)
>>> plt.show()
pyExSi.signals.random_gaussian(N, PSD, fs, rg=None)[source]

Stationary Gaussian realization of random process, characterized by PSD.

Random process is obtained with IFFT of amplitude spectra with random phase [1]. Area under PSD curve represents variance of random process.

Parameters:
  • N (int) – Number of points.

  • PSD (array) – one-sided power spectral density [unit^2].

  • fs (int,float) – sampling frequency [Hz].

  • rg (numpy.random._generator.Generator) – Initialized Generator object

Returns:

stationary Gaussian realization of random process

References

[1] D. E. Newland. An Introduction to Random Vibrations, Spectral & Wavelet Analysis. Dover Publications, 2005

Example

>>> import numpy as np
>>> import matplotlib.pyplot as plt
>>> import pyExSi as es
>>> N = 1000 # number of data points of time signal
>>> fs = 100 # sampling frequency [Hz]
>>> t = np.arange(0,N)/fs # time vector
>>> M = N // 2 + 1 # number of data points in frequency vector
>>> freq = np.arange(0, M, 1) * fs / N # frequency vector
>>> freq_lower = 10 # PSD lower frequency limit  [Hz]
>>> freq_upper = 20 # PSD upper frequency limit [Hz]
>>> PSD = es.get_psd(freq, freq_lower, freq_upper) # one-sided flat-shaped PSD
>>> x = es.random_gaussian(N, PSD, fs)
>>> plt.plot(t,x)
>>> plt.xlabel(t [s])
>>> plt.ylabel(x [unit])
>>> plt.show()
pyExSi.signals.sine_sweep(time, phi=0, freq_start=1, sweep_rate=None, freq_stop=None, mode='linear', phi_end=False, freq_return=False)[source]

Generate a sine sweep signal time series.

Parameters:
  • time – array of shape (N,), time vector.

  • phi – float, initial phase of the sine signal in radians. Defaults to 0.

  • freq_start – float, initial frequency in Hz.

  • sweep_rate – float, the rate of sweep. In Hz/s for a linear sweep, in octaves/minute for a logarithmic sweep. If not given it is calculated from time, freq_start and freq_stop.

  • freq_stop – float, final frequency in Hz.

  • mode – ‘linear’ or ‘logarithmic’, type of sweep, optional. Defaults to ‘linear’.

  • phi_end – If True, return (sweep_sine, phi_end), where phi_end is the end phase which can be used as phi if this function is called for another sweep. Defaults to False.

  • freq_return – If True, return (sweep_sine, phi_end, freq), where freq is the frequency array. Defaults to False.

Returns:

array of shape (N,), the generated sine sweep signal

Example

>>> import numpy as np
>>> import matplotlib.pyplot as plt
>>> import pyExSi as es
>>> t = np.linspace(0,10,1000)
>>> x = es.sine_sweep(time=t, freq_start=0, freq_stop=5)
>>> plt.plot(t, x)
>>> plt.show()
pyExSi.signals.stationary_nongaussian_signal(N, PSD, fs, s_k=0, k_u=3, mean=0, rg=None)[source]

Stationary non-Gaussian realization of random process.

Random process is obtained with IFFT of amplitude spectra with random phase [1]. Non-Gaussianity is obtained by Winterstein polynomials [2].

Parameters:
  • N (int) – number of data points in returned signal

  • PSD (array) – one-sided power spectral density

  • fs (int, float) – sampling frequency

  • s_k (int, float) – skewness of returned signal

  • k_u (int, float) – kurtossis of returned signal

  • mean (int, float) – mean value of returned signal

  • rg (numpy.random._generator.Generator) – Initialized Generator object

Returns:

stationary non-Gaussian realization of random process.

References

[1] D. E. Newland. An Introduction to Random Vibrations, Spectral & Wavelet Analysis. Dover Publications, 2005

[2] Steven R. Winterstein. Nonlinear vibration models for extremes and fatigue. ASCE Journal of Engineering Mechanics, 114:1772–1790, 1988.

Example

>>> import numpy as np
>>> import matplotlib.pyplot as plt
>>> import pyExSi as es
>>> N = 1000 # number of data points of time signal
>>> fs = 100 # sampling frequency [Hz]
>>> t = np.arange(0,N)/fs # time vector
>>> M = N // 2 + 1 # number of data points of frequency vector
>>> freq = np.arange(0, M, 1) * fs / N # frequency vector
>>> freq_lower = 10 # PSD lower frequency limit  [Hz]
>>> freq_upper = 20 # PSD upper frequency limit [Hz]
>>> PSD = es.get_psd(freq, freq_lower, freq_upper) # one-sided flat-shaped PSD
>>> x_gauss = es.random_gaussian(N, PSD, fs)
>>> x_ngauss = es.stationary_nongaussian_signal(N, PSD, fs, k_u = 5)
>>> plt.plot(t, x_gauss, label='gaussian')
>>> plt.plot(t, x_ngauss, label='non-gaussian')
>>> plt.xlabel(t [s])
>>> plt.ylabel(x [unit])
>>> plt.legend()
>>> plt.show()
pyExSi.signals.uniform_random(N, rg=None)[source]

Uniform random distribution

Parameters:
  • N (int) – Number of points.

  • rg (numpy.random._generator.Generator) – Initialized Generator object

Returns:

Random samples from a “uniform” distribution

Example

>>> import numpy as np
>>> import matplotlib.pyplot as plt
>>> import pyExSi as es
>>> N = 100
>>> x = es.uniform_random(N=N)
>>> plt.plot(x)
>>> plt.show()