pygrt.utils
- file:
utils.py
- author:
Zhu Dengda (zhudengda@mail.iggcas.ac.cn)
- date:
2024-07-24
- 该文件包含一些数据处理操作上的补充:
1、剪切源、张裂源、单力源、爆炸源、矩张量源 通过格林函数合成理论地震图的函数
2、Stream类型的时域卷积、微分、积分 (基于numpy和scipy)
3、Stream类型写到本地sac文件,自定义名称
4、读取波数积分和峰谷平均法过程文件
5、其它辅助函数
- pygrt.utils.gen_syn_from_gf_DC(st: Stream | dict, M0: float, strike: float, dip: float, rake: float, az: float = -999, ZNE=False, calc_upar: bool = False, **kwargs)[源代码]
Shear source, the unit of angles is all degrees(°)
- 参数:
st – Green’s functions in a
obspy.Stream(dynamic-case) or a dict (static-case)M0 – scalar seismic moment (dyne*cm)
strike – 0 <= strike <= 360 (north=0, clockwise as positive)
dip – 0 <= dip <= 90
rake – -180 <= rake <= 180 (on the fault plane, counterclockwise as positive)
az – azimuth, 0 <= az <= 360 (not used for static case)
ZNE – whether output in ‘ZNE’-coord, default is ‘ZRT’
calc_upar – whether calculate the spatial derivatives of displacements.
kwargs – For static rsults, you can set “xarr” and “yarr” to define a new XY grid, the program will automatically find the nearest Green’s functions for each node.
- 返回:
stream -
obspy.Stream
- pygrt.utils.gen_syn_from_gf_SF(st: Stream | dict, S: float, fN: float, fE: float, fZ: float, az: float = -999, ZNE=False, calc_upar: bool = False, **kwargs)[源代码]
Single-force source (dyne)
- 参数:
st – Green’s functions in a
obspy.Stream(dynamic-case) or a dict (static-case)S – scaling factor (dyne)
fN – coefficient of Northward force
fE – coefficient of Eastward force
fZ – coefficient of Vertical(Downward) force
az – azimuth, 0 <= az <= 360 (not used for static case)
ZNE – whether output in ‘ZNE’-coord, default is ‘ZRT’
calc_upar – whether calculate the spatial derivatives of displacements.
kwargs – For static rsults, you can set “xarr” and “yarr” to define a new XY grid, the program will automatically find the nearest Green’s functions for each node.
- 返回:
stream -
obspy.Stream
- pygrt.utils.gen_syn_from_gf_EX(st: Stream | dict, M0: float, az: float = -999, ZNE=False, calc_upar: bool = False, **kwargs)[源代码]
Explosion
- 参数:
st – Green’s functions in a
obspy.Stream(dynamic-case) or a dict (static-case)M0 – scalar seismic moment (dyne*cm)
az – azimuth, 0 <= az <= 360 (not used for static case)
ZNE – whether output in ‘ZNE’-coord, default is ‘ZRT’
calc_upar – whether calculate the spatial derivatives of displacements.
kwargs – For static rsults, you can set “xarr” and “yarr” to define a new XY grid, the program will automatically find the nearest Green’s functions for each node.
- 返回:
stream -
obspy.Stream
- pygrt.utils.gen_syn_from_gf_MT(st: Stream | dict, M0: float, MT: _SupportsArray[dtype[Any]] | _NestedSequence[_SupportsArray[dtype[Any]]] | bool | int | float | complex | str | bytes | _NestedSequence[bool | int | float | complex | str | bytes], az: float = -999, ZNE=False, calc_upar: bool = False, **kwargs)[源代码]
Moment tensor
- 参数:
st – Green’s functions in a
obspy.Stream(dynamic-case) or a dict (static-case)M0 – scalar seismic moment (dyne*cm)
MT – coefficient of Moment tensor (M11, M12, M13, M22, M23, M33), subscripts 1,2,3 denote Northward,Eastward,Downward
az – azimuth, 0 <= az <= 360 (not used for static case)
ZNE – whether output in ‘ZNE’-coord, default is ‘ZRT’
calc_upar – whether calculate the spatial derivatives of displacements.
kwargs – For static rsults, you can set “xarr” and “yarr” to define a new XY grid, the program will automatically find the nearest Green’s functions for each node.
- 返回:
stream -
obspy.Stream
- pygrt.utils.compute_strain(st: Stream | dict)[源代码]
Compute dynamic/static strain tensor from synthetic spatial derivatives.
- 参数:
st – synthetic spatial derivatives
obspy.Streamclass for dynamic case, dict class for static case.- 返回:
stres - dynamic/static strain tensor, in
obspy.Streamclass or dict class.
- pygrt.utils.compute_rotation(st: Stream | dict)[源代码]
Compute dynamic/static rotation tensor from synthetic spatial derivatives.
- 参数:
st – synthetic spatial derivatives
obspy.Streamclass for dynamic case, dict class for static case.- 返回:
stres - dynamic/static rotation tensor, in
obspy.Streamclass or dict class.
- pygrt.utils.compute_stress(st: Stream | dict)[源代码]
Compute dynamic/static stress tensor from synthetic spatial derivatives.
- 参数:
st – synthetic spatial derivatives
obspy.Streamclass for dynamic case, dict class for static case.- 返回:
stres - dynamic/static stress tensor (unit: dyne/cm^2 = 0.1 Pa), in
obspy.Streamclass or dict class.
- pygrt.utils.stream_convolve(st0: Stream, signal0: ndarray, inplace=True)[源代码]
convolve each trace with a signal
- 参数:
st0 –
obspy.Streamsignal0 – convolution signal
inplace – whether change in-place
- 返回:
stream - convolution result,
obspy.Stream
- pygrt.utils.stream_integral(st0: Stream, inplace=True)[源代码]
Perform integration on each trace
- 参数:
st0 –
obspy.Streaminplace – whether change in-place
- 返回:
stream - integration result,
obspy.Stream
- pygrt.utils.stream_diff(st0: Stream, inplace=True)[源代码]
Perform central difference on each trace
- 参数:
st0 –
obspy.Streaminplace – whether change in-place
- 返回:
stream - difference result,
obspy.Stream
- pygrt.utils.stream_write_sac(st: Stream, dir: str)[源代码]
save each trace to “dir/{channel}.sac”
- 参数:
st –
obspy.Streamdir – saving directory
- pygrt.utils.read_kernels_freqs(statsdir: str, vels: ndarray | None = None, ktypes: List[str] | None = None)[源代码]
read all statsfiles in statsdir (except that of 0 frequency). If record wavenumber, interpolate to the phase velocity.
- 参数:
statsdir – directory path
vels – When a positive-order vels (km/s) is specified, files starting with K_ are read and linear interpolation from wavenumber to phase velocity is performed. Otherwise read the files starting with C_
ktype – Specify the return of a series of kernel function names, such as EX_q, DS_w, etc. By default, all are returned
- 返回:
kerDct - kernel functions in a dict
- pygrt.utils.read_statsfile(statsfile: str)[源代码]
read a statsfile
- 参数:
statsfile – File path (Wildcards can be used to simplify input)
- 返回:
data - numpy.ndarray custom type array
- pygrt.utils.read_statsfile_ptam(statsfile: str)[源代码]
read a statsfile from PTAM process
- 参数:
statsfile – File path (Wildcards can be used to simplify input)
- 返回:
data1 - numpy.ndarray custom type array, during DCM or (SA)FIM
data2 - numpy.ndarray custom type array, during PTAM
ptam_data - numpy.ndarray custom type array, record the peak/trough from PTAM
dist - epicentral distance from the filename (km)
- pygrt.utils.plot_statsdata(statsdata: ndarray, dist: float, srctype: str, ptype: str, RorI: bool | int = True, fig: Figure | None = None, axs: Axes | None = None)[源代码]
Based on the data read by the
read_statsfilefunction, plot the kernel function \(F(k,\omega)\), the integrand \(F(k,\omega)J_m(kr)k\), and calculate the cumulative integral \(\sum F(k,\omega)J_m(kr)k\) .备注
Not every source type corresponds to every order and every integration type, see 格林函数分类 for details.
- 参数:
statsdata – return value of
read_statsfilefunctiondist – epicentral distance (km)
srctype – abbreviation of source type, including EX, VF, HF, DD, DS, SS
ptype – integration type (0,1,2,3)
RorI – whether to plot real or imaginary part, default is real part, pass 2 to plot both
fig – user-defined matplotlib.Figure object, default is None
axs – user-defined matplotlib.Axes object array (three elements), default is None
- 返回:
fig - matplotlib.Figure object
(ax1,ax2,ax3) - matplotlib.Axes object array
- pygrt.utils.plot_statsdata_ptam(statsdata1: ndarray, statsdata2: ndarray, statsdata_ptam: ndarray, dist: float, srctype: str, ptype: str, RorI: bool | int = True, fig: Figure | None = None, axs: Axes | None = None)[源代码]
Based on data read by the
read_statsfile_ptamfunction, simply calculate and plot the cumulative integral as well as the peak/trough positions used by PTAM.备注
Not every source type corresponds to every order and every integration type, see 格林函数分类 for details.
- 参数:
statsdata1 – integral process data during DWM or FIM
statsdata2 – integral process data during PTAM
statsdata_ptam – peak/trough positions and amplitudes from PTAM
dist – epicentral distance (km)
srctype – abbreviation of source type, including EX, VF, HF, DD, DS, SS
ptype – integration type (0, 1, 2, 3)
RorI – whether to plot real or imaginary part, default is real part, pass 2 to plot both
fig – user-defined matplotlib.Figure object, default is None
axs – user-defined matplotlib.Axes object array (three elements), default is None
- 返回:
fig - matplotlib.Figure object
(ax1, ax2, ax3) - matplotlib.Axes object array
- pygrt.utils.solve_lamb1(nu: float, ts: ndarray, azimuth: float)[源代码]
solve the first-kind Lamb’s problem using the generalized closed-form solution, see:
张海明, 冯禧 著. 2024. 地震学中的 Lamb 问题(下). 科学出版社
- 参数:
nu – possion ratio in (0, 0.5)
ts – dimensionless time \(\bar{t}\) ,\(\bar{t}=\dfrac{t}{r/\beta}\)
azimuth – azimuth in degree
- 返回:
Normalized solution with shape (nt, 3, 3). To get the physical solution, divide by \(\pi^2 \mu r\)