obspy.realtime.rttrace.RtTrace

class RtTrace(max_length=None, *args, **kwargs)[source]

Bases: Trace

An object containing data of a continuous series constructed dynamically from sequential data packets.

New data packets may be periodically appended. Registered time-domain processes can be applied to the new data and the resulting trace will be left trimmed to maintain a specified maximum trace length.

Parameters:

max_length (int, optional) – maximum trace length in seconds

Example

RtTrace has been built to handle real time processing of periodically append data packets, such as adding and processing data requested from an SeedLink server. See obspy.clients.seedlink for further information.

For the sake of simplicity we will just split data of an existing example file into multiple chucks (Trace objects) of about equal size (step 1 + 2) and append those chunks in a simple loop (step 4) into an RtTrace object. Additionally there are two real time processing functions registered to the RtTrace object (step 3) which will automatically process any appended data chunks.

  1. Read first trace of example SAC data file and extract contained time offset and epicentral distance of an earthquake:

    >>> from obspy.realtime import RtTrace
    >>> from obspy import read
    >>> from obspy.realtime.signal import calculate_mwp_mag
    >>> data_trace = read(
    ...     '/path/to/II.TLY.BHZ.SAC', round_sampling_interval=False)[0]
    >>> len(data_trace)
    12684
    >>> ref_time_offset = data_trace.stats.sac.a
    >>> print(ref_time_offset)
    301.506
    >>> epicentral_distance = data_trace.stats.sac.gcarc
    >>> print(epicentral_distance)
    30.0855
    
  2. Split given trace into a list of three sub-traces:

    >>> traces = data_trace / 3
    >>> [len(tr) for tr in traces]
    [4228, 4228, 4228]
    
  3. Assemble real time trace and register two processes:

    >>> rt_trace = RtTrace()
    >>> rt_trace.register_rt_process('integrate')
    1
    >>> rt_trace.register_rt_process('mwpintegral', mem_time=240,
    ...     ref_time=(data_trace.stats.starttime + ref_time_offset),
    ...     max_time=120, gain=1.610210e+09)
    2
    
  4. Append and auto-process packet data into RtTrace:

    >>> for tr in traces:
    ...     processed_trace = rt_trace.append(tr, gap_overlap_check=True)
    ...
    >>> len(rt_trace)
    12684
    
  5. Some post processing to get Mwp:

    >>> peak = np.amax(np.abs(rt_trace.data))
    >>> print(peak)
    0.136404
    >>> mwp = calculate_mwp_mag(peak, epicentral_distance)
    >>> print(mwp)  
    8.78902911791...
    

Attributes

have_appended_data

id

Return a SEED compatible identifier of the trace.

meta

Public Methods

append

Appends a Trace object to this RtTrace.

attach_response

Search for and attach channel response to the trace as obspy.core.trace.Trace.stats.response.

copy

Returns a deepcopy of this RtTrace.

count

Return number of data samples of the current trace.

decimate

Downsample trace data by an integer factor.

detrend

Remove a trend from the trace.

differentiate

Differentiate the trace with respect to time.

filter

Filter the data of the current trace.

get_id

Return a SEED compatible identifier of the trace.

integrate

Integrate the trace with respect to time.

interpolate

Interpolate the data using various interpolation techniques.

max

Returns the value of the absolute maximum amplitude in the trace.

newbyteorder

Change byteorder of the data

normalize

Normalize the trace to its absolute maximum.

plot

Create a simple graph of the current trace.

register_rt_process

Adds real-time processing algorithm to processing list of this RtTrace.

remove_response

Deconvolve instrument response.

remove_sensitivity

Remove instrument sensitivity.

resample

Resample trace data using Fourier method.

rt_process_functions_to_string

Return doc string for all predefined real-time processing functions.

simulate

Correct for instrument response / Simulate new instrument response.

slice

Return a new Trace object with data going from start to end time.

slide

Generator yielding equal length sliding windows of the Trace.

spectrogram

Create a spectrogram plot of the trace.

split

Split Trace object containing gaps using a NumPy masked array into several traces.

std

Method to get the standard deviation of amplitudes in the trace.

taper

Taper the trace.

times

For convenient plotting compute a NumPy array with timing information of all samples in the Trace.

trigger

Run a triggering algorithm on the data of the current trace.

trim

Cut current trace to given start and end time.

verify

Verify current trace object against available meta data.

write

Save current trace into a file.

Private Methods

Warning

Private methods are mainly for internal/developer use and their API might change without notice.

RtTrace._get_response(inventories)

Search for and return channel response for the trace.

Parameters:

inventories (Inventory or Network or a list containing objects of these types or a string with a filename of a StationXML file.) – Station metadata to use in search for response for each trace in the stream.

Returns:

obspy.core.inventory.response.Response object

RtTrace._internal_add_processing_info(info)

Add the given informational string to the processing field in the trace’s Stats object.

RtTrace._ltrim(starttime, pad=False, nearest_sample=True, fill_value=None)

Cut current trace to given start time. For more info see trim().

Example

>>> tr = Trace(data=np.arange(0, 10))
>>> tr.stats.delta = 1.0
>>> tr._ltrim(tr.stats.starttime + 8)  
<...Trace object at 0x...>
>>> tr.data
array([8, 9])
>>> tr.stats.starttime
UTCDateTime(1970, 1, 1, 0, 0, 8)
RtTrace._repr_pretty_(p, cycle)
RtTrace._rtrim(endtime, pad=False, nearest_sample=True, fill_value=None)

Cut current trace to given end time. For more info see trim().

Example

>>> tr = Trace(data=np.arange(0, 10))
>>> tr.stats.delta = 1.0
>>> tr._rtrim(tr.stats.starttime + 2)  
<...Trace object at 0x...>
>>> tr.data
array([0, 1, 2])
>>> tr.stats.endtime
UTCDateTime(1970, 1, 1, 0, 0, 2)

Special Methods

RtTrace.__add__(**kwargs)[source]

Too ambiguous, throw an Error.

RtTrace.__delattr__(name, /)

Implement delattr(self, name).

RtTrace.__dir__()

Default dir() implementation.

RtTrace.__eq__(other)[source]

Implements rich comparison of RtTrace objects for “==” operator.

Traces are the same, if both their data and stats are the same.

RtTrace.__format__(format_spec, /)

Default object formatter.

RtTrace.__ge__(other)

Too ambiguous, throw an Error.

RtTrace.__getattribute__(name, /)

Return getattr(self, name).

RtTrace.__getitem__(index)

__getitem__ method of Trace object.

Return type:

list

Returns:

List of data points

RtTrace.__gt__(other)

Too ambiguous, throw an Error.

RtTrace.__init__(max_length=None, *args, **kwargs)[source]

Initializes an RtTrace.

See obspy.core.trace.Trace for all parameters.

RtTrace.__init_subclass__()

This method is called when a class is subclassed.

The default implementation does nothing. It may be overridden to extend subclasses.

RtTrace.__le__(other)

Too ambiguous, throw an Error.

RtTrace.__len__()

Return number of data samples of the current trace.

Return type:

int

Returns:

Number of data samples.

Example

>>> trace = Trace(data=np.array([1, 2, 3, 4]))
>>> trace.count()
4
>>> len(trace)
4
RtTrace.__lt__(other)

Too ambiguous, throw an Error.

RtTrace.__mod__(num)

Split Trace into new Stream containing Traces with num samples.

Parameters:

num (int) – Number of samples in each trace in returned Stream. Last trace may contain lesser samples.

Returns:

New ObsPy Stream object.

Example

>>> from obspy import read
>>> tr = read()[0]
>>> print(tr)  
BW.RJOB..EHZ | 2009-08-24T00:20:03.000000Z ... | 100.0 Hz, 3000 samples
>>> st = tr % 800
>>> print(st)  
4 Trace(s) in Stream:
BW.RJOB..EHZ | 2009-08-24T00:20:03.000000Z ... | 100.0 Hz, 800 samples
BW.RJOB..EHZ | 2009-08-24T00:20:11.000000Z ... | 100.0 Hz, 800 samples
BW.RJOB..EHZ | 2009-08-24T00:20:19.000000Z ... | 100.0 Hz, 800 samples
BW.RJOB..EHZ | 2009-08-24T00:20:27.000000Z ... | 100.0 Hz, 600 samples
RtTrace.__mul__(num)

Create a new Stream containing num copies of this trace.

Parameters:

num (int) – Number of copies.

Returns:

New ObsPy Stream object.

Example

>>> from obspy import read
>>> tr = read()[0]
>>> st = tr * 5
>>> len(st)
5
RtTrace.__ne__(other)

Implements rich comparison of Trace objects for “!=” operator.

Calls __eq__() and returns the opposite.

RtTrace.__new__(**kwargs)
RtTrace.__nonzero__()

No data means no trace.

RtTrace.__original_str__(id_length=None)

Return short summary string of the current trace.

Return type:

str

Returns:

Short summary string of the current trace containing the SEED identifier, start time, end time, sampling rate and number of points of the current trace.

Example

>>> tr = Trace(header={'station':'FUR', 'network':'GR'})
>>> str(tr)  
'GR.FUR.. | 1970-01-01T00:00:00.000000Z - ... | 1.0 Hz, 0 samples'
RtTrace.__reduce__()

Helper for pickle.

RtTrace.__reduce_ex__(protocol, /)

Helper for pickle.

RtTrace.__repr__()

Return repr(self).

RtTrace.__setattr__(key, value)

__setattr__ method of Trace object.

RtTrace.__sizeof__()

Size of object in memory, in bytes.

RtTrace.__str__(*args, **kwargs)

Monkey patch for the __str__ method of the Trace object. SEGY object do not have network, station, channel codes. It just prints the trace sequence number within the line.

RtTrace.__subclasshook__()

Abstract classes can override this to customize issubclass().

This is invoked early on by abc.ABCMeta.__subclasscheck__(). It should return True, False or NotImplemented. If it returns NotImplemented, the normal algorithm is used. Otherwise, it overrides the normal algorithm (and the outcome is cached).

RtTrace.__truediv__(num)

Split Trace into new Stream containing num Traces of the same size.

Parameters:

num (int) – Number of traces in returned Stream. Last trace may contain lesser samples.

Returns:

New ObsPy Stream object.

Example

>>> from obspy import read
>>> tr = read()[0]
>>> print(tr)  
BW.RJOB..EHZ | 2009-08-24T00:20:03.000000Z ... | 100.0 Hz, 3000 samples
>>> st = tr / 7
>>> print(st)  
7 Trace(s) in Stream:
BW.RJOB..EHZ | 2009-08-24T00:20:03.000000Z ... | 100.0 Hz, 429 samples
BW.RJOB..EHZ | 2009-08-24T00:20:07.290000Z ... | 100.0 Hz, 429 samples
BW.RJOB..EHZ | 2009-08-24T00:20:11.580000Z ... | 100.0 Hz, 429 samples
BW.RJOB..EHZ | 2009-08-24T00:20:15.870000Z ... | 100.0 Hz, 429 samples
BW.RJOB..EHZ | 2009-08-24T00:20:20.160000Z ... | 100.0 Hz, 429 samples
BW.RJOB..EHZ | 2009-08-24T00:20:24.450000Z ... | 100.0 Hz, 429 samples
BW.RJOB..EHZ | 2009-08-24T00:20:28.740000Z ... | 100.0 Hz, 426 samples