obspy.signal.cross_correlation.xcorr_pick_correction

xcorr_pick_correction(pick1, trace1, pick2, trace2, t_before, t_after, cc_maxlag, filter=None, filter_options={}, plot=False, filename=None)[source]

Calculate the correction for the differential pick time determined by cross correlation of the waveforms in narrow windows around the pick times. For details on the fitting procedure refer to [Deichmann1992].

The parameters depend on the epicentral distance and magnitude range. For small local earthquakes (Ml ~0-2, distance ~3-10 km) with consistent manual picks the following can be tried:

t_before=0.05, t_after=0.2, cc_maxlag=0.10,
filter="bandpass", filter_options={'freqmin': 1, 'freqmax': 20}

The appropriate parameter sets can and should be determined/verified visually using the option plot=True on a representative set of picks.

To get the corrected differential pick time calculate: ((pick2 + pick2_corr) - pick1). To get a corrected differential travel time using origin times for both events calculate: ((pick2 + pick2_corr - ot2) - (pick1 - ot1))

Parameters:
  • pick1 (UTCDateTime) – Time of pick for trace1.

  • trace1 (Trace) – Waveform data for pick1. Add some time at front/back. The appropriate part of the trace is used automatically.

  • pick2 (UTCDateTime) – Time of pick for trace2.

  • trace2 (Trace) – Waveform data for pick2. Add some time at front/back. The appropriate part of the trace is used automatically.

  • t_before (float) – Time to start cross correlation window before pick times in seconds.

  • t_after (float) – Time to end cross correlation window after pick times in seconds.

  • cc_maxlag (float) – Maximum lag/shift time tested during cross correlation in seconds.

  • filter (str) – None for no filtering or name of filter type as passed on to filter() if filter should be used. To avoid artifacts in filtering provide sufficiently long time series for trace1 and trace2.

  • filter_options (dict) – Filter options that get passed on to filter() if filtering is used.

  • plot (bool) – If True, a plot window illustrating the alignment of the two traces at best cross correlation will be shown. This can and should be used to verify the used parameters before running automatedly on large data sets.

  • filename (str) – If plot option is selected, specifying a filename here (e.g. ‘myplot.pdf’ or ‘myplot.png’) will output the plot to a file instead of opening a plot window.

Return type:

(float, float)

Returns:

Correction time pick2_corr for pick2 pick time as a float and corresponding correlation coefficient.