6.2.15.1.1. eqcorrscan.utils.trigger.network_trigger¶
- eqcorrscan.utils.trigger.network_trigger(st, parameters, thr_coincidence_sum, moveout, max_trigger_length=60, despike=True, parallel=True)[source]¶
Main function to compute triggers for a network of stations. Computes single-channel characteristic functions using given parameters, then combines these to find network triggers on a number of stations within a set moveout window.
- Parameters:
st (obspy.core.stream.Stream) – Stream to compute detections within
parameters (list) – List of parameter class
thr_coincidence_sum (int) – Minimum number of stations required to raise a network trigger.
moveout (float) – Window to find triggers within in the network detection stage.
max_trigger_length (float) – Maximum trigger length in seconds, used to remove long triggers - can set to False to not use.
despike (bool) – Whether to apply simple despiking routine or not
parallel (bool) – Whether to run in parallel or not
- Returns:
List of triggers
- Return type:
list
Example
>>> from obspy import read >>> from eqcorrscan.utils.trigger import TriggerParameters, network_trigger >>> st = read("https://examples.obspy.org/" + ... "example_local_earthquake_3stations.mseed") >>> parameters = [] >>> for tr in st: ... parameters.append(TriggerParameters({'station': tr.stats.station, ... 'channel': tr.stats.channel, ... 'sta_len': 0.5, ... 'lta_len': 10.0, ... 'thr_on': 10.0, ... 'thr_off': 3.0, ... 'lowcut': 2.0, ... 'highcut': 15.0})) >>> triggers = network_trigger(st=st, parameters=parameters, ... thr_coincidence_sum=5, moveout=30, ... max_trigger_length=60, despike=False) >>> print(len(triggers)) 1