| Type | Visibility | Attributes | Name | Initial | |||
|---|---|---|---|---|---|---|---|
| integer(kind=int32), | public, | parameter | :: | FRF_ACCELERANCE_MODEL | = | 1 |
Defines an accelerance frequency response model. |
| integer(kind=int32), | public, | parameter | :: | FRF_RECEPTANCE_MODEL | = | 2 |
Defines a receptance frequency response model. |
Evaluates the specified accelerance FRF model. The model is of the following form.
| Type | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| real(kind=real64), | intent(in), | dimension(:) | :: | mdl |
The model parameter array. The elements of the array are stored as . |
|
| real(kind=real64), | intent(in) | :: | w |
The frequency value, in rad/s, at which to evaluate the model. |
The resulting frequency response function.
Evaluates the specified accelerance FRF model. The model is of the following form.
| Type | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| real(kind=real64), | intent(in), | dimension(:) | :: | mdl |
The model parameter array. The elements of the array are stored as . |
|
| real(kind=real64), | intent(in), | dimension(:) | :: | w |
The frequency value, in rad/s, at which to evaluate the model. |
The resulting frequency response function.
Evaluates the specified receptance FRF model. The model is of the following form.
| Type | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| real(kind=real64), | intent(in), | dimension(:) | :: | mdl |
The model parameter array. The elements of the array are stored as . |
|
| real(kind=real64), | intent(in) | :: | w |
The frequency value, in rad/s, at which to evaluate the model. |
The resulting frequency response function.
Evaluates the specified receptance FRF model. The model is of the following form.
| Type | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| real(kind=real64), | intent(in), | dimension(:) | :: | mdl |
The model parameter array. The elements of the array are stored as . |
|
| real(kind=real64), | intent(in), | dimension(:) | :: | w |
The frequency value, in rad/s, at which to evaluate the model. |
The resulting frequency response function.
Computes the frequency response functions for a system of ODE's.
Computes the frequency response functions for a multi-degree-of-freedom system that uses proportional damping such that the damping matrix is related to the stiffness an mass matrices by proportional damping coefficients and by .
| Type | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| real(kind=real64), | intent(in), | dimension(:,:) | :: | mass |
The N-by-N mass matrix for the system. This matrix must be symmetric. |
|
| real(kind=real64), | intent(in), | dimension(:,:) | :: | stiff |
The N-by-N stiffness matrix for the system. This matrix must be symmetric. |
|
| real(kind=real64), | intent(in) | :: | alpha |
The mass damping factor, . |
||
| real(kind=real64), | intent(in) | :: | beta |
The stiffness damping factor, . |
||
| real(kind=real64), | intent(in), | dimension(:) | :: | freq |
An M-element array of frequency values at which to evaluate the frequency response functions, in units of rad/s. |
|
| procedure(modal_excite), | intent(in), | pointer | :: | frc |
A pointer to a routine used to compute the modal forcing function. |
|
| real(kind=real64), | intent(out), | optional, | allocatable, dimension(:) | :: | modes |
An optional N-element allocatable array that, if supplied, will be used to retrieve the modal frequencies, in units of rad/s. |
| real(kind=real64), | intent(out), | optional, | allocatable, dimension(:,:) | :: | modeshapes |
An optional N-by-N allocatable matrix that, if supplied, will be used to retrieve the N mode shapes with each vector occupying its own column. |
| class(*), | intent(inout), | optional | :: | args |
An optional argument that can be used to communicate with the outside world. |
|
| class(errors), | intent(inout), | optional, | target | :: | err |
An optional errors-based object that if provided can be used to retrieve information relating to any errors encountered during execution. If not provided, a default implementation of the errors class is used internally to provide error handling. Possible errors and warning messages that may be encountered are as follows.
|
The resulting frequency responses.
Computes the frequency response functions for a multi-degree-of-freedom system that uses proportional damping such that the damping matrix is related to the stiffness an mass matrices by proportional damping coefficients and by .
| Type | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| real(kind=real64), | intent(in), | dimension(:,:) | :: | mass |
The N-by-N mass matrix for the system. This matrix must be symmetric. |
|
| real(kind=real64), | intent(in), | dimension(:,:) | :: | stiff |
The N-by-N stiffness matrix for the system. This matrix must be symmetric. |
|
| real(kind=real64), | intent(in) | :: | alpha |
The mass damping factor, . |
||
| real(kind=real64), | intent(in) | :: | beta |
The stiffness damping factor, . |
||
| integer(kind=int32), | intent(in) | :: | nfreq |
The number of frequency values to analyze. This value must be at least 2. |
||
| real(kind=real64), | intent(in) | :: | freq1 |
The starting frequency, in units of rad/s. |
||
| real(kind=real64), | intent(in) | :: | freq2 |
The ending frequency, in units of rad/s. |
||
| procedure(modal_excite), | intent(in), | pointer | :: | frc |
A pointer to a routine used to compute the modal forcing function. |
|
| real(kind=real64), | intent(out), | optional, | allocatable, dimension(:) | :: | modes |
An optional N-element allocatable array that, if supplied, will be used to retrieve the modal frequencies, in units of rad/s. |
| real(kind=real64), | intent(out), | optional, | allocatable, dimension(:,:) | :: | modeshapes |
An optional N-by-N allocatable matrix that, if supplied, will be used to retrieve the N mode shapes with each vector occupying its own column. |
| class(*), | intent(inout), | optional | :: | args |
An optional argument that can be used to communicate with the outside world. |
|
| class(errors), | intent(inout), | optional, | target | :: | err |
An optional errors-based object that if provided can be used to retrieve information relating to any errors encountered during execution. If not provided, a default implementation of the errors class is used internally to provide error handling. Possible errors and warning messages that may be encountered are as follows.
|
The resulting frequency responses.
Estimates the frequency response of a single-input, single-output (SISO) system.
| Type | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| real(kind=real64), | intent(in), | dimension(:) | :: | x |
An N-element array containing the excitation signal. |
|
| real(kind=real64), | intent(in), | dimension(:) | :: | y |
An N-element array containing the response signal. |
|
| real(kind=real64), | intent(in) | :: | fs |
The sampling frequency, in Hz. |
||
| class(window), | intent(in), | optional, | target | :: | win |
The window to apply to the data. If nothing is supplied, no window is applied. |
| integer(kind=int32), | intent(in), | optional | :: | method |
Enter 1 to utilize an H1 estimator; else, enter 2 to utilize an H2 estimator. The default is an H1 estimator. An H1 estimator is defined as the cross-spectrum of the input and response signals divided by the energy spectral density of the input. An H2 estimator is defined as the energy spectral density of the response divided by the cross-spectrum of the input and response signals.
|
|
| class(errors), | intent(inout), | optional, | target | :: | err |
An optional errors-based object that if provided can be used to retrieve information relating to any errors encountered during execution. If not provided, a default implementation of the errors class is used internally to provide error handling. Possible errors and warning messages that may be encountered are as follows.
|
The resulting frequency response function.
Estimates the frequency responses of a multiple-input, multiple-output (MIMO) system.
| Type | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| real(kind=real64), | intent(in), | dimension(:,:) | :: | x |
An N-by-P array containing the P inputs to the system. |
|
| real(kind=real64), | intent(in), | dimension(:,:) | :: | y |
An N-by-M array containing the M outputs from the system. |
|
| real(kind=real64), | intent(in) | :: | fs |
The sampling frequency, in Hz. |
||
| class(window), | intent(in), | optional, | target | :: | win |
The window to apply to the data. If nothing is supplied, no window is applied. |
| integer(kind=int32), | intent(in), | optional | :: | method |
Enter 1 to utilize an H1 estimator; else, enter 2 to utilize an H2 estimator. The default is an H1 estimator. An H1 estimator is defined as the cross-spectrum of the input and response signals divided by the energy spectral density of the input. An H2 estimator is defined as the energy spectral density of the response divided by the cross-spectrum of the input and response signals.
|
|
| class(errors), | intent(inout), | optional, | target | :: | err |
An optional errors-based object that if provided can be used to retrieve information relating to any errors encountered during execution. If not provided, a default implementation of the errors class is used internally to provide error handling. Possible errors and warning messages that may be encountered are as follows.
|
The resulting frequency response functions.
Computes the frequency response of each equation of a system of harmonically excited ODE's by sweeping through frequency.
| Type | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| procedure(harmonic_ode), | intent(in), | pointer | :: | fcn |
A pointer to the routine containing the ODE's to integrate. |
|
| real(kind=real64), | intent(in), | dimension(:) | :: | freq |
An M-element array containing the frequency points at which the solution should be computed. Notice, whatever units are utilized for this array are also the units of the excitation_frequency property in @p sys. It is recommended that the units be set to Hz. Additionally, this array cannot contain any zero-valued elements as the ODE solution time for each frequency is determined by the period of oscillation and number of cycles. |
|
| real(kind=real64), | intent(in), | dimension(:) | :: | iv |
An N-element array containing the initial conditions for each of the N ODEs. |
|
| class(ode_integrator), | intent(inout), | optional, | target | :: | solver |
An optional differential equation solver. The default solver is the Dormand-Prince Runge-Kutta integrator from the DIFFEQ library. |
| integer(kind=int32), | intent(in), | optional | :: | ncycles |
An optional parameter controlling the number of cycles to analyze when determining the amplitude and phase of the response. The default is 20. |
|
| integer(kind=int32), | intent(in), | optional | :: | ntransient |
An optional parameter controlling how many of the initial "transient" cycles to ignore. The default is 200. |
|
| integer(kind=int32), | intent(in), | optional | :: | points |
An optional parameter controlling how many evenly spaced solution points should be considered per cycle. The default is 1000. Notice, there must be at least 2 points per cycle for the analysis to be effective. The algorithm utilizes a discrete Fourier transform to determine the phase and amplitude, and in order to satisfy Nyquist conditions, the value must be at least 2. |
|
| class(*), | intent(inout), | optional | :: | args |
An optional argument allowing for passing of data in/out of the fcn subroutine. |
|
| class(errors), | intent(inout), | optional, | target | :: | err |
An optional errors-based object that if provided can be used to retrieve information relating to any errors encountered during execution. If not provided, a default implementation of the errors class is used internally to provide error handling. Possible errors and warning messages that may be encountered are as follows.
|
The resulting frequency responses.
Computes the frequency response of each equation of a system of harmonically excited ODE's by sweeping through frequency.
| Type | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| procedure(harmonic_ode), | intent(in), | pointer | :: | fcn |
A pointer to the routine containing the ODE's to integrate. |
|
| integer(kind=int32), | intent(in) | :: | nfreq |
The number of frequency values to analyze. This value must be at least 2. |
||
| real(kind=real64), | intent(in) | :: | freq1 |
The starting frequency. It is recommended that the units be set to Hz. |
||
| real(kind=real64), | intent(in) | :: | freq2 |
The ending frequency. It is recommended that the units be set to Hz. |
||
| real(kind=real64), | intent(in), | dimension(:) | :: | iv |
An N-element array containing the initial conditions for each of the N ODEs. |
|
| class(ode_integrator), | intent(inout), | optional, | target | :: | solver |
An optional differential equation solver. The default solver is the Dormand-Prince Runge-Kutta integrator from the DIFFEQ library. |
| integer(kind=int32), | intent(in), | optional | :: | ncycles |
An optional parameter controlling the number of cycles to analyze when determining the amplitude and phase of the response. The default is 20. |
|
| integer(kind=int32), | intent(in), | optional | :: | ntransient |
An optional parameter controlling how many of the initial "transient" cycles to ignore. The default is 200. |
|
| integer(kind=int32), | intent(in), | optional | :: | points |
An optional parameter controlling how many evenly spaced solution points should be considered per cycle. The default is 1000. Notice, there must be at least 2 points per cycle for the analysis to be effective. The algorithm utilizes a discrete Fourier transform to determine the phase and amplitude, and in order to satisfy Nyquist conditions, the value must be at least 2. |
|
| class(*), | intent(inout), | optional | :: | args |
An optional argument allowing for passing of data in/out of the fcn subroutine. |
|
| class(errors), | intent(inout), | optional, | target | :: | err |
An optional errors-based object that if provided can be used to retrieve information relating to any errors encountered during execution. If not provided, a default implementation of the errors class is used internally to provide error handling. Possible errors and warning messages that may be encountered are as follows.
|
The resulting frequency responses.
Defines a system of ODE's exposed to harmonic excitation.
| Type | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| real(kind=real64), | intent(in) | :: | freq |
The excitation frequency. |
||
| real(kind=real64), | intent(in) | :: | t |
The current time step value. |
||
| real(kind=real64), | intent(in), | dimension(:) | :: | x |
The value of the solution estimate at time t. |
|
| real(kind=real64), | intent(out), | dimension(:) | :: | dxdt |
The derivatives as computed by this routine. |
|
| class(*), | intent(inout), | optional | :: | args |
An optional argument allowing the passing of data in/out of this routine. |
Defines the interface to a routine for defining the forcing function for a modal frequency analysis.
| Type | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| real(kind=real64), | intent(in) | :: | freq |
The excitation frequency. When used as a part of a frequency response calculation, this value will have the same units as the frequency values provided to the frequency response routine. |
||
| complex(kind=real64), | intent(out), | dimension(:) | :: | frc |
An N-element array where the forcing function should be written. |
|
| class(*), | intent(inout), | optional | :: | args |
An optional argument that can be used to communicate with the outside world. |
Defines the interface for a ODE excitation function.
| Type | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| real(kind=real64), | intent(in) | :: | t |
The value of the independent variable at which to evaluate the excitation function. |
The result.
A container for a frequency response function, or series of frequency response functions.
| Type | Visibility | Attributes | Name | Initial | |||
|---|---|---|---|---|---|---|---|
| real(kind=real64), | public, | allocatable, dimension(:) | :: | frequency |
An N-element array containing the frequency values at which the FRF is provided. The units of this array are the same as the units of the frequency values passed to the routine used to compute the frequency response. |
||
| complex(kind=real64), | public, | allocatable, dimension(:,:) | :: | responses |
An N-by-M matrix containing the M frequency response functions evaluated at each of the N frequency points. |
A container for the frequency responses of a system of multiple inputs and multiple outputs (MIMO).
| Type | Visibility | Attributes | Name | Initial | |||
|---|---|---|---|---|---|---|---|
| real(kind=real64), | public, | allocatable, dimension(:) | :: | frequency |
An N-element array containing the frequency values at which the FRF is provided. The units of this array are the same as the units of the frequency values passed to the routine used to compute the frequency response. |
||
| complex(kind=real64), | public, | allocatable, dimension(:,:,:) | :: | responses |
An N-by-M-by-P array containing the frequency response functions for each of the M outputs corresponding to each of the P inputs. |
Evaluates a linear chirp function.
| Type | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| real(kind=real64), | intent(in) | :: | t |
The value of the independent variable at which to evaluate the chirp. |
||
| real(kind=real64), | intent(in) | :: | amp |
The amplitude. |
||
| real(kind=real64), | intent(in) | :: | span |
The duration of the time it takes to sweep from the start frequency to the end frequency. |
||
| real(kind=real64), | intent(in) | :: | f1Hz |
The lower excitation frequency, in Hz. |
||
| real(kind=real64), | intent(in) | :: | f2Hz |
The upper excitation frequency, in Hz. |
The value of the function at t.
Computes the modal damping factors given the proportional damping terms and where , , and is the eigenvalue of the system.
| Type | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| real(kind=real64), | intent(in) | :: | lambda |
The square of the modal frequency - the eigen value. |
||
| real(kind=real64), | intent(in) | :: | alpha |
The mass damping factor, . |
||
| real(kind=real64), | intent(in) | :: | beta |
The stiffness damping factor, . |
The modal damping parameter.
Fits an experimentally obtained frequency response by model for either a receptance model:
| Type | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| integer(kind=int32), | intent(in) | :: | mt |
The excitation method. The options are as follows. |
||
| integer(kind=int32), | intent(in) | :: | n |
The model order (# of resonant modes). |
||
| real(kind=real64), | intent(in), | dimension(:) | :: | freq |
An M-element array containing the excitation frequency values in units of rad/s. |
|
| complex(kind=real64), | intent(in), | dimension(:) | :: | rsp |
An M-element array containing the frequency response to fit. |
|
| real(kind=real64), | intent(in), | optional, | dimension(:) | :: | maxp |
An optional 3*N-element array that can be used as upper limits on the parameter values. If no upper limit is requested for a particular parameter, utilize a very large value. The internal default is to utilize huge() as a value. |
| real(kind=real64), | intent(in), | optional, | dimension(:) | :: | minp |
An optional 3*N-element array that can be used as lower limits on the parameter values. If no lower limit is requested for a particalar parameter, utilize a very large magnitude, but negative, value. The internal default is to utilize -huge() as a value. |
| real(kind=real64), | intent(in), | optional, | dimension(:) | :: | init |
An optional 3N-element array that, if supplied, provides an initial guess for each of the 3N model parameters for the iterative solver. If supplied, this array replaces the peak finding algorithm for estimating an initial guess. |
| type(regression_statistics), | intent(out), | optional, | dimension(:) | :: | stats |
An optional 3*N-element array that, if supplied, will be used to return statistics about the fit for each model parameter. |
| real(kind=real64), | intent(in), | optional | :: | alpha |
The significance level at which to evaluate the confidence intervals. The default value is 0.05 such that a 95% confidence interval is calculated. |
|
| type(iteration_controls), | intent(in), | optional | :: | controls |
An optional input providing custom iteration controls. |
|
| type(lm_solver_options), | intent(in), | optional | :: | settings |
An optional input providing custom settings for the solver. |
|
| type(convergence_info), | intent(out), | optional | :: | info |
An optional output that can be used to gain information about the iterative solution and the nature of the convergence. |
|
| class(errors), | intent(inout), | optional, | target | :: | err |
An optional errors-based object that if provided can be used to retrieve information relating to any errors encountered during execution. If not provided, a default implementation of the errors class is used internally to provide error handling. Possible errors and warning messages that may be encountered are as follows. |
An array containing the model parameters stored as .
Computes the modal frequencies and modes shapes for multi-degree-of-freedom system.
| Type | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| real(kind=real64), | intent(in), | dimension(:,:) | :: | mass |
The N-by-N mass matrix for the system. This matrix must be symmetric. |
|
| real(kind=real64), | intent(in), | dimension(:,:) | :: | stiff |
The N-by-N stiffness matrix for the system. This matrix must be symmetric. |
|
| real(kind=real64), | intent(out), | allocatable, dimension(:) | :: | freqs |
An allocatable N-element array where the modal frequencies will be returned in ascending order with units of rad/s. |
|
| real(kind=real64), | intent(out), | optional, | allocatable, dimension(:,:) | :: | modeshapes |
An optional, allocatable N-by-N matrix where the N mode shapes for the system will be returned. The mode shapes are stored in columns. |
| class(errors), | intent(inout), | optional, | target | :: | err |
Normalizes mode shape vectors such that the largest magnitude value in the vector is one.
| Type | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| real(kind=real64), | intent(inout), | dimension(:,:) | :: | x |
The matrix of mode shape vectors with one vector per column. |