atomscf.scf_hf

Hartree-Fock SCF 实现

本模块实现原子 Hartree-Fock 自洽场计算。

实现方法

采用算子-作用式构造 Fock 矩阵:

  1. 交换算子 K 作用于基向量得到矩阵元

  2. 组合局域哈密顿与交换矩阵

  3. 对角化求解占据态

物理背景

Hartree-Fock 方程(径向中心场形式):

\[\begin{split}\\left[-\\frac{1}{2}\\frac{d^2}{dr^2} + \\frac{\\ell(\\ell+1)}{2r^2} + v_{ext}(r) + v_H(r) + K_\\ell\\right] u_\\ell(r) = \\varepsilon_\\ell u_\\ell(r)\end{split}\]

其中: - v_H(r) 是 Hartree 势(局域) - K_ℓ 是交换算子(非局域)

引用

[HFMethod]

Szabo & Ostlund (1996) "Modern Quantum Chemistry" Dover, Chapter 3

Functions

run_hf_minimal(cfg)

对 H 原子运行最小 HF 自洽:交换与 Hartree 相消,仅剩外势。

run_hf_scf(cfg)

运行通用 HF SCF 计算(支持 s, p 等多 l 通道)。

run_hf_scf_s(cfg)

运行 s 轨道 HF SCF 计算。

Classes

HFConfig(Z, r, w[, tol, maxiter])

HF 最小实现配置(教学版,仅支持 H)。

HFResult(converged, iterations, eps_1s, u_1s)

HF 最小实现结果容器(H 验证)。

HFSCFConfig(Z, r, w[, n_occ, occ_nums, ...])

完整 HF SCF 配置(支持 s 轨道)。

HFSCFGeneralConfig(Z, r, w, occ_by_l, eigs_per_l)

通用 HF SCF 配置(支持 RHF/UHF 和多 l 通道)。

HFSCFGeneralResult(converged, iterations, ...)

通用 HF SCF 结果容器(支持 RHF/UHF)。

HFSCFResult(converged, iterations, ...)

HF SCF 结果容器。

class atomscf.scf_hf.HFConfig(Z, r, w, tol=1e-07, maxiter=50)[源代码]

HF 最小实现配置(教学版,仅支持 H)。

参数:
Z

原子序数(当前仅支持 1)。

Type:

int

r

径向网格。

Type:

numpy.ndarray

w

梯形权重。

Type:

numpy.ndarray

tol

收敛阈值。

Type:

float

maxiter

最大迭代数。

Type:

int

Z: int
r: ndarray
w: ndarray
tol: float = 1e-07
maxiter: int = 50
__init__(Z, r, w, tol=1e-07, maxiter=50)
参数:
返回类型:

None

class atomscf.scf_hf.HFResult(converged, iterations, eps_1s, u_1s)[源代码]

HF 最小实现结果容器(H 验证)。

参数:
converged: bool
iterations: int
eps_1s: float
u_1s: ndarray
__init__(converged, iterations, eps_1s, u_1s)
参数:
返回类型:

None

atomscf.scf_hf.run_hf_minimal(cfg)[源代码]

对 H 原子运行最小 HF 自洽:交换与 Hartree 相消,仅剩外势。

该实现用于快速验证 HF 思想与数值骨架。对多电子(如 C)不适用。

返回类型:

HFResult

参数:

cfg (HFConfig)

class atomscf.scf_hf.HFSCFConfig(Z, r, w, n_occ=1, occ_nums=None, mix_alpha=0.3, tol=1e-06, maxiter=100, initial_guess='hydrogen', delta=None, Rp=None)[源代码]

完整 HF SCF 配置(支持 s 轨道)。

参数:
Z

原子序数

Type:

int

r

径向网格

Type:

np.ndarray

w

积分权重

Type:

np.ndarray

n_occ

占据态数量(例 He: n_occ=1 表示 1s²)

Type:

int

occ_nums

各占据态的占据数(例 He: [2.0] 表示 1s²)

Type:

list[float]

mix_alpha

密度混合参数(0=不混合,1=完全更新)

Type:

float

tol

收敛阈值(波函数 RMS 差)

Type:

float

maxiter

最大 SCF 迭代数

Type:

int

initial_guess

初始猜测方式:"hydrogen" 或 "lsda"

Type:

str

delta

变量变换参数 δ(仅用于指数网格)

Type:

float | None

Rp

变量变换参数 R_p(仅用于指数网格)

Type:

float | None

Z: int
r: ndarray
w: ndarray
n_occ: int = 1
occ_nums: list[float] | None = None
mix_alpha: float = 0.3
tol: float = 1e-06
maxiter: int = 100
initial_guess: str = 'hydrogen'
delta: float | None = None
Rp: float | None = None
__init__(Z, r, w, n_occ=1, occ_nums=None, mix_alpha=0.3, tol=1e-06, maxiter=100, initial_guess='hydrogen', delta=None, Rp=None)
参数:
返回类型:

None

class atomscf.scf_hf.HFSCFResult(converged, iterations, eigenvalues, orbitals, E_total, E_kinetic, E_ext, E_hartree, E_exchange)[源代码]

HF SCF 结果容器。

参数:
converged

是否收敛

Type:

bool

iterations

实际迭代次数

Type:

int

eigenvalues

本征能量(占据态)

Type:

np.ndarray

orbitals

占据态波函数

Type:

list[np.ndarray]

E_total

总能量

Type:

float

E_kinetic

动能

Type:

float

E_ext

外势能

Type:

float

E_hartree

Hartree 能量

Type:

float

E_exchange

交换能量

Type:

float

converged: bool
iterations: int
eigenvalues: ndarray
orbitals: list[ndarray]
E_total: float
E_kinetic: float
E_ext: float
E_hartree: float
E_exchange: float
__init__(converged, iterations, eigenvalues, orbitals, E_total, E_kinetic, E_ext, E_hartree, E_exchange)
参数:
返回类型:

None

atomscf.scf_hf.run_hf_scf_s(cfg)[源代码]

运行 s 轨道 HF SCF 计算。

实现完整的 Hartree-Fock 自洽场循环: 1. 初始化波函数猜测 2. 计算 Hartree 势 3. 构造 Fock 矩阵(包含交换算子) 4. 对角化求解本征态 5. 检查收敛并混合 6. 重复直到收敛

参数:

cfg (HFSCFConfig) -- HF SCF 配置

返回:

收敛的 HF 结果

返回类型:

HFSCFResult

备注

收敛判据: 波函数的 RMS 变化 < tol

密度混合: u_new = α * u_scf + (1-α) * u_old

能量计算:
  • E_kin = Σ_i n_i <u_i| -∇²/2 |u_i>

  • E_ext = Σ_i n_i <u_i| v_ext |u_i>

  • E_H = (1/2) ∫ n(r) v_H(r) 4πr² dr

  • E_x = (1/2) Σ_i n_i <u_i| K |u_i>

  • E_total = E_kin + E_ext + E_H + E_x

示例

氢原子 HF:

>>> from atomscf.grid import radial_grid_linear, trapezoid_weights
>>> r, _ = radial_grid_linear(n=1000, rmin=1e-6, rmax=50.0)
>>> w = trapezoid_weights(r)
>>> cfg = HFSCFConfig(Z=1, r=r, w=w, n_occ=1, occ_nums=[1.0])
>>> res = run_hf_scf_s(cfg)
>>> print(f"E_1s = {res.eigenvalues[0]:.6f} Ha")  # 应约 -0.5 Ha
>>> print(f"E_total = {res.E_total:.6f} Ha")      # 应约 -0.5 Ha

参见

exchange_operator_s

s 轨道交换算子

v_hartree

Hartree 势计算

class atomscf.scf_hf.HFSCFGeneralConfig(Z, r, w, occ_by_l, eigs_per_l, spin_mode='RHF', occ_by_l_spin=None, mix_alpha=0.3, tol=1e-06, maxiter=100, delta=None, Rp=None)[源代码]

通用 HF SCF 配置(支持 RHF/UHF 和多 l 通道)。

参数:
Z

原子序数

Type:

int

r

径向网格

Type:

np.ndarray

w

积分权重

Type:

np.ndarray

occ_by_l

按 l 分组的占据数配置(RHF 模式) 格式:{l: [n_1, n_2, ...]} 表示该 l 通道的各占据态占据数 例:C (1s² 2s² 2p²) → {0: [2.0, 2.0], 1: [2.0]}

Type:

dict[int, list[float]]

eigs_per_l

每个 l 通道求解的本征态数量 例:{0: 2, 1: 1} 表示求解 2 个 s 态和 1 个 p 态

Type:

dict[int, int]

spin_mode

自旋处理模式,默认 'RHF' - 'RHF': 限制性 Hartree-Fock(闭壳层精确,开壳层近似) - 'UHF': 非限制性 Hartree-Fock(自旋极化,适用于开壳层)

Type:

str, optional

occ_by_l_spin

自旋分辨占据数配置(UHF 模式) 格式:{l: {'up': [n_1, ...], 'down': [n_1, ...]}} 例:C ³P 态 → {0: {'up': [1.0, 1.0], 'down': [1.0, 1.0]},

1: {'up': [2.0], 'down': [0.0]}}

若为 None 且 spin_mode='UHF',自动从 occ_by_l 均分生成

Type:

dict[int, dict[str, list[float]]] | None, optional

mix_alpha

密度混合参数

Type:

float

tol

收敛阈值

Type:

float

maxiter

最大迭代数

Type:

int

delta

变量变换参数(可选)

Type:

float | None

Rp

变量变换参数(可选)

Type:

float | None

Z: int
r: ndarray
w: ndarray
occ_by_l: dict[int, list[float]]
eigs_per_l: dict[int, int]
spin_mode: str = 'RHF'
occ_by_l_spin: dict[int, dict[str, list[float]]] | None = None
mix_alpha: float = 0.3
tol: float = 1e-06
maxiter: int = 100
delta: float | None = None
Rp: float | None = None
__init__(Z, r, w, occ_by_l, eigs_per_l, spin_mode='RHF', occ_by_l_spin=None, mix_alpha=0.3, tol=1e-06, maxiter=100, delta=None, Rp=None)
参数:
返回类型:

None

class atomscf.scf_hf.HFSCFGeneralResult(converged, iterations, eigenvalues_by_l, orbitals_by_l, E_total, E_kinetic, E_ext, E_hartree, E_exchange, spin_mode='RHF', eigenvalues_by_l_spin=None, orbitals_by_l_spin=None)[源代码]

通用 HF SCF 结果容器(支持 RHF/UHF)。

参数:
converged

是否收敛

Type:

bool

iterations

实际迭代次数

Type:

int

eigenvalues_by_l

按 l 分组的本征能量(RHF 模式)

Type:

dict[int, np.ndarray]

orbitals_by_l

按 l 分组的占据态波函数(RHF 模式)

Type:

dict[int, list[np.ndarray]]

eigenvalues_by_l_spin

按 (l, spin) 分组的本征能量(UHF 模式)

Type:

dict[tuple[int, str], np.ndarray] | None

orbitals_by_l_spin

按 (l, spin) 分组的占据态波函数(UHF 模式)

Type:

dict[tuple[int, str], list[np.ndarray]] | None

E_total

总能量

Type:

float

E_kinetic

动能

Type:

float

E_ext

外势能

Type:

float

E_hartree

Hartree 能量

Type:

float

E_exchange

交换能量

Type:

float

spin_mode

使用的自旋模式('RHF' 或 'UHF')

Type:

str

converged: bool
iterations: int
eigenvalues_by_l: dict[int, ndarray]
orbitals_by_l: dict[int, list[ndarray]]
E_total: float
E_kinetic: float
E_ext: float
E_hartree: float
E_exchange: float
spin_mode: str = 'RHF'
eigenvalues_by_l_spin: dict[tuple[int, str], ndarray] | None = None
orbitals_by_l_spin: dict[tuple[int, str], list[ndarray]] | None = None
__init__(converged, iterations, eigenvalues_by_l, orbitals_by_l, E_total, E_kinetic, E_ext, E_hartree, E_exchange, spin_mode='RHF', eigenvalues_by_l_spin=None, orbitals_by_l_spin=None)
参数:
返回类型:

None

atomscf.scf_hf.run_hf_scf(cfg)[源代码]

运行通用 HF SCF 计算(支持 s, p 等多 l 通道)。

实现按 l 分组的自洽场循环,每个 l 通道独立求解 Fock 方程, 但通过 Hartree 势和交换算子耦合。

参数:

cfg (HFSCFGeneralConfig) -- 通用 HF SCF 配置

返回:

收敛的 HF 结果

返回类型:

HFSCFGeneralResult

示例

碳原子 HF (1s² 2s² 2p²):

>>> cfg = HFSCFGeneralConfig(
...     Z=6,
...     r=r, w=w,
...     occ_by_l={0: [2.0, 2.0], 1: [2.0]},  # 1s², 2s², 2p²
...     eigs_per_l={0: 2, 1: 1},              # 求解 2 个 s + 1 个 p
... )
>>> res = run_hf_scf(cfg)