pyscf.cc package#
Submodules#
pyscf.cc.addons module#
- pyscf.cc.addons.spatial2spin(tx, orbspin=None)[source]#
Convert T1/T2 of spatial orbital representation to T1/T2 of spin-orbital representation
- pyscf.cc.addons.spatial2spinorb(tx, orbspin=None)#
Convert T1/T2 of spatial orbital representation to T1/T2 of spin-orbital representation
pyscf.cc.bccd module#
Brueckner coupled-cluster doubles (BCCD).
- pyscf.cc.bccd.bccd_kernel_(mycc, u=None, conv_tol_normu=1e-05, max_cycle=20, diis=True, canonicalization=True, verbose=4)[source]#
Brueckner coupled-cluster wrapper, using an outer-loop algorithm.
- Args:
mycc: a converged CCSD object. u: initial transformation matrix. conv_tol_normu: convergence tolerance for u matrix. max_cycle: Maximum number of BCC cycles. diis: whether perform DIIS. canonicalization: whether to semi-canonicalize the Brueckner orbitals. verbose: verbose for CCSD inner iterations.
- Returns:
- mycc: a modified CC object with t1 vanished.
mycc._scf and mycc will be modified.
- pyscf.cc.bccd.transform_l1_to_bo(t1, umat)#
Transform t1 to brueckner orbital basis.
- pyscf.cc.bccd.transform_l2_to_bo(t2, umat, umat_b=None)#
Transform t2 to brueckner orbital basis.
pyscf.cc.ccd module#
Coupled cluster doubles
- class pyscf.cc.ccd.CCD(mf, frozen=None, mo_coeff=None, mo_occ=None)[source]#
Bases:
CCSD
- kernel(t2=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.).
- make_rdm1(t2=None, l2=None, ao_repr=False)[source]#
Un-relaxed 1-particle density matrix in MO space
pyscf.cc.ccsd module#
RCCSD for real integrals 8-fold permutation symmetry has been used (ij|kl) = (ji|kl) = (kl|ij) = …
- class pyscf.cc.ccsd.CCSD(mf, frozen=None, mo_coeff=None, mo_occ=None)[source]#
Bases:
CCSDBase
restricted CCSD
- 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-7.
- conv_tol_normtfloat
converge threshold for norm(t1,t2). Default is 1e-5.
- max_cycleint
max number of iterations. Default is 50.
- diis_spaceint
DIIS space size. Default is 6.
- diis_start_cycleint
The step to start DIIS. Default is 0.
- iterative_dampingfloat
The self consistent damping parameter.
- directbool
AO-direct CCSD. Default is False.
- async_iobool
Allow for asynchronous function execution. Default is True.
- incore_completebool
Avoid all I/O (also for DIIS). Default is False.
- level_shiftfloat
A shift on virtual orbital energies to stabilize the CCSD iteration
- frozenint or list
If integer is given, the inner-most orbitals are frozen from CC amplitudes. Given the orbital indices (0-based) in a list, both occupied and virtual orbitals can be frozen in CC 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 >>> mycc = cc.CCSD(mf).set(frozen = 2).run() >>> # auto-generate the number of core orbitals to be frozen (1 in this case) >>> mycc = cc.CCSD(mf).set_frozen().run() >>> # freeze 2 core orbitals and 3 high lying unoccupied orbitals >>> mycc.set(frozen = [0,1,16,17,18]).run()
- callbackfunction(envs_dict) => None
callback function takes one dict as the argument which is generated by the builtin function
locals()
, so that the callback function can access all local variables in the current environment.
Saved results:
- convergedbool
Whether the CCSD iteration converged
- e_corrfloat
CCSD correlation correction
- e_totfloat
Total CCSD energy (HF + correlation)
- t1, t2 :
T amplitudes t1[i,a], t2[i,j,a,b] (i,j in occ, a,b in virt)
- l1, l2 :
Lambda amplitudes l1[i,a], l2[i,j,a,b] (i,j in occ, a,b in virt)
- cyclesint
The number of iteration cycles performed
- EOMEA(*args, **kwargs)#
- EOMEA_Ta(*args, **kwargs)#
Class for EOM EACCSD(T)*(a) method by Matthews and Stanton.
- EOMEE(*args, **kwargs)#
- EOMEESinglet(*args, **kwargs)#
- EOMEESpinFlip(*args, **kwargs)#
- EOMEETriplet(*args, **kwargs)#
- EOMIP(*args, **kwargs)#
- EOMIP_Ta(*args, **kwargs)#
Class for EOM IPCCSD(T)*(a) method by Matthews and Stanton.
- make_rdm1(t1=None, t2=None, l1=None, l2=None, ao_repr=False, with_frozen=True, with_mf=True)[source]#
Un-relaxed 1-particle density matrix in MO space
- class pyscf.cc.ccsd.CCSDBase(mf, frozen=None, mo_coeff=None, mo_occ=None)[source]#
Bases:
StreamObject
restricted CCSD
- 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-7.
- conv_tol_normtfloat
converge threshold for norm(t1,t2). Default is 1e-5.
- max_cycleint
max number of iterations. Default is 50.
- diis_spaceint
DIIS space size. Default is 6.
- diis_start_cycleint
The step to start DIIS. Default is 0.
- iterative_dampingfloat
The self consistent damping parameter.
- directbool
AO-direct CCSD. Default is False.
- async_iobool
Allow for asynchronous function execution. Default is True.
- incore_completebool
Avoid all I/O (also for DIIS). Default is False.
- level_shiftfloat
A shift on virtual orbital energies to stabilize the CCSD iteration
- frozenint or list
If integer is given, the inner-most orbitals are frozen from CC amplitudes. Given the orbital indices (0-based) in a list, both occupied and virtual orbitals can be frozen in CC 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 >>> mycc = cc.CCSD(mf).set(frozen = 2).run() >>> # auto-generate the number of core orbitals to be frozen (1 in this case) >>> mycc = cc.CCSD(mf).set_frozen().run() >>> # freeze 2 core orbitals and 3 high lying unoccupied orbitals >>> mycc.set(frozen = [0,1,16,17,18]).run()
- callbackfunction(envs_dict) => None
callback function takes one dict as the argument which is generated by the builtin function
locals()
, so that the callback function can access all local variables in the current environment.
Saved results:
- convergedbool
Whether the CCSD iteration converged
- e_corrfloat
CCSD correlation correction
- e_totfloat
Total CCSD energy (HF + correlation)
- t1, t2 :
T amplitudes t1[i,a], t2[i,j,a,b] (i,j in occ, a,b in virt)
- l1, l2 :
Lambda amplitudes l1[i,a], l2[i,j,a,b] (i,j in occ, a,b in virt)
- cyclesint
The number of iteration cycles performed
- as_scanner()#
Generating a scanner/solver for CCSD PES.
The returned solver is a function. This function requires one argument “mol” as input and returns total CCSD energy.
The solver will automatically use the results of last calculation as the initial guess of the new calculation. All parameters assigned in the CCSD 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, cc >>> mol = gto.M(atom='H 0 0 0; F 0 0 1') >>> cc_scanner = cc.CCSD(scf.RHF(mol)).as_scanner() >>> e_tot = cc_scanner(gto.M(atom='H 0 0 0; F 0 0 1.1')) >>> e_tot = cc_scanner(gto.M(atom='H 0 0 0; F 0 0 1.5'))
- async_io = True#
- callback = None#
- cc2 = False#
- conv_tol = 1e-07#
- conv_tol_normt = 1e-05#
- diis = True#
- diis_file = None#
- diis_space = 6#
- diis_start_cycle = 0#
- diis_start_energy_diff = 1000000000.0#
- direct = False#
- property e_tot#
- property ecc#
- energy(t1=None, t2=None, eris=None)#
CCSD correlation energy
- 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_nmo()#
- get_nocc()#
- incore_complete = False#
- iterative_damping = 1.0#
- kernel(t1=None, t2=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.).
- make_rdm1(t1=None, t2=None, l1=None, l2=None, ao_repr=False, with_frozen=True, with_mf=True)[source]#
Un-relaxed 1-particle density matrix in MO space
- make_rdm2(t1=None, t2=None, l1=None, l2=None, ao_repr=False, with_frozen=True, with_dm1=True)[source]#
2-particle density matrix in MO space. The density matrix is stored as
dm2[p,r,q,s] = <p^+ q^+ s r>
- max_cycle = 50#
- property nmo#
- property nocc#
- restore_from_diis_(diis_file, inplace=True)#
Reuse an existed DIIS object in the CCSD calculation.
The CCSD amplitudes will be restored from the DIIS object to generate t1 and t2 amplitudes. The t1/t2 amplitudes of the CCSD object will be overwritten by the generated t1 and t2 amplitudes. The amplitudes vector and error vector will be reused in the CCSD calculation.
- update_amps(t1, t2, eris)#
- class pyscf.cc.ccsd.CCSD_Scanner(cc)[source]#
Bases:
SinglePointScanner
- pyscf.cc.ccsd.as_scanner(cc)[source]#
Generating a scanner/solver for CCSD PES.
The returned solver is a function. This function requires one argument “mol” as input and returns total CCSD energy.
The solver will automatically use the results of last calculation as the initial guess of the new calculation. All parameters assigned in the CCSD 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, cc >>> mol = gto.M(atom='H 0 0 0; F 0 0 1') >>> cc_scanner = cc.CCSD(scf.RHF(mol)).as_scanner() >>> e_tot = cc_scanner(gto.M(atom='H 0 0 0; F 0 0 1.1')) >>> e_tot = cc_scanner(gto.M(atom='H 0 0 0; F 0 0 1.5'))
- pyscf.cc.ccsd.get_d1_diagnostic(t1)[source]#
D1 diagnostic given in
Janssen, et. al Chem. Phys. Lett. 290 (1998) 423
- pyscf.cc.ccsd.get_d2_diagnostic(t2)[source]#
D2 diagnostic given in
Nielsen, et. al Chem. Phys. Lett. 310 (1999) 568
Note: This is currently only defined in the literature for restricted closed-shell systems.
- pyscf.cc.ccsd.get_t1_diagnostic(t1)[source]#
Returns the t1 amplitude norm, normalized by number of correlated electrons.
- pyscf.cc.ccsd.kernel(mycc, eris=None, t1=None, t2=None, max_cycle=50, tol=1e-08, tolnormt=1e-06, verbose=None, callback=None)[source]#
- pyscf.cc.ccsd.restore_from_diis_(mycc, diis_file, inplace=True)[source]#
Reuse an existed DIIS object in the CCSD calculation.
The CCSD amplitudes will be restored from the DIIS object to generate t1 and t2 amplitudes. The t1/t2 amplitudes of the CCSD object will be overwritten by the generated t1 and t2 amplitudes. The amplitudes vector and error vector will be reused in the CCSD calculation.
pyscf.cc.ccsd_lambda module#
Restricted CCSD implementation for real integrals. Permutation symmetry for the 4-index integrals (ij|kl) = (ij|lk) = (ji|kl) are assumed.
Note MO integrals are treated in chemist’s notation
pyscf.cc.ccsd_rdm module#
- pyscf.cc.ccsd_rdm.make_rdm1(mycc, t1, t2, l1, l2, ao_repr=False, with_frozen=True, with_mf=True)[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.cc.ccsd_rdm.make_rdm2(mycc, t1, t2, l1, l2, ao_repr=False, with_frozen=True, with_dm1=True)[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.cc.ccsd_rdm_slow module#
pyscf.cc.ccsd_t module#
RHF-CCSD(T) for real integrals
pyscf.cc.ccsd_t_lambda_slow module#
Spin-free lambda equation of RHF-CCSD(T)
Ref: JCP 147, 044104 (2017); DOI:10.1063/1.4994918
pyscf.cc.ccsd_t_rdm_slow module#
pyscf.cc.ccsd_t_slow module#
pyscf.cc.dfccsd module#
pyscf.cc.eom_gccsd module#
- class pyscf.cc.eom_gccsd.EOMEA(cc)[source]#
Bases:
EOMEA
- static amplitudes_to_vector(r1, r2)#
- ccsd_star_contract(eaccsd_evals, eaccsd_evecs, leaccsd_evecs, imds=None)#
- Returns:
- e_star (list of float):
The EA-CCSD* energy.
- Notes:
See ipccsd_star_contract for description of arguments.
- Reference:
Saeh, Stanton “…energy surfaces of radicals” JCP 111, 8275 (1999); DOI:10.1063/1.480171
- get_diag(imds=None)#
- l_matvec(vector, imds=None, diag=None)#
EA-CCSD left eigenvector equation.
For description of args, see eaccsd_matvec.
- matvec(vector, imds=None, diag=None)#
- vector_to_amplitudes(vector=None, nmo=None, nocc=None)#
- class pyscf.cc.eom_gccsd.EOMEA_Ta(cc)[source]#
Bases:
EOMEA
Class for EOM EACCSD(T)*(a) method by Matthews and Stanton.
- class pyscf.cc.eom_gccsd.EOMEE(cc)[source]#
Bases:
EOMEE
- static amplitudes_to_vector(t1, t2, out=None)#
- eeccsd(nroots=1, koopmans=False, guess=None, eris=None, imds=None)#
Calculate N-electron neutral excitations via EOM-EE-CCSD.
- Kwargs:
- nrootsint
Number of roots (eigenvalues) requested
- koopmansbool
Calculate Koopmans’-like (1p1h) excitations only, targeting via overlap.
- guesslist of ndarray
List of guess vectors to use for targeting via overlap.
- get_diag(imds=None)#
- kernel(nroots=1, koopmans=False, guess=None, eris=None, imds=None)#
Calculate N-electron neutral excitations via EOM-EE-CCSD.
- Kwargs:
- nrootsint
Number of roots (eigenvalues) requested
- koopmansbool
Calculate Koopmans’-like (1p1h) excitations only, targeting via overlap.
- guesslist of ndarray
List of guess vectors to use for targeting via overlap.
- matvec(vector, imds=None, diag=None)#
- vector_to_amplitudes(vector=None, nmo=None, nocc=None)#
- class pyscf.cc.eom_gccsd.EOMIP(cc)[source]#
Bases:
EOMIP
- static amplitudes_to_vector(r1, r2)#
- ccsd_star_contract(ipccsd_evals, ipccsd_evecs, lipccsd_evecs, imds=None)#
- Returns:
- e_star (list of float):
The IP-CCSD* energy.
- Notes:
The user should check to make sure the right and left eigenvalues before running the perturbative correction.
The 2hp right amplitudes are assumed to be of the form s^{a }_{ij}, i.e. the (ia) indices are coupled while the left are assumed to be of the form s^{ b}_{ij}, i.e. the (jb) indices are coupled.
- Reference:
Saeh, Stanton “…energy surfaces of radicals” JCP 111, 8275 (1999); DOI:10.1063/1.480171
- get_diag(imds=None)#
- l_matvec(vector, imds=None, diag=None)#
IP-CCSD left eigenvector equation.
For description of args, see ipccsd_matvec.
- matvec(vector, imds=None, diag=None)#
- vector_to_amplitudes(vector=None, nmo=None, nocc=None)#
- class pyscf.cc.eom_gccsd.EOMIP_Ta(cc)[source]#
Bases:
EOMIP
Class for EOM IPCCSD(T)*(a) method by Matthews and Stanton.
- pyscf.cc.eom_gccsd.eaccsd_star_contract(eom, eaccsd_evals, eaccsd_evecs, leaccsd_evecs, imds=None)[source]#
- Returns:
- e_star (list of float):
The EA-CCSD* energy.
- Notes:
See ipccsd_star_contract for description of arguments.
- Reference:
Saeh, Stanton “…energy surfaces of radicals” JCP 111, 8275 (1999); DOI:10.1063/1.480171
- pyscf.cc.eom_gccsd.eeccsd(eom, nroots=1, koopmans=False, guess=None, eris=None, imds=None)[source]#
Calculate N-electron neutral excitations via EOM-EE-CCSD.
- Kwargs:
- nrootsint
Number of roots (eigenvalues) requested
- koopmansbool
Calculate Koopmans’-like (1p1h) excitations only, targeting via overlap.
- guesslist of ndarray
List of guess vectors to use for targeting via overlap.
- pyscf.cc.eom_gccsd.ipccsd_star_contract(eom, ipccsd_evals, ipccsd_evecs, lipccsd_evecs, imds=None)[source]#
- Returns:
- e_star (list of float):
The IP-CCSD* energy.
- Notes:
The user should check to make sure the right and left eigenvalues before running the perturbative correction.
The 2hp right amplitudes are assumed to be of the form s^{a }_{ij}, i.e. the (ia) indices are coupled while the left are assumed to be of the form s^{ b}_{ij}, i.e. the (jb) indices are coupled.
- Reference:
Saeh, Stanton “…energy surfaces of radicals” JCP 111, 8275 (1999); DOI:10.1063/1.480171
- pyscf.cc.eom_gccsd.leaccsd_matvec(eom, vector, imds=None, diag=None)[source]#
EA-CCSD left eigenvector equation.
For description of args, see eaccsd_matvec.
pyscf.cc.eom_rccsd module#
- class pyscf.cc.eom_rccsd.EOM(cc)[source]#
Bases:
StreamObject
- class pyscf.cc.eom_rccsd.EOMEA(cc)[source]#
Bases:
EOM
- static amplitudes_to_vector(r1, r2)#
- ccsd_star_contract(eaccsd_evals, eaccsd_evecs, leaccsd_evecs, imds=None)#
- eaccsd(nroots=1, left=False, koopmans=False, guess=None, partition=None, eris=None, imds=None)#
Calculate (N+1)-electron charged excitations via EA-EOM-CCSD.
- Args:
See also ipccd()
- eaccsd_star(nroots=1, koopmans=False, right_guess=None, left_guess=None, eris=None, imds=None, **kwargs)#
Calculates CCSD* perturbative correction.
- Args:
See also ipccd_star()
- property eea#
- get_diag(imds=None)#
- kernel(nroots=1, left=False, koopmans=False, guess=None, partition=None, eris=None, imds=None)#
Calculate (N+1)-electron charged excitations via EA-EOM-CCSD.
- Args:
See also ipccd()
- l_matvec(vector, imds=None, diag=None)#
- matvec(vector, imds=None, diag=None)#
- static spatial2spin(rx, orbspin=None)#
- vector_to_amplitudes(vector=None, nmo=None, nocc=None)#
- class pyscf.cc.eom_rccsd.EOMEA_Ta(cc)[source]#
Bases:
EOMEA
Class for EOM EACCSD(T)*(a) method by Matthews and Stanton.
- class pyscf.cc.eom_rccsd.EOMEE(cc)[source]#
Bases:
EOM
- eeccsd(nroots=1, koopmans=False, guess=None, eris=None, imds=None)#
Calculate N-electron neutral excitations via EOM-EE-CCSD.
- Kwargs:
- nrootsint
Number of roots (eigenvalues) requested
- koopmansbool
Calculate Koopmans’-like (1p1h) excitations only, targeting via overlap.
- guesslist of ndarray
List of guess vectors to use for targeting via overlap.
- property eee#
- get_diag(imds=None)#
- kernel(nroots=1, koopmans=False, guess=None, eris=None, imds=None)#
Calculate N-electron neutral excitations via EOM-EE-CCSD.
- Kwargs:
- nrootsint
Number of roots (eigenvalues) requested
- koopmansbool
Calculate Koopmans’-like (1p1h) excitations only, targeting via overlap.
- guesslist of ndarray
List of guess vectors to use for targeting via overlap.
- class pyscf.cc.eom_rccsd.EOMEESinglet(cc)[source]#
Bases:
EOMEE
- static amplitudes_to_vector(t1, t2, out=None)#
- eomee_ccsd_singlet(nroots=1, koopmans=False, guess=None, eris=None, imds=None, diag=None)#
EOM-EE-CCSD singlet
- kernel(nroots=1, koopmans=False, guess=None, eris=None, imds=None, diag=None)#
EOM-EE-CCSD singlet
- matvec(vector, imds=None)#
- static spatial2spin(rx, orbspin=None)#
- vector_to_amplitudes(vector=None, nmo=None, nocc=None)#
- class pyscf.cc.eom_rccsd.EOMEESpinFlip(cc)[source]#
Bases:
EOMEE
- static amplitudes_to_vector(t1, t2, out=None)#
- eomsf_ccsd(nroots=1, koopmans=False, guess=None, eris=None, imds=None, diag=None)#
Spin flip EOM-EE-CCSD
- kernel(nroots=1, koopmans=False, guess=None, eris=None, imds=None, diag=None)#
Spin flip EOM-EE-CCSD
- matvec(vector, imds=None)#
Spin flip EOM-CCSD
- static spatial2spin(rx, orbspin=None)#
- vector_to_amplitudes(vector=None, nmo=None, nocc=None)#
- class pyscf.cc.eom_rccsd.EOMEETriplet(cc)[source]#
Bases:
EOMEE
- static amplitudes_to_vector(t1, t2, out=None)#
- eomee_ccsd_triplet(nroots=1, koopmans=False, guess=None, eris=None, imds=None, diag=None)#
EOM-EE-CCSD triplet
- kernel(nroots=1, koopmans=False, guess=None, eris=None, imds=None, diag=None)#
EOM-EE-CCSD triplet
- matvec(vector, imds=None)#
- static spatial2spin(rx, orbspin=None)#
- vector_to_amplitudes(vector=None, nmo=None, nocc=None)#
- class pyscf.cc.eom_rccsd.EOMIP(cc)[source]#
Bases:
EOM
- static amplitudes_to_vector(r1, r2)#
- ccsd_star_contract(ipccsd_evals, ipccsd_evecs, lipccsd_evecs, imds=None)#
- property eip#
- get_diag(imds=None)#
- ipccsd(nroots=1, left=False, koopmans=False, guess=None, partition=None, eris=None, imds=None)#
Calculate (N-1)-electron charged excitations via IP-EOM-CCSD.
- Kwargs:
- nrootsint
Number of roots (eigenvalues) requested
- partitionbool or str
Use a matrix-partitioning for the doubles-doubles block. Can be None, ‘mp’ (Moller-Plesset, i.e. orbital energies on the diagonal), or ‘full’ (full diagonal elements).
- koopmansbool
Calculate Koopmans’-like (quasiparticle) excitations only, targeting via overlap.
- guesslist of ndarray
List of guess vectors to use for targeting via overlap.
- ipccsd_star(nroots=1, koopmans=False, right_guess=None, left_guess=None, eris=None, imds=None)#
Calculates CCSD* perturbative correction.
Simply calls the relevant kernel() function and perturb_star of the eom class.
- Returns:
- e_t_a_star (list of float):
The IP-CCSD* energy.
- kernel(nroots=1, left=False, koopmans=False, guess=None, partition=None, eris=None, imds=None)#
Calculate (N-1)-electron charged excitations via IP-EOM-CCSD.
- Kwargs:
- nrootsint
Number of roots (eigenvalues) requested
- partitionbool or str
Use a matrix-partitioning for the doubles-doubles block. Can be None, ‘mp’ (Moller-Plesset, i.e. orbital energies on the diagonal), or ‘full’ (full diagonal elements).
- koopmansbool
Calculate Koopmans’-like (quasiparticle) excitations only, targeting via overlap.
- guesslist of ndarray
List of guess vectors to use for targeting via overlap.
- l_matvec(vector, imds=None, diag=None)#
For left eigenvector
- matvec(vector, imds=None, diag=None)#
- static spatial2spin(rx, orbspin=None)#
- vector_to_amplitudes(vector=None, nmo=None, nocc=None)#
- class pyscf.cc.eom_rccsd.EOMIP_Ta(cc)[source]#
Bases:
EOMIP
Class for EOM IPCCSD(T)*(a) method by Matthews and Stanton.
- pyscf.cc.eom_rccsd.eaccsd(eom, nroots=1, left=False, koopmans=False, guess=None, partition=None, eris=None, imds=None)[source]#
Calculate (N+1)-electron charged excitations via EA-EOM-CCSD.
- Args:
See also ipccd()
- pyscf.cc.eom_rccsd.eaccsd_star(eom, nroots=1, koopmans=False, right_guess=None, left_guess=None, eris=None, imds=None, **kwargs)[source]#
Calculates CCSD* perturbative correction.
- Args:
See also ipccd_star()
- pyscf.cc.eom_rccsd.eaccsd_star_contract(eom, eaccsd_evals, eaccsd_evecs, leaccsd_evecs, imds=None)[source]#
- pyscf.cc.eom_rccsd.eeccsd(eom, nroots=1, koopmans=False, guess=None, eris=None, imds=None)[source]#
Calculate N-electron neutral excitations via EOM-EE-CCSD.
- Kwargs:
- nrootsint
Number of roots (eigenvalues) requested
- koopmansbool
Calculate Koopmans’-like (1p1h) excitations only, targeting via overlap.
- guesslist of ndarray
List of guess vectors to use for targeting via overlap.
- pyscf.cc.eom_rccsd.eomee_ccsd_singlet(eom, nroots=1, koopmans=False, guess=None, eris=None, imds=None, diag=None)[source]#
EOM-EE-CCSD singlet
- pyscf.cc.eom_rccsd.eomee_ccsd_triplet(eom, nroots=1, koopmans=False, guess=None, eris=None, imds=None, diag=None)[source]#
EOM-EE-CCSD triplet
- pyscf.cc.eom_rccsd.eomsf_ccsd(eom, nroots=1, koopmans=False, guess=None, eris=None, imds=None, diag=None)[source]#
Spin flip EOM-EE-CCSD
- pyscf.cc.eom_rccsd.ipccsd(eom, nroots=1, left=False, koopmans=False, guess=None, partition=None, eris=None, imds=None)[source]#
Calculate (N-1)-electron charged excitations via IP-EOM-CCSD.
- Kwargs:
- nrootsint
Number of roots (eigenvalues) requested
- partitionbool or str
Use a matrix-partitioning for the doubles-doubles block. Can be None, ‘mp’ (Moller-Plesset, i.e. orbital energies on the diagonal), or ‘full’ (full diagonal elements).
- koopmansbool
Calculate Koopmans’-like (quasiparticle) excitations only, targeting via overlap.
- guesslist of ndarray
List of guess vectors to use for targeting via overlap.
- pyscf.cc.eom_rccsd.ipccsd_star(eom, nroots=1, koopmans=False, right_guess=None, left_guess=None, eris=None, imds=None)[source]#
Calculates CCSD* perturbative correction.
Simply calls the relevant kernel() function and perturb_star of the eom class.
- Returns:
- e_t_a_star (list of float):
The IP-CCSD* energy.
- pyscf.cc.eom_rccsd.ipccsd_star_contract(eom, ipccsd_evals, ipccsd_evecs, lipccsd_evecs, imds=None)[source]#
- pyscf.cc.eom_rccsd.kernel(eom, nroots=1, koopmans=False, guess=None, left=False, eris=None, imds=None, **kwargs)[source]#
pyscf.cc.eom_uccsd module#
- class pyscf.cc.eom_uccsd.EOMEA(cc)[source]#
Bases:
EOMEA
- static amplitudes_to_vector(r1, r2)#
- ccsd_star_contract = None#
- eaccsd_star = None#
- get_diag(imds=None)#
- l_matvec = None#
- matvec(vector, imds=None, diag=None)#
For spin orbitals.
R2 operators of the form s_{ j}^{ab}, i.e. indices jb are coupled.
- static spatial2spin(rx, orbspin=None)#
Convert EOMEA spatial-orbital R1/R2 to spin-orbital R1/R2
- static spin2spatial(rx, orbspin)#
Convert EOMEA spin-orbital R1/R2 to spatial-orbital R1/R2
- vector_to_amplitudes(vector=None, nmo=None, nocc=None)#
- class pyscf.cc.eom_uccsd.EOMEE(cc)[source]#
Bases:
EOMEE
- eeccsd(nroots=1, koopmans=False, guess=None, eris=None, imds=None)#
Calculate N-electron neutral excitations via EOM-EE-CCSD.
- Kwargs:
- nrootsint
Number of roots (eigenvalues) requested
- koopmansbool
Calculate Koopmans’-like (1p1h) excitations only, targeting via overlap.
- guesslist of ndarray
List of guess vectors to use for targeting via overlap.
- get_diag(imds=None)#
- kernel(nroots=1, koopmans=False, guess=None, eris=None, imds=None)#
Calculate N-electron neutral excitations via EOM-EE-CCSD.
- Kwargs:
- nrootsint
Number of roots (eigenvalues) requested
- koopmansbool
Calculate Koopmans’-like (1p1h) excitations only, targeting via overlap.
- guesslist of ndarray
List of guess vectors to use for targeting via overlap.
- class pyscf.cc.eom_uccsd.EOMEESpinFlip(cc)[source]#
Bases:
EOMEE
- static amplitudes_to_vector(t1, t2, out=None)#
- eomsf_ccsd(nroots=1, koopmans=False, guess=None, eris=None, imds=None, diag=None)#
Spin flip EOM-EE-CCSD
- kernel(nroots=1, koopmans=False, guess=None, eris=None, imds=None, diag=None)#
Spin flip EOM-EE-CCSD
- matvec(vector, imds=None)#
Spin flip EOM-CCSD
- static spatial2spin(rx, orbspin=None)#
Convert spin-flip EOMEE spatial-orbital R1/R2 to spin-orbital R1/R2
- static spin2spatial(rx, orbspin)#
Convert EOMEE spin-orbital R1/R2 to spin-flip EOMEE spatial-orbital R1/R2
- vector_to_amplitudes(vector=None, nmo=None, nocc=None)#
- class pyscf.cc.eom_uccsd.EOMEESpinKeep(cc)[source]#
Bases:
EOMEE
- static amplitudes_to_vector(t1, t2, out=None)#
- eomee_ccsd(nroots=1, koopmans=False, guess=None, eris=None, imds=None, diag=None)#
- get_diag(imds=None)#
- kernel(nroots=1, koopmans=False, guess=None, eris=None, imds=None, diag=None)#
Calculate N-electron neutral excitations via EOM-EE-CCSD.
- Kwargs:
- nrootsint
Number of roots (eigenvalues) requested
- koopmansbool
Calculate Koopmans’-like (1p1h) excitations only, targeting via overlap.
- guesslist of ndarray
List of guess vectors to use for targeting via overlap.
- matvec(vector, imds=None)#
- static spatial2spin(tx, orbspin=None)#
Convert T1/T2 of spatial orbital representation to T1/T2 of spin-orbital representation
- static spin2spatial(tx, orbspin)#
Convert T1/T2 in spin-orbital basis to T1/T2 in spatial orbital basis
- vector_to_amplitudes(vector=None, nmo=None, nocc=None)#
- class pyscf.cc.eom_uccsd.EOMIP(cc)[source]#
Bases:
EOMIP
- static amplitudes_to_vector(r1, r2)#
For spin orbitals
- ccsd_star_contract = None#
- get_diag(imds=None)#
- ipccsd_star = None#
- l_matvec = None#
- matvec(vector, imds=None, diag=None)#
For spin orbitals R2 operators of the form s_{ij}^{ b}, i.e. indices jb are coupled.
- static spatial2spin(rx, orbspin=None)#
Convert EOMIP spatial-orbital R1/R2 to spin-orbital R1/R2
- static spin2spatial(rx, orbspin)#
Convert EOMIP spin-orbital R1/R2 to spatial-orbital R1/R2
- vector_to_amplitudes(vector=None, nmo=None, nocc=None)#
For spin orbitals
- pyscf.cc.eom_uccsd.eaccsd_matvec(eom, vector, imds=None, diag=None)[source]#
For spin orbitals.
R2 operators of the form s_{ j}^{ab}, i.e. indices jb are coupled.
- pyscf.cc.eom_uccsd.eeccsd(eom, nroots=1, koopmans=False, guess=None, eris=None, imds=None)[source]#
Calculate N-electron neutral excitations via EOM-EE-CCSD.
- Kwargs:
- nrootsint
Number of roots (eigenvalues) requested
- koopmansbool
Calculate Koopmans’-like (1p1h) excitations only, targeting via overlap.
- guesslist of ndarray
List of guess vectors to use for targeting via overlap.
- pyscf.cc.eom_uccsd.eomee_ccsd(eom, nroots=1, koopmans=False, guess=None, eris=None, imds=None, diag=None)[source]#
- pyscf.cc.eom_uccsd.eomsf_ccsd(eom, nroots=1, koopmans=False, guess=None, eris=None, imds=None, diag=None)[source]#
Spin flip EOM-EE-CCSD
- pyscf.cc.eom_uccsd.ipccsd_matvec(eom, vector, imds=None, diag=None)[source]#
For spin orbitals R2 operators of the form s_{ij}^{ b}, i.e. indices jb are coupled.
- pyscf.cc.eom_uccsd.spatial2spin_ea(rx, orbspin=None)[source]#
Convert EOMEA spatial-orbital R1/R2 to spin-orbital R1/R2
- pyscf.cc.eom_uccsd.spatial2spin_eomsf(rx, orbspin=None)[source]#
Convert spin-flip EOMEE spatial-orbital R1/R2 to spin-orbital R1/R2
- pyscf.cc.eom_uccsd.spatial2spin_ip(rx, orbspin=None)[source]#
Convert EOMIP spatial-orbital R1/R2 to spin-orbital R1/R2
- pyscf.cc.eom_uccsd.spin2spatial_ea(rx, orbspin)[source]#
Convert EOMEA spin-orbital R1/R2 to spatial-orbital R1/R2
- pyscf.cc.eom_uccsd.spin2spatial_eomsf(rx, orbspin)[source]#
Convert EOMEE spin-orbital R1/R2 to spin-flip EOMEE spatial-orbital R1/R2
pyscf.cc.gccsd module#
- class pyscf.cc.gccsd.GCCSD(mf, frozen=None, mo_coeff=None, mo_occ=None)[source]#
Bases:
CCSDBase
- EOMEA(*args, **kwargs)#
- EOMEA_Ta(*args, **kwargs)#
Class for EOM EACCSD(T)*(a) method by Matthews and Stanton.
- EOMEE(*args, **kwargs)#
- EOMIP(*args, **kwargs)#
- EOMIP_Ta(*args, **kwargs)#
Class for EOM IPCCSD(T)*(a) method by Matthews and Stanton.
- ccsd(t1=None, t2=None, eris=None, mbpt2=False)[source]#
Ground-state unrestricted (U)CCSD.
- Kwargs:
- mbpt2bool
Use one-shot MBPT2 approximation to CCSD.
- conv_tol = 1e-07#
- conv_tol_normt = 1e-06#
- energy(t1, t2, eris)#
CCSD correlation energy
- kernel(t1=None, t2=None, eris=None, mbpt2=False)[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.).
- make_rdm1(t1=None, t2=None, l1=None, l2=None, ao_repr=False, with_frozen=True, with_mf=True)[source]#
Un-relaxed 1-particle density matrix in MO space
- make_rdm2(t1=None, t2=None, l1=None, l2=None, ao_repr=False, with_frozen=True, with_dm1=True)[source]#
2-particle density matrix in MO space. The density matrix is stored as
dm2[p,r,q,s] = <p^+ q^+ s r>
- 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.
- update_amps(t1, t2, eris)#
pyscf.cc.gccsd_lambda module#
pyscf.cc.gccsd_rdm module#
- pyscf.cc.gccsd_rdm.make_rdm1(mycc, t1, t2, l1, l2, ao_repr=False, with_frozen=True, with_mf=True)[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.cc.gccsd_rdm.make_rdm2(mycc, t1, t2, l1, l2, ao_repr=False, with_frozen=True, with_dm1=True)[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.cc.gccsd_t module#
GHF-CCSD(T) with spin-orbital integrals
pyscf.cc.gccsd_t_lambda module#
Lambda equation of GHF-CCSD(T) with spin-orbital integrals
Ref: JCP 98, 8718 (1993); DOI:10.1063/1.464480 JCP 147, 044104 (2017); DOI:10.1063/1.4994918
pyscf.cc.gccsd_t_rdm module#
pyscf.cc.gccsd_t_slow module#
GHF-CCSD(T) with spin-orbital integrals
pyscf.cc.gintermediates module#
- pyscf.cc.gintermediates.get_t3p2_imds_slow(cc, t1, t2, eris=None, t3p2_ip_out=None, t3p2_ea_out=None)[source]#
Calculates T1, T2 amplitudes corrected by second-order T3 contribution and intermediates used in IP/EA-CCSD(T)a
- Args:
- cc (
GCCSD
): Object containing coupled-cluster results.
- t1 (
ndarray
): T1 amplitudes.
- t2 (
ndarray
): T2 amplitudes from which the T3[2] amplitudes are formed.
- eris (
_PhysicistsERIs
): Antisymmetrized electron-repulsion integrals in physicist’s notation.
- t3p2_ip_out (
ndarray
): Store results of the intermediate used in IP-EOM-CCSD(T)a.
- t3p2_ea_out (
ndarray
): Store results of the intermediate used in EA-EOM-CCSD(T)a.
- cc (
- Returns:
- delta_ccsd (float):
- Difference of perturbed and unperturbed CCSD ground-state energy,
energy(T1 + T1[2], T2 + T2[2]) - energy(T1, T2)
- pt1 (
ndarray
): Perturbatively corrected T1 amplitudes.
- pt2 (
ndarray
): Perturbatively corrected T2 amplitudes.
- Reference:
Matthews, J. F. Stanton “A new approach to approximate…”
JCP 145, 124102 (2016); DOI:10.1063/1.4962910, Equation 14
- Shavitt and Bartlett “Many-body Methods in Physics and Chemistry”
2009, Equation 10.33
pyscf.cc.momgfccsd module#
GF-CCSD solver via moment constraints.
See reference: Backhouse, Booth, arXiv:2206.13198 (2022).
- class pyscf.cc.momgfccsd.MomGFCCSD(mycc, niter=(2, 2))[source]#
Bases:
StreamObject
Green’s function coupled cluster singles and doubles using the moment-resolved solver.
- Attributes:
- verboseint
Print level. Default value equals to
Mole.verbose
.- nitertuple of (int, int)
Number of block Lanczos iterations for occupied and virtual sectors. If either are None then said sector will not be computed.
- weight_tolfloat
Threshold for weight in the physical space to consider a pole an ionisation or removal event. Default value is 1e-1.
- hermi_momentsbool
Whether to Hermitise the moments, default value is False.
- hermi_solverobol
Whether to use the real-valued, symmetric block Lanczos solver, default value is False.
- Results:
- ehndarray
Energies of the compressed hole Green’s function
- vhtuple of ndarray
Left- and right-hand transition amplitudes of the compressed hole Green’s function
- epndarray
Energies of the compressed particle Green’s function
- vptuple of ndarray
Left- and right-hand transition amplitudes of the compressed particle Green’s function
- build_bra_hole(eom, t1, t2, l1, l2, orb)#
Get the first- and second-order contributions to the left-hand transformed vector for a given orbital for the hole part of the Green’s function.
- build_bra_part(eom, t1, t2, l1, l2, orb)#
Get the first- and second-order contributions to the left-hand transformed vector for a given orbital for the particle part of the Green’s function.
- build_hole_moments(t1=None, t2=None, l1=None, l2=None, imds=None, niter=None)[source]#
Build moments of the hole (IP-EOM-CCSD) Green’s function.
- build_part_moments(t1=None, t2=None, l1=None, l2=None, imds=None, niter=None)[source]#
Build moments of the particle (EA-EOM-CCSD) Green’s function.
- property eomea_method#
- property eomip_method#
- kernel(**kwargs)[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.).
- make_rdm1(ao_repr=False, eris=None, imds=None)[source]#
Build the first-order reduced density matrix at the CCSD level using the zeroth-order moment of the hole part of the CCSD Green’s function.
- property nmo#
- property nocc#
- update(chkfile=None, key='gfccsd')#
- update_from_chk(chkfile=None, key='gfccsd')#
- pyscf.cc.momgfccsd.block_lanczos_nosymm(gfccsd, moments, verbose=None)[source]#
Non-Hermitian block Lanczos solver, returns a set of poles that best reproduce the inputted moments.
- Args:
- gfccsdMomGFCCSD
GF-CCSD object
- momentsndarray (2*niter+2, n, n)
Array of moments with which the resulting poles should be consistent with.
- Kwargs:
- verboseint
Level of verbosity.
- Returns:
- andarray (niter+1, n, n)
On-diagonal blocks of the block tridiagonal Hamiltonian.
- bndarray (niter, n, n)
Upper off-diagonal blocks of the block tridiagonal Hamiltonian.
- cndarray (niter, n, n)
Lower off-diagonal blocks of the block tridiagonal Hamiltonian.
- pyscf.cc.momgfccsd.block_lanczos_symm(gfccsd, moments, verbose=None)[source]#
Hermitian block Lanczos solver, returns a set of poles that best reproduce the inputted moments.
- Args:
- gfccsdMomGFCCSD
GF-CCSD object
- momentsndarray (2*niter+2, n, n)
Array of moments with which the resulting poles should be consistent with.
- Kwargs:
- verboseint
Level of verbosity.
- Returns:
- andarray (niter+1, n, n)
On-diagonal blocks of the block tridiagonal Hamiltonian.
- bndarray (niter, n, n)
Off-diagonal blocks of the block tridiagonal Hamiltonian.
- pyscf.cc.momgfccsd.build_block_tridiagonal(a, b, c=None)[source]#
Construct a block tridiagonal matrix from a list of on-diagonal and off-diagonal blocks.
- pyscf.cc.momgfccsd.build_bra_hole(gfccsd, eom, t1, t2, l1, l2, orb)[source]#
Get the first- and second-order contributions to the left-hand transformed vector for a given orbital for the hole part of the Green’s function.
- pyscf.cc.momgfccsd.build_bra_part(gfccsd, eom, t1, t2, l1, l2, orb)[source]#
Get the first- and second-order contributions to the left-hand transformed vector for a given orbital for the particle part of the Green’s function.
- pyscf.cc.momgfccsd.contract_ket_hole(gfccsd, eom, t1, t2, v, orb)[source]#
Contract a vector with bar{a}^dagger_p |Psi>.
- pyscf.cc.momgfccsd.contract_ket_part(gfccsd, eom, t1, t2, v, orb)[source]#
Contract a vector with bar{a}_p |Psi>.
- pyscf.cc.momgfccsd.eig_block_tridiagonal(gfccsd, a, b, c, orth=None)[source]#
Diagonalise a non-Hermitian block-tridiagonal Hamiltonian and transform its eigenvectors appropriately.
- pyscf.cc.momgfccsd.eigh_block_tridiagonal(gfccsd, a, b, orth=None)[source]#
Diagonalise a Hermitian block-tridiagonal Hamiltonian and transform its eigenvectors appropriately.
- pyscf.cc.momgfccsd.kernel(gfccsd, hole_moments=None, part_moments=None, t1=None, t2=None, l1=None, l2=None, eris=None, imds=None, verbose=None)[source]#
pyscf.cc.qcisd module#
QCISD for real integrals 8-fold permutation symmetry has been used (ij|kl) = (ji|kl) = (kl|ij) = …
- class pyscf.cc.qcisd.QCISD(mf, frozen=None, mo_coeff=None, mo_occ=None)[source]#
Bases:
CCSDBase
restricted QCISD
- as_scanner()#
Generating a scanner/solver for CCSD PES.
The returned solver is a function. This function requires one argument “mol” as input and returns total CCSD energy.
The solver will automatically use the results of last calculation as the initial guess of the new calculation. All parameters assigned in the CCSD 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, cc >>> mol = gto.M(atom='H 0 0 0; F 0 0 1') >>> cc_scanner = cc.CCSD(scf.RHF(mol)).as_scanner() >>> e_tot = cc_scanner(gto.M(atom='H 0 0 0; F 0 0 1.1')) >>> e_tot = cc_scanner(gto.M(atom='H 0 0 0; F 0 0 1.5'))
- ccsd(*args, **kwargs)#
- energy(t1=None, t2=None, eris=None)#
QCISD correlation energy
- kernel(t1=None, t2=None, eris=None)#
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.).
- update_amps(t1, t2, eris)#
pyscf.cc.qcisd_slow module#
Restricted QCISD implementation The 4-index integrals are saved on disk entirely (without using any symmetry).
Note MO integrals are treated in chemist’s notation
Ref:
- class pyscf.cc.qcisd_slow.QCISD(mf, frozen=None, mo_coeff=None, mo_occ=None)[source]#
Bases:
RCCSD
restricted QCISD
- kernel(t1=None, t2=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.).
- update_amps(t1, t2, eris)#
pyscf.cc.qcisd_t module#
RHF-QCISD(T) for real integrals
pyscf.cc.qcisd_t_slow module#
pyscf.cc.rccsd module#
Restricted CCSD implementation which supports both real and complex integrals. The 4-index integrals are saved on disk entirely (without using any symmetry). This code is slower than the pyscf.cc.ccsd implementation.
Note MO integrals are treated in chemist’s notation
Ref: Hirata et al., J. Chem. Phys. 120, 2581 (2004)
- class pyscf.cc.rccsd.RCCSD(mf, frozen=None, mo_coeff=None, mo_occ=None)[source]#
Bases:
CCSD
restricted CCSD with IP-EOM, EA-EOM, EE-EOM, and SF-EOM capabilities
Ground-state CCSD is performed in optimized ccsd.CCSD and EOM is performed here.
- ccsd(t1=None, t2=None, eris=None, mbpt2=False)[source]#
Ground-state CCSD.
- Kwargs:
- mbpt2bool
Use one-shot MBPT2 approximation to CCSD.
- energy(t1=None, t2=None, eris=None)#
RCCSD correlation energy
- kernel(t1=None, t2=None, eris=None, mbpt2=False)[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.).
- update_amps(t1, t2, eris)#
pyscf.cc.rccsd_lambda module#
Restricted CCSD lambda equation solver which supports both real and complex integrals. This code is slower than the pyscf.cc.ccsd_lambda implementation.
Note MO integrals are treated in chemist’s notation
pyscf.cc.rccsd_slow module#
Restricted CCSD
Ref: Stanton et al., J. Chem. Phys. 94, 4334 (1990) Ref: Hirata et al., J. Chem. Phys. 120, 2581 (2004)
- class pyscf.cc.rccsd_slow.RCCSD(mf, frozen=None, mo_coeff=None, mo_occ=None)[source]#
Bases:
CCSD
- ccsd(t1=None, t2=None, eris=None, mbpt2=False, cc2=False)[source]#
Ground-state CCSD.
- Kwargs:
- mbpt2bool
Use one-shot MBPT2 approximation to CCSD.
- cc2bool
Use CC2 approximation to CCSD.
- eaccsd(nroots=1, left=False, koopmans=False, guess=None, partition=None)[source]#
Calculate (N+1)-electron charged excitations via EA-EOM-CCSD.
- Kwargs:
See ipccd()
- eeccsd(nroots=1, koopmans=False, guess=None, partition=None)[source]#
Calculate N-electron neutral excitations via EE-EOM-CCSD.
- Kwargs:
- nrootsint
Number of roots (eigenvalues) requested
- partitionbool or str
Use a matrix-partitioning for the doubles-doubles block. Can be None, ‘mp’ (Moller-Plesset, i.e. orbital energies on the diagonal), or ‘full’ (full diagonal elements).
- koopmansbool
Calculate Koopmans’-like (1p1h) excitations only, targeting via overlap.
- guesslist of ndarray
List of guess vectors to use for targeting via overlap.
- energy(t1, t2, eris)#
CCSD correlation energy
- ipccsd(nroots=1, left=False, koopmans=False, guess=None, partition=None)[source]#
Calculate (N-1)-electron charged excitations via IP-EOM-CCSD.
- Kwargs:
- nrootsint
Number of roots (eigenvalues) requested
- partitionbool or str
Use a matrix-partitioning for the doubles-doubles block. Can be None, ‘mp’ (Moller-Plesset, i.e. orbital energies on the diagonal), or ‘full’ (full diagonal elements).
- koopmansbool
Calculate Koopmans’-like (quasiparticle) excitations only, targeting via overlap.
- guesslist of ndarray
List of guess vectors to use for targeting via overlap.
- kernel(t1=None, t2=None, eris=None, mbpt2=False, cc2=False)[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.).
- update_amps(t1, t2, eris)#
pyscf.cc.rintermediates module#
Intermediates for restricted CCSD. Complex integrals are supported.
pyscf.cc.uccsd module#
UCCSD with spatial integrals
- class pyscf.cc.uccsd.UCCSD(mf, frozen=None, mo_coeff=None, mo_occ=None)[source]#
Bases:
CCSDBase
- EOMEA(*args, **kwargs)#
- EOMEE(*args, **kwargs)#
- EOMEESpinFlip(*args, **kwargs)#
- EOMEESpinKeep(*args, **kwargs)#
- EOMIP(*args, **kwargs)#
- ccsd(t1=None, t2=None, eris=None, mbpt2=False)[source]#
Ground-state unrestricted (U)CCSD.
- Kwargs:
- mbpt2bool
Use one-shot MBPT2 approximation to CCSD.
- conv_tol = 1e-07#
- conv_tol_normt = 1e-06#
- energy(t1=None, t2=None, eris=None)#
UCCSD correlation energy
- 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_nmo()#
- get_nocc()#
- kernel(t1=None, t2=None, eris=None, mbpt2=False)[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.).
- make_rdm1(t1=None, t2=None, l1=None, l2=None, ao_repr=False, with_frozen=True, with_mf=True)[source]#
Un-relaxed 1-particle density matrix in MO space
- Returns:
dm1a, dm1b
- make_rdm2(t1=None, t2=None, l1=None, l2=None, ao_repr=False, with_frozen=True, with_dm1=True)[source]#
2-particle density matrix in spin-orbital basis.
- 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.
- uccsd_t(t1=None, t2=None, eris=None)#
- update_amps(t1, t2, eris)#
pyscf.cc.uccsd_lambda module#
pyscf.cc.uccsd_rdm module#
- pyscf.cc.uccsd_rdm.make_rdm1(mycc, t1, t2, l1, l2, ao_repr=False, with_frozen=True, with_mf=True)[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.cc.uccsd_rdm.make_rdm2(mycc, t1, t2, l1, l2, ao_repr=False, with_frozen=True, with_dm1=True)[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.cc.uccsd_slow module#
- class pyscf.cc.uccsd_slow.UCCSD(mf, frozen=None, mo_coeff=None, mo_occ=None)[source]#
Bases:
CCSD
- ccsd(t1=None, t2=None, eris=None, mbpt2=False)[source]#
Ground-state unrestricted (U)CCSD.
- Kwargs:
- mbpt2bool
Use one-shot MBPT2 approximation to CCSD.
- energy(t1, t2, eris)#
CCSD correlation energy
- 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_nmo()#
- get_nocc()#
- kernel(t1=None, t2=None, eris=None, mbpt2=False)[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.).
- property nmo#
- property nocc#
- pyscf.cc.uccsd_slow.uspatial2spin(cc, moidx, mo_coeff)[source]#
Convert the results of an unrestricted mean-field calculation to spin-orbital form.
Spin-orbital ordering is determined by orbital energy without regard for spin.
- Returns:
- fock(nso,nso) ndarray
The Fock matrix in the basis of spin-orbitals
- so_coeff(nao, nso) ndarray
The matrix of spin-orbital coefficients in the AO basis
- spin(nso,) ndarary
The spin (0 or 1) of each spin-orbital
pyscf.cc.uccsd_t module#
UCCSD(T)
pyscf.cc.uccsd_t_lambda module#
Spin-free lambda equation of UHF-CCSD(T)
pyscf.cc.uccsd_t_rdm module#
pyscf.cc.uccsd_t_slow module#
pyscf.cc.uintermediates module#
pyscf.cc.uintermediates_slow module#
Module contents#
Coupled Cluster#
Simple usage:
>>> from pyscf import gto, scf, cc
>>> mol = gto.M(atom='H 0 0 0; H 0 0 1')
>>> mf = scf.RHF(mol).run()
>>> cc.CCSD(mf).run()
cc.CCSD()
returns an instance of CCSD class. Following are parameters
to control CCSD calculation.
- 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-7.
- conv_tol_normtfloat
converge threshold for norm(t1,t2). Default is 1e-5.
- max_cycleint
max number of iterations. Default is 50.
- diis_spaceint
DIIS space size. Default is 6.
- diis_start_cycleint
The step to start DIIS. Default is 0.
- directbool
AO-direct CCSD. Default is False.
- async_iobool
Allow for asynchronous function execution. Default is True.
- incore_completebool
Avoid all I/O. Default is False.
- frozenint or list
If integer is given, the inner-most orbitals are frozen from CC amplitudes. Given the orbital indices (0-based) in a list, both occupied and virtual orbitals can be frozen in CC calculation.
Saved results
- iterinfocommon.IterationInfo
Information about iteration (see pyscf.common.Iteration in detail)
- e_totfloat
Total CCSD energy (HF + correlation)
- t1, t2 :
t1[i,a], t2[i,j,a,b] (i,j in occ, a,b in virt)
- l1, l2 :
Lambda amplitudes l1[i,a], l2[i,j,a,b] (i,j in occ, a,b in virt)
- pyscf.cc.BCCD(mf, frozen=None, u=None, conv_tol_normu=1e-05, max_cycle=20, diis=True, canonicalization=True)[source]#
- pyscf.cc.CCSD(mf, frozen=None, mo_coeff=None, mo_occ=None)[source]#
restricted CCSD
- 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-7.
- conv_tol_normtfloat
converge threshold for norm(t1,t2). Default is 1e-5.
- max_cycleint
max number of iterations. Default is 50.
- diis_spaceint
DIIS space size. Default is 6.
- diis_start_cycleint
The step to start DIIS. Default is 0.
- iterative_dampingfloat
The self consistent damping parameter.
- directbool
AO-direct CCSD. Default is False.
- async_iobool
Allow for asynchronous function execution. Default is True.
- incore_completebool
Avoid all I/O (also for DIIS). Default is False.
- level_shiftfloat
A shift on virtual orbital energies to stabilize the CCSD iteration
- frozenint or list
If integer is given, the inner-most orbitals are frozen from CC amplitudes. Given the orbital indices (0-based) in a list, both occupied and virtual orbitals can be frozen in CC 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 >>> mycc = cc.CCSD(mf).set(frozen = 2).run() >>> # auto-generate the number of core orbitals to be frozen (1 in this case) >>> mycc = cc.CCSD(mf).set_frozen().run() >>> # freeze 2 core orbitals and 3 high lying unoccupied orbitals >>> mycc.set(frozen = [0,1,16,17,18]).run()
- callbackfunction(envs_dict) => None
callback function takes one dict as the argument which is generated by the builtin function
locals()
, so that the callback function can access all local variables in the current environment.
Saved results:
- convergedbool
Whether the CCSD iteration converged
- e_corrfloat
CCSD correlation correction
- e_totfloat
Total CCSD energy (HF + correlation)
- t1, t2 :
T amplitudes t1[i,a], t2[i,j,a,b] (i,j in occ, a,b in virt)
- l1, l2 :
Lambda amplitudes l1[i,a], l2[i,j,a,b] (i,j in occ, a,b in virt)
- cyclesint
The number of iteration cycles performed
- pyscf.cc.FNOCCSD(mf, thresh=1e-06, pct_occ=None, nvir_act=None, frozen=None)[source]#
Frozen natural orbital CCSD
- Attributes:
- threshfloat
Threshold on NO occupation numbers. Default is 1e-6 (very conservative).
- pct_occfloat
Percentage of total occupation number. Default is None. If present, overrides thresh.
- nvir_actint
Number of virtual NOs to keep. Default is None. If present, overrides thresh and pct_occ.
- pyscf.cc.RCCSD(mf, frozen=None, mo_coeff=None, mo_occ=None)[source]#
restricted CCSD
- 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-7.
- conv_tol_normtfloat
converge threshold for norm(t1,t2). Default is 1e-5.
- max_cycleint
max number of iterations. Default is 50.
- diis_spaceint
DIIS space size. Default is 6.
- diis_start_cycleint
The step to start DIIS. Default is 0.
- iterative_dampingfloat
The self consistent damping parameter.
- directbool
AO-direct CCSD. Default is False.
- async_iobool
Allow for asynchronous function execution. Default is True.
- incore_completebool
Avoid all I/O (also for DIIS). Default is False.
- level_shiftfloat
A shift on virtual orbital energies to stabilize the CCSD iteration
- frozenint or list
If integer is given, the inner-most orbitals are frozen from CC amplitudes. Given the orbital indices (0-based) in a list, both occupied and virtual orbitals can be frozen in CC 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 >>> mycc = cc.CCSD(mf).set(frozen = 2).run() >>> # auto-generate the number of core orbitals to be frozen (1 in this case) >>> mycc = cc.CCSD(mf).set_frozen().run() >>> # freeze 2 core orbitals and 3 high lying unoccupied orbitals >>> mycc.set(frozen = [0,1,16,17,18]).run()
- callbackfunction(envs_dict) => None
callback function takes one dict as the argument which is generated by the builtin function
locals()
, so that the callback function can access all local variables in the current environment.
Saved results:
- convergedbool
Whether the CCSD iteration converged
- e_corrfloat
CCSD correlation correction
- e_totfloat
Total CCSD energy (HF + correlation)
- t1, t2 :
T amplitudes t1[i,a], t2[i,j,a,b] (i,j in occ, a,b in virt)
- l1, l2 :
Lambda amplitudes l1[i,a], l2[i,j,a,b] (i,j in occ, a,b in virt)
- cyclesint
The number of iteration cycles performed