
Classes
splitpy
defines the following base classes:
The class Split
contains attributes
and methods for the analysis of teleseismic shear-wave splitting
from three-component seismograms.
The class PickPlot
defines figure handles
for a picking window showing the seismograms and the predicted teleseismic
shear-wave phase arrivals. This figure is interactive and new picks can
be generated to refine the analysis.
The class DiagPlot
defines figure handles
for a diagnostic figure showing a summary of the splitting results. It can
be called after each application of the split.analyze method to show
the summary of the analysis as a figure. This figure can also be saved as
a .png file.
Split
- class splitpy.classes.Split(sta, zcomp='Z')[source]
A Split object contains dictionary attributes that associate station information with single event (i.e., earthquake) metadata, corresponding raw and rotated seismograms and splitting results.
Note
The object is initialized with the
sta
field only, and other attributes are added to the object as the analysis proceeds.- sta
Object containing station information from a
stdb
database.- Type:
StDbElement
- dataZNE
Object containing raw ZNE trace data
- Type:
Stream
- dataLQT
Object containing rotated trace data
- Type:
Stream
- dataZ12
Object containing (resampled) Z12 trace data (for un-oriented data; optional)
- Type:
Stream
- add_event(event, gacmin=85.0, gacmax=120.0, phase='SKS', returned=False)[source]
Adds event metadata to Split object as Meta object.
- Parameters:
event (
event
) – Event metadata
- download_data(client, new_sr=5.0, dts=120.0, returned=False, verbose=False, remove_response=False)[source]
Downloads seismograms based on event origin time and P phase arrival and adds as object attributes.
- Parameters:
client (
Client
) – Client objectnew_sr (float) – New sampling rate (Hz)
dts (float) – Time duration (sec)
remove_response (bool) – Remove instrument response from seismogram and resitute to true ground velocity (m/s) using obspy.core.trace.Trace.remove_response()
returned (bool) – Whether or not to return the
accept
attributeverbose (bool) – Output diagnostics to screen
- Returns:
accept – Whether or not the object is accepted for further analysis
- Return type:
bool
- rotate(align=None)[source]
Rotates 3-component seismograms from vertical (Z), east (E) and north (N) to longitudinal (L), radial (Q) and tangential (T) components of motion. Note that the method ‘rotate’ from
obspy.core.stream.Stream
is used for the rotation'ZNE->LQT'
.Can also rotate Z12 to ZNE.
- Parameters:
align (str) – Alignment of coordinate system for rotation (‘ZNE’ or ‘LQT’)
- calc_snr(t1=None, dt=30.0, fmin=0.02, fmax=0.5)[source]
Calculates signal-to-noise ratio on either Z, L or P component
- Parameters:
t1 (
UTCDateTime
) – Pick time of arrival (sec)dt (float) – Duration (sec)
fmin (float) – Minimum frequency corner for SNR filter (Hz)
fmax (float) – Maximum frequency corner for SNR filter (Hz)
- snrq
Signal-to-noise ratio on vertical component (dB)
- Type:
float
- snrh
Signal-to-noise ratio on radial component (dB)
- Type:
float
- analyze(t1=None, t2=None, maxdt=4.0, ddt=0.1, dphi=1.0, verbose=False)[source]
Calculates the shear-wave splitting parameters based on two alternative method: the Rotation-Correlation (RC) method and the Silver-Chan (SC) method. Each set of results is stored in a Dictionary as attributes of the split object.
- Parameters:
t1 (
UTCDateTime
) – Start time of picking windowt2 (
UTCDateTime
) – End time of picking windowverbose (bool) – Output diagnostics to screen
- is_null(snrTlim=3.0, verbose=False)[source]
Determines if splitting result is a Null result
- Parameters:
snrTlim (float) – Threshold for snr on T component
verbose (bool) – Output diagnostics to screen
- null
Boolean for Null result
- Type:
bool
- get_quality(verbose=False)[source]
Determines the quality of the estimate (either Null or non-Null) based on ratio of delay times and difference in fast axis directions between Rotation-Correlation and Silver-Chan methods
- Parameters:
verbose (bool) – Output diagnostics to screen
- quality
String representing quality of estimate (‘Good’, ‘Fair’, ‘Poor’)
- Type:
str
- display_results(ds=0)[source]
Prints out best fitting results to screen
- Parameters:
ds (int) – Number of spaces to print out to screen (verbiage)
- display_meta(ds=0)[source]
Prints out content of metadata to screen
- Parameters:
ds (int) – Number of spaces to print out to screen (verbiage)
Meta
- class splitpy.classes.Meta(sta, event, gacmin=85.0, gacmax=180.0, phase='SKS', maxdt=4.0, ddt=0.1, dphi=1.0)[source]
A Meta object contains attributes associated with the station-event data for a single Teleseismic event.
- time
Origin time of earthquake
- Type:
UTCDateTime
- dep
Depth of hypocenter (km)
- Type:
float
- lon
Longitude coordinate of epicenter
- Type:
float
- lat
Latitude coordinate of epicenter
- Type:
float
- mag
Magnitude of earthquake
- Type:
float
- gac
Great arc circle between station and epicenter (degrees)
- Type:
float
- epi_dist
Epicentral distance between station and epicenter (km)
- Type:
float
- baz
Back-azimuth - pointing to earthquake from station (degrees)
- Type:
float
- az
Azimuth - pointing to station from earthquake (degrees)
- Type:
float
- ttime
Predicted arrival time (sec)
- Type:
float
- ph
Phase name
- Type:
str
- slow
Horizontal slowness of phase
- Type:
float
- inc
Incidence angle of phase at surface
- Type:
float
- maxdt
Maximum delay time considered in grid search (sec)
- Type:
float
- ddt
Delay time interval in grid search (sec)
- Type:
float
- dphi
Angular interval in grid search (deg)
- Type:
float
- align
Alignment of coordinate system for rotation (‘ZRT’, ‘LQT’, or ‘PVH’)
- Type:
str
- rotated
Whether or not data have been rotated to
align
coordinate system- Type:
bool
- zcomp
Vertical Component Identifier. Should be a single character. This is different then ‘Z’ only for fully unknown component orientation (i.e., components are 1, 2, 3)
- Type:
str
Result
- class splitpy.classes.Result(Emat, trQ_c, trT_c, trFast, trSlow, phi, dtt, phi_min, edtt, ephi, errc)[source]
A Result object contains attributes associated with the result of a single splitting analysis. These are equally applicable to the RC or SC method - see
analyze()
.- Emat
Error minimization matrix
- Type:
ndarray
- trQ_c
Corrected radial (Q) component
- Type:
Trace
- trT_c
Corrected transverse (T) component
- Type:
Trace
- trFast
Corrected Fast component
- Type:
Trace
- trSlow
Corrected Slow component
- Type:
Trace
- phi
Azimuth of fast axis (deg)
- Type:
float
- dtt
Delay time between fast and slow axes (sec)
- Type:
float
- phi_min
Azimuth used in plotting method
- Type:
float
- ephi
Error on azimuth of fast axis (deg)
- Type:
float
- edtt
Error on delay time between fast and slow axes (sec)
- Type:
float
- errc
Error contours on Emat
- Type:
float
PickPlot
- class splitpy.classes.PickPlot(split)[source]
A PickPlot object contains figure handles and method to plot seismic data for picking/refining the SKS time window. The figure displays the LQT seismograms and predicted arrival times for common shear-wave arrivals in the time window. The figure can be picked to refine the time window for calculating the splitting results.
Note
The object is initialized with a
Split
object, which is temporarily stored as an attribute. All the split attributes are therefore available to thePickPlot
object for plotting.- fp
List of figure handles ([fig, ax1, ax2, ax3])
- Type:
List
- plot_LQT_phases(dts, t1=None, t2=None)[source]
Plots rotated three-components of motion for picking.
- Parameters:
tt (
TauPyModel
) – Taup object containing travel time and phase infodts (float) – Time interval (?)
t1 (
UTCDateTime
) – Start time of picking windowt2 (
UTCDateTime
) – End time of picking window
- ll
List of
axvline
objects for plotting the start and end times of the picking window.- Type:
List
DiagPlot
- class splitpy.classes.DiagPlot(split)[source]
A DiagPlot object contains figure handles and methods to plot the diagnostic figure, which displays the LQT seismograms, the corrected/un-corrected seismograms, the particle motions, the minimization matrix and a text box with a summary of the analysis - for each of the two analysis methods (‘RC’ or ‘SC’)
Note
The object is initialized with a
Split
object, which is temporarily stored as an attribute. All the split attributes are therefore available to theDiagPlot
object for plotting.- fd
List of figure handles ([fig, ax0, axt, axRC1, axRC2, axRC3, axRC4, axSC1, axSC2, axSC3, axSC4])
- Type:
List
GUI classes
Module containing the main GUI widget classes based on PyQt5
:
These classes simply open a Qt5 message box and contain a boolean attribute corresponding the anwer to simple yes/no buttons.
Modules
calc
- splitpy.calc.split_SilverChan(trQ, trT, baz, t1, t2, maxdt, ddt, dphi)[source]
Calculates splitting based on the minimization of energy on the corrected transverse component (Silver and Chan, 1990)
- Parameters:
trQ (
Trace
) – Radial component seismogramtrT (
Trace
) – Tangential component seismogrambaz (float) – Back-azimuth - pointing to earthquake from station (degrees)
t1 (
UTCDateTime
) – Start time of picking windowt2 (
UTCDateTime
) – End time of picking windowmaxdt (float) – Maximum delay time considered in grid search (sec)
ddt (float) – Delay time interval in grid search (sec)
dphi (float) – Angular interval in grid search (deg)
- Returns:
Ematrix (
ndarray
) – Matrix of T component energytrQ_c (
Trace
) – Trace of corrected radial component of motiontrT_c (
Trace
) – Trace of corrected tangential component of motiontrFast (
Trace
) – Trace of corrected fast direction of motiontrSlow (
Trace
) – Trace of corrected slow direction of motion phiSC : float Azimuth of fast axis (deg)dttSC (foat) – Delay time between fast and slow axes (sec)
phi_min (float) – Azimuth used in plotting routine
- splitpy.calc.split_RotCorr(trQ, trT, baz, t1, t2, maxdt, ddt, dphi)[source]
Calculates splitting based on the maximum correlation between corrected radial and tangential components of motion
- Parameters:
trQ (
Trace
) – Radial component seismogramtrT (
Trace
) – Tangential component seismogrambaz (float) – Back-azimuth - pointing to earthquake from station (degrees)
t1 (
UTCDateTime
) – Start time of picking windowt2 (
UTCDateTime
) – End time of picking windowmaxdt (float) – Maximum delay time considered in grid search (sec)
ddt (float) – Delay time interval in grid search (sec)
dphi (float) – Angular interval in grid search (deg)
- Returns:
Ematrix (
ndarray
) – Matrix of T component energytrQ_c (
Trace
) – Trace of corrected radial component of motiontrT_c (
Trace
) – Trace of corrected tangential component of motiontrFast (
Trace
) – Trace of corrected fast direction of motiontrSlow (
Trace
) – Trace of corrected slow direction of motion phiSC : float Azimuth of fast axis (deg)dttSC (foat) – Delay time between fast and slow axes (sec)
phi_min (float) – Azimuth used in plotting routine
- splitpy.calc.tshift(trace, tt)[source]
Shifts a
Trace
object- Parameters:
trace (
Trace
) – Seismogram to apply shifttt (float) – Lag time for shifting
- Returns:
rtrace – Shifted version of trace
- Return type:
Trace
- splitpy.calc.split_dof(tr)[source]
Determines the degrees of freedom to calculate the confidence region of the misfit function. From Walsh, JGR, (2013)
- Parameters:
tr (
Trace
) – Seismogram- Returns:
dof – Degrees of freedom
- Return type:
float
- splitpy.calc.split_errorSC(tr, t1, t2, q, Emat, maxdt, ddt, dphi)[source]
Calculate error bars based on a F-test and a given confidence interval q
- Parameters:
tr (
Trace
) – Seismogramt1 (
UTCDateTime
) – Start time of picking windowt2 (
UTCDateTime
) – End time of picking windowq (float) – Confidence level
Emat (
ndarray
) – Energy minimization matrixmaxdt (float) – Maximum delay time considered in grid search (sec)
ddt (float) – Delay time interval in grid search (sec)
dphi (float) – Angular interval in grid search (deg)
- Returns:
err_dtt (float) – Error in dt estimate (sec)
err_phi (float) – Error in phi estimate (degrees)
err_contour (
ndarray
) – Error contour for plotting
- splitpy.calc.split_errorRC(tr, t1, t2, q, Emat, maxdt, ddt, dphi)[source]
Calculates error bars based on a F-test and a given confidence level q.
Note
This version uses a Fisher transformation for correlation-type misfit.
- Parameters:
tr (
Trace
) – Seismogramt1 (
UTCDateTime
) – Start time of picking windowt2 (
UTCDateTime
) – End time of picking windowq (float) – Confidence level
Emat (
ndarray
) – Energy minimization matrixmaxdt (float) – Maximum delay time considered in grid search (sec)
ddt (float) – Delay time interval in grid search (sec)
dphi (float) – Angular interval in grid search (deg)
- Returns:
err_dtt (float) – Error in dt estimate (sec)
err_phi (float) – Error in phi estimate (degrees)
err_contour (
ndarray
) – Error contour for plotting
- splitpy.calc.split_error_average(q, Emat, maxdt, ddt, dphi, n)[source]
Calculate error bars based on a F-test and a given confidence level q
- Parameters:
q (float) – Confidence level
Emat (
ndarray
) – Error surface matrixmaxdt (float) – Maximum delay time considered in grid search (sec)
ddt (float) – Delay time interval in grid search (sec)
dphi (float) – Angular interval in grid search (deg)
n (int) – Number of measurements
- Returns:
err_dtt (float) – Error in dt estimate (sec)
err_phi (float) – Error in phi estimate (degrees)
err_contour (
ndarray
) – Error contour for plotting
utils
- splitpy.utils.traceshift(trace, tt)[source]
Function to shift traces in time given travel time
- Parameters:
trace (
Trace
) – Seismogramtt (float) – Time delay (sec)
- Returns:
rtrace – Shifted copy of the input trace
- Return type:
Trace
- splitpy.utils.download_data(client=None, sta=None, start=UTCDateTime(2025, 5, 20, 8, 31, 41, 593370), end=UTCDateTime(2025, 5, 20, 8, 31, 41, 593375), new_sr=0.0, verbose=False, remove_response=False, zcomp='Z')[source]
Function to build a stream object for a seismogram in a given time window by getting data from a client object, either from a local SDS archive or from an FDSN web-service. The function performs sanity checks for the start times, sampling rates and window lengths.
- Parameters:
client (
Client
) – Client objectsta (Dict) – Station metadata from
StDb
data basestart (
UTCDateTime
) – Start time for requestend (
UTCDateTime
) – End time for requestnew_sr (float) – New sampling rate (Hz)
verbose (bool) – Whether or not to print messages to screen during run-time
remove_response (bool) – Remove instrument response from seismogram and resitute to true ground velocity (m/s) using obspy.core.trace.Trace.remove_response()
zcomp (str) – Vertical Component Identifier. Should be a single character. This is different then ‘Z’ only for fully unknown component orientation (i.e., components are 1, 2, 3)
- Returns:
err (bool) – Boolean for error handling (False is associated with success)
trN (
Trace
) – Trace of North component of motiontrE (
Trace
) – Trace of East component of motiontrZ (
Trace
) – Trace of Vertical component of motion