Fourier Source Spectrum
The source spectrum can be computed through interaction with the SourceParameters
type, or the higher level FourierParameters
type that holds a SourceParameters
instance as a property.
The source spectrum $E(f; \bm{\theta}E)$ is most commonly written in terms of:
\[ E(f; \bm{\theta}_E) = \mathcal{C} M_0 E_s(f; \bm{\theta}_E)\]
where $\mathcal{C}$ is a constant term, to be defined shortly, $M_0$ is the seismic moment, and $E_s(f; \bm{\theta}_E)$ is the source spectral shape.
The most commonly adopted source spectral shape is the $\omega^2$ model that has the form:
\[ E_s(f) = \frac{1}{1 + \left(\frac{f}{f_c}\right)^2}\]
Functionality
StochasticGroundMotionSimulation.fourier_constant
— Functionfourier_constant(src::SourceParameters)
Define the constant source term for the Fourier Amplitude Spectrum. The value provided corresponds to Fourier displacement units of cm-s. Constant set to permit distances to be passed in km, densities in t/m^3, and velocities in km/s. The reference distance is set to 1.0 km (and interpreted to be a rupture distance).
StochasticGroundMotionSimulation.fourier_source_shape
— Functionfourier_source_shape(f::Float64, fa::T, fb::T, ε::T, src::SourceParameters) where T<:Real
Fourier amplitude spectral shape for displacement defined by corner frequencies.
See also: fourier_source_shape
fourier_source_shape(f::Float64, fa::T, fb::T, ε::T, fas::FourierParameters) where T<:Real
Fourier amplitude spectral shape for displacement defined by corner frequencies.
See also: fourier_source_shape
fourier_source_shape(f::T, m::S, src::SourceParameters) where {S<:Real,T<:Float64}
Source shape of the Fourier Amplitude Spectrum of displacement, without the constant term or seismic moment. This simply includes the source spectral shape.
The nature of the source spectral shape depends upon src.model
:
:Brune
gives the single corner omega-squared spectrum (this is the default):Atkinson_Silva_2000
gives the double corner spectrum of Atkinson & Silva (2000):Beresnev_2019
gives a single-corner spectrum with arbitrary fall off rate related tosrc.n
from Beresnev (2019)
fourier_source_shape(f, m, fas::FourierParameters)
Source shape of the Fourier amplitude spectrum of displacement, without the constant term or seismic moment. Defined using a FourierParameters
instance for the source model.
See also: fourier_source_shape
fourier_source_shape(f::Vector{S}, fa::T, fb::T, ε::T, src::SourceParameters) where {S<:Real,T<:Real}
Fourier amplitude spectral shape for displacement defined by corner frequencies.
See also: fourier_source_shape
fourier_source_shape(f::Vector{T}, m::S, src::SourceParameters) where {S<:Real,T<:Float64}
Source shape of the Fourier Amplitude Spectrum of displacement, without the constant term or seismic moment. This simply includes the source spectral shape.
The nature of the source spectral shape depends upon src.model
:
:Brune
gives the single corner omega-squared spectrum (this is the default):Atkinson_Silva_2000
gives the double corner spectrum of Atkinson & Silva (2000):Beresnev_2019
gives a single-corner spectrum with arbitrary fall off rate related tosrc.n
from Beresnev (2019)
StochasticGroundMotionSimulation.fourier_source
— Functionfourier_source(f::T, m::S, src::SourceParameters) where {S<:Real,T<:Float64}
Source Fourier Amplitude Spectrum of displacement, without the constant term. This simply includes the seismic moment and the source spectral shape.
See also: fourier_source_shape
fourier_source(f, m, fas::FourierParameters)
Source Fourier Amplitude Spectrum of displacement, without the constant term. This simply includes the seismic moment and the source spectral shape. Defined using a FourierParameters
instance for the source model.
See also: fourier_source_shape
StochasticGroundMotionSimulation.corner_frequency
— Functioncorner_frequency(m::U, src::SourceParameters{S,T}) where {S<:Float64, T<:Real, U<:Real}
Computes a 3-tuple of corner frequency components, depending upon source spectrum type. By default the single-corner Brune spectrum is considered, but if src.model
equals :Atkinson_Silva_2000
then the components of the double-corner spectrum are returned. If some other symbol is passed then the Brune model is returned.
Examples
m = 6.0
Δσ = 100.0
κ0 = 0.035
fas = FASParams(Δσ, κ0)
# compute single corner frequency
src.model = :Brune
fc, tmp1, tmp2 = corner_frequency(m, src)
# compute double corner frequencies
src.model = :Atkinson_Silva_2000
fa, fb, ε = corner_frequency(m, src)
StochasticGroundMotionSimulation.magnitude_to_moment
— Functionmagnitude_to_moment(m::T) where T<:Real
Converts moment magnitude to seismic moment (in dyne-cm).
# Examples
m = 6.0
M0 = magnitude_to_moment(m)