pygrt.pymod
- file:
pymod.py
- author:
Zhu Dengda (zhudengda@mail.iggcas.ac.cn)
- date:
2024-07-24
该文件包括 Python端使用的模型 pygrt.c_structures.c_PyModel1D
- class pygrt.pymod.PyModel1D(modarr0: ndarray, depsrc: float, deprcv: float)[源代码]
基类:
object- __init__(modarr0: ndarray, depsrc: float, deprcv: float)[源代码]
将震源和台站插入定义模型的数组,转为
PyModel1D实例的形式- 参数:
modarr0 – 模型数组,每行格式为[thickness(km), Vp(km/s), Vs(km/s), Rho(g/cm^3), Qp, Qs]
depsrc – 震源深度(km)
deprcv – 台站深度(km)
allowLiquid – 是否允许液体层
- compute_travt1d(dist: float)[源代码]
调用C程序,计算初至P波和S波的走时
- 参数:
dist – 震中距
- 返回:
travtP - 初至P波走时(s)
travtS - 初至S波走时(s)
- gen_gf_spectra(*args, **kwargs)[源代码]
Bad function name, has already been removed. Use ‘compute_grn’ instead.
- compute_grn(distarr: ndarray | List[float] | float, nt: int, dt: float, upsampling_n: int = 1, freqband: ndarray | List[float] = [-1, -1], zeta: float = 0.8, keepAllFreq: bool = False, vmin_ref: float = 0.0, keps: float = -1.0, ampk: float = 1.15, k0: float = 5.0, Length: float = 0.0, filonLength: float = 0.0, safilonTol: float = 0.0, filonCut: float = 0.0, delayT0: float = 0.0, delayV0: float = 0.0, calc_upar: bool = False, gf_source=['EX', 'VF', 'HF', 'DC'], statsfile: str | None = None, statsidxs: ndarray | List[int] | None = None, print_runtime: bool = True)[源代码]
调用C库计算格林函数的主函数,以列表的形式返回,其中每个元素为对应震中距的格林函数
obspy.Stream类型。- 参数:
distarr – 多个震中距(km) 的数组, 或单个震中距的浮点数
nt – 时间点数,借助于 SciPy,nt不再要求是2的幂次
dt – 采样间隔(s)
upsampling_n – 升采样倍数
freqband – 频率范围(Hz),以此确定待计算的离散频率点
zeta – 定义虚频率的系数 \(\zeta\) , 虚频率 \(\tilde{\omega} = \omega - j*w_I, w_I = \zeta*\pi/T, T=nt*dt\) , T为时窗长度。 使用离散波数积分时为了避开附加源以及奇点的影响, (Bouchon, 1981) 在频率上添加微小虚部, 更多测试见 (张海明, 2021)
:param keepAllFreq 计算所有频点,不论频率多低 :param vmin_ref: 最小参考速度,默认vmin=max(minimum velocity, 0.1),用于定义波数积分上限,小于0则在达到积分上限后使用峰谷平均法
(默认当震源和场点深度差<=1km时自动使用峰谷平均法)
- 参数:
keps – 波数k积分收敛条件,见 (Yao and Harkrider, 1983) (初稿), 为负数代表不提前判断收敛,按照波数积分上限进行积分
ampk – 影响波数k积分上限的系数,见下方
k0 – 波数k积分的上限 \(\tilde{k_{max}}=\sqrt{(k_{0}*\pi/hs)^2 + (ampk*w/vmin_{ref})^2}\) , 波数k积分循环必须退出, hs=max(震源和台站深度差,1.0)
Length – 定义波数k积分的间隔 dk=2pi / (L*rmax), 选取要求见 (Bouchon, 1981) (张海明, 2021),默认自动选择
filonLength – Filon积分的间隔
safilonTol – 自适应Filon积分采样精度
filonCut – 波数积分和Filon积分的分割点filonCut, k*=<filonCut>/rmax
calc_upar – 是否计算位移u的空间导数
gf_source – 待计算的震源类型
statsfile – 波数k积分(包括Filon积分和峰谷平均法)的过程记录文件,常用于debug或者观察积分过程中 \(F(k,\omega)\) 和 \(F(k,\omega)J_m(kr)k\) 的变化
statsidxs – 仅输出特定频点的过程记录文件,建议给定频点,否则默认所有频率点的记录文件都输出,很占空间
print_runtime – 是否打印运行时间
- 返回:
dataLst - 列表,每个元素为
obspy.Stream类型 )
- compute_static_grn(xarr: ndarray | List[float] | float, yarr: ndarray | List[float] | float, vmin_ref: float = 0.0, keps: float = -1.0, k0: float = 5.0, Length: float = 15.0, filonLength: float = 0.0, safilonTol: float = 0.0, filonCut: float = 0.0, calc_upar: bool = False, statsfile: str | None = None)[源代码]
调用C库计算静态格林函数,以字典的形式返回
- 参数:
xarr – 北向坐标数组,或单个浮点数
yarr – 东向坐标数组,或单个浮点数
vmin_ref – 最小参考速度(具体数值不使用),小于0则在达到积分上限后使用峰谷平均法 (默认当震源和场点深度差<=0.5km时自动使用峰谷平均法)
keps – 波数k积分收敛条件,见 (Yao and Harkrider, 1983) (初稿), 为负数代表不提前判断收敛,按照波数积分上限进行积分
k0 – 波数k积分的上限 \(\tilde{k_{max}}=(k_{0}*\pi/hs)^2\) , 波数k积分循环必须退出, hs=max(震源和台站深度差,1.0)
Length – 定义波数k积分的间隔 dk=2pi / (L*rmax), 默认15;负数表示使用Filon积分
filonLength – Filon积分的间隔
safilonTol – 自适应Filon积分采样精度
filonCut – 波数积分和Filon积分的分割点filonCut, k*=<filonCut>/rmax
calc_upar – 是否计算位移u的空间导数
statsfile – 波数k积分(包括Filon积分和峰谷平均法)的过程记录文件,常用于debug或者观察积分过程中 \(F(k,\omega)\) 和 \(F(k,\omega)J_m(kr)k\) 的变化
- 返回:
dataDct - 字典形式的格林函数