sdynpy.core.sdynpy_data.TransferFunctionArray

class TransferFunctionArray(shape, nelements, buffer=None, offset=0, strides=None, order=None)[source]

Bases: NDDataArray

Data array used to store transfer functions (for example FRFs)

__init__()

Methods

compute_cmif([part, tracking])

Computes a complex mode indicator function from the TransferFunctionArray

compute_mif(mif_type, *mif_args, **mif_kwargs)

Compute a mode indicator functions from the transfer functions

compute_mmif([part, mass_matrix])

Computes a Multi Mode indicator function from the TransferFunctionArray

compute_nmif([part])

Computes a normal mode indicator function from the TransferFunctionArray

delay_response(dt)

Adjusts the FRF phases as if the response had been shifted dt in time

enforce_causality([method, ...])

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_time_data(reference_data, response_data)

Computes a transfer function from reference and response time histories

ifft([norm])

Converts frequency response functions to impulse response functions via an inverse fourier transform.

interpolate_by_zero_pad(irf_padded_length[, ...])

Interpolates a transfer function by zero padding or truncating its impulse response

plot([one_axis, subplots_kwargs, plot_kwargs])

Plot the transfer functions

plot_cond_num([number_retained_values, ...])

Plots the condition number of the FRF matrix

plot_singular_values([rcond, ...])

Plot the singular values of an FRF matrix with a visualization of the rcond tolerance

substructure_by_constraint_matrix(dofs, ...)

Performs frequency based substructuring using the

substructure_by_coordinate(dof_pairs)

Performs frequency based substructuring by constraining pairs of degrees of freedom

svd([full_matrices, compute_uv, as_matrix])

Compute the SVD of the provided FRF matrix

Attributes

function_type

Returns the function type of the data array

compute_cmif(part='both', tracking=None)[source]

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.

Raises

ValueError – Raised if an invalid tracking is specified

Returns

output_array – Complex Mode Indicator Function

Return type

ModeIndicatorFunctionArray

compute_mif(mif_type, *mif_args, **mif_kwargs)[source]

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

Raises

ValueError – If an invalid mif name is provided.

Returns

Mode indicator function

Return type

ModeIndicatorFunctionArray

compute_mmif(part='real', mass_matrix=None)[source]

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

Raises

ValueError – Raised if an invalid part is specified

Returns

output_array – Multi Mode Indicator Function

Return type

ModeIndicatorFunctionArray

compute_nmif(part='real')[source]

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’.

Raises

ValueError – Raised if an invalid part is specified

Returns

output_array – Normal Mode Indicator Function

Return type

ModeIndicatorFunctionArray

delay_response(dt)[source]

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 – A copy of the transfer function with the phase shifted

Return type

TransferFunctionArray

enforce_causality(method='exponential_taper', window_parameter=None, end_of_ringdown=None)[source]

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

The FRF with causality enforced.

Return type

TransferFunctionArray

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://scholarsmine.mst.edu/masters_theses/7490/

static from_time_data(reference_data: TimeHistoryArray, response_data: TimeHistoryArray, samples_per_average: Optional[int] = None, overlap: float = 0.0, method: str = 'H1', window=array([1.]), **timedata2frf_kwargs)[source]

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.

  • 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.

Raises

ValueError – Raised if reference and response functions do not have consistent abscissa

Returns

A transfer function array computed from the specified references and responses.

Return type

TransferFunctionArray

property function_type

Returns the function type of the data array

ifft(norm='backward', **scipy_irfft_kwargs)[source]

Converts frequency response functions to impulse response functions via an inverse fourier transform.

Paramters

normstr, optional

The type of normalization applied to the fft computation.

scipy_irfft_kwargs :

Additional keywords that will be passed to SciPy’s irfft function.

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.

raises 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.

returns

The impulse response function array computed from the transfer function array.

rtype

ImpulseResponseFunctionArray

interpolate_by_zero_pad(irf_padded_length, return_irf=False)[source]

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.

Returns

Transfer function array with appropriately spaced abscissa

Return type

TransferFunctionArray or ImpulseResponseFunctionArray

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(one_axis=True, subplots_kwargs={}, plot_kwargs={})[source]

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.

  • 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 {}.

Returns

axis – On which the data were plotted

Return type

matplotlib axis or array of axes

plot_cond_num(number_retained_values=None, min_frequency=None, max_frequency=None)[source]

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.

Return type

None.

plot_singular_values(rcond=None, condition_number=None, number_retained_values=None, regularization_parameter=None, min_frequency=None, max_frequency=None)[source]

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

substructure_by_constraint_matrix(dofs, constraint_matrix)[source]

Performs frequency based substructuring using the

Parameters
  • dofs (CoordinateArray) – Coordinates to use in the constraints

  • constraint_matrix (np.ndarray) – Constraints to apply to the frequency response functions

Raises

ValueError – If listed degrees of freedom are not found in the function.

Returns

constrained_frfs – Constrained Frequency Response Functions

Return type

TransferFunctionArray

substructure_by_coordinate(dof_pairs)[source]

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.

Returns

Constrained frequency response functions

Return type

TransferFunctionArray

svd(full_matrices=True, compute_uv=True, as_matrix=True)[source]

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.