Defined as:
sdynpy.core.sdynpy_data.TransferFunctionArrayModule:
sdynpy.core.sdynpy_dataSource: GitHub
Parent:
sdynpy.NDDataArrayParent:
sdynpy.SdynpyArrayParent:
numpy.ndarray
Signature¶
class sdynpy.TransferFunctionArray(shape, nelements, buffer=None, offset=0, strides=None, order=None)Data array used to store transfer functions (for example FRFs)
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 reference and response transformations to the transfer functions. |
block_diagonal_frf | Assembles a block diagonal FRF TransferFunctionArray from the supplied FRFs. |
compute_cmif | Computes a complex mode indicator function from the TransferFunctionArray |
compute_mif | Compute a mode indicator functions from the transfer functions |
compute_mmif | Computes a Multi Mode indicator function from the TransferFunctionArray |
compute_nmif | Computes a normal mode indicator function from the TransferFunctionArray |
delay_response | Adjusts the FRF phases as if the response had been shifted dt in time |
enforce_causality | Enforces causality on the frequency response function via a conversion to a impulse response function, applying a cutoff window, then converting back to a frequency response function. |
from_exodus | Reads FRF data from a Sierra/SD ModalFRF Exodus file |
from_time_data | Computes a transfer function from reference and response time histories |
ifft | Converts frequency response functions to impulse response functions via an inverse fourier transform. |
interpolate_by_zero_pad | Interpolates a transfer function by zero padding or truncating its impulse response |
plot | Plot the transfer functions |
plot_cond_num | Plots the condition number of the FRF matrix |
plot_singular_values | Plot the singular values of an FRF matrix with a visualization of the rcond tolerance |
plot_with_coherence | |
substructure_by_constraint_matrix | Performs frequency based substructuring using the supplied constraint matrix with accompanying dof list. |
substructure_by_coordinate | Performs frequency based substructuring by constraining pairs of degrees of freedom |
svd | Compute the SVD of the provided FRF matrix |
apply_transformation¶
Source: GitHub
def sdynpy.TransferFunctionArray.apply_transformation(self, response_transformation=None, reference_transformation=None, invert_response_transformation=False, invert_reference_transformation=True)Applies reference and response transformations to the transfer functions.
Parameters¶
response_transformation : Matrix, optional The response transformation to apply to the (rows of the) transfer functions. It should be a SDynPy matrix object with the “transformed” coordinates on the rows and the “physical” coordinates on the columns. The matrix can be either 2D or 3D (for a frequency dependent transform).
reference_transformation : Matrix, optional The reference transformation to apply to the (columns of the) transfer functions. It should be a SDynPy matrix object with the “transformed” coordinates on the rows and the “physical” coordinates on the columns. The matrix can be either 2D or 3D (for a frequency dependent transform).
invert_response_transformation : bool, optional Whether or not to invert the response transformation when applying it to the transfer functions. The default is false, which is standard practice. The row/column ordering in the reference transformation should be flipped if this is set to true.
invert_reference_transformation : bool, optional Whether or not to invert the reference transformation when applying it to the transfer functions. The default is true, which is standard practice. The row/column ordering in the reference transformation should be flipped if this is set to false.
Returns¶
transformed_transfer_function : TransferFunctionArray The transfer functions with the transformations applied.
Raises¶
ValueError
If the physical degrees of freedom in the transformations don’t
match the transfer functions
Notes¶
This method can be used with just a response transformation, just a reference transformation, or both a response and reference transformation. The transformation will be set to identity if it is not supplied.
References¶
.. [1] M. Van der Seijs, D. van den Bosch, D. Rixen, and D. Klerk, “An improved methodology for the virtual point transformation of measured frequency response functions in dynamic substructuring,” in Proceedings of the 4th International Conference on Computational Methods in Structural Dynamics and Earthquake Engineering, Kos Island, 2013, pp. 4334-4347, doi: 10.7712/120113.4816.C1539.
block_diagonal_frf¶
Source: GitHub
def sdynpy.TransferFunctionArray.block_diagonal_frf(cls, component_frfs: tuple, coordinate_node_offset: int = 0)Assembles a block diagonal FRF TransferFunctionArray from the supplied FRFs.
Parameters¶
component_frfs : iterable of TransferFunctionArrays A set of FRFs to be assembled into a block diagonal FRF matrix.
coordinate_node_offset : int, optional An offset that is applied to the nodes in the supplied FRFs. The default is zero, meaning that an offset is not applied.
Returns¶
TransferFunctionArray
The FRFs organized in block diagonal format.
Raises¶
ValueError
If the objects in component FRFs are not TransferFunctionArrays.
ValueError
If the TransferFunctionArrays do not share the same abscissa.
Notes¶
The coordinate_node_offset is incremented for each system. For example, if the first set of FRFs has nodes [1,2,3,4], the second set of FRFs has nodes [3,4,5,6], and the supplied offset is 10, the node numbers in the returned block diagonal FRFs would be [11,12,13,14,23,24,25,26].
compute_cmif¶
Source: GitHub
def sdynpy.TransferFunctionArray.compute_cmif(self, part='both', tracking=None)Computes a complex mode indicator function from the TransferFunctionArray
Parameters¶
part : str, optional Specifies which part(s) of the transfer functions are used to compute the CMIF. Can be ‘real’, ‘imag’, or ‘both’. The default is ‘both’.
tracking : str or None, optional Specifies if any singular value tracking should be used. Can be ‘left’ or ‘right’. The default is None.
Returns¶
output_array : ModeIndicatorFunctionArray Complex Mode Indicator Function
Raises¶
ValueError
Raised if an invalid tracking is specified
compute_mif¶
Source: GitHub
def sdynpy.TransferFunctionArray.compute_mif(self, mif_type, *mif_args, **mif_kwargs)Compute a mode indicator functions from the transfer functions
Parameters¶
mif_type : str Mode indicator function type, one of ‘cmif’,‘nmif’, or ‘mmif’
*mif_args : list Arguments passed to the compute_*mif function
**mif_kwargs : dict Keyword arguments passed to the compute_*mif function
Returns¶
ModeIndicatorFunctionArray
Mode indicator function
Raises¶
ValueError
If an invalid mif name is provided.
compute_mmif¶
Source: GitHub
def sdynpy.TransferFunctionArray.compute_mmif(self, part='real', mass_matrix=None)Computes a Multi Mode indicator function from the TransferFunctionArray
Parameters¶
part : str, optional Specifies which part(s) of the transfer functions are used to compute the NMIF. Can be ‘real’ or ‘imag’. The default is ‘real’.
mass_matrix : np.ndarray, optional Matrix used to compute the MMIF, Identity is used if not specified
Returns¶
output_array : ModeIndicatorFunctionArray Multi Mode Indicator Function
Raises¶
ValueError
Raised if an invalid part is specified
compute_nmif¶
Source: GitHub
def sdynpy.TransferFunctionArray.compute_nmif(self, part='real')Computes a normal mode indicator function from the TransferFunctionArray
Parameters¶
part : str, optional Specifies which part(s) of the transfer functions are used to compute the NMIF. Can be ‘real’ or ‘imag’. The default is ‘real’.
Returns¶
output_array : ModeIndicatorFunctionArray Normal Mode Indicator Function
Raises¶
ValueError
Raised if an invalid part is specified
delay_response¶
Source: GitHub
def sdynpy.TransferFunctionArray.delay_response(self, dt)Adjusts the FRF phases as if the response had been shifted dt in time
Parameters¶
dt : float Time shift to apply to the responses
Returns¶
shifted_transfer_function : TransferFunctionArray A copy of the transfer function with the phase shifted
enforce_causality¶
Source: GitHub
def sdynpy.TransferFunctionArray.enforce_causality(self, method='exponential_taper', window_parameter=None, end_of_ringdown=None)Enforces causality on the frequency response function via a conversion to a impulse response function, applying a cutoff window, then converting back to a frequency response function.
Parameters¶
method : str The window type that is applied to the data to enforce causality. Note that these options are not necessarily traditional windows (used for data processing). The current options are:
exponential_taper (default) - this applies a exponential taper to the end of a boxcar window on the IRF.
boxcar - this applies a boxcar (uniform) window to the IRF with the cuttoff at a specified sample.
exponential - this applies an exponential window to the IRF with the 40 dB down point (of the window) at a specified sample. Care should be taken when using this window type, since it can lead to erratic behavior.
window_parameter : int, optional This is a parameter that defines the window for the causality enforcement. Methods exist to define this parameter automatically if it isn’t provided. The behaviors for the options are:
boxcar - the window_paramter is the sample after which the IRF is set to zero. It is the same as the end_of_ringdown parameter for this window type.
exponential - the window_parameter is where the 40 dB down point is for the window. It is the same as the end_of_ringdown parameter for this window type.
exponential_taper - the window_parameter is where the end point of the window (where the amplitude is 0.001), as defined by the number of samples after the uniform section of the window.
end_of_ringdown : int, optional This is a parameter that defines the end of the uniform section of the exponetional_taper window. It is not used for either the boxcar or exponential window. Methods exist to define this parameter automatically if it isn’t provided.
Returns¶
TransferFunctionArray
The FRF with causality enforced.
Notes¶
This is a wrapper around the method in the impulse response function class and it may be wiser to use that function instead.
Although optional, it is best practice for the user to supply a parameter for the end_of_ringdown variable if the “exponential_taper” method is being used or a window_parameter if the “exponential” or “boxcar” methods are being used. The code will attempt to find the end of the ring-down in the IRF and use use that as the end_of_ringdown parameter for the “exponential_taper” window or the window_parameter for the exponential and boxcar windows.
It is not suggested that the user provide a window_paramter if the “exponential_taper” method is being used, since the default is likely the most logical choice.
References¶
.. [1] Zvonkin, M. (2015). Methods for checking and enforcing physical quality of linear electrical network models
[Masters Theses, Missouri University of Science and Technology], Missouri S&T Scholars’ Mine, https://
from_exodus¶
Source: GitHub
def sdynpy.TransferFunctionArray.from_exodus(cls, exo, reference_coordinate=None, xval='DispX', xvali='ImagDispX', yval='DispY', yvali='ImagDispY', zval='DispZ', zvali='ImagDispZ')Reads FRF data from a Sierra/SD ModalFRF Exodus file
Parameters¶
exo : Exodus or ExodusInMemory The exodus data from which FRF data will be created.
reference_coord : CoordinateArray The coordinate to be assigned as the reference coordinate, by default 0
xval : str, optional The variable name representing the real part of the X value, by default ‘DispX’
xvali : str, optional The variable name representing the imaginary part of the X value, by default ‘ImagDispX’
yval : str, optional The variable name representing the real part of the Y value, by default ‘DispY’
yvali : str, optional The variable name representing the imaginary part of the Y value, by default ‘ImagDispY’
zval : str, optional The variable name representing the real part of the Z value, by default ‘DispZ’
zvali : str, optional The variable name representing the imaginary part of the Z value, by default ‘ImagDispZ’
Returns¶
TransferFunctionArray
FRF data from the exodus file
from_time_data¶
Source: GitHub
def sdynpy.TransferFunctionArray.from_time_data(reference_data: sdynpy.core.sdynpy_data.TimeHistoryArray, response_data: sdynpy.core.sdynpy_data.TimeHistoryArray, samples_per_average: int = None, overlap: float = 0.0, method: str = 'H1', window=array([1.]), return_model_data=False, **timedata2frf_kwargs)Computes a transfer function from reference and response time histories
Parameters¶
reference_data : TimeHistoryArray Time data to be used as a reference
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.
method : str, optional The method for creating the frequency response function. ‘H1’ is default if not specified. samples_per_average, overlap, and window are not used if method==‘LRM’.
window : np.ndarray or str, optional A 1D ndarray with length samples_per_average that specifies the coefficients of the window. No window is applied if not specified. If a string is specified, then the window will be obtained from scipy.
**timedata2frf_kwargs : various Additional keyword arguments that may be passed into the timedata2frf function in sdynpy.frf. If method==‘LRM’, see also frf_local_model in sdynpy.lrm for more options.
Returns¶
TransferFunctionArray
A transfer function array computed from the specified references and
responses.
Raises¶
ValueError
Raised if reference and response functions do not have consistent
abscissa
ifft¶
Source: GitHub
def sdynpy.TransferFunctionArray.ifft(self, norm='backward', odd_num_samples=False, **scipy_irfft_kwargs)Converts frequency response functions to impulse response functions via an inverse fourier transform.
Returns¶
ImpulseResponseFunctionArray
The impulse response function array computed from the transfer function
array.
Raises¶
Warning
Raised if the transfer function array does not have evenly spaced
frequency data in the 0-maximum frequency range, but appears to have been
high pass filtered.
ValueError
Raised if the transfer function array does not have evenly spaced
frequency data in the 0-maximum frequency range and it does not appear
to have been high pass filtered.
Paramters¶
norm : str, optional The type of normalization applied to the fft computation. odd_num_samples : bool, optional If True, then it is assumed that the output signal has an odd number of samples, meaning the signal will have a length of 2*(m-1)+1 where m is the number of frequency lines. Otherwise, the default value of 2*(m-1) is used, assuming an even signal. This is ignored if num_samples is specified. scipy_irfft_kwargs : Additional keywords that will be passed to SciPy’s irfft function.
interpolate_by_zero_pad¶
Source: GitHub
def sdynpy.TransferFunctionArray.interpolate_by_zero_pad(self, irf_padded_length, return_irf=False, odd_num_samples=False)Interpolates a transfer function by zero padding or truncating its impulse response
Parameters¶
irf_padded_length : int Length of the final zero-padded impulse response function
return_irf : bool, optional If True, the zero-padded impulse response function will be returned. If False, it will be transformed back to a transfer function prior to being returned.
odd_num_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 odd_num_samples = True to make it round-trip equivalent.
Returns¶
TransferFunctionArray or ImpulseResponseFunctionArray Transfer function array with appropriately spaced abscissa
Notes¶
This function will automatically set the last frequency line of the
TransferFunctionArray to zero because it won’t be accurate anyway.
If irf_padded_length is less than the current function’s num_elements,
then it will be truncated instead of zero-padded.
plot¶
Source: GitHub
def sdynpy.TransferFunctionArray.plot(self, one_axis=True, part=None, subplots_kwargs={}, plot_kwargs={}, abscissa_markers=None, abscissa_marker_labels=None, abscissa_marker_type='vline', abscissa_marker_plot_kwargs={})Plot the transfer functions
Parameters¶
one_axis : bool, optional Set to True to plot all data on one axis. Set to False to plot data on multiple subplots. one_axis can also be set to a matplotlib axis to plot data on an existing axis. The default is True.
part : str, optional The part of the FRF to plot. This can be, ‘real’, ‘imag’ or ‘imaginary’, ‘mag’ or ‘magnitude’, or ‘phase’. If not specified, magnitude and phase will be plotted if
one_axisis True, and magnitude will be plotted ifone_axisis False.subplots_kwargs : dict, optional Keywords passed to the matplotlib subplots function to create the figure and axes. The default is {}.
plot_kwargs : dict, optional Keywords passed to the matplotlib plot function. The default is {}.
abscissa_markers : ndarray, optional Array containing abscissa values to mark on the plot to denote significant events.
abscissa_marker_labels : str or ndarray Array of strings to label the abscissa_markers with, or alternatively a format string that accepts index and abscissa inputs (e.g. ‘{index:}: {abscissa:0.2f}’). By default no label will be applied.
abscissa_marker_type : str The type of marker to use. This can either be the string ‘vline’ or a valid matplotlib symbol specifier (e.g. ‘o’, ‘x’, ‘.’).
abscissa_marker_plot_kwargs : dict Additional keyword arguments used when plotting the abscissa label markers.
Returns¶
axis : matplotlib axis or array of axes On which the data were plotted
plot_cond_num¶
Source: GitHub
def sdynpy.TransferFunctionArray.plot_cond_num(self, number_retained_values=None, min_frequency=None, max_frequency=None)Plots the condition number of the FRF matrix
Parameters¶
min_freqency : float, optional Minimum frequency to plot. The default is None.
max_frequency : float, optional Maximum frequency to plot. The default is None.
Returns¶
None.
plot_singular_values¶
Source: GitHub
def sdynpy.TransferFunctionArray.plot_singular_values(self, rcond=None, condition_number=None, number_retained_values=None, regularization_parameter=None, min_frequency=None, max_frequency=None)Plot the singular values of an FRF matrix with a visualization of the rcond tolerance
Parameters¶
rcond : float or ndarray, 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.
condition_number : float or ndarray, optional Condition number threshold for small singular values. The condition number is the reciprocal of rcond. This is to visualize the effect of condition number threshold and is used for display purposes only.
number_retained_values : float or ndarray, optional Cutoff for small singular values a an integer value of number of values to retain. This is to visualize the effect of singular value truncation and is used for display purposes only.
regularization_parameter : float or ndarray, optional Regularization parameter to compute the modified singular values. This is to visualize the effect of Tikhonov regularization and is used for display purposes only.
min_frequency : float, optional Minimum frequency to plot
max_frequency : float, optional Maximum frequency to plot
plot_with_coherence¶
Source: GitHub
def sdynpy.TransferFunctionArray.plot_with_coherence(self, coherence, part=None, subplots_kwargs={}, plot_kwargs={})substructure_by_constraint_matrix¶
Source: GitHub
def sdynpy.TransferFunctionArray.substructure_by_constraint_matrix(self, dofs, constraint_matrix)Performs frequency based substructuring using the supplied constraint matrix with accompanying dof list.
Parameters¶
dofs : CoordinateArray Coordinates to use in the constraints.
constraint_matrix : np.ndarray Constraints to apply to the frequency response functions. Should be sized [number of constraints, number of dofs].
Returns¶
constrained_frfs : TransferFunctionArray Constrained Frequency Response Functions
Raises¶
ValueError
If listed degrees of freedom are not found in the function.
substructure_by_coordinate¶
Source: GitHub
def sdynpy.TransferFunctionArray.substructure_by_coordinate(self, dof_pairs)Performs frequency based substructuring by constraining pairs of degrees of freedom
Parameters¶
dof_pairs : CoordinateArray or None Pairs of coordinates to constrain together. To constain a coordinate to ground (i.e. fix it so it cannot move), the coordinate should be paired with None. This should be size [number of constraints, 2].
Returns¶
TransferFunctionArray
Constrained frequency response functions.
svd¶
Source: GitHub
def sdynpy.TransferFunctionArray.svd(self, full_matrices=True, compute_uv=True, as_matrix=True)Compute the SVD of the provided FRF matrix
Parameters¶
full_matrices : bool, optional This is an optional input for np.linalg.svd, the default for this function is true (which differs from the np.linalg.svd function).
compute_uv : bool, optional This is an optional input for np.linalg.svd, the default for this function is true (which differs from the np.linalg.svd function).
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 or Matrix Left hand singular vectors, sized [..., num_responses, num_responses]. Only returned when compute_uv is True.
s : ndarray or Matrix Singular values, sized [..., num_references]
vh : ndarray or Matrix Right hand singular vectors, sized [..., num_references, num_references]. Only returned when compute_uv is True.