dynamics_vibrations Module



Contents


Functions

public pure function damping_from_fractional_overshoot(x) result(rst)

Employs the method of fractional overshoot to estimate the damping ratio from the response of a system to a step input. This method is useful for cases where the damping ratio is between approximately 0.5 to 0.8. In such range, the logarithmic decrement approach becomes less precise.

Read more…

Arguments

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

The step response of the system.

Return Value real(kind=real64)

The estimated damping ratio.

public pure elemental function damping_from_log_decrement(delta) result(rst)

Computes the damping ratio from the logarithmic decrement . The damping ratio is related to the logarithmic decrement by the following relationship.

Read more…

Arguments

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

The logarithmic decrement.

Return Value real(kind=real64)

The damping ratio.

public pure elemental function estimate_bandwidth(fn, zeta) result(rst)

Estimates the bandwidth of the resonant mode of a vibratory system. The bandwidth is the width of the range of frequencies for which the energy is at least half its peak value and is computed as .

Arguments

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

The resonant frequency. The units are not important; however, the units of the output will be the same as the units of this parameter.

real(kind=real64), intent(in) :: zeta

The damping ratio.

Return Value real(kind=real64)

The bandwidth.

public pure elemental function evaluate_step_response(wn, zeta, xs, t) result(rst)

Evaluates the response of an underdamped single-degree-of-freedom, linear system to a step function of amplitude .

Read more…

Arguments

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

The resonant frequency, in rad/s.

real(kind=real64), intent(in) :: zeta

The damping ratio.

real(kind=real64), intent(in) :: xs

The amplitude of the step input.

real(kind=real64), intent(in) :: t

The point in time at which to evaluate the response (units = s).

Return Value real(kind=real64)

The step response.

public pure function find_settling_amplitude(x) result(rst)

Estimates the settling amplitude for a step response.

Arguments

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

The step response of the system.

Return Value real(kind=real64)

The settling amplitude of the step response.

public pure elemental function logarithmic_decrement(x1, x2, n) result(rst)

Computes the logarithmic decrement given the value of two successive peaks in the time history of the free vibratory response of the system. The logarithmic decrement is calculated as follows.

Read more…

Arguments

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

The amplitude of the first peak.

real(kind=real64), intent(in) :: x2

The amplitude of the second peak that occurs N periods after the first.

integer(kind=int32), intent(in) :: n

The number of periods of oscillation seperating the two peaks.

Return Value real(kind=real64)

The logarithmic decrement .

public pure elemental function q_factor(zeta) result(rst)

Estimates the Q-factor for a vibratory system. The Q-factor is computed .

Arguments

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

The damping ratio.

Return Value real(kind=real64)

The Q-factor.

public pure elemental function rise_time(wn, zeta) result(rst)

Computes the rise time for an underdamped, second-order system. The rise time is the time it takes for the system response to go from 0% to 100% of its final value and is given by the following relationship.

Read more…

Arguments

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

The resonant frequency of the system, in rad/s.

real(kind=real64), intent(in) :: zeta

The damping ratio of the system. This value must be less than 1 as this relationship is only valid for an underdamped system.

Return Value real(kind=real64)

The rise time, in units of seconds.


Subroutines

public subroutine find_free_response_properties(t, x, delta, fn, x1, x2, t1, t2, s, n)

Given a free-response time history, this routine attempts to find the logarithmic decrement and resonant frequency of a vibratory system. The logarithmic decrement is estimated by finding successive peaks by means of peak detection.

Arguments

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

An N-element array containing the values in time

real(kind=real64), intent(in), dimension(:) :: x

An N-element array containing the response sampled at the time points given in t.

real(kind=real64), intent(out) :: delta

The logarithmic decrement estimate. If sufficient peaks cannot be located, the routine returns NaN.

real(kind=real64), intent(out) :: fn

The damped resonant frequency in units of Hz, assuming that the time values are in seconds. If the time units are not in seconds, the units will be cycle/unit time with unit time being the units in which t is supplied. If sufficient peaks cannot be located, the routine returns NaN.

real(kind=real64), intent(out), optional :: x1

An optional parameter that, if provided, allows for the routine to return the amplitude of the first peak. If sufficient peaks cannot be located, the routine returns NaN.

real(kind=real64), intent(out), optional :: x2

An optional parameter that, if provided, allows for the routine to return the amplitude of the second peak. If sufficient peaks cannot be located, the routine returns NaN.

real(kind=real64), intent(out), optional :: t1

An optional parameter that, if provided, allows for the routine to return the time at which the first peak was located. If sufficient peaks cannot be located, the routine returns NaN.

real(kind=real64), intent(out), optional :: t2

An optional parameter that, if provided, allows for the routine to return the time at which the second peak was located. If sufficient peaks cannot be located, the routine returns NaN.

real(kind=real64), intent(in), optional :: s

An optional input that, if provided, allows for control of the sensitivity of the peak detection algorithm. The default is 0.1% of the peak-peak amplitude of the signal.

integer(kind=int32), intent(in), optional :: n

An optional input that, if provided, determines the number of periods to allow between peak selection for the logarithmic decrement calculation. The default is 1.