core.match_filter.party

Functions for network matched-filter detection of seismic data.

Designed to cross-correlate templates generated by template_gen function with data and output the detections.

copyright:

EQcorrscan developers.

license:

GNU Lesser General Public License, Version 3 (https://www.gnu.org/copyleft/lesser.html)


class eqcorrscan.core.match_filter.party.Party(families=None)[source]

Container for multiple Family objects.

Methods

copy()

Returns a copy of the Party.

filter([dates, min_dets])

Return a new Party filtered according to conditions.

plot([plot_grouped, dates, min_dets, rate])

Plot the cumulative detections in time.

decluster(trig_int[, timing, metric, ...])

De-cluster a Party of detections by enforcing a detection separation.

get_catalog()

Get an obspy catalog object from the party.

lag_calc(stream, pre_processed[, shift_len, ...])

Compute picks based on cross-correlation alignment.

min_chans(min_chans)

Remove detections using min_chans or fewer channels.

read([filename, read_detection_catalog, ...])

Read a Party from a file.

rethreshold(new_threshold[, ...])

Remove detections from the Party that are below a new threshold.

select(template_name)

Select a specific family from the party.

sort()

Sort the families by template name.

write(filename[, format, ...])

Write Family out, select output format.

__init__(families=None)[source]

Instantiate the Party object.

copy()[source]

Returns a copy of the Party.

Returns:

Copy of party

Example

>>> from eqcorrscan import Family, Template
>>> party = Party(families=[Family(template=Template(name='a'))])
>>> party_b = party.copy()
>>> party == party_b
True
filter(dates=None, min_dets=1)[source]

Return a new Party filtered according to conditions.

Return a new Party with only detections within a date range and only families with a minimum number of detections.

Parameters:
  • dates (list of obspy.core.UTCDateTime objects) – A start and end date for the new Party

  • min_dets (int) – Minimum number of detections per family

Example

>>> from obspy import UTCDateTime
>>> Party().read().filter(dates=[UTCDateTime(2016, 1, 1),
...                              UTCDateTime(2017, 1, 1)],
...                       min_dets=30) 
plot(plot_grouped=False, dates=None, min_dets=1, rate=False, **kwargs)[source]

Plot the cumulative detections in time.

Parameters:
  • plot_grouped (bool) – Whether to plot all families together (plot_grouped=True), or each as a separate line.

  • dates (list) – list of obspy.core.UTCDateTime objects bounding the plot. The first should be the start date, the last the end date.

  • min_dets (int) – Plot only families with this number of detections or more.

  • rate (bool) – Whether or not to plot the daily rate of detection as opposed to cumulative number. Only works with plot_grouped=True.

  • **kwargs – Any other arguments accepted by eqcorrscan.utils.plotting.cumulative_detections()

Examples

Plot cumulative detections for all templates individually:

>>> Party().read().plot()  

Plot cumulative detections for all templates grouped together:

>>> Party().read().plot(plot_grouped=True) 

Plot the rate of detection for all templates grouped together:

>>> Party().read().plot(plot_grouped=True, rate=True) 

Plot cumulative detections for all templates with more than five detections between June 1st, 2012 and July 31st, 2012:

>>> from obspy import UTCDateTime
>>> Party().read().plot(dates=[UTCDateTime(2012, 6, 1),
...                            UTCDateTime(2012, 7, 31)],
...                     min_dets=5) 
decluster(trig_int, timing='detect', metric='avg_cor', hypocentral_separation=None, min_chans=0, absolute_values=False, num_threads=None)[source]

De-cluster a Party of detections by enforcing a detection separation.

De-clustering occurs between events detected by different (or the same) templates. If multiple detections occur within trig_int then the preferred detection will be determined by the metric argument. This can be either the average single-station correlation coefficient which is calculated as Detection.detect_val / Detection.no_chans, or the raw cross channel correlation sum which is simply Detection.detect_val.

Parameters:
  • trig_int (float) – Minimum detection separation in seconds.

  • metric (str) – What metric to sort peaks by. Either ‘avg_cor’ which takes the single station average correlation, ‘cor_sum’ which takes the total correlation sum across all channels, or ‘thresh_exc’ which takes the factor by how much the detection exceeded the input threshold.

  • timing (str) – Either ‘detect’ or ‘origin’ to decluster based on either the detection time or the origin time.

  • hypocentral_separation (float) – Maximum inter-event separation in km to decluster events within. If an event happens within this distance of another event, and within the trig_int defined time, then they will be declustered.

  • min_chans (int) – Minimum number of channels for a detection to be retained. If a detection was made on very few channels, then the MAD detection statistics become much less robust.

  • absolute_values (bool) – Use the absolute value of the metric to choose the preferred detection.

  • num_threads (int) – Number of threads to use for internal c-funcs if available. Only valid if hypocentral_separation used.

Warning

Works in place on object, if you need to keep the original safe then run this on a copy of the object!

Example

>>> party = Party().read()
>>> len(party)
4
>>> declustered = party.decluster(20)
>>> len(party)
3

Example using epicentral-distance

>>> party = Party().read()
>>> len(party)
4
>>> # Calculate the expected origins based on the template origin.
>>> for family in party:
...     for detection in family:
...         _ = detection._calculate_event(template=family.template)
>>> declustered = party.decluster(20, hypocentral_separation=1)
>>> len(party)
4
>>> declustered = party.decluster(20, hypocentral_separation=100)
>>> len(party)
3
get_catalog()[source]

Get an obspy catalog object from the party.

Returns:

obspy.core.event.Catalog

Example

>>> party = Party().read()
>>> cat = party.get_catalog()
>>> print(len(cat))
4
lag_calc(stream, pre_processed, shift_len=0.2, min_cc=0.4, min_cc_from_mean_cc_factor=None, vertical_chans=['Z'], horizontal_chans=['E', 'N', '1', '2'], cores=1, interpolate=False, plot=False, plotdir=None, parallel=True, process_cores=None, ignore_length=False, ignore_bad_data=False, export_cc=False, cc_dir=None, **kwargs)[source]

Compute picks based on cross-correlation alignment.

Works in-place on events in Party.

Parameters:
  • stream (obspy.core.stream.Stream) – All the data needed to cut from - can be a gappy Stream.

  • pre_processed (bool) – Whether the stream has been pre-processed or not to match the templates. See note below.

  • shift_len (float) – Shift length allowed for the pick in seconds, will be plus/minus this amount - default=0.2

  • min_cc (float) – Minimum cross-correlation value to be considered a pick, default=0.4.

  • min_cc_from_mean_cc_factor (float) – If set to a value other than None, then the minimum cross- correlation value for a trace is set individually for each detection based on: min(detect_val / n_chans * min_cc_from_mean_cc_factor, min_cc).

  • horizontal_chans (list) – List of channel endings for horizontal-channels, on which S-picks will be made.

  • vertical_chans (list) – List of channel endings for vertical-channels, on which P-picks will be made.

  • cores (int) – Number of cores to use in parallel processing, defaults to one.

  • interpolate (bool) – Interpolate the correlation function to achieve sub-sample precision.

  • plot (bool) – To generate a plot for every detection or not, defaults to False

  • plotdir (str) – The path to save plots to. If plotdir=None (default) then the figure will be shown on screen.

  • export_cc (bool) – To generate a binary file in NumPy for every detection or not, defaults to False

  • cc_dir (str) – Path to saving folder, NumPy files will be output here.

  • parallel (bool) – Turn parallel processing on or off.

  • process_cores (int) – Number of processes to use for pre-processing (if different to cores).

  • ignore_length (bool) – If using daylong=True, then dayproc will try check that the data are there for at least 80% of the day, if you don’t want this check (which will raise an error if too much data are missing) then set ignore_length=True. This is not recommended!

  • ignore_bad_data (bool) – If False (default), errors will be raised if data are excessively gappy or are mostly zeros. If True then no error will be raised, but an empty trace will be returned (and not used in detection).

Returns:

Catalog of events with picks. No origin information is included. These events can then be written out via obspy.core.event.Catalog.write(), or to Nordic Sfiles using eqcorrscan.utils.sfile_util.eventtosfile() and located externally.

Return type:

obspy.core.event.Catalog

Note

Note on pre-processing: You can provide a pre-processed stream, which may be beneficial for detections over large time periods (the stream can have gaps, which reduces memory usage). However, in this case the processing steps are not checked, so you must ensure that all the template in the Party have the same sampling rate and filtering as the stream. If pre-processing has not be done then the data will be processed according to the parameters in the templates, in this case templates will be grouped by processing parameters and run with similarly processed data. In this case, all templates do not have to have the same processing parameters.

Note

Picks are corrected for the template pre-pick time.

min_chans(min_chans)[source]

Remove detections using min_chans or fewer channels.

Parameters:

min_chans (int) – Detections using more than this number of channels are maintained. Note that this is a strict if detection.no_chans > min_chans: rather than >=. Maintained for backwards compatability.

Returns:

Party

Note

Works in place on Party.

Example

>>> party = Party().read()
>>> print(len(party))
4
>>> party = party.min_chans(5)
>>> print(len(party))
1
read(filename=None, read_detection_catalog=True, estimate_origin=True)[source]

Read a Party from a file.

Parameters:
  • filename (str) – File to read from - can be a list of files, and can contain wildcards.

  • read_detection_catalog (bool) – Whether to read the detection catalog or not, if False, catalog will be regenerated - for large catalogs this can be faster.

  • estimate_origin (bool) – If True and no catalog is found, or read_detection_catalog is False then new events with origins estimated from the template origin time will be created.

Example

>>> Party().read()
Party of 4 Families.
rethreshold(new_threshold, new_threshold_type='MAD', abs_values=False)[source]

Remove detections from the Party that are below a new threshold.

Note

threshold can only be set higher.

Warning

Works in place on Party.

Parameters:
  • new_threshold (float) – New threshold level

  • new_threshold_type (str) – Either ‘MAD’, ‘absolute’ or ‘av_chan_corr’

  • abs_values (bool) – Whether to compare the absolute value of the detection-value.

Examples

Using the MAD threshold on detections made using the MAD threshold:

>>> party = Party().read()
>>> len(party)
4
>>> party = party.rethreshold(10.0)
>>> len(party)
4
>>> # Note that all detections are self detections

Using the absolute thresholding method on the same Party:

>>> party = Party().read().rethreshold(5.9, 'absolute')
>>> len(party)
1

Using the av_chan_corr method on the same Party:

>>> party = Party().read().rethreshold(0.9, 'av_chan_corr')
>>> len(party)
4
select(template_name)[source]

Select a specific family from the party.

Parameters:

template_name (str) – Template name of Family to select from a party.

Returns:

Family

sort()[source]

Sort the families by template name.

Example

>>> from eqcorrscan import Family, Template
>>> party = Party(families=[Family(template=Template(name='b')),
...                         Family(template=Template(name='a'))])
>>> party[0]
Family of 0 detections from template b
>>> party.sort()[0]
Family of 0 detections from template a
write(filename, format='tar', write_detection_catalog=True, catalog_format='QUAKEML', overwrite=False)[source]

Write Family out, select output format.

Parameters:
  • format (str) – One of either ‘tar’, ‘csv’, or any obspy supported catalog output. See note below on formats

  • filename (str) – Path to write file to.

  • write_detection_catalog (bool) – Whether to write the detection catalog object or not - writing large catalog files can be slow, and catalogs can be reconstructed from the Tribe.

  • catalog_format (str) – What format to write the detection-catalog with. Only Nordic, SC3ML, QUAKEML are supported. Note that not all information is written for all formats (QUAKEML is the most complete, but is slow for IO).

  • overwrite (bool) – Specifies whether detection-files are overwritten if they exist already. By default, no files are overwritten.

Note

csv format will write out detection objects, all other outputs will write the catalog. These cannot be rebuilt into a Family object. The only format that can be read back into Family objects is the ‘tar’ type.

Note

We recommend writing to the ‘tar’ format, which will write out all the template information (wavefiles as miniseed and metadata) alongside the detections and store these in a tar archive. This is readable by other programs and maintains all information required for further study.

Example

>>> party = Party().read()
>>> party.write('test_tar_write', format='tar')
Party of 4 Families.
>>> party.write('test_csv_write.csv', format='csv')
Party of 4 Families.
>>> party.write('test_quakeml.xml', format='quakeml')
Party of 4 Families.

Functions

read_party

Read detections and metadata from a tar archive.