10.14. hci — Heat-bath selected CI

The HCI module is a native implementation of the Heat-bath selected CI algorithm, where selection is done based on strings rather than determinants. For example:

#!/usr/bin/env python
#
# Author: Qiming Sun <osirpt.sun@gmail.com>
#

'''
A simple example to run heat-bath selected CI
'''

from pyscf import gto, scf, ao2mo, hci

mol = gto.M(
    atom = 'H 0 0 0; F 0 0 1.1',
    basis = '6-31g',
)
myhf = scf.RHF(mol).run()

cisolver = hci.SCI(mol)
nmo = myhf.mo_coeff.shape[1]
nelec = mol.nelec
h1 = myhf.mo_coeff.T.dot(myhf.get_hcore()).dot(myhf.mo_coeff)
h2 = ao2mo.full(mol, myhf.mo_coeff)
e, civec = cisolver.kernel(h1, h2, nmo, nelec, verbose=4)

See also cornell_shci.

10.14.1. Examples

Relevant examples examples/hci/00-simple_hci.py examples/hci/10-spin.py

10.14.2. Program reference

10.14.2.1. Interface

Selected CI using Heat-Bath CI algorithm (JCTC 2016, 12, 3674-3680)

Simple usage:

.. py:attribute:: SCI
module

pyscf.hci.hci

alias of pyscf.hci.hci.SelectedCI

class pyscf.hci.hci.SelectedCI(mol=None)[source]
absorb_h1e(eri, *args, **kwargs)[source]

Modify 2e Hamiltonian to include 1e Hamiltonian contribution.

contract_2e(h1_h2, civec, norb, nelec, hdiag=None, **kwargs)[source]

Contract the 4-index tensor eri[pqrs] with a FCI vector

\[ \begin{align}\begin{aligned}\begin{split}|output\rangle = E_{pq} E_{rs} eri_{pq,rs} |CI\rangle \\\end{split}\\\begin{split}E_{pq}E_{rs} = E_{pr,qs} + \delta_{qr} E_{ps} \\\end{split}\\E_{pq} = p^+ q + \bar{p}^+ \bar{q}\\E_{pr,qs} = p^+ r^+ s q + \bar{p}^+ r^+ s \bar{q} + ...\end{aligned}\end{align} \]

\(p,q,...\) means spin-up orbitals and \(\bar{p}, \bar{q}\) means spin-down orbitals.

Note the input argument eri is NOT the 2e hamiltonian tensor. 2e hamiltonian is

\[\begin{split}h2e &= (pq|rs) E_{pr,qs} \\ &= (pq|rs) (E_{pq}E_{rs} - \delta_{qr} E_{ps}) \\ &= eri_{pq,rs} E_{pq}E_{rs} \\\end{split}\]

So the relation between eri and hamiltonian (the 2e-integral tensor) is

\[eri_{pq,rs} = (pq|rs) - (1/Nelec) \sum_q (pq|qs)\]

to restore the symmetry between pq and rs,

\[eri_{pq,rs} = (pq|rs) - (.5/Nelec) [\sum_q (pq|qs) + \sum_p (pq|rp)]\]

See also direct_spin1.absorb_h1e()

make_hdiag(h1e, eri, strs, norb, nelec)[source]

Diagonal Hamiltonian for Davidson preconditioner

make_rdm12s(civec, norb, nelec)[source]

Spin separated 1- and 2-particle density matrices. The return values include two lists, a list of 1-particle density matrices and a list of 2-particle density matrices. The density matrices are: (alpha,alpha), (beta,beta) for 1-particle density matrices; (alpha,alpha,alpha,alpha), (alpha,alpha,beta,beta), (beta,beta,beta,beta) for 2-particle density matrices.

1pdm[p,q] = \(\langle q^\dagger p\rangle\); 2pdm[p,q,r,s] = \(\langle p^\dagger r^\dagger s q\rangle\).

Energy should be computed as E = einsum(‘pq,qp’, h1, 1pdm) + 1/2 * einsum(‘pqrs,pqrs’, eri, 2pdm) where h1[p,q] = <p|h|q> and eri[p,q,r,s] = (pq|rs)

spin_square(civec, norb, nelec)[source]

Spin square for RHF-FCI CI wfn only (obtained from spin-degenerated Hamiltonian)

pyscf.hci.hci.fix_spin(myci, shift=0.2, ss=None, **kwargs)[source]

If Selected CI solver cannot stick on spin eigenfunction, modify the solver by adding a shift on spin square operator

\[(H + shift*S^2) |\Psi\rangle = E |\Psi\rangle\]
Args:

myci : An instance of SelectedCI

Kwargs:
shiftfloat

Level shift for states which have different spin

ssnumber

S^2 expection value == s*(s+1)

Returns

A modified Selected CI object based on myci.

pyscf.hci.hci.make_rdm12s(civec, norb, nelec)[source]

Spin orbital 1- and 2-particle reduced density matrices (aa, bb, aaaa, aabb, bbbb)