core.match_filter.family

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.family.Family(template, detections=None, catalog=None)[source]

Container for Detection objects from a single template.

Parameters:
  • template (eqcorrscan.core.match_filter.Template) – The template used to detect the family

  • detections (list) – list of Detection objects

  • catalog (obspy.core.event.Catalog) – Catalog of detections, with information for the individual detections.

Methods

append(other)

Add another family or detection to the family.

copy()

Returns a copy of the family.

extract_streams(stream, length, prepick[, ...])

Generate a dictionary of cut streams around detections.

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

Compute picks based on cross-correlation alignment.

plot([plot_grouped])

Plot the cumulative number of detections in time.

sort()

Sort by detection time.

write(filename[, format, overwrite])

Write Family out, select output format.

__init__(template, detections=None, catalog=None)[source]

Instantiation of Family object.

append(other)[source]

Add another family or detection to the family.

Example

Append a family to a family

>>> from eqcorrscan import Template, Detection
>>> family_a = Family(template=Template(name='a'))
>>> family_b = Family(template=Template(name='a'))
>>> family_a.append(family_b)
Family of 0 detections from template a

Append a detection to the family

>>> family_a = Family(template=Template(name='a'))
>>> detection = Detection(
...     template_name='a', detect_time=UTCDateTime(), no_chans=5,
...     detect_val=2.5, threshold=1.23, typeofdet='corr',
...     threshold_type='MAD', threshold_input=8.0)
>>> family_a.append(detection)
Family of 1 detections from template a
copy()[source]

Returns a copy of the family.

Returns:

Copy of family

Example

>>> from eqcorrscan import Template, Detection
>>> family = Family(
...     template=Template(name='a'), detections=[
...     Detection(template_name='a', detect_time=UTCDateTime(0),
...               no_chans=8, detect_val=4.2, threshold=1.2,
...               typeofdet='corr', threshold_type='MAD',
...               threshold_input=8.0),
...     Detection(template_name='a', detect_time=UTCDateTime(0) + 10,
...               no_chans=8, detect_val=4.5, threshold=1.2,
...               typeofdet='corr', threshold_type='MAD',
...               threshold_input=8.0)])
>>> family == family.copy()
True
extract_streams(stream, length, prepick, all_vert=False, all_horiz=False, vertical_chans=['Z'], horizontal_chans=['E', 'N', '1', '2'])[source]

Generate a dictionary of cut streams around detections.

Parameters:
  • stream (obspy.core.stream.Stream) – Stream of data to cut from

  • length (float) – Length of data to extract in seconds.

  • prepick (float) – Length before the expected pick on each channel to start the cut stream in seconds.

Return type:

dict

Returns:

Dictionary of cut streams keyed by detection id.

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 the Family

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.

  • 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).

  • 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.

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 the template in the Family has 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 template.

Note

Picks are corrected for the template pre-pick time.

plot(plot_grouped=False)[source]

Plot the cumulative number of detections in time.

Example

>>> from eqcorrscan import Template, Detection
>>> family = Family(
...     template=Template(name='a'), detections=[
...     Detection(template_name='a', detect_time=UTCDateTime(0) + 200,
...               no_chans=8, detect_val=4.2, threshold=1.2,
...               typeofdet='corr', threshold_type='MAD',
...               threshold_input=8.0),
...     Detection(template_name='a', detect_time=UTCDateTime(0),
...               no_chans=8, detect_val=4.5, threshold=1.2,
...               typeofdet='corr', threshold_type='MAD',
...               threshold_input=8.0),
...     Detection(template_name='a', detect_time=UTCDateTime(0) + 10,
...               no_chans=8, detect_val=4.5, threshold=1.2,
...               typeofdet='corr', threshold_type='MAD',
...               threshold_input=8.0)])
>>> family.plot(plot_grouped=True)  

(Source code)

sort()[source]

Sort by detection time.

Example

>>> from eqcorrscan import Template, Detection
>>> family = Family(
...     template=Template(name='a'), detections=[
...     Detection(template_name='a', detect_time=UTCDateTime(0) + 200,
...               no_chans=8, detect_val=4.2, threshold=1.2,
...               typeofdet='corr', threshold_type='MAD',
...               threshold_input=8.0),
...     Detection(template_name='a', detect_time=UTCDateTime(0),
...               no_chans=8, detect_val=4.5, threshold=1.2,
...               typeofdet='corr', threshold_type='MAD',
...               threshold_input=8.0),
...     Detection(template_name='a', detect_time=UTCDateTime(0) + 10,
...               no_chans=8, detect_val=4.5, threshold=1.2,
...               typeofdet='corr', threshold_type='MAD',
...               threshold_input=8.0)])
>>> family[0].detect_time
UTCDateTime(1970, 1, 1, 0, 3, 20)
>>> family.sort()[0].detect_time
UTCDateTime(1970, 1, 1, 0, 0)
write(filename, format='tar', overwrite=False)[source]

Write Family out, select output format.

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

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

  • 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

csv format will append detections to filename, all others will overwrite any existing files.

Example

>>> from eqcorrscan import Template, Detection
>>> from obspy import read
>>> family = Family(
...     template=Template(name='a', st=read()), detections=[
...     Detection(template_name='a', detect_time=UTCDateTime(0) + 200,
...               no_chans=8, detect_val=4.2, threshold=1.2,
...               typeofdet='corr', threshold_type='MAD',
...               threshold_input=8.0),
...     Detection(template_name='a', detect_time=UTCDateTime(0),
...               no_chans=8, detect_val=4.5, threshold=1.2,
...               typeofdet='corr', threshold_type='MAD',
...               threshold_input=8.0),
...     Detection(template_name='a', detect_time=UTCDateTime(0) + 10,
...               no_chans=8, detect_val=4.5, threshold=1.2,
...               typeofdet='corr', threshold_type='MAD',
...               threshold_input=8.0)])
>>> family.write('test_family')