obspy.core.trace.Trace
- class Trace(data=array([], dtype=float64), header=None)[source]
Bases:
object
An object containing data of a continuous series, such as a seismic trace.
- Parameters:
data (
ndarray
orMaskedArray
) – Array of data samplesheader (dict or
Stats
) – Dictionary containing header fields
- Variables:
id – A SEED compatible identifier of the trace.
stats – A container
Stats
for additional header information of the trace.data – Data samples in a
ndarray
orMaskedArray
Note
The
.data
attribute containing the time series samples as anumpy.ndarray
will always be made contiguous in memory. This way it is always safe to use.data
in routines that internally pass the array to C code. On the other hand this might result in some unwanted copying of data in memory. Experts can opt-out by settingTrace._always_contiguous = False
, in this case the user has to make sure themselves that no C operations are performed on potentially incontiguous data.Supported Operations
trace = traceA + traceB
Merges traceA and traceB into one new trace object. See also:
Trace.__add__()
.len(trace)
Returns the number of samples contained in the trace. That is it es equal to
len(trace.data)
. See also:Trace.__len__()
.str(trace)
Returns basic information about the trace object. See also:
Trace.__str__()
.
Attributes
Return a SEED compatible identifier of the trace. |
|
Public Methods
Search for and attach channel response to the trace as |
|
Returns a deepcopy of the trace. |
|
Return number of data samples of the current trace. |
|
Downsample trace data by an integer factor. |
|
Remove a trend from the trace. |
|
Differentiate the trace with respect to time. |
|
Filter the data of the current trace. |
|
Return a SEED compatible identifier of the trace. |
|
Integrate the trace with respect to time. |
|
Interpolate the data using various interpolation techniques. |
|
Returns the value of the absolute maximum amplitude in the trace. |
|
Normalize the trace to its absolute maximum. |
|
Create a simple graph of the current trace. |
|
Deconvolve instrument response. |
|
Remove instrument sensitivity. |
|
Resample trace data using Fourier method. |
|
Correct for instrument response / Simulate new instrument response. |
|
Return a new Trace object with data going from start to end time. |
|
Generator yielding equal length sliding windows of the Trace. |
|
Create a spectrogram plot of the trace. |
|
Split Trace object containing gaps using a NumPy masked array into several traces. |
|
Method to get the standard deviation of amplitudes in the trace. |
|
Taper the trace. |
|
For convenient plotting compute a NumPy array with timing information of all samples in the Trace. |
|
Run a triggering algorithm on the data of the current trace. |
|
Cut current trace to given start and end time. |
|
Verify current trace object against available meta data. |
|
Save current trace into a file. |
Private Methods
Warning
Private methods are mainly for internal/developer use and their API might change without notice.
- Trace._get_response(inventories)[source]
Search for and return channel response for the trace.
- Trace._internal_add_processing_info(info)[source]
Add the given informational string to the processing field in the trace’s
Stats
object.
- Trace._ltrim(starttime, pad=False, nearest_sample=True, fill_value=None)[source]
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)
- Trace._rtrim(endtime, pad=False, nearest_sample=True, fill_value=None)[source]
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
- Trace.__add__(trace, method=0, interpolation_samples=0, fill_value=None, sanity_checks=True)[source]
Add another Trace object to current trace.
- Parameters:
method (int, optional) – Method to handle overlaps of traces. Defaults to
0
. See the Handling Overlaps section below for further details.fill_value (int, float, str or
None
, optional) – Fill value for gaps. Defaults toNone
. Traces will be converted to NumPy masked arrays if no value is given and gaps are present. If the keyword'latest'
is provided it will use the latest value before the gap. If keyword'interpolate'
is provided, missing values are linearly interpolated (not changing the data type e.g. of integer valued traces). See the Handling Gaps section below for further details.interpolation_samples (int, optional) – Used only for
method=1
. It specifies the number of samples which are used to interpolate between overlapping traces. Defaults to0
. If set to-1
all overlapping samples are interpolated.sanity_checks (bool, optional) – Enables some sanity checks before merging traces. Defaults to
True
.
Trace data will be converted into a NumPy masked array data type if any gaps are present. This behavior may be prevented by setting the
fill_value
parameter. Themethod
argument controls the handling of overlapping data values.Sampling rate, data type and trace.id of both traces must match.
Handling Overlaps
Method
Description
0
Discard overlapping data. Overlaps are essentially treated the same way as gaps:
Trace 1: AAAAAAAA Trace 2: FFFFFFFF 1 + 2 : AAAA----FFFF
Contained traces with differing data will be marked as gap:
Trace 1: AAAAAAAAAAAA Trace 2: FF 1 + 2 : AAAA--AAAAAA
Missing data can be merged in from a different trace:
Trace 1: AAAA--AAAAAA (contained trace, missing samples) Trace 2: FF 1 + 2 : AAAAFFAAAAAA
1
Discard data of the previous trace assuming the following trace contains data with a more correct time value. The parameter
interpolation_samples
specifies the number of samples used to linearly interpolate between the two traces in order to prevent steps. Note that if there are gaps inside, the returned array is still a masked array, only iffill_value
is set, the returned array is a normal array and gaps are filled with fill value.No interpolation (
interpolation_samples=0
):Trace 1: AAAAAAAA Trace 2: FFFFFFFF 1 + 2 : AAAAFFFFFFFF
Interpolate first two samples (
interpolation_samples=2
):Trace 1: AAAAAAAA Trace 2: FFFFFFFF 1 + 2 : AAAACDFFFFFF (interpolation_samples=2)
Interpolate all samples (
interpolation_samples=-1
):Trace 1: AAAAAAAA Trace 2: FFFFFFFF 1 + 2 : AAAABCDEFFFF
Any contained traces with different data will be discarded:
Trace 1: AAAAAAAAAAAA (contained trace) Trace 2: FF 1 + 2 : AAAAAAAAAAAA
Missing data can be merged in from a different trace:
Trace 1: AAAA--AAAAAA (contained trace, missing samples) Trace 2: FF 1 + 2 : AAAAFFAAAAAA
Handling gaps
Traces with gaps and
fill_value=None
(default):Trace 1: AAAA Trace 2: FFFF 1 + 2 : AAAA----FFFF
Traces with gaps and given
fill_value=0
:Trace 1: AAAA Trace 2: FFFF 1 + 2 : AAAA0000FFFF
Traces with gaps and given
fill_value='latest'
:Trace 1: ABCD Trace 2: FFFF 1 + 2 : ABCDDDDDFFFF
Traces with gaps and given
fill_value='interpolate'
:Trace 1: AAAA Trace 2: FFFF 1 + 2 : AAAABCDEFFFF
- Trace.__delattr__(name, /)
Implement delattr(self, name).
- Trace.__dir__()
Default dir() implementation.
- Trace.__eq__(other)[source]
Implements rich comparison of Trace objects for “==” operator.
Traces are the same, if both their data and stats are the same.
- Trace.__format__(format_spec, /)
Default object formatter.
- Trace.__getattribute__(name, /)
Return getattr(self, name).
- Trace.__getitem__(index)[source]
__getitem__ method of Trace object.
- Return type:
- Returns:
List of data points
- Trace.__init_subclass__()
This method is called when a class is subclassed.
The default implementation does nothing. It may be overridden to extend subclasses.
- Trace.__len__()[source]
Return number of data samples of the current trace.
- Return type:
- Returns:
Number of data samples.
Example
>>> trace = Trace(data=np.array([1, 2, 3, 4])) >>> trace.count() 4 >>> len(trace) 4
- Trace.__mod__(num)[source]
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
- Trace.__mul__(num)[source]
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
- Trace.__ne__(other)[source]
Implements rich comparison of Trace objects for “!=” operator.
Calls __eq__() and returns the opposite.
- Trace.__new__(**kwargs)
- Trace.__reduce__()
Helper for pickle.
- Trace.__reduce_ex__(protocol, /)
Helper for pickle.
- Trace.__repr__()
Return repr(self).
- Trace.__sizeof__()
Size of object in memory, in bytes.
- Trace.__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.
- Trace.__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).
- Trace.__truediv__(num)[source]
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