seismic.receiver_fn package

Submodules

seismic.receiver_fn.bulk_rf_report module

seismic.receiver_fn.generate_rf module

Generate Receiver Functions (RF) from collection of 3-channel seismic traces.

seismic.receiver_fn.generate_rf.event_waveforms_to_rf(input_file, output_file, config, network_list='*', station_list='*', only_corrections=False)[source]

Main entry point for generating RFs from event traces.

Config file consists of 3 sub-dictionaries. One named “filtering” for input stream filtering settings, one named “processing” for RF processing settings, and one named “correction” for rotating/swapping/negating channel data for one or more named stations with potential orientation discrepancies. Each of these sub-dicts is described below:

"filtering":  # Filtering settings
{
  "resample_rate": float # Resampling rate in Hz
  "taper_limit": float   # Fraction of signal to taper at end, between 0 and 0.5
  "filter_band": (float, float) # Filter pass band (Hz). Not required for freq-domain deconvolution.
  "channel_pattern": # Ordered list of preferred channels, e.g. 'HH*,BH*',
                     # where channel selection is ambiguous.
  "baz_range": (float, float) or [(float, float), ...] # Discrete ranges of source back azimuth to use (degrees).
      # Each value must be between 0 and 360. May be a pair or a list of pairs for multiple ranges.
}

"processing":  # RF processing settings
{
  "custom_preproc":
  {
    "import": 'import custom symbols',  # statement to import required symbols
    "func": 'preproc functor'  # expression to get handle to custom preprocessing functor
    "args": {}  # additional kwargs to pass to func
  }
  "trim_start_time": float # Trace trim start time in sec, relative to onset
  "trim_end_time": float # Trace trim end time in sec, relative to onset
  "rotation_type": str # Choice of ['zrt', 'lqt']. Rotational coordinate system
                       # for aligning ZNE trace components with incident wave direction
  "deconv_domain": str # Choice of ['time', 'freq', 'iter']. Whether to perform deconvolution
                       # in time or freq domain, or iterative technique
  "gauss_width": float # Gaussian freq domain filter width. Only required for freq-domain deconvolution
  "water_level": float # Water-level for freq domain spectrum. Only required for freq-domain deconvolution
  "spiking": float # Spiking factor (noise suppression), only required for time-domain deconvolution
  "normalize": bool # Whether to normalize RF amplitude
}

"correction": # corrections to be applied to data for named stations prior to RF computation
{
  "plot_dir": str # path to folder where plots related to orientation corrections are to be saved
  "swap_ne": list # list of NET.STA.LOC for which N and E channels are to be swapped, e.g ["OA.BL27."],
  "rotate": list # list of NET.STA.LOC that are to be rotated to maximize P-arrival energy on                            the primary RF component, e.g ["OA.BL27."]
  "negate": list # list of NET.STA.LOC.CHA that are to be negated, e.g ["OA.BL27..HHZ"]
}
Parameters:
  • input_file (str or pathlib.Path) – Event waveform source file for seismograms, generated using extract_event_traces.py script

  • output_file (str or pathlib.Path) – Name of hdf5 file to produce containing RFs

  • config (dict) – Dictionary containing job configuration parameters

Returns:

None

seismic.receiver_fn.plot_ccp module

seismic.receiver_fn.plot_ccp_batch module

seismic.receiver_fn.plot_spatial_map module

Use cartopy to plot moho grid and gradient onto a map.

seismic.receiver_fn.plot_spatial_map.from_params(params)[source]

Create plots from WorkflowParameters as part of moho workflow.

seismic.receiver_fn.plot_spatial_map.plot_spatial_map(grid_data, gradient_data, methods_datasets=None, projection_code=None, title=None, feature_label=None, bounds=None, scale=None)[source]

Make spatial plot of point dataset with filled contours overlaid on map.

Parameters:
  • point_dataset – Name of point dataset file. Should be in format produced by script pointsets2grid.py

  • projection_code – EPSG projection code, e.g. 3577 for Australia

  • title – Title string for top of plot

  • feature_label – Label for the color bar of the plotted feature

  • bounds – 4 element tuple of (L, R, B, T) for limiting plot extent

  • scale – 2 element tuple of (vmin, vmax) for limiting colormap scale

Returns:

seismic.receiver_fn.pointsets2grid module

Blend data from multiple methods together onto a common grid.

Multiple methods with per-method settings are passed in using a JSON configuration file.

See config_moho_workflow_example.json and the wiki page at https://github.com/GeoscienceAustralia/hiperseis/wiki/Blending-and-Plotting-Point-Datasets for an explanation of the config.

Bounds are optional. If provided, the interpolation grid will be limited to this extent. If not provided, the interpolation grid is bounded by the maximum extent of the aggregate datasets.

The gridded data will be written to output directory as ‘moho_grid.csv’. If output directory is not provided, the current working directory will be used.

The gridded gradient will be written to output directory as ‘moho_gradient.csv’.

Reference: B. L. N. Kennett 2019, “Areal parameter estimates from multiple datasets”, Proc. R. Soc. A. 475:20190352, http://dx.doi.org/10.1098/rspa.2019.0352

Requires:

  • pyepsg

  • cartopy

seismic.receiver_fn.pointsets2grid.make_grid(params)[source]

Run multi point dataset weighted averaging over Gaussian interpolation functions to produce aggregate dataset. Source data coordinates are assumed to be in lon/lat (WGS84). First 2 rows of output file contain grid dimensions, followed by CSV data.

seismic.receiver_fn.rf_3dmigrate module

seismic.receiver_fn.rf_deconvolution module

Custom receiver function deconvolution methods and support functions.

seismic.receiver_fn.rf_deconvolution.iter_deconv_pulsetrain(numerator, denominator, sampling_rate, time_shift, max_pulses=1000, tol=0.001, gwidth=2.5, only_positive=False, log=None)[source]

Iterative deconvolution of source and response signal to generate seismic receiver function. Adapted to Python by Andrew Medlin, Geoscience Australia (2019), from Chuck Ammon’s (Saint Louis University) iterdeconfd Fortran code, version 1.04.

Note this is not really a frequency-domain deconvolution method, since the deconvolution is based on generating a time-domain pulse train which is filtered and convolved with source signal in time domain in order to try to replicate observation. Results should be the same even if some of the spectral techniques used in functions (such as gauss_filter) were replaced by non-spectral equivalents.

Parameters:
  • numerator (numpy.array(float)) – The observed response signal (e.g. R, Q or T component)

  • denominator (numpy.array(float)) – The source signal (e.g. L or Z component)

  • sampling_rate (float) – The sampling rate in Hz of the numerator and denominator signals

  • time_shift (float) – Time shift (sec) from start of input signals until expected P-wave arrival onset.

  • max_pulses (int) – Maximum number of delta function pulses to synthesize in the unfiltered RF (up to 1000)

  • tol (float) – Convergence tolerance, iteration stops if change in error falls below this value

  • gwidth (float) – Gaussian filter width in the normalized frequency domain.

  • only_positive (bool) – If true, only use positive pulses in the RF.

  • log (logging.Logger) – Log instance to log messages to.

Returns:

RF trace, pulse train, expected response signal, predicted response signal, quality of fit statistic

Return type:

numpy.array(float), numpy.array(float), numpy.array(float), numpy.array(float), float

seismic.receiver_fn.rf_deconvolution.rf_iter_deconv(response_data, source_data, sr, tshift, min_fit_threshold=80.0, normalize=0, **kwargs)[source]

Adapter function to rf library. To use, add arguments deconvolve=’func’, func=rf_iter_deconv to rf.RFStream.rf() function call.

Parameters:
  • response_data (list of numpy.array(float)) – List of response signals for which to compute receiver function

  • source_data (numpy.array(float)) – Source signal to use for computing the receiver functions

  • sampling_rate (float) – Sampling rate of the input signals (Hz)

  • time_shift (float) – Time shift (seconds) from start of signal to onset

  • min_fit_threshold (float) – Minimum percentage of fit to include trace in results, otherwise will be returned as empty numpy array.

  • normalize (int or None) – Component of stream to use for normalization, usually component 0 (the vertical component). Set to None to disable normalization.

Returns:

Receiver functions corresponding to the list of input response signals.

Return type:

list of numpy.array(float)

seismic.receiver_fn.rf_h5_file_event_iterator module

Helper class to iterate over 3-channel event traces in h5 file generated by rf library, without loading all the traces into memory. This is a scalable solution for very large files.

class seismic.receiver_fn.rf_h5_file_event_iterator.IterRfH5FileEvents(h5_filename, memmap=False, channel_pattern=None)[source]

Bases: object

Helper class to iterate over events in h5 file generated by extract_event_traces.py and pass them to RF generator. This class avoids having to load the whole file up front via obspy which is slow and not scalable.

Yields (station_id, event_id, event_time, stream)

Due to the expected hierarchy structure of the input H5 file, yielded event traces are grouped by station ID.

Data yielded per event can easily be hundreds of kB in size, depending on the length of the event traces. rf library defaults to window of (-50, 150) sec about the P wave arrival time.

Initializer

Parameters:
  • h5_filename (str or pathlib.Path) – Name of file containing event seismograms in HDF5 format, indexed in seismic.stream_io.EVENTIO_H5INDEX format

  • memmap (bool) – If True, memmap the open file. Can improve tractability of handling very large files.

  • channel_pattern (str) – [OPTIONAL] Ordered list of preferred channels, e.g. ‘HH*,BH*’. Channel mask to use to select channels returned.

seismic.receiver_fn.rf_h5_file_station_iterator module

Helper class to iterate over station events in h5 file generated by rf library, without loading all the traces into memory. This is a scalable solution for very large files.

class seismic.receiver_fn.rf_h5_file_station_iterator.IterRfH5StationEvents(h5_filename, memmap=False)[source]

Bases: object

Helper class to iterate over stations in h5 file generated by extract_event_traces.py and pass them to RF generator. This class avoids having to load the whole file up front via obspy which is slow and not scalable.

This class yields 3-channel traces per station per event.

Data yielded per station can easily be many MB in size.

Initializer

Parameters:
  • h5_filename (str or pathlib.Path) – Name of file containing event seismograms in HDF5 format, indexed in seismic.stream_io.EVENTIO_H5INDEX format

  • memmap (bool) – If True, memmap the open file. Can improve tractability of handling very large files.

seismic.receiver_fn.rf_handpick_tool module

seismic.receiver_fn.rf_inversion_export module

seismic.receiver_fn.rf_network_dict module

Class encapsulating a collection of receiver functions for stations of one network.

class seismic.receiver_fn.rf_network_dict.NetworkRFDict(rf_stream)[source]

Bases: object

Collection of RFs for a given network indexed by station code, channel code. Note that location codes are not taken into account.

Initialize from rf.RFStream

Parameters:

rf_stream (rf.RFStream) – RFStream data

keys()[source]

Accessor for the top level keys (station codes) of the network in an iterable container.

Returns:

Iterable of top level keys to the dictionary

Return type:

Python iterable

seismic.receiver_fn.rf_plot_utils module

seismic.receiver_fn.rf_plot_vespagram module

Plot vespagrams from receiver function data

seismic.receiver_fn.rf_process_io module

Helper for asynchronous writing of RF data to file.

seismic.receiver_fn.rf_process_io.async_write(rfstream_queue, outfile_name, max_buffered=100, metadata='')[source]

Monitors asynchronous queue for data, removes from queue to buffer, then flushes buffer intermittently and when queue termination signal is put.

When None is received on the queue, this is taken as the signal to terminate monitoring the queue.

Parameters:
  • rfstream_queue (multiprocessing.Manager.Queue) – Queue into which RFStreams are pushed for writing to file.

  • outfile_name (str or Path) – Name of file into which queued RFStream results are periodically written.

  • max_buffered (int, optional) – Maximum number of RFStreams to buffer before flushing to file, defaults to 100

  • metadata (str, optional) – Metadata string to write to root attribute

seismic.receiver_fn.rf_quality_filter module

Filter out invalid RFs, and compute a variety of quality metrics for the remaining RFs. These metrics can be combined in various ways downstream to perform different kinds of filtering. Quality metrics are stored in the stats of each trace.

seismic.receiver_fn.rf_quality_filter.compute_max_coherence(orig_stream, f1, f2)[source]

Finding coherence between two signals in frequency domain f1 and f2 - normalised min and max frequencies, f1 < f2 <= ~0.5 returns array of indexes for coherent traces with median

Suggested minimum level of coherence for good results: 0.6

seismic.receiver_fn.rf_quality_filter.compute_rf_quality_metrics(station_id, station_stream3c, similarity_eps)[source]

Top level function for adding quality metrics to trace metadata.

Parameters:
  • station_id (str) – Station ID

  • station_stream3c (list(rf.RFStream) with 3 components) – 3-channel stream

  • similarity_eps (float) – Distance threshold used for DBSCAN clustering

Returns:

Triplet of RF streams with Z, R or Q, and T components with populated quality metrics. Otherwise return None in case of failure.

seismic.receiver_fn.rf_quality_filter.get_rf_stream_components(stream)[source]

Identify the RF component types and return them.

Parameters:

stream (rf.RFStream) – Stream containing mixed RF components.

Returns:

(RF component type, primary RF component (R or Q), transverse RF component (T), source component (Z or L))

Return type:

(str, rf.RFStream, rf.RFStream, rf.RFStream)

seismic.receiver_fn.rf_quality_filter.rf_group_by_similarity(swipe, similarity_eps)[source]

Cluster waveforms by similarity

Parameters:
  • swipe (numpy.array) – Numpy array of RF rowwise

  • similarity_eps (float) – Tolerance on similarity between traced to be considered in the same group.

Returns:

Index of the group for each trace. -1 if no group is found for a given trace.

Return type:

numpy.array

seismic.receiver_fn.rf_quality_filter.rf_quality_metrics_queue(oqueue, station_id, station_stream3c, similarity_eps, drop_z=True)[source]

Produce RF quality metrics in a stream and queue the QC’d components for downstream processing.

Parameters:
  • oqueue (queue or multiprocessing.Manager.Queue) – Output queue where filtered streams are queued

  • station_id (str) – Station ID

  • station_stream3c (list(rf.RFStream) with 3 components) – 3-channel stream

  • similarity_eps (float) – Distance threshold used for DBSCAN clustering

seismic.receiver_fn.rf_quality_filter.spectral_entropy(stream)[source]

Compute the spectral entropy of a trace

Parameters:

trace (rf.RFTrace) – Single channel seismic trace

Returns:

Spectral entropy of the trace waveform

Return type:

float

seismic.receiver_fn.rf_stacking module

seismic.receiver_fn.rf_stacking.compute_hk_stack(cha_data, Vp=6.5, h_range=None, k_range=None, weights=array([0.5, 0.4, 0.1]), root_order=1, semblance_weighted=True)[source]

Compute H-k stacking array on a dataset of receiver functions.

Parameters:
  • cha_data (Iterable(rf.RFTrace)) – List or iterable of RF traces to use for H-k stacking.

  • Vp (float, optional) – Average crustal Vp for computing H-k stacks

  • h_range (numpy.array [1D], optional) – Range of h values (Moho depth) values to cover, defaults to np.linspace(20.0, 70.0, 251)

  • k_range (numpy.array [1D], optional) – Range of k values to cover, defaults to np.linspace(1.4, 2.0, 301)

  • weights (numpy.array [1D], optional) – numpy array of length 3, containing weights for the three phases (Ps, PpPs and (PpSs + PsPs))

  • root_order (int, optional) – Exponent for nth root stacking as per K.J.Muirhead (1968), defaults to 1

Returns:

k-grid values [2D], h-grid values [2D], H-k stack values [2D]

Return type:

numpy.array [2D], numpy.array [2D], numpy.array [2D]

seismic.receiver_fn.rf_stacking.compute_sediment_hk_stack(cha_data, H_c, k_c, Vp=6.5, h_range=None, k_range=None, root_order=9)[source]

Compute H-k stacking array on a dataset of receiver functions.

Parameters:
  • cha_data (Iterable(rf.RFTrace)) – List or iterable of RF traces to use for H-k stacking.

  • H_c (float, optional) – Crustal thickness estimate from H-k stack

  • k_c (float, optional) – Crustal Vp/Vs ratio estimate from H-k stack

  • Vp (float, optional) – Average crustal Vp for computing H-k stacks

  • h_range (numpy.array [1D], optional) – Range of h values (Moho depth) values to cover, defaults to np.linspace(20.0, 70.0, 251)

  • k_range (numpy.array [1D], optional) – Range of k values to cover, defaults to np.linspace(1.4, 2.0, 301)

  • root_order (int, optional) – Exponent for nth root stacking as per K.J.Muirhead (1968), defaults to 1

Returns:

k-grid values [2D], h-grid values [2D], H-k stack values [2D]

Return type:

numpy.array [2D], numpy.array [2D], numpy.array [2D]

seismic.receiver_fn.rf_stacking.find_global_hk_maximum(k_grid, h_grid, hk_weighted_stack)[source]

Given the weighted stack computed from function compute_weighted_stack and the corresponding k-grid and h-grid, find the location in H-k space of the global maximum.

Parameters:
  • k_grid (Two-dimensional numpy.array) – Grid of k-values

  • h_grid (Two-dimensional numpy.array) – Grid of H-values

  • hk_weighted_stack (Two-dimensional numpy.array) – Grid of stacked RF sample values produced by function rf_stacking.computed_weighted_stack()

Returns:

Location of global maximum on the H-k grid of the maximum stack value.

Return type:

tuple(float, float)

seismic.receiver_fn.rf_stacking.find_local_hk_maxima(k_grid, h_grid, hk_stack_sum, max_number=3)[source]

seismic.receiver_fn.rf_synthetic module

Helper functions for producing synthetic pseudo-Receiver function traces

seismic.receiver_fn.rf_synthetic.convert_inclination_to_distance(inclinations, model='iasp91', nominal_source_depth_km=10.0)[source]

Helper function to convert range of inclinations to teleseismic distance in degrees.

Parameters:
  • inclinations (numpy.array(float)) – Array of inclination angles in degrees

  • model (str, optional) – Name of model to use for ray tracing, defaults to “iasp91”

  • nominal_source_depth_km (float, optional) – Assumed depth of source events, defaults to 10.0

Returns:

Array of teleseismic distances in degrees corresponding to input inclination angles.

Return type:

numpy.array(float)

seismic.receiver_fn.rf_synthetic.generate_synth_rf(arrival_times, arrival_amplitudes, fs_hz=100.0, window_sec=(-10, 30), f_cutoff_hz=2.0)[source]

Simple generator of synthetic R component receiver function with pulses at given arrival times.

Parameters:
  • arrival_times (iterable of float) – Iterable of arrival times as numerical values in seconds

  • arrival_amplitudes (iterable of float) – Iterable of arrival amplitudes

  • fs_hz (float, optional) – Sampling rate (Hz) of output signal, defaults to 100.0

  • window_sec (tuple, optional) – Time window over which to create signal (sec), defaults to (-10, 30)

  • f_cutoff_hz (float, optional) – Cutoff frequency (Hz) for low-pass filtering to generate realistic result, defaults to 2.0

Returns:

Array of times and corresponding signal amplitudes

Return type:

numpy.array, numpy.array

seismic.receiver_fn.rf_synthetic.synthesize_rf_dataset(H, V_p, V_s, inclinations, distances, ds, log=None, include_t3=False, amplitudes=None, baz=0.0)[source]

Synthesize RF R-component data set over range of inclinations and distances and get result as a rf.RFStream instance.

Parameters:
  • H (float) – Moho depth (km)

  • V_p (float) – P body wave velocity in uppermost layer

  • V_s (float) – S body wave velocity in uppermost layer

  • inclinations (numpy.array(float)) – Array of inclinations for which to create RFs

  • distances (numpy.array(float)) – Array of teleseismic distances corresponding to inclinations

  • ds (float) – Final sampling rate (Hz) for the downsampled output signal

  • log (logger, optional) – Logger to send output to, defaults to None

  • include_t3 (bool, optional) – If True, include the third expected multiple PpSs+PsPs

  • amplitudes (list(float), optional) – Custom amplitudes to apply to the multiples

  • baz (float, optional) – Back azimuth for metadata

Returns:

Stream containing synthetic RFs

Return type:

rf.RFStream

seismic.receiver_fn.rf_util module

Utility functions to help with RF processing and analysis.

seismic.receiver_fn.rf_util.choose_rf_source_channel(rf_type, db_station)[source]

Choose source channel for RF analysis.

Parameters:
  • rf_type (str) – The RF rotation type, should be either ‘ZRT’ or ‘LQT’

  • db_station (dict(str, list(rf.RFTrace))) – Dict of traces for a given station keyed by channel code.

Returns:

Channel code of the primary RF source channel

Return type:

str

seismic.receiver_fn.rf_util.compute_extra_rf_stats(stream)[source]

Compute extra statistics for each trace and add it to the RFTrace.stats structure.

Parameters:

stream (rf.RFStream) – RFStream to augment with additional metadata statistics.

seismic.receiver_fn.rf_util.compute_rf_snr(rf_stream)[source]

Compute signal to noise (S/N) ratio of the RF itself about the onset pulse (key ‘snr’). This SNR is a ratio of RMS amplitudes. Stores results in metadata of input stream traces.

In the LQT rotation case when rotation is working ideally, the onset pulse of the rotated transverse signals should be minimal, and a large pulse at t = 0 indicates lack of effective rotation of coordinate system, so for ‘snr’ we use a long time window after onset pulse, deliberately excluding the onset pulse, to maximize contribution to the SNR from the multiples after the onset pulse.

Parameters:

rf_stream (rf.RFStream) – R or Q component of Receiver Function

Returns:

SNR for each trace in the input stream

Return type:

numpy.array

seismic.receiver_fn.rf_util.compute_vertical_snr(src_stream)[source]

Compute the SNR of the Z component (Z before deconvolution) including the onset pulse (key ‘snr_prior’). Stores results in metadata of input stream traces. This SNR is a ratio of max envelopes.

Some authors compute this prior SNR on signal after rotation but before deconvolution, however that doesn’t make sense for LQT rotation where the optimal rotation will result in the least energy in the L component. For simplicity we compute it on Z-component only which is a reasonable estimate for teleseismic events.

Parameters:

src_stream (rf.RFStream or obspy.Stream) – Seismic traces before RF deconvolution of raw stream.

seismic.receiver_fn.rf_util.filter_by_distance(rf_stream, min_dist, max_dist)[source]

Discard RFs that fall outside the distance range(min_dist, max_dist) @param rf_stream: RFStream @param min_dist: minimum angular distance @param max_dist: maximum angular distance @return: trimmed RFStream

seismic.receiver_fn.rf_util.filter_crosscorr_coeff(rf_stream, time_window=(-2, 25), threshold_cc=0.7, min_fraction=0.15, apply_moveout=False)[source]

For each trace in the stream, compute its correlation coefficient with the other traces. Return only traces matching cross correlation coefficient criteria based on C.Sippl (2016) [see http://dx.doi.org/10.1016/j.tecto.2016.03.031]

Parameters:
  • rf_stream (rf.RFStream) – Stream of RF traces to filter, should be for a single component of a single station

  • time_window (tuple, optional) – Time window to filter by, defaults to (-2, 25)

  • threshold_cc (float, optional) – Threshold cross-correlation coefficient, defaults to 0.70. Denoted Xi in Sippl, who used value 0.80.

  • min_fraction (float, optional) – Minimum fraction of coefficients above threshold_cc, defaults to 0.15. Denoted tau in Sippl, who used value 0.15.

  • apply_moveout (bool) – Whether to apply moveout correction to Ps phase prior to computing correlation coefficients.

Returns:

Filtered stream of RF traces

Return type:

rf.RFStream

seismic.receiver_fn.rf_util.filter_invalid_radial_component(rf_stream)[source]

Filter out invalid radial RFs with amplitudes > 1 or troughs around onset time :param rf_stream: Stream of RF traces to filter, should be for a single component of a single station :type rf_stream: rf.RFStream :return: Filtered RF stream :rtype: rf.RFStream

seismic.receiver_fn.rf_util.filter_station_streams(db_station, freq_band=(None, None))[source]

Perform frequency filtering on all channels’ traces for a given station. Returns a copy of db_station with streams containing filtered results.

seismic.receiver_fn.rf_util.filter_station_to_mean_signal(db_station, min_correlation=1.0)[source]

Filter out streams which are not ‘close enough’ to the mean signal, based on simple correlation score. The term “correlation” here really just means a similarity dot product (projection of individual trace onto the mean).

seismic.receiver_fn.rf_util.find_rf_group_ids(stream)[source]

For the given stream, which is expected to have an rf_group attribute in its traces’ metadata, determine the unique set of group ids that the traces contain.

Parameters:

stream (obspy.Stream) – Stream containing traces with rf_group ids associated with them.

Returns:

Set of rf_group ids found in the traces

Return type:

set(int)

seismic.receiver_fn.rf_util.label_rf_quality_simple_amplitude(rf_type, traces, snr_cutoff=2.0, rms_amp_cutoff=0.2, max_amp_cutoff=1.0)[source]

Add RF quality label for a collection of RFs based on simple amplitude criteria computed by quality filter script. Adds quality label in-place.

Parameters:
  • rf_type (str) – The RF rotation type, should be either ‘ZRT’ or ‘LQT’

  • traces (Iterable collection of rf.RFTrace) – Iterable collection of rf.RFTrace

  • snr_cutoff (float, optional) – Minimum signal SNR, defaults to 2.0

  • rms_amp_cutoff (float, optional) – Maximum accepted RMS amplitude of signal, defaults to 0.2

  • max_amp_cutoff (float, optional) – Maximum accepted amplitude of signal, defaults to 1.0

seismic.receiver_fn.rf_util.phase_weights(stream)[source]

Phase weighting takes all the traces in a stream and computes a weighting for each sample in the stream between 0 and 1. The weighting represents how consistent is the phase angle of the signal at the same point in the time series across all streams.

If phase weighting to accentuate higher multiples than Ps, then moveout should be applied first before calling this function.

See https://doi.org/10.1111/j.1365-246X.1997.tb05664.x

Note: this function should not be applied to streams with mixed components.

Parameters:

stream (Iterable container of obspy.Trace) – Stream containing one or more traces from which phase coherence weightings will be generated.

Returns:

Array of normalized weighting factors with same length as traces in stream.

Return type:

numpy.array

seismic.receiver_fn.rf_util.read_h5_rf(src_file, network=None, station=None, loc='', root='/waveforms')[source]

Helper function to load data from hdf5 file generated by rf library or script rf_quality_filter.py. 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 a rf.RFStream container.

Return type:

rf.RFStream

seismic.receiver_fn.rf_util.rf_to_dict(rf_data)[source]

Convert RF data loaded from function read_h5_rf() into a dict format for easier addressing of selected station and channel RF traces.

Parameters:

rf_data (rf.RFStream) – RFStream data

Returns:

Nested dicts to find traces by station then channel code, with attached metadata.

Return type:

seismic.receiver_fn.rf_network_dict.NetworkRFDict

seismic.receiver_fn.rf_util.signed_nth_power(arr, order)[source]

As per DOI https://doi.org/10.1038/217533a0. Muirhead, K.J. “Eliminating False Alarms when detecting Seismic Events Automatically”

Parameters:
  • arr (numpy.array) – Compute n-th power of input array, preserving sign of original data.

  • order (float or int) – Order of the power to compute

Returns:

Input array raised to nth power.

Return type:

numpy.array

seismic.receiver_fn.rf_util.signed_nth_root(arr, order)[source]

As per DOI https://doi.org/10.1038/217533a0. Muirhead, K.J. “Eliminating False Alarms when detecting Seismic Events Automatically”

Parameters:
  • arr (numpy.array) – Compute n-th root of input array, preserving sign of original data.

  • order (float or int) – Order of the root to compute

Returns:

Input array raised to 1/nth power.

Return type:

numpy.array

seismic.receiver_fn.rf_util.trim_hdf_keys(hdf_key_list: [<class 'str'>], networks_string: str, stations_string: str) [<class 'str'>][source]

Trims a list of hdf_keys, filtering out unwanted networks and stations. @param hdf_key_list: @param networks_string: a space-separated list of networks. ‘*’ includes all. @param stations_string: a space-separated list of stations or a text file

with station names in each row, w/wo location codes. ‘*’ includes all.

@return: trimmed list

Module contents