Random Vibration Theory Parameters

Definition of the custom type RandomVibrationParameters to represent the model components/approaches used for the random vibration theory calculations. In particular, the type stores symbols to define the:

  • excitation_duration model to use,
  • rms_duration model to use, and
  • peak_factor model to use.
StochasticGroundMotionSimulation.RandomVibrationParametersType
RandomVibrationParameters

Struct holding parameters/methods for Random Vibration Theory.

  • pf_method is the method used for peak factor computation
- `:DK80` (default) is Der Kiureghian (1980), building on Vanmarcke (1975)
- `:CL56` is Cartwright Longuet-Higgins (1956)
  • dur_ex is the model for excitation duration
- `:BT14` (default) is the Boore & Thompson (2014) model - note that this is adpated to work with `r_ps`
  • dur_rms is the model for rms duration
- `:BT12` is the Boore & Thompson (2012) model
- `:BT15` (default) is the Boore & Thompson (2015) model
  • dur_region is the region specified for the duration model
- `:WNA` (default) is western North America
- `:ENA` is eastern North America
source

Note that the default specification is:

RandomVibrationParameters() = RandomVibrationParameters(:DK80, :BT14, :BT15, :WNA)
RandomVibrationParameters (generic function with 1 method)

However, an alternative constructor exists that takes a pf_method as an argument. For this constructor, the rms_duration model is linked to the peak factor method:

  • DK80 is paired with :BT15 as a default
  • CL56 is paired with :BT12 as a default

As these are currently the only two rms_duration models implemented, the constructor is specified as:

RandomVibrationParameters(pf) = RandomVibrationParameters(pf, :BT14, ((pf == :DK80) ? :BT15 : :BT12), :WNA)
RandomVibrationParameters (generic function with 1 method)

In all cases, the Boore & Thompson (2014) excitation duration model is employed as that is the only model currently implemented.

## Functionality

The overall goal of these random vibration methods is to compute:

\[S_a = \psi \sqrt{ \frac{m_0}{D_{rms}}}\]

where $\psi$ is the peak factor computed from peak_factor, $m_0$ is the zeroth order spectral moment computed from spectral_moment, and $D_{rms}$ is the root-mean-square duration computed from rms_duration.

The main methods used to interact with RandomVibrationParameters are:

StochasticGroundMotionSimulation.spectral_momentFunction
spectral_moment(order::Int, m::S, r_ps::T, fas::FourierParameters, sdof::Oscillator; nodes::Int=31, control_freqs::Vector{Float64}=[1e-3, 1e-1, 1.0, 10.0, 100.0, 300.0] ) where {S<:Real,T<:Real}

Compute spectral moment of a specified order.

Evaluates the expression:

\[ m_k = 2\int_{0}^{\infty} \left(2\pi f\right)^k |H(f;f_n,\zeta_n)|^2 |A(f)|^2 df\]

where $k$ is the order of the moment.

Integration is performed using Gauss-Legendre integration using nodes nodes and weights. The integration domain is partitioned over the control_freqs as well as two inserted frequencies at f_n/1.5 and f_n*1.5 in order to ensure good approximation of the integral around the sdof resonant frequency.

See also: spectral_moments

source
StochasticGroundMotionSimulation.spectral_momentsFunction
spectral_moments(order::Vector{Int}, m::S, r_ps::T, fas::FourierParameters, sdof::Oscillator; nodes::Int=31, control_freqs::Vector{Float64}=[1e-3, 1e-1, 1.0, 10.0, 100.0, 300.0] ) where {S<:Real,T<:Real}

Compute a vector of spectral moments for the specified order.

Evaluates the expression:

\[ m_k = 2\int_{0}^{\infty} \left(2\pi f\right)^k |H(f;f_n,\zeta_n)|^2 |A(f)|^2 df\]

for each order, where $k$ is the order of the moment.

Integration is performed using Gauss-Legendre integration using nodes nodes and weights. The integration domain is partitioned over the control_freqs as well as two inserted frequencies at f_n/1.5 and f_n*1.5 in order to ensure good approximation of the integral around the sdof resonant frequency.

See also: spectral_moment, spectral_moments_gk

source
StochasticGroundMotionSimulation.excitation_durationFunction

excitationduration(m, rps::U, src::SourceParameters{S,T}, rvt::RandomVibrationParameters) where {S<:Float64,T<:Real,U<:Real}

Generic function implementing excitation duration models.

Currently, only the Boore & Thompson (2014) model is implemented.

source
StochasticGroundMotionSimulation.rms_durationFunction
rms_duration(m::S, r_ps::T, src::SourceParameters, sdof::Oscillator, rvt::RandomVibrationParameters) where {S<:Float64,T<:Real}

Returns a 3-tuple of (Drms, Dex, Dratio), using a switch on rvt.dur_rms. Default :BT12 makes use of the :BT14 model for excitation duration, Dex.

  • m is magnitude
  • r_ps is an equivalent point source distance
source
StochasticGroundMotionSimulation.peak_factorFunction
peak_factor(m::S, r_ps::T, fas::FourierParameters, sdof::Oscillator; pf_method::Symbol=:DK80) where {S<:Real,T<:Real}

Peak factor $u_{max} / u_{rms}$ with a switch of pf_method to determine the approach adopted. pf_method can currently be one of: - :CL56 for Cartright Longuet-Higgins (1956) - :DK80 for Der Kiureghian (1980), building on Vanmarcke (1975)

Defaults to :DK80.

source
peak_factor(m::S, r_ps::T, Dex::U, m0::V, fas::FourierParameters, sdof::Oscillator, rvt::RandomVibrationParameters) where {S<:Real,T<:Real,U<:Real,V<:Real}

Peak factor umax / urms with a switch of pf_method to determine the approach adopted. pf_method can currently be one of: - :CL56 for Cartright Longuet-Higgins (1956) - :DK80 for Der Kiureghian (1980), building on Vanmarcke (1975)

Defaults to :DK80.

source
StochasticGroundMotionSimulation.rvt_response_spectral_ordinateFunction
rvt_response_spectral_ordinate(m::S, r_ps::T, fas::FourierParameters, sdof::Oscillator, rvt::RandomVibrationParameters) where {S<:Real,T<:Real}

Response spectral ordinate (units of $g$) for the specified scenario.

The spectral ordinate is computed using the expression:

\[S_a = \psi \sqrt{\frac{m_0}{D_{rms}}}\]

where $\psi$ is the peak factor computed from peak_factor, $m_0$ is the zeroth order spectral moment from spectral_moment, and $D_{rms}$ is the RMS duration computed from rms_duration.

See also: rvt_response_spectral_ordinate, rvt_response_spectrum, rvt_response_spectrum!

source
rvt_response_spectral_ordinate(period::U, m::S, r_ps::T, fas::FourierParameters, rvt::RandomVibrationParameters) where {S<:Real,T<:Real,U<:Float64}

Response spectral ordinate (units of $g$) for the specified scenario.

The spectral ordinate is computed using the expression:

\[S_a = \psi \sqrt{\frac{m_0}{D_{rms}}}\]

where $\psi$ is the peak factor computed from peak_factor, $m_0$ is the zeroth order spectral moment from spectral_moment, and $D_{rms}$ is the RMS duration computed from rms_duration.

See also: rvt_response_spectral_ordinate, rvt_response_spectrum, rvt_response_spectrum!

source
StochasticGroundMotionSimulation.rvt_response_spectrumFunction
rvt_response_spectrum(period::Vector{U}, m::S, r_ps::T, fas::FourierParameters, rvt::RandomVibrationParameters) where {S<:Real,T<:Real,U<:Float64}

Response spectrum (units of $g$) for the vector of periods period and the specified scenario.

Each spectral ordinate is computed using the expression:

\[S_a = \psi \sqrt{\frac{m_0}{D_{rms}}}\]

where $\psi$ is the peak factor computed from peak_factor, $m_0$ is the zeroth order spectral moment from spectral_moment, and $D_{rms}$ is the RMS duration computed from rms_duration. The various terms are all functions of the oscillator period.

See also: rvt_response_spectral_ordinate, rvt_response_spectrum!

source