API reference

Site

class pystrata.site.NonlinearProperty(name='', strains=None, values=None, param=None)[source]

Class for nonlinear property with a method for log-linear interpolation.

Parameters
  • name (str, optional) – used for identification

  • strains (numpy.ndarray, optional) – strains for each of the values [decimal].

  • values (numpy.ndarray, optional) – value of the property corresponding to each strain. Damping should be specified in decimal, e.g., 0.05 for 5%.

  • param (str, optional) –

    type of parameter. Possible values are:

    mod_reduc

    Shear-modulus reduction curve

    damping

    Damping ratio curve [decimal]

property strains

Strains [decimal].

property values

Values of either shear-modulus reduction or damping ratio.

property param

Nonlinear parameter name.

class pystrata.site.SoilType(name='', unit_wt=0.0, mod_reduc=None, damping=None)[source]

Soiltype that combines nonlinear behavior and material properties.

Parameters
  • name (str, optional) – used for identification

  • unit_wt (float) – unit weight of the material in [kN/m³]

  • mod_reduc (NonlinearProperty or None) – shear-modulus reduction curves. If None, linear behavior with no reduction is used

  • damping (NonlinearProperty or float) – damping ratio. [decimal] If float, then linear behavior with constant damping is used.

property density

Density of the soil in kg/m³.

property damping_min

Return the small-strain damping.

property is_nonlinear

If nonlinear properties are specified.

class pystrata.site.ModifiedHyperbolicSoilType(name, unit_wt, damping_min, strains=None)[source]
property masing_scaling

Scaling of the Masing damping component.

property curvature

Curvature of the model.

property strain_ref

Reference strain.

class pystrata.site.TwoParamModifiedHyperbolicCoeffs(a: 'float', b1: 'float', b2: 'float', g1: 'float', g2: 'float', d0: 'float', d1: 'float', d: 'float', c: 'float', gd1: 'float', gd2: 'float')[source]
class pystrata.site.TwoParamModifiedHyperbolicSoilType(name: str = '', unit_wt: float = 0, stress_mean: float = 101.3, strains: Optional[Union[_SupportsArray[dtype], _NestedSequence[_SupportsArray[dtype]], bool, int, float, complex, str, bytes, _NestedSequence[Union[bool, int, float, complex, str, bytes]]]] = None, coeffs: Optional[TwoParamModifiedHyperbolicCoeffs] = None)[source]
class pystrata.site.DarendeliSoilType(unit_wt=0.0, plas_index=0, ocr=1, stress_mean=101.3, freq=1, num_cycles=10, damping_min=None, strains=None)[source]

Darendeli (2001) model for fine grained soils.

Parameters
  • unit_wt (float) – unit weight of the material [kN/m³]

  • plas_index (float, default=0) – plasticity index [percent]

  • ocr (float, default=1) – over-consolidation ratio

  • stress_mean (float, default=101.3) – mean effective stress [kN/m²]

  • freq (float, default=1) – excitation frequency [Hz]

  • num_cycles (float, default=10) – number of cycles of loading

  • strains (array_like, default: np.logspace(-6, -1.5, num=20)) – shear strains levels [decimal]

property masing_scaling

Scaling of the Masing damping component.

property strain_ref

reference strain [decimal]

property curvature

Curvature of the model.

class pystrata.site.MenqSoilType(unit_wt=0.0, uniformity_coeff=10, diam_mean=5, stress_mean=101.3, num_cycles=10, damping_min=None, strains=None)[source]

Menq SoilType for gravelly soils.

Parameters
  • unit_wt (float) – unit weight of the material [kN/m³]

  • uniformity_coeff (float, default=10) – uniformity coeffecient (Cᵤ)

  • diam_mean (float, default=5) – mean diameter (D₅₀) [mm]

  • stress_mean (float, default=101.3) – mean effective stress [kN/m²]

  • num_cycles (float, default=10) – number of cycles of loading

  • strains (array_like, default: np.logspace(-4, 0.5, num=20)) – shear strains levels [decimal]

property strain_ref

reference strain [decimal]

property curvature

Curvature of the model.

class pystrata.site.FixedValues(**kwds)[source]

Utility class to store fixed values

class pystrata.site.KishidaSoilType(name='', unit_wt=None, stress_vert=101.3, organic_content=10, lab_consol_ratio=1, strains=None)[source]

Empirical nonlinear model for highly organic soils.

Parameters
  • name (str, optional) – used for identification

  • unit_wt (float or None, default=None) – unit weight of the material [kN/m³]. If None, then unit weight is computed by the empirical model.

  • stress_vert (float) – vertical effective stress [kN/m²]

  • organic_content (float) – organic_content [percent]

  • lab_consol_ratio (float, default=1) – laboratory consolidation ratio. This parameter is included for completeness, but the default value of 1 should be used for field applications.

  • strains (array_like or None) – shear strains levels. If None, a default of np.logspace(-6, 0.5, num=20) will be used. The first strain should be small such that the shear modulus reduction is equal to 1. [decimal]

class pystrata.site.Layer(soil_type, thickness, shear_vel)[source]

Docstring for Layer

property depth

Depth to the top of the layer [m].

property depth_mid

Depth to the middle of the layer [m].

property depth_base

Depth to the base of the layer [m].

classmethod duplicate(other)[source]

Create a copy of the layer.

property density

Density of soil in [kg/m³].

property damping

Strain-compatible damping.

property initial_shear_mod

Initial (small-strain) shear modulus [kN/m²].

property initial_shear_vel

Initial (small-strain) shear-wave velocity [m/s].

property comp_shear_mod

Strain-compatible complex shear modulus [kN/m²].

property comp_shear_vel

Strain-compatible complex shear-wave velocity [m/s].

property shear_mod

Strain-compatible shear modulus [kN//m²].

property shear_vel

Strain-compatible shear-wave velocity [m/s].

property stress_shear_eff

Effective shear stress at layer midpoint

property stress_shear_max

Maximum shear stress at layer midpoint

property travel_time

Travel time through the layer.

stress_vert(depth_within=0, effective=False)[source]

Vertical stress from the top of the layer [kN//m²].

stress_mean(depth_within=0, effective=False, k0=0.5)[source]

Mean effective stress from the top of the layer [kN//m²].

class pystrata.site.Location(index, layer, wave_field, depth_within=0)[source]

Location within a profile

class pystrata.site.Profile(layers=None, wt_depth=0)[source]

Soil profile with an infinite halfspace at the base.

classmethod from_dataframe(df, wt_depth=0)[source]

Create a profile based on a table with columns: - thickness (m) - vel_shear (m) - unit_wt (kN/m³) - damping (dec)

reset_layers()[source]

Set initial properties from the soil types.

auto_discretize(max_freq: float = 50.0, wave_frac: float = 0.2, nonlinear_only: bool = True) Profile[source]

Subdivide the layers to capture strain variation.

Parameters
  • max_freq (float) – Maximum frequency of interest [Hz].

  • wave_frac (float) – Fraction of wavelength required. Typically 1/3 to 1/5.

  • max_thick (float optional) – If provided, layers are limited to be at most that thick. This is applied to all layers regardless of nonlinearity.

Returns

profile – A new profile with modified layer thicknesses

Return type

Profile

pore_pressure(depth)[source]

Pore pressure at a given depth in [kN//m²].

Parameters

depth

Return type

pore_pressure

lookup_depth(depth: float) Tuple[int, float][source]

Look up the layer and the depth within the layer for a specific depth.

Parameters

depth (float) – Depth corresponding to the location of interest.

Returns

  • index (int) – Layer index

  • depth_within (float) – Depth from the top of the layer to achieve the specific depth.

location(wave_field, depth=None, index=None)[source]

Create a Location for a specific depth.

Parameters
  • wave_field (str) – Wave field. See Location for possible values.

  • depth (float, optional) – Depth corresponding to the :class`Location` of interest. If provided, then index is ignored.

  • index (int, optional) –

    Index corresponding to layer of interest in Profile. If

    provided, then depth is ignored and location is provided a top of layer.

Returns

Corresponding Location object.

Return type

Location

time_average_vel(depth)[source]

Calculate the time-average velocity.

Parameters

depth (float) – Depth over which the average velocity is computed.

Returns

avg_vel – Time averaged velocity.

Return type

float

simplified_rayliegh_vel()[source]

Simplified Rayliegh velocity of the site.

This follows the simplifications proposed by Urzua et al. (2017)

Returns

rayleigh_vel – Equivalent shear-wave velocity.

Return type

float

Variation

class pystrata.variation.TruncatedNorm(limit)[source]

Truncated normal random number generator.

Parameters

limit (float) – Standard normal limits to impose

class pystrata.variation.ToroThicknessVariation(c_1=10.86, c_2=- 0.89, c_3=1.98)[source]

Toro (1995) [T95] thickness variation model.

The recommended values are provided as defaults to this model.

References

T95

Toro, G. R. (1995). Probabilistic models of site velocity profiles for generic and site-specific ground-motion amplification studies. Brookhaven National Laboratory Technical Report: 779574.

Parameters
  • c_1 (float, optional) – \(c_1\) model parameter.

  • c_2 (float, optional) – \(c_2\) model parameter.

  • c_3 (float, optional) – \(c_3\) model parameter.

iter_thickness(depth_total)[source]

Iterate over the varied thicknesses.

The layering is generated using a non-homogenous Poisson process. The following routine is used to generate the layering. The rate function, \(\lambda(t)\), is integrated from 0 to t to generate cumulative rate function, \(\Lambda(t)\). This function is then inverted producing \(\Lambda^{-1}(t)\). Random variables are produced using the a exponential random variation with \(\mu = 1\) and converted to the nonhomogenous variables using the inverted function.

Parameters

depth_total (float) – Total depth generated. Last thickness is truncated to achieve this depth.

Yields

float – Varied thickness.

class pystrata.variation.VelocityVariation(vary_bedrock=False)[source]

Abstract model for varying the velocity.

class pystrata.variation.ToroVelocityVariation(ln_std: float, rho_0: float, delta: float, rho_200: float, h_0: float, b: float, vary_bedrock: bool = False)[source]

Toro (1995) [T95] velocity variation model.

Default values can be selected with generic_model().

Parameters
  • ln_std (float, optional) – \(\sigma_{ln}\) model parameter.

  • rho_0 (float, optional) – \(ρ_0\) model parameter.

  • delta (float, optional) – \(\Delta\) model parameter.

  • rho_200 (float, optional) – \(ρ_200\) model parameter.

  • h_0 (float, optional) – \(h_0\) model parameter.

  • b (float, optional) – \(b\) model parameter.

classmethod generic_model(site_class, **kwds)[source]

Use generic model parameters based on site class.

Parameters

site_class (str) –

Site classification. Possible options are:
  • Geomatrix AB

  • Geomatrix CD

  • USGS AB

  • USGS CD

  • USGS A

  • USGS B

  • USGS C

  • USGS D

See the report for definitions of the Geomatrix site classication. USGS site classification is based on \(V_{s30}\):

Site Class

\(V_{s30}\) (m/s)

A

>750 m/s

B

360 to 750 m/s

C

180 to 360 m/s

D

<180 m/s

Returns

Initialized ToroVelocityVariation with generic parameters.

Return type

ToroVelocityVariation

class pystrata.variation.DepthDependToroVelVariation(depth: Union[_SupportsArray[dtype], _NestedSequence[_SupportsArray[dtype]], bool, int, float, complex, str, bytes, _NestedSequence[Union[bool, int, float, complex, str, bytes]]], ln_std: Union[_SupportsArray[dtype], _NestedSequence[_SupportsArray[dtype]], bool, int, float, complex, str, bytes, _NestedSequence[Union[bool, int, float, complex, str, bytes]]], rho_0: float, delta: float, rho_200: float, h_0: float, b: float, ln_std_map=typing.Optional[typing.Dict[str, float]], vary_bedrock=False)[source]

Toro (1995) [T95] velocity variation model modified for a depth dependent standard deviation that can be overridden by the soil_type name.

Default values can be selected with generic_model().

Parameters
  • depth (array_like, optional) – Depths defining the standard deviation model. Default is [0, 15] following the SPID model.

  • ln_std (array_like, optional) – \(\sigma_{ln}\) model parameter. Default is [0.25, 0.15] following the SPID model.

  • rho_0 (float, optional) – \(ρ_0\) model parameter.

  • delta (float, optional) – \(\Delta\) model parameter.

  • rho_200 (float, optional) – \(ρ_200\) model parameter.

  • h_0 (float, optional) – \(h_0\) model parameter.

  • b (float, optional) – \(b\) model parameter.

  • ln_std_map (dict[str, float], optional) – Mapping between the soil_type and the defined ln_std. Default is None.

classmethod generic_model(site_class, /, *, ln_std_map=None, **kwds)[source]

Use generic model parameters based on site class.

Parameters
  • site_class (str) –

    Site classification. Possible options are:
    • Geomatrix AB

    • Geomatrix CD

    • USGS AB

    • USGS CD

    • USGS A

    • USGS B

    • USGS C

    • USGS D

    See the report for definitions of the Geomatrix site classication. USGS site classification is based on \(V_{s30}\):

    Site Class

    \(V_{s30}\) (m/s)

    A

    >750 m/s

    B

    360 to 750 m/s

    C

    180 to 360 m/s

    D

    <180 m/s

  • ln_std_map (dict[str, float], optional) – Mapping between the soil_type and the defined ln_std. Default is None.

Returns

Initialized DepthDependToroVelVariation with generic parameters.

Return type

DepthAndSoilTypeDependToroVelVariation

class pystrata.variation.SpidVariation(correlation, limits_mod_reduc=[0, 1], limits_damping=[0, 0.15], std_mod_reduc=0.15, std_damping=0.3)[source]

Variation defined by the EPRI SPID (2013) and documented in PNNL (2014).

Motion

Classes used to define input motions.

class pystrata.motion.WaveField(value)[source]

An enumeration.

class pystrata.motion.RvtMotion(freqs, fourier_amps, duration=None, peak_calculator=None, calc_kwds=None)[source]
class pystrata.motion.CompatibleRvtMotion(osc_freqs, osc_accels_target, duration=None, osc_damping=0.05, event_kwds=None, window_len=None, peak_calculator=None, calc_kwds=None)[source]
class pystrata.motion.SourceTheoryRvtMotion(magnitude, distance, region, stress_drop=None, depth=8, peak_calculator=None, calc_kwds=None)[source]

Propagation

class pystrata.propagation.QuarterWaveLenCalculator(site_atten=None)[source]

Compute quarter-wave length site amplification.

No consideration for nolninearity is made by this calculator.

fit(target_type, target, adjust_thickness=False, adjust_site_atten=False, adjust_source_vel=False)[source]

Fit to a target crustal amplification or site term.

The fitting process adjusts the velocity, site attenuation, and layer thickness (if enabled) to fit a target values. The frequency range is specified by the input motion.

Parameters
  • target_type (str) –

    Options are ‘crustal_amp’ to only fit to the crustal amplification,

    or ‘site_term’ to fit both the velocity and the site attenuation parameter.

  • target (array_like) – Target values.

  • adjust_thickness (bool (optional)) – If the thickness of the layers is adjusted as well, default: False.

  • adjust_site_atten (bool (optional)) – If the site attenuation is adjusted as well, default: False.

  • adjust_source_vel (bool (optional)) – If the source velocity should be adjusted, default: False.

Returns

profile – profile optimized to fit a target amplification.

Return type

pyrsa.site.Profile

class pystrata.propagation.LinearElasticCalculator[source]

Class for performing linear elastic site response.

wave_at_location(l)[source]

Compute the wave field at specific location.

Parameters

l (site.Location) – site.Location of the input

Returns

Amplitude and phase of waves

Return type

np.ndarray

calc_accel_tf(lin, lout)[source]

Compute the acceleration transfer function.

Parameters
  • lin (Location) – Location of input

  • lout (Location) – Location of output. Note that this would typically be midheight of the layer.

calc_stress_tf(lin, lout, damped)[source]

Compute the stress transfer function.

Parameters
  • lin (Location) – Location of input

  • lout (Location) – Location of output. Note that this would typically be midheight of the layer.

calc_strain_tf(lin, lout)[source]

Compute the strain transfer function from lout to location_in.

The strain transfer function from the acceleration at layer n (outcrop) to the mid-height of layer m (within) is defined as

Parameters
  • lin (Location) – Location of input

  • lout (Location) – Location of output. Note that this would typically be midheight of the layer.

Returns

strain_tf – Transfer function to be applied to an acceleration FAS.

Return type

numpy.ndarray

class pystrata.propagation.EquivalentLinearCalculator(strain_ratio=0.65, tolerance=0.01, max_iterations=15, strain_limit=0.05)[source]

Class for performing equivalent-linear elastic site response.

classmethod calc_strain_ratio(mag)[source]

Compute the effective strain ratio using Idriss and Sun (1992).

Parameters

mag (float) – Magnitude of the input motion.

Returns

strain_ratio – Effective strain ratio

Return type

float

References

1

Idriss, I. M., & Sun, J. I. (1992). SHAKE91: A computer program for conducting equivalent linear seismic response analyses of horizontally layered soil deposits. Center for Geotechnical Modeling, Department of Civil and Environmental Engineering, University of California, Davis, CA.

class pystrata.propagation.FrequencyDependentEqlCalculator(use_smooth_spectrum=False, strain_ratio=1.0, tolerance=0.01, max_iterations=15)[source]

Class for performing equivalent-linear elastic site response with frequency-dependent modulii and damping.

Parameters
  • use_smooth_spectrum (bool, default=False) – Use the Kausel & Assimaki (2002) smooth spectrum for the strain. Otherwise, the complete Fourier amplitude spectrum is used.

  • strain_ratio (float, default=1.00) – ratio between the maximum strain and effective strain used to compute strain compatible properties. There is not clear guidance the use of the effective strain ratio. However, given the nature of the method, it would make sense not to include the an effective strain ratio.

  • tolerance (float, default=0.01) – tolerance in the iterative properties, which would cause the iterative process to terminate.

  • max_iterations (int, default=15) – maximum number of iterations to perform.

References

1

Kausel, E., & Assimaki, D. (2002). Seismic simulation of inelastic soils via frequency-dependent moduli and damping. Journal of Engineering Mechanics, 128(1), 34-47.

Output

pystrata.output.append_arrays(many, single)[source]

Append an array to another padding with NaNs for constant length.

Parameters
  • many (array_like of rank (j, k)) – values appended to a copy of this array. This may be a 1-D or 2-D array.

  • single (array_like of rank l) – values to append. This should be a 1-D array.

Returns

append – 2-D array with rank (j + 1, max(k, l)) with missing values padded with numpy.nan

Return type

numpy.ndarray

Tools

pystrata.tools.to_str(s)[source]

Parse a string and strip the extra characters.

pystrata.tools.to_float(s)[source]

Try to parse a float.

pystrata.tools.parse_fixed_width(types, lines)[source]

Parse a fixed width line.

pystrata.tools.split_line(line, parsers, sep=' ')[source]

Split a line into pieces and parse the strings.

pystrata.tools.read_nrattle_ctl(fpath)[source]

Read an nrattle control file.