obspy.realtime.rttrace.RtTrace

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

Bases: obspy.core.trace.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')[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

__dict__
__doc__ str(object=’‘) -> str
__hash__
__module__ str(object=’‘) -> str
__weakref__ list of weak references to the object (if defined)
have_appended_data bool(x) -> bool
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
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.
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. Spectra are linearly
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
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
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.

_get_response Search for and return channel response for the trace.
_internal_add_processing_info Add the given informational string to the processing field in the
_ltrim Cut current trace to given start time.
_repr_pretty_
_rtrim Cut current trace to given end time.

Special Methods

__add__ Too ambiguous, throw an Error.
__dir__ default dir() implementation
__div__ Split Trace into new Stream containing num Traces of the same size.
__eq__ Implements rich comparison of RtTrace objects for “==” operator.
__format__ default object formatter
__ge__ Too ambiguous, throw an Error.
__getitem__ __getitem__ method of Trace object.
__gt__ Too ambiguous, throw an Error.
__init__ Initializes an RtTrace.
__le__ Too ambiguous, throw an Error.
__len__ Return number of data samples of the current trace.
__lt__ Too ambiguous, throw an Error.
__mod__ Split Trace into new Stream containing Traces with num samples.
__mul__ Create a new Stream containing num copies of this trace.
__ne__ Implements rich comparison of Trace objects for ”!=” operator.
__new__ Create and return a new object.
__nonzero__ No data means no trace.
__original_str__ Return short summary string of the current trace.
__reduce__ helper for pickle
__reduce_ex__ helper for pickle
__setattr__ __setattr__ method of Trace object.
__sizeof__ size of object in memory, in bytes
__str__ Monkey patch for the __str__ method of the Trace object. SEGY object do not
__subclasshook__ Abstract classes can override this to customize issubclass().
__truediv__ Split Trace into new Stream containing num Traces of the same size.