pyscf.fci package

Submodules

pyscf.fci.addons module

class pyscf.fci.addons.SpinPenaltyFCISolver(fcibase, shift, ss_value)[source]

Bases: object

base_contract_2e(*args, **kwargs)[source]
contract_2e(eri, fcivec, norb, nelec, link_index=None, **kwargs)[source]
undo_fix_spin()[source]
pyscf.fci.addons.cre_a(ci0, norb, neleca_nelecb, ap_id)[source]

Construct (N+1)-electron wavefunction by adding an alpha electron in the N-electron wavefunction.

… math:

|N+1\rangle = \hat{a}^+_p |N\rangle
Args:
ci02D array

CI coefficients, row for alpha strings and column for beta strings.

norbint

Number of orbitals.

(neleca,nelecb)(int,int)

Number of (alpha, beta) electrons of the input CI function

ap_idint

Orbital index (0-based), for the creation operator

Returns:

2D array, row for alpha strings and column for beta strings. Note it has different number of rows to the input CI coefficients.

pyscf.fci.addons.cre_b(ci0, norb, neleca_nelecb, ap_id)[source]

Construct (N+1)-electron wavefunction by adding a beta electron in the N-electron wavefunction.

Args:
ci02D array

CI coefficients, row for alpha strings and column for beta strings.

norbint

Number of orbitals.

(neleca,nelecb)(int,int)

Number of (alpha, beta) electrons of the input CI function

ap_idint

Orbital index (0-based), for the creation operator

Returns:

2D array, row for alpha strings and column for beta strings. Note it has different number of columns to the input CI coefficients.

pyscf.fci.addons.cylindrical_init_guess(mol, norb, nelec, orbsym, wfnsym=0, singlet=True, nroots=1)[source]

FCI initial guess for system of cylindrical symmetry. (In testing)

Examples:

>>> mol = gto.M(atom='O; O 1 1.2', spin=2, symmetry=True)
>>> orbsym = [6,7,2,3]
>>> ci0 = fci.addons.cylindrical_init_guess(mol, 4, (3,3), orbsym, wfnsym=10)[0]
>>> print(ci0.reshape(4,4))
>>> ci0 = fci.addons.cylindrical_init_guess(mol, 4, (3,3), orbsym, wfnsym=10, singlet=False)[0]
>>> print(ci0.reshape(4,4))
pyscf.fci.addons.des_a(ci0, norb, neleca_nelecb, ap_id)[source]

Construct (N-1)-electron wavefunction by removing an alpha electron from the N-electron wavefunction.

… math:

|N-1\rangle = \hat{a}_p |N\rangle
Args:
ci02D array

CI coefficients, row for alpha strings and column for beta strings.

norbint

Number of orbitals.

(neleca,nelecb)(int,int)

Number of (alpha, beta) electrons of the input CI function

ap_idint

Orbital index (0-based), for the annihilation operator

Returns:

2D array, row for alpha strings and column for beta strings. Note it has different number of rows to the input CI coefficients

pyscf.fci.addons.des_b(ci0, norb, neleca_nelecb, ap_id)[source]

Construct (N-1)-electron wavefunction by removing a beta electron from N-electron wavefunction.

Args:
ci02D array

CI coefficients, row for alpha strings and column for beta strings.

norbint

Number of orbitals.

(neleca,nelecb)(int,int)

Number of (alpha, beta) electrons of the input CI function

ap_idint

Orbital index (0-based), for the annihilation operator

Returns:

2D array, row for alpha strings and column for beta strings. Note it has different number of columns to the input CI coefficients.

pyscf.fci.addons.det_overlap(string1, string2, norb, s=None)[source]

Determinants overlap on non-orthogonal one-particle basis

pyscf.fci.addons.fix_spin(fciobj, shift=0.2, ss=None, **kwargs)[source]

If FCI solver cannot stay on spin eigenfunction, this function can add a shift to the states which have wrong spin.

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

fciobj : An instance of FCISolver

Kwargs:
shiftfloat

Level shift for states which have different spin

ssnumber

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

Returns

A modified FCI object based on fciobj.

pyscf.fci.addons.fix_spin_(fciobj, shift=0.1, ss=None)[source]
pyscf.fci.addons.guess_wfnsym(ci, norb, nelec, orbsym)[source]

Guess the wavefunction symmetry based on the non-zero elements in the given CI coefficients.

Args:
ci2D array

CI coefficients, row for alpha strings and column for beta strings.

norbint

Number of orbitals.

nelecint or 2-item list

Number of electrons, or 2-item list for (alpha, beta) electrons

orbsymlist of int

The irrep ID for each orbital.

Returns:

Irrep ID

pyscf.fci.addons.initguess_triplet(norb, nelec, binstring)[source]

Generate a triplet initial guess for FCI solver

pyscf.fci.addons.large_ci(ci, norb, nelec, tol=0.1, return_strs=True)[source]

Search for the largest CI coefficients

pyscf.fci.addons.overlap(bra, ket, norb, nelec, s=None)[source]

Overlap between two CI wavefunctions

Args:
s2D array or a list of 2D array

The overlap matrix of non-orthogonal one-particle basis

pyscf.fci.addons.symm_initguess(norb, nelec, orbsym, wfnsym=0, irrep_nelec=None)[source]

Generate CI wavefunction initial guess which has the given symmetry.

Args:
norbint

Number of orbitals.

nelecint or 2-item list

Number of electrons, or 2-item list for (alpha, beta) electrons

orbsymlist of int

The irrep ID for each orbital.

Kwags:
wfnsymint

The irrep ID of target symmetry

irrep_nelecdict

Freeze occupancy for certain irreps

Returns:

CI coefficients 2D array which has the target symmetry.

pyscf.fci.addons.symmetrize_wfn(ci, norb, nelec, orbsym, wfnsym=0)[source]

Symmetrize the CI wavefunction by zeroing out the determinants which do not have the right symmetry.

Args:
ci2D array

CI coefficients, row for alpha strings and column for beta strings.

norbint

Number of orbitals.

nelecint or 2-item list

Number of electrons, or 2-item list for (alpha, beta) electrons

orbsymlist of int

The irrep ID for each orbital.

Kwags:
wfnsymint

The irrep ID of target symmetry

Returns:

2D array which is the symmetrized CI coefficients

pyscf.fci.addons.transform_ci(ci, nelec, u)[source]

Transform CI coefficients to the representation in new one-particle basis. Solving CI problem for Hamiltonian h1, h2 defined in old basis, CI_old = fci.kernel(h1, h2, …) Given orbital rotation u, the CI problem can be either solved by transforming the Hamiltonian, or transforming the coefficients. CI_new = fci.kernel(u^T*h1*u, …) = transform_ci_for_orbital_rotation(CI_old, u)

Args:
u2D array or a list of 2D array

the orbital rotation to transform the old one-particle basis to new one-particle basis. If u is not a squared matrix, the resultant CI coefficients array may have different shape to the input CI coefficients.

pyscf.fci.addons.transform_ci_for_orbital_rotation(ci, norb, nelec, u)[source]

Transform CI coefficients (dimension conserved) to the representation in new one-particle basis. Solving CI problem for Hamiltonian h1, h2 defined in old basis, CI_old = fci.kernel(h1, h2, …) Given orbital rotation u, the CI problem can be either solved by transforming the Hamiltonian, or transforming the coefficients. CI_new = fci.kernel(u^T*h1*u, …) = transform_ci_for_orbital_rotation(CI_old, u)

Args:
ua squared 2D array or a list of 2D array

the orbital rotation to transform the old one-particle basis to new one-particle basis

pyscf.fci.cistring module

class pyscf.fci.cistring.OIndexList[source]

Bases: ndarray

pyscf.fci.cistring.addr2str(norb, nelec, addr)[source]

Convert CI determinant address to string

pyscf.fci.cistring.addrs2str(norb, nelec, addrs)[source]

Convert a list of CI determinant address to string

pyscf.fci.cistring.cre_des_sign(p, q, string0)[source]
pyscf.fci.cistring.cre_sign(p, string0)[source]
pyscf.fci.cistring.des_sign(p, string0)[source]
pyscf.fci.cistring.gen_cre_str_index(orb_list, nelec)[source]

linkstr_index to map between N electron string to N+1 electron string. It maps the given string to the address of the string which is generated by the creation operator.

For given string str0, index[str0] is nvir x 4 array. Each entry [i(cre),–,str1,sign] means starting from str0, creating i, to get str1.

Returns [[cre, -, target_address, parity], …]

pyscf.fci.cistring.gen_des_str_index(orb_list, nelec)[source]

linkstr_index to map between N electron string to N-1 electron string. It maps the given string to the address of the string which is generated by the annihilation operator.

For given string str0, index[str0] is nvir x 4 array. Each entry [–,i(des),str1,sign] means starting from str0, annihilating i, to get str1.

Returns [[-, des, target_address, parity], …]

pyscf.fci.cistring.gen_linkstr_index(orb_list, nocc, strs=None, tril=False)[source]

Look up table, for the strings relationship in terms of a creation-annihilating operator pair.

For given string str0, index[str0] is (nocc+nocc*nvir) x 4 array. The first nocc rows [i(:occ),i(:occ),str0,sign] are occupied-occupied excitations, which do not change the string. The next nocc*nvir rows [a(:vir),i(:occ),str1,sign] are occupied-virtual exciations, starting from str0, annihilating i, creating a, to get str1.

pyscf.fci.cistring.gen_linkstr_index_o1(orb_list, nelec, strs=None, tril=False)[source]

Look up table, for the strings relationship in terms of a creation-annihilating operator pair.

Similar to gen_linkstr_index. The only difference is the input argument strs, which needs to be a list of OIndexList in this function.

pyscf.fci.cistring.gen_linkstr_index_trilidx(orb_list, nocc, strs=None)[source]

Generate linkstr_index with the assumption that \(p^+ q|0\rangle\) where \(p > q\). So the resultant link_index has the structure [pq, *, str1, sign]. It is identical to a call to reform_linkstr_index(gen_linkstr_index(...)).

pyscf.fci.cistring.gen_occslst(orb_list, nelec)[source]

Generate occupied orbital list for each string.

Returns:

List of lists of int32. Each inner list has length equal to the number of electrons, and contains the occupied orbitals in the corresponding string.

Example:

>>> [bin(x) for x in make_strings((0, 1, 2, 3), 2)]
['0b11', '0b101', '0b110', '0b1001', '0b1010', '0b1100']
>>> gen_occslst((0, 1, 2, 3), 2)
OIndexList([[0, 1],
            [0, 2],
            [1, 2],
            [0, 3],
            [1, 3],
            [2, 3]], dtype=int32)
pyscf.fci.cistring.gen_strings4orblist(orb_list, nelec)

Generate string from the given orbital list.

Returns:

list of int64. One element represents one string in binary format. The binary format takes the convention that the one bit stands for one orbital, bit-1 means occupied and bit-0 means unoccupied. The lowest (right-most) bit corresponds to the lowest orbital in the orb_list.

Examples:

>>> [bin(x) for x in make_strings((0,1,2,3),2)]
[0b11, 0b101, 0b110, 0b1001, 0b1010, 0b1100]
>>> [bin(x) for x in make_strings((3,1,0,2),2)]
[0b1010, 0b1001, 0b11, 0b1100, 0b110, 0b101]
pyscf.fci.cistring.make_strings(orb_list, nelec)[source]

Generate string from the given orbital list.

Returns:

list of int64. One element represents one string in binary format. The binary format takes the convention that the one bit stands for one orbital, bit-1 means occupied and bit-0 means unoccupied. The lowest (right-most) bit corresponds to the lowest orbital in the orb_list.

Examples:

>>> [bin(x) for x in make_strings((0,1,2,3),2)]
[0b11, 0b101, 0b110, 0b1001, 0b1010, 0b1100]
>>> [bin(x) for x in make_strings((3,1,0,2),2)]
[0b1010, 0b1001, 0b11, 0b1100, 0b110, 0b101]
pyscf.fci.cistring.num_strings(n, m)[source]
pyscf.fci.cistring.parity(string0, string1)[source]
pyscf.fci.cistring.reform_linkstr_index(link_index)[source]

Compress the (a, i) pair index in linkstr_index to a lower triangular index. The compressed indices can match the 4-fold symmetry of integrals.

pyscf.fci.cistring.str2addr(norb, nelec, string)[source]

Convert string to CI determinant address

pyscf.fci.cistring.strs2addr(norb, nelec, strings)[source]

Convert a list of string to CI determinant address

pyscf.fci.cistring.sub_addrs(norb, nelec, orbital_indices, sub_nelec=0)[source]

The addresses of the determinants which include the specified orbital indices. The size of the returned addresses is equal to the number of determinants of (norb, nelec) system.

pyscf.fci.cistring.tn_strs(norb, nelec, n)[source]

Generate strings for Tn amplitudes. Eg n=1 (T1) has nvir*nocc strings, n=2 (T2) has nvir*(nvir-1)/2 * nocc*(nocc-1)/2 strings.

pyscf.fci.direct_ep module

Electron phonon coupling

Ref:

arXiv:1602.04195

pyscf.fci.direct_ep.contract_1e(h1e, fcivec, nsite, nelec, nphonon)[source]
pyscf.fci.direct_ep.contract_2e(eri, fcivec, nsite, nelec, nphonon)[source]
pyscf.fci.direct_ep.contract_2e_hubbard(u, fcivec, nsite, nelec, nphonon)[source]
pyscf.fci.direct_ep.contract_all(t, u, g, hpp, ci0, nsite, nelec, nphonon)[source]
pyscf.fci.direct_ep.contract_ep(g, fcivec, nsite, nelec, nphonon)[source]
pyscf.fci.direct_ep.contract_pp(hpp, fcivec, nsite, nelec, nphonon)[source]
pyscf.fci.direct_ep.cre_phonon(fcivec, nsite, nelec, nphonon, site_id)[source]
pyscf.fci.direct_ep.des_phonon(fcivec, nsite, nelec, nphonon, site_id)[source]
pyscf.fci.direct_ep.kernel(t, u, g, hpp, nsite, nelec, nphonon, tol=1e-09, max_cycle=100, verbose=0, ecore=0, **kwargs)[source]
pyscf.fci.direct_ep.make_hdiag(t, u, g, hpp, nsite, nelec, nphonon)[source]
pyscf.fci.direct_ep.make_rdm12e(fcivec, nsite, nelec)[source]

1-electron and 2-electron density matrices dm_pq = <|p^+ q|> dm_{pqrs} = <|p^+ r^+ q s|> (note 2pdm is ordered in chemist notation)

pyscf.fci.direct_ep.make_rdm1e(fcivec, nsite, nelec)[source]

1-electron density matrix dm_pq = <|p^+ q|>

pyscf.fci.direct_ep.make_rdm1p(fcivec, nsite, nelec, nphonon)[source]

1-phonon density matrix dm_pq = <|p^+ q|>

pyscf.fci.direct_ep.make_shape(nsite, nelec, nphonon)[source]
pyscf.fci.direct_ep.slices_for(psite_id, nsite, nphonon)[source]
pyscf.fci.direct_ep.slices_for_cre(psite_id, nsite, nphonon)[source]
pyscf.fci.direct_ep.slices_for_des(psite_id, nsite, nphonon)[source]

pyscf.fci.direct_nosym module

Different FCI solvers are implemented to support different type of symmetry.

Symmetry

File Point group Spin singlet Real hermitian* Alpha/beta degeneracy direct_spin0_symm Yes Yes Yes Yes direct_spin1_symm Yes No Yes Yes direct_spin0 No Yes Yes Yes direct_spin1 No No Yes Yes direct_uhf No No Yes No direct_nosym No No No** Yes

  • Real hermitian Hamiltonian implies (ij|kl) = (ji|kl) = (ij|lk) = (ji|lk)

** Hamiltonian is real but not hermitian, (ij|kl) != (ji|kl) …

pyscf.fci.direct_nosym.FCI

alias of FCISolver

class pyscf.fci.direct_nosym.FCISolver(*args, **kwargs)[source]

Bases: FCISolver

absorb_h1e(h1e, eri, norb, nelec, fac=1)[source]

Modify 2e Hamiltonian to include 1e Hamiltonian contribution.

contract_1e(h1e, fcivec, norb, nelec, link_index=None)[source]

Contract the 1-electron Hamiltonian with a FCI vector to get a new FCI vector.

contract_2e(eri, fcivec, norb, nelec, link_index=None)[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()

eig(op, x0=None, precond=None, **kwargs)[source]
kernel(h1e, eri, norb, nelec, ci0=None, tol=None, lindep=None, max_cycle=None, max_space=None, nroots=None, davidson_only=None, pspace_size=None, orbsym=None, wfnsym=None, ecore=0, **kwargs)[source]
Args:
h1e: ndarray

1-electron Hamiltonian

eri: ndarray

2-electron integrals in chemist’s notation

norb: int

Number of orbitals

nelec: int or (int, int)

Number of electrons of the system

Kwargs:
ci0: ndarray

Initial guess

link_index: ndarray

A lookup table to cache the addresses of CI determinants in wave-function vector

tol: float

Convergence tolerance

lindep: float

Linear dependence threshold

max_cycle: int

Max. iterations for diagonalization

max_space: int

Max. trial vectors to store for sub-space diagonalization method

nroots: int

Number of states to solve

davidson_only: bool

Whether to call subspace diagonlization (davidson solver) or do a full diagonlization (lapack eigh) for small systems

pspace_size: int

Number of determinants as the threshold of “small systems”,

hop: function(c) => array_like_c

Function to use for the Hamiltonian multiplication with trial vector

Note: davidson solver requires more arguments. For the parameters not dispatched, they can be passed to davidson solver via the extra keyword arguments **kwargs

pyscf.fci.direct_nosym.absorb_h1e(h1e, eri, norb, nelec, fac=1)[source]

Modify 2e Hamiltonian to include 1e Hamiltonian contribution.

pyscf.fci.direct_nosym.contract_1e(h1e, fcivec, norb, nelec, link_index=None)[source]
pyscf.fci.direct_nosym.contract_2e(eri, fcivec, norb, nelec, link_index=None)[source]

Contract the 2-electron Hamiltonian with a FCI vector to get a new FCI vector.

Note the input arg eri is NOT the 2e hamiltonian matrix, the 2e hamiltonian is

\[\begin{split}h2e &= eri_{pq,rs} p^+ q r^+ s \\ &= (pq|rs) p^+ r^+ s q - (pq|rs) \delta_{qr} p^+ s\end{split}\]

So eri is defined as

\[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_nosym.absorb_h1e()

pyscf.fci.direct_nosym.energy(h1e, eri, fcivec, norb, nelec, link_index=None)[source]

Compute the FCI electronic energy for given Hamiltonian and FCI vector.

pyscf.fci.direct_spin0 module

FCI solver for Singlet state

Different FCI solvers are implemented to support different type of symmetry.

Symmetry

File Point group Spin singlet Real hermitian* Alpha/beta degeneracy direct_spin0_symm Yes Yes Yes Yes direct_spin1_symm Yes No Yes Yes direct_spin0 No Yes Yes Yes direct_spin1 No No Yes Yes direct_uhf No No Yes No direct_nosym No No No** Yes

  • Real hermitian Hamiltonian implies (ij|kl) = (ji|kl) = (ij|lk) = (ji|lk)

** Hamiltonian is real but not hermitian, (ij|kl) != (ji|kl) …

direct_spin0 solver is specified for singlet state. However, calling this solver sometimes ends up with the error “State not singlet x.xxxxxxe-06” due to numerical issues. Calling direct_spin1 for singlet state is slightly slower but more robust than direct_spin0 especially when combining to energy penalty method (fix_spin_())

pyscf.fci.direct_spin0.FCI

alias of FCISolver

class pyscf.fci.direct_spin0.FCISolver(mol=None)[source]

Bases: FCISolver

contract_1e(f1e, fcivec, norb, nelec, link_index=None, **kwargs)[source]

Contract the 1-electron Hamiltonian with a FCI vector to get a new FCI vector.

contract_2e(eri, fcivec, norb, nelec, link_index=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()

energy(h1e, eri, fcivec, norb, nelec, link_index=None)[source]

Compute the FCI electronic energy for given Hamiltonian and FCI vector.

gen_linkstr(norb, nelec, tril=True, spin=None)[source]
get_init_guess(norb, nelec, nroots, hdiag)[source]

Initial guess is the single Slater determinant

kernel(h1e, eri, norb, nelec, ci0=None, tol=None, lindep=None, max_cycle=None, max_space=None, nroots=None, davidson_only=None, pspace_size=None, orbsym=None, wfnsym=None, ecore=0, **kwargs)[source]
Args:
h1e: ndarray

1-electron Hamiltonian

eri: ndarray

2-electron integrals in chemist’s notation

norb: int

Number of orbitals

nelec: int or (int, int)

Number of electrons of the system

Kwargs:
ci0: ndarray

Initial guess

link_index: ndarray

A lookup table to cache the addresses of CI determinants in wave-function vector

tol: float

Convergence tolerance

lindep: float

Linear dependence threshold

max_cycle: int

Max. iterations for diagonalization

max_space: int

Max. trial vectors to store for sub-space diagonalization method

nroots: int

Number of states to solve

davidson_only: bool

Whether to call subspace diagonlization (davidson solver) or do a full diagonlization (lapack eigh) for small systems

pspace_size: int

Number of determinants as the threshold of “small systems”,

hop: function(c) => array_like_c

Function to use for the Hamiltonian multiplication with trial vector

Note: davidson solver requires more arguments. For the parameters not dispatched, they can be passed to davidson solver via the extra keyword arguments **kwargs

static make_hdiag(h1e, eri, norb, nelec, compress=False)

Diagonal Hamiltonian for Davidson preconditioner

Kwargs:

compress (bool) : whether to remove symmetry forbidden elements

make_rdm1(fcivec, norb, nelec, link_index=None)[source]

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(fcivec, norb, nelec, link_index=None, reorder=True)[source]

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)

make_rdm1s(fcivec, norb, nelec, link_index=None)[source]

Spin separated 1-particle density matrices. The return values include two density matrices: (alpha,alpha), (beta,beta)

dm1[p,q] = <q^dagger p>

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)

trans_rdm1(cibra, ciket, norb, nelec, link_index=None)[source]

Spin traced transition 1-particle transition density matrices.

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

trans_rdm12(cibra, ciket, norb, nelec, link_index=None, reorder=True)[source]

Spin traced transition 1- and 2-particle transition density matrices.

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

trans_rdm1s(cibra, ciket, norb, nelec, link_index=None)[source]

Spin separated transition 1-particle density matrices. The return values include two density matrices: (alpha,alpha), (beta,beta). See also function make_rdm1s()

1pdm[p,q] = \(\langle q^\dagger p \rangle\)

pyscf.fci.direct_spin0.contract_1e(f1e, fcivec, norb, nelec, link_index=None)[source]

Contract the 1-electron Hamiltonian with a FCI vector to get a new FCI vector.

pyscf.fci.direct_spin0.contract_2e(eri, fcivec, norb, nelec, link_index=None)[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()

pyscf.fci.direct_spin0.energy(h1e, eri, fcivec, norb, nelec, link_index=None)[source]
pyscf.fci.direct_spin0.get_init_guess(norb, nelec, nroots, hdiag)[source]
pyscf.fci.direct_spin0.kernel(h1e, eri, norb, nelec, ci0=None, level_shift=0.001, tol=1e-10, lindep=1e-14, max_cycle=50, max_space=12, nroots=1, davidson_only=False, pspace_size=400, orbsym=None, wfnsym=None, ecore=0, **kwargs)[source]
pyscf.fci.direct_spin0.kernel_ms0(fci, h1e, eri, norb, nelec, ci0=None, link_index=None, tol=None, lindep=None, max_cycle=None, max_space=None, nroots=None, davidson_only=None, pspace_size=None, hop=None, max_memory=None, verbose=None, ecore=0, **kwargs)[source]
pyscf.fci.direct_spin0.make_hdiag(h1e, eri, norb, nelec, compress=False)[source]

Diagonal Hamiltonian for Davidson preconditioner

Kwargs:

compress (bool) : whether to remove symmetry forbidden elements

pyscf.fci.direct_spin0.make_rdm1(fcivec, norb, nelec, link_index=None)[source]

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)

pyscf.fci.direct_spin0.make_rdm12(fcivec, norb, nelec, link_index=None, reorder=True)[source]

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)

pyscf.fci.direct_spin0.make_rdm1s(fcivec, norb, nelec, link_index=None)[source]

Spin separated 1-particle density matrices. The return values include two density matrices: (alpha,alpha), (beta,beta)

dm1[p,q] = <q^dagger p>

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)

pyscf.fci.direct_spin0.trans_rdm1(cibra, ciket, norb, nelec, link_index=None)[source]

Spin traced transition 1-particle transition density matrices.

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

pyscf.fci.direct_spin0.trans_rdm12(cibra, ciket, norb, nelec, link_index=None, reorder=True)[source]

Spin traced transition 1- and 2-particle transition density matrices.

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

pyscf.fci.direct_spin0.trans_rdm1s(cibra, ciket, norb, nelec, link_index=None)[source]

Spin separated transition 1-particle density matrices. The return values include two density matrices: (alpha,alpha), (beta,beta). See also function make_rdm1s()

1pdm[p,q] = \(\langle q^\dagger p \rangle\)

pyscf.fci.direct_spin0_symm module

Different FCI solvers are implemented to support different type of symmetry.

Symmetry

File Point group Spin singlet Real hermitian* Alpha/beta degeneracy direct_spin0_symm Yes Yes Yes Yes direct_spin1_symm Yes No Yes Yes direct_spin0 No Yes Yes Yes direct_spin1 No No Yes Yes direct_uhf No No Yes No direct_nosym No No No** Yes

  • Real hermitian Hamiltonian implies (ij|kl) = (ji|kl) = (ij|lk) = (ji|lk)

** Hamiltonian is real but not hermitian, (ij|kl) != (ji|kl) …

pyscf.fci.direct_spin0_symm.FCI

alias of FCISolver

class pyscf.fci.direct_spin0_symm.FCISolver(mol=None, **kwargs)[source]

Bases: FCISolver

absorb_h1e(h1e, eri, norb, nelec, fac=1)

Modify 2e Hamiltonian to include 1e Hamiltonian contribution.

contract_1e(f1e, fcivec, norb, nelec, link_index=None, **kwargs)[source]

Contract the 1-electron Hamiltonian with a FCI vector to get a new FCI vector.

contract_2e(eri, fcivec, norb, nelec, link_index=None, orbsym=None, wfnsym=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()

contract_ss(fcivec, norb, nelec)[source]
davidson_only = True
dump_flags(verbose=None)[source]
get_init_guess(norb, nelec, nroots, hdiag, orbsym=None, wfnsym=None)[source]

Initial guess is the single Slater determinant

guess_wfnsym(norb, nelec, fcivec=None, orbsym=None, wfnsym=None, **kwargs)

Guess point group symmetry of the FCI wavefunction. If fcivec is given, the symmetry of fcivec is used. Otherwise the symmetry is same to the HF determinant.

kernel(h1e, eri, norb, nelec, ci0=None, tol=None, lindep=None, max_cycle=None, max_space=None, nroots=None, davidson_only=None, pspace_size=None, orbsym=None, wfnsym=None, ecore=0, **kwargs)[source]
Args:
h1e: ndarray

1-electron Hamiltonian

eri: ndarray

2-electron integrals in chemist’s notation

norb: int

Number of orbitals

nelec: int or (int, int)

Number of electrons of the system

Kwargs:
ci0: ndarray

Initial guess

link_index: ndarray

A lookup table to cache the addresses of CI determinants in wave-function vector

tol: float

Convergence tolerance

lindep: float

Linear dependence threshold

max_cycle: int

Max. iterations for diagonalization

max_space: int

Max. trial vectors to store for sub-space diagonalization method

nroots: int

Number of states to solve

davidson_only: bool

Whether to call subspace diagonlization (davidson solver) or do a full diagonlization (lapack eigh) for small systems

pspace_size: int

Number of determinants as the threshold of “small systems”,

hop: function(c) => array_like_c

Function to use for the Hamiltonian multiplication with trial vector

Note: davidson solver requires more arguments. For the parameters not dispatched, they can be passed to davidson solver via the extra keyword arguments **kwargs

make_hdiag(h1e, eri, norb, nelec, compress=False)[source]

Diagonal Hamiltonian for Davidson preconditioner

Kwargs:

compress (bool) : whether to remove symmetry forbidden elements

pspace(h1e, eri, norb, nelec, hdiag, np=400)[source]

pspace Hamiltonian to improve Davidson preconditioner. See, CPL, 169, 463

pspace_size = 400
pyscf.fci.direct_spin0_symm.contract_2e(eri, fcivec, norb, nelec, link_index=None, orbsym=None, wfnsym=0)[source]
pyscf.fci.direct_spin0_symm.get_init_guess(norb, nelec, nroots, hdiag, orbsym, wfnsym=0)[source]
pyscf.fci.direct_spin0_symm.get_init_guess_cyl_sym(norb, nelec, nroots, hdiag, orbsym, wfnsym=0)[source]

pyscf.fci.direct_spin1 module

Full CI solver for spin-free Hamiltonian. This solver can be used to compute doublet, triplet,…

The CI wfn are stored as a 2D array [alpha,beta], where each row corresponds to an alpha string. For each row (alpha string), there are total-num-beta-strings of columns. Each column corresponds to a beta string.

Different FCI solvers are implemented to support different type of symmetry.

Symmetry

File Point group Spin singlet Real hermitian* Alpha/beta degeneracy direct_spin0_symm Yes Yes Yes Yes direct_spin1_symm Yes No Yes Yes direct_spin0 No Yes Yes Yes direct_spin1 No No Yes Yes direct_uhf No No Yes No direct_nosym No No No** Yes

  • Real hermitian Hamiltonian implies (ij|kl) = (ji|kl) = (ij|lk) = (ji|lk)

** Hamiltonian is real but not hermitian, (ij|kl) != (ji|kl) …

pyscf.fci.direct_spin1.FCI

alias of FCISolver

class pyscf.fci.direct_spin1.FCIBase(mol=None)[source]

Bases: StreamObject

Full CI solver

Attributes:
verboseint

Print level. Default value equals to Mole.verbose.

max_cycleint

Total number of iterations. Default is 100

max_spacetuple of int

Davidson iteration space size. Default is 14.

conv_tolfloat

Energy convergence tolerance. Default is 1e-10.

level_shiftfloat

Level shift applied in the preconditioner to avoid singularity. Default is 1e-3

davidson_onlybool

By default, the entire Hamiltonian matrix will be constructed and diagonalized if the system is small (see attribute pspace_size). Setting this parameter to True will enforce the eigenvalue problems being solved by Davidson subspace algorithm. This flag should be enabled when initial guess is given or particular spin symmetry or point-group symmetry is required because the initial guess or symmetry are completely ignored in the direct diagonlization.

pspace_sizeint

The dimension of Hamiltonian matrix over which Davidson iteration algorithm will be used for the eigenvalue problem. Default is 400. This is roughly corresponding to a (6e,6o) system.

nrootsint

Number of states to be solved. Default is 1, the ground state.

spinint or None

Spin (2S = nalpha-nbeta) of the system. If this attribute is None, spin will be determined by the argument nelec (number of electrons) of the kernel function.

wfnsymstr or int

Symmetry of wavefunction. It is used only in direct_spin1_symm and direct_spin0_symm solver.

Saved results

ecifloat or a list of float

FCI energy(ies)

cinparray

FCI wfn vector(s)

convergedbool (or a list of bool for multiple roots)

Whether davidson iteration is converged

Examples:

>>> from pyscf import gto, scf, ao2mo, fci
>>> mol = gto.M(atom='Li 0 0 0; Li 0 0 1', basis='sto-3g')
>>> mf = scf.RHF(mol).run()
>>> h1 = mf.mo_coeff.T.dot(mf.get_hcore()).dot(mf.mo_coeff)
>>> eri = ao2mo.kernel(mol, mf.mo_coeff)
>>> cisolver = fci.direct_spin1.FCI(mol)
>>> e, ci = cisolver.kernel(h1, eri, h1.shape[1], mol.nelec, ecore=mol.energy_nuc())
>>> print(e)
-14.4197890826
absorb_h1e(h1e, eri, norb, nelec, fac=1)[source]

Modify 2e Hamiltonian to include 1e Hamiltonian contribution.

contract_1e(f1e, fcivec, norb, nelec, link_index=None, **kwargs)[source]

Contract the 1-electron Hamiltonian with a FCI vector to get a new FCI vector.

contract_2e(eri, fcivec, norb, nelec, link_index=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()

contract_ss(fcivec, norb, nelec)[source]
conv_tol = 1e-10
conv_tol_residual = None
davidson_only = False
dump_flags(verbose=None)[source]
property e_tot
eig(op, x0=None, precond=None, **kwargs)[source]
energy(h1e, eri, fcivec, norb, nelec, link_index=None)[source]

Compute the FCI electronic energy for given Hamiltonian and FCI vector.

gen_linkstr(norb, nelec, tril=True, spin=None)[source]
get_init_guess(norb, nelec, nroots, hdiag)[source]

Initial guess is the single Slater determinant

kernel(h1e, eri, norb, nelec, ci0=None, tol=None, lindep=None, max_cycle=None, max_space=None, nroots=None, davidson_only=None, pspace_size=None, orbsym=None, wfnsym=None, ecore=0, **kwargs)[source]
Args:
h1e: ndarray

1-electron Hamiltonian

eri: ndarray

2-electron integrals in chemist’s notation

norb: int

Number of orbitals

nelec: int or (int, int)

Number of electrons of the system

Kwargs:
ci0: ndarray

Initial guess

link_index: ndarray

A lookup table to cache the addresses of CI determinants in wave-function vector

tol: float

Convergence tolerance

lindep: float

Linear dependence threshold

max_cycle: int

Max. iterations for diagonalization

max_space: int

Max. trial vectors to store for sub-space diagonalization method

nroots: int

Number of states to solve

davidson_only: bool

Whether to call subspace diagonlization (davidson solver) or do a full diagonlization (lapack eigh) for small systems

pspace_size: int

Number of determinants as the threshold of “small systems”,

hop: function(c) => array_like_c

Function to use for the Hamiltonian multiplication with trial vector

Note: davidson solver requires more arguments. For the parameters not dispatched, they can be passed to davidson solver via the extra keyword arguments **kwargs

large_ci(fcivec, norb, nelec, tol=0.1, return_strs=True)[source]
lessio = False
level_shift = 0.001
lindep = 1e-14
make_hdiag(h1e, eri, norb, nelec, compress=False)[source]

Diagonal Hamiltonian for Davidson preconditioner

Kwargs:

compress (bool) : whether to remove symmetry forbidden elements

make_precond(hdiag, pspaceig=None, pspaceci=None, addr=None)[source]
make_rdm1(fcivec, norb, nelec, link_index=None)[source]

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(fcivec, norb, nelec, link_index=None, reorder=True)[source]

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)

make_rdm12s(fcivec, norb, nelec, link_index=None, reorder=True)[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)

make_rdm1s(fcivec, norb, nelec, link_index=None)[source]

Spin separated 1-particle density matrices. The return values include two density matrices: (alpha,alpha), (beta,beta)

dm1[p,q] = <q^dagger p>

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_rdm2(fcivec, norb, nelec, link_index=None, reorder=True)[source]

Spin traced 2-particle density matrice

NOTE the 2pdm is \(\langle p^\dagger q^\dagger s r\rangle\) but stored as [p,r,q,s]

max_cycle = 100
max_space = 12
property nstates
pspace(h1e, eri, norb, nelec, hdiag=None, np=400)[source]

pspace Hamiltonian to improve Davidson preconditioner. See, CPL, 169, 463

pspace_size = 400
spin_square(fcivec, norb, nelec)[source]

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

threads = None
trans_rdm1(cibra, ciket, norb, nelec, link_index=None)[source]

Spin traced transition 1-particle transition density matrices.

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

trans_rdm12(cibra, ciket, norb, nelec, link_index=None, reorder=True)[source]

Spin traced transition 1- and 2-particle transition density matrices.

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

trans_rdm12s(cibra, ciket, norb, nelec, link_index=None, reorder=True)[source]

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

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

trans_rdm1s(cibra, ciket, norb, nelec, link_index=None)[source]

Spin separated transition 1-particle density matrices. The return values include two density matrices: (alpha,alpha), (beta,beta). See also function make_rdm1s()

1pdm[p,q] = \(\langle q^\dagger p \rangle\)

class pyscf.fci.direct_spin1.FCISolver(mol=None)[source]

Bases: FCIBase

transform_ci_for_orbital_rotation(fcivec, norb, nelec, u)[source]
class pyscf.fci.direct_spin1.FCIvector[source]

Bases: ndarray

An 2D np array for FCI coefficients

pyscf.fci.direct_spin1.absorb_h1e(h1e, eri, norb, nelec, fac=1)[source]

Modify 2e Hamiltonian to include 1e Hamiltonian contribution.

pyscf.fci.direct_spin1.contract_1e(f1e, fcivec, norb, nelec, link_index=None)[source]

Contract the 1-electron Hamiltonian with a FCI vector to get a new FCI vector.

pyscf.fci.direct_spin1.contract_2e(eri, fcivec, norb, nelec, link_index=None)[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()

pyscf.fci.direct_spin1.energy(h1e, eri, fcivec, norb, nelec, link_index=None)[source]

Compute the FCI electronic energy for given Hamiltonian and FCI vector.

pyscf.fci.direct_spin1.get_init_guess(norb, nelec, nroots, hdiag)[source]

Initial guess is the single Slater determinant

pyscf.fci.direct_spin1.kernel(h1e, eri, norb, nelec, ci0=None, level_shift=0.001, tol=1e-10, lindep=1e-14, max_cycle=50, max_space=12, nroots=1, davidson_only=False, pspace_size=400, orbsym=None, wfnsym=None, ecore=0, **kwargs)[source]
pyscf.fci.direct_spin1.kernel_ms1(fci, h1e, eri, norb, nelec, ci0=None, link_index=None, tol=None, lindep=None, max_cycle=None, max_space=None, nroots=None, davidson_only=None, pspace_size=None, hop=None, max_memory=None, verbose=None, ecore=0, **kwargs)[source]
Args:
h1e: ndarray

1-electron Hamiltonian

eri: ndarray

2-electron integrals in chemist’s notation

norb: int

Number of orbitals

nelec: int or (int, int)

Number of electrons of the system

Kwargs:
ci0: ndarray

Initial guess

link_index: ndarray

A lookup table to cache the addresses of CI determinants in wave-function vector

tol: float

Convergence tolerance

lindep: float

Linear dependence threshold

max_cycle: int

Max. iterations for diagonalization

max_space: int

Max. trial vectors to store for sub-space diagonalization method

nroots: int

Number of states to solve

davidson_only: bool

Whether to call subspace diagonlization (davidson solver) or do a full diagonlization (lapack eigh) for small systems

pspace_size: int

Number of determinants as the threshold of “small systems”,

hop: function(c) => array_like_c

Function to use for the Hamiltonian multiplication with trial vector

Note: davidson solver requires more arguments. For the parameters not dispatched, they can be passed to davidson solver via the extra keyword arguments **kwargs

pyscf.fci.direct_spin1.make_diag_precond(hdiag, pspaceig, pspaceci, addr, level_shift=0)[source]
pyscf.fci.direct_spin1.make_hdiag(h1e, eri, norb, nelec, compress=False)[source]

Diagonal Hamiltonian for Davidson preconditioner

Kwargs:

compress (bool) : whether to remove symmetry forbidden elements

pyscf.fci.direct_spin1.make_pspace_precond(hdiag, pspaceig, pspaceci, addr, level_shift=0)[source]
pyscf.fci.direct_spin1.make_rdm1(fcivec, norb, nelec, link_index=None)[source]

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)

pyscf.fci.direct_spin1.make_rdm12(fcivec, norb, nelec, link_index=None, reorder=True)[source]

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)

pyscf.fci.direct_spin1.make_rdm12s(fcivec, norb, nelec, link_index=None, reorder=True)[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)

pyscf.fci.direct_spin1.make_rdm1s(fcivec, norb, nelec, link_index=None)[source]

Spin separated 1-particle density matrices. The return values include two density matrices: (alpha,alpha), (beta,beta)

dm1[p,q] = <q^dagger p>

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)

pyscf.fci.direct_spin1.pspace(h1e, eri, norb, nelec, hdiag=None, np=400)[source]

pspace Hamiltonian to improve Davidson preconditioner. See, CPL, 169, 463

pyscf.fci.direct_spin1.trans_rdm1(cibra, ciket, norb, nelec, link_index=None)[source]

Spin traced transition 1-particle transition density matrices.

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

pyscf.fci.direct_spin1.trans_rdm12(cibra, ciket, norb, nelec, link_index=None, reorder=True)[source]

Spin traced transition 1- and 2-particle transition density matrices.

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

pyscf.fci.direct_spin1.trans_rdm12s(cibra, ciket, norb, nelec, link_index=None, reorder=True)[source]

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

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

pyscf.fci.direct_spin1.trans_rdm1s(cibra, ciket, norb, nelec, link_index=None)[source]

Spin separated transition 1-particle density matrices. The return values include two density matrices: (alpha,alpha), (beta,beta). See also function make_rdm1s()

1pdm[p,q] = \(\langle q^\dagger p \rangle\)

pyscf.fci.direct_spin1_cyl_sym module

Cylindrical symmetry

This module is much slower than direct_spin1_symm.

In this implementation, complex orbitals is used to construct the Hamiltonian. FCI wavefunction (called complex wavefunction here) is solved using the complex Hamiltonian. For 2D irreps, the real part and the imaginary part of the complex FCI wavefunction are identical to the Ex and Ey wavefunction obtained from direct_spin1_symm module. However, any observables from the complex FCI wavefunction should have an identical one from either Ex or Ey wavefunction of direct_spin1_symm.

pyscf.fci.direct_spin1_cyl_sym.FCI

alias of FCISolver

class pyscf.fci.direct_spin1_cyl_sym.FCISolver(mol=None, **kwargs)[source]

Bases: FCISolver

absorb_h1e(h1e, eri, norb, nelec, fac=1)[source]

Modify 2e Hamiltonian to include 1e Hamiltonian contribution.

contract_1e(f1e, fcivec, norb, nelec, link_index=None, **kwargs)[source]

Contract the 1-electron Hamiltonian with a FCI vector to get a new FCI vector.

contract_2e(eri, fcivec, norb, nelec, link_index=None, orbsym=None, wfnsym=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()

get_init_guess(norb, nelec, nroots, hdiag, orbsym=None, wfnsym=None)[source]

Initial guess is the single Slater determinant

guess_wfnsym(norb, nelec, fcivec=None, orbsym=None, wfnsym=None, **kwargs)

Guess point group symmetry of the FCI wavefunction. If fcivec is given, the symmetry of fcivec is used. Otherwise the symmetry is same to the HF determinant.

kernel(h1e, eri, norb, nelec, ci0=None, tol=None, lindep=None, max_cycle=None, max_space=None, nroots=None, davidson_only=None, pspace_size=None, orbsym=None, wfnsym=None, ecore=0, **kwargs)[source]
Args:
h1e: ndarray

1-electron Hamiltonian

eri: ndarray

2-electron integrals in chemist’s notation

norb: int

Number of orbitals

nelec: int or (int, int)

Number of electrons of the system

Kwargs:
ci0: ndarray

Initial guess

link_index: ndarray

A lookup table to cache the addresses of CI determinants in wave-function vector

tol: float

Convergence tolerance

lindep: float

Linear dependence threshold

max_cycle: int

Max. iterations for diagonalization

max_space: int

Max. trial vectors to store for sub-space diagonalization method

nroots: int

Number of states to solve

davidson_only: bool

Whether to call subspace diagonlization (davidson solver) or do a full diagonlization (lapack eigh) for small systems

pspace_size: int

Number of determinants as the threshold of “small systems”,

hop: function(c) => array_like_c

Function to use for the Hamiltonian multiplication with trial vector

Note: davidson solver requires more arguments. For the parameters not dispatched, they can be passed to davidson solver via the extra keyword arguments **kwargs

make_hdiag(h1e, eri, norb, nelec, compress=False)[source]

Diagonal Hamiltonian for Davidson preconditioner

Kwargs:

compress (bool) : whether to remove symmetry forbidden elements

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

make_rdm12s(norb, nelec, link_index=None, reorder=True)

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)

make_rdm1s(norb, nelec, link_index=None)

Spin separated 1-particle density matrices. The return values include two density matrices: (alpha,alpha), (beta,beta)

dm1[p,q] = <q^dagger p>

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)

pspace(h1e, eri, norb, nelec, hdiag, np=400)[source]

pspace Hamiltonian to improve Davidson preconditioner. See, CPL, 169, 463

trans_rdm1(ciket, norb, nelec, link_index=None)

Spin traced transition 1-particle transition density matrices.

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

trans_rdm12(ciket, norb, nelec, link_index=None, reorder=True)

Spin traced transition 1- and 2-particle transition density matrices.

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

trans_rdm12s(ciket, norb, nelec, link_index=None, reorder=True)

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

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

trans_rdm1s(ciket, norb, nelec, link_index=None)

Spin separated transition 1-particle density matrices. The return values include two density matrices: (alpha,alpha), (beta,beta). See also function make_rdm1s()

1pdm[p,q] = \(\langle q^\dagger p \rangle\)

pyscf.fci.direct_spin1_cyl_sym.argsort_strs_by_irrep(strs, orbsym, max_momentum, max_gerades)[source]
pyscf.fci.direct_spin1_cyl_sym.contract_2e(eri, fcivec, norb, nelec, link_index=None, orbsym=None, wfnsym=0)[source]
pyscf.fci.direct_spin1_cyl_sym.gen_str_irrep(strs, orbsym, link_index, rank_eri, irrep_eri, max_momentum, max_gerades)[source]
pyscf.fci.direct_spin1_cyl_sym.get_init_guess(norb, nelec, nroots, hdiag, orbsym, wfnsym=0, sym_allowed_idx=None)[source]
pyscf.fci.direct_spin1_cyl_sym.guess_wfnsym(solver, norb, nelec, fcivec=None, orbsym=None, wfnsym=None, **kwargs)[source]

Guess point group symmetry of the FCI wavefunction. If fcivec is given, the symmetry of fcivec is used. Otherwise the symmetry is same to the HF determinant.

pyscf.fci.direct_spin1_cyl_sym.reorder_eri(eri, norb, orbsym, max_momentum, max_gerades)[source]
pyscf.fci.direct_spin1_cyl_sym.sym_allowed_indices(nelec, orbsym, wfnsym)[source]

Indices of symmetry allowed determinants for each irrep

pyscf.fci.direct_spin1_symm module

Different FCI solvers are implemented to support different type of symmetry.

Symmetry

File Point group Spin singlet Real hermitian* Alpha/beta degeneracy direct_spin0_symm Yes Yes Yes Yes direct_spin1_symm Yes No Yes Yes direct_spin0 No Yes Yes Yes direct_spin1 No No Yes Yes direct_uhf No No Yes No direct_nosym No No No** Yes

  • Real hermitian Hamiltonian implies (ij|kl) = (ji|kl) = (ij|lk) = (ji|lk)

** Hamiltonian is real but not hermitian, (ij|kl) != (ji|kl) …

pyscf.fci.direct_spin1_symm.FCI

alias of FCISolver

class pyscf.fci.direct_spin1_symm.FCISolver(mol=None, **kwargs)[source]

Bases: FCISolver

absorb_h1e(h1e, eri, norb, nelec, fac=1)

Modify 2e Hamiltonian to include 1e Hamiltonian contribution.

contract_1e(f1e, fcivec, norb, nelec, link_index=None, **kwargs)[source]

Contract the 1-electron Hamiltonian with a FCI vector to get a new FCI vector.

contract_2e(eri, fcivec, norb, nelec, link_index=None, orbsym=None, wfnsym=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()

contract_ss(fcivec, norb, nelec)[source]
dump_flags(verbose=None)[source]
get_init_guess(norb, nelec, nroots, hdiag, orbsym=None, wfnsym=None)[source]

Initial guess is the single Slater determinant

guess_wfnsym(norb, nelec, fcivec=None, orbsym=None, wfnsym=None, **kwargs)

Guess point group symmetry of the FCI wavefunction. If fcivec is given, the symmetry of fcivec is used. Otherwise the symmetry is same to the HF determinant.

kernel(h1e, eri, norb, nelec, ci0=None, tol=None, lindep=None, max_cycle=None, max_space=None, nroots=None, davidson_only=None, pspace_size=None, orbsym=None, wfnsym=None, ecore=0, **kwargs)[source]
Args:
h1e: ndarray

1-electron Hamiltonian

eri: ndarray

2-electron integrals in chemist’s notation

norb: int

Number of orbitals

nelec: int or (int, int)

Number of electrons of the system

Kwargs:
ci0: ndarray

Initial guess

link_index: ndarray

A lookup table to cache the addresses of CI determinants in wave-function vector

tol: float

Convergence tolerance

lindep: float

Linear dependence threshold

max_cycle: int

Max. iterations for diagonalization

max_space: int

Max. trial vectors to store for sub-space diagonalization method

nroots: int

Number of states to solve

davidson_only: bool

Whether to call subspace diagonlization (davidson solver) or do a full diagonlization (lapack eigh) for small systems

pspace_size: int

Number of determinants as the threshold of “small systems”,

hop: function(c) => array_like_c

Function to use for the Hamiltonian multiplication with trial vector

Note: davidson solver requires more arguments. For the parameters not dispatched, they can be passed to davidson solver via the extra keyword arguments **kwargs

make_hdiag(h1e, eri, norb, nelec, compress=False)[source]

Diagonal Hamiltonian for Davidson preconditioner

Kwargs:

compress (bool) : whether to remove symmetry forbidden elements

pspace(h1e, eri, norb, nelec, hdiag, np=400)[source]

pspace Hamiltonian to improve Davidson preconditioner. See, CPL, 169, 463

pspace_size = 400
pyscf.fci.direct_spin1_symm.argsort_strs_by_irrep(strs, orbsym)[source]
pyscf.fci.direct_spin1_symm.contract_2e(eri, fcivec, norb, nelec, link_index=None, orbsym=None, wfnsym=0)[source]
pyscf.fci.direct_spin1_symm.energy(h1e, eri, fcivec, norb, nelec, link_index=None, orbsym=None, wfnsym=0)[source]
pyscf.fci.direct_spin1_symm.gen_str_irrep(strs, orbsym, link_index, rank_eri, irrep_eri)[source]
pyscf.fci.direct_spin1_symm.get_init_guess(norb, nelec, nroots, hdiag, orbsym, wfnsym=0)[source]
pyscf.fci.direct_spin1_symm.get_init_guess_cyl_sym(norb, nelec, nroots, hdiag, orbsym, wfnsym=0)[source]
pyscf.fci.direct_spin1_symm.guess_wfnsym(solver, norb, nelec, fcivec=None, orbsym=None, wfnsym=None, **kwargs)[source]

Guess point group symmetry of the FCI wavefunction. If fcivec is given, the symmetry of fcivec is used. Otherwise the symmetry is same to the HF determinant.

pyscf.fci.direct_spin1_symm.kernel(h1e, eri, norb, nelec, ci0=None, level_shift=0.001, tol=1e-10, lindep=1e-14, max_cycle=50, max_space=12, nroots=1, davidson_only=False, pspace_size=400, orbsym=None, wfnsym=None, ecore=0, **kwargs)[source]
pyscf.fci.direct_spin1_symm.reorder_eri(eri, norb, orbsym)[source]
pyscf.fci.direct_spin1_symm.sym_allowed_indices(nelec, orbsym, wfnsym)[source]

Indices of symmetry allowed determinants for each irrep

pyscf.fci.direct_uhf module

Different FCI solvers are implemented to support different type of symmetry.

Symmetry

File Point group Spin singlet Real hermitian* Alpha/beta degeneracy direct_spin0_symm Yes Yes Yes Yes direct_spin1_symm Yes No Yes Yes direct_spin0 No Yes Yes Yes direct_spin1 No No Yes Yes direct_uhf No No Yes No direct_nosym No No No** Yes

  • Real hermitian Hamiltonian implies (ij|kl) = (ji|kl) = (ij|lk) = (ji|lk)

** Hamiltonian is real but not hermitian, (ij|kl) != (ji|kl) …

pyscf.fci.direct_uhf.FCI

alias of FCISolver

class pyscf.fci.direct_uhf.FCISolver(mol=None)[source]

Bases: FCISolver

static absorb_h1e(h1e, eri, norb, nelec, fac=1)

Modify 2e Hamiltonian to include 1e Hamiltonian contribution.

contract_1e(f1e=None, fcivec=None, norb=None, nelec=None, link_index=None)

Contract the 1-electron Hamiltonian with a FCI vector to get a new FCI vector.

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

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

static make_hdiag(h1e, eri, norb, nelec, compress=False)

Diagonal Hamiltonian for Davidson preconditioner

Kwargs:

compress (bool) : whether to remove symmetry forbidden elements

make_rdm1(fcivec=None, norb=None, nelec=None, 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)

static pspace(h1e, eri, norb, nelec, hdiag=None, np=400)

pspace Hamiltonian to improve Davidson preconditioner. See, CPL, 169, 463

spin_square(fcivec=None, norb=None, nelec=None, mo_coeff=None, ovlp=1, fcisolver=None)

General spin square operator.

… math:

<CI|S_+*S_-|CI> &= n_\alpha + \delta_{ik}\delta_{jl}Gamma_{i\alpha k\beta ,j\beta l\alpha } \\
<CI|S_-*S_+|CI> &= n_\beta + \delta_{ik}\delta_{jl}Gamma_{i\beta k\alpha ,j\alpha l\beta } \\
<CI|S_z*S_z|CI> &= \delta_{ik}\delta_{jl}(Gamma_{i\alpha k\alpha ,j\alpha l\alpha }
                 - Gamma_{i\alpha k\alpha ,j\beta l\beta }
                 - Gamma_{i\beta k\beta ,j\alpha l\alpha}
                 + Gamma_{i\beta k\beta ,j\beta l\beta})
                 + (n_\alpha+n_\beta)/4

Given the overlap between non-degenerate alpha and beta orbitals, this function can compute the expectation value spin square operator for UHF-FCI wavefunction

trans_rdm1(cibra=None, ciket=None, norb=None, nelec=None, link_index=None)

Spin traced transition 1-particle transition density matrices.

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

pyscf.fci.direct_uhf.absorb_h1e(h1e, eri, norb, nelec, fac=1)[source]
pyscf.fci.direct_uhf.contract_1e(f1e, fcivec, norb, nelec, link_index=None)[source]
pyscf.fci.direct_uhf.contract_2e(eri, fcivec, norb, nelec, link_index=None)[source]
pyscf.fci.direct_uhf.contract_2e_hubbard(u, fcivec, norb, nelec, opt=None)[source]
pyscf.fci.direct_uhf.energy(h1e, eri, fcivec, norb, nelec, link_index=None)[source]
pyscf.fci.direct_uhf.kernel(h1e, eri, norb, nelec, ci0=None, level_shift=0.001, tol=1e-10, lindep=1e-14, max_cycle=50, max_space=12, nroots=1, davidson_only=False, pspace_size=400, orbsym=None, wfnsym=None, ecore=0, **kwargs)[source]
pyscf.fci.direct_uhf.make_hdiag(h1e, eri, norb, nelec, compress=False)[source]
pyscf.fci.direct_uhf.make_rdm1(fcivec, norb, nelec, link_index=None)[source]
pyscf.fci.direct_uhf.pspace(h1e, eri, norb, nelec, hdiag=None, np=400)[source]
pyscf.fci.direct_uhf.trans_rdm1(cibra, ciket, norb, nelec, link_index=None)[source]

pyscf.fci.fci_dhf_slow module

pyscf.fci.fci_dhf_slow.FCI

alias of FCISolver

class pyscf.fci.fci_dhf_slow.FCISolver(mol=None)[source]

Bases: FCISolver

absorb_h1e(h1e, eri, norb, nelec, fac=1)[source]

Modify 2e Hamiltonian to include 1e Hamiltonian contribution.

contract_1e(f1e, fcivec, norb, nelec, link_index=None, **kwargs)[source]

Contract the 1-electron Hamiltonian with a FCI vector to get a new FCI vector.

contract_2e(eri, fcivec, norb, nelec, link_index=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()

energy(h1e, eri, fcivec, norb, nelec, link_index=None)[source]

Compute the FCI electronic energy for given Hamiltonian and FCI vector.

kernel(h1e, eri, norb, nelec, ci0=None, tol=None, lindep=None, max_cycle=None, max_space=None, nroots=None, davidson_only=None, pspace_size=None, orbsym=None, wfnsym=None, ecore=0, **kwargs)[source]
Args:
h1e: ndarray

1-electron Hamiltonian

eri: ndarray

2-electron integrals in chemist’s notation

norb: int

Number of orbitals

nelec: int or (int, int)

Number of electrons of the system

Kwargs:
ci0: ndarray

Initial guess

link_index: ndarray

A lookup table to cache the addresses of CI determinants in wave-function vector

tol: float

Convergence tolerance

lindep: float

Linear dependence threshold

max_cycle: int

Max. iterations for diagonalization

max_space: int

Max. trial vectors to store for sub-space diagonalization method

nroots: int

Number of states to solve

davidson_only: bool

Whether to call subspace diagonlization (davidson solver) or do a full diagonlization (lapack eigh) for small systems

pspace_size: int

Number of determinants as the threshold of “small systems”,

hop: function(c) => array_like_c

Function to use for the Hamiltonian multiplication with trial vector

Note: davidson solver requires more arguments. For the parameters not dispatched, they can be passed to davidson solver via the extra keyword arguments **kwargs

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

Diagonal Hamiltonian for Davidson preconditioner

Kwargs:

compress (bool) : whether to remove symmetry forbidden elements

make_rdm1(fcivec, norb, nelec, link_index=None)[source]

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(fcivec, norb, nelec, link_index=None, reorder=True)[source]

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)

make_rdm2(fcivec, norb, nelec, link_index=None, reorder=True)[source]

Spin traced 2-particle density matrice

NOTE the 2pdm is \(\langle p^\dagger q^\dagger s r\rangle\) but stored as [p,r,q,s]

pspace(h1e, eri, norb, nelec, hdiag, np=400)[source]

pspace Hamiltonian to improve Davidson preconditioner. See, CPL, 169, 463

spin_square(fcivec, norb, nelec)[source]

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

pyscf.fci.fci_dhf_slow.absorb_h1e(h1e, eri, norb, nelec, fac=1)[source]

Modify 2e Hamiltonian to include 1e Hamiltonian contribution.

pyscf.fci.fci_dhf_slow.contract_2e(eri, fcivec, norb, nelec, opt=None)[source]

Compute a^dagger_p a_q a^dagger_r a_s |CI>

pyscf.fci.fci_dhf_slow.kernel(h1e, eri, norb, nelec, ecore=0, nroots=1, verbose=3)[source]
pyscf.fci.fci_dhf_slow.kernel_dhf(fci, h1e, eri, norb, nelec, ci0=None, link_index=None, tol=None, lindep=None, max_cycle=None, max_space=None, nroots=None, davidson_only=None, pspace_size=None, max_memory=None, verbose=None, ecore=0, **kwargs)[source]
pyscf.fci.fci_dhf_slow.make_hdiag(h1e, eri, norb, nelec, opt=None)[source]
pyscf.fci.fci_dhf_slow.make_rdm1(fcivec, norb, nelec, link_index=None)[source]
pyscf.fci.fci_dhf_slow.make_rdm12(fcivec, norb, nelec, link_index=None, reorder=True)[source]
pyscf.fci.fci_dhf_slow.reorder_rdm(rdm1, rdm2, inplace=True)[source]

pyscf.fci.fci_slow module

pyscf.fci.fci_slow.absorb_h1e(h1e, eri, norb, nelec, fac=1)[source]

Modify 2e Hamiltonian to include 1e Hamiltonian contribution.

pyscf.fci.fci_slow.contract_1e(f1e, fcivec, norb, nelec)[source]
pyscf.fci.fci_slow.contract_2e(eri, fcivec, norb, nelec, link_index=None)[source]

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

pyscf.fci.fci_slow.contract_2e_hubbard(u, fcivec, norb, nelec, opt=None)[source]
pyscf.fci.fci_slow.kernel(h1e, eri, norb, nelec, ecore=0)[source]
pyscf.fci.fci_slow.make_hdiag(h1e, eri, norb, nelec, opt=None)[source]
pyscf.fci.fci_slow.make_rdm1(fcivec, norb, nelec, opt=None)[source]
pyscf.fci.fci_slow.make_rdm12(fcivec, norb, nelec, opt=None)[source]
pyscf.fci.fci_slow.reorder_rdm(rdm1, rdm2, inplace=True)[source]

reorder from rdm2(pq,rs) = <E^p_q E^r_s> to rdm2(pq,rs) = <e^{pr}_{qs}>. Although the “reoredered rdm2” is still in Mulliken order (rdm2[e1,e1,e2,e2]), it is the true 2e DM (dotting it with int2e gives the energy of 2e parts)

pyscf.fci.rdm module

FCI 1, 2, 3, 4-particle density matrices.

Note the 1-particle density matrix has the same convention as the mean-field 1-particle density matrix (see McWeeney’s book Eq 5.4.20), which is

dm[p,q] = < q^+ p >

The contraction between 1-particle Hamiltonian and 1-pdm is

E = einsum(‘pq,qp’, h1, 1pdm)

Different conventions are used in the high order density matrices:

dm[p,q,r,s,…] = < p^+ r^+ … s q >

pyscf.fci.rdm.make_dm123(fname, cibra, ciket, norb, nelec)[source]

Spin traced 1, 2 and 3-particle density matrices.

Note

In this function, 2pdm[p,q,r,s] is \(\langle p^\dagger q r^\dagger s\rangle\); 3pdm[p,q,r,s,t,u] is \(\langle p^\dagger q r^\dagger s t^\dagger u\rangle\).

After calling reorder_dm123, the 2pdm and 3pdm are transformed to the normal density matrices: 2pdm[p,r,q,s] = \(\langle p^\dagger q^\dagger s r\rangle\) 3pdm[p,s,q,t,r,u] = \(\langle p^\dagger q^\dagger r^\dagger u t s\rangle\).

pyscf.fci.rdm.make_dm1234(fname, cibra, ciket, norb, nelec)[source]

Spin traced 1, 2, 3 and 4-particle density matrices.

Note

In this function, 2pdm[p,q,r,s] is \(\langle p^\dagger q r^\dagger s\rangle\); 3pdm[p,q,r,s,t,u] is \(\langle p^\dagger q r^\dagger s t^\dagger u\rangle\); 4pdm[p,q,r,s,t,u,v,w] is \(\langle p^\dagger q r^\dagger s t^\dagger u v^\dagger w\rangle\).

After calling reorder_dm123, the 2pdm and 3pdm are transformed to the normal density matrices: 2pdm[p,r,q,s] = \(\langle p^\dagger q^\dagger s r\rangle\) 3pdm[p,s,q,t,r,u] = \(\langle p^\dagger q^\dagger r^\dagger u t s\rangle\). 4pdm[p,t,q,u,r,v,s,w] = \(\langle p^\dagger q^\dagger r^\dagger s^dagger w v u t\rangle\).

pyscf.fci.rdm.make_rdm1(fname, cibra, ciket, norb, nelec, link_index=None)
pyscf.fci.rdm.make_rdm12(fname, cibra, ciket, norb, nelec, link_index=None, symm=0)
pyscf.fci.rdm.make_rdm12_ms0(fname, cibra, ciket, norb, nelec, link_index=None, symm=0)[source]
pyscf.fci.rdm.make_rdm12_spin1(fname, cibra, ciket, norb, nelec, link_index=None, symm=0)[source]
pyscf.fci.rdm.make_rdm1_ms0(fname, cibra, ciket, norb, nelec, link_index=None)[source]
pyscf.fci.rdm.make_rdm1_spin1(fname, cibra, ciket, norb, nelec, link_index=None)[source]
pyscf.fci.rdm.reorder_dm12(rdm1, rdm2, inplace=True)[source]
pyscf.fci.rdm.reorder_dm123(rdm1, rdm2, rdm3, inplace=True)[source]
pyscf.fci.rdm.reorder_dm1234(rdm1, rdm2, rdm3, rdm4, inplace=True)[source]
pyscf.fci.rdm.reorder_rdm(rdm1, rdm2, inplace=False)[source]

pyscf.fci.selected_ci module

Selected CI.

This is an inefficient dialect of Selected CI using the same structure as determinant based FCI algorithm. For the efficient Selected CI programs, Dice program (https://github.com/sanshar/Dice.git) is a good candidate.

Simple usage:

>>> from pyscf import gto, scf, ao2mo, fci
>>> mol = gto.M(atom='C 0 0 0; C 0 0 1')
>>> mf = scf.RHF(mol).run()
>>> h1 = mf.mo_coeff.T.dot(mf.get_hcore()).dot(mf.mo_coeff)
>>> h2 = ao2mo.kernel(mol, mf.mo_coeff)
>>> e = fci.selected_ci.kernel(h1, h2, mf.mo_coeff.shape[1], mol.nelectron)[0]
pyscf.fci.selected_ci.SCI

alias of SelectedCI

class pyscf.fci.selected_ci.SCIvector[source]

Bases: ndarray

An 2D np array for selected CI coefficients

class pyscf.fci.selected_ci.SelectedCI(mol=None)[source]

Bases: FCISolver

ci_coeff_cutoff = 0.0005
contract_2e(eri, civec_strs, norb, nelec, link_index=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()

contract_ss(fcivec, norb, nelec)[source]
conv_tol = 1e-09
dump_flags(verbose=None)[source]
enlarge_space(civec_strs, eri, norb, nelec)
gen_linkstr(norb, nelec, tril=True, spin=None, ci_strs=None)[source]
get_init_guess(ci_strs, norb, nelec, nroots, hdiag)[source]

Initial guess is the single Slater determinant

kernel(h1e, eri, norb, nelec, ci0=None, tol=None, lindep=None, max_cycle=None, max_space=None, nroots=None, davidson_only=None, max_memory=None, verbose=None, ecore=0, **kwargs)
Args:
h1e: ndarray

1-electron Hamiltonian

eri: ndarray

2-electron integrals in chemist’s notation

norb: int

Number of orbitals

nelec: int or (int, int)

Number of electrons of the system

Kwargs:
ci0: ndarray

Initial guess

link_index: ndarray

A lookup table to cache the addresses of CI determinants in wave-function vector

tol: float

Convergence tolerance

lindep: float

Linear dependence threshold

max_cycle: int

Max. iterations for diagonalization

max_space: int

Max. trial vectors to store for sub-space diagonalization method

nroots: int

Number of states to solve

davidson_only: bool

Whether to call subspace diagonlization (davidson solver) or do a full diagonlization (lapack eigh) for small systems

pspace_size: int

Number of determinants as the threshold of “small systems”,

hop: function(c) => array_like_c

Function to use for the Hamiltonian multiplication with trial vector

Note: davidson solver requires more arguments. For the parameters not dispatched, they can be passed to davidson solver via the extra keyword arguments **kwargs

kernel_fixed_space(h1e, eri, norb, nelec, ci_strs, ci0=None, tol=None, lindep=None, max_cycle=None, max_space=None, nroots=None, davidson_only=None, max_memory=None, verbose=None, ecore=0, **kwargs)
large_ci(civec_strs, norb, nelec, tol=0.1, return_strs=True)[source]
static make_hdiag(h1e, eri, ci_strs, norb, nelec, compress=False)

Diagonal Hamiltonian for Davidson preconditioner

Kwargs:

compress (bool) : whether to remove symmetry forbidden elements

make_rdm1(civec_strs, norb, nelec, link_index=None)[source]

Spin-traced 1-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_strs, norb, nelec, link_index=None, **kwargs)[source]

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)

make_rdm12s(civec_strs, norb, nelec, link_index=None, **kwargs)[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)

make_rdm1s(civec_strs, norb, nelec, link_index=None)[source]

Spin separated 1-particle density matrices. The return values include two density matrices: (alpha,alpha), (beta,beta)

dm1[p,q] = <q^dagger p>

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_rdm2(civec_strs, norb, nelec, link_index=None, **kwargs)[source]

Spin-traced two-particle density matrix.

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`.

make_rdm2s(civec_strs, norb, nelec, link_index=None, **kwargs)[source]

Spin separated 2-particle density matrices. The return values include three density matrices: (alpha,alpha,alpha,alpha), (alpha,alpha,beta,beta), (beta,beta,beta,beta)

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

select_cutoff = 0.0005
spin_square(civec_strs, norb, nelec)[source]

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

start_tol = 0.0003
tol_decay_rate = 0.3
trans_rdm1(cibra, ciket, norb, nelec, link_index=None)[source]

Spin traced transition 1-particle density matrices. See also function make_rdm1()

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

trans_rdm1s(cibra, ciket, norb, nelec, link_index=None)[source]

Spin separated transition 1-particle density matrices. See also function make_rdm1s()

1pdm[p,q] = \(\langle q^\dagger p \rangle\)

pyscf.fci.selected_ci.contract_2e(eri, civec_strs, norb, nelec, link_index=None)[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()

pyscf.fci.selected_ci.contract_ss(civec_strs, norb, nelec)[source]

S^2 |Psirangle

pyscf.fci.selected_ci.cre_des_linkstr(strs, norb, nelec, tril=False)[source]

Given intermediates, the link table to generate input strs

pyscf.fci.selected_ci.cre_des_linkstr_tril(strs, norb, nelec)[source]

Given intermediates, the link table to generate input strs

pyscf.fci.selected_ci.des_des_linkstr(strs, norb, nelec, tril=False)[source]

Given intermediates, the link table to generate input strs

pyscf.fci.selected_ci.des_des_linkstr_tril(strs, norb, nelec)[source]

Given intermediates, the link table to generate input strs

pyscf.fci.selected_ci.enlarge_space(myci, civec_strs, eri, norb, nelec)[source]
pyscf.fci.selected_ci.from_fci(fcivec, ci_strs, norb, nelec)[source]
pyscf.fci.selected_ci.gen_cre_linkstr(strs, norb, nelec)[source]

Given intermediates, the link table to generate input strs

pyscf.fci.selected_ci.gen_des_linkstr(strs, norb, nelec)[source]

Given intermediates, the link table to generate input strs

pyscf.fci.selected_ci.kernel(h1e, eri, norb, nelec, ci0=None, level_shift=0.001, tol=1e-10, lindep=1e-14, max_cycle=50, max_space=12, nroots=1, davidson_only=False, pspace_size=400, orbsym=None, wfnsym=None, select_cutoff=0.001, ci_coeff_cutoff=0.001, ecore=0, **kwargs)[source]
pyscf.fci.selected_ci.kernel_fixed_space(myci, h1e, eri, norb, nelec, ci_strs, ci0=None, tol=None, lindep=None, max_cycle=None, max_space=None, nroots=None, davidson_only=None, max_memory=None, verbose=None, ecore=0, **kwargs)[source]
pyscf.fci.selected_ci.kernel_float_space(myci, h1e, eri, norb, nelec, ci0=None, tol=None, lindep=None, max_cycle=None, max_space=None, nroots=None, davidson_only=None, max_memory=None, verbose=None, ecore=0, **kwargs)[source]
pyscf.fci.selected_ci.make_hdiag(h1e, eri, ci_strs, norb, nelec, compress=False)[source]
pyscf.fci.selected_ci.make_rdm1(civec_strs, norb, nelec, link_index=None)[source]

Spin-traced 1-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)

pyscf.fci.selected_ci.make_rdm1s(civec_strs, norb, nelec, link_index=None)[source]

Spin separated 1-particle density matrices. The return values include two density matrices: (alpha,alpha), (beta,beta)

dm1[p,q] = <q^dagger p>

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)

pyscf.fci.selected_ci.make_rdm2(civec_strs, norb, nelec, link_index=None, **kwargs)[source]

Spin-traced two-particle density matrix.

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`.

pyscf.fci.selected_ci.make_rdm2s(civec_strs, norb, nelec, link_index=None, **kwargs)[source]

Spin separated 2-particle density matrices. The return values include three density matrices: (alpha,alpha,alpha,alpha), (alpha,alpha,beta,beta), (beta,beta,beta,beta)

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

pyscf.fci.selected_ci.select_strs(myci, eri, eri_pq_max, civec_max, strs, norb, nelec)[source]
pyscf.fci.selected_ci.spin_square(civec_strs, norb, nelec)[source]

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

pyscf.fci.selected_ci.to_fci(civec_strs, norb, nelec)[source]
pyscf.fci.selected_ci.trans_rdm1(cibra_strs, ciket_strs, norb, nelec, link_index=None)[source]

Spin traced transition 1-particle density matrices. See also function make_rdm1()

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

pyscf.fci.selected_ci.trans_rdm1s(cibra_strs, ciket_strs, norb, nelec, link_index=None)[source]

Spin separated transition 1-particle density matrices. See also function make_rdm1s()

1pdm[p,q] = \(\langle q^\dagger p \rangle\)

pyscf.fci.selected_ci_slow module

class pyscf.fci.selected_ci_slow.SelectedCI[source]

Bases: object

pyscf.fci.selected_ci_slow.contract_2e(eri, civec_strs, norb, nelec, link_index=None)[source]

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

pyscf.fci.selected_ci_slow.cre_des_linkstr(strs, norb, nelec)[source]
pyscf.fci.selected_ci_slow.des_des_linkstr(strs, norb, nelec)[source]
pyscf.fci.selected_ci_slow.enlarge_space(myci, civec_strs, eri, norb, nelec)[source]
pyscf.fci.selected_ci_slow.kernel(h1e, eri, norb, nelec, ecore=0, verbose=3)[source]
pyscf.fci.selected_ci_slow.make_hdiag(h1e, g2e, ci_strs, norb, nelec)[source]
pyscf.fci.selected_ci_slow.make_rdm1(civec_strs, norb, nelec)[source]
pyscf.fci.selected_ci_slow.make_rdm2(civec_strs, norb, nelec)[source]

pyscf.fci.selected_ci_spin0 module

pyscf.fci.selected_ci_spin0.SCI

alias of SelectedCI

class pyscf.fci.selected_ci_spin0.SelectedCI(mol=None)[source]

Bases: SelectedCI

contract_2e(eri, civec_strs, norb, nelec, link_index=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()

enlarge_space(civec_strs, eri, norb, nelec)
static make_hdiag(h1e, eri, ci_strs, norb, nelec, compress=False)

Diagonal Hamiltonian for Davidson preconditioner

Kwargs:

compress (bool) : whether to remove symmetry forbidden elements

pyscf.fci.selected_ci_spin0.contract_2e(eri, civec_strs, norb, nelec, link_index=None)[source]
pyscf.fci.selected_ci_spin0.enlarge_space(myci, civec_strs, eri, norb, nelec)[source]
pyscf.fci.selected_ci_spin0.kernel(h1e, eri, norb, nelec, ci0=None, level_shift=0.001, tol=1e-10, lindep=1e-14, max_cycle=50, max_space=12, nroots=1, davidson_only=False, pspace_size=400, orbsym=None, wfnsym=None, select_cutoff=0.001, ci_coeff_cutoff=0.001, ecore=0, **kwargs)[source]
pyscf.fci.selected_ci_spin0.make_hdiag(h1e, eri, ci_strs, norb, nelec, compress=False)[source]

pyscf.fci.selected_ci_spin0_symm module

pyscf.fci.selected_ci_spin0_symm.SCI

alias of SelectedCI

class pyscf.fci.selected_ci_spin0_symm.SelectedCI(mol=None)[source]

Bases: SelectedCI

contract_2e(eri, civec_strs, norb, nelec, link_index=None, orbsym=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()

enlarge_space(civec_strs, eri, norb, nelec)
static make_hdiag(h1e, eri, ci_strs, norb, nelec, compress=False)

Diagonal Hamiltonian for Davidson preconditioner

Kwargs:

compress (bool) : whether to remove symmetry forbidden elements

pyscf.fci.selected_ci_spin0_symm.contract_2e(eri, civec_strs, norb, nelec, link_index=None, orbsym=None)[source]
pyscf.fci.selected_ci_spin0_symm.kernel(h1e, eri, norb, nelec, ci0=None, level_shift=0.001, tol=1e-10, lindep=1e-14, max_cycle=50, max_space=12, nroots=1, davidson_only=False, pspace_size=400, orbsym=None, wfnsym=None, select_cutoff=0.001, ci_coeff_cutoff=0.001, ecore=0, **kwargs)[source]

pyscf.fci.selected_ci_symm module

pyscf.fci.selected_ci_symm.SCI

alias of SelectedCI

class pyscf.fci.selected_ci_symm.SelectedCI(mol=None)[source]

Bases: SelectedCI

contract_2e(eri, civec_strs, norb, nelec, link_index=None, orbsym=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()

get_init_guess(ci_strs, norb, nelec, nroots, hdiag)[source]

Initial guess is the single Slater determinant

guess_wfnsym(norb, nelec, fcivec=None, orbsym=None, wfnsym=None, **kwargs)[source]
kernel(h1e, eri, norb, nelec, ci0=None, tol=None, lindep=None, max_cycle=None, max_space=None, nroots=None, davidson_only=None, pspace_size=None, orbsym=None, wfnsym=None, ecore=0, **kwargs)[source]
Args:
h1e: ndarray

1-electron Hamiltonian

eri: ndarray

2-electron integrals in chemist’s notation

norb: int

Number of orbitals

nelec: int or (int, int)

Number of electrons of the system

Kwargs:
ci0: ndarray

Initial guess

link_index: ndarray

A lookup table to cache the addresses of CI determinants in wave-function vector

tol: float

Convergence tolerance

lindep: float

Linear dependence threshold

max_cycle: int

Max. iterations for diagonalization

max_space: int

Max. trial vectors to store for sub-space diagonalization method

nroots: int

Number of states to solve

davidson_only: bool

Whether to call subspace diagonlization (davidson solver) or do a full diagonlization (lapack eigh) for small systems

pspace_size: int

Number of determinants as the threshold of “small systems”,

hop: function(c) => array_like_c

Function to use for the Hamiltonian multiplication with trial vector

Note: davidson solver requires more arguments. For the parameters not dispatched, they can be passed to davidson solver via the extra keyword arguments **kwargs

pyscf.fci.selected_ci_symm.contract_2e(eri, civec_strs, norb, nelec, link_index=None, orbsym=None)[source]
pyscf.fci.selected_ci_symm.kernel(h1e, eri, norb, nelec, ci0=None, level_shift=0.001, tol=1e-10, lindep=1e-14, max_cycle=50, max_space=12, nroots=1, davidson_only=False, pspace_size=400, orbsym=None, wfnsym=None, select_cutoff=0.001, ci_coeff_cutoff=0.001, ecore=0, **kwargs)[source]
pyscf.fci.selected_ci_symm.reorder4irrep(eri, norb, link_index, orbsym, offdiag=0)[source]

pyscf.fci.spin_op module

pyscf.fci.spin_op.contract_ss(fcivec, norb, nelec)[source]

Contract spin square operator with FCI wavefunction \(S^2 |CI>\)

pyscf.fci.spin_op.local_spin(fcivec, norb, nelec, mo_coeff=None, ovlp=1, aolst=[])[source]

Local spin expectation value, which is defined as

<CI|(local S^2)|CI>

The local S^2 operator only couples the orbitals specified in aolst. The cross term which involves the interaction between the local part (in aolst) and non-local part (not in aolst) is not included. As a result, the value of local_spin is not additive. In other words, if local_spin is computed twice with the complementary aolst in the two runs, the summation does not equal to the S^2 of the entire system.

For a complete list of AOs, the value of local_spin is equivalent to <CI|S^2|CI>

pyscf.fci.spin_op.make_rdm2_abba(fcivec, norb, nelec)[source]
pyscf.fci.spin_op.make_rdm2_baab(fcivec, norb, nelec)[source]
pyscf.fci.spin_op.spin_square(fcivec, norb, nelec, mo_coeff=None, ovlp=1, fcisolver=None)[source]

General spin square operator.

… math:

<CI|S_+*S_-|CI> &= n_\alpha + \delta_{ik}\delta_{jl}Gamma_{i\alpha k\beta ,j\beta l\alpha } \\
<CI|S_-*S_+|CI> &= n_\beta + \delta_{ik}\delta_{jl}Gamma_{i\beta k\alpha ,j\alpha l\beta } \\
<CI|S_z*S_z|CI> &= \delta_{ik}\delta_{jl}(Gamma_{i\alpha k\alpha ,j\alpha l\alpha }
                 - Gamma_{i\alpha k\alpha ,j\beta l\beta }
                 - Gamma_{i\beta k\beta ,j\alpha l\alpha}
                 + Gamma_{i\beta k\beta ,j\beta l\beta})
                 + (n_\alpha+n_\beta)/4

Given the overlap between non-degenerate alpha and beta orbitals, this function can compute the expectation value spin square operator for UHF-FCI wavefunction

pyscf.fci.spin_op.spin_square0(fcivec, norb, nelec)[source]

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

pyscf.fci.spin_op.spin_square_general(dm1a, dm1b, dm2aa, dm2ab, dm2bb, mo_coeff, ovlp=1)[source]

General spin square operator.

… math:

<CI|S_+*S_-|CI> &= n_\alpha + \delta_{ik}\delta_{jl}Gamma_{i\alpha k\beta ,j\beta l\alpha } \\
<CI|S_-*S_+|CI> &= n_\beta + \delta_{ik}\delta_{jl}Gamma_{i\beta k\alpha ,j\alpha l\beta } \\
<CI|S_z*S_z|CI> &= \delta_{ik}\delta_{jl}(Gamma_{i\alpha k\alpha ,j\alpha l\alpha }
                 - Gamma_{i\alpha k\alpha ,j\beta l\beta }
                 - Gamma_{i\beta k\beta ,j\alpha l\alpha}
                 + Gamma_{i\beta k\beta ,j\beta l\beta})
                 + (n_\alpha+n_\beta)/4

Given the overlap between non-degenerate alpha and beta orbitals, this function can compute the expectation value spin square operator for UHF-FCI wavefunction

Module contents

Different FCI solvers are implemented to support different type of symmetry.

Symmetry

File Point group Spin singlet Real hermitian* Alpha/beta degeneracy direct_spin0_symm Yes Yes Yes Yes direct_spin1_symm Yes No Yes Yes direct_spin0 No Yes Yes Yes direct_spin1 No No Yes Yes direct_uhf No No Yes No direct_nosym No No No** Yes fci_dhf_slow No No No*** No

  • Real hermitian Hamiltonian implies (ij|kl) = (ji|kl) = (ij|lk) = (ji|lk)

** Hamiltonian is real but not hermitian, (ij|kl) != (ji|kl) … *** Hamiltonian is complex hermitian (DHF case) or real hermitian (GHF case)

pyscf.fci.FCI(mol_or_mf, mo=None, singlet=False)[source]

FCI solver

Args:
mol_or_mf :

A Mole object or an SCF object

Kwargs:
mo :

Molecular orbital coefficients

singlet :

Whether to enable spin symmetry for S=0 RHF-based FCI solver.

Returns:

A FCI object

pyscf.fci.solver(mol=None, singlet=False, symm=None)[source]