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)[源代码]
剪切源,角度单位均为度
- 参数:
st – 计算好的时域格林函数,
obspy.Stream类型,或者静态格林函数(字典类型)M0 – 标量地震矩, 单位dyne*cm
strike – 走向,以北顺时针为正,0<=strike<=360
dip – 倾角,以水平面往下为正,0<=dip<=90
rake – 滑动角,在断层面相对于走向方向逆时针为正,-180<=rake<=180
az – 台站方位角,以北顺时针为正,0<=az<=360(静态情况不需要)
ZNE – 是否以ZNE分量输出?
calc_upar – 是否计算位移u的空间导数
- 返回:
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)[源代码]
单力源,力分量单位均为dyne
- 参数:
st – 计算好的时域格林函数,
obspy.Stream类型,或者静态格林函数(字典类型)S – 力的放大系数
fN – 北向力,向北为正
fE – 东向力,向东为正
fZ – 垂向力,向下为正
az – 台站方位角,以北顺时针为正,0<=az<=360 (静态情况不需要)
ZNE – 是否以ZNE分量输出?
calc_upar – 是否计算位移u的空间导数
- 返回:
stream - 三分量地震图,
obspy.Stream类型
- pygrt.utils.gen_syn_from_gf_EX(st: Stream | dict, M0: float, az: float = -999, ZNE=False, calc_upar: bool = False)[源代码]
爆炸源
- 参数:
st – 计算好的时域格林函数,
obspy.Stream类型,或者静态格林函数(字典类型)M0 – 标量地震矩, 单位dyne*cm
az – 台站方位角,以北顺时针为正,0<=az<=360 [不用于计算] (静态情况不需要)
ZNE – 是否以ZNE分量输出?
calc_upar – 是否计算位移u的空间导数
- 返回:
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)[源代码]
矩张量源,单位为dyne*cm
- 参数:
st – 计算好的时域格林函数,
obspy.Stream类型,或者静态格林函数(字典类型)M0 – 标量地震矩
MT – 矩张量 (M11, M12, M13, M22, M23, M33),下标1,2,3分别代表北向,东向,垂直向下
az – 台站方位角,以北顺时针为正,0<=az<=360 (静态情况不需要)
ZNE – 是否以ZNE分量输出?
calc_upar – 是否计算位移u的空间导数
- 返回:
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)[源代码]
对stream中每一道信号做线性卷积
- 参数:
st0 – 记录多个Trace的
obspy.Stream类型signal0 – 卷积信号
inplace – 是否做原地更改
- 返回:
stream - 处理后的结果,
obspy.Stream类型
- pygrt.utils.stream_integral(st0: Stream, inplace=True)[源代码]
对stream中每一道信号做梯形积分
- 参数:
st0 – 记录多个Trace的
obspy.Stream类型inplace – 是否做原地更改
- 返回:
stream - 处理后的结果,
obspy.Stream类型
- pygrt.utils.stream_diff(st0: Stream, inplace=True)[源代码]
对stream中每一道信号做中心差分
- 参数:
st0 – 记录多个Trace的
obspy.Stream类型inplace – 是否做原地更改
- 返回:
stream - 处理后的结果,
obspy.Stream类型
- pygrt.utils.stream_write_sac(st: Stream, dir: str)[源代码]
将一系列Trace以SAC形式保存到本地,以发震时刻作为参考0时刻
- 参数:
st – 记录多个Trace的
obspy.Stream类型dir – 保存目录
- pygrt.utils.read_kernels_freqs(statsdir: str, vels: ndarray, ktypes: List[str] | None = None)[源代码]
读取statsdir目录下所有频率(除了零频)的积分过程文件(K_开头),再将核函数线性插值到指定的速度数组vels
- 参数:
statsdir – 存储积分过程文件的目录
vels – 待插值的速度数组(km/s),必须正序
ktype – 指定返回一系列的核函数名称,如EX_q,DS_w等,默认返回全部
- 返回:
kerDct - 字典格式的核函数插值结果
- pygrt.utils.read_statsfile(statsfile: str)[源代码]
读取单个频率下波数积分(或Filon积分)的记录文件
- 参数:
statsfile – 文件路径(可使用通配符简化输入)
- 返回:
data - numpy.ndarray 自定义类型数组
- pygrt.utils.read_statsfile_ptam(statsfile: str)[源代码]
读取单个频率下峰谷平均法的记录文件
- 参数:
statsfile – PTAM文件路径(可使用通配符简化输入)
- 返回:
data1 - numpy.ndarray 自定义类型数组,DWM或FIM过程中的积分过程数据
data2 - numpy.ndarray 自定义类型数组,PTAM过程中的积分过程数据
ptam_data - numpy.ndarray 自定义类型数组,PTAM的峰谷位置及幅值
dist - 文件对应的震中距(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)[源代码]
根据
read_statsfile函数函数读取的数据, 绘制核函数 \(F(k,\omega)\)、被积函数 \(F(k,\omega)J_m(kr)k\),以及简单计算累积积分 \(\sum F(k,\omega)J_m(kr)k\) 并绘制。备注
并不是每个震源类型对应的每一阶每种积分类型都存在,详见 格林函数分类。
- 参数:
statsdata –
read_statsfile函数返回值dist – 震中距(km)
srctype – 震源类型的缩写,包括EX、VF、HF、DD、DS、SS
ptype – 积分类型(0,1,2,3)
RorI – 绘制实部还是虚部,默认实部,传入2表示实部虚部都绘制
fig – 传入自定义的matplotlib.Figure对象,默认为None
axs – 传入自定义的matplotlib.Axes对象数组(三个),默认为None
- 返回:
fig - matplotlib.Figure对象
(ax1,ax2,ax3) - matplotlib.Axes对象数组
- 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)[源代码]
根据
read_statsfile_ptam函数读取的数据, 简单计算并绘制累积积分以及峰谷平均法使用的波峰波谷的位置。备注
并不是每个震源类型对应的每一阶每种积分类型都存在,详见 格林函数分类。
- 参数:
statsdata1 – DWM或FIM过程中的积分过程数据
statsdata2 – PTAM过程中的积分过程数据
statsdata_ptam – PTAM的峰谷位置及幅值
srctype – 震源类型的缩写,包括EX、VF、HF、DD、DS、SS
ptype – 积分类型(0,1,2,3)
RorI – 绘制实部还是虚部,默认实部,传入2表示实部虚部都绘制
fig – 传入自定义的matplotlib.Figure对象,默认为None
axs – 传入自定义的matplotlib.Axes对象数组(三个),默认为None
- 返回:
fig - matplotlib.Figure对象
(ax1,ax2,ax3) - matplotlib.Axes对象数组