Network stream#

In order to ease the analysis of network seismic data, we offer a set of classes and methods designed for reading, pre-processing, and analyzing seismic data from these networks. The workflow involves working on traces that have been synchronized and pre-processed using consistent methods. To facilitate this, we introduce the concept of a NetworkStream, a subclass of the ObsPy Stream object, which retains all the methods of the Stream object while adding specialized methods for the pre-processing and synchronization of traces.

class covseisnet.stream.NetworkStream(*args, **kwargs)[source]#

Bases: Stream

Subclass of the ObsPy Stream tailored for managing continuous data from seismic networks. The class is designed to handle multiple traces from different stations, and provides additional methods for pre-processing and synchronization of the traces. It also provide network-wide methods such as the calculation of the common time vector of the traces.

Note

The following list of methods and attributes are the ones that are specific to the NetworkStream object. The full list of methods and attributes of the ObsPy Stream object are available in the ObsPy documentation.

Boolean attributes

  • synced — traces share the same time vector.

  • equal_rates — traces have the same sampling rate.

  • equal_length — traces have the same number of points.

Numeric attributes

  • sampling_rate — common traces sampling rate.

  • npts — common traces number of samples.

  • coordinates — geographical coordinates of the stations.

Methods

  • stats() — stats dictionaries of traces.

  • read() — read seismic waveforms files into a NetworkStream object.

  • times() — common traces time vector.

  • time_extent() — minimal time extent of traces in a stream.

  • cut() — trim stream between given start and end times with string format.

  • synchronize() — synchronize the traces into the same times with interpolation.

  • process() — process the traces with a dictionary of processing steps.

  • whiten() — whiten traces in the spectral domain.

  • time_normalize() — normalize the traces in the temporal domain.

  • assign_coordinates() — assign the geographical coordinates to the traces.

  • download_inventory() — get the inventory of each trace.

Tip

There are three main ways to create a NetworkStream object:

  1. Use the NetworkStream.read() method to read seismic waveforms files into a NetworkStream object. This method is the core method to read seismic data into the package.

  2. Use the covseisnet.stream.read() function, which is a wrapper to the first method. Note that this function is directly available from the package root.

  3. Pass any ObsPy Stream object to the NetworkStream constructor. This is the case if you have a special reader for your data.

Other standard ObsPy methods are available for instanciating a NetworkStream object, directly documented in the ObsPy documentation.

Examples

  1. Read seismic waveforms with the NetworkStream.read() method.

>>> import covseisnet as csn
>>> stream = csn.NetworkStream.read()
>>> print(type(stream))
<class 'covseisnet.stream.NetworkStream'>
>>> print(stream)
NetworkStream of 3 traces from 1 station(s) (synced):
BW.RJOB..EHZ | 2009-08-24T00:20:03.000000Z - 2009-08-24T00:20:32.990000Z | 100.0 Hz, 3000 samples
BW.RJOB..EHN | 2009-08-24T00:20:03.000000Z - 2009-08-24T00:20:32.990000Z | 100.0 Hz, 3000 samples
BW.RJOB..EHE | 2009-08-24T00:20:03.000000Z - 2009-08-24T00:20:32.990000Z | 100.0 Hz, 3000 samples
  1. Read seismic waveforms with the read() function.

>>> import covseisnet as csn
>>> stream = csn.read()
>>> print(type(stream))
<class 'covseisnet.stream.NetworkStream'>
  1. Create a NetworkStream object from an ObsPy Stream object.

>>> import obspy
>>> import covseisnet as csn
>>> stream = obspy.read()
>>> network_stream = csn.NetworkStream(stream)
>>> print(type(network_stream))
<class 'covseisnet.stream.NetworkStream'>
stats(index: int = 0, key: str | None = None) Any | Stats[source]#

Stats dictionary of one of the traces in the stream.

The purpose of this method is to extract the stats dictionary of one of the traces (by default the first trace) in the stream when the traces are synchronized. This is not to be confused with the stats() attribute available at the trace level.

Parameters:
  • index (int, optional) -- The index of the trace in the stream. Default to 0.

  • key (str, optional) -- The key of the stats dictionary. If set, the method returns the value of the key in the stats dictionary. If not set, the method returns the stats dictionary itself.

Returns:

Any or Stats -- The stats dictionary of the trace. If the key is set, the method returns the value of the key in the stats dictionary.

Example

>>> stream = csn.read()
>>> stream.stats()
         network: BW
         station: RJOB
        location:
         channel: EHZ
       starttime: 2009-08-24T00:20:03.000000Z
         endtime: 2009-08-24T00:20:32.990000Z
   sampling_rate: 100.0
           delta: 0.01
            npts: 3000
           calib: 1.0
    back_azimuth: 100.0
     inclination: 30.0
        response: Channel Response
        From M/S (Velocity in Meters Per Second) to COUNTS (Digital Counts)
        Overall Sensitivity: 2.5168e+09 defined at 0.020 Hz
        4 stages:
                Stage 1: PolesZerosResponseStage from M/S to V, gain: 1500
                Stage 2: CoefficientsTypeResponseStage from V to COUNTS, gain: 1.67785e+06
                Stage 3: FIRResponseStage from COUNTS to COUNTS, gain: 1
                Stage 4: FIRResponseStage from COUNTS to COUNTS, gain: 1
classmethod read(pathname_or_url: str | BytesIO | Path | None = None, **kwargs) NetworkStream[source]#

Read seismic waveforms files into a NetworkStream object.

This function uses the obspy.core.stream.read() function to read the streams. A detailed list of arguments and options are available in the documentation. This function opens either one or multiple waveform files given via file name or URL using the pathname_or_url attribute. The format of the waveform file will be automatically detected if not given. See the Supported Formats section in the obspy.core.stream.read() function.

This function returns an NetworkStream object which directly inherits from the obspy.core.stream.Stream object.

Parameters:
  • pathname_or_url (str or io.BytesIO or None) -- String containing a file name or a URL or a open file-like object. Wildcards are allowed for a file name. If this attribute is omitted, an example NetworkStream object will be returned.

  • **kwargs (dict, optional) -- Other parameters are passed to the obspy.core.stream.read() directly.

Returns:

NetworkStream -- The seismic waveforms.

Example

>>> import covseisnet as csn
>>> stream = csn.NetworkStream.read()
>>> print(stream)
Network Stream of 3 traces from 1 stations (synced):
BW.RJOB..EHZ | 2009-08-24T00:20:03.000000Z... | 100.0 Hz, 3000 samples
BW.RJOB..EHN | 2009-08-24T00:20:03.000000Z... | 100.0 Hz, 3000 samples
BW.RJOB..EHE | 2009-08-24T00:20:03.000000Z... | 100.0 Hz, 3000 samples

See also

read()

times(*args, **kwargs) ndarray[source]#

Common time vector.

Because the NetworkStream handles traces sampled on the same time vector, this function only returns the times of the first trace with the Trace method times() if the traces are synchronized.

Parameters:
  • *args (tuple) -- Arguments passed to the Trace method times(). For instance, passing "matplotlib" allows to recover matplotlib timestamps instead of seconds from the start of the trace (default).

  • **kwargs (dict, optional) -- Arguments passed to the method times(). For instance, passing type="matplotlib" allows to recover matplotlib timestamps instead of seconds from the start of the trace (default).

Returns:

numpy.ndarray -- The timestamps.

Raises:

AssertionError -- If the traces are not synchronized.

Tip

By default, the method returns the times in seconds since the start of the trace. In order to extract times in matplotlib format, you can set the type parameter of the times() method such as

>>> import covseisnet as csn
>>> stream = csn.read()
>>> stream.times("matplotlib")
array([14480.01392361, 14480.01392373, 14480.01392384, ...,
14480.01427049, 14480.0142706 , 14480.01427072])
time_extent() tuple[UTCDateTime, UTCDateTime][source]#

Get the minimal time extent of traces in a stream.

This function returns the minimal start and end times of the traces in the stream. This is useful when synchronizing the traces to the same time vector. The start time is defined as the maximum start time of the traces, and the end time is defined as the minimum end time of the traces.

Parameters:
  • stream (NetworkStream or)

  • :class:`~obspy.core.stream.Stream` -- The stream object.

Returns:

tuple of UTCDateTime -- The minimal start and end times of the traces.

Example

>>> import covseisnet as csn
>>> stream = csn.read()
>>> stream.time_extent()
(UTCDateTime(2009, 8, 24, 0, 20, 3), UTCDateTime(2009, 8, 24, 0, 20, 32, 990000))
cut(starttime: str | UTCDateTime, endtime: str | UTCDateTime | None = None, duration: float | None = None, **kwargs) None[source]#

Trim traces between start and end date times.

This function is a wrapper to the ObsPy trim() method, but supports string format for the start and end times, enabling a more user-friendly interface. The function uses the ObsPy UTCDateTime function in order to parse the start and end times.

Parameters:
  • starttime (str or UTCDateTime) -- The start date time.

  • endtime (str or UTCDateTime) -- The end date time.

  • duration (float, optional) -- The duration of the trace in seconds. If set, the end time is calculated as starttime + duration. This parameter is ignored if the endtime parameter is set.

  • **kwargs (dict, optional) -- Arguments passed to the trim() method.

Example

This example shows how to cut a stream between two given times. The stream is first read from the example data, and then cut between two given times.

>>> import covseisnet as csn
>>> stream = csn.read()
>>> stream.cut("2009-08-24 00:20:05", "2009-08-24 00:20:12")
>>> print(stream)
NetworkStream of 3 traces from 1 station(s) (synced):
BW.RJOB..EHZ | 2009-08-24T00:20:05.000000Z - 2009-08-24T00:20:12.000000Z | 100.0 Hz, 701 samples
BW.RJOB..EHN | 2009-08-24T00:20:05.000000Z - 2009-08-24T00:20:12.000000Z | 100.0 Hz, 701 samples
BW.RJOB..EHE | 2009-08-24T00:20:05.000000Z - 2009-08-24T00:20:12.000000Z | 100.0 Hz, 701 samples
synchronize(interpolation_method: str = 'linear', sampling_rate: float | None = None, starttime: UTCDateTime | None = None, endtime: UTCDateTime | None = None, npts: int | None = None, **kwargs) None[source]#

Synchronize seismic traces with interpolation.

This method synchronizes the seismic traces in the stream by interpolating the traces to a common time vector. The method uses the largest start time and the smallest end time of the traces to interpolate all traces to the same time vector with the ObsPy method interpolate().

Parameters:
  • method (str, default) -- Interpolation method. Default to "linear".

  • sampling_rate (float, optional) -- The sampling rate of the traces. If not set and if the traces have all the same sampling rate, the method uses the sampling rate of the first trace. If not set and if the traces have different sampling rates, the method raises a ValueError.

  • **kwargs (dict, optional) -- Additional keyword arguments passed to the interpolate() method. Check the ObsPy documentation for more details on the available options.

Raises:

ValueError -- If the traces have different sampling rates and the sampling_rate parameter is not set.

process(processing: dict[str, Any]) None[source]#

Process the seismic traces in the stream.

This method processes the seismic traces in the stream with a dictionary of processing steps. The dictionary must contain the processing steps as keys and the parameters as values. The method applies the processing steps to each trace in the stream.

Parameters:

processing (dict) -- The dictionary of processing steps.

Example

>>> stream = csn.NetworkStream.read()
>>> processing = {
...     "detrend": "linear",
...     "taper": 0.05,
...     "filter": {"type": "bandpass", "freqmin": 1, "freqmax": 10},
...     "time_normalize": {"method": "onebit"}
}
>>> stream.process(processing)
whiten(window_duration: float, smooth_length: int = 0, smooth_order: int = 1, epsilon: float = 1e-10, **kwargs: Any) None[source]#

Whiten traces in the spectral domain.

The action of whitening a seismic trace is to normalize the trace in the spectral domain. Typically, the spectrum becomes flat after whitening, resembling white noise. This strategy is often used to remove the influence of time-localized signal and diminish the site effects from a seismic station to another. Any local source is also drastically reduced thanks to the whitening process.

The following description is applied to every trace in the stream. For the sake of simplicity, we consider a single trace \(x(t)\). Note that the method is applied in every window of a short-time Fourier transform of the trace, namely \(s(t, \omega)\) before applying the inverse short-time Fourier transform to obtain the whitened seismogram \(\hat x(t)\). We here nore \(x(\omega)\) the spectrum of the trace within a given window. For more information on the short-time Fourier transform, see the ShortTimeFourierTransform class documentation.

We define the whitening process as

\[\hat x(\omega) = \frac{x(\omega)}{\mathcal{S}x(\omega) + \epsilon},\]

where \(\mathcal{S}\) is a smoothing operator applied to the spectrum \(x(\omega)\), and \(\epsilon\) is a regularization parameter to avoid division by zero. The smoothing operator is defined by the smooth_length parameter. We distinguish two cases:

  • If the smooth_length parameter is set to 0, the operator \(\mathcal{S}\) is defined as \(\mathcal{S}x(\omega) = |x(\omega)|\), and therefore

    \[\hat x(\omega) = \frac{x(\omega)}{|x(\omega)| + \epsilon} \approx e^{i\phi}.\]

    In this case, the method calls the modulus_division().

  • If the smooth_length parameter is set to a value greater than 0, the operator \(\mathcal{S}\) is defined as a Savitzky-Golay filter with parameters set by smooth_length and smooth_order. This allows to introduce less artifacts in the whitening process. In this case, the method calls the smooth_modulus_division().

Parameters:
  • window_duration (float) -- The duration of the window for the short-time Fourier transform.

  • smooth_length (int, optional) -- The length of the Savitzky-Golay filter for smoothing the spectrum. If set to 0, the spectrum is not smoothed (default).

  • smooth_order (int, optional) -- The order of the Savitzky-Golay filter for smoothing the spectrum. This parameter is only used if smooth_length is greater than 0.

  • epsilon (float, optional) -- Regularization parameter in division, set to 1e-10 by default.

  • **kwargs (dict, optional) -- Additional keyword arguments passed to the covseisnet ShortTimeFourierTransform() class constructor.

time_normalize(method: str = 'onebit', smooth_length: int = 11, smooth_order: int = 1, epsilon: float = 1e-10) None[source]#

Normalize the seismic traces in temporal domain.

Considering the seismic trace \(x(t)\), the normalized trace \(\hat x(t)\) is obtained with

\[\hat x(t) = \frac{x(t)}{\mathcal{A}x(t) + \epsilon}\]

where \(A\) is an operator applied to the trace \(x(t)\), and \(\epsilon > 0\) is a regularization value to avoid division by 0. The operator \(\mathcal{A}\) is defined by the method parameter. We distinguish two cases:

  • If the method parameter is set to "onebit", the operator \(\mathcal{A}\) is defined as \(\mathcal{A}x(t) = |x(t)|\), and therefore

    \[\hat x(t) = \frac{x(t)}{|x(t)| + \epsilon} \approx \text{sign}(x(t)).\]

    In this case, the method calls the modulus_division().

  • If the method parameter is set to "smooth", the operator \(\mathcal{A}\) is defined as a Savitzky-Golay filter applied to the Hilbert envelope of the trace. The Savitzky-Golay filter is defined by the smooth_length and smooth_order parameters. This allows to introduce less artifacts in the normalization process. In this case, the method calls the smooth_envelope_division().

Parameters:
  • method (str, optional) -- Must be one of "onebit" (default) or "smooth".

    • "onebit": compress the seismic trace into a series of -1 and 1.

    • "smooth": normalize each trace by a smooth version of its envelope.

  • smooth_length (int, optional) -- If the method keyword argument is set to "smooth", the normalization is performed with the smoothed trace envelopes, calculated over a sliding window of smooth_length samples.

  • smooth_order (int, optional) -- If the method keyword argument is set to "smooth", the normalization is performed with the smoothed trace envelopes. The smoothing order is set by the smooth_order parameter.

  • epsilon (float, optional) -- Regularization parameter in division, set to 1e-10 by default.

assign_coordinates(inventory: str | Inventory) None[source]#

Assign the geographical coordinates to the traces.

This method assigns the geographical coordinates of the stations to the traces in the stream from Inventory or the inventory file containing the inventory. The method uses the ObsPy read_inventory() method to extract the coordinates of each trace in the stream if a filepath is provided.

The method adds a coordinate dictionary to each trace in the stream. The coordinate dictionary contains the keys latitude, longitude, and elevation.

Parameters:

inventory (str or Inventory) -- The inventory file or the inventory object containing the station coordinates.

Example

>>> import covseisnet as csn
>>> stream = csn.read()
>>> stream.assign_coordinates("inventory.xml")

See also

read_inventory()

download_inventory(datacenter: str = 'IRIS') Inventory[source]#

Get the inventory of each trace in the stream.

This method uses the ObsPy get_stations() method to extract the inventory of each trace in the stream.

Parameters:

datacenter (str, optional) -- The datacenter to use for retrieving the station coordinates. The default is "IRIS".

Example

>>> import covseisnet as csn
>>> stream = csn.read()
>>> stream.download_inventory(datacenter="LMU")
    Inventory created at 2024-08-20T09:10:33.403152Z
        Created by: ObsPy 1.4.1
                    https://www.obspy.org
        Sending institution: ObsPy 1.4.1,SeisComP (LMU)
        Contains:
                Networks (3):
                        BW (3x)
                Stations (3):
                        BW.RJOB (Jochberg, Bavaria, BW-Net) (3x)
                Channels (0):
property synced: bool#

Check if traces are sampled on the same time vector.

This method checks if all traces are sampled on the same time vector. This is useful to ensure that the traces are synchronized before performing any operation that requires the traces to be sampled on the same time vector.

Returns:

bool -- True if all traces are sampled on the same time vector, False otherwise.

property equal_rates: bool#

Check if all traces have the same sampling rate.

This method checks if all traces have the same sampling rate. This is useful to ensure that the traces are synchronized before performing any operation that requires the traces to be sampled on the same time vector.

Returns:

bool -- True if all traces have the same sampling rate, False otherwise.

property equal_length: bool#

Check if all traces have the same number of samples.

This method checks if all traces have the same number of samples. This is useful to ensure that the traces are synchronized before performing any operation that requires the traces to be sampled on the same time vector.

Returns:

bool -- True if all traces have the same number of samples, False otherwise.

property sampling_rate: float#

Common sampling rate in Hertz.

This property asserts that all traces have the same sampling rate and returns the sampling rate of the first trace via the NetworkStream stats() custom method. This is equivalent to getting stream[0].stats.sampling_rate directly, minus the assertion.

Example

>>> stream = csn.read()
>>> stream.sampling_rate
100.0
property npts: int#

Common number of samples.

This property asserts that all traces have the same number of samples and returns the number of samples of the first trace via the NetworkStream stats() custom method. This is equivalent to getting stream[0].stats.npts directly, minus the assertion.

Example

>>> stream = csn.read()
>>> stream.npts
3000
property coordinates: list[tuple[float, float, float]]#

Geographical coordinates of the traces.

This property is also available directly from looping over the traces and accessing the stats.coordinates attribute of each trace.

covseisnet.stream.read(pathname_or_url: str | BytesIO | Path | None = None, **kwargs: dict) NetworkStream[source]#

Read seismic waveforms files into an NetworkStream object.

This function uses the obspy.core.stream.read() function to read the streams. A detailed list of arguments and options are available in the documentation. This function opens either one or multiple waveform files given via file name or URL using the pathname_or_url attribute. The format of the waveform file will be automatically detected if not given. See the Supported Formats section in the obspy.core.stream.read() function.

This function returns an NetworkStream object which directly inherits from the obspy.core.stream.Stream object.

Parameters:
  • pathname_or_url (str or io.BytesIO or None) -- String containing a file name or a URL or a open file-like object. Wildcards are allowed for a file name. If this attribute is omitted, an example NetworkStream object will be returned.

  • **kwargs (dict, optional) -- Other parameters are passed to the obspy.core.stream.read() directly.

Returns:

NetworkStream -- The seismic waveforms.

Example

For a quick start you may omit all arguments and ObsPy will load a basic example seismogram. Further usages of this function can be seen in the ObsPy documentation.

>>> import covseisnet as csn
>>> stream = csn.read()
>>> print(stream)
Network Stream of 3 traces from 1 stations (synced):
BW.RJOB..EHZ | 2009-08-24T00:20:03.000000Z... | 100.0 Hz, 3000 samples
BW.RJOB..EHN | 2009-08-24T00:20:03.000000Z... | 100.0 Hz, 3000 samples
BW.RJOB..EHE | 2009-08-24T00:20:03.000000Z... | 100.0 Hz, 3000 samples

See also

read()