Skip to content

Latest commit

 

History

History
163 lines (126 loc) · 6.63 KB

scf_developer.rst

File metadata and controls

163 lines (126 loc) · 6.63 KB

Self-consistent field (SCF) methods

Modules: :mod:`scf`, :mod:`pbc.scf`, :mod:`soscf`

Overview

The class :class:`pyscf.scf.hf.SCF` is the base class for all mean-field methods, which are implemented as derived classes. These include

:class:`pyscf.scf.hf.RHF` RHF
:class:`pyscf.scf.uhf.UHF` UHF
:class:`pyscf.scf.rohf.ROHF` ROHF
:class:`pyscf.scf.ghf.GHF` GHF
:class:`pyscf.scf.dhf.DHF` Dirac-HF
:class:`pyscf.scf.dhf.RDHF` restricted Dirac-HF
:class:`pyscf.dft.rks.RKS` RKS-DFT
:class:`pyscf.dft.uks.UKS` UKS-DFT
:class:`pyscf.dft.roks.ROKS` ROKS-DFT
:class:`pyscf.dft.gks.GKS` GKS-DFT
:class:`pyscf.dft.dks.DKS` Dirac-KS-DFT
:class:`pyscf.dft.dks.RDKS` restricted Dirac-KS-DFT
:class:`pyscf.pbc.scf.hf.SCF` base class for SCF methods with PBC at \Gamma point
:class:`pyscf.pbc.scf.hf.RHF` \Gamma-point RHF
:class:`pyscf.pbc.scf.uhf.UHF` \Gamma-point UHF
:class:`pyscf.pbc.scf.rohf.ROHF` \Gamma-point ROHF
:class:`pyscf.pbc.scf.ghf.GHF` \Gamma-point GHF
:class:`pyscf.pbc.dft.rks.RKS` \Gamma-point RKS-DFT
:class:`pyscf.pbc.dft.uks.UKS` \Gamma-point UKS-DFT
:class:`pyscf.pbc.dft.roks.ROKS` \Gamma-point ROKS-DFT
:class:`pyscf.pbc.scf.khf.KSCF` base class for SCF methods with PBC and with k-point sampling
:class:`pyscf.pbc.scf.khf.KRHF` KRHF
:class:`pyscf.pbc.scf.kuhf.KUHF` KUHF
:class:`pyscf.pbc.scf.krohf.KROHF` KROHF
:class:`pyscf.pbc.scf.kghf.KGHF` KGHF
:class:`pyscf.pbc.dft.krks.KRKS` KRKS-DFT
:class:`pyscf.pbc.dft.kuks.KUKS` KUKS-DFT
:class:`pyscf.pbc.dft.kroks.KROKS` KROKS-DFT

The key attributes of the SCF class include

:attr:`init_guess` initial guess method
:attr:`DIIS` DIIS method, can be :class:`pyscf.scf.diis.DIIS`, :class:`pyscf.scf.diis.EDIIS`, or :class:`pyscf.scf.diis.ADIIS`
:attr:`mo_coeff` saved MO coefficients
:attr:`mo_energy` saved MO energies
:attr:`mo_occ` saved MO occupations

The SCF iterative loop is the :func:`pyscf.scf.hf.kernel` function which takes a mean-field object. The :func:`kernel` function carries out the following steps:

Internally, different methods reuse this kernel by overwriting the mean-field methods.

Integrals

The molecular HF operator relies on the following underlying AO integral functions:

Custom Hamiltonians

The HF approximation SCF procedure can be used for arbitrary Hamiltonians. This is useful when considering model Hamiltonians. The key functions to supply are the initial guess (which can also be generated by supplying an initial DM), :meth:`get_ovlp`, :meth:`get_hcore`, and the two electron integrals (attribute :attr:`._eri` of the mean-field object). The following shows an example of HF with a Hubbard model Hamiltonian:

import numpy
from pyscf import gto, scf, ao2mo, mcscf

mol = gto.M()

# incore_anyway=True ensures the customized Hamiltonian (the _eri attribute)
# is used.  Without this parameter, the MO integral transformation used in
# subsequent post-HF calculations may
# ignore the customized Hamiltonian if there is not enough memory.
mol.incore_anyway = True

n = 10
mol.nelectron = n

h1 = numpy.zeros((n,n))
for i in range(n-1):
    h1[i,i+1] = h1[i+1,i] = -1.0
h1[n-1,0] = h1[0,n-1] = -1.0  # PBC
eri = numpy.zeros((n,n,n,n))
for i in range(n):
    eri[i,i,i,i] = 4.0

mf = scf.RHF(mol)
mf.get_hcore = lambda *args: h1
mf.get_ovlp = lambda *args: numpy.eye(n)
# ao2mo.restore(8, eri, n) to get 8-fold permutation symmetry of the integrals
# ._eri only supports two-electron integrals with 4-fold or 8-fold symmetry.
mf._eri = ao2mo.restore(8, eri, n)
mf.init_guess = '1e'
mf.kernel()

# post-SCF calculation
mycas = mcscf.CASSCF(mf, 4, 4)
mycas.kernel()