pysep.utils.process

Utilities for processing waveforms

Module Contents

Functions

format_streams_for_rotation(st)

ObsPy requires specific channel naming to get stream rotation working.

append_back_azimuth_to_stats(st, inv, event)

ObsPy rotation requires back azimuth information which is accessed through

_append_null_trace(st, component[, cmpaz, cmpinc])

Create a copy of a trace with zeroed out data, used for rotation

_check_component_availability(st)

Quick logic check to see if stream is missing horizontals or vertical

rotate_to_uvw(st)

UVW orthogonal frame rotation

trim_start_end_times(st[, starttime, endtime, fill_value])

Trim all traces in a Stream to a uniform start and end times. If no

merge_gapped_data(st[, fill_value, gap_fraction])

Merges traces with the same ID, filling gappy data with values if requested,

resample_data(st, resample_freq[, method])

Apply a lowpass filter and taper data before resampling it using

zerophase_chebychev_lowpass_filter(tr, freqmax)

Custom Chebychev type two zerophase lowpass filter useful for decimation

estimate_prefilter_corners(tr)

Estimates corners for pre-deconvolution (instrument response removal)

pysep.utils.process.format_streams_for_rotation(st)[source]

ObsPy requires specific channel naming to get stream rotation working. Comb through the stream and check that this is correct before passing stream to rotation

Parameters:

st (obspy.core.stream.Stream) – Stream to format for rotation

Rtype st:

obspy.core.stream.Stream

Return st:

Stream that has been formatted for rotation

pysep.utils.process.append_back_azimuth_to_stats(st, inv, event)[source]

ObsPy rotation requires back azimuth information which is accessed through trace.stats.back_azimuth. This function appends back azimuth values to the stream by cacluating using metadata from the inventory and event objects.

Note

This was written and then I realized that _append_sac_headers_trace already does this, so it is unused, but left incase it’s useful later.

Parameters:
  • st (obspy.core.stream.Stream) – Stream with missing components

  • inv (obspy.core.inventory.Inventory) – optional user-provided inventory object which will force a skip over StationXML/inventory searching

  • event (obspy.core.event.Event) – optional user-provided event object which will force a skip over QuakeML/event searching

Return type:

obspy.core.stream.Stream

Returns:

stream with back azimuth values appended

pysep.utils.process._append_null_trace(st, component, cmpaz=None, cmpinc=None)[source]

Create a copy of a trace with zeroed out data, used for rotation

Parameters:
  • st (obspy.core.stream.Stream) – Stream with missing components

  • component (str) – component of data to turn into a null trace

  • cmpaz (float) – component azimuth for SAC header

  • cmpinc (float) – component inclination for SAC header

Rtype st:

obspy.core.stream.Stream

Return st:

Stream that has had null traces appended

pysep.utils.process._check_component_availability(st)[source]

Quick logic check to see if stream is missing horizontals or vertical

Parameters:

st (obspy.core.stream.Stream) – Steam to check

Return type:

tuple of bool

Returns:

(is the vertical component missing?, is one or both horizontal components missing?)

pysep.utils.process.rotate_to_uvw(st)[source]

UVW orthogonal frame rotation

In Symmetric Triaxial Seismometers, the sensing elements are also arranged to be mutually orthogonal, but instead of one axis being vertical, all three are inclined upwards from the horizontal at precisely the same angle, as if they were aligned with the edges of a cube balanced on a corner.

TODO test this, not sure how it’s supposed to work

See reference: http://link.springer.com/referenceworkentry/10.1007/978-3-642-36197-5_194-1

Rotation matrix reference: http://link.springer.com/referenceworkentry/10.1007/

978-3-642-36197-5_194-1#page-1

Parameters:

st (obspy.core.stream.Stream) – three-component stream to rotate the UVW

Return type:

obspy.core.stream.Stream

Returns:

Stream that has been rotated to UVW

pysep.utils.process.trim_start_end_times(st, starttime=None, endtime=None, fill_value=None)[source]

Trim all traces in a Stream to a uniform start and end times. If no starttime or endtime are provided, they are selected as the outer bounds of the Streams time bounds.

Parameters:
  • st (obspy.core.stream.Stream) – stream to merge and trim start and end times for

  • starttime (obspy.core.UTCDateTime) – user-defined starttime to trim stream. If None, will get the maximum starttime in entire stream

  • endtime (obspy.core.UTCDateTime) – user-defined endtime to trim stream. If None, will get the minimum endtime in entire stream

  • fill_value (str or int or float) –

    How to deal with data gaps (missing sections of waveform before trace start or trace end w.r.t requested start and end times). NoneType by default, which means data with gaps are removed completely.

    Options include:

    • ’mean’: fill with the mean of all data values in the gappy data

    • <int or float>: fill with a constant, user-defined value, e.g.,

    0 or 1.23 or 9.999 - None: do not fill data gaps

pysep.utils.process.merge_gapped_data(st, fill_value=None, gap_fraction=1.0)[source]

Merges traces with the same ID, filling gappy data with values if requested, otherwise removes traces that have gapped data.

Note

Be careful about fill_value data types, as there are no checks that the fill value matches the internal data types. This may cause unexpected errors.

Parameters:
  • st (obspy.core.stream.Stream) – stream to merge and trim start and end times for

  • fill_value (str or int or float) –

    How to deal with data gaps (missing sections of waveform over a continuous time span). NoneType by default, which means data with gaps are removed completely. Users who want access to data with gaps must choose how gaps are filled. See API for ObsPy.core.stream.Stream.merge() for how merge is handled:

    Options include:

    • ’mean’: fill with the mean of all data values in the gappy data

    • <int or float>: fill with a constant, user-defined value, e.g.,

    0 or 1.23 or 9.999 - ‘interpolate’: linearly interpolate from the last value pre-gap to the first value post-gap - ‘latest’: fill with the last value of pre-gap data - False: do not fill data gaps, which will lead to stations w/ data gaps being removed.

  • gap_fraction (float) – if fill_data_gaps is not None, determines the maximum allowable fraction (percentage) of data that gaps can comprise. For example, a value of 0.3 means that 30% of the data (in samples) can be gaps that will be filled by fill_data_gaps. Traces with gap fractions that exceed this value will be removed. Defaults to 1. (100%) of data can be gaps.

Return type:

obspy.core.stream.Stream

Returns:

Stream that has been trimmed

pysep.utils.process.resample_data(st, resample_freq, method='interpolate')[source]

Apply a lowpass filter and taper data before resampling it using one of a choice of resampling strategies

TODO did not refactor resample_cut which doesn’t seem to do anything but

is included in old code

Parameters:
  • st (obspy.core.stream.Stream) – Stream to resample data for

  • resample_freq (float) – frequency to resample for

  • method (str) –

    choice of ObsPy resampling strategy: - interpolate: interpolate data, requires anti-aliasing lowpass filter

    before interpolating

    • resample: reample data with fourier methods

Return type:

obspy.core.stream.Stream

Returns:

Stream that has been resampled and maybe filtered

pysep.utils.process.zerophase_chebychev_lowpass_filter(tr, freqmax)[source]

Custom Chebychev type two zerophase lowpass filter useful for decimation filtering. This filter is stable up to a reduction in frequency with a factor of 10. If more reduction is desired, simply decimate in steps. Partly based on a filter in ObsPy.

Filter parameters:

rp: maximum ripple of passband rs: attenuation of stopband ws: stop band frequency wp: pass band frequency

Note

Should be replaced once ObsPy has a proper decimation filter.

Note

This function affects the trace in place!

Parameters:
  • tr (obspy.core.trace.Trace) – The trace to be filtered.

  • freqmax (float) – The desired lowpass frequency.

pysep.utils.process.estimate_prefilter_corners(tr)[source]

Estimates corners for pre-deconvolution (instrument response removal) Essentially band limits the data based on the minimum and maximum allowable periods based on the sampling rate and overall waveform length.

See also: https://ds.iris.edu/files/sac-manual/commands/transfer.html

Replaces the old get_pre_filt function

Parameters:

tr (obspy.core.trace.Trace) – trace from which stats are taken

Return type:

tuple of float

Returns:

frequency corners (f0, f1, f2, f3)