pyscf.doci package

Submodules

pyscf.doci.doci_mcscf module

class pyscf.doci.doci_mcscf.CASCI(mf_or_mol, ncas, nelecas, ncore=None)

Bases: pyscf.mcscf.casci.CASCI

DOCI-CASCI

class pyscf.doci.doci_mcscf.CASSCF(mf_or_mol, ncas, nelecas, ncore=None, frozen=None)

Bases: pyscf.mcscf.mc1step.CASSCF

DOCI-CASSCF

pyscf.doci.doci_slow module

Doubly occupied configuration interaction (DOCI)

class pyscf.doci.doci_slow.DOCI(*args, **kwargs)

Bases: pyscf.fci.direct_spin1.FCISolver

contract_2e(eri, fcivec, norb, nelec, link_index=None, **kwargs)

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()

gen_linkstr(norb, nelec)
get_init_guess(norb, nelec, nroots, hdiag)

Initial guess is the single Slater determinant

kernel(h1e, eri, norb, nelec, ci0=None, ecore=0, **kwargs)

Kernel function is the main driver of a method. Every method should define the kernel function as the entry of the calculation. Note the return value of kernel function is not strictly defined. It can be anything related to the method (such as the energy, the wave-function, the DFT mesh grids etc.).

make_hdiag(h1e, eri, norb, nelec, *args, **kwargs)

Diagonal Hamiltonian for Davidson preconditioner

make_rdm1(civec, norb, nelec, link_index=None)

Spin-traced one-particle density matrix

dm1[p,q] = <q_alpha^dagger p_alpha> + <q_beta^dagger p_beta>

The convention is based on McWeeney’s book, Eq (5.4.20) The contraction between 1-particle Hamiltonian and rdm1 is E = einsum(‘pq,qp’, h1, rdm1)

make_rdm12(civec, norb, nelec, link_index=None, reorder=True)

Spin traced 1- and 2-particle density matrices.

1pdm[p,q] = :math:`langle q_alpha^dagger p_alpha rangle +

langle q_beta^dagger p_beta rangle`;

2pdm[p,q,r,s] = :math:`langle p_alpha^dagger r_alpha^dagger s_alpha q_alpharangle +

langle p_beta^dagger r_alpha^dagger s_alpha q_betarangle + langle p_alpha^dagger r_beta^dagger s_beta q_alpharangle + langle p_beta^dagger r_beta^dagger s_beta q_betarangle`.

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 = None
pyscf.doci.doci_slow.contract_2e(eri, civec, norb, nelec, link_index=None)

Compute E_{pq}E_{rs}|CI>

pyscf.doci.doci_slow.kernel(h1e, eri, norb, nelec, ecore=0)
pyscf.doci.doci_slow.make_hdiag(h1e, eri, norb, nelec, opt=None)
pyscf.doci.doci_slow.make_rdm1(civec, norb, nelec, link_index=None)
pyscf.doci.doci_slow.make_rdm12(civec, norb, nelec, link_index=None, reorder=True)

Module contents

DOCI and DOCI based MCSCF