SourceSpec API
SourceSpec has a modular structure. Each module corresponds to a specific function or class of functions.
Main modules
SourceSpec main modules are presented below, following the logical order on which
they’re used within source_spec.py
.
source_spec.py
Earthquake source parameters from inversion of P- or S-wave spectra.
- copyright:
2012 Claudio Satriano <satriano@ipgp.fr>
- 2013-2014 Claudio Satriano <satriano@ipgp.fr>,
Emanuela Matrullo <matrullo@geologie.ens.fr>, Agnes Chounet <chounet@ipgp.fr>
2015-2024 Claudio Satriano <satriano@ipgp.fr>
- license:
CeCILL Free Software License Agreement v2.1 (http://www.cecill.info/licences.en.html)
ssp_parse_arguments
Argument parser for sourcespec.
- copyright:
2021-2024 Claudio Satriano <satriano@ipgp.fr>
- license:
CeCILL Free Software License Agreement v2.1 (http://www.cecill.info/licences.en.html)
ssp_setup
Setup functions for sourcespec.
- copyright:
2012 Claudio Satriano <satriano@ipgp.fr>
- 2013-2014 Claudio Satriano <satriano@ipgp.fr>,
Emanuela Matrullo <matrullo@geologie.ens.fr>, Agnes Chounet <chounet@ipgp.fr>
2015-2024 Claudio Satriano <satriano@ipgp.fr>
- license:
CeCILL Free Software License Agreement v2.1 (http://www.cecill.info/licences.en.html)
- ssp_setup.configure(options, progname, config_overrides=None)[source]
Parse command line arguments and read config file.
- ssp_setup.move_outdir(config)[source]
Move outdir to a new dir named from evid (and optional run_id).
- ssp_setup.setup_logging(config, basename=None, progname='source_spec')[source]
Set up the logging infrastructure.
This function is typically called twice: the first time without basename and a second time with a basename (typically the eventid). When called the second time, the previous logfile is renamed using the given basename.
ssp_read_traces
Read traces in multiple formats of data and metadata.
- copyright:
2012 Claudio Satriano <satriano@ipgp.fr>
- 2013-2014 Claudio Satriano <satriano@ipgp.fr>,
Emanuela Matrullo <matrullo@geologie.ens.fr>
- 2015-2024 Claudio Satriano <satriano@ipgp.fr>,
Sophie Lambotte <sophie.lambotte@unistra.fr>
- license:
CeCILL Free Software License Agreement v2.1 (http://www.cecill.info/licences.en.html)
ssp_process_traces
Trace processing for sourcespec.
- copyright:
2012 Claudio Satriano <satriano@ipgp.fr>
- 2013-2014 Claudio Satriano <satriano@ipgp.fr>,
Emanuela Matrullo <matrullo@geologie.ens.fr>
2015-2024 Claudio Satriano <satriano@ipgp.fr>
- license:
CeCILL Free Software License Agreement v2.1 (http://www.cecill.info/licences.en.html)
ssp_build_spectra
Build spectral objects.
- copyright:
2012 Claudio Satriano <satriano@ipgp.fr>
- 2013-2014 Claudio Satriano <satriano@ipgp.fr>,
Emanuela Matrullo <matrullo@geologie.ens.fr>, Agnes Chounet <chounet@ipgp.fr>
2015-2024 Claudio Satriano <satriano@ipgp.fr>
- license:
CeCILL Free Software License Agreement v2.1 (http://www.cecill.info/licences.en.html)
ssp_plot_traces
Trace plotting routine.
- copyright:
2015-2024 Claudio Satriano <satriano@ipgp.fr>
- license:
CeCILL Free Software License Agreement v2.1 (http://www.cecill.info/licences.en.html)
ssp_inversion
Spectral inversion routines for sourcespec.
- copyright:
2012 Claudio Satriano <satriano@ipgp.fr>
- 2013-2014 Claudio Satriano <satriano@ipgp.fr>,
Emanuela Matrullo <matrullo@geologie.ens.fr>, Agnes Chounet <chounet@ipgp.fr>
2015-2024 Claudio Satriano <satriano@ipgp.fr>
- license:
CeCILL Free Software License Agreement v2.1 (http://www.cecill.info/licences.en.html)
ssp_radiated_energy
Compute radiated energy from spectral integration.
- copyright:
2012 Claudio Satriano <satriano@ipgp.fr>
- 2013-2014 Claudio Satriano <satriano@ipgp.fr>,
Emanuela Matrullo <matrullo@geologie.ens.fr>, Agnes Chounet <chounet@ipgp.fr>
2015-2024 Claudio Satriano <satriano@ipgp.fr>
- license:
CeCILL Free Software License Agreement v2.1 (http://www.cecill.info/licences.en.html)
- ssp_radiated_energy.radiated_energy_and_apparent_stress(config, spec_st, specnoise_st, sspec_output)[source]
Compute radiated energy (in N.m) and apparent stress (in MPa).
- Parameters:
config – Config object
type (sspec_output) –
sourcespec.config.Config
spec_st – Stream of spectra
type –
obspy.core.stream.Stream
specnoise_st – Stream of noise spectra
type –
obspy.core.stream.Stream
sspec_output – Output of spectral inversion
type –
sourcespec.ssp_data_types.SourceSpecOutput
ssp_local_magnitude
Local magnitude calculation for sourcespec.
- copyright:
2012 Claudio Satriano <satriano@ipgp.fr>
- 2013-2014 Claudio Satriano <satriano@ipgp.fr>,
Emanuela Matrullo <matrullo@geologie.ens.fr>
2015-2024 Claudio Satriano <satriano@ipgp.fr>
- license:
CeCILL Free Software License Agreement v2.1 (http://www.cecill.info/licences.en.html)
ssp_summary_statistics
Post processing of station source parameters.
- copyright:
2012-2024 Claudio Satriano <satriano@ipgp.fr>
- license:
CeCILL Free Software License Agreement v2.1 (http://www.cecill.info/licences.en.html)
ssp_output
Output functions for source_spec.
- copyright:
2012 Claudio Satriano <satriano@ipgp.fr>
- 2013-2014 Claudio Satriano <satriano@ipgp.fr>,
Emanuela Matrullo <matrullo@geologie.ens.fr>
2015-2024 Claudio Satriano <satriano@ipgp.fr>
- license:
CeCILL Free Software License Agreement v2.1 (http://www.cecill.info/licences.en.html)
ssp_residuals
Spectral residual routine for sourcespec.
- copyright:
- 2013-2014 Claudio Satriano <satriano@ipgp.fr>,
Agnes Chounet <chounet@ipgp.fr>
2015-2024 Claudio Satriano <satriano@ipgp.fr>
- license:
CeCILL Free Software License Agreement v2.1 (http://www.cecill.info/licences.en.html)
- ssp_residuals.spectral_residuals(config, spec_st, sspec_output)[source]
Compute spectral residuals with respect to an average spectral model. Saves a stream of residuals to disk in HDF5 format.
- Parameters:
config (
Config
) – Configuration objectspec_st (
SpectrumStream
) – Stream of spectrasspec_output (
SourceSpecOutput
) – Output of the source spectral parameter estimation
ssp_plot_spectra
Spectral plotting routine.
- copyright:
2012 Claudio Satriano <satriano@ipgp.fr>
- 2013-2014 Claudio Satriano <satriano@ipgp.fr>,
Emanuela Matrullo <matrullo@geologie.ens.fr>
2015-2024 Claudio Satriano <satriano@ipgp.fr>
- license:
CeCILL Free Software License Agreement v2.1 (http://www.cecill.info/licences.en.html)
ssp_plot_stacked_spectra
Plot stacked spectra, along with summary inverted spectrum.
- copyright:
2023-2024 Claudio Satriano <satriano@ipgp.fr>
- license:
CeCILL Free Software License Agreement v2.1 (http://www.cecill.info/licences.en.html)
ssp_plot_params_stats
Plot parameter statistics.
- copyright:
2022-2024 Claudio Satriano <satriano@ipgp.fr>
- license:
CeCILL Free Software License Agreement v2.1 (http://www.cecill.info/licences.en.html)
ssp_plot_stations
Station plotting routine.
- copyright:
2018-2024 Claudio Satriano <satriano@ipgp.fr>
- license:
CeCILL Free Software License Agreement v2.1 (http://www.cecill.info/licences.en.html)
- ssp_plot_stations.plot_stations(config, sspec_output)[source]
Plot station map, color coded by magnitude or fc.
- Parameters:
config (config.Config) – Configuration object
sspec_output (ssp_data_types.SourceSpecOutput) – SourceSpecOutput object
ssp_html_report
Generate an HTML report for source_spec.
- copyright:
2021-2024 Claudio Satriano <satriano@ipgp.fr>
- license:
CeCILL Free Software License Agreement v2.1 (http://www.cecill.info/licences.en.html)
Other modules
These modules, in alphabetical order, are used by the main modules.
ssp_correction
Spectral station correction calculated from ssp_residuals.
- copyright:
- 2013-2014 Claudio Satriano <satriano@ipgp.fr>,
Agnes Chounet <chounet@ipgp.fr>
2015-2024 Claudio Satriano <satriano@ipgp.fr>
- license:
CeCILL Free Software License Agreement v2.1 (http://www.cecill.info/licences.en.html)
ssp_data_types
Classes for spectral inversion routines.
- copyright:
2017-2024 Claudio Satriano <satriano@ipgp.fr>
- license:
CeCILL Free Software License Agreement v2.1 (http://www.cecill.info/licences.en.html)
- class ssp_data_types.Bounds(config, spec, initial_values)[source]
Bounds for bounded spectral inversion.
- property bounds
Get bounds for minimize() as sequence of (min, max) pairs.
- class ssp_data_types.InitialValues(Mw_0=None, fc_0=None, t_star_0=None)[source]
Initial values for spectral inversion.
- class ssp_data_types.OrderedAttribDict[source]
An ordered dictionary whose values can be accessed as classattributes.
- class ssp_data_types.SourceSpecOutput[source]
The output of SourceSpec.
- error_array(key, filter_outliers=False)[source]
Return an array of errors (two columns) for the given key.
- find_outliers(key, n)[source]
Find outliers using the IQR method.
Q1-n*IQR Q1 median Q3 Q3+n*IQR |-----:-----| o |--------| : |--------| o o |-----:-----| outlier <-----------> outliers IQR
If
n
isNone
, then the above check is skipped.Nan
andinf
values are also marked as outliers.
- percentiles_nobs()[source]
Return a dictionary of number of observations used for computing percentiles.
- reference_summary_parameters()[source]
Return a dictionary of reference summary parameters, each being a SummaryStatistics() object.
- class ssp_data_types.SpectralParameter(param_id, name=None, units=None, value=None, uncertainty=None, lower_uncertainty=None, upper_uncertainty=None, confidence_level=None, format_spec=None)[source]
A spectral parameter measured at one station.
- class ssp_data_types.StationParameters(station_id, instrument_type=None, latitude=None, longitude=None, hypo_dist_in_km=None, epi_dist_in_km=None, azimuth=None)[source]
The parameters describing a given station (e.g., its id and location) and the spectral parameters measured at that station.
Spectral parameters are provided as attributes, using SpectralParameter() objects.
- class ssp_data_types.SummarySpectralParameter(param_id, name=None, units=None, format_spec=None)[source]
A summary spectral parameter comprising one ore more summary statistics.
- class ssp_data_types.SummaryStatistics(stat_type, value=None, uncertainty=None, lower_uncertainty=None, upper_uncertainty=None, confidence_level=None, lower_percentage=None, mid_percentage=None, upper_percentage=None, nobs=None, message=None, format_spec=None)[source]
A summary statistics (e.g., mean, weighted_mean, percentile), along with its uncertainty.
ssp_event
SourceSpec event class and supporting classes.
- copyright:
2023-2024 Claudio Satriano <satriano@ipgp.fr>
- license:
CeCILL Free Software License Agreement v2.1 (http://www.cecill.info/licences.en.html)
- class ssp_event.SSPCoordinate(value=None, units=None, **kwargs)[source]
SourceSpec coordinate class.
Stores a longitude or latitude value. Currently only degrees are supported.
- class ssp_event.SSPDepth(value=None, units=None, **kwargs)[source]
SourceSpec depth class.
- property value_in_km
Return the depth value in km.
- property value_in_m
Return the depth value in m.
- class ssp_event.SSPFocalMechanism(strike=None, dip=None, rake=None, units=None, **kwargs)[source]
SourceSpec focal mechanism class.
Angles are in degrees.
- from_moment_tensor(moment_tensor)[source]
Initialize the focal mechanism object from a moment tensor.
- Parameters:
moment_tensor (
SSPMomentTensor
) – moment tensor object
- class ssp_event.SSPHypocenter(longitude=None, latitude=None, depth=None, origin_time=None, **kwargs)[source]
SourceSpec hypocenter class.
- class ssp_event.SSPMagnitude(value=None, mag_type=None, **kwargs)[source]
SourceSpec magnitude class.
- from_scalar_moment(scalar_moment)[source]
Initialize the magnitude object from a scalar moment.
- Parameters:
scalar_moment (
SSPScalarMoment
) – scalar moment object
- class ssp_event.SSPMomentTensor(units=None, m_rr=None, m_tt=None, m_pp=None, m_rt=None, m_rp=None, m_tp=None, **kwargs)[source]
SourceSpec moment tensor class.
- class ssp_event.SSPScalarMoment(value=None, units=None, **kwargs)[source]
SourceSpec scalar moment class.
- from_moment_tensor(moment_tensor)[source]
Initialize the scalar moment object from a moment tensor.
- Parameters:
moment_tensor (
SSPMomentTensor
) – moment tensor object
ssp_grid_sampling
A class for sampling a parameter space over a grid.
Sampling can be performed by several approaches. The class provides optimal solutions, uncertainties and plotting methods.
- copyright:
2022-2024 Claudio Satriano <satriano@ipgp.fr>
- license:
CeCILL Free Software License Agreement v2.1 (http://www.cecill.info/licences.en.html)
- class ssp_grid_sampling.GridSampling(misfit_func, bounds, nsteps, sampling_mode, params_name, params_unit)[source]
A class for sampling a parameter space over a grid.
Sampling can be performed by several approaches. The class provides optimal solutions, uncertainties and plotting methods.
- property conditional_misfit
Compute conditional misfit along each dimension.
Conditional misfit is computed by fixing the other parameters to their optimal value.
- property conditional_peak_widths
Find width of conditional misfit around its minimum.
- property min_idx
Find the index of the minimum of the misfit function.
- property params_err
Return optimal parameters uncertainties.
- property params_opt
Return optimal parameters.
- property values
Return a meshgrid of parameter values.
- property values_1d
Extract a 1D array of parameter values along one dimension.
ssp_pick
SourceSpec pick class.
- copyright:
2023-2024 Claudio Satriano <satriano@ipgp.fr>
- license:
CeCILL Free Software License Agreement v2.1 (http://www.cecill.info/licences.en.html)
ssp_spectral_model
Spectral model and objective function.
- copyright:
2012 Claudio Satriano <satriano@ipgp.fr>
- 2013-2014 Claudio Satriano <satriano@ipgp.fr>,
Emanuela Matrullo <matrullo@geologie.ens.fr>, Agnes Chounet <chounet@ipgp.fr>
2015-2024 Claudio Satriano <satriano@ipgp.fr>
- license:
CeCILL Free Software License Agreement v2.1 (http://www.cecill.info/licences.en.html)
- ssp_spectral_model.objective_func(xdata, ydata, weight)[source]
Objective function generator for bounded inversion.
- ssp_spectral_model.spectral_model(freq, Mw, fc, t_star, alpha=1.0)[source]
Spectral model.
\[Y_{data} = M_w + \frac{2}{3} \left[ - \log_{10} \left( 1+\left(\frac{f}{f_c}\right)^2 \right) - \pi \, f t^* \log_{10} e \right]\]see Theoretical Background for a detailed derivation of this model.
ssp_qml_output
QuakeML output for source_spec.
- copyright:
2016-2024 Claudio Satriano <satriano@ipgp.fr>
- license:
CeCILL Free Software License Agreement, Version 2.1 (http://www.cecill.info/index.en.html)
ssp_radiation_pattern
Compute radiation pattern.
- copyright:
2021-2024 Claudio Satriano <satriano@ipgp.fr>
- license:
CeCILL Free Software License Agreement v2.1 (http://www.cecill.info/licences.en.html)
- ssp_radiation_pattern.get_radiation_pattern_coefficient(stats, config)[source]
Get radiation pattern coefficient.
ssp_read_sac_header
Read metadata from SAC file headers.
- copyright:
2023-2024 Claudio Satriano <satriano@ipgp.fr>
- license:
CeCILL Free Software License Agreement v2.1 (http://www.cecill.info/licences.en.html)
- ssp_read_sac_header.compute_sensitivity_from_SAC(trace, config)[source]
Compute sensitivity from SAC header fields.
ssp_read_station_metadata
Read station metadata in StationXML, dataless SEED, SEED RESP, PAZ (SAC polezero format).
- copyright:
2012-2024 Claudio Satriano <satriano@ipgp.fr>
- license:
CeCILL Free Software License Agreement v2.1 (http://www.cecill.info/licences.en.html)
ssp_sqlite_output
SQLite output for source_spec.
- copyright:
2013-2024 Claudio Satriano <satriano@ipgp.fr>
- license:
CeCILL Free Software License Agreement v2.1 (http://www.cecill.info/licences.en.html)
- ssp_sqlite_output.write_sqlite(config, sspec_output)[source]
Write SSP output to SQLite database.
- Parameters:
config (config.Config) – SSP configuration object
sspec_output (ssp_data_types.SourceSpecOutput) – SSP output object
ssp_update_db
Update an existing SourceSpec database from a previous version.
- copyright:
2013-2024 Claudio Satriano <satriano@ipgp.fr>
- license:
CeCILL Free Software License Agreement v2.1 (http://www.cecill.info/licences.en.html)
ssp_util
Utility functions for sourcespec.
- copyright:
2012-2024 Claudio Satriano <satriano@ipgp.fr>
- license:
CeCILL Free Software License Agreement v2.1 (http://www.cecill.info/licences.en.html)
- class ssp_util.MediumProperties(lon, lat, depth_in_km, config)[source]
Class to retrieve medium properties from config.
- Parameters:
- get(mproperty, where)[source]
Get medium property (P- or S-wave velocity, density) at a given point from NLL grid, config or taup model.
- get_from_config_param_source(mproperty)[source]
Get medium property at the source from config parameter.
- get_from_config_param_station(mproperty)[source]
Get medium property at the station from config parameter.
- get_from_taup(mproperty)[source]
Get medium property (P- or S-wave velocity, density) at a given depth from taup model.
- ssp_util.geom_spread_boatwright(hypo_dist_in_km, cutoff_dist_in_km, freqs)[source]
” Geometrical spreading coefficient from Boatwright et al. (2002), eq. 8.
Except that we take the square root of eq. 8, since we correct amplitude and not energy.
- Parameters:
hypo_dist_in_km (float) – Hypocentral distance (km).
cutoff_dist_in_km (float) – Cutoff distance (km).
freqs (numpy.ndarray) – Frequencies (Hz).
- Returns:
Geometrical spreading correction (in m)
- Return type:
- ssp_util.geom_spread_r_power_n(hypo_dist_in_km, exponent)[source]
rⁿ geometrical spreading coefficient.
- ssp_util.geom_spread_teleseismic(angular_distance, source_depth_in_km, station_depth_in_km, phase)[source]
Calculate geometrical spreading coefficient for teleseismic body waves.
Implements eq (4) in Okal (1992) for a spherically symmetric Earth. This equations is derived from the conservation of the kinetic energy flux along a ray tube between the source and the receiver.
- ssp_util.quality_factor(travel_time_in_s, t_star_in_s)[source]
Compute quality factor from travel time and t_star.
- ssp_util.remove_instr_response(trace, pre_filt=(0.5, 0.6, 40.0, 45.0))[source]
Remove instrument response from a trace.
Trace is converted to the sensor units (m for a displacement sensor, m/s for a short period or broadband velocity sensor, m/s**2 for a strong motion sensor).
- ssp_util.select_trace(stream, traceid, instrtype)[source]
Select trace from stream using traceid and instrument type.
- ssp_util.smooth(signal, window_len=11, window='hanning')[source]
Smooth the signal using a window with requested size.
- ssp_util.source_radius(fc_in_hz, vs_in_m_per_s, k_coeff=0.3724)[source]
Compute source radius in meters.
Madariaga (2009), doi:10.1007/978-1-4419-7695-6_22, eq. 31 Kaneko and Shearer (2014), doi:10.1093/gji/ggu030, eq. 2, 15, 16
- ssp_util.spec_minmax(amp, freq, amp_minmax=None, freq_minmax=None)[source]
Get minimum and maximum values of spectral amplitude and frequency.
- ssp_util.static_stress_drop(Mo_in_N_m, ra_in_m)[source]
Compute static stress drop in MPa.
Madariaga (2009), doi:10.1007/978-1-4419-7695-6_22, eq. 27
ssp_wave_arrival
Arrival time calculation for sourcespec.
- copyright:
2012-2024 Claudio Satriano <satriano@ipgp.fr>
- license:
CeCILL Free Software License Agreement v2.1 (http://www.cecill.info/licences.en.html)
config
Config class for sourcespec.
- copyright:
2013-2024 Claudio Satriano <satriano@ipgp.fr>
- license:
CeCILL Free Software License Agreement v2.1 (http://www.cecill.info/licences.en.html)
kdtree
Grid importance sampling using a k-d tree.
- copyright:
2022-2024 Claudio Satriano <satriano@ipgp.fr>
- license:
CeCILL Free Software License Agreement v2.1 (http://www.cecill.info/licences.en.html)
- class kdtree.KDTCell(extent, calc_pdf, min_cell_prob=0, ndiv=None, maxdiv=None)[source]
A cell of a k-d tree.
savefig
Save Matplotlib figure. Optimize PNG format using PIL.
- copyright:
2022-2024 Claudio Satriano <satriano@ipgp.fr>
- license:
CeCILL Free Software License Agreement v2.1 (http://www.cecill.info/licences.en.html)
spectrum
Define Spectrum and SpectrumStream classes, similar to ObsPy’s Trace and Stream.
Provides the high-level function read_spectra() to read SpectrumStream objects from HDF5 or TEXT files.
- copyright:
2012-2024 Claudio Satriano <satriano@ipgp.fr>
- license:
CeCILL Free Software License Agreement v2.1 (http://www.cecill.info/licences.en.html)
- class spectrum.AttributeDict(*args, **kwargs)[source]
A dictionary that allows attribute-style access.
- class spectrum.Spectrum(obspy_trace=None)[source]
A class to handle amplitude spectra.
- Parameters:
obspy_trace – An ObsPy Trace object to compute the spectrum from.
- property data
Return the array containing the amplitude spectrum.
- property data_logspaced
Return the array containing the amplitude spectrum in logspaced frequencies.
- property data_mag
Return the array containing the amplitude spectrum in mangitude units.
- property data_mag_logspaced
Return the array containing the amplitude spectrum in logspaced frequencies in magnitude units.
- property freq
Return the frequency axis of the spectrum.
- property freq_logspaced
Return the logspaced frequency axis of the spectrum.
- property id
Return the id of the spectrum.
- slice(fmin, fmax, nearest_sample=True, pad=False, fill_value=None)[source]
Slice the spectrum between fmin and fmax.
- Parameters:
fmin – Minimum frequency.
fmax – Maximum frequency.
nearest_sample – If True, the slice will include the nearest frequency to fmin and fmax.
pad – If True, the slice will be padded with the value of fill_value until fmin and fmax are included.
fill_value – The value to use for padding.
- Note:
Only the linear spaced frequencies, data and data_mag are sliced. If the original spectrum contains logspaced frequencies, data, and data_mag, those are not preserved in the sliced spectrum.
- class spectrum.SpectrumStream(iterable=(), /)[source]
A class to handle a collection of amplitude spectra.
clipping_detection
Check trace for clipping using kernel density estimation of the trace= amplitude values.
- Two methods are available:
clipping_score()
: compute a trace clipping score based on the shape of the kernel density estimation.clipping_peaks()
: check if trace is clipped, based on the number of peaks in the kernel density estimation;
- copyright:
- 2023-2024 Claudio Satriano <satriano@ipgp.fr>,
Kris Vanneste <kris.vanneste@oma.be>
- license:
CeCILL Free Software License Agreement v2.1 (http://www.cecill.info/licences.en.html)
- clipping_detection.clipping_peaks(trace, sensitivity=3, clipping_percentile=10, debug=False)[source]
Check if a trace is clipped, based on the number of peaks in the kernel density estimation of the trace amplitude values.
The algorithm is based on the following steps:
The trace is demeaned.
A kernel density estimation is computed on the trace amplitude values.
The kernel density estimation is weighted by the distance from the zero mean amplitude value, using a parabolic function, between 1 and 5.
Peaks are detected in the weighted kernel density estimation. The sensitivity of the peak detection algorithm is controlled by the
sensitivity
parameter, based on which a minimum prominence threshold is set.The trace is considered clipped if there is at least one peak in the amplitude range corresponding to the
clipping_percentile
parameter.
- Parameters:
trace (
obspy.core.trace.Trace
) – Trace to check.sensitivity (
int
) – Sensitivity level, from 1 (least sensitive) to 5 (most sensitive). (default: 3). The sensitivity level controls the minimum prominence threshold used for peak detection in the kernel density estimation. See thescipy.signal.find_peaks()
documentation for more details.clipping_percentile (
float
,between 0
and100
) – Percentile of trace amplitude range (expressed as percentage) to check for clipping. Default is 10, which means that the 10% highest and lowest values of the trace amplitude will be checked for clipping. A value of 0 means that no clipping check will be performed.debug (
bool
) – If True, plot trace, samples histogram and kernel density.
- Returns:
trace_clipped (
bool
) – True if trace is clipped, False otherwise.properties (
dict
) – Dictionary with the properties used for the clipping detection. The dictionary contains the following keys:npeaks (
int
): number of peaks found in the kernel densitynpeaks_clipped (
int
): number of peaks found in the kernel density that are considered clippedpeaks (
numpy.ndarray
): list of peaks found in the kernel densityprominences (
numpy.ndarray
): list of prominences of the peaks found in the kernel density
- clipping_detection.compute_clipping_score(trace, remove_baseline=False, debug=False)[source]
Compute a trace clipping score based on the shape of the kernel density estimation of the trace amplitude values.
The algorithm is based on the following steps:
The trace is detrended and demeaned. Optionally, the trace baseline can be removed.
A kernel density estimation is computed on the trace amplitude values.
Two weighted kernel density functions are computed:
a full weighted kernel density, where the kernel density is weighted by the distance from the zero mean amplitude value, using a 8th order power function between 1 and 100.
a weighted kernel density without the central peak, where the kernel density is weighted by the distance from the zero mean amplitude value, using a 8th order power function between 0 and 100.
In both cases, the weight gives more importance to samples far from the zero mean value. In the second case, the central peak is ignored.
The score, ranging from 0 to 100, is the sum of the squared weighted kernel density without the central peak, normalized by the sum of the squared full weighted kernel density. The score is 0 if there is no additional peak beyond the central peak.
- Parameters:
trace (
Trace
) – Trace to check for clipping.remove_baseline (
bool
, optional) – If set, remove the trace baseline before computing the score.debug (
bool
, optional) – If set, plot the trace, the samples histogram, the kernel density (unweighted and weighted), the kernel baseline model, and the misfit. Default is False.
- Returns:
clipping_score – Clipping score, in percentage.
- Return type:
Note
Distorted traces (e.g., signals with strong baselines) can also get a high clipping score. To avoid this, the trace baseline can be removed.
A debug mode is available to plot the trace, the samples histogram, the kernel density (unweighted and weighted), the kernel baseline model, and the misfit.
plot_sourcepars
plot_sourcepars.py
1D or 2D plot of source parameters from a sqlite parameter file.
- copyright:
2023-2024 Claudio Satriano <satriano@ipgp.fr>
- license:
CeCILL Free Software License Agreement v2.1 (http://www.cecill.info/licences.en.html)
- class plot_sourcepars.Annot(xdata, ydata, labels, yformat)[source]
Annotate the plot with the evid, Mw and the value of the parameter.
- class plot_sourcepars.Params(args)[source]
Class to handle the parameters from a sqlite file.
- filter(stamin=None, stamax=None, magmin=None, magmax=None, ssdmin=None, ssdmax=None, sigmaamin=None, sigmaamax=None)[source]
Filter the parameters based on one or more conditions.
- plot_Er_mw(hist=False, fit=False, nbins=None)[source]
Plot the logarithm of the radiated energy vs the moment magnitude.
- plot_fc_mw(hist=False, fit=False, slope=False, nbins=None, wave_type='S')[source]
Plot the logarithm of the corner frequency vs the moment magnitude.
- plot_sourcepars.apparent_stress_curve_Er_mw(sigma_a, mu, mw)[source]
Constant apparent stress curve in Er vs Mw.
Madariaga (2009), doi:10.1007/978-1-4419-7695-6_22, eq. 33., page 374.
- plot_sourcepars.mag_to_moment(mag, b=0.5)[source]
Convert magnitude to moment.
The parameter b is used to change the slope of the fc-Mw stress drop curve.
The standard value of b is 0.5, which corresponds to a self-similar model.
- plot_sourcepars.query_event_params_into_numpy(cursor, param, param_type, query_condition)[source]
Query a parameter from the Events table and return it as a numpy array.
- plot_sourcepars.stress_drop_curve_Er_mw(delta_sigma, mu, mw)[source]
Constant stress drop curve in Er vs Mw.
Madariaga (2009), doi:10.1007/978-1-4419-7695-6_22, eq. 33., page 374.
- plot_sourcepars.stress_drop_curve_fc_mw(delta_sigma, vel, mw, k=0.3724, b=-0.5)[source]
Constant stress drop curve in fc vs Mw.
Obtained by combining the equation for stress drop:
delta_sigma = 7/16 * Mo / a^3
with the equation for source radius:
a = k * vel * (delta_sigma / fc)
where k is a coefficient discussed in Kaneko and Shearer (2014).
For the Brune source model, k=0.3724.
- Parameters:
- Returns:
fc – Corner frequency in Hz.
- Return type: