pyatoa.utils.process
Tools for processing obspy.Stream or obspy.Trace objects Used for preprocessing data through filtering and tapering, zero padding etc. Also contains tools for synthetic traces such as source time function convolutions
Module Contents
Functions
|
Choose the appropriate filter depending on the ranges given. |
|
Taper the leading edge of the waveform. If a time offset is given, |
|
Zero pad the data of a stream, change the starttime to reflect the change. |
|
Trim two streams to common start and end times, |
|
Resampling can cause sample number differences which will lead to failure |
|
Normalize amplitudes in traces to the absolute maximum amplitude of the |
|
Check to make sure a stream object has not yet been run through |
|
Convolve function with a Gaussian window source time function. |
- pyatoa.utils.process.apply_filter(st, min_period=None, max_period=None, min_freq=None, max_freq=None, corners=2, zerophase=True, **kwargs)[source]
Choose the appropriate filter depending on the ranges given. Either periods or frequencies can be given. Periods will be prioritized. Uses Butterworth filters by default.
Filters the stream in place. Kwargs passed to filter functions.
- Parameters:
st (obspy.core.stream.Stream) – stream object to be filtered
min_period (float) – minimum filter bound in units of seconds
max_period (float) – maximum filter bound in units of seconds
min_freq (float) – optional minimum filter bound in units of Hz, will be overwritten by max_period if given
max_freq (float) – optional maximum filter bound in units of Hz, will be overwritten by min_period if given
corners (int) – number of filter corners to be passed to ObsPy filter functions
zerophase (bool) – if True, run filter backwards and forwards to avoid any phase shifting
- Return type:
obspy.core.stream.Stream
- Returns:
Filtered stream object
- pyatoa.utils.process.taper_time_offset(st, taper_percentage=0.05, time_offset_sec=0)[source]
Taper the leading edge of the waveform. If a time offset is given, e.g. 20s before the event origin time (T_0), taper all the way up from T=0 to T=T_0, to ensure that there are no impulse-like signals prior to the event origin.
- Parameters:
st (obspy.core.stream.Stream) – Stream object to be tapered
taper_percentage (float) – default taper percentage
time_offset_sec (float) – Any time offset between the start of the stream to the event origin time. All time between these two points will be tapered to reduce any signals prior to the event origin.
- Return type:
obspy.core.stream.Stream
- Returns:
tapered Stream object
- pyatoa.utils.process.zero_pad(st, pad_length_in_seconds, before=True, after=True)[source]
Zero pad the data of a stream, change the starttime to reflect the change. Useful for if e.g. observed data starttime comes in later than synthetic.
- Parameters:
- Rtype st:
obspy.stream.Stream
- Return st:
stream with zero padded data object
- pyatoa.utils.process.trim_streams(st_a, st_b, precision=0.001, force=None)[source]
Trim two streams to common start and end times, Do some basic preprocessing before trimming. Allows user to force one stream to conform to another. Assumes all traces in a stream have the same time. Prechecks make sure that the streams are actually different
- Parameters:
- Return type:
tuple of obspy.stream.Stream
- Returns:
trimmed stream objects in the same order as input
- pyatoa.utils.process.match_npts(st_a, st_b, force=None)[source]
Resampling can cause sample number differences which will lead to failure of some preprocessing or processing steps. This function ensures that npts matches between traces by extending one of the traces with zeros. A small taper is applied to ensure the new values do not cause discontinuities.
Note
its assumed that all traces within a single stream have the same npts
- Parameters:
st_a (obspy.stream.Stream) – one stream to match samples with
st_b (obspy.stream.Stream) – one stream to match samples with
force (str) – choose which stream to use as the default npts, defaults to ‘a’, options: ‘a’, ‘b’
- Return type:
tuple (obspy.stream.Stream, obspy.stream.Stream)
- Returns:
streams that may or may not have adjusted npts, returned in the same order as provided
- pyatoa.utils.process.normalize(st_a, st_b, choice)[source]
Normalize amplitudes in traces to the absolute maximum amplitude of the other, or normalize all traces so that their maximum value is a given value.
- Parameters:
st_a (obspy.stream.Stream) – one stream to match samples with
st_b (obspy.stream.Stream) – one stream to match samples with
choice (str or int or float) – choose which stream to use as the default npts, defaults to ‘a’, options: ‘a’, ‘b’ or a value that you want to set the max amplitude to be equal to, e.g., 1
- Return type:
tuple (obspy.stream.Stream, obspy.stream.Stream)
- Returns:
streams that have been normalized based on choice
- Raises:
NotImplementedError – If incorrect value of choice is provided
- pyatoa.utils.process.is_preprocessed(st, filter_only=True)[source]
Check to make sure a stream object has not yet been run through preprocessing. Assumes that a fresh stream will have no processing attribute in their stats, or if they do, will not have been filtered (getting cut waveforms from FDSN appends a ‘trim’ stat).
- Parameters:
st (obspy.stream.Stream) – stream to check processing on
filter_only (bool) – only check if the stream has been filtered, as other processing steps (e.g., demeaning, rotating) will also lead to a ‘processing’ stat. Usually this is what you want to check as filtering is one of the last steps in the processing chain.
- Return type:
- Returns:
if preprocessing has occurred
- pyatoa.utils.process.stf_convolve(st, half_duration, source_decay=4.0, time_shift=None, time_offset=None)[source]
Convolve function with a Gaussian window source time function. Design follows Specfem3D Cartesian “comp_source_time_function.f90”
hdur given is `hdur`_Gaussian = hdur/SOURCE_DECAY_MIMIC_TRIANGLE with SOURCE_DECAY_MIMIC_TRIANGLE ~ 1.68
This gaussian uses a strong decay rate to avoid non-zero onset times, while still miicking a triangle source time function
- Parameters:
st (obspy.stream.Stream) – stream object to convolve with source time function
half_duration (float) – the half duration of the source time function, usually provided in moment tensor catalogs
source_decay (float) – the decay strength of the source time function, the default value of 4 gives a Gaussian. A value of 1.68 mimics a triangle.
time_shift (float) – Time shift of the source time function in seconds
- Return type:
obspy.stream.Stream
- Returns:
stream object which has been convolved with a source time function