pytomo3d.signal package¶
Submodules¶
pytomo3d.signal.process module¶
Methods that handles signal data processing
copyright: | Wenjie Lei (lei@princeton.edu), 2016 |
---|---|
license: | GNU General Public License, Version 3 (http://www.gnu.org/copyleft/gpl.html) |
- pytomo3d.signal.process.check_array_order(array, order='ascending')[source]¶
Check whether a array is in ascending order or descending order :param array: the input array :param order: “ascending” or “descending” :return:
- pytomo3d.signal.process.filter_stream(st, pre_filt)[source]¶
Filter a stream
Parameters: - st –
- per_filt –
Returns:
- pytomo3d.signal.process.filter_trace(tr, pre_filt)[source]¶
Perform a frequency domain taper mimicing the behavior during the response removal, without a actual response removal.
Parameters: - tr – input trace
- pre_filt (Numpy.array or list) – frequency array(Hz) in ascending order, to define the four corners of filter, for example, [0.01, 0.1, 0.2, 0.5].
Returns: filtered trace
- pytomo3d.signal.process.flex_cut_stream(st, cut_start, cut_end, dynamic_npts=0)[source]¶
Flexible cut stream
Parameters: - st – input stream
- cut_start – cut starttime
- cut_end – cut endtime
- dynamic_npts – the dynamic number of points before cut_start and after cut_end
Returns:
- pytomo3d.signal.process.flex_cut_trace(trace, cut_starttime, cut_endtime, dynamic_npts=0)[source]¶
not cut strictly(but also based on the original trace length)
Parameters: - trace (obspy.Trace) – input trace
- cut_starttime (obspy.UTCDateTime) – starttime of cutting
- cut_endtime (obspy.UTCDateTime) – endtime of cutting
- pytomo3d.signal.process.interpolate_stream(stream, sampling_rate, starttime=None, npts=None)[source]¶
For a fairly large stream, use stream.interpolate() is not a wise choice since if there is one trace fails, then the whole interpolation will stop. So it is better to operate interpolation on the trace level
- pytomo3d.signal.process.process(st, remove_response_flag=False, inventory=None, filter_flag=False, pre_filt=None, starttime=None, endtime=None, resample_flag=False, sampling_rate=1.0, taper_type='hann', taper_percentage=0.05, rotate_flag=False, event_latitude=None, event_longitude=None, **kwargs)[source]¶
Stream processing function defined for general purpose of tomography. The advantage of using Stream, rather than than Trace, is that rotation could be operated if the Stream contains multiple channels. But this function also deals with Trace, but you need to set rotate_flag to False
Parameters: - st (obspy.Stream) – input stream
- remove_response_flag (bool) – flag for remove instrument response. If True, then inv should be specified, and filter_flag would not be taken caren of. If you want just filter the seismogram, please leave this to False and set filter_flag to True.
- inventory (obspy.Inventory) – station inventory information
:param filter_flag:flag for filter the seismogram :type filter_flag: bool :param pre_filt: list of tuple of 4 corner frequency for filter,
in ascending order(unit: Hz)Parameters: - starttime (obspy.UTCDateTime) – starttime of cutting
- endtime (obspy.UTCDateTime) – endtime of cutting
- resample_flag (bool) – flag for data resampling
- sampling_rate (float) – resampling rate(unit: Hz)
- taper_type (str) – taper type, options from obspy taper
- taper_percentage (float) – percentage of taper
- rotate_flag – rotate flag. If true, both inv and event location information should be provided
- event_latitude (float) – latitude of event, for rotation usage
- event_longitude (float) – longitude of event, for rotation usage
Returns: processed stream
pytomo3d.signal.rotate module¶
Methods that handles rotation of seismograms :copyright:
Wenjie Lei (lei@princeton.edu), 2016
license: | GNU General Public License, Version 3 (http://www.gnu.org/copyleft/gpl.html) |
---|
- pytomo3d.signal.rotate.check_orthogonality(azim1, azim2)[source]¶
Check if two azimuth are orthogonal
- pytomo3d.signal.rotate.rotate_12_NE(d1, d2, azim1, azim2)[source]¶
Rotate from any two orthogonal horizontal components to EN components
Parameters: Returns: East and North component of seismogram.
- pytomo3d.signal.rotate.rotate_12_RT(d1, d2, baz, azim1, azim2)[source]¶
Rotate from any two orthogonal horizontal components to RT components
Parameters: Returns: Radial and Transeversal component of seismogram. (None, None) returned if input two components are not orthogonal
- pytomo3d.signal.rotate.rotate_12_RT_func(st, inv, method='12->RT', back_azimuth=None)[source]¶
Rotate 12 component to RT
Parameters: - st – input stream
- inv – station inventory information
- method – rotation method
- back_azimuth – back azimuth(station to event azimuth)
Returns: rotated stream
- pytomo3d.signal.rotate.rotate_NE_12(n, e, azim1, azim2)[source]¶
Rotate from East and North components to give two orghogonal horizontal components. Returned values are (d1, d2) and (d1, d2, Vertical) will form a left-handed coordinate system.
Parameters: Returns: two horizontal orthogonal seismogram after rotation.
- pytomo3d.signal.rotate.rotate_RT_12(r, t, baz, azim1, azim2)[source]¶
Rotate from any two orthogonal horizontal components to RT components
Parameters: Returns: Radial and Transeversal component of seismogram.
- pytomo3d.signal.rotate.rotate_certain_angle(d1, d2, angle)[source]¶
Basic rotating function.
- (d1, d2, Vertical) should form a left-handed coordinate system, i.e.,
- Azimuth_{d2} = Azimuth_{d1} + 90.0
For example, (North, East, Vertical) & (Radial, Transverse, Vertical) are both left-handed coordinate systems. The return value (dnew1, dnew2, vertic) should also form a left-handed coordinate system. The angle is azimuth differnce between d1 and dnew1, i.e.,
angle = Azimuth_{dnew1} - Azimuth_{d1}Parameters: - d1 (ndarray) – Data of one of the two horizontal components
- d2 (ndarray) – Data of one of the two horizontal components
- angle (float) – component azimuth of data2
Returns: two new components after rotation
- pytomo3d.signal.rotate.rotate_stream(st, event_latitude, event_longitude, station_latitude=None, station_longitude=None, inventory=None, mode='ALL')[source]¶
Rotate a stream to radial and transverse components based on the station information and event information
Parameters: - st (obspy.Stream) – input stream
- event_latitude (float) – event latitude
- event_longitude (float) – event longitude
- station_latitude (float) – station latitude. If not provided, extract information from inventory
- station_longitude (float) – station longitude. If not provided, extrace information from inventory
- inv – station inventory information. If you want to rotate
“12” components, you need to provide inventory since only station and station_longitude is not enough. :type inv: obspy.Inventory :param mode: rotation mode, could be one of:
- “NE”: rotate only North and East channel to RT
- “12”: rotate only 1 and 2 channel, like “BH1” and “BH2” to RT
- “all”: rotate all components to RT
Returns: rotated stream(obspy.Stream)