Stochastic Loading Module
Public Member Functions | List of all members
stochastic::VlachosEtAl Class Reference

#include <vlachos_et_al.h>

Inheritance diagram for stochastic::VlachosEtAl:
Inheritance graph
[legend]
Collaboration diagram for stochastic::VlachosEtAl:
Collaboration graph
[legend]

Public Member Functions

 VlachosEtAl ()=default
 
 VlachosEtAl (double moment_magnitude, double rupture_distance, double vs30, double orientation, unsigned int num_spectra, unsigned int num_sims)
 
 VlachosEtAl (double moment_magnitude, double rupture_distance, double vs30, double orientation, unsigned int num_spectra, unsigned int num_sims, int seed_value)
 
virtual ~VlachosEtAl ()
 
 VlachosEtAl (const VlachosEtAl &)=delete
 
VlachosEtAloperator= (const VlachosEtAl &)=delete
 
utilities::JsonObject generate (const std::string &event_name, bool g_units=false) override
 
bool generate (const std::string &event_name, const std::string &output_location, bool g_units=false) override
 
bool time_history_family (std::vector< std::vector< double >> &time_histories, const Eigen::VectorXd &parameters) const
 
void simulate_time_history (std::vector< double > &time_history, const Eigen::MatrixXd &power_spectrum) const
 
bool post_process (std::vector< double > &time_history, const std::vector< double > &filter_imp_resp) const
 
Eigen::VectorXd identify_parameters (const Eigen::VectorXd &initial_params) const
 
std::vector< double > modal_frequencies (const std::vector< double > &parameters, const std::vector< double > &energy) const
 
std::vector< double > energy_accumulation (const std::vector< double > &parameters, const std::vector< double > &times) const
 
std::vector< double > modal_participation_factor (const std::vector< double > &parameters, const std::vector< double > &energy) const
 
std::vector< double > amplitude_modulating_function (double duration, double total_energy, const std::vector< double > &parameters, const std::vector< double > &times) const
 
Eigen::VectorXd kt_2 (const std::vector< double > &parameters, const std::vector< double > &frequencies, const std::vector< double > &highpass_butter) const
 
void rotate_acceleration (const std::vector< double > &acceleration, std::vector< double > &x_accels, std::vector< double > &y_accels, bool g_units) const
 
- Public Member Functions inherited from stochastic::StochasticModel
 StochasticModel ()=default
 
virtual ~StochasticModel ()
 
 StochasticModel (const StochasticModel &)=delete
 
StochasticModeloperator= (const StochasticModel &)=delete
 
std::string model_name () const
 

Additional Inherited Members

- Protected Attributes inherited from stochastic::StochasticModel
std::string model_name_ = "StochasticModel"
 

Detailed Description

Stochastic model for generating scenario specific ground motion time histories. This is based on the paper: Vlachos C., Papakonstantinou K.G., & Deodatis G. (2018). Predictive model for site specific simulation of ground motions based on earthquake scenarios. Earthquake Engineering & Structural Dynamics, 47(1), 195-218.

Definition at line 22 of file vlachos_et_al.h.

Constructor & Destructor Documentation

stochastic::VlachosEtAl::VlachosEtAl ( )
default

Delete default constructor

stochastic::VlachosEtAl::VlachosEtAl ( double  moment_magnitude,
double  rupture_distance,
double  vs30,
double  orientation,
unsigned int  num_spectra,
unsigned int  num_sims 
)

Construct scenario specific ground motion model based on input parameters

Parameters
[in]moment_magnitudeMoment magnitude of earthquake scenario
[in]rupture_distanceClosest-to-site rupture distance in kilometers
[in]vs30Soil shear wave velocity averaged over top 30 meters in meters per second
[in]orientationOrientation of acceleration relative to global coordinates. Represents counter-clockwise angle (in degrees) away from x-axis rotating around z-axis in right-handed coordinate system.
[in]num_spectraNumber of evolutionary power spectra that should be generated.
[in]num_simsNumber of simulated ground motion time histories that should be generated per evolutionary power

Definition at line 26 of file vlachos_et_al.cc.

stochastic::VlachosEtAl::VlachosEtAl ( double  moment_magnitude,
double  rupture_distance,
double  vs30,
double  orientation,
unsigned int  num_spectra,
unsigned int  num_sims,
int  seed_value 
)

Construct scenario specific ground motion model based on input parameters

Parameters
[in]moment_magnitudeMoment magnitude of earthquake scenario
[in]rupture_distanceClosest-to-site rupture distance in kilometers
[in]vs30Soil shear wave velocity averaged over top 30 meters in meters per second
[in]orientationOrientation of acceleration relative to global coordinates. Represents counter-clockwise angle (in degrees) away from x-axis rotating around z-axis in right-handed coordinate system.
[in]num_spectraNumber of evolutionary power spectra that should be generated.
[in]num_simsNumber of simulated ground motion time histories that should be generated per evolutionary power
[in]seed_valueValue to seed random variables with to ensure repeatability

Definition at line 213 of file vlachos_et_al.cc.

virtual stochastic::VlachosEtAl::~VlachosEtAl ( )
inlinevirtual

Virtual destructor

Definition at line 74 of file vlachos_et_al.h.

stochastic::VlachosEtAl::VlachosEtAl ( const VlachosEtAl )
delete

Delete copy constructor

Member Function Documentation

std::vector< double > stochastic::VlachosEtAl::amplitude_modulating_function ( double  duration,
double  total_energy,
const std::vector< double > &  parameters,
const std::vector< double > &  times 
) const

Calculate amplitude modulating function as defined by Eq-7 on page 6

Parameters
[in]durationTotal duration of target seismic record
[in]total_energyTotal energy content of seismic ground acceleration
[in]parametersVector of values containing energy parameters gamma and delta
[in]timesVector containing non-dimensional times at which to calculate energy accumulation
Returns
Vector containing values of amplitude modulating function at input times

Definition at line 849 of file vlachos_et_al.cc.

std::vector< double > stochastic::VlachosEtAl::energy_accumulation ( const std::vector< double > &  parameters,
const std::vector< double > &  times 
) const

Calculate the energy accumulation over the times specified in the inputs as defined by Eq-5 on page 6.

Parameters
[in]parametersVector of values containing energy parameters gamma and delta
[in]timesVector containing non-dimensional times at which to calculate energy accumulation
Returns
Vector containing accumulated energy values at specified input times

Definition at line 816 of file vlachos_et_al.cc.

utilities::JsonObject stochastic::VlachosEtAl::generate ( const std::string &  event_name,
bool  g_units = false 
)
overridevirtual

Generate ground motion time histories based on input parameters and store outputs as JSON object. Throws exception if errors are encountered during time history generation.

Parameters
[in]event_nameName to assign to event
[in]g_unitsIndicates that time histories should be returned in units of g. Defaults to false where time histories are returned in units of m/s^2
Returns
JsonObject containing time histories

Implements stochastic::StochasticModel.

Definition at line 400 of file vlachos_et_al.cc.

bool stochastic::VlachosEtAl::generate ( const std::string &  event_name,
const std::string &  output_location,
bool  g_units = false 
)
overridevirtual

Generate ground motion time histories based on input parameters and write results to file in JSON format. Throws exception if errors are encountered during time history generation.

Parameters
[in]event_nameName to assign to event
[in,out]output_locationLocation to write outputs to
[in]g_unitsIndicates that time histories should be returned in units of g. Defaults to false where time histories are returned in units of m/s^2
Returns
Returns true if successful, false otherwise

Implements stochastic::StochasticModel.

Definition at line 476 of file vlachos_et_al.cc.

Eigen::VectorXd stochastic::VlachosEtAl::identify_parameters ( const Eigen::VectorXd &  initial_params) const

Identifies modal frequency parameters for mode 1 and 2

Parameters
[in]initial_paramsInitial set of parameters
Returns
Vector of identified parameters

Definition at line 715 of file vlachos_et_al.cc.

Eigen::VectorXd stochastic::VlachosEtAl::kt_2 ( const std::vector< double > &  parameters,
const std::vector< double > &  frequencies,
const std::vector< double > &  highpass_butter 
) const

Calculate evolutionary power spectrum at specific time using parametric, bimodal, fully non-stationary Kinai-Tajimi (K-T) model of seismic ground acceleration signal consisting of 2 distinct spectral modes. This is described by Eq-1 on page 4.

Parameters
[in]parametersVector of parameters fg_1, zeta_1, S0_1, fg_2, zeta_2 and S0_2 which define the time varying modal dominant frequency, model apparent damping ratio, and modal participation factor associated with the two modes. NOTE: Values of fg_1, fg_2 and S0_2 must be for the particular time at which the power spectrum is desired
[in]frequenciesVector of frequencies at which to evaluate evolutionary power spectrum
[in]highpass_butterVector of Butterworth filter transfer function energy content at input frequencies
Returns
Vector containing values for model evolutionary power spectrum at time defined by input parameters fg_1, fg_2 and S0_2.

Definition at line 866 of file vlachos_et_al.cc.

std::vector< double > stochastic::VlachosEtAl::modal_frequencies ( const std::vector< double > &  parameters,
const std::vector< double > &  energy 
) const

Calculate the dominant modal frequencies as a function of non-dimensional cumulative energy and the model parameters Q_k, alpha_k, and beta_k where k is the mode (either 1 or 2). This is defined by Eq-8 on page 6.

Parameters
[in]parametersVector of values for alpha_k, beta_k and Q_k at k-th mode
[in]energyVector of non-dimensional energy values at which to calculate frequency
Returns
Vector of frequency values at input energy values

Definition at line 803 of file vlachos_et_al.cc.

std::vector< double > stochastic::VlachosEtAl::modal_participation_factor ( const std::vector< double > &  parameters,
const std::vector< double > &  energy 
) const

Calculate the logarithmic normalized modal participation factor for the second mode as defined by Eq-11 on page 7.

Parameters
[in]parametersVector of values for F(I), mu(I), sigma(I), F(II), mu(II) and sigma(II)
[in]energyVector of non-dimensional energy values at which to calculate frequency
Returns
Vector containing values of logarithmic normalized modal participation factor for second mode at input energy values

Definition at line 830 of file vlachos_et_al.cc.

VlachosEtAl& stochastic::VlachosEtAl::operator= ( const VlachosEtAl )
delete

Delete assignment operator

bool stochastic::VlachosEtAl::post_process ( std::vector< double > &  time_history,
const std::vector< double > &  filter_imp_resp 
) const

Post-process the input time history as described in Vlachos et al. using multiple-window estimation technique after Conte & Peng (1997) and highpass Butterworth filter

Parameters
[in,out]time_historyTime history to post-process. Post-processed results are also stored here.
[in]filter_imp_respImpulse response of Butterworth filter
Returns
Returns true if successful, false otherwise

Definition at line 662 of file vlachos_et_al.cc.

void stochastic::VlachosEtAl::rotate_acceleration ( const std::vector< double > &  acceleration,
std::vector< double > &  x_accels,
std::vector< double > &  y_accels,
bool  g_units 
) const

Rotate acceleration based on orientation angle

Parameters
[in]accelerationAcceleration to rotate
[in,out]x_accelsVector to store x-component of acceleration to
[in,out]y_accelsVector to story y-component of acceleration to
[in]g_unitsIndicates that time histories should be returned in units of g

Definition at line 893 of file vlachos_et_al.cc.

void stochastic::VlachosEtAl::simulate_time_history ( std::vector< double > &  time_history,
const Eigen::MatrixXd &  power_spectrum 
) const

Simulate fully non-stationary ground motion sample realization based on time and frequency discretization and the discretized evolutionary power spectrum. This is described by Eq-19 on page 8.

Parameters
[in,out]time_historyLocation where time history should be stored
[in]power_spectrumMatrix containing values of power spectrum over range of frequencies at specified times.

Definition at line 612 of file vlachos_et_al.cc.

bool stochastic::VlachosEtAl::time_history_family ( std::vector< std::vector< double >> &  time_histories,
const Eigen::VectorXd &  parameters 
) const

Compute a family of time histories for a particular power spectrum

Parameters
[in,out]time_historiesLocation where time histories should be stored
[in]parametersSet of model parameters to use for calculating power specturm and time histories
Returns
Returns true if successful, false otherwise

Definition at line 494 of file vlachos_et_al.cc.


The documentation for this class was generated from the following files: