disba
is a computationally efficient Python library for the modeling of surface wave dispersion that implements a subset of codes from Computer Programs in Seismology (CPS) in Python compiled just-in-time with numba
. Such implementation alleviates the usual prerequisite for a Fortran compiler needed by other libraries also based on CPS (e.g. pysurf96
, srfpython
and PyLayeredModel
) which often leads to further installation troubleshooting, especially on Windows platform.
disba
aims to be lightweight and portable without compromising on the performance. For instance, it yields similar speed compared to CPS's surf96 program compiled with f2py
for Rayleigh-wave but is significantly faster for Love-wave with increasing number of layers. disba
also implements the fast delta matrix algorithm for Rayleigh-wave which, albeit ironically slower, is more robust and handles reversion of phase velocity caused by low velocity zones.
Forward modeling:
- Compute Rayleigh-wave phase or group dispersion curves using Dunkin's matrix or fast delta matrix algorithms,
- Compute Love-wave phase or group dispersion curves using Thomson-Haskell method,
- Compute Rayleigh-wave ellipticity.
Eigenfunctions and sensitivity kernels:
- Compute Rayleigh- and Love- wave eigenfunctions,
- Compute Rayleigh- and Love- wave phase or group sensitivity kernels with respect to layer thickness, P- and S- wave velocities, and density.
The recommended way to install disba
and all its dependencies is through the Python Package Index:
pip install disba --user
Otherwise, clone and extract the package, then run from the package location:
pip install . --user
The following example computes the Rayleigh- and Love- wave phase velocity dispersion curves for the 3 first modes.
import numpy
from disba import PhaseDispersion
# Velocity model
# thickness, Vp, Vs, density
# km, km/s, km/s, g/cm3
velocity_model = numpy.array([
[10.0, 7.00, 3.50, 2.00],
[10.0, 6.80, 3.40, 2.00],
[10.0, 7.00, 3.50, 2.00],
[10.0, 7.60, 3.80, 2.00],
[10.0, 8.40, 4.20, 2.00],
[10.0, 9.00, 4.50, 2.00],
[10.0, 9.40, 4.70, 2.00],
[10.0, 9.60, 4.80, 2.00],
[10.0, 9.50, 4.75, 2.00],
])
# Periods must be sorted starting with low periods
t = numpy.logspace(0.0, 3.0, 100)
# Compute the 3 first Rayleigh- and Love- wave modal dispersion curves
# Fundamental mode corresponds to mode 0
pd = PhaseDispersion(*velocity_model.T)
cpr = [pd(t, mode=i, wave="rayleigh") for i in range(3)]
cpl = [pd(t, mode=i, wave="love") for i in range(3)]
# pd returns a namedtuple (period, velocity, mode, wave, type)
Likewise, GroupDispersion
can be used for group velocity.
disba
's API is consistent across all its classes which are initialized and called in the same fashion. Thus, eigenfunctions are calculated as follow:
from disba import EigenFunction
eigf = EigenFunction(*velocity_model.T)
eigr = eigf(20.0, mode=0, wave="rayleigh")
eigl = eigf(20.0, mode=0, wave="love")
# eigf returns a namedtuple
# - (depth, ur, uz, tz, tr, period, mode) for Rayleigh-wave
# - (depth, uu, tt, period, mode) for Love-wave
And sensitivity kernels:
from disba import PhaseSensitivity
ps = PhaseSensitivity(*velocity_model.T)
parameters = ["thickness", "velocity_p", "velocity_s", "density"]
skr = [ps(20.0, mode=0, wave="rayleigh", parameter=parameter) for parameter in parameters]
skl = [ps(20.0, mode=0, wave="love", parameter=parameter) for parameter in parameters]
# ps returns a namedtuple (depth, kernel, period, velocity, mode, wave, type, parameter)