obspy.realtime.rttrace.RtTrace.interpolate

RtTrace.interpolate(sampling_rate, method='weighted_average_slopes', starttime=None, npts=None, time_shift=0.0, *args, **kwargs)

Interpolate the data using various interpolation techniques.

Be careful when downsampling data and make sure to apply an appropriate anti-aliasing lowpass filter before interpolating in case it’s necessary.

Note

The Trace object has three different methods to change the sampling rate of its data: resample(), decimate(), and interpolate().

Make sure to choose the most appropriate one for the problem at hand.

Note

This operation is performed in place on the actual data arrays. The raw data will no longer be accessible afterwards. To keep your original data, use copy() to create a copy of your Trace object.

Interpolation Methods:

The chosen method is crucial and we will elaborate a bit about the choices here:

  • "lanczos": This offers the highest quality interpolation and should be chosen whenever possible. It is only due to legacy reasons that this is not the default method. The one downside it has is that it can be fairly expensive. See the lanczos_interpolation() function for more details.
  • "weighted_average_slopes": This is the interpolation method used by SAC. Refer to weighted_average_slopes() for more details.
  • "slinear", "quadratic" and "cubic": spline interpolation of first, second or third order.
  • "linear": Linear interpolation.
  • "nearest": Nearest neighbour interpolation.
  • "zero": Last encountered value interpolation.

Parameters:

Parameters:
  • sampling_rate The new sampling rate in Hz.
  • method The kind of interpolation to perform as a string. One of "linear", "nearest", "zero", "slinear", "quadratic", "cubic", "lanczos", or "weighted_average_slopes". Alternatively an integer specifying the order of the spline interpolator to use also works. Defaults to "weighted_average_slopes".
  • starttime (UTCDateTime or int) The start time (or timestamp) for the new interpolated stream. Will be set to current start time of the trace if not given.
  • npts (int) The new number of samples. Will be set to the best fitting number to retain the current end time of the trace if not given.
  • time_shift (float)

    Shift the trace by adding time_shift to the starttime. The time shift is always given in seconds. A positive shift means the data is shifted towards the future, e.g. a positive time delta. Note that this parameter solely affects the metadata. The actual interpolation of the underlaying data is governed by the parameters sampling_rate, starttime and npts.

    Note

    Interpolation can also shift the data with subsample accuracy. Please note that a time shift in the Fourier domain is always more accurate than this. When using Lanczos interpolation with large values of a and away from the boundaries this is nonetheless pretty good.

New in version 0.11:

  • New parameter time_shift.
  • New interpolation method lanczos.

Usage Examples

>>> from obspy import read
>>> tr = read()[0]
>>> print(tr)  
BW.RJOB..EHZ | 2009-08-24T00:20:03... - ... | 100.0 Hz, 3000 samples
>>> tr.interpolate(sampling_rate=111.1)  
<obspy.core.trace.Trace object at 0x...>
>>> print(tr)  
BW.RJOB..EHZ | 2009-08-24T00:20:03... - ... | 111.1 Hz, 3332 samples

Setting starttime and/or npts will interpolate to sampling points with the given start time and/or number of samples. Extrapolation will not be performed.

>>> tr = read()[0]
>>> print(tr)  
BW.RJOB..EHZ | 2009-08-24T00:20:03... - ... | 100.0 Hz, 3000 samples
>>> tr.interpolate(sampling_rate=111.1,
...                starttime=tr.stats.starttime + 10)         
<obspy.core.trace.Trace object at 0x...>
>>> print(tr)  
BW.RJOB..EHZ | 2009-08-24T00:20:13... - ... | 111.1 Hz, 2221 samples