data:image/s3,"s3://crabby-images/33acc/33acc853121be92ac495fe322dc56154775bf781" alt="_images/OrientPy_logo.png"
Classes
orientpy
defines the following base classes:
Orient
- class orientpy.classes.Orient(sta)
A Orient object contains Class attributes that associate station information with a single event (i.e., earthquake) metadata.
Note
The object is initialized with the
sta
field only, and other attributes are added to the object as the analysis proceeds.- Parameters:
sta (object) – Object containing station information - from
stdb
database.meta (
Meta
) – Object of metadata information for single event (initially set to None)data (
Stream
) – Stream object containing the three-component seismograms
- add_event(event, gacmin=None, gacmax=None, depmax=1000.0, returned=False)
Adds event metadata to Orient object.
- Parameters:
event (
event
) – Event XML objectgacmin (float) – Minimum great-arc circle distance to consider (deg)
gacmax (float) – Maximum great-arc circle distance to consider (deg)
depmax (float) – Maximum event depth to consider (km)
returned (bool) – Whether or not to return the
accept
attribute
- Returns:
accept – Whether or not the object is accepted for further analysis
- Return type:
bool
- download_data(client, stdata=[], ndval=nan, new_sr=None, t1=None, t2=None, returned=False, verbose=False)
Downloads seismograms based on a time interval.
- Parameters:
client (
Client
) – Client objectstdata (List) – Station list
ndval (float) – Fill in value for missing data
new_sr (float) – New sampling rate (Hz)
t1 (float) – Start time to consider (sec). Can be float or
UTCDateTime
object.t2 (float) – End time to consider (sec). Can be float or
UTCDateTime
object.returned (bool) – Whether or not to return the
accept
attribute
- Returns:
accept – Whether or not the object is accepted for further analysis
- Return type:
bool
- data
Stream containing
Trace
objects- Type:
Stream
Class Meta
- class orientpy.classes.Meta(sta, event, gacmin=5.0, gacmax=30.0, depmax=1000.0)
A Meta object contains attributes associated with the metadata for a single earthquake event.
- Parameters:
time (
UTCDateTime
) – Origin time of earthquakedep (float) – Depth of hypocenter (km)
lon (float) – Longitude coordinate of epicenter
lat (float) – Latitude coordinate of epicenter
mag (float) – Magnitude of earthquake
gac (float) – Great arc circle between station and epicenter (degrees)
epi_dist (float) – Epicentral distance between station and epicenter (km)
baz (float) – Back-azimuth - pointing to earthquake from station (degrees)
az (float) – Azimuth - pointing to station from earthquake (degrees)
BNG
- class orientpy.classes.BNG(sta)
A BNG object inherits from
Orient
. This object contains a method to estimate station orientation based on P-wave particle motion.Notes
This class is designed after the method developed by Braunmiller, Nabelek and Ghods (2020) [1]. It is slightly more flexible, however, in that it can handle either regional or teleseismic P-wave data and provides additional quality-control mesures (e.g., SNR, 1-R/Z).
References
- Parameters:
sta (object) – Object containing station information - from
stdb
database.meta (
Meta
) – Object of metadata information for single event (initially set to None)data (
Stream
) – Stream object containing the three-component seismograms
- calc(dphi, dts, tt, bp=None, showplot=False)
Method to estimate azimuth of component ?H1 (or ?HN). This method minimizes the energy (RMS) of the transverse component of P-wave data within some bandwidth.
- Parameters:
dphi (float) – Azimuth increment for search (deg)
dts (float) – Length of time window on either side of predicted P-wave arrival time (sec)
tt (list) – List of two floats containing the time picks relative to P-wave time, within which to perform the estimation of station orientation (sec)
bp (list) – List of two floats containing the low- and high-frequency corners of a bandpass filter (Hz)
showplot (bool) – Whether or not to plot waveforms.
- meta.phi
Azimuth of H1 (or HN) component (deg)
- Type:
float
- meta.cc
Cross-correlation coefficient between vertical and radial component
- Type:
float
- meta.snr
Signal-to-noise ratio of P-wave measured on the vertical seismogram
- Type:
float
- meta.TR
Measure of the transverse to radial ratio. In reality this is 1 - T/R
- Type:
float
- meta.RZ
Measure of the radial to vertical ratio. In reality this is 1 - R/Z
- Type:
float
DL
- class orientpy.classes.DL(sta)
A DL object inherits from
Orient
. This object contains a method to estimate station orientation based on Rayleigh-wave particle motion.Notes
This algorithm is heavily based on the code DLOPy by Doran and Laske (2017) [2]. The original code can be found here: https://igppweb.ucsd.edu/~adoran/DLOPy.html
There are also a couple of versions available on GitHub:
References
[2] Doran, A. K., and G. Laske (2017). Ocean-bottom seismometer instrument orientations via automated Rayleigh-wave arrival-angle measurements, Bulletin of the Seismological Society of America, 107(2), doi:10.1785/0120160165.
- Parameters:
sta (object) – Object containing station information - from
stdb
database.meta (
Meta
) – Object of metadata information for single event (initially set to None)data (
Stream
) – Stream object containing the three-component seismograms
- calc(showplot=False)
Method to estimate azimuth of component ?H1 (or ?HN). This method maximizes the normalized covariance between the Hilbert-transformed vertical component and the radial component of Rayleigh-wave data measured at multiple periods. This is done for the shortest (R1) and longest (R2) orbits of fundamental-mode Rayleigh waves. Window selection is done based on average group velocity extracted from a global model of Rayleigh-wave dispersion.
- Parameters:
showplot (bool) – Whether or not to plot waveforms.
- meta.R1phi
List of azimuth of H1 (or HN) component (deg) based on R1, measured at different periods
- Type:
list
- meta.R1cc
Cross-correlation coefficient between hilbert-transformed vertical and radial component for R1, measured at different periods
- Type:
float
- meta.R2phi
List of azimuth of H1 (or HN) component (deg) based on R2, measured at different periods
- Type:
list
- meta.R2cc
Cross-correlation coefficient between hilbert-transformed vertical and radial component for R2, measured at different periods
- Type:
float
Modules
utils
- orientpy.utils.catclean(cat)
This function looks for repeat events in a catalogue of earthquakes.
ADRIAN. K. DORAN and GABI LASKE, DLOPy VERSION 1.0, RELEASED APRIL 2017
- Parameters:
cat (
Catalog
) – Catalogue XML object, container for Event objects.- Returns:
reps – Array of ID of repeat events in catalogue.
- Return type:
ndarray
- orientpy.utils.checklen(st, hrs)
Function to check if there is enough downloaded data to run program
ADRIAN. K. DORAN and GABI LASKE, DLOPy VERSION 1.0, RELEASED APRIL 2017
- Parameters:
st (
Stream
) – Stream containing waveformshrs (float) – Amount of time (hours) that should be available for analysis
- Returns:
Whether or not the check is successful
- Return type:
boolean
- orientpy.utils.resiz(x1, x2, x3)
Function to resize arrays to all identical shapes
ADRIAN. K. DORAN and GABI LASKE, DLOPy VERSION 1.0, RELEASED APRIL 2017
- Parameters:
x* (
ndarray
) – Array for which to check length (here it’s trace.data)- Returns:
x* – Array trimmed to shortest segment
- Return type:
ndarray
- orientpy.utils.getf(freq, A)
Function to extract frequency in array.
ADRIAN. K. DORAN and GABI LASKE, DLOPy VERSION 1.0, RELEASED APRIL 2017
- Parameters:
freq (float) – Frequency of interest (Hz)
A (
ndarray
) – Array of frequencies obtained from global dispersion maps
- Returns:
v – Nearest frequency in A
- Return type:
float
- orientpy.utils.nv(x, v)
Function to find nearest value in an np.array.
ADRIAN. K. DORAN and GABI LASKE, DLOPy VERSION 1.0, RELEASED APRIL 2017
- Parameters:
x (
ndarray
) – Input arrayv (float) – Value to compare within array
- Returns:
x[idx] – Nearest value to v in x
- Return type:
float
- orientpy.utils.rms(x)
Function to calculate root-mean-square of array
ADRIAN. K. DORAN and GABI LASKE, DLOPy VERSION 1.0, RELEASED APRIL 2017
- Parameters:
x (
ndarray
) – Input array- Returns:
rms – Root-Mean-Square value of x
- Return type:
float
- orientpy.utils.mad(x)
Function to calculate Median Absolute Deviation
ADRIAN. K. DORAN and GABI LASKE, DLOPy VERSION 1.0, RELEASED APRIL 2017
- Parameters:
x (
ndarray
) – Input array- Returns:
mad – Median Absolute deviation of x
- Return type:
float
- orientpy.utils.outlier(x, lim=5.0)
Function to remove outliers based on MAD threshold
ADRIAN. K. DORAN and GABI LASKE, DLOPy VERSION 1.0, RELEASED APRIL 2017
- Parameters:
x (
ndarray
) – Input array- Returns:
x – Shortened array where MAD is lower than lim
- Return type:
ndarray
- orientpy.utils.boot(x, bootnum)
Function to calculate directional mean value from bootstrap
ADRIAN. K. DORAN and GABI LASKE, DLOPy VERSION 1.0, RELEASED APRIL 2017
- Parameters:
x (
ndarray
) – Input arraybootnum (int) – Number of bootstrap estimates to produce
- Returns:
m – Array of directional mean values
- Return type:
ndarray
- orientpy.utils.centerat(phi, m=0.0)
Function to re-center azimuth array to the mean value
ADRIAN. K. DORAN and GABI LASKE, DLOPy VERSION 1.0, RELEASED APRIL 2017
- Parameters:
phi (
ndarray
) – Input arraym (float) – Value around which to center the azimuth
- Returns:
phinew – Re-centered azimuth array
- Return type:
ndarray
- orientpy.utils.pathvels(lat1, lon1, lat2, lon2, map10, map15, map20, map25, map30, map35, map40)
Overall function to get path-averaged group velocity.
ADRIAN. K. DORAN and GABI LASKE, DLOPy VERSION 1.0, RELEASED APRIL 2017
- Parameters:
lat1 (float) – Latitude of origin point (deg)
lon1 (float) – Longitude of origin point (deg)
lat2 (float) – Latitude of end point (deg)
lon2 (float) – Longitude of end point (deg)
map* (
ndarray
) – maps of Rayleigh-wave dispersion at various frequencies
- Returns:
R1 (
ndarray
) – R1 velocity pathR2 (
ndarray
) – R2 velocity path
- orientpy.utils.DLcalc(stream, Rf, LPF, HPF, epi, baz, A, winlen=10.0, ptype=0)
DORAN-LASKE calculation for one freq, one orbit of surface wave
ADRIAN. K. DORAN and GABI LASKE, DLOPy VERSION 1.0, RELEASED APRIL 2017
- Parameters:
stream (float) – Latitude of origin point (deg)
lon1 (float) – Longitude of origin point (deg)
lat2 (float) – Latitude of end point (deg)
lon2 (float) – Longitude of end point (deg)
map* (
ndarray
) – maps of Rayleigh-wave dispersion at various frequencies
- Returns:
R1 (
ndarray
) – R1 velocity pathR2 (
ndarray
) – R2 velocity path
- orientpy.utils.estimate(phi, ind)
Function to estimate final azimuth from
ADRIAN. K. DORAN and GABI LASKE, DLOPy VERSION 1.0, RELEASED APRIL 2017
- Parameters:
phi (
ndarray
) – Input array of estimated azimuthind (list) – List of index values that satisfy some QC condition for phi
- Returns:
m (float) – Mean value of robust, bootstrapped estimates
s (float) – 2 Standard deviations of robust, bootstrapped estimates
plotting
- orientpy.plotting.density_estimate(values, x, alpha)
This function estimates the KDE, mode of distribution, and 95% confidence intervals for the
values
variables evaluated at locationsx
- Parameters:
values (
ndarray
) – Array of values for which to estimate KDEx (
ndarray
) – Array of locations at which to evaluate the KDEalpha (float) – Significance level (alpha=0.05 means 95% confidence)
- Returns:
kde (
gaussian_kde
) – KDE object for datavalues
evaluated atx
x_max (float) – Value at which the KDE is maximum
CI_min (float) – Lower bound of confidence interval
CI_max (float) – Upper bound of confidence interval
- orientpy.plotting.plot_bng_waveforms(bng, stream, dts, tt)
This function plots the original and rotated waveforms following the BNG processing for quality control.
- Parameters:
bng (
BNG
) – BNG object for a single eventstream (
Stream
) – Stream of un-rotated and rotated waveformsdts (float) – Length of time window over which to plot waveforms
tt (list) – Zoomed in time window over which processing is done
- Returns:
plt – Handle to final plot
- Return type:
pyplot
- orientpy.plotting.plot_bng_conditions(stkey, snr, cc, TR, RZ, ind)
This function plots all parameters and the threshold condition for contributing to the final estimate.
- Parameters:
stkey (str) – Station key
snr (
ndarray
) – Array of signal-to-noise ratio values for each earthquakecc (
ndarray
) – Array of cross-correlation values between rotated radial and vertical componentsTR (
ndarray
) – Array of transverse-to-radial component ratios for rotated componentsRZ (
ndarray
) – Array of radial-to-vertical component ratios for rotated componentsind (
ndarray
) – Array of boolean (index) values where conditions on previous parameters are examined
- Returns:
plt – Handle to final plot
- Return type:
pyplot
- orientpy.plotting.plot_bng_results(stkey, phi, snr, cc, TR, RZ, baz, mag, ind, val, err, alpha=0.05)
This function plots the results of all BNG estimates with final estimates from those that pass the conditions.
- Parameters:
stkey (str) – Station key
phi (
ndarray
) – Array of estimated azimuth for each earthquakesnr (
ndarray
) – Array of signal-to-noise ratio values for each earthquakecc (
ndarray
) – Array of cross-correlation values between rotated radial and vertical componentsTR (
ndarray
) – Array of transverse-to-radial component ratios for rotated componentsRZ (
ndarray
) – Array of radial-to-vertical component ratios for rotated componentsbaz (
ndarray
) – Array of back-azimuth values corresponding to each earthquakemag (
ndarray
) – Array of magnitude values corresponding to each earthquakeind (
ndarray
) – Array of boolean values where conditions on previous parameters are examinedval (float) – Final estimated azimuth for station orientation
err (float) – Final error estimate on azimuth
alpha (float) – Significance level (alpha=0.05 means 95% confidence)
- Returns:
plt – Handle to final plot
- Return type:
pyplot
- orientpy.plotting.plot_dl_results(stkey, R1phi, R1cc, R2phi, R2cc, ind, val, err, phi, cc, cc0, alpha=0.05)
This function plots the results of all DL estimates with final estimates from those that pass the conditions.
- Parameters:
stkey (str) – Station key
R1phi (
ndarray
) – Array of azimuth values from each estimate for direct Rayleigh-wave passR1cc (
ndarray
) – Array of cross-correlation values between rotated radial and vertical components for direct Rayleigh-wave passR2phi (
ndarray
) – Array of azimuth values from each estimate for complementary Rayleigh-wave passR2cc (
ndarray
) – Array of cross-correlation values between rotated radial and vertical components for complementary Rayleigh-wave passind (
ndarray
) – Array of boolean values where conditions on previous parameters are examinedval (float) – Final estimated azimuth for station orientation
err (float) – Final error estimate on azimuth
phi (
ndarray
) – All azimuth estimates from both R1 (direct pass) and R2 (complementary pass)cc (
ndarray
) – All cross-correlation estimates from both R1 (direct pass) and R2 (complementary pass)cc0 (float) – Threshold cross-correlation value
alpha (float) – Significance level (alpha=0.05 means 95% confidence)
- Returns:
plt – Handle to final plot
- Return type:
pyplot