pyscf.mp package¶
Submodules¶
pyscf.mp.dfmp2 module¶
density fitting MP2, 3-center integrals incore.
-
class
pyscf.mp.dfmp2.
DFMP2
(mf, frozen=None, mo_coeff=None, mo_occ=None)¶ Bases:
pyscf.mp.mp2.MP2
-
ao2mo
(mo_coeff=None)¶
-
init_amps
(mo_energy=None, mo_coeff=None, eris=None, with_t2=True)¶
-
loop_ao2mo
(mo_coeff, nocc)¶
-
make_rdm1
(t2=None, ao_repr=False)¶ Spin-traced one-particle density matrix. The occupied-virtual orbital response is 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)
- Kwargs:
- ao_reprboolean
Whether to transfrom 1-particle density matrix to AO representation.
-
make_rdm2
(t2=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)
-
nuc_grad_method
()¶
-
reset
(mol=None)¶
-
update_amps
(t2, eris)¶ Update non-canonical MP2 amplitudes
-
-
pyscf.mp.dfmp2.
MP2
¶ alias of
pyscf.mp.dfmp2.DFMP2
-
pyscf.mp.dfmp2.
kernel
(mp, mo_energy=None, mo_coeff=None, eris=None, with_t2=True, verbose=3)¶
pyscf.mp.gmp2 module¶
GMP2 in spin-orbital form E(MP2) = 1/4 <ij||ab><ab||ij>/(ei+ej-ea-eb)
-
class
pyscf.mp.gmp2.
GMP2
(mf, frozen=None, mo_coeff=None, mo_occ=None)¶ Bases:
pyscf.mp.mp2.MP2
-
ao2mo
(mo_coeff=None)¶
-
energy
(t2, eris)¶ MP2 energy
-
init_amps
(mo_energy=None, mo_coeff=None, eris=None, with_t2=True)¶
-
make_rdm1
(t2=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
(t2=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)
-
nuc_grad_method
()¶
-
update_amps
(t2, eris)¶ Update non-canonical MP2 amplitudes
-
-
pyscf.mp.gmp2.
MP2
¶ alias of
pyscf.mp.gmp2.GMP2
-
pyscf.mp.gmp2.
energy
(mp, t2, eris)¶ MP2 energy
-
pyscf.mp.gmp2.
kernel
(mp, mo_energy=None, mo_coeff=None, eris=None, with_t2=True, verbose=None)¶
-
pyscf.mp.gmp2.
make_rdm1
(mp, t2=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)
-
pyscf.mp.gmp2.
make_rdm2
(mp, t2=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)
-
pyscf.mp.gmp2.
update_amps
(mp, t2, eris)¶ Update non-canonical MP2 amplitudes
pyscf.mp.mp2 module¶
RMP2
-
class
pyscf.mp.mp2.
MP2
(mf, frozen=None, mo_coeff=None, mo_occ=None)¶ Bases:
pyscf.lib.misc.StreamObject
restricted MP2 with canonical HF and non-canonical HF reference
- 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
For non-canonical MP2, converge threshold for MP2 correlation energy. Default value is 1e-7.
- conv_tol_normtfloat
For non-canonical MP2, converge threshold for norm(t2). Default value is 1e-5.
- max_cycleint
For non-canonical MP2, max number of MP2 iterations. Default value is 50.
- diis_spaceint
For non-canonical MP2, DIIS space size in MP2 iterations. Default is 6.
- level_shiftfloat
A shift on virtual orbital energies to stablize the MP2 iterations.
- frozenint or list
If integer is given, the inner-most orbitals are excluded from MP2 amplitudes. Given the orbital indices (0-based) in a list, both occupied and virtual orbitals can be frozen in MP2 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 >>> pt = mp.MP2(mf).set(frozen = 2).run() >>> # freeze 2 core orbitals and 3 high lying unoccupied orbitals >>> pt.set(frozen = [0,1,16,17,18]).run()
Saved results
- e_corrfloat
MP2 correlation correction
- e_totfloat
Total MP2 energy (HF + correlation)
- t2 :
T amplitudes t2[i,j,a,b] (i,j in occ, a,b in virt)
-
Gradients
(*args, **kwargs)¶
-
ao2mo
(mo_coeff=None)¶
-
as_scanner
()¶ Generating a scanner/solver for MP2 PES.
The returned solver is a function. This function requires one argument “mol” as input and returns total MP2 energy.
The solver will automatically use the results of last calculation as the initial guess of the new calculation. All parameters assigned in the MP2 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, mp >>> mol = gto.M(atom='H 0 0 0; F 0 0 1') >>> mp2_scanner = mp.MP2(scf.RHF(mol)).as_scanner() >>> e_tot = mp2_scanner(gto.M(atom='H 0 0 0; F 0 0 1.1')) >>> e_tot = mp2_scanner(gto.M(atom='H 0 0 0; F 0 0 1.5'))
-
conv_tol
= 1e-07¶
-
conv_tol_normt
= 1e-05¶
-
density_fit
(auxbasis=None, with_df=None)¶
-
dump_flags
(verbose=None)¶
-
property
e_tot
¶
-
property
emp2
¶
-
energy
(t2, eris)¶ MP2 energy
-
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 corresonds to the frozen orbital.
-
get_nmo
()¶
-
get_nocc
()¶
-
init_amps
(mo_energy=None, mo_coeff=None, eris=None, with_t2=True)¶
-
kernel
(mo_energy=None, mo_coeff=None, eris=None, with_t2=True)¶ - Args:
- with_t2bool
Whether to generate and hold t2 amplitudes in memory.
-
make_fno
(thresh=1e-06, pct_occ=None, nvir_act=None, t2=None)¶ Frozen natural orbitals
- Returns:
- frozenlist or ndarray
List of orbitals to freeze
- no_coeffndarray
Semicanonical NO coefficients in the AO basis
-
make_rdm1
(t2=None, eris=None, ao_repr=False)¶ Spin-traced one-particle density matrix. The occupied-virtual orbital response is 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)
- Kwargs:
- ao_reprboolean
Whether to transfrom 1-particle density matrix to AO representation.
-
make_rdm2
(t2=None, eris=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¶
-
property
nmo
¶
-
property
nocc
¶
-
nuc_grad_method
()¶
-
reset
(mol=None)¶
-
update_amps
(t2, eris)¶ Update non-canonical MP2 amplitudes
-
pyscf.mp.mp2.
RMP2
¶ alias of
pyscf.mp.mp2.MP2
-
pyscf.mp.mp2.
as_scanner
(mp)¶ Generating a scanner/solver for MP2 PES.
The returned solver is a function. This function requires one argument “mol” as input and returns total MP2 energy.
The solver will automatically use the results of last calculation as the initial guess of the new calculation. All parameters assigned in the MP2 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, mp >>> mol = gto.M(atom='H 0 0 0; F 0 0 1') >>> mp2_scanner = mp.MP2(scf.RHF(mol)).as_scanner() >>> e_tot = mp2_scanner(gto.M(atom='H 0 0 0; F 0 0 1.1')) >>> e_tot = mp2_scanner(gto.M(atom='H 0 0 0; F 0 0 1.5'))
-
pyscf.mp.mp2.
energy
(mp, t2, eris)¶ MP2 energy
-
pyscf.mp.mp2.
get_frozen_mask
(mp)¶ Get boolean mask for the restricted reference orbitals.
In the returned boolean (mask) array of frozen orbital indices, the element is False if it corresonds to the frozen orbital.
-
pyscf.mp.mp2.
get_nmo
(mp)¶
-
pyscf.mp.mp2.
get_nocc
(mp)¶
-
pyscf.mp.mp2.
kernel
(mp, mo_energy=None, mo_coeff=None, eris=None, with_t2=True, verbose=None)¶
-
pyscf.mp.mp2.
make_fno
(mp, thresh=1e-06, pct_occ=None, nvir_act=None, t2=None)¶ Frozen natural orbitals
- Returns:
- frozenlist or ndarray
List of orbitals to freeze
- no_coeffndarray
Semicanonical NO coefficients in the AO basis
-
pyscf.mp.mp2.
make_rdm1
(mp, t2=None, eris=None, ao_repr=False)¶ Spin-traced one-particle density matrix. The occupied-virtual orbital response is 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)
- Kwargs:
- ao_reprboolean
Whether to transfrom 1-particle density matrix to AO representation.
-
pyscf.mp.mp2.
make_rdm2
(mp, t2=None, eris=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)
-
pyscf.mp.mp2.
update_amps
(mp, t2, eris)¶ Update non-canonical MP2 amplitudes
pyscf.mp.mp2f12_slow module¶
MP2-F12 (In testing)
Refs: * JCC 32, 2492 (2011); DOI:10.1002/jcc.21825 * JCP 139, 084112 (2013); DOI:10.1063/1.4818753
With strong orthogonalization ansatz 2
-
pyscf.mp.mp2f12_slow.
energy_f12
(mf, auxmol, zeta)¶
-
pyscf.mp.mp2f12_slow.
find_cabs
(mol, auxmol, lindep=1e-08)¶
-
pyscf.mp.mp2f12_slow.
trans
(eri, mos)¶
pyscf.mp.ump2 module¶
UMP2 with spatial integals
-
pyscf.mp.ump2.
MP2
¶ alias of
pyscf.mp.ump2.UMP2
-
class
pyscf.mp.ump2.
UMP2
(mf, frozen=None, mo_coeff=None, mo_occ=None)¶ Bases:
pyscf.mp.mp2.MP2
-
Gradients
(*args, **kwargs)¶
-
ao2mo
(mo_coeff=None)¶
-
energy
(t2, eris)¶ MP2 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 corresonds to the frozen orbital.
-
get_nmo
()¶
-
get_nocc
()¶
-
init_amps
(mo_energy=None, mo_coeff=None, eris=None, with_t2=True)¶
-
make_fno
(thresh=1e-06, pct_occ=None, t2=None, eris=None)¶ Frozen natural orbitals
- Returns:
- frozenlist or ndarray
Length-2 list of orbitals to freeze
- no_coeffndarray
Length-2 list of semicanonical NO coefficients in the AO basis
-
make_rdm1
(t2=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
(t2=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 )
-
nuc_grad_method
()¶
-
update_amps
(t2, eris)¶ Update non-canonical MP2 amplitudes
-
-
pyscf.mp.ump2.
energy
(mp, t2, eris)¶ MP2 energy
-
pyscf.mp.ump2.
get_frozen_mask
(mp)¶ Get boolean mask for the unrestricted reference orbitals.
In the returned boolean (mask) array of frozen orbital indices, the element is False if it corresonds to the frozen orbital.
-
pyscf.mp.ump2.
get_nmo
(mp)¶
-
pyscf.mp.ump2.
get_nocc
(mp)¶
-
pyscf.mp.ump2.
kernel
(mp, mo_energy=None, mo_coeff=None, eris=None, with_t2=True, verbose=None)¶
-
pyscf.mp.ump2.
make_fno
(mp, thresh=1e-06, pct_occ=None, t2=None, eris=None)¶ Frozen natural orbitals
- Returns:
- frozenlist or ndarray
Length-2 list of orbitals to freeze
- no_coeffndarray
Length-2 list of semicanonical NO coefficients in the AO basis
-
pyscf.mp.ump2.
make_rdm1
(mp, t2=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).
-
pyscf.mp.ump2.
make_rdm2
(mp, t2=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 )
-
pyscf.mp.ump2.
update_amps
(mp, t2, eris)¶ Update non-canonical MP2 amplitudes
Module contents¶
Moller-Plesset perturbation theory
-
pyscf.mp.
GMP2
(mf, frozen=None, mo_coeff=None, mo_occ=None)¶
-
pyscf.mp.
MP2
(mf, frozen=None, mo_coeff=None, mo_occ=None)¶ restricted MP2 with canonical HF and non-canonical HF reference
- 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
For non-canonical MP2, converge threshold for MP2 correlation energy. Default value is 1e-7.
- conv_tol_normtfloat
For non-canonical MP2, converge threshold for norm(t2). Default value is 1e-5.
- max_cycleint
For non-canonical MP2, max number of MP2 iterations. Default value is 50.
- diis_spaceint
For non-canonical MP2, DIIS space size in MP2 iterations. Default is 6.
- level_shiftfloat
A shift on virtual orbital energies to stablize the MP2 iterations.
- frozenint or list
If integer is given, the inner-most orbitals are excluded from MP2 amplitudes. Given the orbital indices (0-based) in a list, both occupied and virtual orbitals can be frozen in MP2 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 >>> pt = mp.MP2(mf).set(frozen = 2).run() >>> # freeze 2 core orbitals and 3 high lying unoccupied orbitals >>> pt.set(frozen = [0,1,16,17,18]).run()
Saved results
- e_corrfloat
MP2 correlation correction
- e_totfloat
Total MP2 energy (HF + correlation)
- t2 :
T amplitudes t2[i,j,a,b] (i,j in occ, a,b in virt)
-
pyscf.mp.
RMP2
(mf, frozen=None, mo_coeff=None, mo_occ=None)¶ restricted MP2 with canonical HF and non-canonical HF reference
- 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
For non-canonical MP2, converge threshold for MP2 correlation energy. Default value is 1e-7.
- conv_tol_normtfloat
For non-canonical MP2, converge threshold for norm(t2). Default value is 1e-5.
- max_cycleint
For non-canonical MP2, max number of MP2 iterations. Default value is 50.
- diis_spaceint
For non-canonical MP2, DIIS space size in MP2 iterations. Default is 6.
- level_shiftfloat
A shift on virtual orbital energies to stablize the MP2 iterations.
- frozenint or list
If integer is given, the inner-most orbitals are excluded from MP2 amplitudes. Given the orbital indices (0-based) in a list, both occupied and virtual orbitals can be frozen in MP2 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 >>> pt = mp.MP2(mf).set(frozen = 2).run() >>> # freeze 2 core orbitals and 3 high lying unoccupied orbitals >>> pt.set(frozen = [0,1,16,17,18]).run()
Saved results
- e_corrfloat
MP2 correlation correction
- e_totfloat
Total MP2 energy (HF + correlation)
- t2 :
T amplitudes t2[i,j,a,b] (i,j in occ, a,b in virt)
-
pyscf.mp.
UMP2
(mf, frozen=None, mo_coeff=None, mo_occ=None)¶