core.match_filter.template

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.template.Template(name=None, st=None, lowcut=None, highcut=None, samp_rate=None, filt_order=None, process_length=None, prepick=None, event=None)[source]

Template object holder.

Contains waveform data and metadata parameters used to generate the template.

Parameters:
  • name (str) – A name associated with the template as a parsable string (no spaces allowed)

  • st (obspy.core.stream.Stream) – The Stream object containing the Template waveform data

  • lowcut (float) – The low-cut filter used to achieve st (float, Hz)

  • highcut (float) – The high-cut filter used to achieve st (float, Hz)

  • samp_rate (float) – The Sampling-rate of the template (float, Hz) - note that this should be the same as the sampling-rate of all traces in st.

  • filt_order (int) – The order of the filter applied to achieve st (int), see pre-processing functions for a more complete description

  • process_length (float) – The length of data (in seconds) processed to achieve st. e.g. if you processed a day of data then cut the template st from this, then process_length would be 86400.0 (float)

  • prepick (float) – The time before picks that waveforms were cut to make st (seconds, float)

  • event (obspy.core.event.Event) – The Event associated with the template (obspy.core.event.Event)

Methods

construct(method, name, lowcut, highcut, ...)

Construct a template using a given method.

copy()

Returns a copy of the template.

detect(stream, threshold, threshold_type, ...)

Detect using a single template within a continuous stream.

read(filename)

Read template from tar format with metadata.

same_processing(other)

Check is the templates are processed the same.

write(filename[, format])

Write template.

__init__(name=None, st=None, lowcut=None, highcut=None, samp_rate=None, filt_order=None, process_length=None, prepick=None, event=None)[source]
construct(method, name, lowcut, highcut, samp_rate, filt_order, length, prepick, swin='all', process_len=86400, all_horiz=False, delayed=True, plot=False, plotdir=None, min_snr=None, parallel=False, num_cores=False, skip_short_chans=False, **kwargs)[source]

Construct a template using a given method.

Parameters:
  • method (str) – Method to make the template, the only available method is: from_sac. For all other methods (from_client and from_meta_file) use Tribe.construct().

  • name (str) – Name for the template

  • lowcut (float) – Low cut (Hz), if set to None will not apply a lowcut

  • highcut (float) – High cut (Hz), if set to None will not apply a highcut.

  • samp_rate (float) – New sampling rate in Hz.

  • filt_order (int) – Filter level (number of corners).

  • length (float) – Length of template waveform in seconds.

  • prepick (float) – Pre-pick time in seconds

  • swin (str) – P, S, P_all, S_all or all, defaults to all: see note in eqcorrscan.core.template_gen.template_gen()

  • process_len (int) – Length of data in seconds to download and process.

  • all_horiz (bool) – To use both horizontal channels even if there is only a pick on one of them. Defaults to False.

  • delayed (bool) – If True, each channel will begin relative to it’s own pick-time, if set to False, each channel will begin at the same time.

  • plot (bool) – Plot templates or not.

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

  • min_snr (float) – Minimum signal-to-noise ratio for a channel to be included in the template, where signal-to-noise ratio is calculated as the ratio of the maximum amplitude in the template window to the rms amplitude in the whole window given.

  • parallel (bool) – Whether to process data in parallel or not.

  • num_cores (int) – Number of cores to try and use, if False and parallel=True, will use either all your cores, or as many traces as in the data (whichever is smaller).

  • skip_short_chans (bool) – Whether to ignore channels that have insufficient length data or not. Useful when the quality of data is not known, e.g. when downloading old, possibly triggered data from a datacentre

Note

method=from_sac requires the following kwarg(s): :param list sac_files:

osbpy.core.stream.Stream of sac waveforms, or list of paths to sac waveforms.

Note

See eqcorrscan.utils.sac_util.sactoevent for details on how pick information is collected.

Example

>>> # Get the path to the test data
>>> import eqcorrscan
>>> import os, glob
>>> TEST_PATH = (
...     os.path.dirname(eqcorrscan.__file__) + '/tests/test_data')
>>> sac_files = glob.glob(TEST_PATH + '/SAC/2014p611252/*')
>>> template = Template().construct(
...     method='from_sac', name='test', lowcut=2.0, highcut=8.0,
...     samp_rate=20.0, filt_order=4, prepick=0.1, swin='all',
...     length=2.0, sac_files=sac_files)
>>> print(template) 
Template test:
 12 channels;
 lowcut: 2.0 Hz;
 highcut: 8.0 Hz;
 sampling rate 20.0 Hz;
 filter order: 4;
 process length: 300.0 s

This will raise an error if the method is unsupported:

>>> template = Template().construct(
...     method='from_meta_file', name='test', lowcut=2.0, highcut=8.0,
...     samp_rate=20.0, filt_order=4, prepick=0.1, swin='all',
...     length=2.0) 
Traceback (most recent call last):
NotImplementedError: Method is not supported, use         Tribe.construct instead.
copy()[source]

Returns a copy of the template.

Returns:

Copy of template

Example

>>> from obspy import read
>>> template_a = Template(
...     name='a', st=read(), lowcut=2.0, highcut=8.0, samp_rate=100,
...     filt_order=4, process_length=3600, prepick=0.5)
>>> template_b = template_a.copy()
>>> template_a == template_b
True
detect(stream, threshold, threshold_type, trig_int, plot=False, plotdir=None, pre_processed=False, daylong=False, parallel_process=True, xcorr_func=None, concurrency=None, cores=None, ignore_length=False, overlap='calculate', full_peaks=False, **kwargs)[source]

Detect using a single template within a continuous stream.

Parameters:
  • stream (obspy.core.stream.Stream) – Continuous data to detect within using the Template.

  • threshold (float) – Threshold level, if using threshold_type=’MAD’ then this will be the multiple of the median absolute deviation.

  • threshold_type (str) – The type of threshold to be used, can be MAD, absolute or av_chan_corr. See Note on thresholding below.

  • trig_int (float) – Minimum gap between detections in seconds. If multiple detections occur within trig_int of one-another, the one with the highest cross-correlation sum will be selected.

  • plot (bool) – Turn plotting on or off.

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

  • pre_processed (bool) – Set to True if stream has already undergone processing, in this case eqcorrscan will only check that the sampling rate is correct. Defaults to False, which will use the eqcorrscan.utils.pre_processing routines to resample and filter the continuous data.

  • daylong (bool) – Set to True to use the eqcorrscan.utils.pre_processing.dayproc() routine, which preforms additional checks and is more efficient for day-long data over other methods.

  • parallel_process (bool) –

  • xcorr_func (str or callable) – A str of a registered xcorr function or a callable for implementing a custom xcorr function. For more details see eqcorrscan.utils.correlate.register_array_xcorr().

  • concurrency (str) – The type of concurrency to apply to the xcorr function. Options are ‘multithread’, ‘multiprocess’, ‘concurrent’. For more details see eqcorrscan.utils.correlate.get_stream_xcorr().

  • cores (int) – Number of workers for processing and detection.

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

  • overlap (float) – Either None, “calculate” or a float of number of seconds to overlap detection streams by. This is to counter the effects of the delay-and-stack in calculating cross-correlation sums. Setting overlap = “calculate” will work out the appropriate overlap based on the maximum lags within templates.

  • full_peaks (bool) – See eqcorrscan.utils.findpeaks.find_peaks2_short()

Returns:

Family of detections.

Note

When using the “fftw” correlation backend the length of the fft can be set. See eqcorrscan.utils.correlate for more info.

Note

stream must not be pre-processed. If your data contain gaps you should NOT fill those gaps before using this method. The pre-process functions (called within) will fill the gaps internally prior to processing, process the data, then re-fill the gaps with zeros to ensure correlations are not incorrectly calculated within gaps. If your data have gaps you should pass a merged stream without the fill_value argument (e.g.: stream = stream.merge()).

Note

Data overlap:

Internally this routine shifts and trims the data according to the offsets in the template (e.g. if trace 2 starts 2 seconds after trace 1 in the template then the continuous data will be shifted by 2 seconds to align peak correlations prior to summing). Because of this, detections at the start and end of continuous data streams may be missed. The maximum time-period that might be missing detections is the maximum offset in the template.

To work around this, if you are conducting matched-filter detections through long-duration continuous data, we suggest using some overlap (a few seconds, on the order of the maximum offset in the templates) in the continous data. You will then need to post-process the detections (which should be done anyway to remove duplicates).

Note

Thresholding:

MAD threshold is calculated as the:

\[threshold {\times} (median(abs(cccsum)))\]

where \(cccsum\) is the cross-correlation sum for a given template.

absolute threshold is a true absolute threshold based on the cccsum value.

av_chan_corr is based on the mean values of single-channel cross-correlations assuming all data are present as required for the template, e.g:

\[av\_chan\_corr\_thresh=threshold \times (cccsum / len(template))\]

where \(template\) is a single template from the input and the length is the number of channels within this template.

Note

See tutorials for example.

read(filename)[source]

Read template from tar format with metadata.

Parameters:

filename (str) – Filename to read template from.

Example

>>> from obspy import read
>>> template_a = Template(
...     name='a', st=read(), lowcut=2.0, highcut=8.0, samp_rate=100,
...     filt_order=4, process_length=3600, prepick=0.5)
>>> template_a.write(
...     'test_template_read') 
Template a:
 3 channels;
 lowcut: 2.0 Hz;
 highcut: 8.0 Hz;
 sampling rate 100 Hz;
 filter order: 4;
 process length: 3600 s
>>> template_b = Template().read('test_template_read.tgz')
>>> template_a == template_b
True
same_processing(other)[source]

Check is the templates are processed the same.

Example

>>> from obspy import read
>>> template_a = Template(
...     name='a', st=read(), lowcut=2.0, highcut=8.0, samp_rate=100,
...     filt_order=4, process_length=3600, prepick=0.5)
>>> template_b = template_a.copy()
>>> template_a.same_processing(template_b)
True
>>> template_b.lowcut = 5.0
>>> template_a.same_processing(template_b)
False
write(filename, format='tar')[source]

Write template.

Parameters:
  • filename (str) – Filename to write to, if it already exists it will be opened and appended to, otherwise it will be created.

  • format (str) – Format to write to, either ‘tar’ (to retain metadata), or any obspy supported waveform format to just extract the waveform.

Example

>>> from obspy import read
>>> template_a = Template(
...     name='a', st=read(), lowcut=2.0, highcut=8.0, samp_rate=100,
...     filt_order=4, process_length=3600, prepick=0.5)
>>> template_a.write('test_template') 
Template a:
 3 channels;
 lowcut: 2.0 Hz;
 highcut: 8.0 Hz;
 sampling rate 100 Hz;
 filter order: 4;
 process length: 3600 s
>>> template_a.write('test_waveform.ms',
...                  format='MSEED') 
Template a:
 3 channels;
 lowcut: 2.0 Hz;
 highcut: 8.0 Hz;
 sampling rate 100 Hz;
 filter order: 4;
 process length: 3600 s

Functions

eqcorrscan.core.match_filter.template.read_template(fname)[source]

Read a Template object from a tar archive.

Parameters:

fname (str) – Filename to read from

Returns:

eqcorrscan.core.match_filter.Template