pyscf.ci package#
Submodules#
pyscf.ci.addons module#
pyscf.ci.cisd module#
Solve CISD equation H C = C e where e = E_HF + E_CORR
- class pyscf.ci.cisd.CISD(mf, frozen=None, mo_coeff=None, mo_occ=None)[source]#
Bases:
StreamObject
restricted CISD
- Attributes:
- verboseint
Print level. Default value equals to
Mole.verbose
- max_memoryfloat or int
Allowed memory in MB. Default value equals to
Mole.max_memory
- conv_tolfloat
converge threshold. Default is 1e-9.
- max_cycleint
max number of iterations. Default is 50.
- max_spaceint
Davidson diagonalization space size. Default is 12.
- directbool
AO-direct CISD. Default is False.
- async_iobool
Allow for asynchronous function execution. Default is True.
- frozenint or list
If integer is given, the inner-most orbitals are frozen from CI amplitudes. Given the orbital indices (0-based) in a list, both occupied and virtual orbitals can be frozen in CI calculation.
>>> mol = gto.M(atom = 'H 0 0 0; F 0 0 1.1', basis = 'ccpvdz') >>> mf = scf.RHF(mol).run() >>> # freeze 2 core orbitals >>> myci = ci.CISD(mf).set(frozen = 2).run() >>> # freeze 2 core orbitals and 3 high lying unoccupied orbitals >>> myci.set(frozen = [0,1,16,17,18]).run()
Saved results
- convergedbool
CISD converged or not
- e_corrfloat
CISD correlation correction
- e_totfloat
Total CCSD energy (HF + correlation)
- ci :
CI wavefunction coefficients
- as_scanner()#
Generating a scanner/solver for CISD PES.
The returned solver is a function. This function requires one argument “mol” as input and returns total CISD energy.
The solver will automatically use the results of last calculation as the initial guess of the new calculation. All parameters assigned in the CISD and the underlying SCF objects (conv_tol, max_memory etc) are automatically applied in the solver.
Note scanner has side effects. It may change many underlying objects (_scf, with_df, with_x2c, …) during calculation.
Examples:
>>> from pyscf import gto, scf, ci >>> mol = gto.M(atom='H 0 0 0; F 0 0 1') >>> ci_scanner = ci.CISD(scf.RHF(mol)).as_scanner() >>> e_tot = ci_scanner(gto.M(atom='H 0 0 0; F 0 0 1.1')) >>> e_tot = ci_scanner(gto.M(atom='H 0 0 0; F 0 0 1.5'))
- async_io = True#
- contract(civec, eris)#
Application of CISD hamiltonian onto civec.
- Args:
myci : CISD (inheriting) object civec : numpy array, same length as a CI vector. eris : ccsd._ChemistsERIs (inheriting) object (poss diff for df)
Contains the various (pq|rs) integrals needed.
- Returns:
numpy array, same length as a CI vector.
- conv_tol = 1e-09#
- direct = False#
- property e_tot#
- get_e_hf(mo_coeff=None)#
- get_frozen_mask()#
Get boolean mask for the restricted reference orbitals.
In the returned boolean (mask) array of frozen orbital indices, the element is False if it corresponds to the frozen orbital.
- get_init_guess(eris=None, nroots=1, diag=None)[source]#
MP2 energy and MP2 initial guess(es) for CISD coefficients.
- Kwargs:
- erisccsd._ChemistsERIs (inheriting) object (poss diff for df)
Contains the various (pq|rs) integrals needed.
- nrootsinteger
Number of CISD solutions to be found.
- diagnumpy array (1D)
e.g. CISD Hamiltonian diagonal in Slater determinant space with HF energy subtracted.
- Returns:
Tuple of float and numpy array or tuple of float and list of numpy arrays (if nroots > 1) MP2 energy and initial guess(es) for CISD coefficients.
- get_nmo()#
- get_nocc()#
- kernel(ci0=None, eris=None)[source]#
Kernel function is the main driver of a method. Every method should define the kernel function as the entry of the calculation. Note the return value of kernel function is not strictly defined. It can be anything related to the method (such as the energy, the wave-function, the DFT mesh grids etc.).
- level_shift = 0.001#
- lindep = 1e-14#
- make_diagonal(eris)#
Return diagonal of CISD hamiltonian in Slater determinant basis.
Note that a constant has been substracted of all elements. The first element is the HF energy (minus the constant), the next elements are the diagonal elements with singly excited determinants (<D_i^a|H|D_i^a> within the constant), then doubly excited determinants (<D_ij^ab|H|D_ij^ab> within the constant).
- Args:
myci : CISD (inheriting) object eris : ccsd._ChemistsERIs (inheriting) object (poss diff for df)
Contains the various (pq|rs) integrals needed.
- Returns:
- numpy array (size: (1, 1 + #single excitations from HF det
#double excitations from HF det))
Diagonal elements of hamiltonian matrix within a constant, see above.
- make_rdm1(civec=None, nmo=None, nocc=None, ao_repr=False)#
Spin-traced one-particle density matrix in MO basis (the occupied-virtual blocks from the orbital response contribution are not included).
dm1[p,q] = <q_alpha^dagger p_alpha> + <q_beta^dagger p_beta>
The convention of 1-pdm 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=None, nmo=None, nocc=None, ao_repr=False)#
Spin-traced two-particle density matrix in MO basis
dm2[p,q,r,s] = sum_{sigma,tau} <p_sigma^dagger r_tau^dagger s_tau q_sigma>
Note the contraction between ERIs (in Chemist’s notation) and rdm2 is E = einsum(‘pqrs,pqrs’, eri, rdm2)
- max_cycle = 50#
- max_space = 12#
- property nmo#
- property nocc#
- property nstates#
- to_gpu(out=None)#
Convert a method to its corresponding GPU variant, and recursively converts all attributes of a method to cupy objects or gpu4pyscf objects.
- trans_rdm1(cibra, ciket, nmo=None, nocc=None)#
Spin-traced one-particle transition density matrix in MO basis.
dm1[p,q] = <q_alpha^dagger p_alpha> + <q_beta^dagger p_beta>
The convention of 1-pdm 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)
- vector_size()[source]#
The size of the vector which was returned from
amplitudes_to_cisdvec()
- class pyscf.ci.cisd.CISD_Scanner(ci)[source]#
Bases:
SinglePointScanner
- pyscf.ci.cisd.as_scanner(ci)[source]#
Generating a scanner/solver for CISD PES.
The returned solver is a function. This function requires one argument “mol” as input and returns total CISD energy.
The solver will automatically use the results of last calculation as the initial guess of the new calculation. All parameters assigned in the CISD and the underlying SCF objects (conv_tol, max_memory etc) are automatically applied in the solver.
Note scanner has side effects. It may change many underlying objects (_scf, with_df, with_x2c, …) during calculation.
Examples:
>>> from pyscf import gto, scf, ci >>> mol = gto.M(atom='H 0 0 0; F 0 0 1') >>> ci_scanner = ci.CISD(scf.RHF(mol)).as_scanner() >>> e_tot = ci_scanner(gto.M(atom='H 0 0 0; F 0 0 1.1')) >>> e_tot = ci_scanner(gto.M(atom='H 0 0 0; F 0 0 1.5'))
- pyscf.ci.cisd.contract(myci, civec, eris)[source]#
Application of CISD hamiltonian onto civec.
- Args:
myci : CISD (inheriting) object civec : numpy array, same length as a CI vector. eris : ccsd._ChemistsERIs (inheriting) object (poss diff for df)
Contains the various (pq|rs) integrals needed.
- Returns:
numpy array, same length as a CI vector.
- pyscf.ci.cisd.from_fcivec(ci0, norb, nelec, frozen=None)[source]#
Extract CISD coefficients from FCI coefficients
- pyscf.ci.cisd.kernel(myci, eris, ci0=None, max_cycle=50, tol=1e-08, verbose=4)[source]#
Run CISD calculation.
- Args:
myci : CISD (inheriting) object eris : ccsd._ChemistsERIs (inheriting) object (poss diff for df)
Contains the various (pq|rs) integrals needed.
- Kwargs:
- ci0(List of) numpy array(s) (if None it will set)
Initial guess for CISD coeffs.
- max_cycleinteger
Maximum number of iterations to converge to CISD solution. If not converged before, calculation stops without having converged.
- tolfloat
Convergence tolerance.
- verboseinteger
Level of output (roughly: the higher, the more output).
- Returns:
- convbool
Is it converged?
- ecisdList of floats or float
The lowest
myci.nroots
eigenvalues.- ciList of 1D arrays or 1D array
The lowest
myci.nroots
eigenvectors.
- pyscf.ci.cisd.make_diagonal(myci, eris)[source]#
Return diagonal of CISD hamiltonian in Slater determinant basis.
Note that a constant has been substracted of all elements. The first element is the HF energy (minus the constant), the next elements are the diagonal elements with singly excited determinants (<D_i^a|H|D_i^a> within the constant), then doubly excited determinants (<D_ij^ab|H|D_ij^ab> within the constant).
- Args:
myci : CISD (inheriting) object eris : ccsd._ChemistsERIs (inheriting) object (poss diff for df)
Contains the various (pq|rs) integrals needed.
- Returns:
- numpy array (size: (1, 1 + #single excitations from HF det
#double excitations from HF det))
Diagonal elements of hamiltonian matrix within a constant, see above.
- pyscf.ci.cisd.make_rdm1(myci, civec=None, nmo=None, nocc=None, ao_repr=False)[source]#
Spin-traced one-particle density matrix in MO basis (the occupied-virtual blocks from the orbital response contribution are not included).
dm1[p,q] = <q_alpha^dagger p_alpha> + <q_beta^dagger p_beta>
The convention of 1-pdm 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.ci.cisd.make_rdm2(myci, civec=None, nmo=None, nocc=None, ao_repr=False)[source]#
Spin-traced two-particle density matrix in MO basis
dm2[p,q,r,s] = sum_{sigma,tau} <p_sigma^dagger r_tau^dagger s_tau q_sigma>
Note the contraction between ERIs (in Chemist’s notation) and rdm2 is E = einsum(‘pqrs,pqrs’, eri, rdm2)
- pyscf.ci.cisd.overlap(cibra, ciket, nmo, nocc, s=None)[source]#
Overlap between two CISD wavefunctions.
- Args:
- s2D array
The overlap matrix of non-orthogonal one-particle basis
- pyscf.ci.cisd.t1strs(norb, nelec)[source]#
Compute the FCI strings (address) for CIS single-excitation amplitudes and the signs of the coefficients when transferring the reference from physics vacuum to HF vacuum.
- pyscf.ci.cisd.tn_addrs_signs(norb, nelec, n_excite)[source]#
Compute the FCI strings (address) for CIS n-excitation amplitudes and the signs of the coefficients when transferring the reference from physics vacuum to HF vacuum.
If the excitation level is not compatible with the number of electrons and holes, empty lists are returned for the addresses and signs.
- pyscf.ci.cisd.to_fcivec(cisdvec, norb, nelec, frozen=None)[source]#
Convert CISD coefficients to FCI coefficients
- pyscf.ci.cisd.trans_rdm1(myci, cibra, ciket, nmo=None, nocc=None)[source]#
Spin-traced one-particle transition density matrix in MO basis.
dm1[p,q] = <q_alpha^dagger p_alpha> + <q_beta^dagger p_beta>
The convention of 1-pdm 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.ci.gcisd module#
General spin-orbital CISD
- class pyscf.ci.gcisd.GCISD(mf, frozen=None, mo_coeff=None, mo_occ=None)[source]#
Bases:
CISD
- contract(civec, eris)#
Application of CISD hamiltonian onto civec.
- Args:
myci : CISD (inheriting) object civec : numpy array, same length as a CI vector. eris : ccsd._ChemistsERIs (inheriting) object (poss diff for df)
Contains the various (pq|rs) integrals needed.
- Returns:
numpy array, same length as a CI vector.
- from_rcisdvec(civec, nocc=None, orbspin=None)#
Convert the (spin-separated) CISD coefficient vector to GCISD coefficient vector
- from_ucisdvec(civec, nocc=None, orbspin=None)[source]#
Convert the (spin-separated) CISD coefficient vector to GCISD coefficient vector
- get_init_guess(eris=None, nroots=1, diag=None)[source]#
MP2 energy and MP2 initial guess(es) for CISD coefficients.
- Kwargs:
- erisccsd._ChemistsERIs (inheriting) object (poss diff for df)
Contains the various (pq|rs) integrals needed.
- nrootsinteger
Number of CISD solutions to be found.
- diagnumpy array (1D)
e.g. CISD Hamiltonian diagonal in Slater determinant space with HF energy subtracted.
- Returns:
Tuple of float and numpy array or tuple of float and list of numpy arrays (if nroots > 1) MP2 energy and initial guess(es) for CISD coefficients.
- make_diagonal(eris)#
Return diagonal of CISD hamiltonian in Slater determinant basis.
Note that a constant has been substracted of all elements. The first element is the HF energy (minus the constant), the next elements are the diagonal elements with singly excited determinants (<D_i^a|H|D_i^a> within the constant), then doubly excited determinants (<D_ij^ab|H|D_ij^ab> within the constant).
- Args:
myci : CISD (inheriting) object eris : ccsd._ChemistsERIs (inheriting) object (poss diff for df)
Contains the various (pq|rs) integrals needed.
- Returns:
- numpy array (size: (1, 1 + #single excitations from HF det
#double excitations from HF det))
Diagonal elements of hamiltonian matrix within a constant, see above.
- make_rdm1(civec=None, nmo=None, nocc=None, ao_repr=False)#
One-particle density matrix in the molecular spin-orbital representation (the occupied-virtual blocks from the orbital response contribution are not included).
dm1[p,q] = <q^dagger p> (p,q are spin-orbitals)
The convention of 1-pdm 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=None, nmo=None, nocc=None, ao_repr=False)#
Two-particle density matrix in the molecular spin-orbital representation
dm2[p,q,r,s] = <p^dagger r^dagger s q>
where p,q,r,s are spin-orbitals. p,q correspond to one particle and r,s correspond to another particle. The contraction between ERIs (in Chemist’s notation) and rdm2 is E = einsum(‘pqrs,pqrs’, eri, rdm2)
- to_ucisdvec(civec, orbspin=None)[source]#
Convert the GCISD coefficient vector to UCISD coefficient vector
- trans_rdm1(cibra, ciket, nmo=None, nocc=None)#
One-particle transition density matrix in the molecular spin-orbital representation.
dm1[p,q] = <q^dagger p> (p,q are spin-orbitals)
The convention of 1-pdm 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)
- vector_size()[source]#
The size of the vector which was returned from
amplitudes_to_cisdvec()
- pyscf.ci.gcisd.from_rcisdvec(civec, nocc, orbspin)#
Convert the (spin-separated) CISD coefficient vector to GCISD coefficient vector
- pyscf.ci.gcisd.from_ucisdvec(civec, nocc, orbspin)[source]#
Convert the (spin-separated) CISD coefficient vector to GCISD coefficient vector
- pyscf.ci.gcisd.make_rdm1(myci, civec=None, nmo=None, nocc=None, ao_repr=False)[source]#
One-particle density matrix in the molecular spin-orbital representation (the occupied-virtual blocks from the orbital response contribution are not included).
dm1[p,q] = <q^dagger p> (p,q are spin-orbitals)
The convention of 1-pdm 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.ci.gcisd.make_rdm2(myci, civec=None, nmo=None, nocc=None, ao_repr=False)[source]#
Two-particle density matrix in the molecular spin-orbital representation
dm2[p,q,r,s] = <p^dagger r^dagger s q>
where p,q,r,s are spin-orbitals. p,q correspond to one particle and r,s correspond to another particle. The contraction between ERIs (in Chemist’s notation) and rdm2 is E = einsum(‘pqrs,pqrs’, eri, rdm2)
- pyscf.ci.gcisd.to_ucisdvec(civec, nmo, nocc, orbspin)[source]#
Convert the GCISD coefficient vector to UCISD coefficient vector
- pyscf.ci.gcisd.trans_rdm1(myci, cibra, ciket, nmo=None, nocc=None)[source]#
One-particle transition density matrix in the molecular spin-orbital representation.
dm1[p,q] = <q^dagger p> (p,q are spin-orbitals)
The convention of 1-pdm 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.ci.ucisd module#
Unrestricted CISD
- class pyscf.ci.ucisd.UCISD(mf, frozen=None, mo_coeff=None, mo_occ=None)[source]#
Bases:
CISD
- contract(civec, eris)#
Application of CISD hamiltonian onto civec.
- Args:
myci : CISD (inheriting) object civec : numpy array, same length as a CI vector. eris : ccsd._ChemistsERIs (inheriting) object (poss diff for df)
Contains the various (pq|rs) integrals needed.
- Returns:
numpy array, same length as a CI vector.
- get_frozen_mask()#
Get boolean mask for the unrestricted reference orbitals.
In the returned boolean (mask) array of frozen orbital indices, the element is False if it corresponds to the frozen orbital.
- get_init_guess(eris=None, nroots=1, diag=None)[source]#
MP2 energy and MP2 initial guess(es) for CISD coefficients.
- Kwargs:
- erisccsd._ChemistsERIs (inheriting) object (poss diff for df)
Contains the various (pq|rs) integrals needed.
- nrootsinteger
Number of CISD solutions to be found.
- diagnumpy array (1D)
e.g. CISD Hamiltonian diagonal in Slater determinant space with HF energy subtracted.
- Returns:
Tuple of float and numpy array or tuple of float and list of numpy arrays (if nroots > 1) MP2 energy and initial guess(es) for CISD coefficients.
- get_nmo()#
- get_nocc()#
- make_diagonal(eris)#
Return diagonal of CISD hamiltonian in Slater determinant basis.
Note that a constant has been substracted of all elements. The first element is the HF energy (minus the constant), the next elements are the diagonal elements with singly excited determinants (<D_i^a|H|D_i^a> within the constant), then doubly excited determinants (<D_ij^ab|H|D_ij^ab> within the constant).
- Args:
myci : CISD (inheriting) object eris : ccsd._ChemistsERIs (inheriting) object (poss diff for df)
Contains the various (pq|rs) integrals needed.
- Returns:
- numpy array (size: (1, 1 + #single excitations from HF det
#double excitations from HF det))
Diagonal elements of hamiltonian matrix within a constant, see above.
- make_rdm1(civec=None, nmo=None, nocc=None, ao_repr=False)#
One-particle spin density matrices dm1a, dm1b in MO basis (the occupied-virtual blocks due to the orbital response contribution are not included).
dm1a[p,q] = <q_alpha^dagger p_alpha> dm1b[p,q] = <q_beta^dagger p_beta>
The convention of 1-pdm is based on McWeeney’s book, Eq (5.4.20).
- make_rdm2(civec=None, nmo=None, nocc=None, ao_repr=False)#
Two-particle spin density matrices dm2aa, dm2ab, dm2bb in MO basis
dm2aa[p,q,r,s] = <q_alpha^dagger s_alpha^dagger r_alpha p_alpha> dm2ab[p,q,r,s] = <q_alpha^dagger s_beta^dagger r_beta p_alpha> dm2bb[p,q,r,s] = <q_beta^dagger s_beta^dagger r_beta p_beta>
(p,q correspond to one particle and r,s correspond to another particle) Two-particle density matrix should be contracted to integrals with the pattern below to compute energy
E = numpy.einsum(‘pqrs,pqrs’, eri_aa, dm2_aa) E+= numpy.einsum(‘pqrs,pqrs’, eri_ab, dm2_ab) E+= numpy.einsum(‘pqrs,rspq’, eri_ba, dm2_ab) E+= numpy.einsum(‘pqrs,pqrs’, eri_bb, dm2_bb)
where eri_aa[p,q,r,s] = (p_alpha q_alpha | r_alpha s_alpha ) eri_ab[p,q,r,s] = ( p_alpha q_alpha | r_beta s_beta ) eri_ba[p,q,r,s] = ( p_beta q_beta | r_alpha s_alpha ) eri_bb[p,q,r,s] = ( p_beta q_beta | r_beta s_beta )
- trans_rdm1(cibra, ciket, nmo=None, nocc=None)#
One-particle spin density matrices dm1a, dm1b in MO basis (the occupied-virtual blocks due to the orbital response contribution are not included).
dm1a[p,q] = <q_alpha^dagger p_alpha> dm1b[p,q] = <q_beta^dagger p_beta>
The convention of 1-pdm is based on McWeeney’s book, Eq (5.4.20).
- vector_size()[source]#
The size of the vector which was returned from
amplitudes_to_cisdvec()
- pyscf.ci.ucisd.from_fcivec(ci0, norb, nelec, frozen=None)[source]#
Extract CISD coefficients from FCI coefficients
- pyscf.ci.ucisd.make_rdm1(myci, civec=None, nmo=None, nocc=None, ao_repr=False)[source]#
One-particle spin density matrices dm1a, dm1b in MO basis (the occupied-virtual blocks due to the orbital response contribution are not included).
dm1a[p,q] = <q_alpha^dagger p_alpha> dm1b[p,q] = <q_beta^dagger p_beta>
The convention of 1-pdm is based on McWeeney’s book, Eq (5.4.20).
- pyscf.ci.ucisd.make_rdm2(myci, civec=None, nmo=None, nocc=None, ao_repr=False)[source]#
Two-particle spin density matrices dm2aa, dm2ab, dm2bb in MO basis
dm2aa[p,q,r,s] = <q_alpha^dagger s_alpha^dagger r_alpha p_alpha> dm2ab[p,q,r,s] = <q_alpha^dagger s_beta^dagger r_beta p_alpha> dm2bb[p,q,r,s] = <q_beta^dagger s_beta^dagger r_beta p_beta>
(p,q correspond to one particle and r,s correspond to another particle) Two-particle density matrix should be contracted to integrals with the pattern below to compute energy
E = numpy.einsum(‘pqrs,pqrs’, eri_aa, dm2_aa) E+= numpy.einsum(‘pqrs,pqrs’, eri_ab, dm2_ab) E+= numpy.einsum(‘pqrs,rspq’, eri_ba, dm2_ab) E+= numpy.einsum(‘pqrs,pqrs’, eri_bb, dm2_bb)
where eri_aa[p,q,r,s] = (p_alpha q_alpha | r_alpha s_alpha ) eri_ab[p,q,r,s] = ( p_alpha q_alpha | r_beta s_beta ) eri_ba[p,q,r,s] = ( p_beta q_beta | r_alpha s_alpha ) eri_bb[p,q,r,s] = ( p_beta q_beta | r_beta s_beta )
- pyscf.ci.ucisd.overlap(cibra, ciket, nmo, nocc, s=None)[source]#
Overlap between two CISD wavefunctions.
- Args:
- sa list of 2D arrays
The overlap matrix of non-orthogonal one-particle basis
- pyscf.ci.ucisd.to_fcivec(cisdvec, norb, nelec, frozen=None)[source]#
Convert CISD coefficients to FCI coefficients
- pyscf.ci.ucisd.trans_rdm1(myci, cibra, ciket, nmo=None, nocc=None)[source]#
One-particle spin density matrices dm1a, dm1b in MO basis (the occupied-virtual blocks due to the orbital response contribution are not included).
dm1a[p,q] = <q_alpha^dagger p_alpha> dm1b[p,q] = <q_beta^dagger p_beta>
The convention of 1-pdm is based on McWeeney’s book, Eq (5.4.20).
Module contents#
- pyscf.ci.CISD(mf, frozen=None, mo_coeff=None, mo_occ=None)[source]#
restricted CISD
- Attributes:
- verboseint
Print level. Default value equals to
Mole.verbose
- max_memoryfloat or int
Allowed memory in MB. Default value equals to
Mole.max_memory
- conv_tolfloat
converge threshold. Default is 1e-9.
- max_cycleint
max number of iterations. Default is 50.
- max_spaceint
Davidson diagonalization space size. Default is 12.
- directbool
AO-direct CISD. Default is False.
- async_iobool
Allow for asynchronous function execution. Default is True.
- frozenint or list
If integer is given, the inner-most orbitals are frozen from CI amplitudes. Given the orbital indices (0-based) in a list, both occupied and virtual orbitals can be frozen in CI calculation.
>>> mol = gto.M(atom = 'H 0 0 0; F 0 0 1.1', basis = 'ccpvdz') >>> mf = scf.RHF(mol).run() >>> # freeze 2 core orbitals >>> myci = ci.CISD(mf).set(frozen = 2).run() >>> # freeze 2 core orbitals and 3 high lying unoccupied orbitals >>> myci.set(frozen = [0,1,16,17,18]).run()
Saved results
- convergedbool
CISD converged or not
- e_corrfloat
CISD correlation correction
- e_totfloat
Total CCSD energy (HF + correlation)
- ci :
CI wavefunction coefficients