sdynpy.core.sdynpy_data.PowerSpectralDensityArray
- class PowerSpectralDensityArray(shape, nelements, buffer=None, offset=0, strides=None, order=None)[source]
Bases:
NDDataArray
Data array used to store power spectral density arrays
- __init__()
Methods
angle
()Computes the angle of a PSD matrix
bandwidth_average
(band_lb, band_ub)Integrates the PSD over frequency to get the power spectrum for each frequency bin (line)
Computes the coherence of a PSD matrix
compare_asds
([figure_kwargs, linewidth])Plot the diagonals of the CPSD matrix, as well as the level
error_summary
([figure_kwargs, linewidth, ...])Plots an error summary compared to the current array
eye
(frequencies, coordinates[, rms, ...])Computes a diagonal CPSD matrix
from_time_data
(response_data[, ...])Computes a PSD matrix from reference and response time histories
generate_time_history
([time_length, ...])Generates a time history from a CPSD matrix
get_asd
()Get functions where the response coordinate is equal to the reference coordinate
mimo_forward
(transfer_function)Compute the forward MIMO problem Gxx = Hxv@Gvv@Hxv*
mimo_inverse
(transfer_function[, method, ...])Computes input estimation for MIMO random vibration problems
plot_asds
([figure_kwargs, linewidth])Plots the magnitude, coherence, and phase of a CPSD matrix.
plot_singular_values
([rcond, min_freqency, ...])Plot the singular values of an FRF matrix with a visualization of the rcond tolerance
rms
()Compute RMSs of the PSDs using the diagonals
set_coherence_phase
(coherence_array, angle_array)Sets the coherence and phase of a PSD matrix while maintaining the ASDs
svd
([full_matrices, compute_uv, as_matrix])Compute the SVD of the provided CPSD matrix
to_rattlesnake_specification
([filename, ...])Attributes
Returns the function type of the data array
- angle()[source]
Computes the angle of a PSD matrix
- Returns
Data array consisting of the angle of each function at each frequency line
- Return type
- bandwidth_average(band_lb, band_ub)[source]
Integrates the PSD over frequency to get the power spectrum for each frequency bin (line)
- Parameters
band_lb (ndarray) – (n_bands,1) array of bandwidth lower bounds
band_ub (ndarray) – (n_bands,1) array of bandwidth upper bounds
- Returns
PowerSpectralDensityArray with abscissa given by the mean of band_lb
and band_ub
Notes
Determines which freq bins (lines) contribute to each band. Contribute means the freq bin is at least partially within the band limits
The portion of the bin which contributes to the band is computed based multiplied by the fraction of the contributing frequency to get how much bin PS adds to the band PS
- coherence()[source]
Computes the coherence of a PSD matrix
- Raises
ValueError – If abscissa are not consistent.
- Returns
CoherenceArray containing the values of coherence for each function.
- Return type
- static compare_asds(figure_kwargs={}, linewidth=1, **cpsd_matrices)[source]
Plot the diagonals of the CPSD matrix, as well as the level
- Parameters
figure_kwargs (dict, optional) – Optional arguments to use when creating the figure. The default is {}.
linewidth (float, optional) – Width of plotted lines. The default is 1.
**cpsd_matrices (PowerSpectralDensityArray) – PSDs to plot. Only gets plotted if response and reference are identical. The key will be used as the label with _ replaced by a space.
- Raises
ValueError – If degrees of freedom are not consistent between PSDs
- Return type
None.
- error_summary(figure_kwargs={}, linewidth=1, plot_kwargs={}, **cpsd_matrices)[source]
Plots an error summary compared to the current array
- Parameters
figure_kwargs (dict, optional) – Arguments to use when creating the figure. The default is {}.
linewidth (float, optional) – Widths of the lines on the plot. The default is 1.
plot_kwargs (dict, optional) – Arguments to use when plotting the lines. The default is {}.
**cpsd_matrices (PowerSpectralDensityArray) – Data to compare against the current CPSD matrix. The keys will be used as labels with _ replaced with a space.
- Raises
ValueError – If CPSD abscissa do not match
- Returns
A tuple of dictionaries of error metrics
- Return type
Error Metrics
- classmethod eye(frequencies, coordinates, rms=None, full_matrix=False, breakpoint_frequencies=None, breakpoint_levels=None, breakpoint_interpolation='lin', min_frequency=0.0, max_frequency=None)[source]
Computes a diagonal CPSD matrix
- Parameters
frequencies (ndarray) – Frequencies at which the CPSD should be constructed
coordinates (CoordinateArray) – CoordinateArray to use to set the CPSD values
rms (ndarray, optional) – Value to scale the RMS of each CPSD to
full_matrix (bool, optional) – If True, a full, square CPSD matrix will be computed. If False, only the ASDs will be computed. The default is False.
breakpoint_frequencies (iterable, optional) – A list of frequencies that breakpoints are defined at.
breakpoint_levels (iterable, optional) – A list of levels that breakpoints are defined at
breakpoint_interpolation (str, optional) – ‘lin’ or ‘log’ to specify the type of interpolation. The default is ‘lin’.
min_frequency (float, optional) – Low frequency cutoff for the CPSD. Frequency lines below this value will be set to zero.
max_frequency (float, optional) – High frequency cutoff for the CPSD. Frequency lines above this value will be set to zero.
- Raises
ValueError – If invalid interpolation is specified, or if RMS is specified with inconsistent frequency spacing.
- Returns
A set of PSDs.
- Return type
- static from_time_data(response_data: TimeHistoryArray, samples_per_average: Optional[int] = None, overlap: float = 0.0, window=array([1.]), reference_data: Optional[TimeHistoryArray] = None, only_asds=False)[source]
Computes a PSD matrix from reference and response time histories
- Parameters
response_data (TimeHistoryArray) – Time data to be used as responses
samples_per_average (int, optional) – Number of samples used to split up the signals into averages. The default is None, meaning the data is treated as a single measurement frame.
overlap (float, optional) – The overlap as a fraction of the frame (e.g. 0.5 specifies 50% overlap). The default is 0.0, meaning no overlap is used.
window (np.ndarray or str, optional) – A 1D ndarray with length samples_per_average that specifies the coefficients of the window. A Hann window is applied if not specified. If a string is specified, then the window will be obtained from scipy.
reference_data (TimeHistoryArray) – Time data to be used as reference. If not specified, the response data will be used as references, resulting in a square CPSD matrix.
- Raises
ValueError – Raised if reference and response functions do not have consistent abscissa
- Returns
A PSD array computed from the specified reference and response signals.
- Return type
- property function_type
Returns the function type of the data array
- generate_time_history(time_length=None, output_oversample=1)[source]
Generates a time history from a CPSD matrix
- Parameters
time_length (float, optional) – The length (in time, not samples) of the signal. If not specified, the signal length will be based on the frequency spacing and nyquist frequency of the CPSD matrix. If specified, a signal will be constructed using constant overlap and add techniques. A whole number of realizations will be constructed, so the output signal can be longer than the time_length specified.
output_oversample (int, optional) – Oversample factor applied to the output signal. The default is 1.
- Raises
ValueError – If the entries in the CPSD matrix do not have consistent abscissa or equally spaced frequency bins.
- Returns
time_history – A time history satisfying the properties of the CPSD matrix.
- Return type
- get_asd()[source]
Get functions where the response coordinate is equal to the reference coordinate
- Returns
PowerSpectralDensityArrays where the response is equal to the reference
- Return type
- mimo_forward(transfer_function)[source]
Compute the forward MIMO problem Gxx = Hxv@Gvv@Hxv*
- Parameters
transfer_function (TransferFunctionArray) – Transfer function used to transform the input matrix to the response matrix
- Raises
ValueError – If abscissa do not match between self and transfer function
- Returns
Response CPSD matrix
- Return type
- mimo_inverse(transfer_function, 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)[source]
Computes input estimation for MIMO random vibration problems
- Parameters
transfer_function (TransferFunctionArray) – System transfer functions used to estimate the input from the given response matrix
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, optional) – Diagonal matrix used to weight response degrees of freedom (to solve the problem as a weight least squares) by multiplying the rows of the FRF matrix by a scalar weights. This matrix can also be a 3D matrix such that the the weights are different for each frequency line. The matrix should be sized [number of lines, number of references, number of references], where the number of lines either be one (the same weights at all frequencies) or the length of the abscissa (for the case where a 3D matrix is supplied).
reference_weighting_matrix (sdpy.Matrix, optional) – Diagonal matrix used to weight reference degrees of freedom (generally for normalization) by multiplying the columns of the FRF matrix by a scalar weights. This matrix can also be a 3D matrix such that the the weights are different for each frequency line. The matrix should be sized [number of lines, number of references, number of references], where the number of lines either be one (the same weights at all frequencies) or the length of the abscissa (for the case where a 3D matrix is supplied).
regularization_weighting_matrix (sdpy.Matrix, optional) – Matrix used to weight input degrees of freedom via Tikhonov regularization. This matrix can also be a 3D matrix such that the the weights are different for each frequency line. The matrix should be sized [number of lines, number of references, number of references], where the number of lines either be one (the same weights at all frequencies) or the length of the abscissa (for the case where a 3D matrix is supplied).
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 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 abscissa in this case.
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 abscissa in this case.
- Raises
ValueError – If Abscissa are not consistent
- Returns
Input CPSD matrix
- Return type
Notes
This function solves the MIMO problem Gxx = Hxv@Gvv@Hxv^* using the pseudoinverse. Gvv = Hxv^+@Gxx@Hxv^+^*, where Gvv is the source.
References
- 1
Wikipedia, “Moore-Penrose inverse”. https://en.wikipedia.org/wiki/Moore%E2%80%93Penrose_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, https://doi.org/10.1016/S0022-460X(02)01203-8.
- 3
Wikipedia, “Ridge regression”. https://en.wikipedia.org/wiki/Ridge_regression
- plot_magnitude_coherence_phase(compare_data=None, plot_axes=False, sharex=True, sharey=True, logx=False, logy=True, magnitude_plot_kwargs={}, coherence_plot_kwargs={}, angle_plot_kwargs={}, figure_kwargs={})[source]
Plots the magnitude, coherence, and phase of a CPSD matrix.
Coherence is plotted on the upper triangle, phase on the lower triangle, and magnitude on the diagonal.
- Parameters
compare_data (PowerSpectralDensityArray, optional) – An optional dataset to compare against. The default is None.
plot_axes (bool, optional) – If True, axes tick labels will be plotted. If false, the plots will be pushed right against one another without room for labels. The default is False.
sharex (bool, optional) – If True, all plots will share the same range on the X axis. The default is True.
sharey (bool, optional) – If true, all plots of the same type will share the same range on the Y axis. The default is True.
logx (bool, optional) – If true, the x-axis will be logarithmic. The default is False.
logy (bool, optional) – If true, the y-axis on magnitude plots will be logrithmic. The default is True.
magnitude_plot_kwargs (dict, optional) – Optional keywards to use when plotting magnitude. The default is {}.
coherence_plot_kwargs (dict, optional) – Optional keywards to use when plotting coherence. The default is {}.
angle_plot_kwargs (dict, optional) – Optional keywards to use when plotting phase. The default is {}.
figure_kwargs (dict, optional) – Optional keywards to use when creating the figure. The default is {}.
- Return type
None.
- plot_singular_values(rcond=None, min_freqency=None, max_frequency=None)[source]
Plot the singular values of an FRF matrix with a visualization of the rcond tolerance
- Parameters
rcond (value of float, optional) – Cutoff for small singular values. Implemented such that the cutoff is rcond* largest_singular_value (the same as np.linalg.pinv). This is to visualize the effect of rcond and is used for display purposes only.
min_frequency (float, optional) – Minimum frequency to plot
max_frequency (float, optional) – Maximum frequency to plot
- rms()[source]
Compute RMSs of the PSDs using the diagonals
- Returns
RMS values for the ASDS
- Return type
ndarray
- set_coherence_phase(coherence_array, angle_array)[source]
Sets the coherence and phase of a PSD matrix while maintaining the ASDs
- Parameters
coherence_array (CoherenceArray) – Coherence to which the PSD will be set
angle_array (NDDataArray) – Phase to which the PSD will be set
- Returns
output – PSD with coherence and phase matching that of the input argument
- Return type
- svd(full_matrices=True, compute_uv=True, as_matrix=True)[source]
Compute the SVD of the provided CPSD matrix
- Parameters
full_matrices (bool, optional) – This is an optional input for np.linalg.svd
compute_uv (bool, optional) – This is an optional input for np.linalg.svd
as_matrix (bool, optional) – If True, matrices are returned as a SDynPy Matrix class with named rows and columns. Otherwise, a simple numpy array is returned
- Returns
u (ndarray) – Left hand singular vectors, sized […, num_responses, num_responses]. Only returned when compute_uv is True.
s (ndarray) – Singular values, sized […, num_references]
vh (ndarray) – Right hand singular vectors, sized […, num_references, num_references]. Only returned when compute_uv is True.