dynamics_frequency_response Module


Uses


Contents


Variables

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.


Interfaces

  • private pure function evaluate_accelerance_frf_model_scalar(mdl, w) result(rst)

    Evaluates the specified accelerance FRF model. The model is of the following form.

    Arguments

    Type IntentOptional 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.

    Return Value complex(kind=real64)

    The resulting frequency response function.

  • private pure function evaluate_accelerance_frf_model_array(mdl, w) result(rst)

    Evaluates the specified accelerance FRF model. The model is of the following form.

    Arguments

    Type IntentOptional 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.

    Return Value complex(kind=real64), allocatable, dimension(:)

    The resulting frequency response function.

  • private pure function evaluate_receptance_frf_model_scalar(mdl, w) result(rst)

    Evaluates the specified receptance FRF model. The model is of the following form.

    Arguments

    Type IntentOptional 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.

    Return Value complex(kind=real64)

    The resulting frequency response function.

  • private pure function evaluate_receptance_frf_model_array(mdl, w) result(rst)

    Evaluates the specified receptance FRF model. The model is of the following form.

    Arguments

    Type IntentOptional 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.

    Return Value complex(kind=real64), allocatable, dimension(:)

    The resulting frequency response function.

public interface frequency_response

Computes the frequency response functions for a system of ODE's.

  • private function frf_modal_prop_damp(mass, stiff, alpha, beta, freq, frc, modes, modeshapes, args, err) result(rst)

    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 .

    Arguments

    Type IntentOptional 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.

    • DYN_MEMORY_ERROR: Occurs if there are issues allocating memory.
    • DYN_MATRIX_SIZE_ERROR: Occurs if the mass or stiffness matrices are not square, or if the mass and stiffness matrices are different sized.
    • DYN_NULL_POINTER_ERROR: Occurs if the forcing function pointer is undefined.

    Return Value type(frf)

    The resulting frequency responses.

  • private function frf_modal_prop_damp_2(mass, stiff, alpha, beta, nfreq, freq1, freq2, frc, modes, modeshapes, args, err) result(rst)

    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 .

    Arguments

    Type IntentOptional 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.

    • DYN_MEMORY_ERROR: Occurs if there are issues allocating memory.
    • DYN_MATRIX_SIZE_ERROR: Occurs if the mass or stiffness matrices are not square, or if the mass and stiffness matrices are different sized.
    • DYN_NULL_POINTER_ERROR: Occurs if the forcing function pointer is undefined.

    Return Value type(frf)

    The resulting frequency responses.

  • private function siso_freqres(x, y, fs, win, method, err) result(rst)

    Estimates the frequency response of a single-input, single-output (SISO) system.

    Arguments

    Type IntentOptional 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.

    • DYN_MEMORY_ERROR: Occurs if there are issues allocating memory.

    • DYN_ARRAY_SIZE_ERROR: Occurs if x and y are not the same size.

    Return Value type(frf)

    The resulting frequency response function.

  • private function mimo_freqres(x, y, fs, win, method, err) result(rst)

    Estimates the frequency responses of a multiple-input, multiple-output (MIMO) system.

    Arguments

    Type IntentOptional 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.

    • DYN_MEMORY_ERROR: Occurs if there are issues allocating memory.

    • DYN_ARRAY_SIZE_ERROR: Occurs if x and y do not have the same number of rows.

    Return Value type(mimo_frf)

    The resulting frequency response functions.

public interface frequency_sweep

  • private function frf_sweep_1(fcn, freq, iv, solver, ncycles, ntransient, points, args, err) result(rst)

    Computes the frequency response of each equation of a system of harmonically excited ODE's by sweeping through frequency.

    Arguments

    Type IntentOptional 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.

    • DYN_MEMORY_ERROR: Occurs if there are issues allocating memory.
    • DYN_NULL_POINTER_ERROR: Occurs if a null pointer is supplied.
    • DYN_INVALID_INPUT_ERROR: Occurs if an invalid parameter is given.
    • DYN_ZERO_VALUED_FREQUENCY_ERROR: Occurs if a zero-valued frequency was supplied.

    Return Value type(frf)

    The resulting frequency responses.

  • private function frf_sweep_2(fcn, nfreq, freq1, freq2, iv, solver, ncycles, ntransient, points, args, err) result(rst)

    Computes the frequency response of each equation of a system of harmonically excited ODE's by sweeping through frequency.

    Arguments

    Type IntentOptional 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.

    • DYN_MEMORY_ERROR: Occurs if there are issues allocating memory.
    • DYN_NULL_POINTER_ERROR: Occurs if a null pointer is supplied.
    • DYN_INVALID_INPUT_ERROR: Occurs if an invalid parameter is given.
    • DYN_ZERO_VALUED_FREQUENCY_ERROR: Occurs if a zero-valued frequency was supplied.

    Return Value type(frf)

    The resulting frequency responses.

interface

  • public subroutine harmonic_ode(freq, t, x, dxdt, args)

    Defines a system of ODE's exposed to harmonic excitation.

    Arguments

    Type IntentOptional 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.

interface

  • public subroutine modal_excite(freq, frc, args)

    Defines the interface to a routine for defining the forcing function for a modal frequency analysis.

    Arguments

    Type IntentOptional 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.

interface

  • public function ode_excite(t) result(rst)

    Defines the interface for a ODE excitation function.

    Arguments

    Type IntentOptional Attributes Name
    real(kind=real64), intent(in) :: t

    The value of the independent variable at which to evaluate the excitation function.

    Return Value real(kind=real64)

    The result.


Derived Types

type, public ::  frf

A container for a frequency response function, or series of frequency response functions.

Components

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.

type, public ::  mimo_frf

A container for the frequency responses of a system of multiple inputs and multiple outputs (MIMO).

Components

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.


Functions

public pure elemental function chirp(t, amp, span, f1Hz, f2Hz) result(rst)

Evaluates a linear chirp function.

Arguments

Type IntentOptional 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.

Return Value real(kind=real64)

The value of the function at t.

public pure elemental function compute_modal_damping(lambda, alpha, beta) result(rst)

Computes the modal damping factors given the proportional damping terms and where , , and is the eigenvalue of the system.

Arguments

Type IntentOptional 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, .

Return Value real(kind=real64)

The modal damping parameter.

public function fit_frf(mt, n, freq, rsp, maxp, minp, init, stats, alpha, controls, settings, info, err) result(rst)

Fits an experimentally obtained frequency response by model for either a receptance model:

Read more…

Arguments

Type IntentOptional Attributes Name
integer(kind=int32), intent(in) :: mt

The excitation method. The options are as follows.

Read more…
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.

Read more…

Return Value real(kind=real64), allocatable, dimension(:)

An array containing the model parameters stored as .


Subroutines

public subroutine modal_response(mass, stiff, freqs, modeshapes, err)

Computes the modal frequencies and modes shapes for multi-degree-of-freedom system.

Arguments

Type IntentOptional 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

public subroutine normalize_mode_shapes(x)

Normalizes mode shape vectors such that the largest magnitude value in the vector is one.

Arguments

Type IntentOptional Attributes Name
real(kind=real64), intent(inout), dimension(:,:) :: x

The matrix of mode shape vectors with one vector per column.