Source code for openquake.hazardlib.gsim.kotha_2020
# -*- coding: utf-8 -*-# vim: tabstop=4 shiftwidth=4 softtabstop=4## Copyright (C) 2014-2025 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:`KothaEtAl2020`, :class:`KothaEtAl2020Site`, :class:`KothaEtAl2020Slope`, :class:`KothaEtAl2020ESHM20`, :class:`KothaEtAl2020ESHM20SlopeGeology` :class:`KothaEtAl2020regional`"""importosimportnumpyasnpfromscipy.constantsimportgfromshapely.geometryimportPoint,shapefromshapely.preparedimportprepfromopenquake.baselib.generalimportCallableDictfromopenquake.hazardlib.geo.packagerimportfionafromopenquake.hazardlib.gsim.baseimportGMPE,CoeffsTable,add_aliasfromopenquake.hazardlibimportconstfromopenquake.hazardlib.imtimportPGA,PGV,SA,from_stringfromopenquake.hazardlib.gsim.nga_eastimport(get_tau_at_quantile,ITPL,TAU_EXECUTION,TAU_SETUP)DATA_FOLDER=os.path.join(os.path.dirname(__file__),'Kotha_2020')CONSTANTS={"Mref":4.5,"Rref":30.,"Mh":5.7,"h_D10":4.0,"h_10D20":8.0,"h_D20":12.0}# The large-magnitude statistical standard deviation values are taken from data# supplied by Kotha et al. (2020)SIGMA_MU_COEFFS=CoeffsTable(sa_damping=5,table="""\ imt sigma_mu_m8_shallow sigma_mu_m8_intermediate sigma_mu_m8_deep sigma_mu_m7p4_shallow sigma_mu_m7p4_intermediate sigma_mu_m7p4_deep pgv 0.2865 0.2829 0.2814 0.2108 0.2072 0.2057 pga 0.3040 0.3003 0.2986 0.2250 0.2213 0.2197 0.010 0.3039 0.3002 0.2986 0.2250 0.2213 0.2197 0.025 0.3026 0.2988 0.2972 0.2243 0.2205 0.2189 0.040 0.3010 0.2972 0.2955 0.2241 0.2203 0.2186 0.050 0.3053 0.3014 0.2997 0.2278 0.2239 0.2222 0.070 0.3133 0.3093 0.3076 0.2340 0.2301 0.2284 0.100 0.3219 0.3179 0.3162 0.2403 0.2364 0.2346 0.150 0.3199 0.3159 0.3141 0.2377 0.2337 0.2319 0.200 0.3174 0.3134 0.3117 0.2343 0.2303 0.2285 0.250 0.3118 0.3078 0.3061 0.2297 0.2257 0.2240 0.300 0.3094 0.3055 0.3038 0.2275 0.2236 0.2220 0.350 0.3038 0.2999 0.2982 0.2230 0.2191 0.2174 0.400 0.2989 0.2950 0.2933 0.2197 0.2157 0.2140 0.450 0.2964 0.2926 0.2909 0.2180 0.2142 0.2125 0.500 0.2916 0.2878 0.2861 0.2145 0.2106 0.2090 0.600 0.2897 0.2860 0.2844 0.2131 0.2094 0.2078 0.700 0.2888 0.2852 0.2836 0.2124 0.2088 0.2072 0.750 0.2902 0.2867 0.2851 0.2134 0.2098 0.2083 0.800 0.2923 0.2888 0.2873 0.2147 0.2112 0.2097 0.900 0.2948 0.2915 0.2900 0.2165 0.2132 0.2117 1.000 0.2964 0.2932 0.2918 0.2175 0.2142 0.2128 1.200 0.2961 0.2930 0.2917 0.2170 0.2139 0.2126 1.400 0.3019 0.2990 0.2977 0.2211 0.2182 0.2169 1.600 0.3041 0.3013 0.3000 0.2225 0.2197 0.2184 1.800 0.3060 0.3032 0.3020 0.2235 0.2207 0.2195 2.000 0.3094 0.3067 0.3055 0.2258 0.2231 0.2219 2.500 0.3121 0.3095 0.3083 0.2275 0.2249 0.2237 3.000 0.3279 0.3254 0.3243 0.2392 0.2366 0.2355 3.500 0.3256 0.3230 0.3219 0.2378 0.2351 0.2339 4.000 0.3269 0.3243 0.3232 0.2386 0.2359 0.2348 4.500 0.3483 0.3456 0.3444 0.2537 0.2510 0.2498 5.000 0.3525 0.3498 0.3486 0.2567 0.2539 0.2527 6.000 0.3458 0.3422 0.3406 0.2514 0.2478 0.2462 7.000 0.3453 0.3417 0.3402 0.2513 0.2477 0.2461 8.000 0.3428 0.3392 0.3376 0.2497 0.2460 0.2444 """)AVGSA_SIGMA_MU_COEFFS=CoeffsTable(sa_damping=5,table="""\ imt sigma_mu_m8_shallow sigma_mu_m8_intermediate sigma_mu_m8_deep sigma_mu_m7p4_shallow sigma_mu_m7p4_intermediate sigma_mu_m7p4_deep AvgSA(0.050) 0.36237117 0.36194128 0.36153721 0.26814685 0.26770184 0.26728207 AvgSA(0.100) 0.37100348 0.37058454 0.37019017 0.27464169 0.27420650 0.27379506 AvgSA(0.150) 0.37113960 0.37073150 0.37034714 0.27399830 0.27357286 0.27317036 AvgSA(0.200) 0.36853324 0.36813709 0.36776392 0.27144161 0.27102774 0.27063616 AvgSA(0.250) 0.36438013 0.36399264 0.36362762 0.26790615 0.26750087 0.26711746 AvgSA(0.300) 0.36004639 0.35966739 0.35931043 0.26444337 0.26404673 0.26367159 AvgSA(0.400) 0.34960757 0.34924208 0.34889809 0.25636041 0.25597756 0.25561573 AvgSA(0.500) 0.34196320 0.34160761 0.34127316 0.25046920 0.25009644 0.24974445 AvgSA(0.600) 0.33628994 0.33594377 0.33561851 0.24613837 0.24577525 0.24543272 AvgSA(0.700) 0.33265779 0.33232082 0.33200448 0.24327786 0.24292413 0.24259077 AvgSA(0.800) 0.33080518 0.33047781 0.33017078 0.24172055 0.24137672 0.24105304 AvgSA(0.900) 0.33028960 0.32997162 0.32967376 0.24123027 0.24089610 0.24058191 AvgSA(1.000) 0.33118486 0.33087499 0.33058504 0.24174675 0.24142106 0.24111520 AvgSA(1.250) 0.33684608 0.33655872 0.33629100 0.24552957 0.24522646 0.24494296 AvgSA(1.500) 0.34131151 0.34103635 0.34078037 0.24857633 0.24828558 0.24801406 AvgSA(1.750) 0.34468757 0.34442002 0.34417138 0.25099193 0.25070866 0.25044432 AvgSA(2.000) 0.34962212 0.34936230 0.34912098 0.25441729 0.25414208 0.25388544 AvgSA(2.500) 0.35772559 0.35748828 0.35726965 0.25993600 0.25968178 0.25944620 AvgSA(3.000) 0.38024588 0.38006170 0.37989860 0.27645421 0.27624724 0.27606079 AvgSA(3.500) 0.37303807 0.37283697 0.37265703 0.27159333 0.27137037 0.27116809 AvgSA(4.000) 0.37757969 0.37738041 0.37720210 0.27464125 0.27442015 0.27421954 AvgSA(4.500) 0.39949828 0.39927571 0.39907300 0.29012395 0.28988040 0.28965658 AvgSA(5.000) 0.40186133 0.40163856 0.40143559 0.29173110 0.29148734 0.29126324 """)def_get_h(C,hypo_depth):""" Returns the depth-specific coefficient """returnnp.where(hypo_depth<=10.,CONSTANTS["h_D10"],np.where(hypo_depth>20.,CONSTANTS["h_D20"],CONSTANTS["h_10D20"]))get_distance_coefficients=CallableDict()
[docs]@get_distance_coefficients.add("base","site","slope","avgsa_base")defget_distance_coefficients_1(kind,c3,c3_epsilon,C,imt,sctx):""" Returns either the directly specified c3 value or the c3 from the existing tau_c3 distribution """ifc3:# Use the c3 that has been defined on inputreturnc3else:# Define the c3 as a number of standard deviation multiplied# by tau_c3returnC["c3"]+(c3_epsilon*C["tau_c3"])
[docs]@get_distance_coefficients.add("ESHM20","geology","avgsa_ESHM20","avgsa_ESHM20_geology","avgsa_ESHM20_homoskedastic")defget_distance_coefficients_2(kind,c3,c3_epsilon,C,imt,sctx):""" Returns the c3 term. If c3 was input directly into the GMPE then this over-rides the c3 regionalisation. Otherwise the c3 and tau_c3 are determined according to the region to which each site is assigned. Note that no regionalisation is defined for PGV and hence the default values from Kotha et al. (2020) are taken unless defined otherwise in the input c3 """ifc3:# If c3 is input then this over-rides the regionalisation# assumed within this modelreturnc3[imt]["c3"]*np.ones(sctx.region.shape)# Default c3 and tau values to the original GMPE c3 and tauc3_=C["c3"]+np.zeros(sctx.region.shape)tau_c3=C["tau_c3"]+np.zeros(sctx.region.shape)ifnotnp.any(sctx.region)or("PGV"instr(imt)):# No regionalisation - take the default C3 and multiply tau_c3# by the original epsilonreturn(c3_+c3_epsilon*tau_c3)+np.zeros(sctx.region.shape)# Some ctx belong to the calibrated regions - loop through themC3_R=C3_REGIONS_AVGSA[imt]ifkind.startswith("avgsa")elseC3_REGIONS[imt]foriinrange(1,6):idx=sctx.region==ic3_[idx]=C3_R["region_{:s}".format(str(i))]tau_c3[idx]=C3_R["tau_region_{:s}".format(str(i))]returnc3_+c3_epsilon*tau_c3
[docs]defget_distance_coefficients_3(att,delta_c3_epsilon,C,imt,sctx):""" Return site-specific coefficient 'C3'. The function retrieves the value of delta_c3 and the standard error of delta_c3 from the 'att' geojson file depending on the location of site. This delta_c3 is added to the generic coefficient 'c3' from the GMPE. A delta_c3_epsilon value of +/- 1.6 gives the 95% confidence interval for delta_c3. """s=[(Point(lon,lat))forlon,latinzip(sctx.lon,sctx.lat)]delta_c3=np.zeros((len(sctx.lat),2),dtype=float)fori,featureinenumerate(att):prepared_polygon=prep(shape(feature['geometry']))contained=list(filter(prepared_polygon.contains,s))ifcontained:ll=np.concatenate([np.where((sctx['lon']==p.x)&(sctx['lat']==p.y))[0]forpincontained])delta_c3[ll,0]=feature['properties'][str(imt)]delta_c3[ll,1]=feature['properties'][str(imt)+'_se']returnC["c3"]+delta_c3[:,0]+delta_c3_epsilon*delta_c3[:,1]
[docs]defget_distance_term(kind,c3,c3_epsilon,C,ctx,imt):""" Returns the distance attenuation factor """h=_get_h(C,ctx.hypo_depth)rval=np.sqrt(ctx.rjb**2.+h**2.)rref_val=np.sqrt(CONSTANTS["Rref"]**2.+h**2.)ifkind!='regional':c3=get_distance_coefficients(kind,c3,c3_epsilon,C,imt,ctx)f_r=(C["c1"]+C["c2"]*(ctx.mag-CONSTANTS["Mref"]))*\
np.log(rval/rref_val)+(c3*(rval-rref_val)/100.)returnf_r
[docs]defget_magnitude_scaling(C,mag):""" Returns the magnitude scaling term """d_m=mag-CONSTANTS["Mh"]returnnp.where(mag<=CONSTANTS["Mh"],C["e1"]+C["b1"]*d_m+C["b2"]*d_m**2.0,C["e1"]+C["b3"]*d_m)
[docs]defget_dl2l(tec,ctx,imt,delta_l2l_epsilon):""" Returns rupture source specific delta_l2l values. The method retrieves the delta_l2l and standard error of delta_l2l values. if delta_l2l_epsilon is provided, standard error of delta_c3 will be included. A delta_l2l_epsilon value of +/- 1.6 gives the 95% confidence interval for delta_l2l. """f=[(Point(lon,lat))forlon,latinzip(ctx.hypo_lon,ctx.hypo_lat)]dl2l=np.zeros((len(ctx.hypo_lon),2),dtype=float)fori,featureinenumerate(tec):prepared_polygon=prep(shape(feature['geometry']))contained=list(filter(prepared_polygon.contains,f))ifcontained:ll=np.concatenate([np.where((ctx['hypo_lon']==p.x)&(ctx['hypo_lat']==p.y))[0]forpincontained])dl2l[ll,0]=feature['properties'][str(imt)]dl2l[ll,1]=feature['properties'][str(imt)+'_se']returndl2l[:,0]+delta_l2l_epsilon*dl2l[:,1]
[docs]defget_sigma_mu_adjustment(kind,C,imt,ctx):""" Returns the sigma_mu adjusment factor, which is taken as the maximum of tau_L2L and the sigma_mu. For M < 7.4 the sigma statistical does not exceed tau_L2L at any period or distance. For M > 7.4, sigma_mu is approximately linear up to M 8.0 so we interpolate between the two values and cap sigma statistical at M 8.0 """C_SIG_MU=AVGSA_SIGMA_MU_COEFFS[imt]ifkind.startswith("avgsa")else\
SIGMA_MU_COEFFS[imt]uf=np.full_like(ctx.mag,C_SIG_MU["sigma_mu_m8_intermediate"])lf=np.full_like(ctx.mag,C_SIG_MU["sigma_mu_m7p4_intermediate"])idx=ctx.hypo_depth<10.0uf[idx]=C_SIG_MU["sigma_mu_m8_shallow"]lf[idx]=C_SIG_MU["sigma_mu_m7p4_shallow"]idx=ctx.hypo_depth>=20.0uf[idx]=C_SIG_MU["sigma_mu_m8_deep"]lf[idx]=C_SIG_MU["sigma_mu_m7p4_deep"]adj=np.maximum(C["tau_l2l"],ITPL(ctx.mag,uf,lf,7.4,0.6))# Below M 7.4 tau_L2L is always larger than sigma muadj[ctx.mag<7.4]=C["tau_l2l"]# Cap the sigma mu as the value for M 8.0adj[ctx.mag>=8.0]=np.maximum(C["tau_l2l"],uf[ctx.mag>=8.0])returnadj
[docs]defget_site_amplification(kind,extra,C,ctx,imt):""" Apply the correct site amplification depending on the kind of GMPE """ifkindin{"base","avgsa_base"}:# no site amplificationampl=0.elifkindin{"site","regional"}:# Render with respect to 800 m/s reference Vs30sref=np.log(ctx.vs30/800.)ampl=(C["g0_vs30"]+C["g1_vs30"]*sref+C["g2_vs30"]*(sref**2.))elifkind=="slope":# Render with respect to 0.1 m/m reference slopesref=np.log(ctx.slope/0.1)ampl=(C["g0_slope"]+C["g1_slope"]*sref+C["g2_slope"]*(sref**2.))elifkindin{"ESHM20","avgsa_ESHM20","avgsa_ESHM20_homoskedastic"}:vs30=np.copy(ctx.vs30)vs30[vs30>1100.]=1100.ampl=np.zeros(vs30.shape)# For observed vs30 ctxampl[ctx.vs30measured]=(C["d0_obs"]+C["d1_obs"]*np.log(vs30[ctx.vs30measured]))# For inferred Vs30 ctxidx=np.logical_not(ctx.vs30measured)ampl[idx]=(C["d0_inf"]+C["d1_inf"]*np.log(vs30[idx]))elifkindin{"geology","avgsa_ESHM20_geology"}:C_AMP_FIXED=extra['COEFFS_FIXED'][imt]C_AMP_RAND_INT=extra['COEFFS_RANDOM_INT'][imt]C_AMP_RAND_GRAD=extra['COEFFS_RANDOM_GRAD'][imt]ampl=np.zeros(ctx.slope.shape)geol_units=np.unique(ctx.geology)t_slope=np.copy(ctx.slope)t_slope[t_slope>0.3]=0.3# Slope lower than 0.0005 m/m takes value for 0.0005 m/mt_slope[t_slope<0.0005]=0.0005forgeol_unitingeol_units:idx=ctx.geology==geol_unitifgeol_unitinextra['GEOLOGICAL_UNITS']:unit=geol_unit.decode()# Supported geological unit, use the random effects modelv1=C_AMP_FIXED["V1"]+C_AMP_RAND_INT[unit]v2=C_AMP_FIXED["V2"]+C_AMP_RAND_GRAD[unit]else:# Unrecognised geological unit, use the fixed effects modelv1=C_AMP_FIXED["V1"]v2=C_AMP_FIXED["V2"]ampl[idx]=v1+v2*np.log(t_slope[idx])returnampl
[docs]defget_stddevs(kind,ergodic,phi_s2s,C,ctx,imt):""" Returns the homoskedastic standard deviation model """mag=ctx.magifkindin{"ESHM20","geology"}:# Get the heteroskedastic tau and phi0tau=get_tau(imt,mag)phi=get_phi_ss(imt,mag)elifkindin{'avgsa_ESHM20',"avgsa_ESHM20_geology"}:# Get the heteroskedastic tau and phi0 for AvgSAtau,phi=get_heteroskedastic_tau_phi0_avgsa(imt,ctx.mag)else:# Get the homoskedastic tau and phi0tau=C["tau_event_0"]phi=C["phi_0"]ifergodic:ifkindin{'ESHM20',"geology","avgsa_ESHM20_geology","avgsa_ESHM20","avgsa_ESHM20_homoskedastic"}:# phi_s2s retrieved in the compute() function of the GMMphi=np.sqrt(phi**2.+phi_s2s**2.)elifkindin{"site","regional"}:phi=np.sqrt(phi**2.0+C["phi_s2s_vs30"]**2.)elifkind=='slope':phi=np.sqrt(phi**2.0+C["phi_s2s_slope"]**2.)else:phi=np.sqrt(phi**2.+C["phis2s"]**2.)return[np.sqrt(tau**2.+phi**2.),tau,phi]
[docs]classKothaEtAl2020(GMPE):""" Implements the first complete version of the newly derived GMPE for Shallow Crustal regions using the Engineering Strong Motion Flatfile. Kotha, S. R., Weatherill, G., Bindi, D., Cotton F. (2020) "A regionally- adaptable ground-motion model for shallow crustal earthquakes in Europe. Bulletin of Earthquake Engineering, 18:4091-4125 The GMPE is desiged for regional adaptation within a logic-tree framework, and as such contains several parameters that can be calibrated on input: 1) Source-region scaling, a simple scalar factor that defines how much to increase or decrease the "regional average" ground motion in the region. This value is taken as the maximum of the source-region variability term (tau_l2l) and the statistical uncertainty (sigma_mu). The latter defines the within-model uncertainty owing to the data set from which the model is derived and only exceeds the former at large magnitudes 2) Residual attenuation scaling "c3", a factor that controls the residual attenuation part of the model to make the ground motion decay more or less rapidly with distance than the regional average. Both factors are period dependent. The two adaptable factors can be controlled either by direct specification at input (in the form of an imt-dependent dictionary) or by a number of standard deviations multiplying the existing variance terms. The two approaches are mutually exclusive, with the directly specified parameters always being used if defined on input. In the core form of the GMPE no site term is included. This is added in the subclasses. :param float sigma_mu_epsilon: Parameter to control the source-region scaling as a number of standard deviations by which to multiply the source-region to source- region variance, max(tau_l2l, sigma_mu) :param float c3_epsilon: Parameter to control the residual attenuation scaling as a number of standard deviations by which to multiply the attenuation-region variance, tau_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 bool ergodic: Use the ergodic standard deviation (True) or non-ergodic standard deviation (False) :param dict dl2l: If specifying the source-region scaling directly, defines the increase or decrease of the ground motion in the form of an imt- dependent dictionary of delta L2L factors :param dict c3: If specifying the residual attenuation scaling directly, defines the apparent anelastic attenuation term, c3, as an imt-dependent dictionary """kind="base"#: 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={PGA,PGV,SA}#: Supported intensity measure component is the geometric mean of two#: horizontal componentsDEFINED_FOR_INTENSITY_MEASURE_COMPONENT=const.IMC.RotD50#: Supported standard deviation types are inter-event, intra-event#: and totalDEFINED_FOR_STANDARD_DEVIATION_TYPES={const.StdDev.TOTAL,const.StdDev.INTER_EVENT,const.StdDev.INTRA_EVENT}#: Required site parameter is not setREQUIRES_SITES_PARAMETERS=set()#: Required rupture parameters are magnitude and hypocentral depthREQUIRES_RUPTURE_PARAMETERS={'mag','hypo_depth'}#: Required distance measure is Rjb (eq. 1).REQUIRES_DISTANCES={'rjb'}def__init__(self,sigma_mu_epsilon=0.0,c3_epsilon=0.0,ergodic=True,dl2l=None,c3=None):""" Instantiate setting the sigma_mu_epsilon and c3 terms """self.sigma_mu_epsilon=sigma_mu_epsilonself.c3_epsilon=c3_epsilonself.ergodic=ergodicifdl2l:# Check that the input is a dictionary and pifnotisinstance(dl2l,dict):raiseIOError("For Kotha et al. (2020) GMM, source-region ""scaling parameter (dl2l) must be input in the ""form of a dictionary, if specified")self.dl2l={}forkeyindl2l:self.dl2l[from_string(key)]={"dl2l":dl2l[key]}self.dl2l=CoeffsTable.fromdict(self.dl2l)else:self.dl2l=Noneifc3:ifnotisinstance(c3,dict):raiseIOError("For Kotha et al. (2020) GMM, residual ""attenuation scaling (c3) must be input in the ""form of a dictionary, if specified")self.c3={}forkeyinc3:self.c3[from_string(key)]={"c3":c3[key]}self.c3=CoeffsTable.fromdict(self.c3)else:self.c3=None
[docs]defcompute(self,ctx:np.recarray,imts,mean,sig,tau,phi):""" See :meth:`superclass method <.base.GroundShakingIntensityModel.compute>` for spec of input and result values. """form,imtinenumerate(imts):C=self.COEFFS[imt]extra={}ifself.kindin{'ESHM20',"avgsa_ESHM20","avgsa_ESHM20_homoskedastic"}:phi_s2s=np.zeros(ctx.vs30measured.shape,dtype=float)phi_s2s[ctx.vs30measured]+=C["phi_s2s_obs"]phi_s2s[np.logical_not(ctx.vs30measured)]+=C["phi_s2s_inf"]elifself.kindin{'geology',"avgsa_ESHM20_geology"}:phi_s2s=self.COEFFS_FIXED[imt]["phi_s2s"]extra['COEFFS_FIXED']=self.COEFFS_FIXEDextra['COEFFS_RANDOM_INT']=self.COEFFS_RANDOM_INTextra['COEFFS_RANDOM_GRAD']=self.COEFFS_RANDOM_GRADextra['GEOLOGICAL_UNITS']=self.GEOLOGICAL_UNITSelse:phi_s2s=Noneifself.kind=='regional':c3=get_distance_coefficients_3(self.att,self.delta_c3_epsilon,C,imt,ctx)else:c3=self.c3fp=get_distance_term(self.kind,c3,self.c3_epsilon,C,ctx,imt)mean[m]=(get_magnitude_scaling(C,ctx.mag)+fp+get_site_amplification(self.kind,extra,C,ctx,imt))# GMPE originally in cm/s/s - convert to gifimt.string.startswith(('PGA','SA','AvgSA')):mean[m]-=np.log(100.0*g)sig[m],tau[m],phi[m]=get_stddevs(self.kind,self.ergodic,phi_s2s,C,ctx,imt)ifself.dl2l:# The source-region parameter is specified explicitymean[m]+=self.dl2l[imt]["dl2l"]elifself.kind=='regional':dl2l=get_dl2l(self.tec,ctx,imt,self.delta_l2l_epsilon)mean[m]+=dl2lelifself.sigma_mu_epsilon:# epistemic uncertainty factor (sigma_mu) multiplied by# the number of standard deviationssigma_mu=get_sigma_mu_adjustment(self.kind,C,imt,ctx)mean[m]+=self.sigma_mu_epsilon*sigma_mu
[docs]classKothaEtAl2020regional(KothaEtAl2020):""" Adaptation of the Kotha et al. (2020) GMPE using the source and site specific adjustments. """experimental=True#: Required rupture parameters are magnitude, hypocentral locationREQUIRES_RUPTURE_PARAMETERS={'mag','hypo_lat','hypo_lon','hypo_depth'}#: Required site parameter are vs30, lat and lon of the siteREQUIRES_SITES_PARAMETERS={'vs30','lat','lon'}kind="regional"def__init__(self,delta_l2l_epsilon=0.0,delta_c3_epsilon=0.0,ergodic=True,c3=None,dl2l=None):""" Instantiate setting the dl2l and c3 terms. """super().__init__()# importantself.delta_l2l_epsilon=delta_l2l_epsilonself.delta_c3_epsilon=delta_c3_epsilonself.ergodic=ergodicattenuation_file=os.path.join(DATA_FOLDER,'kotha_attenuation_regions.geojson')self.att=list(fiona.open(attenuation_file))tectonic_file=os.path.join(DATA_FOLDER,'kotha_tectonic_regions.geojson')self.tec=list(fiona.open(tectonic_file))
[docs]classKothaEtAl2020Site(KothaEtAl2020):""" Preliminary adaptation of the Kotha et al. (2020) GMPE using a polynomial site amplification function dependent on Vs30 (m/s) """#: Required site parameter is not setREQUIRES_SITES_PARAMETERS={"vs30"}kind="site"
[docs]classKothaEtAl2020Slope(KothaEtAl2020):""" Preliminary adaptation of the Kotha et al. (2020) GMPE using a polynomial site amplification function dependent on slope (m/m) """#: Required site parameter is not setREQUIRES_SITES_PARAMETERS={"slope"}kind="slope"
[docs]defget_tau(imt,mag):""" Heteroskedastic Tau model adopts the "global" model from Al Atik (2015) """tau_model=TAU_SETUP["global"]tau=get_tau_at_quantile(tau_model["MEAN"],tau_model["STD"],None)returnTAU_EXECUTION["global"](imt,mag,tau)
[docs]defget_phi_ss(imt,mag):""" Returns the single station phi (or it's variance) for a given magnitude and intensity measure type according to equation 5.14 of Al Atik (2015) with coefficients calibrated on the ESM data set and Kotha et al. (2020) GMPE """C=HETERO_PHI0[imt]phi=C["a"]+(mag-5.0)*((C["b"]-C["a"])/1.5)phi[mag<=5.0]=C["a"]phi[mag>6.5]=C["b"]returnphi
[docs]classKothaEtAl2020ESHM20(KothaEtAl2020):""" Adaptation of the Kotha et al. (2020) GMPE for application to the 2020 European Seismic Hazard Model, as described in Weatherill et al. (2020) Weatherill, G., Kotha, S. R. and Cotton, F. (2020) "A regionally-adaptable 'scaled-backbone' ground motion logic tree for shallow seismicity in Europe: application to the 2020 European seismic hazard model". Bulletin of Earthquake Engineering, 18:5087 - 5117 There are three key adaptations of the original Kotha et al. (2020) GMM: 1) The use of the residual attenuation regions, which represent the five main sub-regions of Europe with similar attenuation characteristics. The assignment to a particular group is now a site-dependent property, requiring the definition of the "eshm20_region", an integer value between 0 and 5 indicating the residual attenuation region to which the site belongs (1 - 5) or else the default values (0). For each region an expected c3 and variance, tau_c3, are defined from which the resulting c3 is taken as a multiple of the number of standard deviations of tau_c3. 2) The site amplification is defined using a two-segment piecewise linear 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. 3) A magnitude-dependent heteroskedastic aleatroy uncertainty model is used for the region-corrected between-event residuals and the site- corrected within event residuals. The former taken from the "global" tau model of Al Atik (2015), while the later is adapted from the "global" phi0 model of Al Atik (2015) adapted to the distribution of site-corrected within-event residuals determined by the original regression of Kotha et al. (2020). Al Atik, L. (2015) NGA-East: Ground-Motion Standard Deviation Models for Central and Eastern North America, PEER Technical Report, No 2015/07 """#: Required site parameters are vs30, vs30measured and the eshm20_regionREQUIRES_SITES_PARAMETERS=set(("region","vs30","vs30measured"))kind="ESHM20"COEFFS=CoeffsTable(sa_damping=5,table="""\ imt e1 b1 b2 b3 c1 c2 c3 tau_c3 phi_s2s tau_event_0 tau_l2l phi_0 d0_obs d1_obs phi_s2s_obs d0_inf d1_inf phi_s2s_inf pgv 1.11912161648479 2.55771078860152 0.353267224391297 0.879839839344054 -1.41931258132547 0.2706807258213520 -0.304426142175370 0.178233997535235 0.560627759977840 0.422935885699239 0.258560350227890 0.446525247049620 3.30975201 -0.53326451 0.36257068 2.78401517 -0.43790954 0.42677529 pga 3.93782347219377 2.06573167101440 0.304988012209292 0.444773874960317 -1.49787542346412 0.2812414746313380 -0.609876182476899 0.253818777234181 0.606771946180224 0.441761487685862 0.355279206886721 0.467151252053241 2.65261454 -0.43301831 0.38806156 1.88258216 -0.29656277 0.51606938 0.010 3.94038760011295 2.06441772899445 0.305294151898347 0.444352974827805 -1.50006146971318 0.2816120431678390 -0.608869451197394 0.253797652143759 0.607030265833062 0.441635449735044 0.356047209347534 0.467206938011971 2.56961762 -0.41981270 0.40044760 1.82057082 -0.28687880 0.51867018 0.025 3.97499686979384 2.04519749120013 0.308841647142436 0.439374383710060 -1.54376149680542 0.2830031280602480 -0.573207556417252 0.252734624432000 0.610030865927204 0.437676505154608 0.368398604288111 0.468698397037258 2.52820436 -0.41328371 0.40623719 1.79206766 -0.28244435 0.52160624 0.040 4.08702279605872 1.99149766561616 0.319673428428720 0.418531185104657 -1.63671359040283 0.2984823762486280 -0.535139204130152 0.244894143623498 0.626413180170373 0.429637401735540 0.412921240156940 0.473730661220076 2.42784360 -0.39762162 0.41977221 1.72300482 -0.27169228 0.53093819 0.050 4.18397570399970 1.96912968528742 0.328982074841989 0.389853296189063 -1.66358950776148 0.3121928913488560 -0.555191107011420 0.260330694464557 0.638967955474841 0.433639923327438 0.444324049044753 0.479898166019243 2.30956730 -0.37937894 0.43465421 1.64224336 -0.25906654 0.54404664 0.070 4.38176649786342 1.92450788134500 0.321182873495225 0.379581373255289 -1.64352914575492 0.3138101953091510 -0.641089475725666 0.286976037026550 0.661064599433347 0.444338223383705 0.470938801038256 0.487060899687138 2.21859665 -0.36551691 0.44921838 1.56920377 -0.24754055 0.55532276 0.100 4.60722959404894 1.90125096928647 0.298805051330753 0.393002352641809 -1.54339428982169 0.2849395739776680 -0.744270750619733 0.321927482439715 0.663309669119995 0.458382304191096 0.478737965504940 0.496152397155402 2.22143266 -0.36624939 0.46432610 1.53915732 -0.24268225 0.56118134 0.150 4.78583314367062 1.92620172077838 0.249893333649662 0.435396192976506 -1.38136438628699 0.2254113422224680 -0.815688997995934 0.322145126407981 0.655406109737959 0.459702777214781 0.414046169030935 0.497805936702476 2.35118737 -0.38662423 0.47703588 1.59963888 -0.25206957 0.55911690 0.200 4.81847463780069 1.97006598187863 0.218722883323200 0.469713318293785 -1.30697558633587 0.1826533194804230 -0.773372802995208 0.301795870071949 0.643585009231006 0.464006126996261 0.321975745683642 0.494075956910651 2.55240529 -0.41806691 0.48025344 1.75423282 -0.27634242 0.54824186 0.250 4.75134747347049 2.01097445156370 0.195062831156806 0.532210412551561 -1.26259484078950 0.1551575007473110 -0.722012122448262 0.274998157533509 0.623240061418664 0.457687642192569 0.293329526713994 0.488950837091220 2.74904047 -0.44882046 0.46891833 1.96527860 -0.30954933 0.53109975 0.300 4.65252285968525 2.09278551802016 0.194929941231544 0.557034893811231 -1.24071282395616 0.1370008066985060 -0.660466290850886 0.260774631679394 0.609748615552919 0.457514283978959 0.266836791529257 0.482157450259502 2.93212957 -0.47759683 0.44983953 2.19913556 -0.34634476 0.51454301 0.350 4.53350897671045 2.14179725762371 0.189511462582876 0.609892595327716 -1.21514531872583 0.1247122464559250 -0.618593385936676 0.254261888951322 0.609506191611413 0.450960093750492 0.231614185359720 0.480254056040507 3.12993498 -0.50873128 0.43569377 2.44212272 -0.38459154 0.50459028 0.400 4.44193244811952 2.22862498827440 0.200305171692326 0.614767001033243 -1.18897228839914 0.1156387616270450 -0.591574546068960 0.243643375298288 0.615477199296824 0.441122908694716 0.240825814626397 0.475193646646757 3.33033435 -0.54013326 0.43045602 2.67707249 -0.42163058 0.50107926 0.450 4.33697728548038 2.29103572171716 0.209573442606565 0.634252522127606 -1.18013993982454 0.1100834686500940 -0.555234498707119 0.245883260391068 0.619384591074073 0.436294164198843 0.249245758570064 0.469672671050266 3.50290267 -0.56696060 0.43223316 2.88578405 -0.45456492 0.50146998 0.500 4.23507897753587 2.35399193121686 0.218088423514177 0.658541873692286 -1.17726165949601 0.1026978146186720 -0.519413341065942 0.238559829231160 0.624993564560933 0.428500398327627 0.243778652813106 0.463165027132890 3.65227902 -0.58990263 0.43887979 3.06576841 -0.48290522 0.50314566 0.600 4.02306439391925 2.42753387249929 0.218787915039312 0.754615594874153 -1.16678688970027 0.0940582863096094 -0.454043559543982 0.216855298090451 0.635090711921061 0.426296731581312 0.246117069779268 0.451206692163190 3.78937389 -0.61070144 0.44724118 3.20894580 -0.50535303 0.50313816 0.700 3.83201580121827 2.51268432884949 0.225024841305000 0.765438564882833 -1.16236278470164 0.0865917976706938 -0.397781532595396 0.215716276719833 0.633635835573626 0.425379430268476 0.246750734502549 0.446704739768374 3.90172707 -0.62754331 0.45268279 3.29999705 -0.51955858 0.50200072 0.750 3.74614211993052 2.55840246083607 0.231604957273506 0.793480645885641 -1.15333203234665 0.0824927940948198 -0.376630503031279 0.209593410875067 0.637877956868669 0.428563811859323 0.245166749142241 0.444311331912854 3.97560847 -0.63847685 0.45583313 3.34616641 -0.52673049 0.50236259 0.800 3.65168809980226 2.59467404437385 0.237334498546207 0.828241777740572 -1.14645090256437 0.0837439530041729 -0.363246464853852 0.192106714053294 0.638753820813416 0.433880652259324 0.240072953116796 0.439300059540554 4.01969394 -0.64478309 0.46384687 3.37966751 -0.53196741 0.50266660 0.900 3.51228638217709 2.68810225072750 0.251716558693382 0.845561170244942 -1.13599614124436 0.0834018259445213 -0.333908265367165 0.177456610405390 0.640328521929993 0.438913972406961 0.247662698012904 0.433043490235851 4.05410191 -0.64939631 0.47448247 3.42678904 -0.53940883 0.49912472 1.000 3.36982044793917 2.74249776483975 0.256784133033388 0.896648260528882 -1.12443352348542 0.0854384622609198 -0.317465939881623 0.171997778367260 0.638429444564638 0.444086895369946 0.238111905941701 0.426703815544157 4.07365692 -0.65153510 0.48134887 3.49473194 -0.55015995 0.49404787 1.200 3.10224418952824 2.82683484364226 0.262683442221073 0.982921357727718 -1.12116148624672 0.0973231293288241 -0.275616235541070 0.160445653296358 0.640086303643832 0.446121165446841 0.226825215617356 0.416539877732589 4.05048971 -0.64704214 0.48708350 3.57270165 -0.56244631 0.49375397 1.400 2.84933745949861 2.89911332547612 0.272065572034688 1.040000637056720 -1.12848926976065 0.1002887249133400 -0.234977212668109 0.150949141990859 0.649359928046388 0.457011583377380 0.231922092201736 0.409641113489270 3.99349305 -0.63756820 0.49596280 3.64615783 -0.57391983 0.49885402 1.600 2.63503429015231 2.98365736561984 0.289670716036571 1.073002118658300 -1.14064711059980 0.1100788214866130 -0.198050139347725 0.148738498099927 0.650540540696659 0.462781403376806 0.223897549097876 0.404985162254916 3.94048869 -0.62914699 0.50237219 3.70614492 -0.58319956 0.50427003 1.800 2.43032254290751 3.06358840071518 0.316828766785138 1.109809835991900 -1.15419967841818 0.1131278831612640 -0.167123738873435 0.156141593013035 0.656949311785981 0.468432106332010 0.205207971335941 0.399057812399511 3.90126474 -0.62332928 0.49599967 3.73733460 -0.58797931 0.50406486 2.000 2.24716354703519 3.11067747935049 0.326774527695550 1.132479221218060 -1.16620971948721 0.1162990300931710 -0.140731664063789 0.155054491423268 0.647763389017009 0.476577198889343 0.196850466599025 0.396502973620567 3.84084468 -0.61459972 0.47661567 3.71781492 -0.58487198 0.49679447 2.500 1.83108464781202 3.23289020747997 0.374214285707986 1.226390493979360 -1.17531326311999 0.1395412164588280 -0.120745041347963 0.176744551716694 0.629481669044830 0.479859874942997 0.190867925368865 0.393288023064441 3.71684077 -0.59605682 0.44991701 3.63149526 -0.57133201 0.48588889 3.000 1.58259215964414 3.44640772476285 0.454951810817816 1.313954219909490 -1.15664484431459 0.1494902905791280 -0.149050671035371 0.174876785480317 0.616446588503561 0.488309107285476 0.220914253465451 0.390859427279163 3.54176439 -0.56936072 0.42220113 3.49013277 -0.54916732 0.47625314 3.500 1.32153652077149 3.56445182133655 0.518610571029448 1.394984393379380 -1.16368470057735 0.1543445278711660 -0.142873831246493 0.193619214137258 0.600202108018105 0.479187019962682 0.237281350236338 0.388102875218375 3.34546112 -0.53906501 0.39951709 3.34520093 -0.52645323 0.47012445 4.000 1.10607064193676 3.64336885536264 0.555331865800278 1.418144933323620 -1.17757508691221 0.1730832048262120 -0.142053716741244 0.193571789393738 0.593046407283143 0.482524831704549 0.233827536969510 0.386956009422453 3.13392178 -0.50620694 0.38303088 3.23169516 -0.50870031 0.46555128 4.500 1.05987610378773 3.82152567982841 0.666476453600402 1.430548279466630 -1.17323633891422 0.1936210609543320 -0.156076448842833 0.152553585766189 0.581331910387036 0.456765160173852 0.196697785051230 0.372827866334900 2.90740942 -0.47082887 0.36840706 3.13020974 -0.49278809 0.46035806 5.000 0.82373381739570 3.84747968562771 0.684665144355361 1.496536314224210 -1.20969230916539 0.2213041109459350 -0.126052481240424 0.137919529808920 0.558954997903623 0.464229101930025 0.195572800413952 0.377458812369736 2.68344324 -0.43562070 0.35254196 2.99932475 -0.47213713 0.45347349 6.000 0.50685354955206 3.80040950285788 0.700805222359295 1.625591116375650 -1.22440411739130 0.2292764533844400 -0.113766839623945 0.141669390606605 0.538973145096788 0.439059204276786 0.190680023411634 0.384862538848542 2.50354874 -0.40714992 0.33854229 2.83412987 -0.44598168 0.44328149 7.000 0.19675504234642 3.78431011962409 0.716569352050671 1.696310364814470 -1.28517895409644 0.2596896867469380 -0.070585399916418 0.146488759166368 0.523331606096182 0.434396029381517 0.208231539543981 0.385850838707000 2.39499327 -0.38989994 0.33074643 2.69365804 -0.42370171 0.43214765 8.000 -0.08979569600589 3.74815514351616 0.726493405776986 1.695347146909250 -1.32882937608962 0.2849197966362740 -0.051296439369391 0.150981191615944 0.508537123776905 0.429104860654150 0.216201318346277 0.387633769846605 2.35979253 -0.38432385 0.32874669 2.64017872 -0.41521615 0.42722298 """)
[docs]defget_heteroskedastic_tau_phi0_avgsa(imt,mag):"""Returns the heteroskedastic between-event and single-station within- event variability for AvgSa """C=COEFFS_AVGSA_TAU_PHI[imt]tau=np.full_like(mag,C["tau1"])tau[mag>6.5]=C["tau4"]idx=(mag>5.5)&(mag<=6.5)tau[idx]=ITPL(mag[idx],C["tau4"],C["tau3"],5.5,1.0)idx=(mag>5.0)&(mag<=5.5)tau[idx]=ITPL(mag[idx],C["tau3"],C["tau2"],5.0,0.5)idx=(mag>4.5)&(mag<=5.0)tau[idx]=ITPL(mag[idx],C["tau2"],C["tau1"],4.5,0.5)phi0=C["phi0_a"]+(mag-5.0)*((C["phi0_b"]-C["phi0_a"])/1.5)phi0[mag<=5.0]=C["phi0_a"]phi0[mag>6.5]=C["phi0_b"]returntau,phi0