势反演模块 (invert)

概述

从 Troullier-Martins 伪轨道反演半局域势 \(V_l(r)\)

核心功能

atomppgen.invert.invert_semilocal_potential(tm_result, r, u_ps=None, node_tol=1e-10, V_max_clip=1000.0, smooth_rc=False, smooth_width=0.1)[源代码]

从 TM 伪轨道反演半局域势

使用径向 Schrödinger 方程反演:

V_l(r) = ε + (1/2) · u''(r)/u(r) - l(l+1)/(2r²)

内区(r ≤ rc)使用 TM 解析导数,外区使用样条法导数。

参数:
  • tm_result (TMResult) -- TM 伪化结果(包含 a_coeff, rc, eps, l)

  • r (np.ndarray) -- 径向网格(Bohr)

  • u_ps (np.ndarray, optional) -- 伪轨道,若不提供则从 tm_result.u_ps 获取

  • node_tol (float, default=1e-10) -- 节点检测阈值(|u| < node_tol 视为节点)

  • V_max_clip (float, default=1000.0) -- 势的裁剪上限(防止除零发散)

  • smooth_rc (bool, default=False) -- 是否在 rc 处应用平滑过渡(内外区势值混合)

  • smooth_width (float, default=0.1) -- 平滑区域半宽度(Bohr),仅在 smooth_rc=True 时有效 平滑范围为 [rc-smooth_width, rc+smooth_width]

返回:

包含反演势和诊断信息

返回类型:

InvertResult

抛出:

ValueError -- 如果网格与 tm_result 不匹配

备注

  1. 内区(r ≤ rc):使用 TM 解析导数,从 _eval_tm_at_rc() 获取 u, u', u''

  2. 外区(r > rc):使用 eval_derivatives_at() 样条法

  3. 节点保护:在 |u| < node_tol 附近使用线性插值填充势值

  4. 离心势项 l(l+1)/(2r²) 在原点附近奇异,需特殊处理

  5. **平滑过渡**(可选):若 smooth_rc=True,在 rc±smooth_width 区域使用 三次样条平滑内外区势值,适用于特殊 rc/网格组合导致跳变较大的情况。 默认关闭(测试显示大多数情况 rc 处相对跳变 < 1%)。

示例

>>> from atomppgen import solve_ae_atom, tm_pseudize
>>> ae = solve_ae_atom(Z=13, spin_mode='LDA', lmax=0)
>>> tm = tm_pseudize(r=ae.r, w=ae.w, u_ae=ae.u_by_l[0][2],
...                   eps=ae.eps_by_l[0][2], l=0, rc=2.0)
>>> inv = invert_semilocal_potential(tm, ae.r)
>>> print(inv.diagnostics['V_at_rc'])  # rc 处势值

数据类

class atomppgen.invert.InvertResult(V_l, r, l, rc, eps, diagnostics)[源代码]

势反演结果

参数:
V_l

半局域势 V_l(r)(Hartree)

Type:

np.ndarray

r

径向网格(Bohr)

Type:

np.ndarray

l

角动量量子数

Type:

int

rc

伪化半径(Bohr)

Type:

float

eps

伪化能量(Hartree)

Type:

float

diagnostics

诊断信息: - n_nodes : 节点数 - V_max : 最大势值 - V_min : 最小势值 - V_at_rc : rc 处势值 - method_inner : 内区方法 ('analytical') - method_outer : 外区方法 ('spline')

Type:

Dict

V_l: ndarray
r: ndarray
l: int
rc: float
eps: float
diagnostics: Dict
__init__(V_l, r, l, rc, eps, diagnostics)
参数:
返回类型:

None

内部函数

这些函数通常不需要直接调用,由 invert_semilocal_potential 内部使用。

atomppgen.invert._invert_inner_analytical(r, u, a_coeff, l, eps, node_tol, V_max_clip)[源代码]

内区使用 TM 解析导数反演势

对于 TM 形式 u = r^{l+1} exp(p(r)),可以解析计算:

u' = [(l+1)/r + p'] u u'' = [l(l+1)/r² + 2(l+1)p'/r + p'² + p''] u

因此:

u''/u = l(l+1)/r² + 2(l+1)p'/r + p'² + p''

代入:
V_l = ε + (1/2) u''/u - l(l+1)/(2r²)

= ε + (1/2)[l(l+1)/r² + 2(l+1)p'/r + p'² + p''] - l(l+1)/(2r²) = ε + (l+1)p'/r + (1/2)(p'² + p'')

参数:
  • r (np.ndarray) -- 内区网格

  • u (np.ndarray) -- 内区伪轨道

  • a_coeff (np.ndarray) -- TM 系数 [a_0, a_2, a_4, ...]

  • l (int) -- 角动量

  • eps (float) -- 伪化能量

  • node_tol (float) -- 节点容忍度

  • V_max_clip (float) -- 势裁剪上限

返回:

内区势 V_l(r)

返回类型:

np.ndarray

atomppgen.invert._invert_outer_spline(r_full, u_full, i_start, l, eps, node_tol, V_max_clip)[源代码]

外区使用样条法导数反演势

V_l(r) = ε + (1/2) u''/u - l(l+1)/(2r²)

参数:
  • r_full (np.ndarray) -- 完整网格(包括内外区)

  • u_full (np.ndarray) -- 完整伪轨道(包括内外区)

  • i_start (int) -- 外区起始索引

  • l (int) -- 角动量

  • eps (float) -- 伪化能量

  • node_tol (float) -- 节点容忍度

  • V_max_clip (float) -- 势裁剪上限

返回:

外区势 V_l(r)

返回类型:

np.ndarray

atomppgen.invert._interpolate_nodes(r, u, V_l, node_tol)[源代码]

在节点附近插值势值

参数:
  • r (np.ndarray) -- 网格

  • u (np.ndarray) -- 轨道

  • V_l (np.ndarray) -- 势(节点处为占位值)

  • node_tol (float) -- 节点容忍度

返回:

插值后的势

返回类型:

np.ndarray

atomppgen.invert._count_nodes(u, node_tol)[源代码]

计算轨道的节点数(符号变化次数)

参数:
  • u (np.ndarray) -- 轨道

  • node_tol (float) -- 容忍度

返回:

节点数

返回类型:

int