# -*- coding: utf-8 -*-
# vim: tabstop=4 shiftwidth=4 softtabstop=4
#
# Copyright (C) 2012-2023 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:`BooreEtAl2014`,
               :class:`BooreEtAl2014HighQ`,
               :class:`BooreEtAl2014LowQ`
"""
import numpy as np
from openquake.baselib.general import CallableDict
from openquake.hazardlib.gsim.base import GMPE, CoeffsTable, add_alias
from openquake.hazardlib.gsim.abrahamson_2014 import get_epistemic_sigma
from openquake.hazardlib import const
from openquake.hazardlib.imt import PGA, PGV, SA
#: Equation constants that are IMT-independent
CONSTS = {
    "Mref": 4.5,
    "Rref": 1.0,
    "Vref": 760.0,
    "f1": 0.0,
    "f3": 0.1,
    "v1": 225.0,
    "v2": 300.0}
def _get_basin_depth_term(region, C, ctx, period):
    """
    In the case of the base model the basin depth term is switched off.
    Therefore we return an array of zeros.
    """
    if region == "nobasin" or period < 0.65:  # switched off
        return np.zeros(len(ctx.vs30), dtype=float)
    bmodel = (japan_basin_model(ctx.vs30) if region == "JPN"
              else california_basin_model(ctx.vs30))
    f_dz1 = C["f7"] + np.zeros(len(ctx.vs30), dtype=float)
    f_ratio = C["f7"] / C["f6"]
    dz1 = (ctx.z1pt0 / 1000.0) - bmodel
    idx = dz1 <= f_ratio
    f_dz1[idx] = C["f6"] * dz1[idx]
    return f_dz1
def _get_inter_event_tau(C, mag):
    """
    Returns the inter-event standard deviation (tau), which is dependent
    on magnitude
    """
    tau = np.full_like(mag, C["tau1"] + (C["tau2"] - C["tau1"]) * (mag - 4.5))
    tau[mag <= 4.5] = C["tau1"]
    tau[mag >= 5.5] = C["tau2"]
    return tau
_get_intra_event_phi = CallableDict()
@_get_intra_event_phi.add("base")
def _get_intra_event_phi_base(kind, C, mag, rjb, vs30):
    """
    Returns the intra-event standard deviation (phi), dependent on
    magnitude, distance and vs30
    """
    # Magnitude Dependent phi (Equation 17)
    base_vals = np.full_like(rjb, C["f1"] + (C["f2"] - C["f1"]) * (mag - 4.5))
    base_vals[mag <= 4.5] = C["f1"]
    base_vals[mag >= 5.5] = C["f2"]
    # Distance dependent phi (Equation 16)
    idx1 = rjb > C["R2"]
    base_vals[idx1] += C["DfR"]
    idx2 = np.logical_and(rjb > C["R1"], rjb <= C["R2"])
    base_vals[idx2] += C["DfR"] * (np.log(rjb[idx2] / C["R1"]) /
                                   np.log(C["R2"] / C["R1"]))
    # Site-dependent phi (Equation 15)
    idx1 = vs30 <= CONSTS["v1"]
    base_vals[idx1] -= C["DfV"]
    idx2 = np.logical_and(vs30 >= CONSTS["v1"], vs30 <= CONSTS["v2"])
    base_vals[idx2] -= C["DfV"] * (np.log(CONSTS["v2"] / vs30[idx2]) /
                                   np.log(CONSTS["v2"] / CONSTS["v1"]))
    return base_vals
@_get_intra_event_phi.add("stewart")
def _get_intra_event_phi_stewart(kind, C, mag, rjb, vs30):
    """
    Returns the intra-event standard deviation (phi)
    In SBSA15, the intra-event standard deviation only depends on magnitude
    """
    phi = np.full_like(mag, C["phi1"] + (C["phi2"] - C["phi1"]) * (mag - 4.5))
    phi[mag <= 4.5] = C["phi1"]
    phi[mag >= 5.5] = C["phi2"]
    return phi
def _get_linear_site_term(C, vs30):
    """
    Returns the linear site scaling term (equation 6)
    """
    flin = vs30 / CONSTS["Vref"]
    flin[vs30 > C["Vc"]] = C["Vc"] / CONSTS["Vref"]
    return C["c"] * np.log(flin)
def _get_magnitude_scaling_term(sof, C, ctx):
    """
    Returns the magnitude scling term defined in equation (2)
    """
    mag_term = np.zeros_like(ctx.mag)
    below = ctx.mag <= C["Mh"]
    dmag = ctx.mag - C["Mh"]
    mag_term[below] = C["e4"] * dmag[below] + C["e5"] * dmag[below] ** 2.0
    mag_term[~below] = C["e6"] * dmag[~below]
    if not sof:  # unspecified style-of-faulting
        return C["e0"] + mag_term
    return _get_style_of_faulting_term(C, ctx) + mag_term
def _get_nonlinear_site_term(C, vs30, pga_rock):
    """
    Returns the nonlinear site scaling term (equation 7)
    """
    v_s = np.copy(vs30)
    v_s[vs30 > 760.] = 760.
    # Nonlinear controlling parameter (equation 8)
    f_2 = C["f4"] * (np.exp(C["f5"] * (v_s - 360.)) - np.exp(C["f5"] * 400.))
    fnl = CONSTS["f1"] + f_2 * np.log((pga_rock + CONSTS["f3"]) / CONSTS["f3"])
    return fnl
_get_path_scaling = CallableDict()
@_get_path_scaling.add("base")
def _get_path_scaling_base(kind, region, C, ctx):
    """
    Returns the path scaling term given by equation (3)
    """
    rval = np.sqrt(ctx.rjb ** 2.0 + C["h"] ** 2.0)
    scaling = (C["c1"] + C["c2"] * (ctx.mag - CONSTS["Mref"])) * np.log(
        rval / CONSTS["Rref"])
    return scaling + ((C["c3"] + C["Dc3"]) * (rval - CONSTS["Rref"]))
@_get_path_scaling.add("stewart")
def _get_path_scaling_stewart(kind, region, C, ctx):
    """
    Returns the path scaling term defined in Equation 3
    """
    # Calculate R in Equation 4
    rval = np.sqrt(ctx.rjb ** 2.0 + C["h"] ** 2.0)
    if region == "CAL":
        delta_c3 = 0
    elif region == "CHN":
        delta_c3 = C['Dc3CH']
    elif region == "JPN":
        delta_c3 = C['Dc3JP']
    else:
        raise ValueError("region=%s" % region)
    # Calculate geometric spreading component of path scaling term
    fp_geom = (C["c1"] + C["c2"] * (ctx.mag - CONSTS["Mref"])) * np.log(
        rval / CONSTS["Rref"])
    # Calculate anelastic attenuation component of path scaling term, with
    # delta c3 accounting for regional effects
    fp_atten = (C["c3"] + delta_c3) * (rval - CONSTS["Rref"])
    return fp_geom + fp_atten
# imported by other GMPEs
def _get_pga_on_rock(kind, region, sof, C, ctx):
    """
    Returns the median PGA on rock, which is a sum of the
    magnitude and distance scaling
    """
    return np.exp(_get_magnitude_scaling_term(sof, C, ctx) +
                  _get_path_scaling(kind, region, C, ctx))
def _get_site_scaling(kind, region, C, pga_rock, ctx, period, rjb):
    """
    Returns the site-scaling term (equation 5), broken down into a
    linear scaling, a nonlinear scaling and a basin scaling term
    """
    flin = _get_linear_site_term(C, ctx.vs30)
    fnl = _get_nonlinear_site_term(C, ctx.vs30, pga_rock)
    if kind == 'stewart':
        fbd = 0.  # there is no basin term in the Stewart models
    else:
        fbd = _get_basin_depth_term(region, C, ctx, period)
    return flin + fnl + fbd
def _get_stddevs(kind, C, ctx):
    """
    Returns the aleatory uncertainty terms described in equations (13) to
    (17)
    """
    tau = _get_inter_event_tau(C, ctx.mag)
    phi = _get_intra_event_phi(kind, C, ctx.mag, ctx.rjb, ctx.vs30)
    return [np.sqrt(tau**2 + phi**2), tau, phi]
def _get_style_of_faulting_term(C, ctx):
    """
    Get fault type dummy variables
    Fault type (Strike-slip, Normal, Thrust/reverse) is
    derived from rake angle.
    Rakes angles within 30 of horizontal are strike-slip,
    angles from 30 to 150 are reverse, and angles from
    -30 to -150 are normal.
    Note that the 'Unspecified' case is not considered here as
    rake is always given.
    """
    # normal
    res = np.full_like(ctx.rake, C["e2"])
    # strike-slip
    res[(np.abs(ctx.rake) <= 30.) |
        ((180. - np.abs(ctx.rake)) <= 30.)] = C["e1"]
    # reverse
    res[(ctx.rake > 30.) & (ctx.rake < 150.)] = C["e3"]
    return res
[docs]class BooreEtAl2014(GMPE):
    """
    Implements GMPE developed by David M. Boore, Jonathan P. Stewart,
    Emel Seyhan and Gail Atkinson, and published as "NGA-West2 Equations for
    Predicting PGA, PGV, nd 5 % Damped PGA for Shallow Crustal Earthquakes
    (2014, Earthquake Spectra, Volume 30, No. 3, pages 1057 - 1085).
    """
    #: Supported tectonic region type is active shallow crust
    DEFINED_FOR_TECTONIC_REGION_TYPE = const.TRT.ACTIVE_SHALLOW_CRUST
    #: Supported intensity measure types are spectral acceleration,
    #: peak ground velocity and peak ground acceleration
    DEFINED_FOR_INTENSITY_MEASURE_TYPES = {PGA, PGV, SA}
    #: Supported intensity measure component is orientation-independent
    #: measure :attr:`~openquake.hazardlib.const.IMC.RotD50`
    DEFINED_FOR_INTENSITY_MEASURE_COMPONENT = const.IMC.RotD50
    #: Supported standard deviation types are inter-event, intra-event
    #: and total, see equation 2, pag 106.
    DEFINED_FOR_STANDARD_DEVIATION_TYPES = {
        const.StdDev.TOTAL, const.StdDev.INTER_EVENT, const.StdDev.INTRA_EVENT}
    #: Required site parameters is Vs30
    REQUIRES_SITES_PARAMETERS = {'vs30'}
    #: Required rupture parameters are magnitude, and rake.
    REQUIRES_RUPTURE_PARAMETERS = {'mag', 'rake'}
    #: Required distance measure is Rjb
    REQUIRES_DISTANCES = {'rjb'}
    kind = "base"
    def __init__(self, region='nobasin', sof=True, sigma_mu_epsilon=0.0, **kwargs):
        super().__init__(sigma_mu_epsilon = sigma_mu_epsilon, **kwargs)
        self.region = region
        self.sof = sof
        self.sigma_mu_epsilon = sigma_mu_epsilon
        if region != "nobasin":  # z1pt0 is used if period >= 0.65
            self.REQUIRES_SITES_PARAMETERS |= {'z1pt0'}
[docs]    def compute(self, ctx: np.recarray, imts, mean, sig, tau, phi):
        """
        See :meth:`superclass method
        <.base.GroundShakingIntensityModel.compute>`
        for spec of input and result values.
        """
        C_PGA = self.COEFFS[PGA()]
        for m, imt in enumerate(imts):
            C = self.COEFFS[imt]
            pga_rock = _get_pga_on_rock(
                self.kind, self.region, self.sof, C_PGA, ctx)
            mean[m] = (
                _get_magnitude_scaling_term(self.sof, C, ctx) +
                _get_path_scaling(self.kind, self.region, C, ctx) +
                _get_site_scaling(self.kind, self.region,
                                  C, pga_rock, ctx, imt.period, ctx.rjb))
            if self.sigma_mu_epsilon:
                mean[m] += (self.sigma_mu_epsilon*get_epistemic_sigma(ctx))
            sig[m], tau[m], phi[m] = _get_stddevs(self.kind, C, ctx) 
    COEFFS = CoeffsTable(sa_damping=5, table="""\
    IMT            e0          e1          e2          e3         e4          e5          e6         Mh          c1         c2          c3          h        Dc3           c            Vc          f4          f5          f6          f7           R1           R2        DfR        DfV         f1         f2         tau1         tau2
    pgv      5.037000    5.078000    4.849000    5.033000   1.073000   -0.153600    0.225200   6.200000   -1.243000   0.148900   -0.003440   5.300000   0.000000   -0.840000   1300.000000   -0.100000   -0.008440   -9.900000   -9.900000   105.000000   272.000000   0.082000   0.080000   0.644000   0.552000   0.401000   0.346000
    pga      0.447300    0.485600    0.245900    0.453900   1.431000    0.050530   -0.166200   5.500000   -1.134000   0.191700   -0.008088   4.500000   0.000000   -0.600000   1500.000000   -0.150000   -0.007010   -9.900000   -9.900000   110.000000   270.000000   0.100000   0.070000   0.695000   0.495000   0.398000   0.348000
    0.010    0.453400    0.491600    0.251900    0.459900   1.421000    0.049320   -0.165900   5.500000   -1.134000   0.191600   -0.008088   4.500000   0.000000   -0.603720   1500.200000   -0.148330   -0.007010   -9.900000   -9.900000   111.670000   270.000000   0.096000   0.070000   0.698000   0.499000   0.402000   0.345000
    0.020    0.485980    0.523590    0.297070    0.488750   1.433100    0.053388   -0.165610   5.500000   -1.139400   0.189620   -0.008074   4.500000   0.000000   -0.573880   1500.360000   -0.147100   -0.007280   -9.900000   -9.900000   113.100000   270.000000   0.092000   0.030000   0.702000   0.502000   0.409000   0.346000
    0.022    0.498660    0.536470    0.313470    0.499730   1.433600    0.054888   -0.165200   5.500000   -1.140500   0.189240   -0.008095   4.500000   0.000000   -0.566750   1500.680000   -0.148010   -0.007320   -9.900000   -9.900000   113.370000   270.000000   0.088000   0.027000   0.707000   0.505000   0.418000   0.349000
    0.025    0.522830    0.561300    0.344260    0.519990   1.432800    0.057529   -0.164990   5.500000   -1.141900   0.188750   -0.008153   4.500000   0.000000   -0.555200   1501.040000   -0.150150   -0.007360   -9.900000   -9.900000   113.070000   270.000000   0.086000   0.026000   0.711000   0.508000   0.427000   0.354000
    0.029    0.559490    0.599230    0.391460    0.549950   1.427900    0.060732   -0.166320   5.500000   -1.142300   0.188440   -0.008290   4.500000   0.000000   -0.538500   1501.260000   -0.153870   -0.007370   -9.900000   -9.900000   112.360000   270.000000   0.084000   0.028000   0.716000   0.510000   0.436000   0.359000
    0.030    0.569160    0.609200    0.403910    0.557830   1.426100    0.061444   -0.166900   5.500000   -1.142100   0.188420   -0.008336   4.490000   0.000000   -0.534140   1502.950000   -0.154850   -0.007350   -9.900000   -9.900000   112.130000   270.000000   0.081000   0.029000   0.721000   0.514000   0.445000   0.364000
    0.032    0.588020    0.628750    0.427880    0.573300   1.422700    0.062806   -0.168130   5.500000   -1.141200   0.188400   -0.008445   4.450000   0.000000   -0.525290   1503.120000   -0.156850   -0.007310   -9.900000   -9.900000   111.650000   270.000000   0.078000   0.030000   0.726000   0.516000   0.454000   0.369000
    0.035    0.616360    0.658180    0.462520    0.597040   1.417400    0.064559   -0.170150   5.500000   -1.138800   0.188390   -0.008642   4.400000   0.000000   -0.511920   1503.240000   -0.160160   -0.007210   -9.900000   -9.900000   110.640000   270.000000   0.077000   0.031000   0.730000   0.518000   0.462000   0.374000
    0.036    0.625540    0.667720    0.473380    0.604960   1.415800    0.065028   -0.170830   5.500000   -1.137800   0.188370   -0.008715   4.380000   0.000000   -0.507520   1503.320000   -0.161420   -0.007170   -9.900000   -9.900000   109.530000   270.000000   0.075000   0.031000   0.734000   0.520000   0.470000   0.379000
    0.040    0.662810    0.706040    0.515320    0.638280   1.409000    0.066183   -0.173570   5.500000   -1.132400   0.188160   -0.009030   4.320000   0.000000   -0.490650   1503.350000   -0.167770   -0.006980   -9.900000   -9.900000   108.280000   270.000000   0.073000   0.032000   0.738000   0.521000   0.478000   0.384000
    0.042    0.680870    0.724430    0.534450    0.655050   1.405900    0.066438   -0.174850   5.500000   -1.129200   0.187970   -0.009195   4.290000   0.000000   -0.482900   1503.340000   -0.171930   -0.006870   -9.900000   -9.900000   106.990000   270.000000   0.072000   0.032000   0.742000   0.523000   0.484000   0.390000
    0.044    0.698820    0.742770    0.552820    0.672250   1.403300    0.066663   -0.176190   5.500000   -1.125900   0.187750   -0.009360   4.270000   0.000000   -0.475720   1503.130000   -0.176640   -0.006770   -9.900000   -9.900000   105.410000   270.000000   0.070000   0.031000   0.745000   0.525000   0.490000   0.397000
    0.045    0.708220    0.752320    0.562220    0.681390   1.402100    0.066774   -0.176930   5.500000   -1.124200   0.187640   -0.009441   4.250000   0.000000   -0.472360   1502.840000   -0.179140   -0.006720   -9.900000   -9.900000   103.610000   270.000000   0.069000   0.031000   0.748000   0.527000   0.496000   0.405000
    0.046    0.717790    0.762020    0.571660    0.690760   1.400900    0.066891   -0.177690   5.500000   -1.122400   0.187520   -0.009521   4.240000   0.000000   -0.469150   1502.470000   -0.181700   -0.006670   -9.900000   -9.900000   101.700000   270.000000   0.067000   0.031000   0.750000   0.529000   0.499000   0.412000
    0.048    0.735740    0.780150    0.588880    0.708540   1.399100    0.067127   -0.179200   5.500000   -1.119200   0.187300   -0.009676   4.220000   0.000000   -0.463210   1502.010000   -0.186880   -0.006560   -9.900000   -9.900000    99.760000   270.000000   0.065000   0.031000   0.752000   0.530000   0.502000   0.419000
    0.050    0.754360    0.799050    0.606520    0.727260   1.397400    0.067357   -0.180820   5.500000   -1.115900   0.187090   -0.009819   4.200000   0.000000   -0.457950   1501.420000   -0.192000   -0.006470   -9.900000   -9.900000    97.930000   270.000000   0.063000   0.030000   0.753000   0.532000   0.503000   0.426000
    0.055    0.799600    0.844500    0.647700    0.773700   1.394700    0.067797   -0.184800   5.500000   -1.108200   0.186550   -0.010120   4.150000   0.000000   -0.447870   1500.710000   -0.203690   -0.006250   -9.900000   -9.900000    96.030000   270.000000   0.062000   0.029000   0.753000   0.534000   0.502000   0.434000
    0.060    0.843940    0.888840    0.685620    0.820670   1.395400    0.068591   -0.188580   5.500000   -1.100900   0.185820   -0.010330   4.110000   0.000000   -0.441860   1499.830000   -0.213740   -0.006070   -9.900000   -9.900000    94.100000   270.010000   0.061000   0.027000   0.753000   0.536000   0.499000   0.441000
    0.065    0.886550    0.931160    0.719410    0.867240   1.400400    0.070127   -0.191760   5.500000   -1.094200   0.184850   -0.010480   4.080000   0.000000   -0.439510   1498.740000   -0.222250   -0.005930   -9.900000   -9.900000    92.080000   270.020000   0.061000   0.025000   0.752000   0.538000   0.495000   0.448000
    0.067    0.902700    0.947110    0.731710    0.885260   1.403200    0.070895   -0.192910   5.500000   -1.091800   0.184420   -0.010520   4.070000   0.000000   -0.439500   1497.420000   -0.225240   -0.005880   -9.900000   -9.900000    90.010000   270.020000   0.061000   0.025000   0.750000   0.540000   0.489000   0.455000
    0.070    0.926520    0.970570    0.749400    0.912270   1.408200    0.072075   -0.194510   5.500000   -1.088400   0.183690   -0.010560   4.060000   0.000000   -0.440400   1495.850000   -0.229310   -0.005820   -9.900000   -9.900000    87.970000   270.030000   0.062000   0.024000   0.748000   0.541000   0.483000   0.461000
    0.075    0.964470    1.007700    0.776780    0.956300   1.417400    0.073549   -0.196650   5.500000   -1.083100   0.182250   -0.010580   4.040000   0.000000   -0.444110   1494.000000   -0.235000   -0.005730   -9.900000   -9.900000    85.990000   270.040000   0.064000   0.022000   0.745000   0.542000   0.474000   0.466000
    0.080    1.000300    1.042600    0.801610    0.998180   1.426100    0.073735   -0.198160   5.500000   -1.078500   0.180520   -0.010560   4.020000   0.000000   -0.450200   1491.820000   -0.239440   -0.005670   -9.900000   -9.900000    84.230000   270.050000   0.067000   0.020000   0.741000   0.543000   0.464000   0.468000
    0.085    1.034000    1.075500    0.824230    1.037900   1.432200    0.071940   -0.199020   5.510000   -1.074500   0.178560   -0.010510   4.030000   0.000000   -0.458130   1489.290000   -0.242850   -0.005630   -9.900000   -9.900000    82.740000   270.060000   0.072000   0.019000   0.737000   0.543000   0.452000   0.468000
    0.090    1.066600    1.107600    0.845910    1.076200   1.435000    0.068097   -0.199290   5.520000   -1.070900   0.176430   -0.010420   4.070000   0.000000   -0.467320   1486.360000   -0.245440   -0.005610   -9.900000   -9.900000    81.540000   270.070000   0.076000   0.017000   0.734000   0.542000   0.440000   0.466000
    0.095    1.098100    1.138500    0.867030    1.112700   1.433900    0.062327   -0.199000   5.530000   -1.067800   0.174200   -0.010320   4.100000   0.000000   -0.477210   1482.980000   -0.247470   -0.005600   -9.900000   -9.900000    80.460000   270.080000   0.082000   0.016000   0.731000   0.542000   0.428000   0.464000
    0.100    1.126800    1.166900    0.887100    1.145400   1.429300    0.055231   -0.198380   5.540000   -1.065200   0.172030   -0.010200   4.130000   0.000000   -0.487240   1479.120000   -0.249160   -0.005600   -9.900000   -9.900000    79.590000   270.090000   0.087000   0.014000   0.728000   0.541000   0.415000   0.458000
    0.110    1.178500    1.217900    0.927020    1.203000   1.411000    0.037389   -0.196010   5.570000   -1.060700   0.167700   -0.009964   4.190000   0.000000   -0.506320   1474.740000   -0.252130   -0.005620   -9.900000   -9.900000    79.050000   270.110000   0.093000   0.012000   0.726000   0.540000   0.403000   0.451000
    0.120    1.223000    1.262100    0.966160    1.250200   1.383100    0.016373   -0.192650   5.620000   -1.057200   0.163520   -0.009722   4.240000   0.000000   -0.524380   1469.750000   -0.254550   -0.005670   -9.900000   -9.900000    78.850000   270.130000   0.099000   0.011000   0.724000   0.539000   0.392000   0.441000
    0.130    1.259600    1.298600    1.003100    1.286900   1.349700   -0.005158   -0.188980   5.660000   -1.054900   0.159820   -0.009476   4.290000   0.000000   -0.542140   1464.090000   -0.256280   -0.005720   -9.900000   -9.900000    78.990000   270.150000   0.104000   0.011000   0.723000   0.538000   0.381000   0.430000
    0.133    1.269200    1.308200    1.013500    1.296100   1.339500   -0.011354   -0.187920   5.670000   -1.054500   0.158820   -0.009402   4.300000   0.000000   -0.547520   1457.760000   -0.256650   -0.005740   -9.900000   -9.900000    79.470000   270.150000   0.110000   0.011000   0.722000   0.538000   0.371000   0.417000
    0.140    1.288300    1.327000    1.036000    1.313700   1.316200   -0.024711   -0.185660   5.700000   -1.053700   0.156720   -0.009228   4.340000   0.000000   -0.560320   1450.710000   -0.257190   -0.005780   -9.900000   -9.900000    80.260000   270.160000   0.115000   0.012000   0.721000   0.537000   0.362000   0.403000
    0.150    1.309500    1.348100    1.064800    1.332400   1.284400   -0.042065   -0.182340   5.740000   -1.053200   0.154010   -0.008977   4.390000   0.000000   -0.579620   1442.850000   -0.257130   -0.005850   -9.900000   -9.900000    81.330000   270.160000   0.120000   0.015000   0.720000   0.537000   0.354000   0.388000
    0.160    1.323500    1.361500    1.087600    1.343700   1.254100   -0.057593   -0.178530   5.780000   -1.053300   0.151580   -0.008725   4.440000   0.000000   -0.600520   1434.220000   -0.256040   -0.005910   -9.900000   -9.900000    82.860000   270.160000   0.125000   0.020000   0.720000   0.536000   0.349000   0.372000
    0.170    1.330600    1.367900    1.104000    1.348700   1.224400   -0.071861   -0.174210   5.820000   -1.054100   0.149480   -0.008472   4.490000   0.000000   -0.622520   1424.850000   -0.254140   -0.005970   -9.900000   -9.900000    84.720000   270.140000   0.128000   0.026000   0.718000   0.536000   0.346000   0.357000
    0.180    1.332700    1.368900    1.114900    1.349200   1.194100   -0.085640   -0.169390   5.850000   -1.055600   0.147680   -0.008219   4.530000   0.000000   -0.644860   1414.770000   -0.251730   -0.006020   -9.900000   -9.900000    86.670000   270.110000   0.131000   0.033000   0.717000   0.536000   0.344000   0.341000
    0.190    1.330700    1.365600    1.120800    1.346300   1.163500   -0.098884   -0.164040   5.890000   -1.057900   0.146160   -0.007967   4.570000   0.000000   -0.666810   1403.990000   -0.249110   -0.006080   -9.900000   -9.900000    88.730000   270.060000   0.134000   0.039000   0.714000   0.537000   0.343000   0.324000
    0.200    1.325500    1.359000    1.122000    1.341400   1.134900   -0.110960   -0.158520   5.920000   -1.060700   0.144890   -0.007717   4.610000   0.000000   -0.687620   1392.610000   -0.246580   -0.006140   -9.900000   -9.900000    90.910000   270.000000   0.136000   0.045000   0.711000   0.539000   0.344000   0.309000
    0.220    1.309100    1.339400    1.113300    1.328100   1.082300   -0.133000   -0.147040   5.970000   -1.067000   0.142630   -0.007224   4.680000   0.000000   -0.724310   1380.720000   -0.242350   -0.006260   -9.900000   -9.900000    93.040000   269.830000   0.138000   0.052000   0.708000   0.541000   0.345000   0.294000
    0.240    1.288100    1.315000    1.094500    1.313200   1.036600   -0.152990   -0.134450   6.030000   -1.073700   0.140350   -0.006747   4.750000   0.000000   -0.756460   1368.510000   -0.238230   -0.006380   -9.900000   -9.900000    95.080000   269.590000   0.140000   0.055000   0.703000   0.544000   0.347000   0.280000
    0.250    1.276600    1.301700    1.082800    1.305200   1.016600   -0.162130   -0.127840   6.050000   -1.077300   0.139250   -0.006517   4.780000   0.000000   -0.771770   1356.210000   -0.235740   -0.006440   -9.900000   -9.900000    97.040000   269.450000   0.141000   0.055000   0.698000   0.547000   0.350000   0.266000
    0.260    1.265100    1.288600    1.071000    1.297200   0.999320   -0.170410   -0.121150   6.070000   -1.080800   0.138180   -0.006293   4.820000   0.000000   -0.786970   1343.890000   -0.232800   -0.006500   -9.900000   -9.900000    98.870000   269.300000   0.141000   0.055000   0.693000   0.550000   0.353000   0.255000
    0.280    1.242900    1.263500    1.047600    1.281500   0.972820   -0.184630   -0.107140   6.110000   -1.087900   0.136040   -0.005866   4.880000   0.000000   -0.816130   1331.670000   -0.226010   -0.006600   -9.900000   -9.900000   100.530000   268.960000   0.140000   0.053000   0.687000   0.554000   0.357000   0.244000
    0.290    1.232400    1.251700    1.036300    1.273600   0.963480   -0.190570   -0.100110   6.120000   -1.091300   0.134990   -0.005666   4.900000   0.000000   -0.829500   1319.830000   -0.222500   -0.006650   -9.900000   -9.900000   102.010000   268.780000   0.139000   0.051000   0.681000   0.557000   0.360000   0.236000
    0.300    1.221700    1.240100    1.024600    1.265300   0.956760   -0.195900   -0.092855   6.140000   -1.094800   0.133880   -0.005475   4.930000   0.000000   -0.841650   1308.470000   -0.219120   -0.006700   -9.900000   -9.900000   103.150000   268.590000   0.138000   0.050000   0.675000   0.561000   0.363000   0.229000
    0.320    1.200700    1.217700    1.001100    1.247900   0.950040   -0.204540   -0.078923   6.160000   -1.101300   0.131790   -0.005122   4.980000   0.000000   -0.861750   1297.650000   -0.213180   -0.006800   -9.900000   -9.900000   104.000000   268.200000   0.135000   0.048000   0.670000   0.566000   0.366000   0.223000
    0.340    1.179000    1.195500    0.976770    1.228600   0.949560   -0.211340   -0.065134   6.180000   -1.107400   0.129840   -0.004808   5.030000   0.000000   -0.877260   1287.500000   -0.208160   -0.006890   -9.900000   -9.900000   104.700000   267.790000   0.133000   0.047000   0.664000   0.570000   0.369000   0.218000
    0.350    1.167400    1.183600    0.963800    1.217700   0.950770   -0.214460   -0.057921   6.180000   -1.110500   0.128900   -0.004663   5.060000   0.000000   -0.883750   1278.060000   -0.205920   -0.006930   -9.900000   -9.900000   105.260000   267.580000   0.130000   0.047000   0.658000   0.573000   0.372000   0.215000
    0.360    1.155800    1.172000    0.951200    1.206600   0.952780   -0.217160   -0.051040   6.190000   -1.113300   0.128060   -0.004527   5.080000   0.000000   -0.889650   1269.190000   -0.203790   -0.006970   -9.900000   -9.900000   105.610000   267.370000   0.128000   0.047000   0.653000   0.576000   0.375000   0.212000
    0.380    1.130500    1.146800    0.924400    1.181600   0.958990   -0.222140   -0.036755   6.190000   -1.119000   0.126470   -0.004276   5.120000   0.000000   -0.900380   1260.740000   -0.199780   -0.007050   -9.900000   -9.900000   105.870000   266.950000   0.125000   0.048000   0.648000   0.578000   0.378000   0.210000
    0.400    1.104600    1.121400    0.897650    1.155200   0.967660   -0.226080   -0.023189   6.200000   -1.124300   0.125120   -0.004053   5.160000   0.000000   -0.910920   1252.660000   -0.195820   -0.007130   -9.900000   -9.900000   106.020000   266.540000   0.122000   0.049000   0.643000   0.580000   0.381000   0.210000
    0.420    1.078200    1.095500    0.870670    1.127600   0.978620   -0.229240   -0.010417   6.200000   -1.129100   0.123890   -0.003853   5.200000   0.000000   -0.922410   1244.800000   -0.191710   -0.007190   -9.900000   -9.900000   106.030000   266.160000   0.120000   0.051000   0.638000   0.583000   0.384000   0.210000
    0.440    1.051500    1.069700    0.843550    1.099500   0.991440   -0.231660    0.001168   6.200000   -1.133700   0.122780   -0.003673   5.240000   0.000000   -0.934590   1237.030000   -0.187470   -0.007260   -9.900000   -9.900000   105.920000   265.800000   0.117000   0.053000   0.634000   0.585000   0.388000   0.211000
    0.450    1.037600    1.056200    0.829410    1.084700   0.998760   -0.232630    0.006589   6.200000   -1.135900   0.122270   -0.003590   5.250000   0.000000   -0.940750   1229.230000   -0.185340   -0.007290   -9.900000   -9.900000   105.790000   265.640000   0.115000   0.054000   0.629000   0.589000   0.393000   0.213000
    0.460    1.023400    1.042600    0.815090    1.069600   1.006400   -0.233500    0.011871   6.200000   -1.138100   0.121770   -0.003510   5.270000   0.000000   -0.946860   1221.160000   -0.183210   -0.007320   -9.900000   -9.900000   105.690000   265.480000   0.113000   0.055000   0.624000   0.592000   0.398000   0.216000
    0.480    0.997190    1.017200    0.788600    1.041500   1.021500   -0.234640    0.020767   6.200000   -1.142000   0.120930   -0.003360   5.300000   0.000000   -0.958630   1212.740000   -0.179020   -0.007380   -9.900000   -9.900000   105.590000   265.210000   0.111000   0.057000   0.619000   0.595000   0.404000   0.219000
    0.500    0.969910    0.991060    0.761500    1.012000   1.038400   -0.235220    0.029119   6.200000   -1.145900   0.120150   -0.003220   5.340000   0.000000   -0.969300   1203.910000   -0.175000   -0.007440   -9.900000   -9.900000   105.540000   265.000000   0.109000   0.060000   0.615000   0.599000   0.410000   0.224000
    0.550    0.904800    0.928300    0.698400    0.941700   1.083300   -0.234490    0.046932   6.200000   -1.154300   0.118470   -0.002897   5.410000   0.000000   -0.989250   1194.590000   -0.166010   -0.007580   -9.900000   -9.900000   105.610000   264.740000   0.108000   0.066000   0.610000   0.603000   0.417000   0.229000
    0.600    0.841650    0.867150    0.638750    0.873510   1.133600   -0.231280    0.062667   6.200000   -1.161500   0.116710   -0.002610   5.480000   0.000000   -1.001200   1184.930000   -0.158300   -0.007730   -9.900000   -9.900000   105.830000   264.830000   0.106000   0.071000   0.605000   0.607000   0.424000   0.235000
    0.650    0.781810    0.808760    0.582310    0.809480   1.186100   -0.226660    0.077997   6.200000   -1.167600   0.114650   -0.002356   5.530000   0.000000   -1.007800   1175.190000   -0.151440   -0.007870    0.005829    0.003762   106.200000   265.200000   0.105000   0.073000   0.599000   0.611000   0.431000   0.243000
    0.667    0.762620    0.789940    0.564220    0.789160   1.203500   -0.224970    0.083058   6.200000   -1.169400   0.113940   -0.002276   5.540000   0.000000   -1.009300   1165.690000   -0.149230   -0.007920    0.026446    0.017100   106.750000   265.380000   0.103000   0.074000   0.593000   0.615000   0.440000   0.250000
    0.700    0.725130    0.753020    0.528780    0.749850   1.237500   -0.221430    0.093185   6.200000   -1.172800   0.112530   -0.002131   5.560000   0.000000   -1.011700   1156.460000   -0.145030   -0.008000    0.055495    0.035728   107.480000   265.780000   0.102000   0.073000   0.587000   0.619000   0.448000   0.258000
    0.750    0.669030    0.697370    0.475230    0.691730   1.287100   -0.215910    0.108290   6.200000   -1.177700   0.110540   -0.001931   5.600000   0.000000   -1.015400   1147.590000   -0.138660   -0.008120    0.092259    0.059024   108.390000   266.510000   0.100000   0.070000   0.581000   0.622000   0.457000   0.266000
    0.800    0.613460    0.641960    0.421730    0.635190   1.334100   -0.210470    0.122560   6.200000   -1.181900   0.108730   -0.001754   5.630000   0.000000   -1.021000   1139.210000   -0.132010   -0.008220    0.140400    0.088481   109.620000   267.320000   0.099000   0.063000   0.576000   0.624000   0.466000   0.274000
    0.850    0.558530    0.586980    0.368130    0.579690   1.378000   -0.205280    0.136080   6.200000   -1.185400   0.107090   -0.001597   5.660000   0.000000   -1.028200   1131.340000   -0.125170   -0.008300    0.194950    0.120290   111.080000   268.140000   0.099000   0.053000   0.570000   0.625000   0.475000   0.281000
    0.900    0.502960    0.531360    0.313760    0.523610   1.420800   -0.200110    0.149830   6.200000   -1.188400   0.105480   -0.001456   5.690000   0.000000   -1.036000   1123.910000   -0.118310   -0.008360    0.252030    0.151590   112.710000   268.900000   0.098000   0.042000   0.564000   0.626000   0.483000   0.288000
    0.950    0.447010    0.475410    0.259190    0.467060   1.462300   -0.194900    0.164320   6.200000   -1.190900   0.103890   -0.001328   5.720000   0.000000   -1.043600   1116.830000   -0.111600   -0.008410    0.309400    0.180770   114.500000   269.550000   0.098000   0.030000   0.558000   0.626000   0.491000   0.294000
    1.000    0.393200    0.421800    0.207000    0.412400   1.500400   -0.189830    0.178950   6.200000   -1.193000   0.102480   -0.001210   5.740000   0.000000   -1.050000   1109.950000   -0.105210   -0.008440    0.366950    0.207890   116.390000   270.000000   0.098000   0.020000   0.553000   0.625000   0.498000   0.298000
    1.100    0.284840    0.313740    0.101820    0.302090   1.569000   -0.180010    0.210420   6.200000   -1.196600   0.100160   -0.000994   5.820000   0.000000   -1.057300   1103.070000   -0.093837   -0.008470    0.424520    0.233100   118.300000   270.180000   0.099000   0.007000   0.548000   0.624000   0.505000   0.302000
    1.200    0.173400    0.202590   -0.006195    0.188660   1.628200   -0.170900    0.244100   6.200000   -1.199600   0.098482   -0.000803   5.920000   0.000000   -1.058400   1096.040000   -0.084144   -0.008420    0.481160    0.256060   120.190000   269.420000   0.100000   0.002000   0.543000   0.623000   0.511000   0.306000
    1.300    0.061520    0.091060   -0.113450    0.074330   1.679400   -0.162330    0.277990   6.200000   -1.201800   0.097375   -0.000635   6.010000   0.000000   -1.055400   1088.670000   -0.075819   -0.008290    0.535970    0.276480   122.010000   267.820000   0.101000   0.003000   0.539000   0.622000   0.516000   0.309000
    1.400   -0.045750   -0.015700   -0.215500   -0.036070   1.723900   -0.154130    0.309560   6.200000   -1.203900   0.096743   -0.000490   6.100000   0.000000   -1.050400   1080.770000   -0.068543   -0.008060    0.588320    0.294240   123.750000   265.450000   0.102000   0.006000   0.535000   0.620000   0.521000   0.312000
    1.500   -0.149540   -0.118660   -0.313800   -0.143700   1.762200   -0.146700    0.338960   6.200000   -1.206300   0.096445   -0.000365   6.180000   0.000000   -1.045400   1072.390000   -0.062000   -0.007710    0.637890    0.309440   125.380000   262.410000   0.104000   0.010000   0.532000   0.619000   0.525000   0.315000
    1.600   -0.248600   -0.216720   -0.406820   -0.247080   1.795500   -0.139970    0.366160   6.200000   -1.208600   0.096338   -0.000259   6.260000   0.000000   -1.042100   1061.770000   -0.055927   -0.007230    0.688890    0.324490   126.900000   258.780000   0.105000   0.012000   0.529000   0.618000   0.528000   0.318000
    1.700   -0.341450   -0.308400   -0.492950   -0.344650   1.825900   -0.133610    0.390650   6.200000   -1.210600   0.096254   -0.000171   6.330000   0.000000   -1.040400   1049.290000   -0.050286   -0.006660    0.735780    0.337200   128.140000   254.660000   0.106000   0.012000   0.527000   0.618000   0.530000   0.321000
    1.800   -0.429750   -0.395580   -0.573880   -0.438180   1.856400   -0.126860    0.412440   6.200000   -1.212300   0.096207   -0.000099   6.400000   0.000000   -1.039700   1036.420000   -0.045096   -0.006030    0.779920    0.349630   129.110000   250.110000   0.106000   0.012000   0.526000   0.618000   0.531000   0.323000
    1.900   -0.512760   -0.477310   -0.648990   -0.526820   1.886800   -0.119590    0.431510   6.200000   -1.214100   0.096255   -0.000042   6.480000   0.000000   -1.039500   1023.140000   -0.040373   -0.005400    0.824270    0.364130   129.860000   245.250000   0.106000   0.010000   0.526000   0.618000   0.532000   0.326000
    2.000   -0.586690   -0.550030   -0.714660   -0.606580   1.915200   -0.112370    0.447880   6.200000   -1.215900   0.096361    0.000000   6.540000   0.000000   -1.039200   1009.490000   -0.036136   -0.004790    0.871380    0.382450   130.370000   240.140000   0.105000   0.008000   0.526000   0.618000   0.532000   0.329000
    2.200   -0.721430   -0.682200   -0.830030   -0.754020   1.968100   -0.098017    0.480240   6.200000   -1.219000   0.096497    0.000000   6.660000   0.000000   -1.036800    995.520000   -0.029105   -0.003780    0.920330    0.404020   130.670000   229.550000   0.103000   0.005000   0.527000   0.619000   0.533000   0.332000
    2.400   -0.848100   -0.806900   -0.932600   -0.894100   2.017000   -0.083765    0.518730   6.200000   -1.220200   0.096198    0.000000   6.730000   0.000000   -1.032300    981.330000   -0.023710   -0.003020    0.969310    0.427420   130.810000   219.050000   0.100000   0.003000   0.528000   0.619000   0.533000   0.335000
    2.500   -0.909660   -0.867650   -0.982280   -0.961870   2.040600   -0.076308    0.538830   6.200000   -1.220100   0.096106    0.000000   6.770000   0.000000   -1.029400    966.940000   -0.021509   -0.002720    1.016700    0.451420   130.810000   214.040000   0.097000   0.002000   0.530000   0.619000   0.534000   0.337000
    2.600   -0.968630   -0.925770   -1.031300   -1.026600   2.062800   -0.068925    0.558100   6.200000   -1.219800   0.096136    0.000000   6.810000   0.000000   -1.026200    952.340000   -0.019576   -0.002460    1.060100    0.474220   130.720000   209.320000   0.094000   0.001000   0.531000   0.620000   0.535000   0.340000
    2.800   -1.081700   -1.036700   -1.130100   -1.149500   2.101400   -0.055229    0.593940   6.200000   -1.218900   0.096667    0.000000   6.870000   0.000000   -1.019000    937.520000   -0.016324   -0.002080    1.099000    0.495360   130.570000   201.080000   0.091000   0.000000   0.532000   0.619000   0.536000   0.342000
    3.000   -1.189800   -1.142000   -1.230000   -1.266400   2.132300   -0.043320    0.626940   6.200000   -1.217900   0.097638    0.000000   6.930000   0.000000   -1.011200    922.430000   -0.013577   -0.001830    1.134800    0.515850   130.360000   195.000000   0.088000   0.000000   0.534000   0.619000   0.537000   0.344000
    3.200   -1.291400   -1.240600   -1.325500   -1.376000   2.154500   -0.034440    0.658110   6.200000   -1.216900   0.098649   -0.000023   6.990000   0.000000   -1.003200    908.790000   -0.011030   -0.001670    1.163800    0.534110   130.130000   191.610000   0.084000   0.000000   0.535000   0.618000   0.538000   0.345000
    3.400   -1.386000   -1.332200   -1.415000   -1.478600   2.170400   -0.027889    0.687550   6.200000   -1.216000   0.099553   -0.000040   7.080000   0.000000   -0.995120    896.150000   -0.008665   -0.001580    1.188200    0.551450   129.900000   190.730000   0.081000   0.000000   0.535000   0.618000   0.540000   0.346000
    3.500   -1.433200   -1.377800   -1.459900   -1.529700   2.177500   -0.024997    0.702160   6.200000   -1.215600   0.099989   -0.000045   7.120000   0.000000   -0.991000    883.160000   -0.007568   -0.001550    1.211400    0.569810   129.710000   191.110000   0.078000   0.000000   0.536000   0.617000   0.541000   0.347000
    3.600   -1.476200   -1.419300   -1.501400   -1.576400   2.183400   -0.022575    0.715230   6.200000   -1.215600   0.100430   -0.000049   7.160000   0.000000   -0.986820    870.050000   -0.006537   -0.001540    1.233600    0.589320   129.560000   191.980000   0.075000   0.000000   0.536000   0.616000   0.542000   0.348000
    3.800   -1.561700   -1.501400   -1.586500   -1.668500   2.193800   -0.018362    0.740280   6.200000   -1.215800   0.101420   -0.000053   7.240000   0.000000   -0.978260    857.070000   -0.004702   -0.001520    1.253400    0.609180   129.490000   195.010000   0.072000   0.000000   0.536000   0.616000   0.543000   0.349000
    4.000   -1.638800   -1.574800   -1.667300   -1.751600   2.204000   -0.014642    0.763030   6.200000   -1.216200   0.102180   -0.000052   7.320000   0.000000   -0.969380    844.480000   -0.003212   -0.001520    1.271100    0.629390   129.490000   199.450000   0.070000   0.000000   0.536000   0.616000   0.543000   0.349000
    4.200   -1.711600   -1.643900   -1.745100   -1.829000   2.212300   -0.012248    0.785520   6.200000   -1.216500   0.102690   -0.000047   7.390000   0.000000   -0.960120    832.450000   -0.002103   -0.001520    1.287000    0.651570   129.570000   204.930000   0.068000   0.000000   0.535000   0.616000   0.542000   0.349000
    4.400   -1.779800   -1.708900   -1.819200   -1.901100   2.218100   -0.011459    0.807920   6.200000   -1.216900   0.103040   -0.000039   7.460000   0.000000   -0.950490    821.180000   -0.001324   -0.001500    1.300400    0.674210   129.710000   211.090000   0.066000   0.000000   0.534000   0.617000   0.540000   0.347000
    4.600   -1.846900   -1.773100   -1.892300   -1.971200   2.223000   -0.011760    0.831260   6.200000   -1.217500   0.103240   -0.000027   7.520000   0.000000   -0.940500    810.790000   -0.000804   -0.001480    1.312500    0.696980   129.870000   217.560000   0.064000   0.000000   0.533000   0.619000   0.538000   0.345000
    4.800   -1.906300   -1.830300   -1.957300   -2.032600   2.226800   -0.012879    0.852400   6.200000   -1.218200   0.103370   -0.000014   7.640000   0.000000   -0.930180    801.410000   -0.000471   -0.001460    1.322500    0.718880   130.050000   223.990000   0.063000   0.000000   0.531000   0.621000   0.535000   0.341000
    5.000   -1.966000   -1.888200   -2.024500   -2.092800   2.229900   -0.014855    0.873140   6.200000   -1.218900   0.103530    0.000000   7.780000   0.000000   -0.919540    793.130000   -0.000255   -0.001440    1.328900    0.738060   130.220000   230.000000   0.061000   0.000000   0.528000   0.622000   0.532000   0.335000
    5.500   -2.105100   -2.023200   -2.190800   -2.228800   2.238900   -0.019502    0.914660   6.200000   -1.220400   0.104600    0.000000   8.070000   0.000000   -0.891760    785.730000    0.000072   -0.001400    1.345000    0.777900   130.390000   241.860000   0.060000   0.000000   0.526000   0.624000   0.528000   0.329000
    6.000   -2.242100   -2.156300   -2.365900   -2.357900   2.237700   -0.026383    0.948700   6.200000   -1.223200   0.107500    0.000000   8.480000   0.000000   -0.862860    779.910000    0.000188   -0.001380    1.349600    0.802780   130.530000   249.340000   0.059000   0.000000   0.524000   0.625000   0.524000   0.321000
    6.500   -2.368600   -2.278500   -2.532200   -2.477000   2.215000   -0.039505    0.976430   6.200000   -1.229900   0.112310    0.000000   8.900000   0.000000   -0.833550    775.600000    0.000159   -0.001370    1.348900    0.814800   130.630000   252.940000   0.059000   0.000000   0.520000   0.634000   0.517000   0.312000
    7.000   -2.482700   -2.388100   -2.681800   -2.585400   2.172000   -0.059140    0.997570   6.200000   -1.240800   0.118530    0.000000   9.200000   0.000000   -0.804570    772.680000    0.000056   -0.001370    1.342200    0.816150   130.700000   253.120000   0.059000   0.000000   0.515000   0.636000   0.514000   0.302000
    7.500   -2.586500   -2.487400   -2.817600   -2.685400   2.118700   -0.081606    1.012100   6.200000   -1.254300   0.125070    0.000000   9.480000   0.000000   -0.776650    771.010000   -0.000055   -0.001370    1.328800    0.809000   130.720000   250.390000   0.058000   0.000000   0.512000   0.634000   0.511000   0.270000
    8.000   -2.686100   -2.582900   -2.943800   -2.782300   2.061300   -0.103820    1.023200   6.200000   -1.268800   0.131460    0.000000   9.570000   0.000000   -0.750330    760.810000   -0.000117   -0.001370    1.308400    0.795380   130.870000   245.230000   0.059000   0.000000   0.510000   0.630000   0.507000   0.278000
    8.500   -2.782000   -2.675200   -3.059700   -2.877600   2.008400   -0.121140    1.033500   6.200000   -1.283900   0.137420    0.000000   9.620000   0.000000   -0.725440    764.500000   -0.000131   -0.001370    1.282200    0.776710   130.710000   238.130000   0.059000   0.000000   0.509000   0.622000   0.503000   0.265000
    9.000   -2.879200   -2.768700   -3.171300   -2.975900   1.960500   -0.134070    1.045300   6.200000   -1.298900   0.142940    0.000000   9.660000   0.000000   -0.701610    768.070000   -0.000108   -0.001370    1.251500    0.754250   130.500000   229.560000   0.060000   0.000000   0.509000   0.613000   0.498000   0.252000
    9.500   -2.976900   -2.863400   -3.278500   -3.076000   1.918900   -0.143640    1.056700   6.200000   -1.313000   0.147810    0.000000   9.660000   0.000000   -0.678500    771.550000   -0.000060   -0.001360    1.217900    0.729250   130.260000   220.020000   0.060000   0.000000   0.509000   0.604000   0.492000   0.239000
    10.00   -3.070200   -2.953700   -3.377600   -3.172600   1.883700   -0.150960    1.065100   6.200000   -1.325300   0.151830    0.000000   9.660000   0.000000   -0.655750    775.000000    0.000000   -0.001360    1.182900    0.703000   130.000000   210.000000   0.060000   0.000000   0.510000   0.604000   0.487000   0.239000
    """) 
[docs]class BooreEtAl2014HighQ(BooreEtAl2014):
    """
    This class implements the Boore et al. (2014) model considering the
    correction to the path scaling term for High Q regions (e.g. China and
    Turkey)
    The modification is made to the "Dc3" coefficient
    """
    COEFFS = CoeffsTable(sa_damping=5, table="""\
    IMT            e0          e1          e2          e3         e4          e5          e6         Mh          c1         c2          c3          h        Dc3           c            Vc          f4          f5          f6          f7           R1           R2        DfR        DfV         f1         f2         tau1         tau2
    pgv      5.037000    5.078000    4.849000    5.033000   1.073000   -0.153600    0.225200   6.200000   -1.243000   0.148900   -0.003440   5.300000   0.004345   -0.840000   1300.000000   -0.100000   -0.008440   -9.900000   -9.900000   105.000000   272.000000   0.082000   0.080000   0.644000   0.552000   0.401000   0.346000
    pga      0.447300    0.485600    0.245900    0.453900   1.431000    0.050530   -0.166200   5.500000   -1.134000   0.191700   -0.008088   4.500000   0.002858   -0.600000   1500.000000   -0.150000   -0.007010   -9.900000   -9.900000   110.000000   270.000000   0.100000   0.070000   0.695000   0.495000   0.398000   0.348000
    0.010    0.453400    0.491600    0.251900    0.459900   1.421000    0.049320   -0.165900   5.500000   -1.134000   0.191600   -0.008088   4.500000   0.002816   -0.603720   1500.200000   -0.148330   -0.007010   -9.900000   -9.900000   111.670000   270.000000   0.096000   0.070000   0.698000   0.499000   0.402000   0.345000
    0.020    0.485980    0.523590    0.297070    0.488750   1.433100    0.053388   -0.165610   5.500000   -1.139400   0.189620   -0.008074   4.500000   0.002780   -0.573880   1500.360000   -0.147100   -0.007280   -9.900000   -9.900000   113.100000   270.000000   0.092000   0.030000   0.702000   0.502000   0.409000   0.346000
    0.022    0.498660    0.536470    0.313470    0.499730   1.433600    0.054888   -0.165200   5.500000   -1.140500   0.189240   -0.008095   4.500000   0.002760   -0.566750   1500.680000   -0.148010   -0.007320   -9.900000   -9.900000   113.370000   270.000000   0.088000   0.027000   0.707000   0.505000   0.418000   0.349000
    0.025    0.522830    0.561300    0.344260    0.519990   1.432800    0.057529   -0.164990   5.500000   -1.141900   0.188750   -0.008153   4.500000   0.002754   -0.555200   1501.040000   -0.150150   -0.007360   -9.900000   -9.900000   113.070000   270.000000   0.086000   0.026000   0.711000   0.508000   0.427000   0.354000
    0.029    0.559490    0.599230    0.391460    0.549950   1.427900    0.060732   -0.166320   5.500000   -1.142300   0.188440   -0.008290   4.500000   0.002762   -0.538500   1501.260000   -0.153870   -0.007370   -9.900000   -9.900000   112.360000   270.000000   0.084000   0.028000   0.716000   0.510000   0.436000   0.359000
    0.030    0.569160    0.609200    0.403910    0.557830   1.426100    0.061444   -0.166900   5.500000   -1.142100   0.188420   -0.008336   4.490000   0.002765   -0.534140   1502.950000   -0.154850   -0.007350   -9.900000   -9.900000   112.130000   270.000000   0.081000   0.029000   0.721000   0.514000   0.445000   0.364000
    0.032    0.588020    0.628750    0.427880    0.573300   1.422700    0.062806   -0.168130   5.500000   -1.141200   0.188400   -0.008445   4.450000   0.002769   -0.525290   1503.120000   -0.156850   -0.007310   -9.900000   -9.900000   111.650000   270.000000   0.078000   0.030000   0.726000   0.516000   0.454000   0.369000
    0.035    0.616360    0.658180    0.462520    0.597040   1.417400    0.064559   -0.170150   5.500000   -1.138800   0.188390   -0.008642   4.400000   0.002784   -0.511920   1503.240000   -0.160160   -0.007210   -9.900000   -9.900000   110.640000   270.000000   0.077000   0.031000   0.730000   0.518000   0.462000   0.374000
    0.036    0.625540    0.667720    0.473380    0.604960   1.415800    0.065028   -0.170830   5.500000   -1.137800   0.188370   -0.008715   4.380000   0.002801   -0.507520   1503.320000   -0.161420   -0.007170   -9.900000   -9.900000   109.530000   270.000000   0.075000   0.031000   0.734000   0.520000   0.470000   0.379000
    0.040    0.662810    0.706040    0.515320    0.638280   1.409000    0.066183   -0.173570   5.500000   -1.132400   0.188160   -0.009030   4.320000   0.002822   -0.490650   1503.350000   -0.167770   -0.006980   -9.900000   -9.900000   108.280000   270.000000   0.073000   0.032000   0.738000   0.521000   0.478000   0.384000
    0.042    0.680870    0.724430    0.534450    0.655050   1.405900    0.066438   -0.174850   5.500000   -1.129200   0.187970   -0.009195   4.290000   0.002847   -0.482900   1503.340000   -0.171930   -0.006870   -9.900000   -9.900000   106.990000   270.000000   0.072000   0.032000   0.742000   0.523000   0.484000   0.390000
    0.044    0.698820    0.742770    0.552820    0.672250   1.403300    0.066663   -0.176190   5.500000   -1.125900   0.187750   -0.009360   4.270000   0.002874   -0.475720   1503.130000   -0.176640   -0.006770   -9.900000   -9.900000   105.410000   270.000000   0.070000   0.031000   0.745000   0.525000   0.490000   0.397000
    0.045    0.708220    0.752320    0.562220    0.681390   1.402100    0.066774   -0.176930   5.500000   -1.124200   0.187640   -0.009441   4.250000   0.002900   -0.472360   1502.840000   -0.179140   -0.006720   -9.900000   -9.900000   103.610000   270.000000   0.069000   0.031000   0.748000   0.527000   0.496000   0.405000
    0.046    0.717790    0.762020    0.571660    0.690760   1.400900    0.066891   -0.177690   5.500000   -1.122400   0.187520   -0.009521   4.240000   0.002923   -0.469150   1502.470000   -0.181700   -0.006670   -9.900000   -9.900000   101.700000   270.000000   0.067000   0.031000   0.750000   0.529000   0.499000   0.412000
    0.048    0.735740    0.780150    0.588880    0.708540   1.399100    0.067127   -0.179200   5.500000   -1.119200   0.187300   -0.009676   4.220000   0.002944   -0.463210   1502.010000   -0.186880   -0.006560   -9.900000   -9.900000    99.760000   270.000000   0.065000   0.031000   0.752000   0.530000   0.502000   0.419000
    0.050    0.754360    0.799050    0.606520    0.727260   1.397400    0.067357   -0.180820   5.500000   -1.115900   0.187090   -0.009819   4.200000   0.002957   -0.457950   1501.420000   -0.192000   -0.006470   -9.900000   -9.900000    97.930000   270.000000   0.063000   0.030000   0.753000   0.532000   0.503000   0.426000
    0.055    0.799600    0.844500    0.647700    0.773700   1.394700    0.067797   -0.184800   5.500000   -1.108200   0.186550   -0.010120   4.150000   0.002964   -0.447870   1500.710000   -0.203690   -0.006250   -9.900000   -9.900000    96.030000   270.000000   0.062000   0.029000   0.753000   0.534000   0.502000   0.434000
    0.060    0.843940    0.888840    0.685620    0.820670   1.395400    0.068591   -0.188580   5.500000   -1.100900   0.185820   -0.010330   4.110000   0.002968   -0.441860   1499.830000   -0.213740   -0.006070   -9.900000   -9.900000    94.100000   270.010000   0.061000   0.027000   0.753000   0.536000   0.499000   0.441000
    0.065    0.886550    0.931160    0.719410    0.867240   1.400400    0.070127   -0.191760   5.500000   -1.094200   0.184850   -0.010480   4.080000   0.002968   -0.439510   1498.740000   -0.222250   -0.005930   -9.900000   -9.900000    92.080000   270.020000   0.061000   0.025000   0.752000   0.538000   0.495000   0.448000
    0.067    0.902700    0.947110    0.731710    0.885260   1.403200    0.070895   -0.192910   5.500000   -1.091800   0.184420   -0.010520   4.070000   0.002968   -0.439500   1497.420000   -0.225240   -0.005880   -9.900000   -9.900000    90.010000   270.020000   0.061000   0.025000   0.750000   0.540000   0.489000   0.455000
    0.070    0.926520    0.970570    0.749400    0.912270   1.408200    0.072075   -0.194510   5.500000   -1.088400   0.183690   -0.010560   4.060000   0.002964   -0.440400   1495.850000   -0.229310   -0.005820   -9.900000   -9.900000    87.970000   270.030000   0.062000   0.024000   0.748000   0.541000   0.483000   0.461000
    0.075    0.964470    1.007700    0.776780    0.956300   1.417400    0.073549   -0.196650   5.500000   -1.083100   0.182250   -0.010580   4.040000   0.002957   -0.444110   1494.000000   -0.235000   -0.005730   -9.900000   -9.900000    85.990000   270.040000   0.064000   0.022000   0.745000   0.542000   0.474000   0.466000
    0.080    1.000300    1.042600    0.801610    0.998180   1.426100    0.073735   -0.198160   5.500000   -1.078500   0.180520   -0.010560   4.020000   0.002945   -0.450200   1491.820000   -0.239440   -0.005670   -9.900000   -9.900000    84.230000   270.050000   0.067000   0.020000   0.741000   0.543000   0.464000   0.468000
    0.085    1.034000    1.075500    0.824230    1.037900   1.432200    0.071940   -0.199020   5.510000   -1.074500   0.178560   -0.010510   4.030000   0.002930   -0.458130   1489.290000   -0.242850   -0.005630   -9.900000   -9.900000    82.740000   270.060000   0.072000   0.019000   0.737000   0.543000   0.452000   0.468000
    0.090    1.066600    1.107600    0.845910    1.076200   1.435000    0.068097   -0.199290   5.520000   -1.070900   0.176430   -0.010420   4.070000   0.002914   -0.467320   1486.360000   -0.245440   -0.005610   -9.900000   -9.900000    81.540000   270.070000   0.076000   0.017000   0.734000   0.542000   0.440000   0.466000
    0.095    1.098100    1.138500    0.867030    1.112700   1.433900    0.062327   -0.199000   5.530000   -1.067800   0.174200   -0.010320   4.100000   0.002897   -0.477210   1482.980000   -0.247470   -0.005600   -9.900000   -9.900000    80.460000   270.080000   0.082000   0.016000   0.731000   0.542000   0.428000   0.464000
    0.100    1.126800    1.166900    0.887100    1.145400   1.429300    0.055231   -0.198380   5.540000   -1.065200   0.172030   -0.010200   4.130000   0.002879   -0.487240   1479.120000   -0.249160   -0.005600   -9.900000   -9.900000    79.590000   270.090000   0.087000   0.014000   0.728000   0.541000   0.415000   0.458000
    0.110    1.178500    1.217900    0.927020    1.203000   1.411000    0.037389   -0.196010   5.570000   -1.060700   0.167700   -0.009964   4.190000   0.002863   -0.506320   1474.740000   -0.252130   -0.005620   -9.900000   -9.900000    79.050000   270.110000   0.093000   0.012000   0.726000   0.540000   0.403000   0.451000
    0.120    1.223000    1.262100    0.966160    1.250200   1.383100    0.016373   -0.192650   5.620000   -1.057200   0.163520   -0.009722   4.240000   0.002847   -0.524380   1469.750000   -0.254550   -0.005670   -9.900000   -9.900000    78.850000   270.130000   0.099000   0.011000   0.724000   0.539000   0.392000   0.441000
    0.130    1.259600    1.298600    1.003100    1.286900   1.349700   -0.005158   -0.188980   5.660000   -1.054900   0.159820   -0.009476   4.290000   0.002831   -0.542140   1464.090000   -0.256280   -0.005720   -9.900000   -9.900000    78.990000   270.150000   0.104000   0.011000   0.723000   0.538000   0.381000   0.430000
    0.133    1.269200    1.308200    1.013500    1.296100   1.339500   -0.011354   -0.187920   5.670000   -1.054500   0.158820   -0.009402   4.300000   0.002817   -0.547520   1457.760000   -0.256650   -0.005740   -9.900000   -9.900000    79.470000   270.150000   0.110000   0.011000   0.722000   0.538000   0.371000   0.417000
    0.140    1.288300    1.327000    1.036000    1.313700   1.316200   -0.024711   -0.185660   5.700000   -1.053700   0.156720   -0.009228   4.340000   0.002802   -0.560320   1450.710000   -0.257190   -0.005780   -9.900000   -9.900000    80.260000   270.160000   0.115000   0.012000   0.721000   0.537000   0.362000   0.403000
    0.150    1.309500    1.348100    1.064800    1.332400   1.284400   -0.042065   -0.182340   5.740000   -1.053200   0.154010   -0.008977   4.390000   0.002786   -0.579620   1442.850000   -0.257130   -0.005850   -9.900000   -9.900000    81.330000   270.160000   0.120000   0.015000   0.720000   0.537000   0.354000   0.388000
    0.160    1.323500    1.361500    1.087600    1.343700   1.254100   -0.057593   -0.178530   5.780000   -1.053300   0.151580   -0.008725   4.440000   0.002765   -0.600520   1434.220000   -0.256040   -0.005910   -9.900000   -9.900000    82.860000   270.160000   0.125000   0.020000   0.720000   0.536000   0.349000   0.372000
    0.170    1.330600    1.367900    1.104000    1.348700   1.224400   -0.071861   -0.174210   5.820000   -1.054100   0.149480   -0.008472   4.490000   0.002734   -0.622520   1424.850000   -0.254140   -0.005970   -9.900000   -9.900000    84.720000   270.140000   0.128000   0.026000   0.718000   0.536000   0.346000   0.357000
    0.180    1.332700    1.368900    1.114900    1.349200   1.194100   -0.085640   -0.169390   5.850000   -1.055600   0.147680   -0.008219   4.530000   0.002699   -0.644860   1414.770000   -0.251730   -0.006020   -9.900000   -9.900000    86.670000   270.110000   0.131000   0.033000   0.717000   0.536000   0.344000   0.341000
    0.190    1.330700    1.365600    1.120800    1.346300   1.163500   -0.098884   -0.164040   5.890000   -1.057900   0.146160   -0.007967   4.570000   0.002660   -0.666810   1403.990000   -0.249110   -0.006080   -9.900000   -9.900000    88.730000   270.060000   0.134000   0.039000   0.714000   0.537000   0.343000   0.324000
    0.200    1.325500    1.359000    1.122000    1.341400   1.134900   -0.110960   -0.158520   5.920000   -1.060700   0.144890   -0.007717   4.610000   0.002612   -0.687620   1392.610000   -0.246580   -0.006140   -9.900000   -9.900000    90.910000   270.000000   0.136000   0.045000   0.711000   0.539000   0.344000   0.309000
    0.220    1.309100    1.339400    1.113300    1.328100   1.082300   -0.133000   -0.147040   5.970000   -1.067000   0.142630   -0.007224   4.680000   0.002560   -0.724310   1380.720000   -0.242350   -0.006260   -9.900000   -9.900000    93.040000   269.830000   0.138000   0.052000   0.708000   0.541000   0.345000   0.294000
    0.240    1.288100    1.315000    1.094500    1.313200   1.036600   -0.152990   -0.134450   6.030000   -1.073700   0.140350   -0.006747   4.750000   0.002505   -0.756460   1368.510000   -0.238230   -0.006380   -9.900000   -9.900000    95.080000   269.590000   0.140000   0.055000   0.703000   0.544000   0.347000   0.280000
    0.250    1.276600    1.301700    1.082800    1.305200   1.016600   -0.162130   -0.127840   6.050000   -1.077300   0.139250   -0.006517   4.780000   0.002444   -0.771770   1356.210000   -0.235740   -0.006440   -9.900000   -9.900000    97.040000   269.450000   0.141000   0.055000   0.698000   0.547000   0.350000   0.266000
    0.260    1.265100    1.288600    1.071000    1.297200   0.999320   -0.170410   -0.121150   6.070000   -1.080800   0.138180   -0.006293   4.820000   0.002380   -0.786970   1343.890000   -0.232800   -0.006500   -9.900000   -9.900000    98.870000   269.300000   0.141000   0.055000   0.693000   0.550000   0.353000   0.255000
    0.280    1.242900    1.263500    1.047600    1.281500   0.972820   -0.184630   -0.107140   6.110000   -1.087900   0.136040   -0.005866   4.880000   0.002314   -0.816130   1331.670000   -0.226010   -0.006600   -9.900000   -9.900000   100.530000   268.960000   0.140000   0.053000   0.687000   0.554000   0.357000   0.244000
    0.290    1.232400    1.251700    1.036300    1.273600   0.963480   -0.190570   -0.100110   6.120000   -1.091300   0.134990   -0.005666   4.900000   0.002250   -0.829500   1319.830000   -0.222500   -0.006650   -9.900000   -9.900000   102.010000   268.780000   0.139000   0.051000   0.681000   0.557000   0.360000   0.236000
    0.300    1.221700    1.240100    1.024600    1.265300   0.956760   -0.195900   -0.092855   6.140000   -1.094800   0.133880   -0.005475   4.930000   0.002196   -0.841650   1308.470000   -0.219120   -0.006700   -9.900000   -9.900000   103.150000   268.590000   0.138000   0.050000   0.675000   0.561000   0.363000   0.229000
    0.320    1.200700    1.217700    1.001100    1.247900   0.950040   -0.204540   -0.078923   6.160000   -1.101300   0.131790   -0.005122   4.980000   0.002155   -0.861750   1297.650000   -0.213180   -0.006800   -9.900000   -9.900000   104.000000   268.200000   0.135000   0.048000   0.670000   0.566000   0.366000   0.223000
    0.340    1.179000    1.195500    0.976770    1.228600   0.949560   -0.211340   -0.065134   6.180000   -1.107400   0.129840   -0.004808   5.030000   0.002124   -0.877260   1287.500000   -0.208160   -0.006890   -9.900000   -9.900000   104.700000   267.790000   0.133000   0.047000   0.664000   0.570000   0.369000   0.218000
    0.350    1.167400    1.183600    0.963800    1.217700   0.950770   -0.214460   -0.057921   6.180000   -1.110500   0.128900   -0.004663   5.060000   0.002104   -0.883750   1278.060000   -0.205920   -0.006930   -9.900000   -9.900000   105.260000   267.580000   0.130000   0.047000   0.658000   0.573000   0.372000   0.215000
    0.360    1.155800    1.172000    0.951200    1.206600   0.952780   -0.217160   -0.051040   6.190000   -1.113300   0.128060   -0.004527   5.080000   0.002097   -0.889650   1269.190000   -0.203790   -0.006970   -9.900000   -9.900000   105.610000   267.370000   0.128000   0.047000   0.653000   0.576000   0.375000   0.212000
    0.380    1.130500    1.146800    0.924400    1.181600   0.958990   -0.222140   -0.036755   6.190000   -1.119000   0.126470   -0.004276   5.120000   0.002098   -0.900380   1260.740000   -0.199780   -0.007050   -9.900000   -9.900000   105.870000   266.950000   0.125000   0.048000   0.648000   0.578000   0.378000   0.210000
    0.400    1.104600    1.121400    0.897650    1.155200   0.967660   -0.226080   -0.023189   6.200000   -1.124300   0.125120   -0.004053   5.160000   0.002107   -0.910920   1252.660000   -0.195820   -0.007130   -9.900000   -9.900000   106.020000   266.540000   0.122000   0.049000   0.643000   0.580000   0.381000   0.210000
    0.420    1.078200    1.095500    0.870670    1.127600   0.978620   -0.229240   -0.010417   6.200000   -1.129100   0.123890   -0.003853   5.200000   0.002127   -0.922410   1244.800000   -0.191710   -0.007190   -9.900000   -9.900000   106.030000   266.160000   0.120000   0.051000   0.638000   0.583000   0.384000   0.210000
    0.440    1.051500    1.069700    0.843550    1.099500   0.991440   -0.231660    0.001168   6.200000   -1.133700   0.122780   -0.003673   5.240000   0.002160   -0.934590   1237.030000   -0.187470   -0.007260   -9.900000   -9.900000   105.920000   265.800000   0.117000   0.053000   0.634000   0.585000   0.388000   0.211000
    0.450    1.037600    1.056200    0.829410    1.084700   0.998760   -0.232630    0.006589   6.200000   -1.135900   0.122270   -0.003590   5.250000   0.002202   -0.940750   1229.230000   -0.185340   -0.007290   -9.900000   -9.900000   105.790000   265.640000   0.115000   0.054000   0.629000   0.589000   0.393000   0.213000
    0.460    1.023400    1.042600    0.815090    1.069600   1.006400   -0.233500    0.011871   6.200000   -1.138100   0.121770   -0.003510   5.270000   0.002250   -0.946860   1221.160000   -0.183210   -0.007320   -9.900000   -9.900000   105.690000   265.480000   0.113000   0.055000   0.624000   0.592000   0.398000   0.216000
    0.480    0.997190    1.017200    0.788600    1.041500   1.021500   -0.234640    0.020767   6.200000   -1.142000   0.120930   -0.003360   5.300000   0.002298   -0.958630   1212.740000   -0.179020   -0.007380   -9.900000   -9.900000   105.590000   265.210000   0.111000   0.057000   0.619000   0.595000   0.404000   0.219000
    0.500    0.969910    0.991060    0.761500    1.012000   1.038400   -0.235220    0.029119   6.200000   -1.145900   0.120150   -0.003220   5.340000   0.002348   -0.969300   1203.910000   -0.175000   -0.007440   -9.900000   -9.900000   105.540000   265.000000   0.109000   0.060000   0.615000   0.599000   0.410000   0.224000
    0.550    0.904800    0.928300    0.698400    0.941700   1.083300   -0.234490    0.046932   6.200000   -1.154300   0.118470   -0.002897   5.410000   0.002399   -0.989250   1194.590000   -0.166010   -0.007580   -9.900000   -9.900000   105.610000   264.740000   0.108000   0.066000   0.610000   0.603000   0.417000   0.229000
    0.600    0.841650    0.867150    0.638750    0.873510   1.133600   -0.231280    0.062667   6.200000   -1.161500   0.116710   -0.002610   5.480000   0.002452   -1.001200   1184.930000   -0.158300   -0.007730   -9.900000   -9.900000   105.830000   264.830000   0.106000   0.071000   0.605000   0.607000   0.424000   0.235000
    0.650    0.781810    0.808760    0.582310    0.809480   1.186100   -0.226660    0.077997   6.200000   -1.167600   0.114650   -0.002356   5.530000   0.002508   -1.007800   1175.190000   -0.151440   -0.007870    0.005829    0.003762   106.200000   265.200000   0.105000   0.073000   0.599000   0.611000   0.431000   0.243000
    0.667    0.762620    0.789940    0.564220    0.789160   1.203500   -0.224970    0.083058   6.200000   -1.169400   0.113940   -0.002276   5.540000   0.002568   -1.009300   1165.690000   -0.149230   -0.007920    0.026446    0.017100   106.750000   265.380000   0.103000   0.074000   0.593000   0.615000   0.440000   0.250000
    0.700    0.725130    0.753020    0.528780    0.749850   1.237500   -0.221430    0.093185   6.200000   -1.172800   0.112530   -0.002131   5.560000   0.002629   -1.011700   1156.460000   -0.145030   -0.008000    0.055495    0.035728   107.480000   265.780000   0.102000   0.073000   0.587000   0.619000   0.448000   0.258000
    0.750    0.669030    0.697370    0.475230    0.691730   1.287100   -0.215910    0.108290   6.200000   -1.177700   0.110540   -0.001931   5.600000   0.002690   -1.015400   1147.590000   -0.138660   -0.008120    0.092259    0.059024   108.390000   266.510000   0.100000   0.070000   0.581000   0.622000   0.457000   0.266000
    0.800    0.613460    0.641960    0.421730    0.635190   1.334100   -0.210470    0.122560   6.200000   -1.181900   0.108730   -0.001754   5.630000   0.002746   -1.021000   1139.210000   -0.132010   -0.008220    0.140400    0.088481   109.620000   267.320000   0.099000   0.063000   0.576000   0.624000   0.466000   0.274000
    0.850    0.558530    0.586980    0.368130    0.579690   1.378000   -0.205280    0.136080   6.200000   -1.185400   0.107090   -0.001597   5.660000   0.002797   -1.028200   1131.340000   -0.125170   -0.008300    0.194950    0.120290   111.080000   268.140000   0.099000   0.053000   0.570000   0.625000   0.475000   0.281000
    0.900    0.502960    0.531360    0.313760    0.523610   1.420800   -0.200110    0.149830   6.200000   -1.188400   0.105480   -0.001456   5.690000   0.002843   -1.036000   1123.910000   -0.118310   -0.008360    0.252030    0.151590   112.710000   268.900000   0.098000   0.042000   0.564000   0.626000   0.483000   0.288000
    0.950    0.447010    0.475410    0.259190    0.467060   1.462300   -0.194900    0.164320   6.200000   -1.190900   0.103890   -0.001328   5.720000   0.002885   -1.043600   1116.830000   -0.111600   -0.008410    0.309400    0.180770   114.500000   269.550000   0.098000   0.030000   0.558000   0.626000   0.491000   0.294000
    1.000    0.393200    0.421800    0.207000    0.412400   1.500400   -0.189830    0.178950   6.200000   -1.193000   0.102480   -0.001210   5.740000   0.002921   -1.050000   1109.950000   -0.105210   -0.008440    0.366950    0.207890   116.390000   270.000000   0.098000   0.020000   0.553000   0.625000   0.498000   0.298000
    1.100    0.284840    0.313740    0.101820    0.302090   1.569000   -0.180010    0.210420   6.200000   -1.196600   0.100160   -0.000994   5.820000   0.002954   -1.057300   1103.070000   -0.093837   -0.008470    0.424520    0.233100   118.300000   270.180000   0.099000   0.007000   0.548000   0.624000   0.505000   0.302000
    1.200    0.173400    0.202590   -0.006195    0.188660   1.628200   -0.170900    0.244100   6.200000   -1.199600   0.098482   -0.000803   5.920000   0.002982   -1.058400   1096.040000   -0.084144   -0.008420    0.481160    0.256060   120.190000   269.420000   0.100000   0.002000   0.543000   0.623000   0.511000   0.306000
    1.300    0.061520    0.091060   -0.113450    0.074330   1.679400   -0.162330    0.277990   6.200000   -1.201800   0.097375   -0.000635   6.010000   0.003007   -1.055400   1088.670000   -0.075819   -0.008290    0.535970    0.276480   122.010000   267.820000   0.101000   0.003000   0.539000   0.622000   0.516000   0.309000
    1.400   -0.045750   -0.015700   -0.215500   -0.036070   1.723900   -0.154130    0.309560   6.200000   -1.203900   0.096743   -0.000490   6.100000   0.003026   -1.050400   1080.770000   -0.068543   -0.008060    0.588320    0.294240   123.750000   265.450000   0.102000   0.006000   0.535000   0.620000   0.521000   0.312000
    1.500   -0.149540   -0.118660   -0.313800   -0.143700   1.762200   -0.146700    0.338960   6.200000   -1.206300   0.096445   -0.000365   6.180000   0.003039   -1.045400   1072.390000   -0.062000   -0.007710    0.637890    0.309440   125.380000   262.410000   0.104000   0.010000   0.532000   0.619000   0.525000   0.315000
    1.600   -0.248600   -0.216720   -0.406820   -0.247080   1.795500   -0.139970    0.366160   6.200000   -1.208600   0.096338   -0.000259   6.260000   0.003041   -1.042100   1061.770000   -0.055927   -0.007230    0.688890    0.324490   126.900000   258.780000   0.105000   0.012000   0.529000   0.618000   0.528000   0.318000
    1.700   -0.341450   -0.308400   -0.492950   -0.344650   1.825900   -0.133610    0.390650   6.200000   -1.210600   0.096254   -0.000171   6.330000   0.003026   -1.040400   1049.290000   -0.050286   -0.006660    0.735780    0.337200   128.140000   254.660000   0.106000   0.012000   0.527000   0.618000   0.530000   0.321000
    1.800   -0.429750   -0.395580   -0.573880   -0.438180   1.856400   -0.126860    0.412440   6.200000   -1.212300   0.096207   -0.000099   6.400000   0.003000   -1.039700   1036.420000   -0.045096   -0.006030    0.779920    0.349630   129.110000   250.110000   0.106000   0.012000   0.526000   0.618000   0.531000   0.323000
    1.900   -0.512760   -0.477310   -0.648990   -0.526820   1.886800   -0.119590    0.431510   6.200000   -1.214100   0.096255   -0.000042   6.480000   0.002968   -1.039500   1023.140000   -0.040373   -0.005400    0.824270    0.364130   129.860000   245.250000   0.106000   0.010000   0.526000   0.618000   0.532000   0.326000
    2.000   -0.586690   -0.550030   -0.714660   -0.606580   1.915200   -0.112370    0.447880   6.200000   -1.215900   0.096361    0.000000   6.540000   0.002923   -1.039200   1009.490000   -0.036136   -0.004790    0.871380    0.382450   130.370000   240.140000   0.105000   0.008000   0.526000   0.618000   0.532000   0.329000
    2.200   -0.721430   -0.682200   -0.830030   -0.754020   1.968100   -0.098017    0.480240   6.200000   -1.219000   0.096497    0.000000   6.660000   0.002870   -1.036800    995.520000   -0.029105   -0.003780    0.920330    0.404020   130.670000   229.550000   0.103000   0.005000   0.527000   0.619000   0.533000   0.332000
    2.400   -0.848100   -0.806900   -0.932600   -0.894100   2.017000   -0.083765    0.518730   6.200000   -1.220200   0.096198    0.000000   6.730000   0.002815   -1.032300    981.330000   -0.023710   -0.003020    0.969310    0.427420   130.810000   219.050000   0.100000   0.003000   0.528000   0.619000   0.533000   0.335000
    2.500   -0.909660   -0.867650   -0.982280   -0.961870   2.040600   -0.076308    0.538830   6.200000   -1.220100   0.096106    0.000000   6.770000   0.002761   -1.029400    966.940000   -0.021509   -0.002720    1.016700    0.451420   130.810000   214.040000   0.097000   0.002000   0.530000   0.619000   0.534000   0.337000
    2.600   -0.968630   -0.925770   -1.031300   -1.026600   2.062800   -0.068925    0.558100   6.200000   -1.219800   0.096136    0.000000   6.810000   0.002708   -1.026200    952.340000   -0.019576   -0.002460    1.060100    0.474220   130.720000   209.320000   0.094000   0.001000   0.531000   0.620000   0.535000   0.340000
    2.800   -1.081700   -1.036700   -1.130100   -1.149500   2.101400   -0.055229    0.593940   6.200000   -1.218900   0.096667    0.000000   6.870000   0.002658   -1.019000    937.520000   -0.016324   -0.002080    1.099000    0.495360   130.570000   201.080000   0.091000   0.000000   0.532000   0.619000   0.536000   0.342000
    3.000   -1.189800   -1.142000   -1.230000   -1.266400   2.132300   -0.043320    0.626940   6.200000   -1.217900   0.097638    0.000000   6.930000   0.002616   -1.011200    922.430000   -0.013577   -0.001830    1.134800    0.515850   130.360000   195.000000   0.088000   0.000000   0.534000   0.619000   0.537000   0.344000
    3.200   -1.291400   -1.240600   -1.325500   -1.376000   2.154500   -0.034440    0.658110   6.200000   -1.216900   0.098649   -0.000023   6.990000   0.002588   -1.003200    908.790000   -0.011030   -0.001670    1.163800    0.534110   130.130000   191.610000   0.084000   0.000000   0.535000   0.618000   0.538000   0.345000
    3.400   -1.386000   -1.332200   -1.415000   -1.478600   2.170400   -0.027889    0.687550   6.200000   -1.216000   0.099553   -0.000040   7.080000   0.002575   -0.995120    896.150000   -0.008665   -0.001580    1.188200    0.551450   129.900000   190.730000   0.081000   0.000000   0.535000   0.618000   0.540000   0.346000
    3.500   -1.433200   -1.377800   -1.459900   -1.529700   2.177500   -0.024997    0.702160   6.200000   -1.215600   0.099989   -0.000045   7.120000   0.002572   -0.991000    883.160000   -0.007568   -0.001550    1.211400    0.569810   129.710000   191.110000   0.078000   0.000000   0.536000   0.617000   0.541000   0.347000
    3.600   -1.476200   -1.419300   -1.501400   -1.576400   2.183400   -0.022575    0.715230   6.200000   -1.215600   0.100430   -0.000049   7.160000   0.002574   -0.986820    870.050000   -0.006537   -0.001540    1.233600    0.589320   129.560000   191.980000   0.075000   0.000000   0.536000   0.616000   0.542000   0.348000
    3.800   -1.561700   -1.501400   -1.586500   -1.668500   2.193800   -0.018362    0.740280   6.200000   -1.215800   0.101420   -0.000053   7.240000   0.002586   -0.978260    857.070000   -0.004702   -0.001520    1.253400    0.609180   129.490000   195.010000   0.072000   0.000000   0.536000   0.616000   0.543000   0.349000
    4.000   -1.638800   -1.574800   -1.667300   -1.751600   2.204000   -0.014642    0.763030   6.200000   -1.216200   0.102180   -0.000052   7.320000   0.002605   -0.969380    844.480000   -0.003212   -0.001520    1.271100    0.629390   129.490000   199.450000   0.070000   0.000000   0.536000   0.616000   0.543000   0.349000
    4.200   -1.711600   -1.643900   -1.745100   -1.829000   2.212300   -0.012248    0.785520   6.200000   -1.216500   0.102690   -0.000047   7.390000   0.002619   -0.960120    832.450000   -0.002103   -0.001520    1.287000    0.651570   129.570000   204.930000   0.068000   0.000000   0.535000   0.616000   0.542000   0.349000
    4.400   -1.779800   -1.708900   -1.819200   -1.901100   2.218100   -0.011459    0.807920   6.200000   -1.216900   0.103040   -0.000039   7.460000   0.002624   -0.950490    821.180000   -0.001324   -0.001500    1.300400    0.674210   129.710000   211.090000   0.066000   0.000000   0.534000   0.617000   0.540000   0.347000
    4.600   -1.846900   -1.773100   -1.892300   -1.971200   2.223000   -0.011760    0.831260   6.200000   -1.217500   0.103240   -0.000027   7.520000   0.002622   -0.940500    810.790000   -0.000804   -0.001480    1.312500    0.696980   129.870000   217.560000   0.064000   0.000000   0.533000   0.619000   0.538000   0.345000
    4.800   -1.906300   -1.830300   -1.957300   -2.032600   2.226800   -0.012879    0.852400   6.200000   -1.218200   0.103370   -0.000014   7.640000   0.002615   -0.930180    801.410000   -0.000471   -0.001460    1.322500    0.718880   130.050000   223.990000   0.063000   0.000000   0.531000   0.621000   0.535000   0.341000
    5.000   -1.966000   -1.888200   -2.024500   -2.092800   2.229900   -0.014855    0.873140   6.200000   -1.218900   0.103530    0.000000   7.780000   0.002604   -0.919540    793.130000   -0.000255   -0.001440    1.328900    0.738060   130.220000   230.000000   0.061000   0.000000   0.528000   0.622000   0.532000   0.335000
    5.500   -2.105100   -2.023200   -2.190800   -2.228800   2.238900   -0.019502    0.914660   6.200000   -1.220400   0.104600    0.000000   8.070000   0.002591   -0.891760    785.730000    0.000072   -0.001400    1.345000    0.777900   130.390000   241.860000   0.060000   0.000000   0.526000   0.624000   0.528000   0.329000
    6.000   -2.242100   -2.156300   -2.365900   -2.357900   2.237700   -0.026383    0.948700   6.200000   -1.223200   0.107500    0.000000   8.480000   0.002584   -0.862860    779.910000    0.000188   -0.001380    1.349600    0.802780   130.530000   249.340000   0.059000   0.000000   0.524000   0.625000   0.524000   0.321000
    6.500   -2.368600   -2.278500   -2.532200   -2.477000   2.215000   -0.039505    0.976430   6.200000   -1.229900   0.112310    0.000000   8.900000   0.002586   -0.833550    775.600000    0.000159   -0.001370    1.348900    0.814800   130.630000   252.940000   0.059000   0.000000   0.520000   0.634000   0.517000   0.312000
    7.000   -2.482700   -2.388100   -2.681800   -2.585400   2.172000   -0.059140    0.997570   6.200000   -1.240800   0.118530    0.000000   9.200000   0.002603   -0.804570    772.680000    0.000056   -0.001370    1.342200    0.816150   130.700000   253.120000   0.059000   0.000000   0.515000   0.636000   0.514000   0.302000
    7.500   -2.586500   -2.487400   -2.817600   -2.685400   2.118700   -0.081606    1.012100   6.200000   -1.254300   0.125070    0.000000   9.480000   0.002600   -0.776650    771.010000   -0.000055   -0.001370    1.328800    0.809000   130.720000   250.390000   0.058000   0.000000   0.512000   0.634000   0.511000   0.270000
    8.000   -2.686100   -2.582900   -2.943800   -2.782300   2.061300   -0.103820    1.023200   6.200000   -1.268800   0.131460    0.000000   9.570000   0.002630   -0.750330    760.810000   -0.000117   -0.001370    1.308400    0.795380   130.870000   245.230000   0.059000   0.000000   0.510000   0.630000   0.507000   0.278000
    8.500   -2.782000   -2.675200   -3.059700   -2.877600   2.008400   -0.121140    1.033500   6.200000   -1.283900   0.137420    0.000000   9.620000   0.002670   -0.725440    764.500000   -0.000131   -0.001370    1.282200    0.776710   130.710000   238.130000   0.059000   0.000000   0.509000   0.622000   0.503000   0.265000
    9.000   -2.879200   -2.768700   -3.171300   -2.975900   1.960500   -0.134070    1.045300   6.200000   -1.298900   0.142940    0.000000   9.660000   0.002760   -0.701610    768.070000   -0.000108   -0.001370    1.251500    0.754250   130.500000   229.560000   0.060000   0.000000   0.509000   0.613000   0.498000   0.252000
    9.500   -2.976900   -2.863400   -3.278500   -3.076000   1.918900   -0.143640    1.056700   6.200000   -1.313000   0.147810    0.000000   9.660000   0.002890   -0.678500    771.550000   -0.000060   -0.001360    1.217900    0.729250   130.260000   220.020000   0.060000   0.000000   0.509000   0.604000   0.492000   0.239000
    10.00   -3.070200   -2.953700   -3.377600   -3.172600   1.883700   -0.150960    1.065100   6.200000   -1.325300   0.151830    0.000000   9.660000   0.003030   -0.655750    775.000000    0.000000   -0.001360    1.182900    0.703000   130.000000   210.000000   0.060000   0.000000   0.510000   0.604000   0.487000   0.239000
    """) 
[docs]class BooreEtAl2014LowQ(BooreEtAl2014):
    """
    This class implements the Boore et al. (2014) model considering the
    correction to the path scaling term for Low Q regions (e.g. Italy and
    Japan)
    The modification is made to the "Dc3" coefficient
    """
    COEFFS = CoeffsTable(sa_damping=5, table="""\
    IMT            e0          e1          e2          e3         e4          e5          e6     Mh          c1         c2          c3      h          Dc3           c        Vc          f4          f5       f6       f7        R1        R2     DfR     DfV      f1      f2      tau1      tau2
    pgv      5.037000    5.078000    4.849000    5.033000   1.073000   -0.153600    0.225200   6.20   -1.243000   0.148900   -0.003440   5.30   -0.0003300   -0.840000   1300.00   -0.100000   -0.008440   -9.900   -9.900   105.000   272.000   0.082   0.080   0.644   0.552   0.401   0.346
    pga      0.447300    0.485600    0.245900    0.453900   1.431000    0.050530   -0.166200   5.50   -1.134000   0.191700   -0.008088   4.50   -0.0025500   -0.600000   1500.00   -0.150000   -0.007010   -9.900   -9.900   110.000   270.000   0.100   0.070   0.695   0.495   0.398   0.348
    0.010    0.453400    0.491600    0.251900    0.459900   1.421000    0.049320   -0.165900   5.50   -1.134000   0.191600   -0.008088   4.50   -0.0024367   -0.603720   1500.20   -0.148330   -0.007010   -9.900   -9.900   111.670   270.000   0.096   0.070   0.698   0.499   0.402   0.345
    0.020    0.485980    0.523590    0.297070    0.488750   1.433100    0.053388   -0.165610   5.50   -1.139400   0.189620   -0.008074   4.50   -0.0023400   -0.573880   1500.36   -0.147100   -0.007280   -9.900   -9.900   113.100   270.000   0.092   0.030   0.702   0.502   0.409   0.346
    0.022    0.498660    0.536470    0.313470    0.499730   1.433600    0.054888   -0.165200   5.50   -1.140500   0.189240   -0.008095   4.50   -0.0022886   -0.566750   1500.68   -0.148010   -0.007320   -9.900   -9.900   113.370   270.000   0.088   0.027   0.707   0.505   0.418   0.349
    0.025    0.522830    0.561300    0.344260    0.519990   1.432800    0.057529   -0.164990   5.50   -1.141900   0.188750   -0.008153   4.50   -0.0022500   -0.555200   1501.04   -0.150150   -0.007360   -9.900   -9.900   113.070   270.000   0.086   0.026   0.711   0.508   0.427   0.354
    0.029    0.559490    0.599230    0.391460    0.549950   1.427900    0.060732   -0.166320   5.50   -1.142300   0.188440   -0.008290   4.50   -0.0022145   -0.538500   1501.26   -0.153870   -0.007370   -9.900   -9.900   112.360   270.000   0.084   0.028   0.716   0.510   0.436   0.359
    0.030    0.569160    0.609200    0.403910    0.557830   1.426100    0.061444   -0.166900   5.50   -1.142100   0.188420   -0.008336   4.49   -0.0021676   -0.534140   1502.95   -0.154850   -0.007350   -9.900   -9.900   112.130   270.000   0.081   0.029   0.721   0.514   0.445   0.364
    0.032    0.588020    0.628750    0.427880    0.573300   1.422700    0.062806   -0.168130   5.50   -1.141200   0.188400   -0.008445   4.45   -0.0021167   -0.525290   1503.12   -0.156850   -0.007310   -9.900   -9.900   111.650   270.000   0.078   0.030   0.726   0.516   0.454   0.369
    0.035    0.616360    0.658180    0.462520    0.597040   1.417400    0.064559   -0.170150   5.50   -1.138800   0.188390   -0.008642   4.40   -0.0020957   -0.511920   1503.24   -0.160160   -0.007210   -9.900   -9.900   110.640   270.000   0.077   0.031   0.730   0.518   0.462   0.374
    0.036    0.625540    0.667720    0.473380    0.604960   1.415800    0.065028   -0.170830   5.50   -1.137800   0.188370   -0.008715   4.38   -0.0020748   -0.507520   1503.32   -0.161420   -0.007170   -9.900   -9.900   109.530   270.000   0.075   0.031   0.734   0.520   0.470   0.379
    0.040    0.662810    0.706040    0.515320    0.638280   1.409000    0.066183   -0.173570   5.50   -1.132400   0.188160   -0.009030   4.32   -0.0020539   -0.490650   1503.35   -0.167770   -0.006980   -9.900   -9.900   108.280   270.000   0.073   0.032   0.738   0.521   0.478   0.384
    0.042    0.680870    0.724430    0.534450    0.655050   1.405900    0.066438   -0.174850   5.50   -1.129200   0.187970   -0.009195   4.29   -0.0020348   -0.482900   1503.34   -0.171930   -0.006870   -9.900   -9.900   106.990   270.000   0.072   0.032   0.742   0.523   0.484   0.390
    0.044    0.698820    0.742770    0.552820    0.672250   1.403300    0.066663   -0.176190   5.50   -1.125900   0.187750   -0.009360   4.27   -0.0020157   -0.475720   1503.13   -0.176640   -0.006770   -9.900   -9.900   105.410   270.000   0.070   0.031   0.745   0.525   0.490   0.397
    0.045    0.708220    0.752320    0.562220    0.681390   1.402100    0.066774   -0.176930   5.50   -1.124200   0.187640   -0.009441   4.25   -0.0020011   -0.472360   1502.84   -0.179140   -0.006720   -9.900   -9.900   103.610   270.000   0.069   0.031   0.748   0.527   0.496   0.405
    0.046    0.717790    0.762020    0.571660    0.690760   1.400900    0.066891   -0.177690   5.50   -1.122400   0.187520   -0.009521   4.24   -0.0019929   -0.469150   1502.47   -0.181700   -0.006670   -9.900   -9.900   101.700   270.000   0.067   0.031   0.750   0.529   0.499   0.412
    0.048    0.735740    0.780150    0.588880    0.708540   1.399100    0.067127   -0.179200   5.50   -1.119200   0.187300   -0.009676   4.22   -0.0019884   -0.463210   1502.01   -0.186880   -0.006560   -9.900   -9.900    99.760   270.000   0.065   0.031   0.752   0.530   0.502   0.419
    0.050    0.754360    0.799050    0.606520    0.727260   1.397400    0.067357   -0.180820   5.50   -1.115900   0.187090   -0.009819   4.20   -0.0019911   -0.457950   1501.42   -0.192000   -0.006470   -9.900   -9.900    97.930   270.000   0.063   0.030   0.753   0.532   0.503   0.426
    0.055    0.799600    0.844500    0.647700    0.773700   1.394700    0.067797   -0.184800   5.50   -1.108200   0.186550   -0.010120   4.15   -0.0020011   -0.447870   1500.71   -0.203690   -0.006250   -9.900   -9.900    96.030   270.000   0.062   0.029   0.753   0.534   0.502   0.434
    0.060    0.843940    0.888840    0.685620    0.820670   1.395400    0.068591   -0.188580   5.50   -1.100900   0.185820   -0.010330   4.11   -0.0020193   -0.441860   1499.83   -0.213740   -0.006070   -9.900   -9.900    94.100   270.010   0.061   0.027   0.753   0.536   0.499   0.441
    0.065    0.886550    0.931160    0.719410    0.867240   1.400400    0.070127   -0.191760   5.50   -1.094200   0.184850   -0.010480   4.08   -0.0020448   -0.439510   1498.74   -0.222250   -0.005930   -9.900   -9.900    92.080   270.020   0.061   0.025   0.752   0.538   0.495   0.448
    0.067    0.902700    0.947110    0.731710    0.885260   1.403200    0.070895   -0.192910   5.50   -1.091800   0.184420   -0.010520   4.07   -0.0020766   -0.439500   1497.42   -0.225240   -0.005880   -9.900   -9.900    90.010   270.020   0.061   0.025   0.750   0.540   0.489   0.455
    0.070    0.926520    0.970570    0.749400    0.912270   1.408200    0.072075   -0.194510   5.50   -1.088400   0.183690   -0.010560   4.06   -0.0021148   -0.440400   1495.85   -0.229310   -0.005820   -9.900   -9.900    87.970   270.030   0.062   0.024   0.748   0.541   0.483   0.461
    0.075    0.964470    1.007700    0.776780    0.956300   1.417400    0.073549   -0.196650   5.50   -1.083100   0.182250   -0.010580   4.04   -0.0021594   -0.444110   1494.00   -0.235000   -0.005730   -9.900   -9.900    85.990   270.040   0.064   0.022   0.745   0.542   0.474   0.466
    0.080    1.000300    1.042600    0.801610    0.998180   1.426100    0.073735   -0.198160   5.50   -1.078500   0.180520   -0.010560   4.02   -0.0022122   -0.450200   1491.82   -0.239440   -0.005670   -9.900   -9.900    84.230   270.050   0.067   0.020   0.741   0.543   0.464   0.468
    0.085    1.034000    1.075500    0.824230    1.037900   1.432200    0.071940   -0.199020   5.51   -1.074500   0.178560   -0.010510   4.03   -0.0022695   -0.458130   1489.29   -0.242850   -0.005630   -9.900   -9.900    82.740   270.060   0.072   0.019   0.737   0.543   0.452   0.468
    0.090    1.066600    1.107600    0.845910    1.076200   1.435000    0.068097   -0.199290   5.52   -1.070900   0.176430   -0.010420   4.07   -0.0023269   -0.467320   1486.36   -0.245440   -0.005610   -9.900   -9.900    81.540   270.070   0.076   0.017   0.734   0.542   0.440   0.466
    0.095    1.098100    1.138500    0.867030    1.112700   1.433900    0.062327   -0.199000   5.53   -1.067800   0.174200   -0.010320   4.10   -0.0023833   -0.477210   1482.98   -0.247470   -0.005600   -9.900   -9.900    80.460   270.080   0.082   0.016   0.731   0.542   0.428   0.464
    0.100    1.126800    1.166900    0.887100    1.145400   1.429300    0.055231   -0.198380   5.54   -1.065200   0.172030   -0.010200   4.13   -0.0024388   -0.487240   1479.12   -0.249160   -0.005600   -9.900   -9.900    79.590   270.090   0.087   0.014   0.728   0.541   0.415   0.458
    0.110    1.178500    1.217900    0.927020    1.203000   1.411000    0.037389   -0.196010   5.57   -1.060700   0.167700   -0.009964   4.19   -0.0024916   -0.506320   1474.74   -0.252130   -0.005620   -9.900   -9.900    79.050   270.110   0.093   0.012   0.726   0.540   0.403   0.451
    0.120    1.223000    1.262100    0.966160    1.250200   1.383100    0.016373   -0.192650   5.62   -1.057200   0.163520   -0.009722   4.24   -0.0025398   -0.524380   1469.75   -0.254550   -0.005670   -9.900   -9.900    78.850   270.130   0.099   0.011   0.724   0.539   0.392   0.441
    0.130    1.259600    1.298600    1.003100    1.286900   1.349700   -0.005158   -0.188980   5.66   -1.054900   0.159820   -0.009476   4.29   -0.0025844   -0.542140   1464.09   -0.256280   -0.005720   -9.900   -9.900    78.990   270.150   0.104   0.011   0.723   0.538   0.381   0.430
    0.133    1.269200    1.308200    1.013500    1.296100   1.339500   -0.011354   -0.187920   5.67   -1.054500   0.158820   -0.009402   4.30   -0.0026263   -0.547520   1457.76   -0.256650   -0.005740   -9.900   -9.900    79.470   270.150   0.110   0.011   0.722   0.538   0.371   0.417
    0.140    1.288300    1.327000    1.036000    1.313700   1.316200   -0.024711   -0.185660   5.70   -1.053700   0.156720   -0.009228   4.34   -0.0026663   -0.560320   1450.71   -0.257190   -0.005780   -9.900   -9.900    80.260   270.160   0.115   0.012   0.721   0.537   0.362   0.403
    0.150    1.309500    1.348100    1.064800    1.332400   1.284400   -0.042065   -0.182340   5.74   -1.053200   0.154010   -0.008977   4.39   -0.0027063   -0.579620   1442.85   -0.257130   -0.005850   -9.900   -9.900    81.330   270.160   0.120   0.015   0.720   0.537   0.354   0.388
    0.160    1.323500    1.361500    1.087600    1.343700   1.254100   -0.057593   -0.178530   5.78   -1.053300   0.151580   -0.008725   4.44   -0.0027509   -0.600520   1434.22   -0.256040   -0.005910   -9.900   -9.900    82.860   270.160   0.125   0.020   0.720   0.536   0.349   0.372
    0.170    1.330600    1.367900    1.104000    1.348700   1.224400   -0.071861   -0.174210   5.82   -1.054100   0.149480   -0.008472   4.49   -0.0028019   -0.622520   1424.85   -0.254140   -0.005970   -9.900   -9.900    84.720   270.140   0.128   0.026   0.718   0.536   0.346   0.357
    0.180    1.332700    1.368900    1.114900    1.349200   1.194100   -0.085640   -0.169390   5.85   -1.055600   0.147680   -0.008219   4.53   -0.0028547   -0.644860   1414.77   -0.251730   -0.006020   -9.900   -9.900    86.670   270.110   0.131   0.033   0.717   0.536   0.344   0.341
    0.190    1.330700    1.365600    1.120800    1.346300   1.163500   -0.098884   -0.164040   5.89   -1.057900   0.146160   -0.007967   4.57   -0.0029111   -0.666810   1403.99   -0.249110   -0.006080   -9.900   -9.900    88.730   270.060   0.134   0.039   0.714   0.537   0.343   0.324
    0.200    1.325500    1.359000    1.122000    1.341400   1.134900   -0.110960   -0.158520   5.92   -1.060700   0.144890   -0.007717   4.61   -0.0029702   -0.687620   1392.61   -0.246580   -0.006140   -9.900   -9.900    90.910   270.000   0.136   0.045   0.711   0.539   0.344   0.309
    0.220    1.309100    1.339400    1.113300    1.328100   1.082300   -0.133000   -0.147040   5.97   -1.067000   0.142630   -0.007224   4.68   -0.0030285   -0.724310   1380.72   -0.242350   -0.006260   -9.900   -9.900    93.040   269.830   0.138   0.052   0.708   0.541   0.345   0.294
    0.240    1.288100    1.315000    1.094500    1.313200   1.036600   -0.152990   -0.134450   6.03   -1.073700   0.140350   -0.006747   4.75   -0.0030849   -0.756460   1368.51   -0.238230   -0.006380   -9.900   -9.900    95.080   269.590   0.140   0.055   0.703   0.544   0.347   0.280
    0.250    1.276600    1.301700    1.082800    1.305200   1.016600   -0.162130   -0.127840   6.05   -1.077300   0.139250   -0.006517   4.78   -0.0031395   -0.771770   1356.21   -0.235740   -0.006440   -9.900   -9.900    97.040   269.450   0.141   0.055   0.698   0.547   0.350   0.266
    0.260    1.265100    1.288600    1.071000    1.297200   0.999320   -0.170410   -0.121150   6.07   -1.080800   0.138180   -0.006293   4.82   -0.0031896   -0.786970   1343.89   -0.232800   -0.006500   -9.900   -9.900    98.870   269.300   0.141   0.055   0.693   0.550   0.353   0.255
    0.280    1.242900    1.263500    1.047600    1.281500   0.972820   -0.184630   -0.107140   6.11   -1.087900   0.136040   -0.005866   4.88   -0.0032351   -0.816130   1331.67   -0.226010   -0.006600   -9.900   -9.900   100.530   268.960   0.140   0.053   0.687   0.554   0.357   0.244
    0.290    1.232400    1.251700    1.036300    1.273600   0.963480   -0.190570   -0.100110   6.12   -1.091300   0.134990   -0.005666   4.90   -0.0032742   -0.829500   1319.83   -0.222500   -0.006650   -9.900   -9.900   102.010   268.780   0.139   0.051   0.681   0.557   0.360   0.236
    0.300    1.221700    1.240100    1.024600    1.265300   0.956760   -0.195900   -0.092855   6.14   -1.094800   0.133880   -0.005475   4.93   -0.0032969   -0.841650   1308.47   -0.219120   -0.006700   -9.900   -9.900   103.150   268.590   0.138   0.050   0.675   0.561   0.363   0.229
    0.320    1.200700    1.217700    1.001100    1.247900   0.950040   -0.204540   -0.078923   6.16   -1.101300   0.131790   -0.005122   4.98   -0.0033033   -0.861750   1297.65   -0.213180   -0.006800   -9.900   -9.900   104.000   268.200   0.135   0.048   0.670   0.566   0.366   0.223
    0.340    1.179000    1.195500    0.976770    1.228600   0.949560   -0.211340   -0.065134   6.18   -1.107400   0.129840   -0.004808   5.03   -0.0032997   -0.877260   1287.50   -0.208160   -0.006890   -9.900   -9.900   104.700   267.790   0.133   0.047   0.664   0.570   0.369   0.218
    0.350    1.167400    1.183600    0.963800    1.217700   0.950770   -0.214460   -0.057921   6.18   -1.110500   0.128900   -0.004663   5.06   -0.0032869   -0.883750   1278.06   -0.205920   -0.006930   -9.900   -9.900   105.260   267.580   0.130   0.047   0.658   0.573   0.372   0.215
    0.360    1.155800    1.172000    0.951200    1.206600   0.952780   -0.217160   -0.051040   6.19   -1.113300   0.128060   -0.004527   5.08   -0.0032669   -0.889650   1269.19   -0.203790   -0.006970   -9.900   -9.900   105.610   267.370   0.128   0.047   0.653   0.576   0.375   0.212
    0.380    1.130500    1.146800    0.924400    1.181600   0.958990   -0.222140   -0.036755   6.19   -1.119000   0.126470   -0.004276   5.12   -0.0032423   -0.900380   1260.74   -0.199780   -0.007050   -9.900   -9.900   105.870   266.950   0.125   0.048   0.648   0.578   0.378   0.210
    0.400    1.104600    1.121400    0.897650    1.155200   0.967660   -0.226080   -0.023189   6.20   -1.124300   0.125120   -0.004053   5.16   -0.0032123   -0.910920   1252.66   -0.195820   -0.007130   -9.900   -9.900   106.020   266.540   0.122   0.049   0.643   0.580   0.381   0.210
    0.420    1.078200    1.095500    0.870670    1.127600   0.978620   -0.229240   -0.010417   6.20   -1.129100   0.123890   -0.003853   5.20   -0.0031768   -0.922410   1244.80   -0.191710   -0.007190   -9.900   -9.900   106.030   266.160   0.120   0.051   0.638   0.583   0.384   0.210
    0.440    1.051500    1.069700    0.843550    1.099500   0.991440   -0.231660    0.001168   6.20   -1.133700   0.122780   -0.003673   5.24   -0.0031322   -0.934590   1237.03   -0.187470   -0.007260   -9.900   -9.900   105.920   265.800   0.117   0.053   0.634   0.585   0.388   0.211
    0.450    1.037600    1.056200    0.829410    1.084700   0.998760   -0.232630    0.006589   6.20   -1.135900   0.122270   -0.003590   5.25   -0.0030794   -0.940750   1229.23   -0.185340   -0.007290   -9.900   -9.900   105.790   265.640   0.115   0.054   0.629   0.589   0.393   0.213
    0.460    1.023400    1.042600    0.815090    1.069600   1.006400   -0.233500    0.011871   6.20   -1.138100   0.121770   -0.003510   5.27   -0.0030212   -0.946860   1221.16   -0.183210   -0.007320   -9.900   -9.900   105.690   265.480   0.113   0.055   0.624   0.592   0.398   0.216
    0.480    0.997190    1.017200    0.788600    1.041500   1.021500   -0.234640    0.020767   6.20   -1.142000   0.120930   -0.003360   5.30   -0.0029639   -0.958630   1212.74   -0.179020   -0.007380   -9.900   -9.900   105.590   265.210   0.111   0.057   0.619   0.595   0.404   0.219
    0.500    0.969910    0.991060    0.761500    1.012000   1.038400   -0.235220    0.029119   6.20   -1.145900   0.120150   -0.003220   5.34   -0.0029065   -0.969300   1203.91   -0.175000   -0.007440   -9.900   -9.900   105.540   265.000   0.109   0.060   0.615   0.599   0.410   0.224
    0.550    0.904800    0.928300    0.698400    0.941700   1.083300   -0.234490    0.046932   6.20   -1.154300   0.118470   -0.002897   5.41   -0.0028483   -0.989250   1194.59   -0.166010   -0.007580   -9.900   -9.900   105.610   264.740   0.108   0.066   0.610   0.603   0.417   0.229
    0.600    0.841650    0.867150    0.638750    0.873510   1.133600   -0.231280    0.062667   6.20   -1.161500   0.116710   -0.002610   5.48   -0.0027892   -1.001200   1184.93   -0.158300   -0.007730   -9.900   -9.900   105.830   264.830   0.106   0.071   0.605   0.607   0.424   0.235
    0.650    0.781810    0.808760    0.582310    0.809480   1.186100   -0.226660    0.077997   6.20   -1.167600   0.114650   -0.002356   5.53   -0.0027273   -1.007800   1175.19   -0.151440   -0.007870    0.006    0.004   106.200   265.200   0.105   0.073   0.599   0.611   0.431   0.243
    0.667    0.762620    0.789940    0.564220    0.789160   1.203500   -0.224970    0.083058   6.20   -1.169400   0.113940   -0.002276   5.54   -0.0026618   -1.009300   1165.69   -0.149230   -0.007920    0.026    0.017   106.750   265.380   0.103   0.074   0.593   0.615   0.440   0.250
    0.700    0.725130    0.753020    0.528780    0.749850   1.237500   -0.221430    0.093185   6.20   -1.172800   0.112530   -0.002131   5.56   -0.0025953   -1.011700   1156.46   -0.145030   -0.008000    0.055    0.036   107.480   265.780   0.102   0.073   0.587   0.619   0.448   0.258
    0.750    0.669030    0.697370    0.475230    0.691730   1.287100   -0.215910    0.108290   6.20   -1.177700   0.110540   -0.001931   5.60   -0.0025271   -1.015400   1147.59   -0.138660   -0.008120    0.092    0.059   108.390   266.510   0.100   0.070   0.581   0.622   0.457   0.266
    0.800    0.613460    0.641960    0.421730    0.635190   1.334100   -0.210470    0.122560   6.20   -1.181900   0.108730   -0.001754   5.63   -0.0024561   -1.021000   1139.21   -0.132010   -0.008220    0.140    0.088   109.620   267.320   0.099   0.063   0.576   0.624   0.466   0.274
    0.850    0.558530    0.586980    0.368130    0.579690   1.378000   -0.205280    0.136080   6.20   -1.185400   0.107090   -0.001597   5.66   -0.0023787   -1.028200   1131.34   -0.125170   -0.008300    0.195    0.120   111.080   268.140   0.099   0.053   0.570   0.625   0.475   0.281
    0.900    0.502960    0.531360    0.313760    0.523610   1.420800   -0.200110    0.149830   6.20   -1.188400   0.105480   -0.001456   5.69   -0.0022932   -1.036000   1123.91   -0.118310   -0.008360    0.252    0.152   112.710   268.900   0.098   0.042   0.564   0.626   0.483   0.288
    0.950    0.447010    0.475410    0.259190    0.467060   1.462300   -0.194900    0.164320   6.20   -1.190900   0.103890   -0.001328   5.72   -0.0021958   -1.043600   1116.83   -0.111600   -0.008410    0.309    0.181   114.500   269.550   0.098   0.030   0.558   0.626   0.491   0.294
    1.000    0.393200    0.421800    0.207000    0.412400   1.500400   -0.189830    0.178950   6.20   -1.193000   0.102480   -0.001210   5.74   -0.0020894   -1.050000   1109.95   -0.105210   -0.008440    0.367    0.208   116.390   270.000   0.098   0.020   0.553   0.625   0.498   0.298
    1.100    0.284840    0.313740    0.101820    0.302090   1.569000   -0.180010    0.210420   6.20   -1.196600   0.100160   -0.000994   5.82   -0.0019774   -1.057300   1103.07   -0.093837   -0.008470    0.425    0.233   118.300   270.180   0.099   0.007   0.548   0.624   0.505   0.302
    1.200    0.173400    0.202590   -0.006195    0.188660   1.628200   -0.170900    0.244100   6.20   -1.199600   0.098482   -0.000803   5.92   -0.0018619   -1.058400   1096.04   -0.084144   -0.008420    0.481    0.256   120.190   269.420   0.100   0.002   0.543   0.623   0.511   0.306
    1.300    0.061520    0.091060   -0.113450    0.074330   1.679400   -0.162330    0.277990   6.20   -1.201800   0.097375   -0.000635   6.01   -0.0017454   -1.055400   1088.67   -0.075819   -0.008290    0.536    0.276   122.010   267.820   0.101   0.003   0.539   0.622   0.516   0.309
    1.400   -0.045750   -0.015700   -0.215500   -0.036070   1.723900   -0.154130    0.309560   6.20   -1.203900   0.096743   -0.000490   6.10   -0.0016298   -1.050400   1080.77   -0.068543   -0.008060    0.588    0.294   123.750   265.450   0.102   0.006   0.535   0.620   0.521   0.312
    1.500   -0.149540   -0.118660   -0.313800   -0.143700   1.762200   -0.146700    0.338960   6.20   -1.206300   0.096445   -0.000365   6.18   -0.0015179   -1.045400   1072.39   -0.062000   -0.007710    0.638    0.309   125.380   262.410   0.104   0.010   0.532   0.619   0.525   0.315
    1.600   -0.248600   -0.216720   -0.406820   -0.247080   1.795500   -0.139970    0.366160   6.20   -1.208600   0.096338   -0.000259   6.26   -0.0014114   -1.042100   1061.77   -0.055927   -0.007230    0.689    0.324   126.900   258.780   0.105   0.012   0.529   0.618   0.528   0.318
    1.700   -0.341450   -0.308400   -0.492950   -0.344650   1.825900   -0.133610    0.390650   6.20   -1.210600   0.096254   -0.000171   6.33   -0.0013231   -1.040400   1049.29   -0.050286   -0.006660    0.736    0.337   128.140   254.660   0.106   0.012   0.527   0.618   0.530   0.321
    1.800   -0.429750   -0.395580   -0.573880   -0.438180   1.856400   -0.126860    0.412440   6.20   -1.212300   0.096207   -0.000099   6.40   -0.0012531   -1.039700   1036.42   -0.045096   -0.006030    0.780    0.350   129.110   250.110   0.106   0.012   0.526   0.618   0.531   0.323
    1.900   -0.512760   -0.477310   -0.648990   -0.526820   1.886800   -0.119590    0.431510   6.20   -1.214100   0.096255   -0.000042   6.48   -0.0012012   -1.039500   1023.14   -0.040373   -0.005400    0.824    0.364   129.860   245.250   0.106   0.010   0.526   0.618   0.532   0.326
    2.000   -0.586690   -0.550030   -0.714660   -0.606580   1.915200   -0.112370    0.447880   6.20   -1.215900   0.096361    0.000000   6.54   -0.0011703   -1.039200   1009.49   -0.036136   -0.004790    0.871    0.382   130.370   240.140   0.105   0.008   0.526   0.618   0.532   0.329
    2.200   -0.721430   -0.682200   -0.830030   -0.754020   1.968100   -0.098017    0.480240   6.20   -1.219000   0.096497    0.000000   6.66   -0.0011566   -1.036800    995.52   -0.029105   -0.003780    0.920    0.404   130.670   229.550   0.103   0.005   0.527   0.619   0.533   0.332
    2.400   -0.848100   -0.806900   -0.932600   -0.894100   2.017000   -0.083765    0.518730   6.20   -1.220200   0.096198    0.000000   6.73   -0.0011548   -1.032300    981.33   -0.023710   -0.003020    0.969    0.427   130.810   219.050   0.100   0.003   0.528   0.619   0.533   0.335
    2.500   -0.909660   -0.867650   -0.982280   -0.961870   2.040600   -0.076308    0.538830   6.20   -1.220100   0.096106    0.000000   6.77   -0.0011593   -1.029400    966.94   -0.021509   -0.002720    1.017    0.451   130.810   214.040   0.097   0.002   0.530   0.619   0.534   0.337
    2.600   -0.968630   -0.925770   -1.031300   -1.026600   2.062800   -0.068925    0.558100   6.20   -1.219800   0.096136    0.000000   6.81   -0.0011684   -1.026200    952.34   -0.019576   -0.002460    1.060    0.474   130.720   209.320   0.094   0.001   0.531   0.620   0.535   0.340
    2.800   -1.081700   -1.036700   -1.130100   -1.149500   2.101400   -0.055229    0.593940   6.20   -1.218900   0.096667    0.000000   6.87   -0.0011803   -1.019000    937.52   -0.016324   -0.002080    1.099    0.495   130.570   201.080   0.091   0.000   0.532   0.619   0.536   0.342
    3.000   -1.189800   -1.142000   -1.230000   -1.266400   2.132300   -0.043320    0.626940   6.20   -1.217900   0.097638    0.000000   6.93   -0.0011885   -1.011200    922.43   -0.013577   -0.001830    1.135    0.516   130.360   195.000   0.088   0.000   0.534   0.619   0.537   0.344
    3.200   -1.291400   -1.240600   -1.325500   -1.376000   2.154500   -0.034440    0.658110   6.20   -1.216900   0.098649   -0.000023   6.99   -0.0011921   -1.003200    908.79   -0.011030   -0.001670    1.164    0.534   130.130   191.610   0.084   0.000   0.535   0.618   0.538   0.345
    3.400   -1.386000   -1.332200   -1.415000   -1.478600   2.170400   -0.027889    0.687550   6.20   -1.216000   0.099553   -0.000040   7.08   -0.0011866   -0.995120    896.15   -0.008665   -0.001580    1.188    0.551   129.900   190.730   0.081   0.000   0.535   0.618   0.540   0.346
    3.500   -1.433200   -1.377800   -1.459900   -1.529700   2.177500   -0.024997    0.702160   6.20   -1.215600   0.099989   -0.000045   7.12   -0.0011739   -0.991000    883.16   -0.007568   -0.001550    1.211    0.570   129.710   191.110   0.078   0.000   0.536   0.617   0.541   0.347
    3.600   -1.476200   -1.419300   -1.501400   -1.576400   2.183400   -0.022575    0.715230   6.20   -1.215600   0.100430   -0.000049   7.16   -0.0011539   -0.986820    870.05   -0.006537   -0.001540    1.234    0.589   129.560   191.980   0.075   0.000   0.536   0.616   0.542   0.348
    3.800   -1.561700   -1.501400   -1.586500   -1.668500   2.193800   -0.018362    0.740280   6.20   -1.215800   0.101420   -0.000053   7.24   -0.0011229   -0.978260    857.07   -0.004702   -0.001520    1.253    0.609   129.490   195.010   0.072   0.000   0.536   0.616   0.543   0.349
    4.000   -1.638800   -1.574800   -1.667300   -1.751600   2.204000   -0.014642    0.763030   6.20   -1.216200   0.102180   -0.000052   7.32   -0.0010829   -0.969380    844.48   -0.003212   -0.001520    1.271    0.629   129.490   199.450   0.070   0.000   0.536   0.616   0.543   0.349
    4.200   -1.711600   -1.643900   -1.745100   -1.829000   2.212300   -0.012248    0.785520   6.20   -1.216500   0.102690   -0.000047   7.39   -0.0010247   -0.960120    832.45   -0.002103   -0.001520    1.287    0.652   129.570   204.930   0.068   0.000   0.535   0.616   0.542   0.349
    4.400   -1.779800   -1.708900   -1.819200   -1.901100   2.218100   -0.011459    0.807920   6.20   -1.216900   0.103040   -0.000039   7.46   -0.0009464   -0.950490    821.18   -0.001324   -0.001500    1.300    0.674   129.710   211.090   0.066   0.000   0.534   0.617   0.540   0.347
    4.600   -1.846900   -1.773100   -1.892300   -1.971200   2.223000   -0.011760    0.831260   6.20   -1.217500   0.103240   -0.000027   7.52   -0.0008445   -0.940500    810.79   -0.000804   -0.001480    1.313    0.697   129.870   217.560   0.064   0.000   0.533   0.619   0.538   0.345
    4.800   -1.906300   -1.830300   -1.957300   -2.032600   2.226800   -0.012879    0.852400   6.20   -1.218200   0.103370   -0.000014   7.64   -0.0007180   -0.930180    801.41   -0.000471   -0.001460    1.323    0.719   130.050   223.990   0.063   0.000   0.531   0.621   0.535   0.341
    5.000   -1.966000   -1.888200   -2.024500   -2.092800   2.229900   -0.014855    0.873140   6.20   -1.218900   0.103530    0.000000   7.78   -0.0005715   -0.919540    793.13   -0.000255   -0.001440    1.329    0.738   130.220   230.000   0.061   0.000   0.528   0.622   0.532   0.335
    5.500   -2.105100   -2.023200   -2.190800   -2.228800   2.238900   -0.019502    0.914660   6.20   -1.220400   0.104600    0.000000   8.07   -0.0004077   -0.891760    785.73    0.000072   -0.001400    1.345    0.778   130.390   241.860   0.060   0.000   0.526   0.624   0.528   0.329
    6.000   -2.242100   -2.156300   -2.365900   -2.357900   2.237700   -0.026383    0.948700   6.20   -1.223200   0.107500    0.000000   8.48   -0.0002284   -0.862860    779.91    0.000188   -0.001380    1.350    0.803   130.530   249.340   0.059   0.000   0.524   0.625   0.524   0.321
    6.500   -2.368600   -2.278500   -2.532200   -2.477000   2.215000   -0.039505    0.976430   6.20   -1.229900   0.112310    0.000000   8.90   -0.0000364   -0.833550    775.60    0.000159   -0.001370    1.349    0.815   130.630   252.940   0.059   0.000   0.520   0.634   0.517   0.312
    7.000   -2.482700   -2.388100   -2.681800   -2.585400   2.172000   -0.059140    0.997570   6.20   -1.240800   0.118530    0.000000   9.20    0.0001684   -0.804570    772.68    0.000056   -0.001370    1.342    0.816   130.700   253.120   0.059   0.000   0.515   0.636   0.514   0.302
    7.500   -2.586500   -2.487400   -2.817600   -2.685400   2.118700   -0.081606    1.012100   6.20   -1.254300   0.125070    0.000000   9.48    0.0003849   -0.776650    771.01   -0.000055   -0.001370    1.329    0.809   130.720   250.390   0.058   0.000   0.512   0.634   0.511   0.270
    8.000   -2.686100   -2.582900   -2.943800   -2.782300   2.061300   -0.103820    1.023200   6.20   -1.268800   0.131460    0.000000   9.57    0.0007200   -0.750330    760.81   -0.000117   -0.001370    1.308    0.795   130.870   245.230   0.059   0.000   0.510   0.630   0.507   0.278
    8.500   -2.782000   -2.675200   -3.059700   -2.877600   2.008400   -0.121140    1.033500   6.20   -1.283900   0.137420    0.000000   9.62    0.0009400   -0.725440    764.50   -0.000131   -0.001370    1.282    0.777   130.710   238.130   0.059   0.000   0.509   0.622   0.503   0.265
    9.000   -2.879200   -2.768700   -3.171300   -2.975900   1.960500   -0.134070    1.045300   6.20   -1.298900   0.142940    0.000000   9.66    0.0011300   -0.701610    768.07   -0.000108   -0.001370    1.252    0.754   130.500   229.560   0.060   0.000   0.509   0.613   0.498   0.252
    9.500   -2.976900   -2.863400   -3.278500   -3.076000   1.918900   -0.143640    1.056700   6.20   -1.313000   0.147810    0.000000   9.66    0.0013100   -0.678500    771.55   -0.000060   -0.001360    1.218    0.729   130.260   220.020   0.060   0.000   0.509   0.604   0.492   0.239
    10.00   -3.070200   -2.953700   -3.377600   -3.172600   1.883700   -0.150960    1.065100   6.20   -1.325300   0.151830    0.000000   9.66    0.0014900   -0.655750    775.00    0.000000   -0.001360    1.183    0.703   130.000   210.000   0.060   0.000   0.510   0.604   0.487   0.239
    """) 
[docs]def california_basin_model(vs30):
    """
    Returns the centred z1.0 (mu_z1) based on the California model
    (equation 11)
    """
    coeff = 570.94 ** 4.0
    model = -7.15 / 4.0 * np.log(
        (vs30 ** 4.0 + coeff) / (1360.0 ** 4.0 + coeff)) - np.log(1000.)
    return np.exp(model) 
[docs]def japan_basin_model(vs30):
    """
    Returns the centred z1.0 (mu_z1) based on the Japan model
    (equation 12)
    """
    coeff = 412.39 ** 2.0
    model = -5.23 / 2.0 * np.log(
        ((vs30 ** 2.0) + coeff) / ((1360.0 ** 2.0) + coeff)) - np.log(1000.)
    return np.exp(model) 
for cls in (BooreEtAl2014LowQ, BooreEtAl2014HighQ, BooreEtAl2014):
    add_alias(f'{cls.__name__}CaliforniaBasin', cls, region="CAL")
    add_alias(f'{cls.__name__}CaliforniaBasinNoSOF', cls, region="CAL", sof=0)
    add_alias(f'{cls.__name__}JapanBasin', cls, region="JPN")
    add_alias(f'{cls.__name__}JapanBasinNoSOF', cls, region="JPN", sof=0)
    add_alias(f'{cls.__name__}NoSOF', cls, sof=0)