gmprocess.waveform_processing.adjust_highpass.adjust_highpass_corner(st, step_factor=1.5, maximum_freq=0.5, max_final_displacement=0.025, max_displacement_ratio=0.2, config=None)[source]

Adjust high pass corner frequency.

Options for further refinement of the highpass corner. Currently, this includes criteria employed by:

Dawood, H.M., Rodriguez-Marek, A., Bayless, J., Goulet, C. and Thompson, E. (2016). A flatfile for the KiK-net database processed using an automated protocol. Earthquake Spectra, 32(2), pp.1281-1302.

This algorithm begins with an initial corner frequency that was selected as configured in the get_corner_frequencies step. It then checks the criteria described below; if the criteria are not met then the high pass corner is increased the multiplicative step factor until the criteria are met.

Parameters:
  • st (StationStream) – Stream of data.

  • step_factor (float) – Multiplicative factor for incrementing high pass corner.

  • maximum_freq (float) – Limit (maximum) frequency on the trial corner frequencies to search across to pass displacement checks.

  • max_final_displacement (float) – Maximum allowable value (in cm) for final displacement.

  • max_displacement_ratio (float) – Maximum allowable ratio of final displacement to max displacement.

  • config (dict) – Configuration dictionary (or None). See get_config().

Returns:

Stream with highpass corner frequency adjusted using Dawood method.

Return type:

StationStream

gmprocess.waveform_processing.filtering.bandpass_filter(st, frequency_domain=True, filter_order=5, number_of_passes=1, config=None)[source]

Apply the bandpass filter.

Parameters:
  • st (StationStream) – Stream of data.

  • frequency_domain (bool) – If true, use gmprocess frequency domain implementation; if false, use ObsPy filters.

  • filter_order (int) – Filter order.

  • number_of_passes (int) – Number of passes.

  • config (dict) – Configuration dictionary (or None). See get_config().

Returns:

Stream with bandpass filter applied.

Return type:

StationStream

gmprocess.waveform_processing.clipping.clipping_check.check_clipping(st, event, threshold=0.2, config=None)[source]

Apply clipping check.

Lower thresholds will pass fewer streams but will give less false negatives (i.e., streams in which clipping actually occurred but were missed).

Parameters:
  • st (gmprocess.core.stationstream.StationStream) – Trace of data.

  • event (gmprocess.utils.scalar_event.ScalarEvent) – ScalarEvent object.

  • threshold (float) – Threshold probability.

  • config (dict) – Configuration dictionary (or None). See get_config().

Returns:

With clipping check applied.

Return type:

StationStream

gmprocess.waveform_processing.pretesting.check_free_field(st, reject_non_free_field=True, config=None)[source]

Checks free field status of stream.

Parameters:
  • st (gmprocess.core.stationstream.StationStream) – Stream of data.

  • reject_non_free_field (bool) – Should non free-field stations be failed?

  • config (dict) – Configuration dictionary (or None). See get_config().

Returns:

Stream that has been checked for free field status.

Return type:

StationStream

gmprocess.waveform_processing.pretesting.check_instrument(st, n_max=3, n_min=2, require_two_horiz=True, config=None)[source]

Test the channels of the station.

The purpose of the maximum limit is to skip over stations with muliple strong motion instruments, which can occur with downhole or structural arrays since our code currently is not able to reliably group by location within an array.

The purpose of the minimum and require_two_horiz checks are to ensure the channels are required for subsequent intensity measures such as ROTD.

Parameters:
  • st (gmprocess.core.stationstream.StationStream) – Stream of data.

  • n_max (int) – Maximum allowed number of streams; default to 3.

  • n_min (int) – Minimum allowed number of streams; default to 1.

  • require_two_horiz (bool) – Require two horizontal components; default to False.

  • config (dict) – Configuration dictionary (or None). See get_config().

Returns:

Stream with instrument criteria applied.

Return type:

StationStream

gmprocess.waveform_processing.pretesting.check_max_amplitude(st, min=5, max=2000000.0, config=None)[source]

Check the maximum amplitude of the traces.

Checks that the maximum amplitude of the traces in the stream are within a defined range. Only applied to counts/raw data. This is a simple way to screen for clipped records, but we now recommend/prefer the check_clipping method.

Parameters:
  • st (gmprocess.core.stationstream.StationStream) – Stream of data.

  • min (float) – Minimum amplitude for the acceptable range. Default is 5.

  • max (float) – Maximum amplitude for the acceptable range. Default is 2e6.

  • config (dict) – Configuration dictionary (or None). See get_config().

Returns:

Stream that has been checked for maximum amplitude criteria.

Return type:

StationStream

gmprocess.waveform_processing.pretesting.check_sta_lta(st, sta_length=1.0, lta_length=20.0, threshold=5.0, config=None)[source]

Apply STA/LTA ratio criteria.

Checks that the maximum STA/LTA ratio of the stream’s traces is above a threshold.

Parameters:
  • st (gmprocess.core.stationstream.StationStream) – Stream of data.

  • sta_length (float) – Length of time window for STA (seconds).

  • lta_length (float) – Length of time window for LTA (seconds).

  • threshold (float) – Required maximum STA/LTA ratio to pass the test.

  • config (dict) – Configuration dictionary (or None). See get_config().

Returns:

Stream that has been checked for sta/lta requirements.

Return type:

StationStream

gmprocess.waveform_processing.sanity_checks.check_tail(st, duration=5.0, max_vel_ratio=0.3, max_dis_ratio=0.9, config=None)[source]

Check for abnormally large values in the tail of the stream.

This QA check looks for the presence of abnormally large values in the tail velocity and displacement traces. This can occur due to baseline shifts or long period noise that has not been properly filtered out that manifest as long-period drifts in the velocity and/or displacement traces.

Note that an additional problem that this check should eliminate is records in which the time window has not captured the full duration of shaking.

Parameters:
  • st (StationStream) – StationStream object.

  • duration (float) – Duration of tail.

  • max_vel_ratio (float) – Trace is labeled as failed if the max absolute velocity in the tail is greater than max_vel_ratio times the max absolute velocity of the whole trace.

  • max_dis_ratio (float) – Trace is labeled as failed if the max absolute displacement in the tail is greater than max_disp_ratio times the max absolute displacement of the whole trace.

  • config (dict) – Configuration dictionary (or None). See get_config().

Returns:

With QA checks applied.

Return type:

StationStream

gmprocess.waveform_processing.zero_crossings.check_zero_crossings(st, min_crossings=0.1, config=None)[source]

Requires a minimum zero crossing rate.

This is intended to screen out instrumental failures or resetting. Value determined empirically from observations on the GeoNet network by R Lee.

Parameters:
  • st (StationStream) – StationStream object.

  • min_crossings (float) – Minimum average number of zero crossings per second for the full trace.

  • config (dict) – Configuration dictionary (or None). See get_config().

Returns:

With zero crossing rate criteria applied.

Return type:

StationStream

gmprocess.waveform_processing.snr.compute_snr(st, event, smoothing_parameter=20.0, config=None)[source]

Compute SNR dictionaries for a stream, looping over all traces.

Parameters:
  • st (StationStream) – Trace of data.

  • event (ScalarEvent) – ScalarEvent object.

  • smoothing_parameter (float) – Konno-Omachi smoothing bandwidth parameter.

  • config (dict) – Configuration dictionary (or None). See get_config().

Returns:

With SNR dictionaries added as trace parameters.

Return type:

StationStream

gmprocess.waveform_processing.convert_units.convert_to_acceleration(st, taper=True, taper_type='hann', taper_width=0.05, taper_side='both', config=None)[source]

Convert stream to acceleration if it isn’t acceleration.

Parameters:
  • st (StationStream) – Stream of data.

  • config (dict) – Configuration dictionary (or None). See get_config().

Returns:

Stream in units of acceleration.

Return type:

StationStream

gmprocess.waveform_processing.windows.cut(st, sec_before_split=5.0, config=None)[source]

Cut/trim the record.

This method minimally requires that the windows.signal_end method has been run, in which case the record is trimmed to the end of the signal that was estimated by that method.

To trim the beginning of the record, the sec_before_split must be specified, which uses the noise/signal split time that was estimated by the windows.signal_split method.

Parameters:
  • st (StationStream) – Stream of data.

  • sec_before_split (float) – Seconds to trim before split. If None, then the beginning of the record will be unchanged.

  • config (dict) – Configuration dictionary (or None). See get_config().

Returns:

With the cut applied.

Return type:

StationStream

gmprocess.waveform_processing.detrend.detrend(st, detrending_method=None, config=None)[source]

Detrend stream.

Parameters:
  • st (StationStream) – Stream of data.

  • detrending_method (str) –

    Method to detrend; valid options include the ‘type’ options supported by obspy.core.trace.Trace.detrend as well as:

    • ’baseline_sixth_order’, which is for a baseline correction

      method that fits a sixth-order polynomial to the displacement time series, and sets the zeroth- and first-order terms to be zero. The second derivative of the fit polynomial is then removed from the acceleration time series.

    • ’pre’, for removing the mean of the pre-event noise window.

  • config (dict) – Configuration dictionary (or None). See get_config().

Returns:

Stream with that has been detrended.

Return type:

StationStream

gmprocess.waveform_processing.spectrum.fit_spectra(st, event, kappa=0.035, RP=0.55, VHC=0.7071068, FSE=2.0, density=2.8, shear_vel=3.7, R0=1.0, moment_factor=100, min_stress=0.1, max_stress=10000, config=None)[source]

Fit spectra varying stress_drop and moment.

Parameters:
  • st (gmprocess.core.stationstream.StationStream) – Stream of data.

  • event (gmprocess.utils.scalar_event.ScalarEvent) – ScalarEvent object.

  • kappa (float) – Site diminution factor (sec). Typical value for active crustal regions is about 0.03-0.04, and stable continental regions is about 0.006.

  • RP (float) – Partition of shear-wave energy into horizontal components.

  • VHC (float) – Partition of shear-wave energy into horizontal components 1 / np.sqrt(2.0).

  • FSE (float) – Free surface effect.

  • density (float) – Density at source (gm/cc).

  • shear_vel (float) – Shear-wave velocity at source (km/s).

  • R0 (float) – Reference distance (km).

  • moment_factor (float) – Multiplicative factor for setting bounds on moment, where the moment (from the catalog moment magnitude) is multiplied and divided by moment_factor to set the bounds for the spectral optimization.

  • min_stress (float) – Min stress for fit search (bars).

  • max_stress (float) – Max stress for fit search (bars).

  • config (dict) – Configuration dictionary (or None). See get_config().

Returns:

StationStream with fitted spectra parameters.

gmprocess.waveform_processing.corner_frequencies.get_corner_frequencies(st, event, method='snr', constant={'highpass': 0.08, 'lowpass': 20.0}, snr={'same_horiz': True}, magnitude={'highpass': [0.5, 0.3, 0.1], 'lowpass': [25.0, 35.0, 40.0], 'minmag': [-999.0, 3.5, 5.5]}, config=None)[source]

Select corner frequencies.

Note that this step only selects the highpass and lowpass corners. The results can be modifed by other steps (such as lowpass_max_frequency) and then the filters are applied with the lowpass_filter and highpass_filter steps.

Parameters:
  • st (gmprocess.core.stationstream.StationStream) – Stream of data.

  • event (gmprocess.utils.scalar_event.ScalarEvent) – ScalarEvent object.

  • method (str) – Which method to use; currently allowed “snr” or “constant”.

  • constant (dict) – Dictionary of constant method config options.

  • snr (dict) – Dictionary of snr method config options.

  • magnitude (dict) – Dictionary of magnitude method config options.

  • config (dict) – Configuration dictionary (or None). See get_config().

Returns:

Stream with selected corner frequencies added.

Return type:

StationStream

gmprocess.waveform_processing.filtering.highpass_filter(st, frequency_domain=True, filter_order=5, number_of_passes=1, config=None)[source]

Apply the highpass filter.

Parameters:
  • st (StationStream) – Stream of data.

  • frequency_domain (bool) – If true, use gmprocess frequency domain implementation; if false, use ObsPy filters.

  • filter_order (int) – Filter order.

  • number_of_passes (int) – Number of passes.

  • config (dict) – Configuration dictionary (or None). See get_config().

Returns:

Stream with highpass filter applied.

Return type:

StationStream

gmprocess.waveform_processing.filtering.lowpass_filter(st, frequency_domain=True, filter_order=5, number_of_passes=1, config=None)[source]

Apply the lowpass filter.

Parameters:
  • st (StationStream) – Stream of data.

  • frequency_domain (bool) – If true, use gmprocess frequency domain implementation; if false, use ObsPy filters.

  • filter_order (int) – Filter order.

  • number_of_passes (int) – Number of passes.

  • config (dict) – Configuration dictionary (or None). See get_config().

Returns:

Filtered streams.

Return type:

StationStream

gmprocess.waveform_processing.corner_frequencies.lowpass_max_frequency(st, fn_fac=0.75, lp_max=40.0, config=None)[source]

Cap lowpass corner frequency.

Options on this include a constant maximum, or as a fraction of the Nyquist.

Parameters:
  • st (gmprocess.core.stationstream.StationStream) – Stream of data.

  • fn_fac (float) – Factor to be multiplied by the Nyquist to cap the lowpass filter.

  • lp_max (float) – Maximum lowpass corner frequency (Hz).

  • config (dict) – Configuration dictionary (or None). See get_config().

Returns:

Stream with lowpass frequency adjustment applied.

Return type:

StationStream

gmprocess.waveform_processing.pretesting.max_traces(st, n_max=3, config=None)[source]

Reject a stream if it has more than n_max traces.

The purpose of this is to skip over stations with muliple strong motion instruments, which can occur with downhole or structural arrays since our code currently is not able to reliably group by location within an array.

Parameters:
  • st (gmprocess.core.stationstream.StationStream) – Stream of data.

  • n_max (int) – Maximum allowed number of streams; default to 3.

  • config (dict) – Configuration dictionary (or None). See get_config().

Returns:

Stream with max number of traces criteria applied.

Return type:

StationStream

gmprocess.waveform_processing.pretesting.min_sample_rate(st, min_sps=20.0, config=None)[source]

Require a minimum sample rate.

Parameters:
  • st (gmprocess.core.stationstream.StationStream) – Stream of data.

  • min_sps (float) – Minimum samples per second.

  • config (dict) – Configuration dictionary (or None). See get_config().

Returns:

Stream checked for sample rate criteria.

Return type:

StationStream

gmprocess.waveform_processing.nn_quality_assurance.nnet_qa(st, acceptance_threshold, model_name, config=None)[source]

Apply the neural network QA algorithm by Bellagamba et al. (2019),

Assess the quality of a stream by analyzing its two horizontal components as described in Bellagamba et al. (2019). Performs three steps: 1) Compute the quality metrics (see paper for more info) 2) Preprocess the quality metrics (deskew, standardize and decorrelate) 3) Evaluate the quality using a neural network-based model Two models are available: ‘Cant’ and ‘CantWell’. To minimize the number of low quality ground motion included, the natural acceptance threshold 0.5 can be raised (up to an extreme value of 0.95). Recommended parameters are: - acceptance_threshold = 0.5 or 0.6 - model_name = ‘CantWell’

Parameters:
  • st (StationStream) – The ground motion record to analyze. Should contain at least 2 orthogonal horizontal traces.

  • acceptance_threshold (float) – Threshold from which GM records are considered acceptable.

  • model_name (string) – name of the used model (‘Cant’ or ‘CantWell’)

  • config (dict) – Configuration dictionary (or None). See get_config().

Returns:

With QA analysis applied.

Return type:

StationStream

gmprocess.waveform_processing.instrument_response.pre_filter(st, lp_factor=0.8, hp_freq=0.05, corners=20)[source]

Filters to be applied prior to instrument response removal.

Using this is an alternative to the “pre_filt” option in ObsPy’s “remove_response” function. The prefiltering they apply is a cosine filter that is somewhat unusual and provide this option as an alternative.

One of the motivations for this is that sometimes the instrument response can radically amplify high-frequency noise near the Nyquist and it is essential to severely suppress the amplitudes within about 20% of the Nyquist prior to removing the response, which is why the default number of corners is so large.

Note that this is only

Parameters:
  • st (StationStream) – Stream of data.

  • lp_factor (float) – Multiplicative factor to get the lowpass corner (flp) from the Nyquist (Fn). Computed as flp = lp_factor*Fn.

  • hp_freq (float) – Highpass corner frequency (Hz).

  • corners (int) – Number of poles for the bandpass filter.

gmprocess.waveform_processing.prism_adaptive_baseline.prism_adaptive_baseline(st, n_iter=100, lead_vel_threshold=0.1, trail_vel_threshold=0.1, config=None)[source]

Apply the PRISM adaptive baseline correction step.

Based on the description in Jones et al. (2017; doi: 10.1785/0220160200).

Please note that the PRISM algorithm is presented as a full entire workflow and that the adaptive baseline correction step is difficult to separate into a singular step. More specifically, other processing steps that we treat as independent in gmprocess, like filtering and QA checks, are intertwined with the adaptive baseline correction (see Figs 2 and 5 in Jones et al., 2017).

The first QA step (applied to the leading and tailing velocity window) is used to decide if the adaptive baseline correction is required. We include this QA step in this implementation because we do not want to apply the baseline correction unnecessarily. This first QA check is meant to be applied prior to any filtering. Thus, this processing step should not be configured to occur after filtering steps. If this QA step passes, then the input stream is not modified. If this QA step fails, then the record is modified using the adaptive baseline correction only, i.e., subtracting the derivative of the trend of the fit to velocity from the acceleration data. This last clarification is necessary because additional filtering is described as part of the adaptive baseline correction in Jones et al., and we prefer to have that occur as a separate processing step.

The final QA step, which checks the mean velocity and displacement in the leading and trailing windows, is described to apply after bandpass filtering. Since we keep the filtering as a separate step, we also need to keep the final QA as a separate step. Thus, the stream is never failed failed by this adaptive baseline correction step.

Parameters:
  • st (StationStream) – Stream of data.

  • n_iter (int) – Number of iterations to select time2.

  • lead_vel_threshold (float) – Threshold for leading window velocity mean (cm/s).

  • trail_vel_threshold (float) – Threshold for trailing window velocity mean (cm/s).

  • config (dict) – Configuration dictionary (or None). See get_config().

Returns:

Stream with the adaptive baseline correction applied.

Return type:

StationStream

gmprocess.waveform_processing.prism_adaptive_baseline.prism_quality_criteria(st, lead_vel_threshold=0.1, trail_vel_threshold=0.1, trail_dis_threshold=0.1, config=None)[source]

Apply the PRISM quality criteria.

Based on the description in Jones et al. (2017; doi: 10.1785/0220160200).

See discussion in the prism_adaptive_baseline docstring for the reasoning for why we have made this a separate processing step and how it relates to the adaptive baseline correction algorithm.

Parameters:
  • st (StationStream) – Stream of data.

  • lead_vel_threshold (float) – Threshold for leading window velocity mean (cm/s).

  • trail_vel_threshold (float) – Threshold for trailing window velocity mean (cm/s).

  • trail_dis_threshold (float) – Threshold for trailing window displacement mean (cm/s).

Returns:

Stream with the PRISM quality criteria applied.

Return type:

StationStream

gmprocess.waveform_processing.instrument_response.remove_response(st, sensitivity_threshold=10.0, output_as_acceleration=True, pre_filt=True, f1=0.001, f2=0.005, f3=None, f4=None, water_level=60, inv=None, custom_pre_filt=None, config=None)[source]

Perform the instrument response correction.

If the response information is not already attached to the stream, then an inventory object must be provided. See the “instrument response” section of the manual for more details on how the instrument response is removed.

The f1, f2, f3, f4, and water_level arguments are all handed off to ObsPy’s remove_response function. If f3 is Null it will be set to 0.9*fn, if f4 is Null it will be set to fn.

The custom_pre_filt argument is an alternative to using the pre_filt argument in ObsPy’s remove_response function in which the prefiltering is a cosine filter that is somewhat unusual. If custom_pre_filt is specified, then pre_filt is overwritten to be None, and the prefiltering is applied with the options specified in the custom_pre_filt dictionary using a Butterworth filter prior to calling remove_response. One of the motivations for this is that sometimes the instrument response can radically amplify high-frequency noise near the Nyquist and it is essential to severely suppress the amplitudes within about 20% of the Nyquist prior to removing the response, which is why the default number of corners is so large.

The required keys for the custom_pre_filt dictionary are: - lp_factor: Multiplicative factor to get the lowpass corner (flp) from the

Nyquist (Fn). Computed as flp = lp_factor*Fn.

  • hp_freq: Highpass corner frequency (Hz).

  • corners: Number of poles for the bandpass filter.

Parameters:
  • st (StationStream) – Stream of data.

  • sensitivity_threshold (float) – Traces will be failed if the overall reported sensitivity differs by more than this threshold (units are percent) when compared to the sensitivity as computed by combining the reported gains of each response stage.

  • output_as_acceleration (bool) – For velocity instruments, differentiate the waveform after the instrument response is removed to convert the units to acceleration.

  • pre_filt (bool) – Apply a bandpass filter in frequency domain to the data before deconvolution?

  • f1 (float) – Frequency 1 for pre-filter.

  • f2 (float) – Frequency 2 for pre-filter.

  • f3 (float) – Frequency 3 for pre-filter.

  • f4 (float) – Frequency 4 for pre-filter.

  • water_level (float) – Water level for deconvolution.

  • inv (obspy.core.inventory.inventory) – Obspy inventory object containing response information.

  • custom_pre_filt (dict, None) – If not None, a dictionary with the following options for pre-filtering.

  • config (dict) – Configuration dictionary (or None). See get_config().

Returns:

With instrument response correction applied.

Return type:

StationStream

gmprocess.waveform_processing.resample.resample(st, new_sampling_rate=None, method=None, a=None, config=None)[source]

Resample stream.

Parameters:
  • st (StationStream) – Stream of data.

  • new_sampling_rate (float) – New sampling rate, in Hz.

  • method (str) – Method for interpolation. Currently only supports ‘lanczos’.

  • a (int) – Width of the Lanczos window, in number of samples.

  • config (dict) – Configuration dictionary (or None). See get_config().

Returns:

Stream with resampling applied.

Return type:

StationStream

gmprocess.waveform_processing.adjust_highpass_ridder.ridder_fchp(st, target=0.02, tol=0.001, maxiter=30, maxfc=0.5, config=None)[source]

Search for highpass corner using Ridder’s method.

Search such that the criterion that the ratio between the maximum of a third order polynomial fit to the displacement time series and the maximum of the displacement timeseries is a target % within a tolerance.

This algorithm searches between a low initial corner frequency a maximum fc.

Method developed originally by Scott Brandenberg

Parameters:
  • st (StationStream) – Stream of data.

  • target (float) – target percentage for ratio between max polynomial value and max displacement.

  • tol (float) – tolerance for matching the ratio target

  • maxiter (float) – maximum number of allowed iterations in Ridder’s method

  • maxfc (float) – Maximum allowable value of the highpass corner freq.

  • int_method (string) – method used to perform integration between acceleration, velocity, and dispacement. Options are “frequency_domain”, “time_domain_zero_init” or “time_domain_zero_mean”

  • config (dict) – Configuration dictionary (or None). See get_config().

Returns:

Stream with the highpass corner adjusted using Ridder’s method.

Return type:

StationStream

gmprocess.waveform_processing.snr.snr_check(st, mag, threshold=2.0, min_freq='f0', max_freq=5.0, f0_options={'ceiling': 2.0, 'floor': 0.1, 'shear_vel': 3.7, 'stress_drop': 10}, config=None)[source]

Check signal-to-noise ratio.

Requires noise/signal windowing to have succeeded.

Parameters:
  • st (StationStream) – Trace of data.

  • mag (float) – Earthquake magnitude.

  • threshold (float) – Threshold SNR value.

  • min_freq (float or str) – Minimum frequency for threshold to be exeeded. If ‘f0’, then the Brune corner frequency will be used.

  • max_freq (float) – Maximum frequency for threshold to be exeeded.

  • smoothing_parameter (float) – Konno-Omachi smoothing bandwidth parameter.

  • f0_options (dict) – Dictionary of f0 options (see config file).

  • config (dict) – Configuration dictionary (or None). See get_config().

Returns:

Trace with SNR check.

Return type:

trace

gmprocess.waveform_processing.zero_pad.strip_zero_pad(st, config=None)[source]

Remove zero pads from streams.

Parameters:
  • st (StationStream) – Stream of data.

  • config (dict) – Configuration dictionary (or None). See get_config().

Returns:

Streams with the zero padding removed.

Return type:

StationStream

gmprocess.waveform_processing.taper.taper(st, type='hann', width=0.05, side='both', config=None)[source]

Taper streams.

Parameters:
  • st (StationStream) – Stream of data.

  • type (str) – Taper type.

  • width (float) – Taper width as percentage of trace length.

  • side (str) – Valid options: “both”, “left”, “right”.

  • config (dict) – Configuration dictionary (or None). See get_config().

Returns:

Stream with the taper applied.

Return type:

StationStream

gmprocess.waveform_processing.zero_pad.zero_pad(st, length='fhp', config=None)[source]

Add zero pads to streams.

Parameters:
  • st (StationStream) – Stream of data.

  • length (float, str) – The length (in sec) to pad with zeros before and after the trace, or “fhp” to compute the zero pad length from the filter order and the highpass corner as 1.5 * filter_order / flc.

  • config (dict) – Configuration dictionary (or None). See get_config().

Returns:

Streams with the zero padding applied.

Return type:

StationStream