Source code for openquake.hazardlib.gsim.kotha_2019

# -*- coding: utf-8 -*-
# vim: tabstop=4 shiftwidth=4 softtabstop=4
#
# Copyright (C) 2014-2019 GEM Foundation
#
# OpenQuake is free software: you can redistribute it and/or modify it
# under the terms of the GNU Affero General Public License as published
# by the Free Software Foundation, either version 3 of the License, or
# (at your option) any later version.
#
# OpenQuake is distributed in the hope that it will be useful,
# but WITHOUT ANY WARRANTY; without even the implied warranty of
# MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
# GNU Affero General Public License for more details.
#
# You should have received a copy of the GNU Affero General Public License
# along with OpenQuake. If not, see <http://www.gnu.org/licenses/>.

"""
Module exports :class:`KothaEtAl2019`,
               :class:`KothaEtAl2019SERA`,
               :class:`KothaEtAl2019Site`,
               :class:`KothaEtAl2019Slope`,
               :class:`KothaEtAl2019SERASlopeGeology`
"""
import os
import h5py
import numpy as np
from scipy.constants import g
from scipy.interpolate import interp1d
from openquake.hazardlib.gsim.base import GMPE, CoeffsTable
from openquake.hazardlib import const
from openquake.hazardlib.imt import PGA, PGV, SA, from_string


BASE_PATH = os.path.join(os.path.dirname(__file__), "kotha_2019_tables")


[docs]class KothaEtAl2019(GMPE): """ Implements the first complete version of the newly derived GMPE for Shallow Crustal regions using the Engineering Strong Motion Flatfile. (Working Title) Kotha, S. R., Weatherill, G., Bindi, D., Cotton F. (2019) A Revised Ground Motion Prediction Equation for Shallow Crustal Earthquakes in Europe and the Middle-East The GMPE is desiged for calibration of the stress parameter term (a multiple of the fault-to-fault variability, tau_f) an attenuation scaling term (c3) and a statistical uncertainty term (sigma_mu). The statistical uncertainty is a scalar factor dependent on period, magnitude and distance. These are read in from hdf5 upon instantiation and interpolated to the necessary values. In the core form of the GMPE no site term is included. This will be added in the subclasses. :param c3: User supplied table for the coefficient c3 controlling the anelastic attenuation as an instance of :class: `openquake.hazardlib.gsim.base.CoeffsTable`. If absent, the value is taken from the normal coefficients table. :param sigma_mu_epsilon: The number by which to multiply the epistemic uncertainty (sigma_mu) for the adjustment of the mean ground motion. """ experimental = True #: Supported tectonic region type is 'active shallow crust' DEFINED_FOR_TECTONIC_REGION_TYPE = const.TRT.ACTIVE_SHALLOW_CRUST #: Set of :mod:`intensity measure types <openquake.hazardlib.imt>` #: this GSIM can calculate. A set should contain classes from module #: :mod:`openquake.hazardlib.imt`. DEFINED_FOR_INTENSITY_MEASURE_TYPES = set([ PGA, PGV, SA ]) #: Supported intensity measure component is the geometric mean of two #: horizontal components DEFINED_FOR_INTENSITY_MEASURE_COMPONENT = const.IMC.RotD50 #: Supported standard deviation types are inter-event, intra-event #: and total DEFINED_FOR_STANDARD_DEVIATION_TYPES = set([ const.StdDev.TOTAL, const.StdDev.INTER_EVENT, const.StdDev.INTRA_EVENT ]) #: Required site parameter is not set REQUIRES_SITES_PARAMETERS = set() #: Required rupture parameters are magnitude and hypocentral depth REQUIRES_RUPTURE_PARAMETERS = set(('mag', 'hypo_depth')) #: Required distance measure is Rjb (eq. 1). REQUIRES_DISTANCES = set(('rjb', )) def __init__(self, sigma_mu_epsilon=0.0, c3=None, ergodic=True): """ Instantiate setting the sigma_mu_epsilon and c3 terms """ super().__init__() if isinstance(c3, dict): # Inputing c3 as a dictionary sorted by the string representation # of the IMT c3in = {} for c3key in c3: c3in[from_string(c3key)] = {"c3": c3[c3key]} self.c3 = CoeffsTable(sa_damping=5, table=c3in) else: self.c3 = c3 self.sigma_mu_epsilon = sigma_mu_epsilon self.ergodic = ergodic if self.sigma_mu_epsilon: # Connect to hdf5 and load tables into memory self.retreive_sigma_mu_data() else: # No adjustments, so skip this step self.mags = None self.dists = None self.s_a = None self.pga = None self.pgv = None self.periods = None
[docs] def retreive_sigma_mu_data(self): """ For the general form of the GMPE this retrieves the sigma mu values from the hdf5 file using the "general" model, i.e. sigma mu factors that are independent of the choice of region or depth """ fle = h5py.File(os.path.join(BASE_PATH, "KothaEtAl2019_SigmaMu_Fixed.hdf5"), "r") self.mags = fle["M"][:] self.dists = fle["R"][:] self.periods = fle["T"][:] self.pga = fle["PGA"][:] self.pgv = fle["PGV"][:] self.s_a = fle["SA"][:] fle.close()
[docs] def get_mean_and_stddevs(self, sites, rup, dists, imt, stddev_types): """ See :meth:`superclass method <.base.GroundShakingIntensityModel.get_mean_and_stddevs>` for spec of input and result values. """ # extracting dictionary of coefficients specific to required # intensity measure type. C = self.COEFFS[imt] mean = (self.get_magnitude_scaling(C, rup.mag) + self.get_distance_term(C, rup, dists.rjb, imt) + self.get_site_amplification(C, sites, imt)) # GMPE originally in cm/s/s - convert to g if imt.name in "PGA SA": mean -= np.log(100.0 * g) stddevs = self.get_stddevs(C, dists.rjb.shape, stddev_types, sites) if self.sigma_mu_epsilon: # Apply the epistemic uncertainty factor (sigma_mu) multiplied by # the number of standard deviations sigma_mu = self.get_sigma_mu_adjustment(C, imt, rup, dists) # Cap sigma_mu at 0.5 ln units sigma_mu[sigma_mu > 0.5] = 0.5 # Sigma mu should not be less than the standard deviation of the # fault-to-fault variability sigma_mu[sigma_mu < C["tau_fault"]] = C["tau_fault"] mean += (self.sigma_mu_epsilon * sigma_mu) return mean, stddevs
[docs] def get_magnitude_scaling(self, C, mag): """ Returns the magnitude scaling term """ d_m = mag - self.CONSTANTS["Mh"] if mag <= self.CONSTANTS["Mh"]: return C["e1"] + C["b1"] * d_m + C["b2"] * (d_m ** 2.0) else: return C["e1"] + C["b3"] * d_m
[docs] def get_distance_term(self, C, rup, rjb, imt): """ Returns the distance attenuation factor """ h = self._get_h(C, rup.hypo_depth) rval = np.sqrt(rjb ** 2. + h ** 2.) rref_val = np.sqrt(self.CONSTANTS["Rref"] ** 2. + h ** 2.) c3 = self.get_distance_coefficients(C, imt) f_r = (C["c1"] + C["c2"] * (rup.mag - self.CONSTANTS["Mref"])) *\ np.log(rval / rref_val) + (c3 * (rval - rref_val) / 100.) return f_r
def _get_h(self, C, hypo_depth): """ Returns the depth-specific coefficient """ if hypo_depth <= 10.0: return self.CONSTANTS["h_D10"] elif hypo_depth > 20.0: return self.CONSTANTS["h_D20"] else: return self.CONSTANTS["h_10D20"]
[docs] def get_distance_coefficients(self, C, imt): """ Returns the c3 term """ c3 = self.c3[imt]["c3"] if self.c3 else C["c3"] return c3
[docs] def get_site_amplification(self, C, sites, imt): """ In base model no site amplification is used """ return 0.0
[docs] def get_sigma_mu_adjustment(self, C, imt, rup, dists): """ Returns the sigma mu adjustment factor """ if imt.name in "PGA PGV": # PGA and PGV are 2D arrays of dimension [nmags, ndists] sigma_mu = getattr(self, imt.name.lower()) if rup.mag <= self.mags[0]: sigma_mu_m = sigma_mu[0, :] elif rup.mag >= self.mags[-1]: sigma_mu_m = sigma_mu[-1, :] else: intpl1 = interp1d(self.mags, sigma_mu, axis=0) sigma_mu_m = intpl1(rup.mag) # Linear interpolation with distance intpl2 = interp1d(self.dists, sigma_mu_m, bounds_error=False, fill_value=(sigma_mu_m[0], sigma_mu_m[-1])) return intpl2(dists.rjb) # In the case of SA the array is of dimension [nmags, ndists, nperiods] # Get values for given magnitude if rup.mag <= self.mags[0]: sigma_mu_m = self.s_a[0, :, :] elif rup.mag >= self.mags[-1]: sigma_mu_m = self.s_a[-1, :, :] else: intpl1 = interp1d(self.mags, self.s_a, axis=0) sigma_mu_m = intpl1(rup.mag) # Get values for period - N.B. ln T, linear sigma mu interpolation if imt.period <= self.periods[0]: sigma_mu_t = sigma_mu_m[:, 0] elif imt.period >= self.periods[-1]: sigma_mu_t = sigma_mu_m[:, -1] else: intpl2 = interp1d(np.log(self.periods), sigma_mu_m, axis=1) sigma_mu_t = intpl2(np.log(imt.period)) intpl3 = interp1d(self.dists, sigma_mu_t, bounds_error=False, fill_value=(sigma_mu_t[0], sigma_mu_t[-1])) return intpl3(dists.rjb)
[docs] def get_stddevs(self, C, stddev_shape, stddev_types, sites): """ Returns the standard deviations """ stddevs = [] tau = C["tau_event"] phi = C["phi0"] if self.ergodic: phi = np.sqrt(phi ** 2. + C["phis2s"] ** 2.) for stddev_type in stddev_types: assert stddev_type in self.DEFINED_FOR_STANDARD_DEVIATION_TYPES if stddev_type == const.StdDev.TOTAL: stddevs.append(np.sqrt(tau ** 2. + phi ** 2.) + np.zeros(stddev_shape)) elif stddev_type == const.StdDev.INTRA_EVENT: stddevs.append(phi + np.zeros(stddev_shape)) elif stddev_type == const.StdDev.INTER_EVENT: stddevs.append(tau + np.zeros(stddev_shape)) return stddevs
COEFFS = CoeffsTable(sa_damping=5, table="""\ imt e1 b1 b2 b3 c1 c2 c3 tau_c3 phis2s tau_event tau_fault phi0 g0_vs30 g1_vs30 g2_vs30 phis2s_vs30 g0_slope g1_slope g2_slope phis2s_slope pgv 2.025713592693850 2.406038866414110 0.215611393306395 0.291297120242859 -1.403387613083740 0.281415036690207 -0.358475654018642 0.158122937530244 0.547281437240368 0.403913422762195 0.299427366002708 0.447957473198010 -0.182438591933237 -0.497305420073801 0.013111800993100 0.336555436174304 -0.032073084806331 -0.162025214998054 -0.011307802652687 0.415124591549937 pga 4.687702649605810 1.953983695221440 0.185885068933981 -0.241738900149521 -1.478079598457430 0.288086630948598 -0.678379247430291 0.263916141885138 0.626209850856568 0.402258312134260 0.416094806703880 0.466947277758691 -0.176814652331689 -0.577689600440376 -0.138853673082031 0.390464689737078 -0.013277224195799 -0.129210051415153 -0.018181356488862 0.522940223323660 0.010 4.689953753124180 1.953153318951010 0.186226399234793 -0.242289907409910 -1.480356424151250 0.288439060230840 -0.677333232895566 0.263885400397128 0.626583207016585 0.402108659213375 0.416770119149264 0.466978946978010 -0.176590552121014 -0.577276296453196 -0.138997049798249 0.390027772318435 -0.013304099558285 -0.129063012511688 -0.018181472980112 0.524762067903723 0.025 4.713377198166330 1.936551887385350 0.189245873356868 -0.242251953373212 -1.527209359618410 0.289740058354006 -0.638666791582353 0.262164920363914 0.631328171215036 0.398276492676289 0.420866218374451 0.468574643779387 -0.173431568326523 -0.564867479776383 -0.136825809696313 0.387036285288661 -0.012581551255242 -0.126800131880132 -0.018534777394440 0.526867826025094 0.040 4.795848182057760 1.898151028508490 0.201336642438332 -0.237871501230899 -1.624031787587650 0.304689714377880 -0.597052092823935 0.251662690566386 0.653172181191298 0.388148609537982 0.442362727924819 0.473434298534627 -0.155802568231936 -0.539764415936919 -0.144284241970030 0.414375380025667 -0.011470899481726 -0.123849772957656 -0.019533153675162 0.537006725941913 0.050 4.891963784319470 1.888058117695470 0.209278189591062 -0.269811787471139 -1.650049238074590 0.318075099135608 -0.621868044857147 0.265892288510560 0.667742830430293 0.385327490491663 0.468784102021719 0.479389514380334 -0.154750389172492 -0.527932614967297 -0.143759700115643 0.412391276624274 -0.012554601043835 -0.120775767836190 -0.018659349384142 0.552523848710959 0.070 5.106825186036690 1.874242004704420 0.212124686343632 -0.320099420113234 -1.620479311265060 0.320046491923063 -0.724343238043036 0.291387067989580 0.693461972699411 0.392248296873683 0.489423169198560 0.485345337361231 -0.162237529464212 -0.518783615288977 -0.152448665014284 0.435631365330907 -0.009349344181752 -0.117606720279786 -0.019815999895356 0.587882707304527 0.100 5.347455916482890 1.868122162855800 0.199938780857780 -0.356072391432107 -1.517710013159360 0.289475493518785 -0.829503523244518 0.332348696555227 0.695072073756251 0.402824216442605 0.512338505510147 0.493702647513133 -0.152133628557392 -0.546449196888745 -0.189976803645333 0.402824194957262 -0.009005626438202 -0.106760536180492 -0.016109590193024 0.582026423711337 0.150 5.549588796784530 1.848063879486740 0.152771408659783 -0.348229196634948 -1.358583560220900 0.230015236938864 -0.889955230608062 0.338152565231841 0.682916099464874 0.408594346059644 0.479823906245635 0.496027258004371 -0.174546733673003 -0.615418367954405 -0.210112041223272 0.446909590011407 -0.005021668511935 -0.109778372621890 -0.018811242281261 0.591287386920046 0.200 5.618412125074440 1.885071394426640 0.129823775224324 -0.292676175821898 -1.293132464925350 0.185865789481763 -0.834473063837879 0.308867853425155 0.669218276129109 0.429871752908997 0.417131887788201 0.494836307784100 -0.189909912690980 -0.695953862273611 -0.219447643781414 0.438438738054802 -0.004873785069078 -0.128463642498102 -0.022525331953886 0.577959445522595 0.250 5.544966205863690 1.889841953164010 0.105763631540550 -0.185382496096154 -1.250626677741170 0.158041062691118 -0.776587060423500 0.277393830052251 0.646372161163196 0.425139089741650 0.389555883970782 0.488371422580152 -0.190776384198166 -0.675387322023739 -0.197048540176423 0.453984574385616 -0.012909132753988 -0.136516396256906 -0.020223378666953 0.536361261247963 0.300 5.450353086526580 1.923037058761700 0.093388646106402 -0.133518495442793 -1.226366053882270 0.141154449809013 -0.718787095472858 0.275157669080336 0.625211078118202 0.423146216277508 0.378059931938700 0.482704664228660 -0.191198552619570 -0.664576488567174 -0.175749756009659 0.425411122163987 -0.016962208092821 -0.144010286645154 -0.019080440896514 0.512960425931692 0.350 5.334813762532070 1.955199629413510 0.086785841616639 -0.028679427285532 -1.204150467410360 0.130192388780628 -0.673786484345681 0.259619616389265 0.620539435099090 0.420137510023815 0.339003245617876 0.481586998560825 -0.201002954926354 -0.669087444350307 -0.142578290289139 0.380441367195110 -0.022875317564253 -0.157635709265738 -0.018165563045804 0.511808957079240 0.400 5.271854940995910 2.029644803968490 0.091116974902809 -0.047769201478982 -1.178699297309130 0.119724051296491 -0.642239681195656 0.246279727018558 0.625819372379235 0.414328187301189 0.334544067635360 0.475978739636533 -0.214631481381591 -0.686924941969383 -0.119932391592713 0.396506175978199 -0.025411983908599 -0.170701008403919 -0.019422521300558 0.499218491920244 0.450 5.178258137551410 2.069576283726890 0.090216666792283 -0.039988422268797 -1.167332473566050 0.114113990185286 -0.603205730169031 0.247585522515150 0.629409711420667 0.411924382020225 0.333145123580464 0.469784548562159 -0.215626912041199 -0.664913647056549 -0.080933694394528 0.428435020434419 -0.032902127419588 -0.178844859528253 -0.017472789991119 0.487209738771248 0.500 5.082700791920750 2.099155311093840 0.086780566158664 0.009747490727000 -1.167604250112560 0.107042637521344 -0.563257403363693 0.239584431802383 0.630142062016198 0.404800890388222 0.327182920069769 0.462936819509800 -0.216242812282768 -0.653744666745149 -0.054735637081128 0.411697056400221 -0.037662523031403 -0.184328658869514 -0.015930481499202 0.483211504332540 0.600 4.893472350204460 2.164821202949550 0.084969812782495 0.098861244388268 -1.157406472584810 0.100755617814698 -0.497303864164013 0.219211803845992 0.630252473395296 0.405655568634848 0.315141786853113 0.451108093286210 -0.212080108734536 -0.606862157074996 0.004643862201540 0.387963928568662 -0.041150172334837 -0.190596930116924 -0.013925137557033 0.491713186187139 0.700 4.705054382687820 2.212434833698710 0.081637278906465 0.100462763700203 -1.152016011593760 0.093631933132715 -0.440618892293181 0.214461662068525 0.627465836385995 0.401597002148329 0.333218982022375 0.447706507848849 -0.211612134937899 -0.577781379445103 0.059708795638049 0.401715451807114 -0.041484083678375 -0.190773062795295 -0.012575628562317 0.498751240066827 0.750 4.616038943804310 2.235375898851390 0.079021945207681 0.154139215140058 -1.142973389858090 0.088852396569925 -0.419724044373000 0.207438203391352 0.628069780963942 0.407735734666991 0.326474473374952 0.445225236252333 -0.208196130269188 -0.564953996039896 0.079855677797207 0.392778956892387 -0.041053305954862 -0.191671551621522 -0.012839677141417 0.488875554454694 0.800 4.527849452198010 2.264397979442970 0.081390475773282 0.202229318059039 -1.136648599451540 0.090459890537105 -0.403124991422463 0.189995354821562 0.625851513946276 0.415349466205398 0.317303066311374 0.440278408506843 -0.207821311493284 -0.544622358152350 0.097916912087933 0.395150348296333 -0.040756368302139 -0.192366021205877 -0.012770161740153 0.496597703262078 0.900 4.391239646181240 2.325276988761640 0.083821056747981 0.228778153653814 -1.127292034561540 0.089501921986240 -0.374788234612032 0.167063633965528 0.622882319152969 0.419993397797855 0.323260336648660 0.434351634884768 -0.201326062730655 -0.506240273665515 0.131764693787857 0.427279270677769 -0.039070225366038 -0.196384036683547 -0.014885815281926 0.487436451683401 1.000 4.257838894688260 2.367725674978040 0.084432350092012 0.341875244589753 -1.116475252259780 0.090216898231893 -0.353209184896043 0.152524622009451 0.618724312927540 0.426363619664496 0.306907261128117 0.427075711161247 -0.195585911487788 -0.471587459636263 0.170765965273430 0.428747852389176 -0.042322588650302 -0.196373316070861 -0.013577948619500 0.474097282940156 1.200 4.084000623932080 2.535516003279140 0.113291072791050 0.341105428564444 -1.113281899461570 0.101198267860225 -0.307609138051120 0.136322298762975 0.619667399305834 0.432812736959479 0.288816783423553 0.416736878326100 -0.197149999245979 -0.454881557729585 0.191345940624623 0.411794472026739 -0.048443328633651 -0.200482529897246 -0.011090618891798 0.473036686625438 1.400 3.858720456780630 2.616885295578470 0.122899590123645 0.433898170105913 -1.123145841250970 0.104462233472800 -0.263665247407548 0.126451619959459 0.624067599637613 0.447550641437210 0.264214688492163 0.408709837205184 -0.192315162660632 -0.428772115769810 0.191177201860379 0.424012365377324 -0.048688377977240 -0.204173690423109 -0.011113848967168 0.494918868371392 1.600 3.662163462447190 2.693561483172600 0.134723407599728 0.414339375834985 -1.138808268506760 0.115930381901348 -0.218847134352892 0.131025001289405 0.621514663534016 0.453656810401114 0.252582646468883 0.404631540526549 -0.182820341877323 -0.408004013819723 0.191459322848169 0.421307232691460 -0.052400028644053 -0.201581753853561 -0.008071819851122 0.479369207289037 1.800 3.504335982115610 2.804194382503690 0.162346194743908 0.350562313097486 -1.152524904581590 0.120074987873089 -0.185702153948541 0.137842414980722 0.621624707587364 0.466001615468431 0.225283551221296 0.398741991979155 -0.197524528430911 -0.416910215242165 0.187853421320560 0.420068790275599 -0.056968393800159 -0.197938713020462 -0.002566124187343 0.486073235669474 2.000 3.352901148637760 2.871065334242860 0.175304278490051 0.301922671452104 -1.163940125236620 0.122626640563880 -0.159248285466522 0.140849131837099 0.614398039046221 0.474830104830787 0.210499929051966 0.396429621856767 -0.195284524507871 -0.415793696081323 0.167469518698839 0.399526349489727 -0.057012355893419 -0.198068134894897 -0.002649466846245 0.477257284963676 2.500 3.022521552262530 3.060144403282750 0.230493278721300 0.299896341050341 -1.172681966940500 0.146359630952262 -0.136859516271201 0.166320529763324 0.592031938560107 0.477869423474038 0.197657986712972 0.392547829352281 -0.200209246687191 -0.368243515228489 0.191269363305152 0.377508193666366 -0.056246216722837 -0.191863014055096 -0.002790191019170 0.463031276921082 3.000 2.783676662198770 3.235592094294130 0.279996868952325 0.327515809704831 -1.151782848584300 0.155626264164908 -0.169153068756301 0.162858122205747 0.587593509864436 0.484152326322702 0.203965885466933 0.389782912231145 -0.195979560707066 -0.343935003431511 0.183734350750094 0.349965866471484 -0.055125397660279 -0.188661961972340 -0.001705071124028 0.459864765313606 3.500 2.550184139271760 3.327248776655720 0.313799508552235 0.427593125870714 -1.157788818464150 0.159502683400810 -0.166147766003538 0.172286178709976 0.573188890082076 0.457777235057253 0.216074174034918 0.387317520863446 -0.196868004602194 -0.298233432015600 0.193122399737634 0.344353816098595 -0.057093151610144 -0.174713255554241 0.001639289848024 0.453901292150480 4.000 2.341796330718510 3.391100460884400 0.334750578303044 0.507655543762405 -1.174149262226720 0.179644555378667 -0.162312356740869 0.174567193920684 0.566535043868093 0.463984870088957 0.202192053366099 0.386005461021883 -0.193109038618532 -0.270710175119568 0.206241030566559 0.343394941733746 -0.054623545271567 -0.168496860827395 0.001461376635318 0.456824250621315 4.500 2.242498493808660 3.476962562487750 0.380318565940050 0.666872895571886 -1.172337191947140 0.199618659323079 -0.172605613836776 0.145926141801792 0.577217595569856 0.446117692446200 0.203666933838210 0.371586179805776 -0.187797880786474 -0.229625458712285 0.205255565254180 0.354032312869895 -0.049478507179597 -0.180378011984364 -0.002198413451819 0.462112739563011 5.000 2.009611059381530 3.508361778014590 0.395686780409462 0.861356359764995 -1.208963806141110 0.228641145313669 -0.140697393026636 0.136942909437525 0.559951442166458 0.448040327766279 0.209658464567696 0.377238739334535 -0.179857347925864 -0.209509661770966 0.196942074093841 0.334943922025071 -0.046278579432035 -0.171052829621680 -0.002650996842423 0.458091860533063 6.000 1.688862627964000 3.516737645416790 0.420264280758498 1.309501218469920 -1.214089541959290 0.237065702688919 -0.133212839440133 0.141314902795929 0.536431141593929 0.478894088472203 0.212035440719214 0.384111550789384 -0.170768904943291 -0.178067139919822 0.202603258158714 0.318505170939545 -0.041013415286417 -0.158353028710237 -0.003054369451784 0.427211301123401 7.000 1.341373165722010 3.465765382316050 0.416036203766797 1.376730174631370 -1.277754681712860 0.267078137770598 -0.088623284188977 0.145375015287528 0.522666557851543 0.427622455537343 0.214045145419707 0.384907261550421 -0.176688339526200 -0.184550689338629 0.199593069787839 0.300934929804920 -0.039851734566937 -0.147823734420266 -0.001648593255893 0.422037425667615 8.000 1.063213165919580 3.450720620500530 0.427534350267617 1.454576290963340 -1.323042041332830 0.291715072522946 -0.064115203294262 0.151582553207261 0.507211779986322 0.422147964378804 0.234785059413908 0.386735725945679 -0.166955752594165 -0.176106622570061 0.195250506134964 0.315321534521396 -0.036026205230018 -0.145317983302120 -0.002646907312481 0.410746028675711 """) CONSTANTS = {"Mref": 4.5, "Rref": 30., "Mh": 6.2, "h_D10": 4.0, "h_10D20": 8.0, "h_D20": 12.0}
[docs]class KothaEtAl2019Site(KothaEtAl2019): """ Preliminary adaptation of the Kotha et al. (2019) GMPE using a polynomial site amplification function dependent on Vs30 (m/s) """ #: Required site parameter is not set REQUIRES_SITES_PARAMETERS = set(("vs30",))
[docs] def get_site_amplification(self, C, sites, imt): """ Defines a second order polynomial site amplification model """ # Render with respect to 800 m/s reference Vs30 sref = np.log(sites.vs30 / 800.) return C["g0_vs30"] + C["g1_vs30"] * sref + C["g2_vs30"] * (sref ** 2.)
[docs] def get_stddevs(self, C, stddev_shape, stddev_types, sites): """ Returns the standard deviations """ stddevs = [] tau = C["tau_event"] phi = np.sqrt(C["phi0"] ** 2.0 + C["phis2s_vs30"] ** 2.) for stddev_type in stddev_types: assert stddev_type in self.DEFINED_FOR_STANDARD_DEVIATION_TYPES if stddev_type == const.StdDev.TOTAL: stddevs.append(np.sqrt(tau ** 2. + phi ** 2.) + np.zeros(stddev_shape)) elif stddev_type == const.StdDev.INTRA_EVENT: stddevs.append(phi + np.zeros(stddev_shape)) elif stddev_type == const.StdDev.INTER_EVENT: stddevs.append(tau + np.zeros(stddev_shape)) return stddevs
[docs]class KothaEtAl2019Slope(KothaEtAl2019): """ Preliminary adaptation of the Kotha et al. (2019) GMPE using a polynomial site amplification function dependent on slope (m/m) """ #: Required site parameter is not set REQUIRES_SITES_PARAMETERS = set(("slope",))
[docs] def get_site_amplification(self, C, sites, imt): """ Defines a second order polynomial site amplification model """ # Render with respect to 0.1 m/m reference slope sref = np.log(sites.slope / 0.1) return C["g0_slope"] + C["g1_slope"] * sref +\ C["g2_slope"] * (sref ** 2.)
[docs] def get_stddevs(self, C, stddev_shape, stddev_types, sites): """ Returns the standard deviations """ stddevs = [] tau = C["tau_event"] phi = C["phi0"] if self.ergodic: phi = np.sqrt(phi ** 2. + C["phis2s_slope"] ** 2.) for stddev_type in stddev_types: assert stddev_type in self.DEFINED_FOR_STANDARD_DEVIATION_TYPES if stddev_type == const.StdDev.TOTAL: stddevs.append(np.sqrt(tau ** 2. + phi ** 2.) + np.zeros(stddev_shape)) elif stddev_type == const.StdDev.INTRA_EVENT: stddevs.append(phi + np.zeros(stddev_shape)) elif stddev_type == const.StdDev.INTER_EVENT: stddevs.append(tau + np.zeros(stddev_shape)) return stddevs
[docs]class KothaEtAl2019SERA(KothaEtAl2019): """ Implementation of the Kotha et al. (2019) GMPE with the site amplification components included using a two-segment piecewise linear function. This form of the GMPE defines the site in terms of a measured or inferred Vs30, with the total aleatory variability adjusted accordingly. """ #: Required site parameter is not set REQUIRES_SITES_PARAMETERS = set(("vs30", "vs30measured"))
[docs] def get_site_amplification(self, C, sites, imt): """ Returns the linear site amplification term depending on whether the Vs30 is observed of inferred """ vs30 = np.copy(sites.vs30) vs30[vs30 > 1100.] = 1100. ampl = np.zeros(vs30.shape) # For observed vs30 sites ampl[sites.vs30measured] = (C["d0_obs"] + C["d1_obs"] * np.log(vs30[sites.vs30measured])) # For inferred Vs30 sites idx = np.logical_not(sites.vs30measured) ampl[idx] = (C["d0_inf"] + C["d1_inf"] * np.log(vs30[idx])) return ampl
[docs] def get_stddevs(self, C, stddev_shape, stddev_types, sites): """ Returns the standard deviations, with different site standard deviation for inferred vs. observed vs30 sites. """ stddevs = [] tau = C["tau_event"] phi = C["phi0"] if self.ergodic: phi_s2s = np.zeros(sites.vs30measured.shape, dtype=float) phi_s2s[sites.vs30measured] += C["phi_s2s_obs"] phi_s2s[np.logical_not(sites.vs30measured)] += C["phi_s2s_inf"] phi = np.sqrt(phi ** 2. + phi_s2s ** 2.) for stddev_type in stddev_types: assert stddev_type in self.DEFINED_FOR_STANDARD_DEVIATION_TYPES if stddev_type == const.StdDev.TOTAL: stddevs.append(np.sqrt(tau ** 2. + phi ** 2.) + np.zeros(stddev_shape)) elif stddev_type == const.StdDev.INTRA_EVENT: stddevs.append(phi + np.zeros(stddev_shape)) elif stddev_type == const.StdDev.INTER_EVENT: stddevs.append(tau + np.zeros(stddev_shape)) return stddevs
COEFFS = CoeffsTable(sa_damping=5, table="""\ imt e1 b1 b2 b3 c1 c2 c3 tau_c3 phis2s tau_event tau_fault phi0 d0_obs d1_obs phi_s2s_obs d0_inf d1_inf phi_s2s_inf pgv 2.025713592693850 2.406038866414110 0.215611393306395 0.291297120242859 -1.403387613083740 0.281415036690207 -0.358475654018642 0.158122937530244 0.547281437240368 0.403913422762195 0.299427366002708 0.447957473198010 3.39409070 -0.53805299 0.34110665 2.83554791 -0.44449383 0.41721486 pga 4.687702649605810 1.953983695221440 0.185885068933981 -0.241738900149521 -1.478079598457430 0.288086630948598 -0.678379247430291 0.263916141885138 0.626209850856568 0.402258312134260 0.416094806703880 0.466947277758691 2.85338089 -0.45683808 0.40353112 2.09171530 -0.32830709 0.52402952 0.010 4.689953753124180 1.953153318951010 0.186226399234793 -0.242289907409910 -1.480356424151250 0.288439060230840 -0.677333232895566 0.263885400397128 0.626583207016585 0.402108659213375 0.416770119149264 0.466978946978010 2.75935019 -0.44173920 0.40693700 2.04328331 -0.32080733 0.52603004 0.025 4.713377198166330 1.936551887385350 0.189245873356868 -0.242251953373212 -1.527209359618410 0.289740058354006 -0.638666791582353 0.262164920363914 0.631328171215036 0.398276492676289 0.420866218374451 0.468574643779387 2.71223448 -0.43418290 0.41023752 2.02257539 -0.31761109 0.52878573 0.040 4.795848182057760 1.898151028508490 0.201336642438332 -0.237871501230899 -1.624031787587650 0.304689714377880 -0.597052092823935 0.251662690566386 0.653172181191298 0.388148609537982 0.442362727924819 0.473434298534627 2.59943543 -0.41635020 0.41678193 1.97105951 -0.30963894 0.53861744 0.050 4.891963784319470 1.888058117695470 0.209278189591062 -0.269811787471139 -1.650049238074590 0.318075099135608 -0.621868044857147 0.265892288510560 0.667742830430293 0.385327490491663 0.468784102021719 0.479389514380334 2.47698919 -0.39744555 0.42368611 1.90316633 -0.29904349 0.55532505 0.070 5.106825186036690 1.874242004704420 0.212124686343632 -0.320099420113234 -1.620479311265060 0.320046491923063 -0.724343238043036 0.291387067989580 0.693461972699411 0.392248296873683 0.489423169198560 0.485345337361231 2.41434842 -0.38832264 0.43249266 1.82839463 -0.28724310 0.57276719 0.100 5.347455916482890 1.868122162855800 0.199938780857780 -0.356072391432107 -1.517710013159360 0.289475493518785 -0.829503523244518 0.332348696555227 0.695072073756251 0.402824216442605 0.512338505510147 0.493702647513133 2.49026634 -0.40091734 0.44418283 1.78554773 -0.28036580 0.58350938 0.150 5.549588796784530 1.848063879486740 0.152771408659783 -0.348229196634948 -1.358583560220900 0.230015236938864 -0.889955230608062 0.338152565231841 0.682916099464874 0.408594346059644 0.479823906245635 0.496027258004371 2.70826227 -0.43537774 0.45576797 1.83102995 -0.28737956 0.58238166 0.200 5.618412125074440 1.885071394426640 0.129823775224324 -0.292676175821898 -1.293132464925350 0.185865789481763 -0.834473063837879 0.308867853425155 0.669218276129109 0.429871752908997 0.417131887788201 0.494836307784100 2.96015539 -0.47471852 0.46074058 1.97209995 -0.30950033 0.56820594 0.250 5.544966205863690 1.889841953164010 0.105763631540550 -0.185382496096154 -1.250626677741170 0.158041062691118 -0.776587060423500 0.277393830052251 0.646372161163196 0.425139089741650 0.389555883970782 0.488371422580152 3.15218029 -0.50435711 0.45479246 2.16603829 -0.33998731 0.54559779 0.300 5.450353086526580 1.923037058761700 0.093388646106402 -0.133518495442793 -1.226366053882270 0.141154449809013 -0.718787095472858 0.275157669080336 0.625211078118202 0.423146216277508 0.378059931938700 0.482704664228660 3.31253402 -0.52897680 0.44012160 2.37526546 -0.37287433 0.52300268 0.350 5.334813762532070 1.955199629413510 0.086785841616639 -0.028679427285532 -1.204150467410360 0.130192388780628 -0.673786484345681 0.259619616389265 0.620539435099090 0.420137510023815 0.339003245617876 0.481586998560825 3.50716680 -0.55908225 0.42696893 2.59119010 -0.40680160 0.50748852 0.400 5.271854940995910 2.029644803968490 0.091116974902809 -0.047769201478982 -1.178699297309130 0.119724051296491 -0.642239681195656 0.246279727018558 0.625819372379235 0.414328187301189 0.334544067635360 0.475978739636533 3.71718652 -0.59162782 0.42367531 2.80329824 -0.44015273 0.49827621 0.450 5.178258137551410 2.069576283726890 0.090216666792283 -0.039988422268797 -1.167332473566050 0.114113990185286 -0.603205730169031 0.247585522515150 0.629409711420667 0.411924382020225 0.333145123580464 0.469784548562159 3.88263065 -0.61701074 0.42473421 2.99144916 -0.46973419 0.49214988 0.500 5.082700791920750 2.099155311093840 0.086780566158664 0.009747490727000 -1.167604250112560 0.107042637521344 -0.563257403363693 0.239584431802383 0.630142062016198 0.404800890388222 0.327182920069769 0.462936819509800 3.99872441 -0.63424929 0.42061370 3.14615994 -0.49395948 0.49003692 0.600 4.893472350204460 2.164821202949550 0.084969812782495 0.098861244388268 -1.157406472584810 0.100755617814698 -0.497303864164013 0.219211803845992 0.630252473395296 0.405655568634848 0.315141786853113 0.451108093286210 4.09744148 -0.64840075 0.41111777 3.25898021 -0.51148327 0.49106444 0.700 4.705054382687820 2.212434833698710 0.081637278906465 0.100462763700203 -1.152016011593760 0.093631933132715 -0.440618892293181 0.214461662068525 0.627465836385995 0.401597002148329 0.333218982022375 0.447706507848849 4.18377374 -0.66059839 0.40312443 3.32389861 -0.52144948 0.49277756 0.750 4.616038943804310 2.235375898851390 0.079021945207681 0.154139215140058 -1.142973389858090 0.088852396569925 -0.419724044373000 0.207438203391352 0.628069780963942 0.407735734666991 0.326474473374952 0.445225236252333 4.23347646 -0.66735737 0.40261382 3.35343875 -0.52592501 0.49296771 0.800 4.527849452198010 2.264397979442970 0.081390475773282 0.202229318059039 -1.136648599451540 0.090459890537105 -0.403124991422463 0.189995354821562 0.625851513946276 0.415349466205398 0.317303066311374 0.440278408506843 4.23278318 -0.66641099 0.40997990 3.37310303 -0.52894692 0.49030209 0.900 4.391239646181240 2.325276988761640 0.083821056747981 0.228778153653814 -1.127292034561540 0.089501921986240 -0.374788234612032 0.167063633965528 0.622882319152969 0.419993397797855 0.323260336648660 0.434351634884768 4.20215961 -0.66071038 0.41900130 3.40498076 -0.53397430 0.48472128 1.000 4.257838894688260 2.367725674978040 0.084432350092012 0.341875244589753 -1.116475252259780 0.090216898231893 -0.353209184896043 0.152524622009451 0.618724312927540 0.426363619664496 0.306907261128117 0.427075711161247 4.16049170 -0.65331837 0.42199891 3.45895434 -0.54252573 0.47906156 1.200 4.084000623932080 2.535516003279140 0.113291072791050 0.341105428564444 -1.113281899461570 0.101198267860225 -0.307609138051120 0.136322298762975 0.619667399305834 0.432812736959479 0.288816783423553 0.416736878326100 4.09710461 -0.64281410 0.42099900 3.52525368 -0.55299963 0.47780748 1.400 3.858720456780630 2.616885295578470 0.122899590123645 0.433898170105913 -1.123145841250970 0.104462233472800 -0.263665247407548 0.126451619959459 0.624067599637613 0.447550641437210 0.264214688492163 0.408709837205184 4.01062132 -0.62915316 0.42342568 3.58063100 -0.56168209 0.48047546 1.600 3.662163462447190 2.693561483172600 0.134723407599728 0.414339375834985 -1.138808268506760 0.115930381901348 -0.218847134352892 0.131025001289405 0.621514663534016 0.453656810401114 0.252582646468883 0.404631540526549 3.92048752 -0.61543771 0.42442677 3.60901149 -0.56604864 0.48225350 1.800 3.504335982115610 2.804194382503690 0.162346194743908 0.350562313097486 -1.152524904581590 0.120074987873089 -0.185702153948541 0.137842414980722 0.621624707587364 0.466001615468431 0.225283551221296 0.398741991979155 3.82918567 -0.60202524 0.41506055 3.60635912 -0.56555023 0.47843511 2.000 3.352901148637760 2.871065334242860 0.175304278490051 0.301922671452104 -1.163940125236620 0.122626640563880 -0.159248285466522 0.140849131837099 0.614398039046221 0.474830104830787 0.210499929051966 0.396429621856767 3.71030084 -0.58461546 0.39755819 3.57408720 -0.56044621 0.47045811 2.500 3.022521552262530 3.060144403282750 0.230493278721300 0.299896341050341 -1.172681966940500 0.146359630952262 -0.136859516271201 0.166320529763324 0.592031938560107 0.477869423474038 0.197657986712972 0.392547829352281 3.53939609 -0.55908639 0.38202997 3.50749336 -0.54996344 0.46322831 3.000 2.783676662198770 3.235592094294130 0.279996868952325 0.327515809704831 -1.151782848584300 0.155626264164908 -0.169153068756301 0.162858122205747 0.587593509864436 0.484152326322702 0.203965885466933 0.389782912231145 3.33681969 -0.52846169 0.37247458 3.40452938 -0.53375590 0.45984712 3.500 2.550184139271760 3.327248776655720 0.313799508552235 0.427593125870714 -1.157788818464150 0.159502683400810 -0.166147766003538 0.172286178709976 0.573188890082076 0.457777235057253 0.216074174034918 0.387317520863446 3.12496030 -0.49612353 0.36497222 3.30189344 -0.51756844 0.45986473 4.000 2.341796330718510 3.391100460884400 0.334750578303044 0.507655543762405 -1.174149262226720 0.179644555378667 -0.162312356740869 0.174567193920684 0.566535043868093 0.463984870088957 0.202192053366099 0.386005461021883 2.89941316 -0.46132625 0.35761356 3.23185300 -0.50647913 0.46077326 4.500 2.242498493808660 3.476962562487750 0.380318565940050 0.666872895571886 -1.172337191947140 0.199618659323079 -0.172605613836776 0.145926141801792 0.577217595569856 0.446117692446200 0.203666933838210 0.371586179805776 2.67694161 -0.42670283 0.34944100 3.16855635 -0.49644663 0.45945400 5.000 2.009611059381530 3.508361778014590 0.395686780409462 0.861356359764995 -1.208963806141110 0.228641145313669 -0.140697393026636 0.136942909437525 0.559951442166458 0.448040327766279 0.209658464567696 0.377238739334535 2.50349509 -0.39960740 0.34047678 3.05922169 -0.47921741 0.45274837 6.000 1.688862627964000 3.516737645416790 0.420264280758498 1.309501218469920 -1.214089541959290 0.237065702688919 -0.133212839440133 0.141314902795929 0.536431141593929 0.478894088472203 0.212035440719214 0.384111550789384 2.41685314 -0.38612724 0.33122324 2.89854557 -0.45396539 0.44059217 7.000 1.341373165722010 3.465765382316050 0.416036203766797 1.376730174631370 -1.277754681712860 0.267078137770598 -0.088623284188977 0.145375015287528 0.522666557851543 0.427622455537343 0.214045145419707 0.384907261550421 2.39868751 -0.38342544 0.32324116 2.75868843 -0.43197944 0.42858057 8.000 1.063213165919580 3.450720620500530 0.427534350267617 1.454576290963340 -1.323042041332830 0.291715072522946 -0.064115203294262 0.151582553207261 0.507211779986322 0.422147964378804 0.234785059413908 0.386735725945679 2.39997563 -0.38372495 0.32003214 2.70660001 -0.42378341 0.42381284 """)
[docs]class KothaEtAl2019SERASlopeGeology(KothaEtAl2019SERA): """ Adaptation of the Kotha et al. (2019) GMPE for use with slope and geology in place of inferred/measured Vs30. """ experimental = True #: Required site parameter is not set REQUIRES_SITES_PARAMETERS = set(("slope", "geology")) #: Geological Units GEOLOGICAL_UNITS = [b"CENOZOIC", b"HOLOCENE", b"MESOZOIC", b"PALEOZOIC", b"PLEISTOCENE", b"PRECAMBRIAN"]
[docs] def get_site_amplification(self, C, sites, imt): """ Returns the site amplification term depending on whether the Vs30 is observed of inferred """ C_AMP_FIXED = self.COEFFS_FIXED[imt] C_AMP_RAND_INT = self.COEFFS_RANDOM_INT[imt] C_AMP_RAND_GRAD = self.COEFFS_RANDOM_GRAD[imt] ampl = np.zeros(sites.slope.shape) geol_units = np.unique(sites.geology) t_slope = np.copy(sites.slope) t_slope[t_slope > 0.1] = 0.1 for geol_unit in geol_units: idx = sites.geology == geol_unit if geol_unit in self.GEOLOGICAL_UNITS: # Supported geological unit - use the random effects model v1 = C_AMP_FIXED["V1"] + C_AMP_RAND_INT[geol_unit.decode()] v2 = C_AMP_FIXED["V2"] + C_AMP_RAND_GRAD[geol_unit.decode()] else: # Unrecognised geological unit - use the fixed effects model v1 = C_AMP_FIXED["V1"] v2 = C_AMP_FIXED["V2"] ampl[idx] = v1 + v2 * np.log(t_slope[idx]) return ampl
[docs] def get_stddevs(self, C, stddev_shape, stddev_types, sites): """ Returns the ergodic standard deviation with phi_s2s_inf based on that of the inferred Vs30 """ stddevs = [] tau = C["tau_event"] phi = C["phi0"] if self.ergodic: phi = np.sqrt(phi ** 2. + C["phi_s2s_inf"] ** 2.) for stddev_type in stddev_types: assert stddev_type in self.DEFINED_FOR_STANDARD_DEVIATION_TYPES if stddev_type == const.StdDev.TOTAL: stddevs.append(np.sqrt(tau ** 2. + phi ** 2.) + np.zeros(stddev_shape)) elif stddev_type == const.StdDev.INTRA_EVENT: stddevs.append(phi + np.zeros(stddev_shape)) elif stddev_type == const.StdDev.INTER_EVENT: stddevs.append(tau + np.zeros(stddev_shape)) return stddevs
COEFFS_FIXED = CoeffsTable(sa_damping=5, table="""\ imt V1 V2 pgv -0.39341866 -0.13481396 pga -0.29452412 -0.10139168 0.0100 -0.29669217 -0.10286746 0.0250 -0.30159629 -0.10517674 0.0400 -0.31078421 -0.10955506 0.0500 -0.31350486 -0.11143807 0.0700 -0.30101410 -0.10762801 0.1000 -0.27677509 -0.09933628 0.1500 -0.25910646 -0.09270073 0.2000 -0.26216715 -0.09236768 0.2500 -0.28252048 -0.09712923 0.3000 -0.30905284 -0.10368510 0.3500 -0.33949889 -0.11183386 0.4000 -0.37663849 -0.12258381 0.4500 -0.41091170 -0.13210332 0.5000 -0.42752568 -0.13473636 0.6000 -0.42540127 -0.13071676 0.7000 -0.41564728 -0.12480806 0.7500 -0.40812831 -0.12069992 0.8000 -0.40365318 -0.11835608 0.9000 -0.40125176 -0.11666171 1.0000 -0.40229742 -0.11549082 1.2000 -0.42326943 -0.12205893 1.4000 -0.47826523 -0.14350172 1.6000 -0.53714148 -0.16790990 1.8000 -0.56669999 -0.18088943 2.0000 -0.57078751 -0.18356059 2.5000 -0.56937347 -0.18485605 3.0000 -0.56348106 -0.18577307 3.5000 -0.53863711 -0.17950416 4.0000 -0.49270393 -0.16369559 5.0000 -0.42753793 -0.13854620 6.0000 -0.35391502 -0.10915034 7.0000 -0.29887783 -0.08733112 8.0000 -0.27997634 -0.07999381 """) COEFFS_RANDOM_INT = CoeffsTable(sa_damping=5, table="""\ imt PRECAMBRIAN PALEOZOIC MESOZOIC CENOZOIC PLEISTOCENE HOLOCENE pgv -0.04511728 -0.13136056 -0.18001284 -0.03707882 0.15836140 0.23520810 pga 0.01219580 -0.04970639 -0.09890170 -0.03764331 0.07670686 0.09734874 0.0100 0.01533802 -0.05013079 -0.10570863 -0.05550472 0.09391737 0.10208875 0.0250 0.01484921 -0.05600987 -0.11752524 -0.07579728 0.11901269 0.11547049 0.0400 0.01481762 -0.06618308 -0.14016495 -0.11760419 0.17022449 0.13891010 0.0500 0.01869939 -0.07164723 -0.15637671 -0.15427718 0.21552290 0.14807883 0.0700 0.02715475 -0.06977797 -0.16327819 -0.16985317 0.23757146 0.13818312 0.1000 0.03442503 -0.06172966 -0.16251306 -0.15668822 0.22653014 0.11997577 0.1500 0.03354460 -0.05549254 -0.15500108 -0.12214645 0.19114828 0.10794719 0.2000 0.02307337 -0.06041465 -0.14445088 -0.08176772 0.15392561 0.10963428 0.2500 0.00709515 -0.07696600 -0.13396619 -0.05060453 0.13421273 0.12022884 0.3000 -0.01198412 -0.10029723 -0.12913823 -0.03015621 0.13412885 0.13744694 0.3500 -0.03816362 -0.12906927 -0.13326702 -0.01385106 0.14836072 0.16599025 0.4000 -0.07310134 -0.16266592 -0.14379833 -0.00138047 0.17515971 0.20578634 0.4500 -0.10398615 -0.18487344 -0.14805115 0.00729355 0.19560013 0.23401706 0.5000 -0.11344841 -0.17621413 -0.13147360 0.01364461 0.18566877 0.22182277 0.6000 -0.10263361 -0.14653758 -0.10424944 0.01724100 0.15518854 0.18099109 0.7000 -0.09003508 -0.12090367 -0.08458995 0.01925331 0.13071759 0.14555781 0.7500 -0.08754214 -0.10968099 -0.07707821 0.02109986 0.12244465 0.13075683 0.8000 -0.08813501 -0.10247871 -0.07351865 0.02111362 0.11868867 0.12433008 0.9000 -0.08205339 -0.08915012 -0.06531793 0.01822001 0.10724245 0.11105899 1.0000 -0.07090768 -0.07322311 -0.05334984 0.01413841 0.09081701 0.09252520 1.2000 -0.08884488 -0.08551099 -0.05780292 0.01058145 0.11286028 0.10871707 1.4000 -0.16629049 -0.15266807 -0.09785700 0.00735734 0.21188351 0.19757471 1.6000 -0.25949238 -0.23487485 -0.14999422 0.00431025 0.33078413 0.30926707 1.8000 -0.31116827 -0.27808781 -0.17967320 0.00087488 0.39538274 0.37267165 2.0000 -0.32219642 -0.28044321 -0.18468110 -0.00496891 0.40718948 0.38510016 2.5000 -0.32729903 -0.27808434 -0.18869352 -0.01434389 0.41388831 0.39453247 3.0000 -0.33054632 -0.28187680 -0.20225424 -0.02705065 0.42448364 0.41724437 3.5000 -0.30366378 -0.27259313 -0.21193373 -0.03747213 0.40662104 0.41904173 4.0000 -0.24211383 -0.23957031 -0.19941948 -0.03754530 0.34734392 0.37130500 5.0000 -0.16158511 -0.17845904 -0.15328860 -0.02599092 0.24995201 0.26937166 6.0000 -0.08403767 -0.10359760 -0.08914993 -0.01036742 0.13917923 0.14797339 7.0000 -0.03574464 -0.04677124 -0.04075787 -0.00053391 0.06073989 0.06306776 8.0000 -0.02197328 -0.02715535 -0.02481612 0.00194432 0.03569322 0.03630720 """) COEFFS_RANDOM_GRAD = CoeffsTable(sa_damping=5, table="""\ imt PRECAMBRIAN PALEOZOIC MESOZOIC CENOZOIC PLEISTOCENE HOLOCENE pgv -0.00377539 -0.01099219 -0.01506339 -0.00310274 0.01325160 0.01968210 pga 0.00082083 -0.00334544 -0.00665649 -0.00253355 0.00516268 0.00655197 0.0100 0.00143450 -0.00559292 -0.01163716 -0.00789349 0.01214384 0.01154522 0.0250 0.00175995 -0.00840252 -0.01742227 -0.01404592 0.02053739 0.01757337 0.0400 0.00268444 -0.01374949 -0.02902116 -0.02692220 0.03794890 0.02905951 0.0500 0.00446022 -0.01736141 -0.03792810 -0.03833820 0.05327756 0.03588992 0.0700 0.00679173 -0.01756870 -0.04090881 -0.04285452 0.05982109 0.03471921 0.1000 0.00791817 -0.01465675 -0.03749064 -0.03779000 0.05394811 0.02807112 0.1500 0.00651200 -0.01074805 -0.02920993 -0.02589674 0.03891131 0.02043141 0.2000 0.00355278 -0.00791426 -0.01977059 -0.01332566 0.02284373 0.01461400 0.2500 0.00134682 -0.00576354 -0.01173464 -0.00560782 0.01190128 0.00985790 0.3000 0.00028407 -0.00396232 -0.00662483 -0.00261403 0.00662457 0.00629255 0.3500 -0.00153070 -0.00451714 -0.00510216 -0.00095266 0.00570660 0.00639605 0.4000 -0.00460162 -0.00838084 -0.00702825 0.00010529 0.00900041 0.01090500 0.4500 -0.00536110 -0.00966305 -0.00784643 0.00025691 0.01021350 0.01240016 0.5000 0.00049697 -0.00070336 -0.00110620 -0.00060095 0.00058357 0.00132997 0.6000 0.01162399 0.01531601 0.01052193 -0.00270846 -0.01668135 -0.01807213 0.7000 0.02210851 0.02892367 0.02012612 -0.00520967 -0.03167407 -0.03427456 0.7500 0.02868897 0.03559221 0.02503883 -0.00699066 -0.03990509 -0.04242425 0.8000 0.03374459 0.03850612 0.02783394 -0.00791762 -0.04506538 -0.04710166 0.9000 0.04014491 0.04274323 0.03144282 -0.00852165 -0.05211671 -0.05369260 1.0000 0.04789522 0.04905975 0.03524228 -0.00907716 -0.06125321 -0.06186688 1.2000 0.04530281 0.04587299 0.03176464 -0.00861789 -0.05720943 -0.05711311 1.4000 0.02053609 0.02171626 0.01472133 -0.00633121 -0.02502446 -0.02561802 1.6000 -0.01119364 -0.00911637 -0.00576127 -0.00322652 0.01558554 0.01371227 1.8000 -0.02933279 -0.02626479 -0.01694104 -0.00119070 0.03817246 0.03555685 2.0000 -0.03385865 -0.02937003 -0.01934717 -0.00125874 0.04331988 0.04051471 2.5000 -0.03793958 -0.03189031 -0.02177657 -0.00246749 0.04834548 0.04572848 3.0000 -0.04414596 -0.03742074 -0.02685797 -0.00388275 0.05677537 0.05553205 3.5000 -0.04187817 -0.03738244 -0.02887518 -0.00502456 0.05588741 0.05727294 4.0000 -0.02831360 -0.02663928 -0.02209343 -0.00540328 0.03960300 0.04284660 5.0000 -0.00650593 -0.00439066 -0.00351352 -0.00504094 0.00820062 0.01125043 6.0000 0.01744442 0.02257143 0.02058306 -0.00426916 -0.02892987 -0.02739987 7.0000 0.03349592 0.04154954 0.03789116 -0.00364860 -0.05460270 -0.05468531 8.0000 0.03828146 0.04740530 0.04337585 -0.00342909 -0.06242891 -0.06320461 """)