Custom Models
The package now supports custom FAS and duration models for maximum flexibility.
Quick Start
using StochasticGroundMotionSimulation
# Define a custom FAS function
my_fas = (f, m, r, params) -> begin
stress_drop, Q0, kappa = params
# ... your model implementation
end
# Create a functional model
fas_model = FunctionalFASModel(my_fas, [100.0, 200.0, 0.04])
# Use it in RVT calculations
Sa = rvt_response_spectral_ordinate(1.0, 6.0, 10.0, fas_model, duration_model)Custom Type Models
For more complex models, define your own types:
struct MyFASModel <: AbstractFASModel
# your parameters
end
function compute_fas(model::MyFASModel, freq, mag, dist)
# your implementation
endForwardDiff Compatibility
All custom models must be ForwardDiff-compatible. Use the validation functions:
validate_fas_model(my_model)
validate_duration_model(my_model)API Reference
Abstract Types
StochasticGroundMotionSimulation.AbstractFASModel — Type
AbstractFASModelAbstract base type for Fourier Amplitude Spectrum models. All custom FAS models should subtype this.
sourceStochasticGroundMotionSimulation.AbstractDurationModel — Type
AbstractDurationModelAbstract base type for duration models. All custom duration models should subtype this.
sourceFAS Model Types
StochasticGroundMotionSimulation.CustomModels.FunctionalFASModel — Type
FunctionalFASModel{F}A wrapper that allows users to pass a function as a FAS model. Useful for quick prototyping without defining a new type.
Example
# Define a simple FAS function
my_fas = (f, m, r, params) -> begin
stress_drop, Q0, kappa = params
corner_freq = 4.9e6 * 1e-6 * (stress_drop/1e6)^(1/3) * 10^(-m/3)
source = 1.0 / (1.0 + (f/corner_freq)^2)
path = exp(-π * f * r / (Q0 * 30.0)) # simplified
site = exp(-π * kappa * f)
return source * path * site
end
# Create model with parameters
model = FunctionalFASModel(my_fas, [100.0, 200.0, 0.04])
# Use in computation
fas_value = compute_fas(model, 1.0, 6.0, 10.0)sourceStochasticGroundMotionSimulation.CustomModels.CustomFASModel — Type
CustomFASModelExample of a custom FAS model with structured parameters. This approach is preferred for more complex models.
sourceStochasticGroundMotionSimulation.CustomModels.HybridFASModel — Type
HybridFASModelCombines multiple FAS models with a custom combining function. Useful for mixing different model components or creating ensemble models.
Fields
models: Vector of AbstractFASModel instances to combinecombine_func: Function that combines results: (results, freq, mag, dist, params) -> combined_resultparams: Additional parameters for the combining function
Example
# Create component models
source_model = FunctionalFASModel((f, m, r, p) -> p[1] * f^(-2), [1e20])
path_model = FunctionalFASModel((f, m, r, p) -> exp(-π*f*r/(p[1]*3.5)), [200.0])
# Combine by multiplication
combine = (results, f, m, r, p) -> results[1] * results[2]
hybrid = HybridFASModel([source_model, path_model], combine, Float64[])
# Use like any other FAS model
fas = compute_fas(hybrid, 1.0, 6.0, 10.0)sourceDuration Model Types
StochasticGroundMotionSimulation.CustomModels.FunctionalDurationModel — Type
FunctionalDurationModel{F}A wrapper that allows users to pass a function as a duration model.
Example
# Define a simple duration function
my_duration = (m, r, params) -> begin
a, b, c = params
return a + b * m + c * log10(r)
end
# Create model with parameters
model = FunctionalDurationModel(my_duration, [1.0, 0.5, 0.3])
# Use in computation
duration = compute_duration(model, 6.0, 10.0)sourceCore Interface Functions
StochasticGroundMotionSimulation.CustomModels.compute_fas — Function
compute_fas(model::AbstractFASModel, freq, mag, dist)Compute the Fourier Amplitude Spectrum at a given frequency.
Arguments
model: FAS model instance containing model parametersfreq: Frequency in Hz (can be Real or ForwardDiff.Dual)mag: Earthquake magnitude (can be Real or ForwardDiff.Dual)dist: Distance in km (can be Real or ForwardDiff.Dual)
Returns
- FAS value at the given frequency
Notes
- This function must be differentiable for ForwardDiff.jl compatibility
- Avoid branching logic that depends on the values of freq, mag, or dist
- Use smooth approximations for discontinuous functions
StochasticGroundMotionSimulation.CustomModels.compute_duration — Function
compute_duration(model::AbstractDurationModel, mag, dist)Compute the duration metric.
Arguments
model: Duration model instance containing model parametersmag: Earthquake magnitude (can be Real or ForwardDiff.Dual)dist: Distance in km (can be Real or ForwardDiff.Dual)
Returns
- Duration value in seconds
Notes
- This function must be differentiable for ForwardDiff.jl compatibility
Validation Functions
StochasticGroundMotionSimulation.CustomModels.validate_fas_model — Function
validate_fas_model(model::AbstractFASModel; freq_range=(0.01, 100), mag=6.0, dist=10.0)Validate that a FAS model is properly implemented and ForwardDiff-compatible.
sourceStochasticGroundMotionSimulation.CustomModels.validate_duration_model — Function
validate_duration_model(model::AbstractDurationModel; mag=6.0, dist=10.0)Validate that a duration model is properly implemented.
sourceRVT Computation with Custom Models
StochasticGroundMotionSimulation.CustomModels.rvt_response_spectral_ordinate_custom — Function
rvt_response_spectral_ordinate_custom(
T::Real,
mag::Real,
dist::Real,
fas_model::AbstractFASModel,
duration_model::AbstractDurationModel,
damping::Real=0.05,
glxi=nothing,
glwi=nothing
)Compute response spectral ordinate using custom FAS and duration models.
Arguments
T: Oscillator period (seconds)mag: Earthquake magnitudedist: Distance (km)fas_model: Custom FAS modelduration_model: Custom duration modeldamping: Damping ratio (default: 0.05)glxi: Gauss-Legendre quadrature points (optional, for compatibility)glwi: Gauss-Legendre quadrature weights (optional, for compatibility)
Keyword Arguments (for backward compatibility)
freq_range: Frequency range for integration (Hz)nfreq: Number of frequency points for integration
Returns
- Response spectral acceleration (g)
Notes
- The
glxiandglwiparameters are included for API compatibility with the traditional RVT interface but are not currently used in this implementation, which uses log-spaced frequency points instead.
Utilities
StochasticGroundMotionSimulation.CustomModels.FourierParametersWrapper — Type
FourierParametersWrapper <: AbstractFASModelWrapper to make existing FourierParameters compatible with AbstractFASModel interface. This allows using traditional parameter-based models with the custom model interface.
Example
# Create traditional parameters
src = SourceParameters(100.0)
geo = GeometricSpreadingParameters([1.0, 50.0, Inf], [1.0, 0.5])
ane = AnelasticAttenuationParameters(200.0, 0.4)
sat = NearSourceSaturationParameters(:BT15)
path = PathParameters(geo, sat, ane)
site = SiteParameters(0.039)
fourier_params = FourierParameters(src, path, site)
# Wrap for use with custom interface
wrapper = FourierParametersWrapper(fourier_params)
# Now can use compute_fas
fas = compute_fas(wrapper, freq, mag, dist)sourceStochasticGroundMotionSimulation.CustomModels.ExistingDurationWrapper — Type
ExistingDurationWrapper <: AbstractDurationModelWrapper to make existing RandomVibrationParameters compatible with AbstractDurationModel interface.
Since excitation_duration requires both FourierParameters and RandomVibrationParameters, this wrapper stores both to enable proper delegation to the existing duration calculation.
Example
# Create traditional parameters
fourier_params = FourierParameters(...)
rvt_params = RandomVibrationParameters(:BT15)
# Wrap for use with custom interface
# Note: Must provide both parameters
wrapper = ExistingDurationWrapper(rvt_params, fourier_params)
# Now can use compute_duration
duration = compute_duration(wrapper, mag, dist)Important
This wrapper is only used when BOTH FourierParameters and RandomVibrationParameters are provided (the fully concrete case). If you're using a custom FAS model, you must also use a custom duration model.
source