atomppgen.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

Functions

check_ghost_states(inv_result, r, w, ...[, ...])

幽灵态检测(径向级别)

check_log_derivative(V_AE, V_PS, r, l, r_test)

对数导数匹配验证

check_norm_conservation(tm_result[, tolerance])

检验范数守恒条件

is_covalent(element)

判断元素是否以共价价键表现为主,影响阈值选择

run_full_validation(ae_result, tm_dict, inv_dict)

完整验证流程

Classes

GhostStateResult(method, l, eigenvalues, ...)

幽灵态检测结果

LogDerivativeResult(l, r_test, energies, ...)

对数导数匹配验证结果

NormConservationResult(l, norm_error, ...)

范数守恒检验结果

ValidationReport(norm_results, ...)

完整验证报告

class atomppgen.validate.NormConservationResult(l, norm_error, passed, rc, tolerance, diagnostics)[源代码]

范数守恒检验结果

参数:
l

角动量量子数

Type:

int

norm_error

范数误差:∫₀^rc |u_PS|² - ∫₀^rc |u_AE|²

Type:

float

passed

是否通过(abs(norm_error) < tolerance

Type:

bool

rc

截断半径(Bohr)

Type:

float

tolerance

容许误差阈值

Type:

float

diagnostics

诊断信息: - norm_ae : 全电子内区范数 - norm_ps : 伪轨道内区范数

Type:

Dict

l: int
norm_error: float
passed: bool
rc: float
tolerance: float
diagnostics: Dict
__init__(l, norm_error, passed, rc, tolerance, diagnostics)
参数:
返回类型:

None

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)[源代码]

对数导数匹配验证结果

参数:
l

角动量量子数

Type:

int

r_test

测试半径(Bohr)

Type:

float

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

zero_crossing_rms

零点均方根偏差(Hartree)

Type:

float

curve_rms_valence

价区曲线均方根差异(-0.05 ~ +0.05 Ha)

Type:

float

curve_rms_full

全能量窗口曲线 RMS(告警用)

Type:

float

passed

是否通过验证(基于 zero_crossing_rms 和 curve_rms_valence)

Type:

bool

diagnostics

诊断信息

Type:

Dict

l: int
r_test: float
energies: ndarray
L_AE: ndarray
L_PS: ndarray
zero_crossings_AE: ndarray
zero_crossings_PS: ndarray
zero_crossing_rms: float
curve_rms_valence: float
curve_rms_full: float
passed: bool
diagnostics: 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)
参数:
返回类型:

None

class atomppgen.validate.GhostStateResult(method, l, eigenvalues, known_valence, ghost_states, box_states, n_ghosts, n_box_states, passed, tail_ratios, diagnostics)[源代码]

幽灵态检测结果

参数:
method

检测方法('radial' 或 'plane_wave')

Type:

str

l

角动量量子数(径向方法)或 -1(平面波方法)

Type:

int

eigenvalues

检测到的本征值(Hartree)

Type:

np.ndarray

known_valence

已知价电子能级

Type:

np.ndarray

ghost_states

幽灵态能级(不含盒态)

Type:

np.ndarray

box_states

盒态能级(波函数尾部未衰减)

Type:

np.ndarray

n_ghosts

真幽灵态数量(不含盒态)

Type:

int

n_box_states

盒态数量

Type:

int

passed

是否通过(n_ghosts ≤ 10)

Type:

bool

tail_ratios

各本征态的尾部比例(\(|\psi(R_{\max})| / \max_r |\psi(r)|\)

Type:

np.ndarray

diagnostics

诊断信息

Type:

Dict

method: str
l: int
eigenvalues: ndarray
known_valence: ndarray
ghost_states: ndarray
box_states: ndarray
n_ghosts: int
n_box_states: int
passed: bool
tail_ratios: ndarray
diagnostics: Dict
__init__(method, l, eigenvalues, known_valence, ghost_states, box_states, n_ghosts, n_box_states, passed, tail_ratios, diagnostics)
参数:
返回类型:

None

class atomppgen.validate.ValidationReport(norm_results, log_deriv_results, ghost_result, overall_passed, diagnostics)[源代码]

完整验证报告

参数:
norm_results

各通道范数守恒结果

Type:

Dict[int, NormConservationResult]

log_deriv_results

各通道对数导数结果

Type:

Dict[int, LogDerivativeResult]

ghost_result

幽灵态检测结果

Type:

GhostStateResult

overall_passed

整体是否通过

Type:

bool

diagnostics

汇总诊断信息

Type:

Dict

norm_results: Dict[int, NormConservationResult]
log_deriv_results: Dict[int, LogDerivativeResult]
ghost_result: Optional[GhostStateResult]
overall_passed: bool
diagnostics: Dict
to_dict()[源代码]

转换为字典(用于 JSON 序列化)

返回类型:

Dict

summary()[源代码]

生成 Markdown 验证摘要,便于在报告或导出文件中引用

返回类型:

str

__init__(norm_results, log_deriv_results, ghost_result, overall_passed, diagnostics)
参数:
返回类型:

None

atomppgen.validate.check_norm_conservation(tm_result, tolerance=1e-06)[源代码]

检验范数守恒条件

验证伪轨道在截断半径内的归一化是否与全电子一致:

∫₀^rc |u_PS(r)|² dr = ∫₀^rc |u_AE(r)|² dr

参数:
  • tm_result (TMResult) -- TM 伪化结果(包含 norm_error 字段)

  • tolerance (float, default=1e-6) -- 容许误差阈值

返回:

范数守恒检验结果

返回类型:

NormConservationResult

备注

本函数实际上是对 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,默认按金属元素阈值处理。

返回:

对数导数匹配结果

返回类型:

LogDerivativeResult

备注

  1. 径向薛定谔方程: 包含离心项 l(l+1)/(2r²),在求解器内部添加

  2. 边界条件: ψ(0)=0,ψ(∞)=0 或外向波

  3. 数值方法: Numerov 或五点有限差分

  4. 能量单位: 统一使用 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)

返回:

幽灵态检测结果

返回类型:

GhostStateResult

备注

  1. **径向方法**(A 级): 对角化径向哈密顿量,查找异常束缚态

  2. **平面波方法**(B 级,可选): 小球平面波基组,包含非局域势

示例

>>> 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) -- 全电子原子解

  • tm_dict (Dict[int, TMResult]) -- 各通道 TM 伪化结果

  • 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)

返回:

完整验证报告

返回类型:

ValidationReport

示例

>>> report = run_full_validation(ae, tm_dict, inv_dict)
>>> print(report.overall_passed)
>>> print(report.to_dict())