validate - 赝势可转移性验证
本模块提供赝势质量评估工具,包括范数守恒、对数导数匹配和幽灵态检测。
赝势可转移性验证模块
提供三类验证功能: 1. 范数守恒检验:伪轨道与全电子轨道在截断半径内的归一化一致性 2. 对数导数匹配:散射性质的能量依赖性(AE vs PS) 3. 幽灵态检测:赝势哈密顿量的病态束缚态检查
主要函数
check_norm_conservation : 范数守恒检验 check_log_derivative : 对数导数匹配验证 check_ghost_states : 幽灵态检测(径向级别) run_full_validation : 完整验证流程
验证阈值说明
- 范数守恒误差:
norm_error < 1e-6(所有元素统一标准)
- 对数导数曲线 RMS:
金属元素(Al, Na, Mg): curve_rms_valence < 16.0
共价元素(Si, C, N): curve_rms_valence < 3.0
物理依据: 金属元素在远离核区(r ~ r_c)的软势中,对数导数 L(E,r) = r·ψ'/ψ 对能量变化不敏感,全电子与赝势的相位差异在过渡区被放大。这是固有 特性而非赝势质量缺陷。共价元素的波函数节点清晰、相位变化更敏感, 因而使用更严格的阈值。对教学实现而言,数值求解与采样点数会放大 RMS 波动,因此采用 3.0 作为“可复现且具区分度”的实践阈值。
- 幽灵态判定:
仅统计能量高于最高价态 0.1 Ha 以上的幽灵态(n_ghosts_above_valence = 0)
理由: TM 伪化在 r < r_c 创造的浅束缚态(能量低于价态)对基态 DFT 计算 影响可忽略。若需高精度激发态计算,建议使用 RRKJ 方法或增加角动量 通道数以抑制幽灵。
- 零点 RMS:
zero_crossing_rms < 0.025 Ha(所有元素统一标准)
这是对数导数曲线零点位置的 RMS 偏差,反映散射相移的准确性。
参考文献
Troullier & Martins, PRB 43, 1993 (1991) - 范数守恒条件 Gonze et al., Comput. Mater. Sci. 25, 478 (2002) - 对数导数方法 Rappe et al., PRB 41, 1227 (1990) - 幽灵态检测 Hamann, PRB 88, 085117 (2013) - ONCVPSP 验证标准 van Setten et al., Comput. Phys. Commun. 226, 39 (2018) - PseudoDojo
- class atomppgen.validate.NormConservationResult(l, norm_error, passed, rc, tolerance, diagnostics)[源代码]
基类:
object范数守恒检验结果
- diagnostics
诊断信息: - norm_ae : 全电子内区范数 - norm_ps : 伪轨道内区范数
- Type:
Dict
- class atomppgen.validate.LogDerivativeResult(l, r_test, energies, L_AE, L_PS, zero_crossings_AE, zero_crossings_PS, zero_crossing_rms, curve_rms_valence, curve_rms_full, passed, diagnostics)[源代码]
基类:
object对数导数匹配验证结果
- 参数:
- energies
能量网格(Hartree)
- Type:
np.ndarray
- L_AE
全电子对数导数 L(E) = r·ψ'/ψ
- Type:
np.ndarray
- L_PS
伪势对数导数
- Type:
np.ndarray
- zero_crossings_AE
全电子零点能量
- Type:
np.ndarray
- zero_crossings_PS
伪势零点能量
- Type:
np.ndarray
- diagnostics
诊断信息
- Type:
Dict
- __init__(l, r_test, energies, L_AE, L_PS, zero_crossings_AE, zero_crossings_PS, zero_crossing_rms, curve_rms_valence, curve_rms_full, passed, diagnostics)
- class atomppgen.validate.GhostStateResult(method, l, eigenvalues, known_valence, ghost_states, box_states, n_ghosts, n_box_states, passed, tail_ratios, diagnostics)[源代码]
基类:
object幽灵态检测结果
- 参数:
- eigenvalues
检测到的本征值(Hartree)
- Type:
np.ndarray
- known_valence
已知价电子能级
- Type:
np.ndarray
- ghost_states
幽灵态能级(不含盒态)
- Type:
np.ndarray
- box_states
盒态能级(波函数尾部未衰减)
- Type:
np.ndarray
- tail_ratios
各本征态的尾部比例(\(|\psi(R_{\max})| / \max_r |\psi(r)|\))
- Type:
np.ndarray
- diagnostics
诊断信息
- Type:
Dict
- __init__(method, l, eigenvalues, known_valence, ghost_states, box_states, n_ghosts, n_box_states, passed, tail_ratios, diagnostics)
- class atomppgen.validate.ValidationReport(norm_results, log_deriv_results, ghost_result, overall_passed, diagnostics)[源代码]
基类:
object完整验证报告
- 参数:
norm_results (Dict[int, NormConservationResult])
log_deriv_results (Dict[int, LogDerivativeResult])
ghost_result (GhostStateResult | None)
overall_passed (bool)
diagnostics (Dict)
- norm_results
各通道范数守恒结果
- Type:
Dict[int, NormConservationResult]
- log_deriv_results
各通道对数导数结果
- Type:
Dict[int, LogDerivativeResult]
- ghost_result
幽灵态检测结果
- Type:
- diagnostics
汇总诊断信息
- Type:
Dict
-
norm_results:
Dict[int,NormConservationResult]
-
log_deriv_results:
Dict[int,LogDerivativeResult]
-
ghost_result:
Optional[GhostStateResult]
- __init__(norm_results, log_deriv_results, ghost_result, overall_passed, diagnostics)
- 参数:
norm_results (Dict[int, NormConservationResult])
log_deriv_results (Dict[int, LogDerivativeResult])
ghost_result (GhostStateResult | None)
overall_passed (bool)
diagnostics (Dict)
- 返回类型:
None
- atomppgen.validate.check_norm_conservation(tm_result, tolerance=1e-06)[源代码]
检验范数守恒条件
验证伪轨道在截断半径内的归一化是否与全电子一致:
- 参数:
- 返回:
范数守恒检验结果
- 返回类型:
备注
本函数实际上是对 TMResult.norm_error 的包装,因为 TM 伪化过程 已经通过非线性方程组保证了范数守恒。此处仅验证误差是否在容许范围内。
示例
>>> from atomppgen import solve_ae_atom, tm_pseudize, check_norm_conservation >>> ae = solve_ae_atom(Z=13, spin_mode='LDA', lmax=0) >>> tm = tm_pseudize(ae.r, ae.w, ae.u_by_l[0][-1], ae.eps_by_l[0][-1], l=0, rc=2.0) >>> result = check_norm_conservation(tm) >>> print(result.passed) # True >>> print(result.norm_error) # < 1e-6
- atomppgen.validate.check_log_derivative(V_AE, V_PS, r, l, r_test, E_range_Ha=(-0.25, 0.25), E_step_Ha=0.025, element=None)[源代码]
对数导数匹配验证
在测试半径 r_test 处,扫描能量窗口,比较全电子和伪势的对数导数:
L(E, r) = r · d ln ψ(r) / dr = r · ψ'(r) / ψ(r)
评价指标: 1. 零点能量均方根偏差:ΔE_RMS < 0.025 Ha (≈0.05 Ry) 2. 价区曲线均方根差异:金属元素 < 16.0,共价元素 < 3.0
- 参数:
V_AE (np.ndarray) -- 全电子半局域势(Hartree),不含离心项
V_PS (np.ndarray) -- 伪势半局域势(Hartree),不含离心项
r (np.ndarray) -- 径向网格(Bohr)
l (int) -- 角动量量子数
r_test (float) -- 测试半径(Bohr),建议 max(rc_l) + 0.5
E_range_Ha (tuple, default=(-0.25, 0.25)) -- 能量扫描范围(Hartree),对应 Ry 的 (-0.5, 0.5)
E_step_Ha (float, default=0.025) -- 能量步长(Hartree),对应 Ry 的 0.05
element (str, optional) -- 元素符号(如
'Al'、'Si'),用于自动选择金属/共价元素的曲线 RMS 阈值。 若为None,默认按金属元素阈值处理。
- 返回:
对数导数匹配结果
- 返回类型:
备注
径向薛定谔方程: 包含离心项 l(l+1)/(2r²),在求解器内部添加
边界条件: ψ(0)=0,ψ(∞)=0 或外向波
数值方法: Numerov 或五点有限差分
能量单位: 统一使用 Hartree
示例
>>> result = check_log_derivative(V_AE, V_PS, r, l=0, r_test=3.0) >>> print(result.zero_crossing_rms) # < 0.025 Ha >>> print(result.passed) # True
- atomppgen.validate.check_ghost_states(inv_result, r, w, valence_energy, E_window_Ha=(-0.25, 0.25), method='radial', radial_grid_n=None)[源代码]
幽灵态检测(径向级别)
检查赝势径向哈密顿量 H_l = T + V_PS(r) + l(l+1)/(2r²) 在能量窗口内 是否有额外的病态束缚态。
- 参数:
inv_result (InvertResult) -- 半局域势反演结果
r (np.ndarray) -- 径向网格(Bohr)
w (np.ndarray) -- 积分权重
valence_energy (float) -- 已知价电子能级(Hartree)
E_window_Ha (tuple, default=(-0.25, 0.25)) -- 能量窗口(Hartree)
method (str, default='radial') -- 检测方法('radial' 或 'plane_wave')
radial_grid_n (int | None)
- 返回:
幽灵态检测结果
- 返回类型:
示例
>>> result = check_ghost_states(inv_result, r, w, valence_energy=-0.5) >>> print(result.n_ghosts) # 0 >>> print(result.passed) # True
- atomppgen.validate.run_full_validation(ae_result, tm_dict, inv_dict, r_test=3.0, E_range_Ry=(-0.5, 0.5), E_step_Ry=0.05, ghost_radial_grid_n=None)[源代码]
完整验证流程
对所有通道执行范数守恒、对数导数匹配和幽灵态检测。
- 参数:
ae_result (AEAtomResult) -- 全电子原子解
inv_dict (Dict[int, InvertResult]) -- 各通道势反演结果
r_test (float, default=3.0) -- 对数导数测试半径(Bohr)
E_range_Ry (tuple, default=(-0.5, 0.5)) -- 能量窗口(Rydberg)
E_step_Ry (float, default=0.05) -- 能量步长(Rydberg)
ghost_radial_grid_n (int | None)
- 返回:
完整验证报告
- 返回类型:
示例
>>> report = run_full_validation(ae, tm_dict, inv_dict) >>> print(report.overall_passed) >>> print(report.to_dict())
主要函数
范数守恒检验
- atomppgen.validate.check_norm_conservation(tm_result, tolerance=1e-06)[源代码]
检验范数守恒条件
验证伪轨道在截断半径内的归一化是否与全电子一致:
- 参数:
- 返回:
范数守恒检验结果
- 返回类型:
备注
本函数实际上是对 TMResult.norm_error 的包装,因为 TM 伪化过程 已经通过非线性方程组保证了范数守恒。此处仅验证误差是否在容许范围内。
示例
>>> from atomppgen import solve_ae_atom, tm_pseudize, check_norm_conservation >>> ae = solve_ae_atom(Z=13, spin_mode='LDA', lmax=0) >>> tm = tm_pseudize(ae.r, ae.w, ae.u_by_l[0][-1], ae.eps_by_l[0][-1], l=0, rc=2.0) >>> result = check_norm_conservation(tm) >>> print(result.passed) # True >>> print(result.norm_error) # < 1e-6
对数导数匹配
- atomppgen.validate.check_log_derivative(V_AE, V_PS, r, l, r_test, E_range_Ha=(-0.25, 0.25), E_step_Ha=0.025, element=None)[源代码]
对数导数匹配验证
在测试半径 r_test 处,扫描能量窗口,比较全电子和伪势的对数导数:
L(E, r) = r · d ln ψ(r) / dr = r · ψ'(r) / ψ(r)
评价指标: 1. 零点能量均方根偏差:ΔE_RMS < 0.025 Ha (≈0.05 Ry) 2. 价区曲线均方根差异:金属元素 < 16.0,共价元素 < 3.0
- 参数:
V_AE (np.ndarray) -- 全电子半局域势(Hartree),不含离心项
V_PS (np.ndarray) -- 伪势半局域势(Hartree),不含离心项
r (np.ndarray) -- 径向网格(Bohr)
l (int) -- 角动量量子数
r_test (float) -- 测试半径(Bohr),建议 max(rc_l) + 0.5
E_range_Ha (tuple, default=(-0.25, 0.25)) -- 能量扫描范围(Hartree),对应 Ry 的 (-0.5, 0.5)
E_step_Ha (float, default=0.025) -- 能量步长(Hartree),对应 Ry 的 0.05
element (str, optional) -- 元素符号(如
'Al'、'Si'),用于自动选择金属/共价元素的曲线 RMS 阈值。 若为None,默认按金属元素阈值处理。
- 返回:
对数导数匹配结果
- 返回类型:
备注
径向薛定谔方程: 包含离心项 l(l+1)/(2r²),在求解器内部添加
边界条件: ψ(0)=0,ψ(∞)=0 或外向波
数值方法: Numerov 或五点有限差分
能量单位: 统一使用 Hartree
示例
>>> result = check_log_derivative(V_AE, V_PS, r, l=0, r_test=3.0) >>> print(result.zero_crossing_rms) # < 0.025 Ha >>> print(result.passed) # True
幽灵态检测
- atomppgen.validate.check_ghost_states(inv_result, r, w, valence_energy, E_window_Ha=(-0.25, 0.25), method='radial', radial_grid_n=None)[源代码]
幽灵态检测(径向级别)
检查赝势径向哈密顿量 H_l = T + V_PS(r) + l(l+1)/(2r²) 在能量窗口内 是否有额外的病态束缚态。
- 参数:
inv_result (InvertResult) -- 半局域势反演结果
r (np.ndarray) -- 径向网格(Bohr)
w (np.ndarray) -- 积分权重
valence_energy (float) -- 已知价电子能级(Hartree)
E_window_Ha (tuple, default=(-0.25, 0.25)) -- 能量窗口(Hartree)
method (str, default='radial') -- 检测方法('radial' 或 'plane_wave')
radial_grid_n (int | None)
- 返回:
幽灵态检测结果
- 返回类型:
示例
>>> result = check_ghost_states(inv_result, r, w, valence_energy=-0.5) >>> print(result.n_ghosts) # 0 >>> print(result.passed) # True
完整验证流程
- atomppgen.validate.run_full_validation(ae_result, tm_dict, inv_dict, r_test=3.0, E_range_Ry=(-0.5, 0.5), E_step_Ry=0.05, ghost_radial_grid_n=None)[源代码]
完整验证流程
对所有通道执行范数守恒、对数导数匹配和幽灵态检测。
- 参数:
ae_result (AEAtomResult) -- 全电子原子解
inv_dict (Dict[int, InvertResult]) -- 各通道势反演结果
r_test (float, default=3.0) -- 对数导数测试半径(Bohr)
E_range_Ry (tuple, default=(-0.5, 0.5)) -- 能量窗口(Rydberg)
E_step_Ry (float, default=0.05) -- 能量步长(Rydberg)
ghost_radial_grid_n (int | None)
- 返回:
完整验证报告
- 返回类型:
示例
>>> report = run_full_validation(ae, tm_dict, inv_dict) >>> print(report.overall_passed) >>> print(report.to_dict())
数据类
- class atomppgen.validate.NormConservationResult(l, norm_error, passed, rc, tolerance, diagnostics)[源代码]
范数守恒检验结果
- diagnostics
诊断信息: - norm_ae : 全电子内区范数 - norm_ps : 伪轨道内区范数
- Type:
Dict
- class atomppgen.validate.LogDerivativeResult(l, r_test, energies, L_AE, L_PS, zero_crossings_AE, zero_crossings_PS, zero_crossing_rms, curve_rms_valence, curve_rms_full, passed, diagnostics)[源代码]
对数导数匹配验证结果
- 参数:
- energies
能量网格(Hartree)
- Type:
np.ndarray
- L_AE
全电子对数导数 L(E) = r·ψ'/ψ
- Type:
np.ndarray
- L_PS
伪势对数导数
- Type:
np.ndarray
- zero_crossings_AE
全电子零点能量
- Type:
np.ndarray
- zero_crossings_PS
伪势零点能量
- Type:
np.ndarray
- diagnostics
诊断信息
- Type:
Dict
- __init__(l, r_test, energies, L_AE, L_PS, zero_crossings_AE, zero_crossings_PS, zero_crossing_rms, curve_rms_valence, curve_rms_full, passed, diagnostics)
- class atomppgen.validate.GhostStateResult(method, l, eigenvalues, known_valence, ghost_states, box_states, n_ghosts, n_box_states, passed, tail_ratios, diagnostics)[源代码]
幽灵态检测结果
- 参数:
- eigenvalues
检测到的本征值(Hartree)
- Type:
np.ndarray
- known_valence
已知价电子能级
- Type:
np.ndarray
- ghost_states
幽灵态能级(不含盒态)
- Type:
np.ndarray
- box_states
盒态能级(波函数尾部未衰减)
- Type:
np.ndarray
- tail_ratios
各本征态的尾部比例(\(|\psi(R_{\max})| / \max_r |\psi(r)|\))
- Type:
np.ndarray
- diagnostics
诊断信息
- Type:
Dict
- __init__(method, l, eigenvalues, known_valence, ghost_states, box_states, n_ghosts, n_box_states, passed, tail_ratios, diagnostics)
- class atomppgen.validate.ValidationReport(norm_results, log_deriv_results, ghost_result, overall_passed, diagnostics)[源代码]
完整验证报告
- 参数:
norm_results (Dict[int, NormConservationResult])
log_deriv_results (Dict[int, LogDerivativeResult])
ghost_result (GhostStateResult | None)
overall_passed (bool)
diagnostics (Dict)
- norm_results
各通道范数守恒结果
- Type:
Dict[int, NormConservationResult]
- log_deriv_results
各通道对数导数结果
- Type:
Dict[int, LogDerivativeResult]
- ghost_result
幽灵态检测结果
- Type:
- diagnostics
汇总诊断信息
- Type:
Dict
-
norm_results:
Dict[int,NormConservationResult]
-
log_deriv_results:
Dict[int,LogDerivativeResult]
-
ghost_result:
Optional[GhostStateResult]
- __init__(norm_results, log_deriv_results, ghost_result, overall_passed, diagnostics)
- 参数:
norm_results (Dict[int, NormConservationResult])
log_deriv_results (Dict[int, LogDerivativeResult])
ghost_result (GhostStateResult | None)
overall_passed (bool)
diagnostics (Dict)
- 返回类型:
None
使用示例
基础验证流程
from atomppgen import solve_ae_atom, tm_pseudize, invert_semilocal_potential
from atomppgen.validate import run_full_validation
# 1. 生成 Al 的全电子解和伪化
ae = solve_ae_atom(Z=13, spin_mode='LDA', lmax=0, grid_params={'n': 600})
tm_s = tm_pseudize(ae.r, ae.w, ae.u_by_l[0][-1], ae.eps_by_l[0][-1],
l=0, rc=2.0)
inv_s = invert_semilocal_potential(tm_s, ae.r)
# 2. 构建通道字典
tm_dict = {0: tm_s}
inv_dict = {0: inv_s}
# 3. 完整验证
report = run_full_validation(
ae, tm_dict, inv_dict,
r_test=3.0, # 测试半径 (Bohr)
E_range_Ry=(-0.5, 0.5), # 能量窗口 (Rydberg)
E_step_Ry=0.05 # 能量步长 (Rydberg)
)
# 4. 查看结果
print(f"整体通过: {report.overall_passed}")
print(f"范数守恒 (s): {report.norm_results[0].passed}")
print(f"对数导数 (s): {report.log_deriv_results[0].passed}")
if report.ghost_result:
print(f"幽灵态数量: {report.ghost_result.n_ghosts}")
# 5. 导出 JSON 报告
import json
with open('validation_report.json', 'w') as f:
json.dump(report.to_dict(), f, indent=2)
范数守恒单项检验
from atomppgen import solve_ae_atom, tm_pseudize
from atomppgen.validate import check_norm_conservation
ae = solve_ae_atom(Z=13, spin_mode='LDA', lmax=0)
tm = tm_pseudize(ae.r, ae.w, ae.u_by_l[0][-1], ae.eps_by_l[0][-1],
l=0, rc=2.0)
# 范数守恒检验
result = check_norm_conservation(tm, tolerance=1e-6)
print(f"范数误差: {result.norm_error:.3e}")
print(f"通过检验: {result.passed}")
print(f"截断半径: {result.rc} Bohr")
print(f"求解器收敛: {result.diagnostics['solver_converged']}")
对数导数匹配单项检验
from atomppgen.validate import check_log_derivative, _extract_ks_potential
# 提取 KS 有效势和伪势
V_AE = _extract_ks_potential(ae)
V_PS = inv.V_l
# 对数导数匹配(缩小参数加速示例)
result = check_log_derivative(
V_AE, V_PS, ae.r,
l=0, # s 轨道
r_test=3.0, # 测试半径 (Bohr)
E_range_Ha=(-0.2, 0.2), # 能量窗口 (Hartree)
E_step_Ha=0.05 # 能量步长 (Hartree)
)
print(f"零点 RMS: {result.zero_crossing_rms:.6f} Ha")
print(f"曲线 RMS: {result.curve_rms:.6f}")
print(f"通过检验: {result.passed}")
# 绘制对数导数曲线
import matplotlib.pyplot as plt
plt.plot(result.energies, result.L_AE, label='AE')
plt.plot(result.energies, result.L_PS, label='PS', linestyle='--')
plt.axhline(0, color='k', linewidth=0.5)
plt.xlabel('Energy (Ha)')
plt.ylabel('L(E, r_test)')
plt.legend()
plt.savefig('log_derivative.png')
幽灵态检测单项检验
from atomppgen.validate import check_ghost_states
# 幽灵态检测 (径向方法)
result = check_ghost_states(
inv, ae.r, ae.w,
valence_energy=tm.eps,
E_window_Ha=(-0.25, 0.25),
method='radial'
)
print(f"检测到幽灵态数量: {result.n_ghosts}")
print(f"能量窗口内束缚态数: {result.diagnostics['n_bound_states_in_window']}")
print(f"通过检验: {result.passed}")
# 查看所有束缚态
if len(result.eigenvalues) > 0:
print("窗口内本征值 (Ha):")
for E in result.eigenvalues:
print(f" {E:.6f}")
多通道验证
# 生成 s, p, d 三通道
ae = solve_ae_atom(Z=13, spin_mode='LDA', lmax=2)
# s 通道
tm_s = tm_pseudize(ae.r, ae.w, ae.u_by_l[0][-1], ae.eps_by_l[0][-1],
l=0, rc=2.0)
inv_s = invert_semilocal_potential(tm_s, ae.r)
# p 通道 (需调参)
tm_p = tm_pseudize(ae.r, ae.w, ae.u_by_l[1][-1], ae.eps_by_l[1][-1],
l=1, rc=1.9, continuity_orders=2)
inv_p = invert_semilocal_potential(tm_p, ae.r)
# d 通道
tm_d = tm_pseudize(ae.r, ae.w, ae.u_by_l[2][-1], ae.eps_by_l[2][-1],
l=2, rc=2.4)
inv_d = invert_semilocal_potential(tm_d, ae.r)
# 完整验证
tm_dict = {0: tm_s, 1: tm_p, 2: tm_d}
inv_dict = {0: inv_s, 1: inv_p, 2: inv_d}
report = run_full_validation(ae, tm_dict, inv_dict)
# 逐通道查看结果
for l in [0, 1, 2]:
label = ['s', 'p', 'd'][l]
norm_pass = report.norm_results[l].passed
ld_pass = report.log_deriv_results[l].passed
print(f"{label} 通道: 范数={norm_pass}, 对数导数={ld_pass}")
技术细节
能量单位转换
所有内部计算使用 Hartree (Ha) 原子单位。外部 API 支持 Rydberg (Ry) 输入, 自动转换:
# 用户指定 Rydberg
E_range_Ry = (-0.5, 0.5)
E_step_Ry = 0.05
# 内部转换为 Hartree
E_range_Ha = (E_range_Ry[0] / 2, E_range_Ry[1] / 2) # (-0.25, 0.25) Ha
E_step_Ha = E_step_Ry / 2 # 0.025 Ha
Numerov 求解器
对于非均匀网格(如 exp_transformed),自动重采样到等距网格:
# 检查网格均匀性
dr = np.diff(r)
is_uniform = np.allclose(dr, dr[0], rtol=1e-6)
if not is_uniform:
# 重采样到均匀网格
n_uniform = len(r)
r_uniform = np.linspace(r[0], r[-1], n_uniform)
V_uniform = interp1d(r, V, kind='cubic')(r_uniform)
# 在均匀网格上求解...
# 结果插值回原网格
对数导数零点匹配
零点通过符号变化检测,线性插值精确化:
def find_zero_crossings(E, L):
zeros = []
for i in range(len(L) - 1):
if L[i] * L[i+1] < 0: # 符号变化
# 线性插值
E_zero = E[i] - L[i] * (E[i+1] - E[i]) / (L[i+1] - L[i])
zeros.append(E_zero)
return zeros
幽灵态判定准则
束缚态筛选条件:
能量 \(E < E_{\\text{max}}\) (默认 0 Ha)
边界衰减:\(|\\psi(r_{\\text{max}})| < 0.1 \\cdot \\max|\\psi|\)
幽灵态判定:
在能量窗口 \([-0.25, +0.25]\) Ha 内
距已知价态 \(>0.1\) Ha 的额外束缚态
性能优化建议
对数导数扫描:粗扫 (0.1 Ry) 定位零点,细扫 (0.025 Ry) 精确化
幽灵态检测:限制网格点数 (≤300) 加速对角化
多通道并行:各通道独立,可并行计算(当前串行)
参考算法文档
详细数学推导和原理见 赝势可转移性验证方法。