qstack.spahm.guesses

Initial guess Hamiltonian methods for SPAHM.

Implements various guess methods: Hcore, Hückel, GWH, SAD, SAP, LB2020.

Provides:

guesses_dict: Dictionary mapping guess names to functions.

qstack.spahm.guesses.GWH(mol, *_)[source]

Compute guess Hamiltonian using Generalized Wolfsberg-Helmholtz (GWH) method.

Uses the formula: H_ij = 0.5 * K * (H_ii + H_jj) * S_ij with K = 1.75.

Reference:

M. Wolfsberg, L. Helmholtz, “The spectra and electronic structure of the tetrahedral ions MnO4-, CrO4–, and ClO4-“, J. Chem. Phys. 20 837-843 (1952), doi:10.1063/1.1700580

Parameters:
  • mol (pyscf Mole) – pyscf Mole object.

  • *_ – Unused positional arguments (for interface compatibility).

Returns:

2D GWH Hamiltonian matrix in AO basis.

Return type:

numpy ndarray

qstack.spahm.guesses.LB(mol, *_)[source]

Compute guess Hamiltonian using Laikov-Briling 2020 model with HF parameters.

Reference:

D. N. Laikov, K. R. Briling, “Atomic effective potentials for starting molecular electronic structure calculations”, Theor. Chem. Acc. 139, 17 (2020), doi:10.1007/s00214-019-2521-3

Parameters:
  • mol (pyscf Mole) – pyscf Mole object.

  • *_ – Unused positional arguments (for interface compatibility).

Returns:

2D effective Hamiltonian matrix from LB2020 model in AO basis.

Return type:

numpy ndarray

qstack.spahm.guesses.LB_HFS(mol, *_)[source]

Compute guess Hamiltonian using Laikov-Briling 2020 model with HFS parameters.

Reference:

D. N. Laikov, K. R. Briling, “Atomic effective potentials for starting molecular electronic structure calculations”, Theor. Chem. Acc. 139, 17 (2020), doi:10.1007/s00214-019-2521-3

Parameters:
  • mol (pyscf Mole) – pyscf Mole object.

  • *_ – Unused positional arguments (for interface compatibility).

Returns:

2D effective Hamiltonian matrix from LB2020-HFS model in AO basis.

Return type:

numpy ndarray

qstack.spahm.guesses.LB_grad(mf)[source]

Return Laikov-Briling Hamiltonian gradient generator function.

Combines core Hamiltonian gradient with LB2020 model gradient.

Parameters:

mf – Mean-field object with hcore_generator method.

Returns:

Function that returns total Hamiltonian gradient for a given atom.

Return type:

callable

qstack.spahm.guesses.SAD(mol, xc)[source]

Compute guess Hamiltonian using Superposition of Atomic Densities (SAD).

Constructs the Fock matrix from atomic Hartree-Fock density matrices summed together as an initial guess for molecular calculations.

References

J. Almlöf, K. Faegri Jr, K. Korsell, “Principles for a direct SCF approach to LICAO–MO ab-initio calculations”, J. Comput. Chem. 3, 385–399 (1982), doi:10.1002/jcc.540030314

L. Amat, R. Carbó-Dorca, “Use of promolecular ASA density functions as a general algorithm to obtain starting MO in SCF calculations”, Int. J. Quantum Chem. 87, 59–67 (2001), doi:10.1002/qua.10068

J. H. Van Lenthe, R. Zwaans, H. J. J. Van Dam and M. F. Guest, “Starting SCF calculations by superposition of atomic densities”, J. Comput. Chem. 27, 926–932 (2006), doi:10.1002/jcc.20393

Parameters:
  • mol (pyscf Mole) – pyscf Mole object.

  • xc (str) – Exchange-correlation functional.

Returns:

2D Fock matrix in AO basis computed from SAD.

Return type:

numpy ndarray

Warns:

RuntimeWarning – If alpha and beta effective potentials differ for the functional.

qstack.spahm.guesses.SAP(mol, *_)[source]

Compute guess Hamiltonian using Superposition of Atomic Potentials (SAP).

Constructs initial Hamiltonian from kinetic energy plus summed atomic potentials.

Reference:

S. Lehtola, “Assessment of initial guesses for self-consistent field calculations.

Superposition of atomic potentials: Simple yet efficient”,

  1. Chem. Theory Comput. 15, 1593 (2019), doi:10.1021/acs.jctc.8b01089

Parameters:
  • mol (pyscf Mole) – pyscf Mole object.

  • *_ – Unused positional arguments (for interface compatibility).

Returns:

2D Hamiltonian matrix (T + V_SAP) in AO basis.

Return type:

numpy ndarray

qstack.spahm.guesses.check_nelec(nelec, nao)[source]

Validate that the number of electrons can be accommodated by available orbitals.

Parameters:
  • nelec (tuple or int) – Number of electrons (alpha, beta) or total.

  • nao (int) – Number of atomic orbitals.

Raises:

RuntimeError – If there are more electrons than available orbitals.

Warns:

RuntimeWarning – If all orbitals are filled (complete shell warning).

qstack.spahm.guesses.eigenvalue_grad(mol, e, c, s1, h1)[source]

Compute nuclear gradients of orbital eigenvalues from generalized eigenvalue problem HC = eSC.

Uses the Hellmann-Feynman theorem for eigenvalue derivatives.

Parameters:
  • mol (pyscf Mole) – pyscf Mole object.

  • e (numpy ndarray) – 1D array (nao,) of orbital eigenvalues.

  • c (numpy ndarray) – 2D array (nao, nao) of MO coefficients (eigenvectors).

  • s1 (numpy ndarray) – 3D array (3, nao, nao) - gradient of overlap matrix.

  • h1 (callable) – Function returning dH/dr[iat] - Hamiltonian gradient for atom iat.

Returns:

3D array (nao, natm, 3) of eigenvalue gradients in Eh/bohr.

Return type:

numpy ndarray

qstack.spahm.guesses.get_dm(v, nelec, spin)[source]

Construct density matrix from occupied molecular orbitals.

Parameters:
  • v (numpy ndarray) – 2D array of MO coefficients (eigenvectors), columns are MOs.

  • nelec (tuple) – Number of (alpha, beta) electrons.

  • spin (int or None) – Spin multiplicity. If None, assumes closed-shell (RHF).

Returns:

Density matrix in AO basis. - Closed-shell: 2D array (nao, nao) - Open-shell: 3D array (2, nao, nao) for alpha and beta

Return type:

numpy ndarray

qstack.spahm.guesses.get_guess(arg)[source]

Return guess Hamiltonian function by name.

Parameters:

arg (str) – Guess method name. Available options: - ‘core’: Core Hamiltonian (H_core). - ‘sad’: Superposition of Atomic Densities. - ‘sap’: Superposition of Atomic Potentials. - ‘gwh’: Generalized Wolfsberg-Helmholtz. - ‘lb’: Laikov-Briling 2020 (HF parameters). - ‘lb-hfs’: Laikov-Briling 2020 (HFS parameters). - ‘huckel’: Extended Hückel method.

Returns:

Guess Hamiltonian function with signature f(mol, xc) -> numpy.ndarray.

Return type:

callable

Raises:

RuntimeError – If the specified guess method is not available.

qstack.spahm.guesses.get_guess_g(arg)[source]

Return both guess Hamiltonian function and its gradient generator.

Parameters:

arg (str) – Guess method name. Available: ‘core’, ‘lb’.

Returns:

(hamiltonian_function, gradient_function) pair.

Return type:

tuple

Raises:

RuntimeError – If the specified guess method is not available for gradients.

qstack.spahm.guesses.get_occ(e, nelec, spin)[source]

Extract occupied orbital eigenvalues/energies.

Parameters:
  • e (numpy ndarray) – Full array of orbital eigenvalues (1D) or possibly arrays of larger dimensionality.

  • nelec (tuple) – Number of (alpha, beta) electrons.

  • spin (int or None) – Spin multiplicity. If None, assumes closed-shell.

Returns:

Occupied eigenvalues. Shape depends on spin: - Closed-shell (spin=None): 1D array of occupied eigenvalues - Open-shell: 2D array (2, nocc) for alpha and beta separately

Return type:

numpy ndarray

qstack.spahm.guesses.hcore(mol, *_)[source]

Compute guess Hamiltonian from core contributions (kinetic + nuclear + ECP).

Parameters:
  • mol (pyscf Mole) – pyscf Mole object.

  • *_ – Unused positional arguments (for interface compatibility).

Returns:

2D array containing the core Hamiltonian matrix in AO basis.

Return type:

numpy ndarray

qstack.spahm.guesses.hcore_grad(mf)[source]

Return core Hamiltonian gradient generator function.

Parameters:

mf – PySCF mean-field object.

Returns:

Function that returns core Hamiltonian gradient for a given atom.

Return type:

callable

qstack.spahm.guesses.solveF(mol, fock)[source]

Solves generalized eigenvalue problem FC = SCε for the Fock/Hamiltonian matrix.

Parameters:
  • mol (pyscf Mole) – pyscf Mole object.

  • fock (numpy ndarray) – 2D Fock or Hamiltonian matrix in AO basis.

Returns:

(eigenvalues, eigenvectors) where: - eigenvalues: 1D array of orbital energies - eigenvectors: 2D array of MO coefficients (columns are MOs)

Return type:

tuple