Defined as:
sdynpy.core.sdynpy_data.TimeHistoryArrayModule:
sdynpy.core.sdynpy_dataSource: GitHub
Parent:
sdynpy.NDDataArrayParent:
sdynpy.SdynpyArrayParent:
numpy.ndarray
Signature¶
class sdynpy.TimeHistoryArray(shape, nelements, buffer=None, offset=0, strides=None, order=None)Data array used to store time history data
Attributes¶
| Name | Summary |
|---|---|
function_type | Returns the function type of the data array |
function_type¶
Returns the function type of the data array
Methods¶
| Name | Summary |
|---|---|
apply_transformation | Applies a transformations to the time traces. |
burst_random_signal | Generates a burst random signal with the specified parameters |
chirp_signal | Creates a chirp (sine sweep) signal with the specified parameters |
cpsd | Computes a CPSD matrix from the time histories |
digital_tracking_filter | |
fft | Computes the frequency spectra of the time signal |
filter | Filter the signal using a butterworth filter of the specified order |
find_signal_shift | Computes the shift between two sets of time signals |
find_zero_crossings | Finds zero crossings in the time history |
from_exodus | Reads time data from displacements in an Exodus file |
haversine_signal | Creates a haversine pulse with the specified parameters |
mimo_forward | Performs the forward mimo calculation via convolution. |
mimo_inverse | Performs the inverse source estimation for time domain (transient) problems using Fourier deconvolution. The response nodes used in the inverse source estimation are the ones contained in the supplied FRF matrix. |
overlap_and_add | Creates a time history by overlapping and adding other time histories. |
pseudorandom_signal | Generates a pseudorandom signal at the specified coordinates |
pulse_signal | Creates a pulse using a cosine function raised to a specified exponent |
random_signal | Generates a random signal with the specified parameters |
remove_rigid_body_motion | Removes rigid body displacements from time data. |
resample | Uses Scipy.signal.resample to resample the time history array |
rms | |
shift_signal | Shift a signal in time by a specified amount. |
sine_signal | Creates a sinusoidal signal with the specified parameters |
sine_sweep_signal | |
split_into_frames | Splits a time history into measurement frames with a given overlap and window function applied. |
srs | Compute a shock response spectrum (SRS) from the time history |
stft | Computes a Short-Time Fourier Transform (STFT) |
to_rattlesnake_specification | |
upsample | Upsamples a time history using frequency domain zero padding. |
vold_kalman_filter | |
zefft | Creates a FFTs from time histories with early times zeroed out |
zero_early_time | Creates a time histories with early times zeroed out |
apply_transformation¶
Source: GitHub
def sdynpy.TimeHistoryArray.apply_transformation(self, transformation, invert_transformation=False)Applies a transformations to the time traces.
Parameters¶
transformation : Matrix The transformation to apply to the time traces. It should be a SDynPy matrix object with the “transformed” coordinates on the rows and the “physical” coordinates on the columns. The matrix can only be be 2D.
invert_reference_transformation : bool, optional Whether or not to invert the transformation when applying it to the time traces. The default is False, which is standard practice. The row/column ordering in the transformation should be flipped if this is set to true.
Returns¶
transformed_data : TimeHistoryArray The time traces with the transformations applied.
Raises¶
ValueError
If the transformation array has more than two dimensions.
ValueError
If the physical degrees of freedom in the transformation does not
match the spectra.
burst_random_signal¶
Source: GitHub
def sdynpy.TimeHistoryArray.burst_random_signal(cls, dt, signal_length, coordinates, min_frequency=None, max_frequency=None, signal_rms=1, frames=1, frequency_shape=None, on_fraction=0.5, delay_fraction=0.0, ramp_fraction=0.05, comment1='', comment2='', comment3='', comment4='', comment5='')Generates a burst random signal with the specified parameters
Parameters¶
dt : float Abscissa spacing in the final signal.
signal_length : int Number of samples in the signal
coordinates : CoordinateArray Coordinate array used to generate the signal. If the last dimension of coordinates is not shape 1, then a new axis will be added to make it shape 1. The shape of the resulting TimeHistoryArray will be determined by the shape of the input coordinates.
min_frequency : float, optional Minimum frequency content in the signal. The default is the lowest nonzero frequency line.
max_frequency : float, optional Maximum frequency content in the signal. The default is the highest frequency content in the signal, e.g. the Nyquist frequency.
signal_rms : float or np.ndarray, optional RMS value for the generated signals. The default is 1. The shape of this value should be broadcastable with the size of the generated TimeHistoryArray if different RMS values are desired for each signal. Note that the RMS will be computed for the “burst” part of the signal and not include the zero portion of the signal.
frames : int, optional Number of frames to generate. These will essentially be repeats of the first frame for the number of frames specified. The default is 1.
frequency_shape : function, optional An optional function that should accept a frequency value and return an amplitude at that frequency. The default is constant scaling across all frequency lines.
on_fraction : float, optional The fraction of the frame that the signal is active, default is 0.5. This portion includes the ramp_fraction, so an on_fraction of 0.5 with a ramp_fraction of 0.05 will be at full level for 0.5-2*0.05 = 0.4 fraction of the full measurement frame.
delay_fraction : float, optional The fraction of the frame that is empty before the signal starts, default is 0.0
ramp_fraction : float, optional The fraction of the frame that is used to ramp between the off and active signal
comment1 : np.ndarray, optional Comment used to describe the data in the data array. The default is ‘’.
comment2 : np.ndarray, optional Comment used to describe the data in the data array. The default is ‘’.
comment3 : np.ndarray, optional Comment used to describe the data in the data array. The default is ‘’.
comment4 : np.ndarray, optional Comment used to describe the data in the data array. The default is ‘’.
comment5 : np.ndarray, optional Comment used to describe the data in the data array. The default is ‘’.
Returns¶
TimeHistoryArray A time history containing the specified burst random signal
chirp_signal¶
Source: GitHub
def sdynpy.TimeHistoryArray.chirp_signal(cls, dt, signal_length, coordinates, start_frequency=None, end_frequency=None, frames=1, amplitude_function=None, force_integer_cycles=True, comment1='', comment2='', comment3='', comment4='', comment5='')Creates a chirp (sine sweep) signal with the specified parameters
Parameters¶
dt : float Abscissa spacing in the final signal.
signal_length : int Number of samples in the signal
coordinates : CoordinateArray Coordinate array used to generate the signal. If the last dimension of coordinates is not shape 1, then a new axis will be added to make it shape 1. The shape of the resulting TimeHistoryArray will be determined by the shape of the input coordinates.
start_frequency : TYPE, optional Starting frequency content in the signal. The default is the lowest nonzero frequency line.
end_frequency : TYPE, optional Stopping frequency content in the signal. The default is the highest non-nyquist frequency line.
frames : int, optional Number of frames to generate. These will essentially be repeats of the first frame for the number of frames specified. The default is 1.
amplitude_function : function, optional An optional function that should accept a frequency value and return an amplitude at that frequency. The default is constant scaling across all frequencies. Multiple amplitudes can be returned as long as they broadcast with the shape of the final TimeHistoryArray.
force_integer_cycles : bool, optional If True, it will force an integer number of cycles, which will adjust the maximum frequency of the signal. This will ensure the signal is continuous if repeated. If False, the
comment1 : np.ndarray, optional Comment used to describe the data in the data array. The default is ‘’.
comment2 : np.ndarray, optional Comment used to describe the data in the data array. The default is ‘’.
comment3 : np.ndarray, optional Comment used to describe the data in the data array. The default is ‘’.
comment4 : np.ndarray, optional Comment used to describe the data in the data array. The default is ‘’.
comment5 : np.ndarray, optional Comment used to describe the data in the data array. The default is ‘’.
Returns¶
TimeHistoryArray A time history containing the specified chirp signal
cpsd¶
Source: GitHub
def sdynpy.TimeHistoryArray.cpsd(self, samples_per_frame: int, overlap: float, window: str, averages_to_keep: int = None, only_asds=False, rtol=1, atol=1e-08)Computes a CPSD matrix from the time histories
Parameters¶
samples_per_frame : int Number of samples per frame
overlap : float Overlap fraction (not percent, e.g. 0.5 not 50)
window : str Name of a window function in scipy.signal.windows
averages_to_keep : int, optional Optional number of averages to use, otherwise as many as possible will be used.
only_asds : bool, optional If True, only compute autospectral densities (diagonal of the CPSD matrix)
rtol : float, optional Tolerance used to check abscissa spacing. The default is 1.
atol : float, optional Tolerance used to check abscissa spacing. The default is 1e-8.
Returns¶
cpsd_array : PowerSpectralDensityArray Cross Power Spectral Density Array.
Raises¶
ValueError
If time history abscissa are not equally spaced.
digital_tracking_filter¶
Source: GitHub
def sdynpy.TimeHistoryArray.digital_tracking_filter(self, frequencies, arguments, cutoff_frequency_ratio=0.15, filter_order=2, phase_estimate=None, amplitude_estimate=None, block_size=None, plot_results=False)fft¶
Source: GitHub
def sdynpy.TimeHistoryArray.fft(self, samples_per_frame=None, norm='backward', rtol=1, atol=1e-08, **scipy_rfft_kwargs)Computes the frequency spectra of the time signal
Parameters¶
samples_per_frame : int, optional Number of samples per measurement frame. If this is specified, then the signal will be split up into frames and averaged together. Be aware that if the time signal is not periodic, averaging it may have the effect of zeroing out the spectrum (because the average time signal is zero). The default is no averaging, the frame size is the length of the signal.
norm : str, optional The type of normalization applied to the fft computation. The default is “backward”.
rtol : float, optional Relative tolerance used in the abcsissa spacing check. The default is 1e-5.
atol : float, optional Relative tolerance used in the abscissa spacing check. The default is 1e-8.
scipy_rfft_kwargs Additional keywords that will be passed to SciPy’s rfft function.
Returns¶
SpectrumArray
The frequency spectra of the TimeHistoryArray.
Raises¶
ValueError
Raised if the time signal passed to this function does not have
equally spaced abscissa.
filter¶
Source: GitHub
def sdynpy.TimeHistoryArray.filter(self, filter_order, frequency, filter_type='low', filter_method='filtfilt', filter_kwargs=None)Filter the signal using a butterworth filter of the specified order
Parameters¶
filter_order : int Order of the butterworth filter
frequency : array_like The critical frequency or frequencies. For lowpass and highpass filters, this is a scalar; for bandpass and bandstop filters, this is a length-2 sequence (low,high).
filter_type : str, optional Type of filter. Must be one of ‘low’,‘high’,‘bandpass’, or ‘bandstop’. The default is ‘low’.
filter_method : str, optional Method of filter application. Must be one of ‘filtfilt’ or ‘lfilter’. The default is ‘filtfilt’.
filter_kwargs : dict, optional Optional keyword arguments to be passed to filtfilt or lfilter.
Returns¶
TimeHistoryArray
The filtered time history array.
find_signal_shift¶
Source: GitHub
def sdynpy.TimeHistoryArray.find_signal_shift(self, other_signal, compute_subsample_shift=True, good_line_threshold=0.01)Computes the shift between two sets of time signals
This is the amount that other_signal leads self. If the time shift
is positive, it means that features in other_signal occur earlier in
time compared to self. If the time shift is negative, it means that
features in other_signal occur later in time compared to self.
To align two signals, you can take the time shift from this function and
pass it into the shift_signal method of other_signal.
Parameters¶
other_signal : TimeHistoryArray The signal against which this signal should be compared in time. It should have the same coordinate ordering and the same number of abscissa as this signal.
compute_subsample_shift : bool, optional If False, this function will simply align to the nearest sample. If True, this function will attempt to use FFT phases to compute a subsample shift between the signals. Default is True.
good_line_threshold : float, optional Threshold to use to compute “good” frequency lines. This function uses phase to compute subsample shifts. If there are frequency lines without content, they should be ignored. Frequency lines less than
good_line_thresholdtimes the maximum of the spectra are ignored. The default is 0.01.
Returns¶
time_shift : float The time difference between the two signals.
find_zero_crossings¶
Source: GitHub
def sdynpy.TimeHistoryArray.find_zero_crossings(self, return_abscissa=False, include_start=False)Finds zero crossings in the time history
Parameters¶
return_abscissa : bool, optional If True, returns abscissa values associated with the zero crossing. By default it is False.
include_start : bool, optional If True, returns the first index as a zero crossing, default is False
Returns¶
sign_changes_array : np.ndarray An array of indices into the ordinates of the function where sine changes occur.
abscissa_array : np.ndarray An array of the abscissa values of the function where sine changes occur. Only returned if return_abscissa is True.
from_exodus¶
Source: GitHub
def sdynpy.TimeHistoryArray.from_exodus(cls, exo, x_disp='DispX', y_disp='DispY', z_disp='DispZ', x_rot=None, y_rot=None, z_rot=None, timesteps=None)Reads time data from displacements in an Exodus file
Parameters¶
exo : Exodus or ExodusInMemory The exodus data from which geometry will be created.
x_disp : str, optional String denoting the nodal variable in the exodus file from which the X-direction displacement should be read. The default is ‘DispX’.
y_disp : str, optional String denoting the nodal variable in the exodus file from which the Y-direction displacement should be read. The default is ‘DispY’.
z_disp : str, optional String denoting the nodal variable in the exodus file from which the Z-direction displacement should be read. The default is ‘DispZ’.
timesteps : iterable, optional A list of timesteps from which data should be read. The default is None, which reads all timesteps.
Returns¶
TYPE
DESCRIPTION.
haversine_signal¶
Source: GitHub
def sdynpy.TimeHistoryArray.haversine_signal(cls, dt, signal_length, coordinates, pulse_width=None, pulse_time=None, pulse_peak=1, frames=1, comment1='', comment2='', comment3='', comment4='', comment5='')Creates a haversine pulse with the specified parameters
Parameters¶
dt : float Abscissa spacing in the final signal.
signal_length : int Number of samples in the signal
coordinates : CoordinateArray, optional Coordinate array used to generate the signal. If the last dimension of coordinates is not shape 1, then a new axis will be added to make it shape 1. The shape of the resulting TimeHistoryArray will be determined by the shape of the input coordinates.
pulse_width : float, optional With of the pulse in the same units as
dt. The default is 5*dt.pulse_time : float, optional The time of the pulse’s peak occurance in the same units as
dt. The default is 5*dt.pulse_peak : float, optional The peak amplitude of the pulse. The default is 1.
frames : int, optional Number of frames to generate. These will essentially be repeats of the first frame for the number of frames specified. The default is 1.
comment1 : np.ndarray, optional Comment used to describe the data in the data array. The default is ‘’.
comment2 : np.ndarray, optional Comment used to describe the data in the data array. The default is ‘’.
comment3 : np.ndarray, optional Comment used to describe the data in the data array. The default is ‘’.
comment4 : np.ndarray, optional Comment used to describe the data in the data array. The default is ‘’.
comment5 : np.ndarray, optional Comment used to describe the data in the data array. The default is ‘’.
Returns¶
TimeHistoryArray A time history containing the specified haversine pulse signal
mimo_forward¶
Source: GitHub
def sdynpy.TimeHistoryArray.mimo_forward(self, transfer_function)Performs the forward mimo calculation via convolution.
Parameters¶
transfer_function : TransferFunctionArray or ImpulseResponseFunctionArray This is the FRFs that will be used in the forward problem. A matrix of IRFs is prefered, but FRFs can also be used, although the FRFs will be immediately converted to IRFs.
Returns¶
TimeHistoryArray
Response time histories
Raises¶
ValueError
If the sampling rates for the data and IRFs/FRFs don’t match.
ValueError
If the references in the IRFs/FRFs don’t match the supplied input
data.
mimo_inverse¶
Source: GitHub
def sdynpy.TimeHistoryArray.mimo_inverse(self, transfer_function, time_method='single_frame', cola_frame_length=None, cola_window='hann', cola_overlap=None, zero_pad_length=None, inverse_method='standard', response_weighting_matrix=None, reference_weighting_matrix=None, regularization_weighting_matrix=None, regularization_parameter=None, cond_num_threshold=None, num_retained_values=None, transfer_function_odd_samples=False)Performs the inverse source estimation for time domain (transient) problems using Fourier deconvolution. The response nodes used in the inverse source estimation are the ones contained in the supplied FRF matrix.
Parameters¶
transer_function : TransferFunctionArray or ImpulseResponseFunctionArray This is the FRFs that will be used in the inverse source estimation
time_method : str, optional The method to used to handle the time data for the inverse source estimation. The available options are:
single_frame - this method performs the Fourier deconvolution via an FFT on a single frame that encompases the entire time signal.
COLA - this method performs the Fourier deconvolution via a series of FFTs on relatively small frames of the time signal using a “constant overlap and add” method. This method may be faster than the single_frame method.
cola_frame_length : float, optional The frame length (in samples) if the COLA method is being used. The default frame length is Fs/df from the transfer function.
cola_window : str, optional The desired window for the COLA procedure, must exist in the scipy window library. The default is a hann window.
cola_overlap : int, optional The number of overlapping samples between measurement frames in the COLA procedure. If not specified, a default value of half the cola_frame_length is used.
zero_pad_length : int, optional The number of zeros used to pre and post pad the response data, to avoid convolution wrap-around error. The default is to use the “determine_zero_pad_length” function to determine the zero_pad_length.
inverse_method : str, optional The method to be used for the FRF matrix inversions. The available methods are:
standard - basic pseudo-inverse via numpy.linalg.pinv with the default rcond parameter, this is the default method
threshold - pseudo-inverse via numpy.linalg.pinv with a specified condition number threshold
tikhonov - pseudo-inverse using the Tikhonov regularization method
truncation - pseudo-inverse where a fixed number of singular values are retained for the inverse
response_weighting_matrix : sdpy.Matrix or np.ndarray, optional Not currently implemented
reference_weighting_matrix : sdpy.Matrix or np.ndarray, optional Not currently implemented
regularization_weighting_matrix : sdpy.Matrix, optional Matrix used to weight input degrees of freedom via Tikhonov regularization.
regularization_parameter : float or np.ndarray, optional Scaling parameter used on the regularization weighting matrix when the tikhonov method is chosen. A vector of regularization parameters can be provided so the regularization is different at each frequency line. The vector must match the length of the FRF abscissa in this case (either be size [num_lines,] or [num_lines, 1]).
cond_num_threshold : float or np.ndarray, optional Condition number used for SVD truncation when the threshold method is chosen. A vector of condition numbers can be provided so it varies as a function of frequency. The vector must match the length of the FRF abscissa in this case (either be size [num_lines,] or [num_lines, 1]).
num_retained_values : float or np.ndarray, optional Number of singular values to retain in the pseudo-inverse when the truncation method is chosen. A vector of can be provided so the number of retained values can change as a function of frequency. The vector must match the length of the FRF abscissa in this case (either be size [num_lines,] or [num_lines, 1]).
transfer_function_odd_samples : bool, optional If True, then it is assumed that the spectrum has been constructed from a signal with an odd number of samples. Note that this function uses the rfft function from scipy to compute the inverse fast fourier transform. The irfft function is not round-trip equivalent for odd functions, because by default it assumes an even signal length. For an odd signal length, the user must either specify transfer_function_odd_samples = True to make it round-trip equivalent.
Returns¶
TimeHistoryArray
Time history array of the estimated sources
Raises¶
NotImplementedError
If a response weighting matrix is supplied
NotImplementedError
If a reference weighting matrix is supplied
ValueError
If the sampling rates for the data and FRFs don’t match.
ValueError
If the number of responses in the FRFs don’t match the supplied response
data.
Notes¶
This function computes the time domain inputs required to match the target time traces using Fourier deconvolution, which is essentially a frequency domain problem. The general method is to compute the frequency spectrum of the target time traces, then solve the inverse problem in the time domain using the supplied FRFs (H^+ * X). The inverse of the FRF matrix is found using the same methods as the mimo_inverse function for the PowerSpectralDensityArray class. The input spectrum is then converted back to the time domain via a inverse fourier transform.
The 0 Hz component is explicitly rejected from the FRFs, so the estimated forces cannot include a 0 Hz component.
References¶
.. [1] Wikipedia, “Moore-Penrose inverse”. https://en.wikipedia.org/wiki/Moore–Penrose_inverse .. [2] A.N. Tithe, D.J. Thompson, The quantification of structure-borne transmission pathsby inverse methods. Part 2: Use of regularization techniques, Journal of Sound and Vibration, Volume 264, Issue 2, 2003, Pages 433-451, ISSN 0022-460X, Thite & Thompson (2003). .. [3] Wikipedia, “Ridge regression”. Ridge regression .. [4] Wikipedia, “Overlap-add Method”. Overlap-add method
overlap_and_add¶
Source: GitHub
def sdynpy.TimeHistoryArray.overlap_and_add(functions_to_overlap, overlap_samples)Creates a time history by overlapping and adding other time histories.
Parameters¶
functions_to_overlap : TimeHistoryArray or list of TimeHistoryArray A set of TimeHistoryArrays to overlap and add together. If a single TimeHistoryArray is specified, then the first dimension will be used to split the signal into segments. All TimeHistoryArrays must have the same shape and metadata, but need not have the same number of elements.
overlap_samples : int Number of samples to overlap the segments as they are added together
Returns¶
TimeHistoryArray
A TimeHistoryArray consisting of the signals overlapped and added
together.
Notes¶
All metadata is taken from the first signal. No checks are performed to make sure that the subsequent functions have common coordinates or abscissa spacing.
pseudorandom_signal¶
Source: GitHub
def sdynpy.TimeHistoryArray.pseudorandom_signal(cls, dt, signal_length, coordinates, min_frequency=None, max_frequency=None, signal_rms=1, frames=1, frequency_shape=None, different_realizations=False, comment1='', comment2='', comment3='', comment4='', comment5='')Generates a pseudorandom signal at the specified coordinates
Parameters¶
dt : float Abscissa spacing in the final signal.
signal_length : int Number of samples in the signal
coordinates : CoordinateArray Coordinate array used to generate the signal. If the last dimension of coordinates is not shape 1, then a new axis will be added to make it shape 1. The shape of the resulting TimeHistoryArray will be determined by the shape of the input coordinates.
min_frequency : float, optional Minimum frequency content in the signal. The default is the lowest nonzero frequency line.
max_frequency : float, optional Maximum frequency content in the signal. The default is the highest frequency content in the signal, e.g. the Nyquist frequency.
signal_rms : float or np.ndarray, optional RMS value for the generated signals. The default is 1. The shape of this value should be broadcastable with the size of the generated TimeHistoryArray if different RMS values are desired for each signal.
frames : int, optional Number of frames to generate. These will essentially be repeats of the first frame for the number of frames specified. The default is 1.
frequency_shape : function, optional An optional function that should accept a frequency value and return an amplitude at that frequency. The default is constant scaling across all frequency lines.
different_realizations : bool An optional argument that specifies whether or not different functions should have different realizations of the pseudorandom signal, or if they should all be identical.
comment1 : np.ndarray, optional Comment used to describe the data in the data array. The default is ‘’.
comment2 : np.ndarray, optional Comment used to describe the data in the data array. The default is ‘’.
comment3 : np.ndarray, optional Comment used to describe the data in the data array. The default is ‘’.
comment4 : np.ndarray, optional Comment used to describe the data in the data array. The default is ‘’.
comment5 : np.ndarray, optional Comment used to describe the data in the data array. The default is ‘’.
Returns¶
TimeHistoryArray A time history containing the specified pseudorandom signal
pulse_signal¶
Source: GitHub
def sdynpy.TimeHistoryArray.pulse_signal(cls, dt, signal_length, coordinates, pulse_width=None, pulse_time=None, pulse_peak=1, sine_exponent=1, frames=1, comment1='', comment2='', comment3='', comment4='', comment5='')Creates a pulse using a cosine function raised to a specified exponent
Parameters¶
dt : float Abscissa spacing in the final signal.
signal_length : int Number of samples in the signal
coordinates : CoordinateArray Coordinate array used to generate the signal. If the last dimension of coordinates is not shape 1, then a new axis will be added to make it shape 1. The shape of the resulting TimeHistoryArray will be determined by the shape of the input coordinates.
pulse_width : float, optional With of the pulse in the same units as
dt. The default is 5*dt.pulse_time : float, optional The time of the pulse’s occurance in the same units as
dt. The default is 5*dt.pulse_peak : float, optional The peak amplitude of the pulse. The default is 1.
sine_exponent : float, optional The exponent that the cosine function is raised to. The default is 1.
frames : int, optional Number of frames to generate. These will essentially be repeats of the first frame for the number of frames specified. The default is 1.
comment1 : np.ndarray, optional Comment used to describe the data in the data array. The default is ‘’.
comment2 : np.ndarray, optional Comment used to describe the data in the data array. The default is ‘’.
comment3 : np.ndarray, optional Comment used to describe the data in the data array. The default is ‘’.
comment4 : np.ndarray, optional Comment used to describe the data in the data array. The default is ‘’.
comment5 : np.ndarray, optional Comment used to describe the data in the data array. The default is ‘’.
Returns¶
TimeHistoryArray A time history containing the specified pulse signal
random_signal¶
Source: GitHub
def sdynpy.TimeHistoryArray.random_signal(cls, dt, signal_length, coordinates, min_frequency=None, max_frequency=None, signal_rms=1, frames=1, frequency_shape=None, comment1='', comment2='', comment3='', comment4='', comment5='')Generates a random signal with the specified parameters
Parameters¶
dt : float Abscissa spacing in the final signal.
signal_length : int Number of samples in the signal
coordinates : CoordinateArray Coordinate array used to generate the signal. If the last dimension of coordinates is not shape 1, then a new axis will be added to make it shape 1. The shape of the resulting TimeHistoryArray will be determined by the shape of the input coordinates.
min_frequency : float, optional Minimum frequency content in the signal. The default is the lowest nonzero frequency line.
max_frequency : float, optional Maximum frequency content in the signal. The default is the highest frequency content in the signal, e.g. the Nyquist frequency.
signal_rms : float or np.ndarray, optional RMS value for the generated signals. The default is 1. The shape of this value should be broadcastable with the size of the generated TimeHistoryArray if different RMS values are desired for each signal.
frames : int, optional Number of frames to generate. These will essentially be repeats of the first frame for the number of frames specified. The default is 1.
frequency_shape : function, optional An optional function that should accept a frequency value and return an amplitude at that frequency. The default is constant scaling across all frequency lines.
comment1 : np.ndarray, optional Comment used to describe the data in the data array. The default is ‘’.
comment2 : np.ndarray, optional Comment used to describe the data in the data array. The default is ‘’.
comment3 : np.ndarray, optional Comment used to describe the data in the data array. The default is ‘’.
comment4 : np.ndarray, optional Comment used to describe the data in the data array. The default is ‘’.
comment5 : np.ndarray, optional Comment used to describe the data in the data array. The default is ‘’.
Returns¶
TimeHistoryArray A time history containing the specified random signal
remove_rigid_body_motion¶
Source: GitHub
def sdynpy.TimeHistoryArray.remove_rigid_body_motion(self, geometry)Removes rigid body displacements from time data.
This function assumes the current TimeHistoryArray is a displacement signal and adds it to the geometry to create node positions over time, then it fits a rigid coordinate transformation to each time step and subtracts off that portion of the motion from the displacement signal.
Parameters¶
geometry : Geometry Geometry with which the node positions are computed
Returns¶
TimeHistoryArray
A TimeHistoryArray with the rigid body component of motion removed
resample¶
Source: GitHub
def sdynpy.TimeHistoryArray.resample(self, num_samples)Uses Scipy.signal.resample to resample the time history array
Parameters¶
num_samples : int The number of samples in the desired signal.
Returns¶
TimeHistoryArray
A TimeHistoryArray object with resampled abscissa and ordinate
rms¶
Source: GitHub
def sdynpy.TimeHistoryArray.rms(self)shift_signal¶
Source: GitHub
def sdynpy.TimeHistoryArray.shift_signal(self, time_shift)Shift a signal in time by a specified amount.
Utilizes the FFT shift theorem to move a signal in time.
Parameters¶
time_shift : float The time shift to apply to the signal. A negative value will cause features to occur earlier in time. A positive value will cause features to occur later in time.
Returns¶
shifted_signal : TimeHistoryArray A shifted version of the original signal.
sine_signal¶
Source: GitHub
def sdynpy.TimeHistoryArray.sine_signal(cls, dt, signal_length, coordinates, frequency, amplitude=1, phase=0, comment1='', comment2='', comment3='', comment4='', comment5='')Creates a sinusoidal signal with the specified parameters
Parameters¶
dt : float Abscissa spacing in the final signal.
signal_length : int Number of samples in the signal
coordinates : CoordinateArray Coordinate array used to generate the signal. If the last dimension of coordinates is not shape 1, then a new axis will be added to make it shape 1. The shape of the resulting TimeHistoryArray will be determined by the shape of the input coordinates.
frequency : float or np.ndarray Frequency of signal that will be generated. If multiple frequencies are specified, they must broadcast with the final size of the TimeHistoryArray.
amplitude : TYPE, optional Amplitude of signal that will be generated. If multiple amplitudes are specified, they must broadcast with the final size of the TimeHistoryArray. The default is 1.
phase : TYPE, optional Phase of signal that will be generated. If multiple phases are specified, they must broadcast with the final size of the TimeHistoryArray.. The default is 0.
comment1 : np.ndarray, optional Comment used to describe the data in the data array. The default is ‘’.
comment2 : np.ndarray, optional Comment used to describe the data in the data array. The default is ‘’.
comment3 : np.ndarray, optional Comment used to describe the data in the data array. The default is ‘’.
comment4 : np.ndarray, optional Comment used to describe the data in the data array. The default is ‘’.
comment5 : np.ndarray, optional Comment used to describe the data in the data array. The default is ‘’.
Returns¶
TimeHistoryArray A time history containing the specified sine signal
sine_sweep_signal¶
Source: GitHub
def sdynpy.TimeHistoryArray.sine_sweep_signal(cls, dt, coordinates, frequency_breakpoints, sweep_types, sweep_rates, amplitudes=1, phases=1, comment1='', comment2='', comment3='', comment4='', comment5='')split_into_frames¶
Source: GitHub
def sdynpy.TimeHistoryArray.split_into_frames(self, samples_per_frame=None, frame_length=None, overlap=None, overlap_samples=None, window=None, check_cola=False, allow_fractional_frames=False)Splits a time history into measurement frames with a given overlap and window function applied.
Parameters¶
samples_per_frame : int, optional Number of samples in each measurement frame. Either this argument or
frame_lengthmust be specified. If both or neither are specified, aValueErroris raised.frame_length : float, optional Length of each measurement frame in the same units as the
abscissafield (samples_per_frame=frame_length/self.abscissa_spacing). Either this argument orsamples_per_framemust be specified. If both or neither are specified, aValueErroris raised.overlap : float, optional Fraction of the measurement frame to overlap (i.e. 0.25 not 25 to overlap a quarter of the frame). Either this argument or
overlap_samplesmust be specified. If both are specified, aValueErroris raised. If neither are specified, no overlap is used.overlap_samples : int, optional Number of samples in the measurement frame to overlap. Either this argument or
overlap_samplesmust be specified. If both are specified, aValueErroris raised. If neither are specified, no overlap is used.window : str or tuple or array_like, optional Desired window to use. If window is a string or tuple, it is passed to
scipy.signal.get_windowto generate the window values, which are DFT-even by default. Seeget_windowfor a list of windows and required parameters. If window is array_like it will be used directly as the window and its length must besamples_per_frame. If not specified, no window will be applied.check_cola : bool, optional If
True, raise aValueErrorif the specified overlap and window function are not compatible with COLA. The default is False.allow_fractional_frames : bool, optional If
False(default), the signal will be split into a number of full frames, and any remaining fractional frame will be discarded. This will not allow COLA to be satisfied. IfTrue, fractional frames will be retained and zero padded to create a full frames.
Returns¶
TimeHistoryArray
Returns a new TimeHistoryArray with shape [num_frames,...] where
... is the shape of the original array.
srs¶
Source: GitHub
def sdynpy.TimeHistoryArray.srs(self, min_frequency=None, max_frequency=None, frequencies=None, damping=0.03, num_points=None, points_per_octave=12, srs_type='MMAA')Compute a shock response spectrum (SRS) from the time history
Parameters¶
min_frequency : float, optional Minimum frequency to compute the SRS. Either
frequenciesormin_frequencyandmax_frequencymust be specified.max_frequency : float, optional Maximum frequency to compute the SRS. Either
frequenciesormin_frequencyandmax_frequencymust be specified.frequencies : np.ndarray, optional Frequency lines at which to compute the SRS. Either
frequenciesormin_frequencyandmax_frequencymust be specified.damping : float, optional Fraction of critical damping to use in the SRS calculation (e.g. you should specify 0.03 to represent 3%, not 3). The default is 0.03.
num_points : int, optional Number of frequency lines to compute from
min_frequencytomax_frequency, log spaced between these two values. Ifmin_frequencyandmax_frequencyare specified, then eithernum_pointsorpoints_per_octavemust be specified. Iffrequenciesis specified, this argument is ignored.points_per_octave : float, optional Number of frequency lines per octave to compute from
min_frequencytomax_frequency. Ifmin_frequencyandmax_frequencyare specified, then eithernum_pointsorpoints_per_octavemust be specified. Iffrequenciesis specified, this argument is ignored. The default is 12.srs_type : str, optional A string encoding for the type of SRS to be computed. See notes for more information.
Returns¶
ShockResponseSpectrumArray
SRSs representing the current time histories. If
srs_typeis'all', then an extra dimension of 9 will be added to the front ofthe array, and the indices in that dimension will be different SRS
types.
Notes¶
The srs_type argument takes a 4 character string that specifies how
the SRS is computed.
stft¶
Source: GitHub
def sdynpy.TimeHistoryArray.stft(self, samples_per_frame=None, frame_length=None, overlap=None, overlap_samples=None, window=None, check_cola=False, allow_fractional_frames=False, norm='backward')Computes a Short-Time Fourier Transform (STFT)
The time history is split up into frames with specified length and computes the spectra for each frame.
Parameters¶
samples_per_frame : int, optional Number of samples in each measurement frame. Either this argument or
frame_lengthmust be specified. If both or neither are specified, aValueErroris raised.frame_length : float, optional Length of each measurement frame in the same units as the
abscissafield (samples_per_frame=frame_length/self.abscissa_spacing). Either this argument orsamples_per_framemust be specified. If both or neither are specified, aValueErroris raised.overlap : float, optional Fraction of the measurement frame to overlap (i.e. 0.25 not 25 to overlap a quarter of the frame). Either this argument or
overlap_samplesmust be specified. If both are specified, aValueErroris raised. If neither are specified, no overlap is used.overlap_samples : int, optional Number of samples in the measurement frame to overlap. Either this argument or
overlap_samplesmust be specified. If both are specified, aValueErroris raised. If neither are specified, no overlap is used.window : str or tuple or array_like, optional Desired window to use. If window is a string or tuple, it is passed to
scipy.signal.get_windowto generate the window values, which are DFT-even by default. Seeget_windowfor a list of windows and required parameters. If window is array_like it will be used directly as the window and its length must besamples_per_frame. If not specified, no window will be applied.check_cola : bool, optional If
True, raise aValueErrorif the specified overlap and window function are not compatible with COLA. The default is False.allow_fractional_frames : bool, optional If
False(default), the signal will be split into a number of full frames, and any remaining fractional frame will be discarded. This will not allow COLA to be satisfied. IfTrue, fractional frames will be retained and zero padded to create a full frames.
Returns¶
frame_abscissa : np.ndarray The abscissa values at the center of each of the STFT frames
stft : SpectrumArray A spectrum array with the first axis corresponding to the time values in
frame_abscissa.
to_rattlesnake_specification¶
Source: GitHub
def sdynpy.TimeHistoryArray.to_rattlesnake_specification(self, filename=None, coordinate_order=None, min_time=None, max_time=None)upsample¶
Source: GitHub
def sdynpy.TimeHistoryArray.upsample(self, factor)Upsamples a time history using frequency domain zero padding.
Parameters¶
factor : int The upsample factor.
Returns¶
TimeHistoryArray
A time history with a sample rate that is factor larger than the
original signal
vold_kalman_filter¶
Source: GitHub
def sdynpy.TimeHistoryArray.vold_kalman_filter(self, arguments, filter_order=None, bandwidth=None, method=None, block_size=None, overlap=0.5, plot_results=False, verbose=False)zefft¶
Source: GitHub
def sdynpy.TimeHistoryArray.zefft(self, num_steps=10, abscissa_range=None, closest_zero_crossing=True)Creates a FFTs from time histories with early times zeroed out
This computes an FFT time history with the first portion of the original time history zeroed out. This is useful for analyzing ring-downs of nonlinear systems to identify how the system changes over time.
Parameters¶
num_steps : int, optional The number of different levels of zeroing to compute for each signal, by default 10
abscissa_range : iterable, optional A length-2 iterable defining the earliest and latest abscissa to use as a time to zero data until, by default, it will spread
num_stepszero locations over the entire time signal.closest_zero_crossing : bool, optional If True, zeroing will start at a zero crossing in a signal that is closest to the specified zero time. Otherwise, the zeroing will start exactly where specified, by default True
Returns¶
SpectrumArray
A SpectrumArray with an additional dimension for the
num_stepsspecified. These are outputs of thezero_early_timemethod with thefftmethodsubsequently applied.
zero_early_time¶
Source: GitHub
def sdynpy.TimeHistoryArray.zero_early_time(self, num_steps=10, abscissa_range=None, closest_zero_crossing=True)Creates a time histories with early times zeroed out
This computes a time history with the first portion of the original time history zeroed out. This is useful for analyzing ring-downs of nonlinear systems to identify how the system changes over time.
Parameters¶
num_steps : int, optional The number of different levels of zeroing to compute for each signal, by default 10
abscissa_range : iterable, optional A length-2 iterable defining the earliest and latest abscissa to use as a time to zero data until, by default, it will spread
num_stepszero locations over the entire time signal.closest_zero_crossing : bool, optional If True, zeroing will start at a zero crossing in a signal that is closest to the specified zero time. Otherwise, the zeroing will start exactly where specified, by default True
Returns¶
TimeHistoryArray
A TimeHistoryArray with an additional dimension for the
num_stepsspecified. Each of thenum_stepsset ofsignals has a different initial zero point.
- Thite, A. N., & Thompson, D. J. (2003). The quantification of structure-borne transmission paths by inverse methods. Part 2: Use of regularization techniques. Journal of Sound and Vibration, 264(2), 433–451. 10.1016/s0022-460x(02)01203-8