seismic package

Subpackages

Submodules

seismic.analyze_station_orientations module

seismic.extract_event_traces module

seismic.model_properties module

Helper classes to encapsulate model properties.

class seismic.model_properties.LayerProps(vp, vs, rho, thickness)[source]

Bases: object

Helper class to contain layer bulk material properties

Constructor for given properties

Parameters:
  • vp (float) – P-wave body wave velocity

  • vs (float) – S-wave body wave velocity

  • rho (float) – Bulk material density

  • thickness (float) – 1D (vertical) thickness of the layer.

property H

Get layer thickness

property Vp

Get P-wave body wave velocity

property Vs

Get S-wave body wave velocity

property rho

Get bulk material density

seismic.network_event_dataset module

Class encapsulating a collection of event waveforms for stations of one network.

class seismic.network_event_dataset.NetworkEventDataset(stream_src, network=None, station=None, location='', ordering='ZNE', root='/waveforms')[source]

Bases: object

Collection of 3-channel ZNE streams with traces aligned to a fixed time window about seismic P-wave arrival events, for a given network.

Two indexes are provided. One indexes hierarchically by station code and event ID, yielding a 3-channel ZNE stream per event, so that you can easily gather all traces for a given station by iterating over events.

The other index indexes hierarchically by event ID and station code, yielding a 3-channel ZNE stream per station. Using this index you can easily gather all traces for a given event across multiple stations.

Preferably each input trace will already have an ‘event_id’ attribute in its stats. If not, an event ID will be invented based on station identifiers and time window.

Initialize from data source (file or obspy.Stream). Traces are COPIED into the dataset in order to leave input object intact, since many obspy functions mutate traces in-place.

All streams in the input data source stream_src are expected to belong to the same network. This is checked as the data is ingested. A discrepant network code is an error condition.

Parameters:
  • stream_src (str, pathlib.Path or obspy.Stream) – Source of input streams. May be a file name or an Obspy Stream

  • network (str) – Network code of streams to load. If stream_src is an Obspy Stream, the streams will be filtered to match this network code.

  • station (str) – Station code of streams to load. If stream_src is an Obspy Stream, the streams will be filtered to match this station code.

  • location (str) – [OPTIONAL] Location code of streams to load. Leave as default (empty string) if location code is empty in the data source.

  • ordering (str) – Channel ordering to be applied to the data after loading. The channel labelling must be consistent with the requested ordering - rotation to the coordinate system implied by the ordering is NOT applied.

Raises:

AssertionError – If discrepant network code is found in input data

apply(_callable)[source]

Apply a callable across all streams. Use to apply uniform processing steps to the whole dataset.

Parameters:

_callable (Any Callable compatible with the call signature.) – Callable object that takes an obspy Stream as input and applies itself to that Stream. Expect that stream may be mutated in-place by the callable.

Returns:

None

by_event()[source]

Iterate over event sub-dictionaries.

Returns:

Iterable over the discrete events, each element consisting of pair containing (event id, station dict).

Return type:

Iterable(tuple)

by_station()[source]

Iterate over station sub-dictionaries.

Returns:

Iterable over the stations, each element consisting of pair containing (station code, event dict).

Return type:

Iterable(tuple)

curate(curator)[source]

Curate the dataset according to a callable curator. Modifies collection in-place to remove streams that do not satisfy the curation criteria of the callable. Curator call signature must be consitent with:

callable(station_code, event_id, stream) -> bool

The callable returns a boolean indicating whether to keep the Stream or not.

Parameters:

curator (Callable) – Function or callable delegate to adjudicate whether to keep each given stream.

Returns:

None

event(event_id)[source]

Accessor for stations for a given event.

Parameters:

event_id (str) – ID of event to look up

Returns:

Station index for given event, if event ID is found, otherwise None

Return type:

SortedDict or NoneType

num_events()[source]

Get number of events in the dataset.

Returns:

Number of events

Return type:

int

num_stations()[source]

Get number of stations in the dataset.

Returns:

Number of stations

Return type:

int

prune(items, cull=True)[source]

Remove a given sequence of (station, event) pairs from the dataset.

Parameters:
  • items (Iterable(tuple)) – Iterable of (station, event) pairs

  • cull (boolean) – If True, then empty entries in the top level index will be removed.

Returns:

None

station(station_code)[source]

Accessor for events for a given station.

Parameters:

station_code (str) – Station to get

Returns:

Event index for station, if station is found

Return type:

SortedDict

write(output_h5_filename, index_format='event')[source]

Write event dataset back out to HDF5 file.

Parameters:
  • output_h5_filename (str or path) – Output file name

  • index_format (str) – Format to use for index. Must be ‘event’ (default) or ‘standard’ (obspy default)

Returns:

True if file was written

Return type:

boolean

seismic.plot_network_event_dataset module

Bulk plotting helper functions based on NetworkEventDataset

seismic.plot_network_event_dataset.plot_ned_seismograms(ned, output_file, channel_order='ZNE')[source]

Plot seismograms in NetworkEventDataset to PDF file. If dataset is very large, this may take a long time to run.

Parameters:
  • ned – NetworkEventDataset containing waveforms.

  • output_file – Output file name

seismic.stream_io module

Helper functions for seismic stream IO.

seismic.stream_io.get_obspyh5_index(src_file, seeds_only=False, root='/waveforms')[source]

Scrape the index (only) from an obspyh5 file.

Parameters:
  • src_file (str or pathlib.Path) – Name of file to extract index from

  • seeds_only (bool) – If True, only get the seed IDs of the traces. Otherwise (default), get full index.

Returns:

Sorted dictionary with index of waveforms in the file

Return type:

sortedcontainers.SortedDict

seismic.stream_io.iter_h5_stream(src_file, headonly=False, root='/waveforms')[source]

Iterate over hdf5 file containing streams in obspyh5 format.

Parameters:
  • src_file (str or pathlib.Path) – Path to file to read

  • headonly (bool) – Only read trace stats, do not read actual time series data

Yield:

obspy.Stream containing traces for a single seismic event.

seismic.stream_io.read_h5_stream(src_file, network=None, station=None, loc='', root='/waveforms')[source]

Helper function to load stream data from hdf5 file saved by obspyh5 HDF5 file IO. Typically the source file is generated using extract_event_traces.py script. For faster loading time, a particular network and station may be specified.

Parameters:
  • src_file (str or Path) – File from which to load data

  • network (str, optional) – Specific network to load, defaults to None

  • station (str, optional) – Specific station to load, defaults to None

  • root (str, optional) – Root path in hdf5 file where to start looking for data, defaults to ‘/waveforms’

Returns:

All the loaded data in an obspy Stream.

Return type:

obspy.Stream

seismic.stream_io.remove_group(hdf_fn, net_sta_loc, root='/waveforms', logger=None)[source]
seismic.stream_io.sac2hdf5(src_folder, basenames, channels, dest_h5_file, tt_model_id='iasp91')[source]

Convert collection of SAC files from a folder into a single HDF5 stream file.

Parameters:
  • src_folder (str or Path) – Path to folder containing SAC files

  • basenames (list of str) – List of base filenames (file name excluding extension) to load.

  • channels (List of str) – List of channels to load. For each base filename in basenames, there is expected to be a file with each channel as the filename extension.

  • dest_h5_file (str or Path) – Path to output file. Will be created, or overwritten if already exists.

  • tt_model_id (str) – Which travel time earth model to use for synthesizing trace metadata. Must be known to obspy.taup.TauPyModel

Returns:

None

seismic.stream_io.safe_h5_root(src_file, root)[source]
seismic.stream_io.safe_iter_event_data(events, inventory, get_waveforms, use_rfstats=True, phase='P', request_window=None, pad=10, pbar=None, **kwargs)[source]

Return iterator yielding three component streams per station and event.

Parameters:
  • events – list of events or ~obspy.core.event.Catalog instance

  • inventory~obspy.core.inventory.inventory.Inventory instance with station and channel information

  • get_waveforms – Function returning the data. It has to take the arguments network, station, location, channel, starttime, endtime.

  • phase – Considered phase, e.g. ‘P’, ‘S’, ‘PP’

  • request_window (tuple (start, end)) – requested time window around the onset of the phase

  • pad (float) – add specified time in seconds to request window and trim afterwards again

  • pbartqdm instance for displaying a progressbar

  • kwargs – all other kwargs are passed to ~rf.rfstream.rfstats()

Returns:

three component streams with raw data

Example usage with progressbar:

from tqdm import tqdm
from rf.util import iter_event_data
with tqdm() as t:
    for stream3c in iter_event_data(*args, pbar=t):
        do_something(stream3c)
seismic.stream_io.write_h5_event_stream(dest_h5_file, stream, index='waveforms/{network}.{station}.{location}/{event_time.datetime:%Y-%m-%dT%H:%M:%S}/{channel}_{starttime.datetime:%Y-%m-%dT%H:%M:%S}_{endtime.datetime:%Y-%m-%dT%H:%M:%S}', mode='a', ignore=())[source]

Write stream to HDF5 file in event indexed format using obspy.

Parameters:
  • dest_h5_file (str or pathlib.Path) – File in which to write the stream.

  • stream (obspy.Stream) – The stream to write

  • mode (str) – Write mode, such as ‘w’ or ‘a’. Use ‘a’ to iteratively write multiple streams to one file.

  • ignore (Any iterable of str) – List of headers to ignore when writing attributes to group. Passed on directly to obspyh5.writeh5

seismic.stream_processing module

Utility stream processing functions.

seismic.stream_processing.assert_homogenous_stream(stream, funcname)[source]

Verify that the given stream does not contain mixture of stations or channels/components.

Parameters:

stream (obspy.Stream or rf.RFStream) – Stream containing one or more traces

Returns:

None

seismic.stream_processing.back_azimuth_filter(back_azi, back_azi_range)[source]

Check if back azimuth back_azi is within range. Inputs must be in the range [0, 360] degrees.

Parameters:
  • back_azi (int or float) – Value to check

  • back_azi_range (List or array of 2 floats, min and max back azimuth) – Pair of angles in degrees.

Returns:

True if back_azi is within back_azi_range, False otherwise.

Return type:

bool

seismic.stream_processing.correct_back_azimuth(_event_id, stream, baz_correction)[source]

Apply modification to the back azimuth value in the stream stats

Parameters:
  • _event_id – Ignored

  • stream (obspy.Stream or rf.RFSTream) – Stream to which correction is applied

  • baz_correction – Any object with a registered scalarize function for generating an angle correction for a trace in degrees. E.g. could be a numeric value, a dictionary of correction values, or a file produced by script rf_station_orientations.py

Returns:

Stream with modified back azimuth

Return type:

Same as type(stream)

seismic.stream_processing.negate_channel(_event_id, stream, channel)[source]

Negate the data in the given channel of the stream

Parameters:
  • _event_id – Ignored

  • stream (obspy.Stream) – Stream containing channel to flip

  • channel (str) – Single character string indicating which component to flip

Returns:

Stream with channel data negated

Return type:

obspy.Stream

seismic.stream_processing.scalarize(_obj, _stats)[source]
seismic.stream_processing.scalarize(val: Number, _stats)
seismic.stream_processing.scalarize(d: dict, stats)
seismic.stream_processing.scalarize(filename: str, stats)

Fallback scalarize function for non-specialized type

seismic.stream_processing.sinc_resampling(t, y, t_new)[source]

Resample signal y for known times t onto new times t_new. Sampling rates do not need to match and time windows do not need to overlap. t_new should not have a lower sampling rate than t.

Parameters:
  • t (numpy.array) – 1D array of times

  • y (numpy.array) – 1D array of sample values

  • t_new (numpy.array) – 1D array of new times to interpolate onto

Returns:

1D array of new interpolated sample values

Return type:

numpy.array

seismic.stream_processing.swap_ne_channels(_event_id, stream)[source]

Swap N and E channels on a stream. Changes the input stream.

Parameters:
  • _event_id – Ignored

  • stream (obspy.Stream) – Stream whose N and E channels are to be swapped

Returns:

Stream with channel swapping applied

Return type:

obspy.Stream

seismic.stream_processing.zerophase_resample(item, resample_hz)[source]
seismic.stream_processing.zne_order(tr)[source]

Channel ordering sort key function for ZNE ordering

Parameters:

tr (obspy.Trace or rf.RFTrace) – Trace whose ordinal is to be determined.

Returns:

Numeric index indicating ZNE sort order of traces in a stream

Return type:

int

seismic.stream_processing.zrt_order(tr)[source]

Channel ordering sort key function for ZRT ordering

Parameters:

tr (obspy.Trace or rf.RFTrace) – Trace whose ordinal is to be determined.

Returns:

Numeric index indicating ZRT sort order of traces in a stream

Return type:

int

seismic.stream_quality_filter module

Helper functions for curating and quality controlling stream objects.

seismic.stream_quality_filter.curate_seismograms(data_all, curation_opts, logger, rotate_to_zrt=True)[source]

Curation function to remove bad data from streams. Note that this function will modify the input dataset during curation.

Parameters:
Returns:

None, curation operates directly on data_all

seismic.stream_quality_filter.curate_stream3c(ev_id, stream3c, logger=None)[source]

Apply quality curation criteria to a stream. Modifies the stream in-place if required. Traces in stream must be in ZNE order. Each trace in the stream is expected to have metadata with starttime, endtime, channel and inclination.

The following checks are made. If any of these checks fails, the function returns False: * Inclination value is not NaN * The stream has 3 channels * Each trace in the stream has the same number of samples * None of the traces have any NaN values in the time series data * None of the traces have zero variance

The following cleanups are attempted on the stream: * All 3 traces have the same time range

Parameters:
  • ev_id (int or str) – The event id

  • stream3c (obspy.Stream) – Stream with 3 components of trace data

  • logger (logging.Logger) – Logger in which to log messages

Returns:

True if checks pass, False otherwise

Return type:

bool

seismic.units_utils module

Constants and utility functions for unit conversion.