pyscf.pbc.scf package¶
Submodules¶
pyscf.pbc.scf.addons module¶
-
pyscf.pbc.scf.addons.
canonical_occ
(mf, nelec=None)¶ Label the occupancies for each orbital for sampled k-points. This is for KUHF objects. Each k-point has a fixed number of up and down electrons in this, which results in a finite size error for metallic systems but can accelerate convergence.
-
pyscf.pbc.scf.addons.
canonical_occ_
(mf, nelec=None)¶ Label the occupancies for each orbital for sampled k-points. This is for KUHF objects. Each k-point has a fixed number of up and down electrons in this, which results in a finite size error for metallic systems but can accelerate convergence.
-
pyscf.pbc.scf.addons.
convert_to_ghf
(mf, out=None)¶ Convert the given mean-field object to the generalized HF/KS object
- Args:
mf : SCF object
- Returns:
An generalized SCF object
-
pyscf.pbc.scf.addons.
convert_to_khf
(mf, out=None)¶ Convert gamma point SCF object to k-point SCF object
-
pyscf.pbc.scf.addons.
convert_to_rhf
(mf, out=None)¶ Convert the given mean-field object to the corresponding restricted HF/KS object
-
pyscf.pbc.scf.addons.
convert_to_uhf
(mf, out=None)¶ Convert the given mean-field object to the corresponding unrestricted HF/KS object
-
pyscf.pbc.scf.addons.
project_mo_nr2nr
(cell1, mo1, cell2, kpts=None)¶ Project orbital coefficients
\[ \begin{align}\begin{aligned}|\psi1> = |AO1> C1\\|\psi2> = P |\psi1> = |AO2>S^{-1}<AO2| AO1> C1 = |AO2> C2\\C2 = S^{-1}<AO2|AO1> C1\end{aligned}\end{align} \]
-
pyscf.pbc.scf.addons.
smearing_
(mf, sigma=None, method='fermi', mu0=None)¶ Fermi-Dirac or Gaussian smearing
pyscf.pbc.scf.cphf module¶
Restricted coupled pertubed Hartree-Fock solver Modified from pyscf.scf.cphf
-
pyscf.pbc.scf.cphf.
kernel
(fvind, mo_energy, mo_occ, h1, s1=None, max_cycle=20, tol=1e-09, hermi=False, verbose=2)¶ - Args:
- fvindfunction
Given density matrix, compute (ij|kl)D_{lk}*2 - (ij|kl)D_{jk}
- Kwargs:
- hermiboolean
Whether the matrix defined by fvind is Hermitian or not.
-
pyscf.pbc.scf.cphf.
solve
(fvind, mo_energy, mo_occ, h1, s1=None, max_cycle=20, tol=1e-09, hermi=False, verbose=2)¶ - Args:
- fvindfunction
Given density matrix, compute (ij|kl)D_{lk}*2 - (ij|kl)D_{jk}
- Kwargs:
- hermiboolean
Whether the matrix defined by fvind is Hermitian or not.
-
pyscf.pbc.scf.cphf.
solve_nos1
(fvind, mo_energy, mo_occ, h1, max_cycle=20, tol=1e-09, hermi=False, verbose=2)¶ For field independent basis. First order overlap matrix is zero
-
pyscf.pbc.scf.cphf.
solve_withs1
(fvind, mo_energy, mo_occ, h1, s1, max_cycle=20, tol=1e-09, hermi=False, verbose=2)¶ For field dependent basis. First order overlap matrix is non-zero. The first order orbitals are set to C^1_{ij} = -1/2 S1 e1 = h1 - s1*e0 + (e0_j-e0_i)*c1 + vhf[c1]
- Kwargs:
- hermiboolean
Whether the matrix defined by fvind is Hermitian or not.
- Returns:
First order orbital coefficients (in MO basis) and first order orbital energy matrix
pyscf.pbc.scf.ghf module¶
Generalized Hartree-Fock for periodic systems at a single k-point
-
class
pyscf.pbc.scf.ghf.
GHF
(cell, kpt=array([0., 0., 0.]), exxdiv='ewald')¶ Bases:
pyscf.pbc.scf.hf.SCF
,pyscf.scf.ghf.GHF
GHF class for PBCs.
-
CCSD
(*args, **kwargs)¶ 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 stablize 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() >>> # freeze 2 core orbitals and 3 high lying unoccupied orbitals >>> mycc.set(frozen = [0,1,16,17,18]).run()
Saved results
- convergedbool
CCSD converged or not
- 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)
-
CISD
(*args, **kwargs)¶
-
MP2
(*args, **kwargs)¶
-
convert_from_
(mf)¶ Convert given mean-field object to RHF/ROHF
-
get_bands
(kpts_band, cell=None, dm=None, kpt=None)¶ Get energy bands at the given (arbitrary) ‘band’ k-points.
- Returns:
- mo_energy(nmo,) ndarray or a list of (nmo,) ndarray
Bands energies E_n(k)
- mo_coeff(nao, nmo) ndarray or a list of (nao,nmo) ndarray
Band orbitals psi_n(k)
-
get_grad
(mo_coeff, mo_occ, fock=None)¶ RHF orbital gradients
- Args:
- mo_coeff2D ndarray
Obital coefficients
- mo_occ1D ndarray
Orbital occupancy
- fock_ao2D ndarray
Fock matrix in AO representation
- Returns:
Gradients in MO representation. It’s a num_occ*num_vir vector.
-
get_hcore
(cell=None, kpt=None)¶
-
get_init_guess
(cell=None, key='minao')¶
-
get_j
(cell=None, dm=None, hermi=0, kpt=None, kpts_band=None, **kwargs)¶ Compute J matrix for the given density matrix and k-point (kpt). When kpts_band is given, the J matrices on kpts_band are evaluated.
J_{pq} = sum_{rs} (pq|rs) dm[s,r]
where r,s are orbitals on kpt. p and q are orbitals on kpts_band if kpts_band is given otherwise p and q are orbitals on kpt.
-
get_jk
(cell=None, dm=None, hermi=0, kpt=None, kpts_band=None, with_j=True, with_k=True, **kwargs)¶ Get Coulomb (J) and exchange (K) following
scf.hf.RHF.get_jk_()
. for particular k-point (kpt).When kpts_band is given, the J, K matrices on kpts_band are evaluated.
J_{pq} = sum_{rs} (pq|rs) dm[s,r] K_{pq} = sum_{rs} (pr|sq) dm[r,s]
where r,s are orbitals on kpt. p and q are orbitals on kpts_band if kpts_band is given otherwise p and q are orbitals on kpt.
-
get_k
(cell=None, dm=None, hermi=0, kpt=None, kpts_band=None, **kwargs)¶ Compute K matrix for the given density matrix.
-
get_occ
(mo_energy=None, mo_coeff=None)¶ Label the occupancies for each orbital
- Kwargs:
- mo_energy1D ndarray
Obital energies
- mo_coeff2D ndarray
Obital coefficients
Examples:
>>> from pyscf import gto, scf >>> mol = gto.M(atom='H 0 0 0; F 0 0 1.1') >>> mf = scf.hf.SCF(mol) >>> energy = numpy.array([-10., -1., 1, -2., 0, -3]) >>> mf.get_occ(energy) array([2, 2, 0, 2, 2, 2])
-
get_ovlp
(cell=None, kpt=None)¶
-
get_veff
(cell=None, dm=None, dm_last=0, vhf_last=0, hermi=1, kpt=None, kpts_band=None)¶ Hartree-Fock potential matrix for the given density matrix. See
scf.hf.get_veff()
andscf.hf.RHF.get_veff()
-
nuc_grad_method
= None¶
-
stability
= None¶
-
-
pyscf.pbc.scf.ghf.
get_jk
(mf, cell=None, dm=None, hermi=0, kpt=None, kpts_band=None, with_j=True, with_k=True, **kwargs)¶
pyscf.pbc.scf.hf module¶
Hartree-Fock for periodic systems at a single k-point
- See Also:
pyscf.pbc.scf.khf.py : Hartree-Fock for periodic systems with k-point sampling
-
class
pyscf.pbc.scf.hf.
RHF
(cell, kpt=array([0., 0., 0.]), exxdiv='ewald')¶ Bases:
pyscf.pbc.scf.hf.SCF
,pyscf.scf.hf.RHF
-
CCSD
(*args, **kwargs)¶ 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 stablize 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() >>> # freeze 2 core orbitals and 3 high lying unoccupied orbitals >>> mycc.set(frozen = [0,1,16,17,18]).run()
Saved results
- convergedbool
CCSD converged or not
- 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)
-
CISD
(*args, **kwargs)¶
-
MP2
(*args, **kwargs)¶ 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)
-
TDA
(*args, **kwargs)¶ Tamm-Dancoff approximation
- Attributes:
- conv_tolfloat
Diagonalization convergence tolerance. Default is 1e-9.
- nstatesint
Number of TD states to be computed. Default is 3.
Saved results:
- convergedbool
Diagonalization converged or not
- e1D array
excitation energy for each excited state.
- xyA list of two 2D arrays
The two 2D arrays are Excitation coefficients X (shape [nocc,nvir]) and de-excitation coefficients Y (shape [nocc,nvir]) for each excited state. (X,Y) are normalized to 1/2 in RHF/RKS methods and normalized to 1 for UHF/UKS methods. In the TDA calculation, Y = 0.
-
TDHF
(*args, **kwargs)¶ Time-dependent Hartree-Fock
- Attributes:
- conv_tolfloat
Diagonalization convergence tolerance. Default is 1e-9.
- nstatesint
Number of TD states to be computed. Default is 3.
Saved results:
- convergedbool
Diagonalization converged or not
- e1D array
excitation energy for each excited state.
- xyA list of two 2D arrays
The two 2D arrays are Excitation coefficients X (shape [nocc,nvir]) and de-excitation coefficients Y (shape [nocc,nvir]) for each excited state. (X,Y) are normalized to 1/2 in RHF/RKS methods and normalized to 1 for UHF/UKS methods. In the TDA calculation, Y = 0.
-
convert_from_
(mf)¶ Convert given mean-field object to RHF
-
stability
(internal=True, external=False, verbose=None)¶ RHF/RKS stability analysis.
See also pyscf.scf.stability.rhf_stability function.
- Kwargs:
- internalbool
Internal stability, within the RHF optimization space.
- externalbool
External stability. Including the RHF -> UHF and real -> complex stability analysis.
- Returns:
New orbitals that are more close to the stable condition. The return value includes two set of orbitals. The first corresponds to the internal stability and the second corresponds to the external stability.
-
-
class
pyscf.pbc.scf.hf.
SCF
(cell, kpt=array([0., 0., 0.]), exxdiv='ewald')¶ Bases:
pyscf.scf.hf.SCF
SCF base class adapted for PBCs.
- Attributes:
- kpt(3,) ndarray
The AO k-point in Cartesian coordinates, in units of 1/Bohr.
- exxdivstr
Exchange divergence treatment, can be one of
None : ignore G=0 contribution in exchange‘ewald’ : Ewald probe charge correction [JCP 122, 234102 (2005); DOI:10.1063/1.1926272]- with_dfdensity fitting object
Default is the FFT based DF model. For all-electron calculation, MDF model is favored for better accuracy. See also
pyscf.pbc.df
.- direct_scfbool
When this flag is set to true, the J/K matrices will be computed directly through the underlying with_df methods. Otherwise, depending the available memory, the 4-index integrals may be cached and J/K matrices are computed based on the 4-index integrals.
-
build
(cell=None)¶
-
check_sanity
()¶ Check input of class/object attributes, check whether a class method is overwritten. It does not check the attributes which are prefixed with “_”. The return value of method set is the object itself. This allows a series of functions/methods to be executed in pipe.
-
density_fit
(auxbasis=None, with_df=None)¶
-
dip_moment
(cell=None, dm=None, unit='Debye', verbose=3, **kwargs)¶ Dipole moment in the unit cell (is it well defined)?
- Args:
cell : an instance of
Cell
dm (ndarray) : density matrix
- Return:
A list: the dipole moment on x, y and z components
-
direct_scf
= False¶
-
dump_chk
(envs)¶
-
dump_flags
(verbose=None)¶
-
energy_nuc
()¶
-
from_chk
(chk=None, project=None, kpt=None)¶ Read the HF results from checkpoint file, then project it to the basis defined by
mol
- Returns:
Density matrix, 2D ndarray
-
get_bands
(kpts_band, cell=None, dm=None, kpt=None)¶ Get energy bands at the given (arbitrary) ‘band’ k-points.
- Returns:
- mo_energy(nmo,) ndarray or a list of (nmo,) ndarray
Bands energies E_n(k)
- mo_coeff(nao, nmo) ndarray or a list of (nao,nmo) ndarray
Band orbitals psi_n(k)
-
get_hcore
(cell=None, kpt=None)¶
-
get_init_guess
(cell=None, key='minao')¶
-
get_j
(cell=None, dm=None, hermi=1, kpt=None, kpts_band=None, omega=None)¶ Compute J matrix for the given density matrix and k-point (kpt). When kpts_band is given, the J matrices on kpts_band are evaluated.
J_{pq} = sum_{rs} (pq|rs) dm[s,r]
where r,s are orbitals on kpt. p and q are orbitals on kpts_band if kpts_band is given otherwise p and q are orbitals on kpt.
-
get_jk
(cell=None, dm=None, hermi=1, kpt=None, kpts_band=None, with_j=True, with_k=True, omega=None, **kwargs)¶ Get Coulomb (J) and exchange (K) following
scf.hf.RHF.get_jk_()
. for particular k-point (kpt).When kpts_band is given, the J, K matrices on kpts_band are evaluated.
J_{pq} = sum_{rs} (pq|rs) dm[s,r] K_{pq} = sum_{rs} (pr|sq) dm[r,s]
where r,s are orbitals on kpt. p and q are orbitals on kpts_band if kpts_band is given otherwise p and q are orbitals on kpt.
-
get_jk_incore
(cell=None, dm=None, hermi=1, kpt=None, omega=None, **kwargs)¶ Get Coulomb (J) and exchange (K) following
scf.hf.RHF.get_jk_()
.Incore version of Coulomb and exchange build only. Currently RHF always uses PBC AO integrals (unlike RKS), since exchange is currently computed by building PBC AO integrals.
-
get_k
(cell=None, dm=None, hermi=1, kpt=None, kpts_band=None, omega=None)¶ Compute K matrix for the given density matrix.
-
get_ovlp
(cell=None, kpt=None)¶
-
get_rho
(dm=None, grids=None, kpt=None)¶ Compute density in real space
-
get_veff
(cell=None, dm=None, dm_last=0, vhf_last=0, hermi=1, kpt=None, kpts_band=None)¶ Hartree-Fock potential matrix for the given density matrix. See
scf.hf.get_veff()
andscf.hf.RHF.get_veff()
-
init_guess_by_1e
(cell=None)¶ Generate initial guess density matrix from core hamiltonian
- Returns:
Density matrix, 2D ndarray
-
init_guess_by_chkfile
(chk=None, project=None, kpt=None)¶ Read the HF results from checkpoint file, then project it to the basis defined by
mol
- Returns:
Density matrix, 2D ndarray
-
property
kpt
¶
-
mix_density_fit
(auxbasis=None, with_df=None)¶
-
nuc_grad_method
(*args, **kwargs)¶ Hook to create object for analytical nuclear gradients.
-
reset
(cell=None)¶ Reset cell and relevant attributes associated to the old cell object
-
sfx2c1e
()¶
-
to_ghf
(mf)¶ Convert the input mean-field object to a GHF/GKS object
-
to_rhf
(mf)¶ Convert the input mean-field object to a RHF/ROHF/RKS/ROKS object
-
to_uhf
(mf)¶ Convert the input mean-field object to a UHF/UKS object
-
x2c
()¶
-
x2c1e
()¶
-
pyscf.pbc.scf.hf.
dip_moment
(cell, dm, unit='Debye', verbose=3, grids=None, rho=None, kpt=array([0., 0., 0.]), origin=None)¶ Dipole moment in the unit cell (is it well defined)?
- Args:
cell : an instance of
Cell
dm (ndarray) : density matrix
- Return:
A list: the dipole moment on x, y and z components
-
pyscf.pbc.scf.hf.
get_bands
(mf, kpts_band, cell=None, dm=None, kpt=None)¶ Get energy bands at the given (arbitrary) ‘band’ k-points.
- Returns:
- mo_energy(nmo,) ndarray or a list of (nmo,) ndarray
Bands energies E_n(k)
- mo_coeff(nao, nmo) ndarray or a list of (nao,nmo) ndarray
Band orbitals psi_n(k)
-
pyscf.pbc.scf.hf.
get_hcore
(cell, kpt=array([0., 0., 0.]))¶ Get the core Hamiltonian AO matrix.
-
pyscf.pbc.scf.hf.
get_j
(cell, dm, hermi=1, vhfopt=None, kpt=array([0., 0., 0.]), kpts_band=None)¶ Get the Coulomb (J) AO matrix for the given density matrix.
- Args:
- dmndarray or list of ndarrays
A density matrix or a list of density matrices
- Kwargs:
- hermiint
Whether J, K matrix is hermitian | 0 : no hermitian or symmetric | 1 : hermitian | 2 : anti-hermitian
- vhfopt :
A class which holds precomputed quantities to optimize the computation of J, K matrices
- kpt(3,) ndarray
The “inner” dummy k-point at which the DM was evaluated (or sampled).
- kpts_band(3,) ndarray or (*,3) ndarray
An arbitrary “band” k-point at which J is evaluated.
- Returns:
The function returns one J matrix, corresponding to the input density matrix (both order and shape).
-
pyscf.pbc.scf.hf.
get_jk
(mf, cell, dm, hermi=1, vhfopt=None, kpt=array([0., 0., 0.]), kpts_band=None, with_j=True, with_k=True, omega=None, **kwargs)¶ Get the Coulomb (J) and exchange (K) AO matrices for the given density matrix.
- Args:
- dmndarray or list of ndarrays
A density matrix or a list of density matrices
- Kwargs:
- hermiint
Whether J, K matrix is hermitian | 0 : no hermitian or symmetric | 1 : hermitian | 2 : anti-hermitian
- vhfopt :
A class which holds precomputed quantities to optimize the computation of J, K matrices
- kpt(3,) ndarray
The “inner” dummy k-point at which the DM was evaluated (or sampled).
- kpts_band(3,) ndarray or (*,3) ndarray
An arbitrary “band” k-point at which J and K are evaluated.
- Returns:
The function returns one J and one K matrix, corresponding to the input density matrix (both order and shape).
-
pyscf.pbc.scf.hf.
get_nuc
(cell, kpt=array([0., 0., 0.]))¶ Get the bare periodic nuc-el AO matrix, with G=0 removed.
See Martin (12.16)-(12.21).
-
pyscf.pbc.scf.hf.
get_ovlp
(cell, kpt=array([0., 0., 0.]))¶ Get the overlap AO matrix.
-
pyscf.pbc.scf.hf.
get_rho
(mf, dm=None, grids=None, kpt=None)¶ Compute density in real space
-
pyscf.pbc.scf.hf.
get_t
(cell, kpt=array([0., 0., 0.]))¶ Get the kinetic energy AO matrix.
-
pyscf.pbc.scf.hf.
init_guess_by_chkfile
(cell, chkfile_name, project=None, kpt=None)¶ Read the HF results from checkpoint file, then project it to the basis defined by
cell
- Returns:
Density matrix, (nao,nao) ndarray
-
pyscf.pbc.scf.hf.
makov_payne_correction
(mf)¶ Makov-Payne correction (Phys. Rev. B, 51, 4014)
-
pyscf.pbc.scf.hf.
normalize_dm_
(mf, dm)¶ Scale density matrix to make it produce the correct number of electrons.
pyscf.pbc.scf.kghf module¶
Generalized Hartree-Fock for periodic systems with k-point sampling
-
class
pyscf.pbc.scf.kghf.
KGHF
(cell, kpts=array([[0., 0., 0.]]), exxdiv='ewald')¶ Bases:
pyscf.pbc.scf.ghf.GHF
,pyscf.pbc.scf.khf.KSCF
GHF class for PBCs.
-
CCSD
(*args, **kwargs)¶ 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 stablize 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() >>> # freeze 2 core orbitals and 3 high lying unoccupied orbitals >>> mycc.set(frozen = [0,1,16,17,18]).run()
Saved results
- convergedbool
CCSD converged or not
- 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)
-
MP2
= None¶
-
convert_from_
(mf)¶ Convert given mean-field object to RHF/ROHF
-
density_fit
(auxbasis=None, with_df=None)¶
-
energy_elec
(dm_kpts=None, h1e_kpts=None, vhf_kpts=None)¶ Following pyscf.scf.hf.energy_elec()
-
gen_response
(mo_coeff=None, mo_occ=None, with_j=True, hermi=0, max_memory=None)¶ Generate a function to compute the product of KGHF response function and KGHF density matrices.
-
get_bands
(kpts_band, cell=None, dm_kpts=None, kpts=None)¶ Get energy bands at the given (arbitrary) ‘band’ k-points.
- Returns:
- mo_energy(nmo,) ndarray or a list of (nmo,) ndarray
Bands energies E_n(k)
- mo_coeff(nao, nmo) ndarray or a list of (nao,nmo) ndarray
Band orbitals psi_n(k)
-
get_grad
(mo_coeff_kpts, mo_occ_kpts, fock=None)¶ returns 1D array of gradients, like non K-pt version note that occ and virt indices of different k pts now occur in sequential patches of the 1D array
-
get_hcore
(cell=None, kpts=None)¶ Get the core Hamiltonian AO matrices at sampled k-points.
- Args:
kpts : (nkpts, 3) ndarray
- Returns:
hcore : (nkpts, nao, nao) ndarray
-
get_init_guess
(cell=None, key='minao')¶
-
get_j
(cell=None, dm_kpts=None, hermi=0, kpts=None, kpts_band=None)¶ Compute J matrix for the given density matrix and k-point (kpt). When kpts_band is given, the J matrices on kpts_band are evaluated.
J_{pq} = sum_{rs} (pq|rs) dm[s,r]
where r,s are orbitals on kpt. p and q are orbitals on kpts_band if kpts_band is given otherwise p and q are orbitals on kpt.
-
get_jk
(cell=None, dm_kpts=None, hermi=0, kpts=None, kpts_band=None, with_j=True, with_k=True, **kwargs)¶ Get Coulomb (J) and exchange (K) following
scf.hf.RHF.get_jk_()
. for particular k-point (kpt).When kpts_band is given, the J, K matrices on kpts_band are evaluated.
J_{pq} = sum_{rs} (pq|rs) dm[s,r] K_{pq} = sum_{rs} (pr|sq) dm[r,s]
where r,s are orbitals on kpt. p and q are orbitals on kpts_band if kpts_band is given otherwise p and q are orbitals on kpt.
-
get_k
(cell=None, dm_kpts=None, hermi=0, kpts=None, kpts_band=None)¶ Compute K matrix for the given density matrix.
-
get_occ
(mo_energy_kpts=None, mo_coeff_kpts=None)¶ Label the occupancies for each orbital for sampled k-points.
This is a k-point version of scf.hf.SCF.get_occ
-
get_ovlp
(cell=None, kpts=None)¶ Get the overlap AO matrices at sampled k-points.
- Args:
kpts : (nkpts, 3) ndarray
- Returns:
ovlp_kpts : (nkpts, nao, nao) ndarray
-
get_veff
(cell=None, dm_kpts=None, dm_last=0, vhf_last=0, hermi=1, kpts=None, kpts_band=None)¶ Hartree-Fock potential matrix for the given density matrix. See
scf.hf.get_veff()
andscf.hf.RHF.get_veff()
-
newton
()¶
-
nuc_grad_method
= None¶
-
stability
= None¶
-
x2c
= None¶
-
-
pyscf.pbc.scf.kghf.
get_jk
(mf, cell=None, dm_kpts=None, hermi=0, kpts=None, kpts_band=None, with_j=True, with_k=True, **kwargs)¶
-
pyscf.pbc.scf.kghf.
get_occ
(mf, mo_energy_kpts=None, mo_coeff_kpts=None)¶ Label the occupancies for each orbital for sampled k-points.
This is a k-point version of scf.hf.SCF.get_occ
pyscf.pbc.scf.khf module¶
Hartree-Fock for periodic systems with k-point sampling
- See Also:
hf.py : Hartree-Fock for periodic systems at a single k-point
-
class
pyscf.pbc.scf.khf.
KRHF
(cell, kpts=array([[0., 0., 0.]]), exxdiv='ewald')¶ Bases:
pyscf.pbc.scf.khf.KSCF
,pyscf.pbc.scf.hf.RHF
-
CCSD
(*args, **kwargs)¶ 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 stablize 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() >>> # freeze 2 core orbitals and 3 high lying unoccupied orbitals >>> mycc.set(frozen = [0,1,16,17,18]).run()
Saved results
- convergedbool
CCSD converged or not
- 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)
-
MP2
(*args, **kwargs)¶ 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)
-
TDA
(*args, **kwargs)¶ Tamm-Dancoff approximation
- Attributes:
- conv_tolfloat
Diagonalization convergence tolerance. Default is 1e-9.
- nstatesint
Number of TD states to be computed. Default is 3.
Saved results:
- convergedbool
Diagonalization converged or not
- e1D array
excitation energy for each excited state.
- xyA list of two 2D arrays
The two 2D arrays are Excitation coefficients X (shape [nocc,nvir]) and de-excitation coefficients Y (shape [nocc,nvir]) for each excited state. (X,Y) are normalized to 1/2 in RHF/RKS methods and normalized to 1 for UHF/UKS methods. In the TDA calculation, Y = 0.
-
TDHF
(*args, **kwargs)¶ Time-dependent Hartree-Fock
- Attributes:
- conv_tolfloat
Diagonalization convergence tolerance. Default is 1e-9.
- nstatesint
Number of TD states to be computed. Default is 3.
Saved results:
- convergedbool
Diagonalization converged or not
- e1D array
excitation energy for each excited state.
- xyA list of two 2D arrays
The two 2D arrays are Excitation coefficients X (shape [nocc,nvir]) and de-excitation coefficients Y (shape [nocc,nvir]) for each excited state. (X,Y) are normalized to 1/2 in RHF/RKS methods and normalized to 1 for UHF/UKS methods. In the TDA calculation, Y = 0.
-
check_sanity
()¶ Check input of class/object attributes, check whether a class method is overwritten. It does not check the attributes which are prefixed with “_”. The return value of method set is the object itself. This allows a series of functions/methods to be executed in pipe.
-
convert_from_
(mf)¶ Convert given mean-field object to KRHF
-
gen_response
(mo_coeff=None, mo_occ=None, singlet=None, hermi=0, max_memory=None)¶ Generate a function to compute the product of RHF response function and RHF density matrices.
- Kwargs:
- singlet (None or boolean)If singlet is None, response function for
orbital hessian or CPHF will be generated. If singlet is boolean, it is used in TDDFT response kernel.
-
nuc_grad_method
()¶ Hook to create object for analytical nuclear gradients.
-
-
class
pyscf.pbc.scf.khf.
KSCF
(cell, kpts=array([[0., 0., 0.]]), exxdiv='ewald')¶ Bases:
pyscf.pbc.scf.hf.SCF
SCF base class with k-point sampling.
Compared to molecular SCF, some members such as mo_coeff, mo_occ now have an additional first dimension for the k-points, e.g. mo_coeff is (nkpts, nao, nao) ndarray
- Attributes:
- kpts(nks,3) ndarray
The sampling k-points in Cartesian coordinates, in units of 1/Bohr.
-
analyze
(verbose=None, with_meta_lowdin=True, **kwargs)¶ Analyze the given SCF object: print orbital energies, occupancies; print orbital coefficients; Mulliken population analysis; Diople moment.
-
as_scanner
()¶ Generating a scanner/solver for HF PES.
The returned solver is a function. This function requires one argument “mol” as input and returns total HF energy.
The solver will automatically use the results of last calculation as the initial guess of the new calculation. All parameters assigned in the SCF object (DIIS, 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 >>> hf_scanner = scf.RHF(gto.Mole().set(verbose=0)).as_scanner() >>> hf_scanner(gto.M(atom='H 0 0 0; F 0 0 1.1')) -98.552190448277955 >>> hf_scanner(gto.M(atom='H 0 0 0; F 0 0 1.5')) -98.414750424294368
-
build
(cell=None)¶
-
canonicalize
(mo_coeff_kpts, mo_occ_kpts, fock=None)¶ Canonicalization diagonalizes the Fock matrix within occupied, open, virtual subspaces separatedly (without change occupancy).
-
check_sanity
()¶ Check input of class/object attributes, check whether a class method is overwritten. It does not check the attributes which are prefixed with “_”. The return value of method set is the object itself. This allows a series of functions/methods to be executed in pipe.
-
conv_tol_grad
= None¶
-
density_fit
(auxbasis=None, with_df=None)¶
-
dip_moment
(cell=None, dm=None, unit='Debye', verbose=3, **kwargs)¶ Dipole moment in the unit cell (is it well defined)?
- Args:
cell : an instance of
Cell
dm_kpts (a list of ndarrays) : density matrices of k-points
- Return:
A list: the dipole moment on x, y and z components
-
direct_scf
= False¶
-
dump_chk
(envs)¶
-
dump_flags
(verbose=None)¶
-
eig
(h_kpts, s_kpts)¶ Solver for generalized eigenvalue problem
\[HC = SCE\]
-
energy_elec
(dm_kpts=None, h1e_kpts=None, vhf_kpts=None)¶ Following pyscf.scf.hf.energy_elec()
-
from_chk
(chk=None, project=None, kpts=None)¶ Read the HF results from checkpoint file, then project it to the basis defined by
mol
- Returns:
Density matrix, 2D ndarray
-
get_bands
(kpts_band, cell=None, dm_kpts=None, kpts=None)¶ Get energy bands at the given (arbitrary) ‘band’ k-points.
- Returns:
- mo_energy(nmo,) ndarray or a list of (nmo,) ndarray
Bands energies E_n(k)
- mo_coeff(nao, nmo) ndarray or a list of (nao,nmo) ndarray
Band orbitals psi_n(k)
-
get_fermi
(mo_energy_kpts=None, mo_occ_kpts=None)¶ Fermi level
-
get_fock
(h1e=None, s1e=None, vhf=None, dm=None, cycle=- 1, diis=None, diis_start_cycle=None, level_shift_factor=None, damp_factor=None)¶ F = h^{core} + V^{HF}
Special treatment (damping, DIIS, or level shift) will be applied to the Fock matrix if diis and cycle is specified (The two parameters are passed to get_fock function during the SCF iteration)
- Kwargs:
- h1e2D ndarray
Core hamiltonian
- s1e2D ndarray
Overlap matrix, for DIIS
- vhf2D ndarray
HF potential matrix
- dm2D ndarray
Density matrix, for DIIS
- cycleint
Then present SCF iteration step, for DIIS
- diisan object of
SCF.DIIS
class DIIS object to hold intermediate Fock and error vectors
- diis_start_cycleint
The step to start DIIS. Default is 0.
- level_shift_factorfloat or int
Level shift (in AU) for virtual space. Default is 0.
-
get_grad
(mo_coeff_kpts, mo_occ_kpts, fock=None)¶ returns 1D array of gradients, like non K-pt version note that occ and virt indices of different k pts now occur in sequential patches of the 1D array
-
get_hcore
(cell=None, kpts=None)¶ Get the core Hamiltonian AO matrices at sampled k-points.
- Args:
kpts : (nkpts, 3) ndarray
- Returns:
hcore : (nkpts, nao, nao) ndarray
-
get_init_guess
(cell=None, key='minao')¶
-
get_j
(cell=None, dm_kpts=None, hermi=1, kpts=None, kpts_band=None, omega=None)¶ Compute J matrix for the given density matrix and k-point (kpt). When kpts_band is given, the J matrices on kpts_band are evaluated.
J_{pq} = sum_{rs} (pq|rs) dm[s,r]
where r,s are orbitals on kpt. p and q are orbitals on kpts_band if kpts_band is given otherwise p and q are orbitals on kpt.
-
get_jk
(cell=None, dm_kpts=None, hermi=1, kpts=None, kpts_band=None, with_j=True, with_k=True, omega=None, **kwargs)¶ Get Coulomb (J) and exchange (K) following
scf.hf.RHF.get_jk_()
. for particular k-point (kpt).When kpts_band is given, the J, K matrices on kpts_band are evaluated.
J_{pq} = sum_{rs} (pq|rs) dm[s,r] K_{pq} = sum_{rs} (pr|sq) dm[r,s]
where r,s are orbitals on kpt. p and q are orbitals on kpts_band if kpts_band is given otherwise p and q are orbitals on kpt.
-
get_k
(cell=None, dm_kpts=None, hermi=1, kpts=None, kpts_band=None, omega=None)¶ Compute K matrix for the given density matrix.
-
get_occ
(mo_energy_kpts=None, mo_coeff_kpts=None)¶ Label the occupancies for each orbital for sampled k-points.
This is a k-point version of scf.hf.SCF.get_occ
-
get_ovlp
(cell=None, kpts=None)¶ Get the overlap AO matrices at sampled k-points.
- Args:
kpts : (nkpts, 3) ndarray
- Returns:
ovlp_kpts : (nkpts, nao, nao) ndarray
-
get_rho
(dm=None, grids=None, kpts=None)¶ Compute density in real space
-
get_veff
(cell=None, dm_kpts=None, dm_last=0, vhf_last=0, hermi=1, kpts=None, kpts_band=None)¶ Hartree-Fock potential matrix for the given density matrix. See
scf.hf.get_veff()
andscf.hf.RHF.get_veff()
-
init_guess_by_1e
(cell=None)¶ Generate initial guess density matrix from core hamiltonian
- Returns:
Density matrix, 2D ndarray
-
init_guess_by_chkfile
(chk=None, project=None, kpts=None)¶ Read the HF results from checkpoint file, then project it to the basis defined by
mol
- Returns:
Density matrix, 2D ndarray
-
property
kpts
¶
-
make_rdm1
(mo_coeff_kpts=None, mo_occ_kpts=None, **kwargs)¶ One-particle density matrix in AO representation
- Args:
- mo_coeff2D ndarray
Orbital coefficients. Each column is one orbital.
- mo_occ1D ndarray
Occupancy
-
mix_density_fit
(auxbasis=None, with_df=None)¶
-
property
mo_coeff_kpts
¶
-
property
mo_energy_kpts
¶
-
property
mo_occ_kpts
¶
-
mulliken_meta
(cell=None, dm=None, verbose=5, pre_orth_method='ANO', s=None)¶ Mulliken population analysis, based on meta-Lowdin AOs. In the meta-lowdin, the AOs are grouped in three sets: core, valence and Rydberg, the orthogonalization are carreid out within each subsets.
- Args:
mol : an instance of
Mole
- dmndarray or 2-item list of ndarray
Density matrix. ROHF dm is a 2-item list of 2D array
- Kwargs:
verbose : int or instance of
lib.logger.Logger
- pre_orth_methodstr
Pre-orthogonalization, which localized GTOs for each atom. To obtain the occupied and unoccupied atomic shells, there are three methods
‘ano’ : Project GTOs to ANO basis‘minao’ : Project GTOs to MINAO basis‘scf’ : Symmetry-averaged fractional occupation atomic RHF
- Returns:
A list : pop, charges
- popnparray
Mulliken population on each atomic orbitals
- chargesnparray
Mulliken charges
-
mulliken_pop
()¶ Mulliken population analysis
\[M_{ij} = D_{ij} S_{ji}\]Mulliken charges
\[\delta_i = \sum_j M_{ij}\]- Returns:
A list : pop, charges
- popnparray
Mulliken population on each atomic orbitals
- chargesnparray
Mulliken charges
-
newton
()¶
-
sfx2c1e
()¶
-
stability
(internal=True, external=False, verbose=None)¶
-
to_ghf
(mf)¶ Convert the input mean-field object to a KGHF/KGKS object
-
to_rhf
(mf)¶ Convert the input mean-field object to a KRHF/KROHF/KRKS/KROKS object
-
to_uhf
(mf)¶ Convert the input mean-field object to a KUHF/KUKS object
-
x2c
()¶
-
x2c1e
()¶
-
pyscf.pbc.scf.khf.
analyze
(mf, verbose=5, with_meta_lowdin=True, **kwargs)¶ Analyze the given SCF object: print orbital energies, occupancies; print orbital coefficients; Mulliken population analysis; Dipole moment
-
pyscf.pbc.scf.khf.
as_scanner
(mf)¶
-
pyscf.pbc.scf.khf.
canonicalize
(mf, mo_coeff_kpts, mo_occ_kpts, fock=None)¶
-
pyscf.pbc.scf.khf.
dip_moment
(cell, dm_kpts, unit='Debye', verbose=3, grids=None, rho=None, kpts=array([[0., 0., 0.]]))¶ Dipole moment in the unit cell (is it well defined)?
- Args:
cell : an instance of
Cell
dm_kpts (a list of ndarrays) : density matrices of k-points
- Return:
A list: the dipole moment on x, y and z components
-
pyscf.pbc.scf.khf.
energy_elec
(mf, dm_kpts=None, h1e_kpts=None, vhf_kpts=None)¶ Following pyscf.scf.hf.energy_elec()
-
pyscf.pbc.scf.khf.
get_fermi
(mf, mo_energy_kpts=None, mo_occ_kpts=None)¶ Fermi level
-
pyscf.pbc.scf.khf.
get_fock
(mf, h1e=None, s1e=None, vhf=None, dm=None, cycle=- 1, diis=None, diis_start_cycle=None, level_shift_factor=None, damp_factor=None)¶
-
pyscf.pbc.scf.khf.
get_grad
(mo_coeff_kpts, mo_occ_kpts, fock)¶ returns 1D array of gradients, like non K-pt version note that occ and virt indices of different k pts now occur in sequential patches of the 1D array
-
pyscf.pbc.scf.khf.
get_hcore
(mf, cell=None, kpts=None)¶ Get the core Hamiltonian AO matrices at sampled k-points.
- Args:
kpts : (nkpts, 3) ndarray
- Returns:
hcore : (nkpts, nao, nao) ndarray
-
pyscf.pbc.scf.khf.
get_j
(mf, cell, dm_kpts, kpts, kpts_band=None)¶ Get the Coulomb (J) AO matrix at sampled k-points.
- Args:
- dm_kpts(nkpts, nao, nao) ndarray or a list of (nkpts,nao,nao) ndarray
Density matrix at each k-point. If a list of k-point DMs, eg, UHF alpha and beta DM, the alpha and beta DMs are contracted separately. It needs to be Hermitian.
- Kwargs:
- kpts_band(k,3) ndarray
A list of arbitrary “band” k-points at which to evalute the matrix.
- Returns:
vj : (nkpts, nao, nao) ndarray or list of vj if the input dm_kpts is a list of DMs
-
pyscf.pbc.scf.khf.
get_jk
(mf, cell, dm_kpts, kpts, kpts_band=None, with_j=True, with_k=True, omega=None, **kwargs)¶ Get the Coulomb (J) and exchange (K) AO matrices at sampled k-points.
- Args:
- dm_kpts(nkpts, nao, nao) ndarray
Density matrix at each k-point. It needs to be Hermitian.
- Kwargs:
- kpts_band(3,) ndarray
A list of arbitrary “band” k-point at which to evalute the matrix.
- Returns:
vj : (nkpts, nao, nao) ndarray vk : (nkpts, nao, nao) ndarray or list of vj and vk if the input dm_kpts is a list of DMs
-
pyscf.pbc.scf.khf.
get_occ
(mf, mo_energy_kpts=None, mo_coeff_kpts=None)¶ Label the occupancies for each orbital for sampled k-points.
This is a k-point version of scf.hf.SCF.get_occ
-
pyscf.pbc.scf.khf.
get_ovlp
(mf, cell=None, kpts=None)¶ Get the overlap AO matrices at sampled k-points.
- Args:
kpts : (nkpts, 3) ndarray
- Returns:
ovlp_kpts : (nkpts, nao, nao) ndarray
-
pyscf.pbc.scf.khf.
get_rho
(mf, dm=None, grids=None, kpts=None)¶ Compute density in real space
-
pyscf.pbc.scf.khf.
init_guess_by_chkfile
(cell, chkfile_name, project=None, kpts=None)¶ Read the KHF results from checkpoint file, then project it to the basis defined by
cell
- Returns:
Density matrix, 3D ndarray
-
pyscf.pbc.scf.khf.
make_rdm1
(mo_coeff_kpts, mo_occ_kpts, **kwargs)¶ One particle density matrices for all k-points.
- Returns:
dm_kpts : (nkpts, nao, nao) ndarray
-
pyscf.pbc.scf.khf.
mulliken_meta
(cell, dm_ao_kpts, verbose=5, pre_orth_method='ANO', s=None)¶ A modified Mulliken population analysis, based on meta-Lowdin AOs.
Note this function only computes the Mulliken population for the gamma point density matrix.
pyscf.pbc.scf.krohf module¶
Restricted open-shell Hartree-Fock for periodic systems with k-point sampling
-
class
pyscf.pbc.scf.krohf.
KROHF
(cell, kpts=array([[0., 0., 0.]]), exxdiv='ewald')¶ Bases:
pyscf.pbc.scf.khf.KRHF
,pyscf.pbc.scf.rohf.ROHF
UHF class with k-point sampling.
-
CCSD
= None¶
-
MP2
= None¶
-
TDA
= None¶
-
TDHF
= None¶
-
analyze
(verbose=None, with_meta_lowdin=True, **kwargs)¶ Analyze the given SCF object: print orbital energies, occupancies; print orbital coefficients; Mulliken population analysis
-
canonicalize
(mo_coeff_kpts, mo_occ_kpts, fock=None)¶ Canonicalization diagonalizes the ROHF Fock matrix within occupied, virtual subspaces separatedly (without change occupancy).
-
conv_tol
= 1e-07¶
-
conv_tol_grad
= None¶
-
convert_from_
(mf)¶ Convert given mean-field object to KUHF
-
dip_moment
(cell=None, dm=None, unit='Debye', verbose=3, **kwargs)¶ Dipole moment in the unit cell.
- Args:
cell : an instance of
Cell
dm_kpts (two lists of ndarrays) : KUHF density matrices of k-points
- Return:
A list: the dipole moment on x, y and z components
-
direct_scf
= False¶
-
dump_flags
(verbose=None)¶
-
eig
(fock, s)¶ Solver for generalized eigenvalue problem
\[HC = SCE\]
-
energy_elec
(dm_kpts=None, h1e_kpts=None, vhf_kpts=None)¶ Following pyscf.scf.hf.energy_elec()
-
gen_response
(mo_coeff=None, mo_occ=None, with_j=True, hermi=0, max_memory=None)¶ Generate a function to compute the product of UHF response function and UHF density matrices.
-
get_fock
(h1e=None, s1e=None, vhf=None, dm=None, cycle=- 1, diis=None, diis_start_cycle=None, level_shift_factor=None, damp_factor=None)¶ Build fock matrix based on Roothaan’s effective fock. See also
get_roothaan_fock()
-
get_grad
(mo_coeff_kpts, mo_occ_kpts, fock=None)¶ returns 1D array of gradients, like non K-pt version note that occ and virt indices of different k pts now occur in sequential patches of the 1D array
-
get_init_guess
(cell=None, key='minao')¶
-
get_occ
(mo_energy_kpts=None, mo_coeff_kpts=None)¶ Label the occupancies for each orbital for sampled k-points.
This is a k-point version of scf.hf.SCF.get_occ
-
get_rho
(dm=None, grids=None, kpts=None)¶ Compute density in real space
-
get_veff
(cell=None, dm_kpts=None, dm_last=0, vhf_last=0, hermi=1, kpts=None, kpts_band=None)¶ Hartree-Fock potential matrix for the given density matrix. See
scf.hf.get_veff()
andscf.hf.RHF.get_veff()
-
init_guess_by_atom
(mol=None)¶ Generate initial guess density matrix from superposition of atomic HF density matrix. The atomic HF is occupancy averaged RHF
- Returns:
Density matrix, 2D ndarray
-
init_guess_by_chkfile
(chk=None, project=True, kpts=None)¶ Read the HF results from checkpoint file, then project it to the basis defined by
mol
- Returns:
Density matrix, 2D ndarray
-
init_guess_by_huckel
(mol=None)¶ Generate initial guess density matrix from a Huckel calculation based on occupancy averaged atomic RHF calculations, doi:10.1021/acs.jctc.8b01089
- Returns:
Density matrix, 2D ndarray
-
init_guess_by_minao
(mol=None)¶ Generate initial guess density matrix based on ANO basis, then project the density matrix to the basis set defined by
mol
- Returns:
Density matrix, 2D ndarray
Examples:
>>> from pyscf import gto, scf >>> mol = gto.M(atom='H 0 0 0; H 0 0 1.1') >>> scf.hf.init_guess_by_minao(mol) array([[ 0.94758917, 0.09227308], [ 0.09227308, 0.94758917]])
-
make_rdm1
(mo_coeff_kpts=None, mo_occ_kpts=None, **kwargs)¶ One-particle densit matrix. mo_occ is a 1D array, with occupancy 1 or 2.
-
mulliken_meta
(cell=None, dm=None, verbose=5, pre_orth_method='ANO', s=None)¶ Mulliken population analysis, based on meta-Lowdin AOs. In the meta-lowdin, the AOs are grouped in three sets: core, valence and Rydberg, the orthogonalization are carreid out within each subsets.
- Args:
mol : an instance of
Mole
- dmndarray or 2-item list of ndarray
Density matrix. ROHF dm is a 2-item list of 2D array
- Kwargs:
verbose : int or instance of
lib.logger.Logger
- pre_orth_methodstr
Pre-orthogonalization, which localized GTOs for each atom. To obtain the occupied and unoccupied atomic shells, there are three methods
‘ano’ : Project GTOs to ANO basis‘minao’ : Project GTOs to MINAO basis‘scf’ : Symmetry-averaged fractional occupation atomic RHF
- Returns:
A list : pop, charges
- popnparray
Mulliken population on each atomic orbitals
- chargesnparray
Mulliken charges
-
property
nelec
¶
-
spin_square
(mo_coeff=None, s=None)¶ Spin square and multiplicity of RHF determinant
-
stability
(internal=True, external=False, verbose=None)¶ ROHF/ROKS stability analysis.
See also pyscf.scf.stability.rohf_stability function.
- Kwargs:
- internalbool
Internal stability, within the RHF optimization space.
- externalbool
External stability. It is not available in current version.
- Returns:
The return value includes two set of orbitals which are more close to the required stable condition.
-
-
pyscf.pbc.scf.krohf.
canonicalize
(mf, mo_coeff_kpts, mo_occ_kpts, fock=None)¶ Canonicalization diagonalizes the ROHF Fock matrix within occupied, virtual subspaces separatedly (without change occupancy).
-
pyscf.pbc.scf.krohf.
get_fock
(mf, h1e=None, s1e=None, vhf=None, dm=None, cycle=- 1, diis=None, diis_start_cycle=None, level_shift_factor=None, damp_factor=None)¶
-
pyscf.pbc.scf.krohf.
get_occ
(mf, mo_energy_kpts=None, mo_coeff_kpts=None)¶ Label the occupancies for each orbital for sampled k-points.
This is a k-point version of scf.hf.SCF.get_occ
-
pyscf.pbc.scf.krohf.
get_roothaan_fock
(focka_fockb, dma_dmb, s)¶ Roothaan’s effective fock.
space
closed
open
virtual
closed
Fc
Fb
Fc
open
Fb
Fc
Fa
virtual
Fc
Fa
Fc
where Fc = (Fa + Fb) / 2
- Returns:
Roothaan effective Fock matrix
-
pyscf.pbc.scf.krohf.
make_rdm1
(mo_coeff_kpts, mo_occ_kpts, **kwargs)¶ Alpha and beta spin one particle density matrices for all k-points.
- Returns:
dm_kpts : (2, nkpts, nao, nao) ndarray
-
pyscf.pbc.scf.krohf.
mulliken_meta
(cell, dm_ao_kpts, verbose=5, pre_orth_method='ANO', s=None)¶ A modified Mulliken population analysis, based on meta-Lowdin AOs.
Note this function only computes the Mulliken population for the gamma point density matrix.
pyscf.pbc.scf.kuhf module¶
Hartree-Fock for periodic systems with k-point sampling
- See Also:
hf.py : Hartree-Fock for periodic systems at a single k-point
-
class
pyscf.pbc.scf.kuhf.
KUHF
(cell, kpts=array([[0., 0., 0.]]), exxdiv='ewald')¶ Bases:
pyscf.pbc.scf.khf.KSCF
,pyscf.pbc.scf.uhf.UHF
UHF class with k-point sampling.
-
CCSD
(*args, **kwargs)¶ 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 stablize 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() >>> # freeze 2 core orbitals and 3 high lying unoccupied orbitals >>> mycc.set(frozen = [0,1,16,17,18]).run()
Saved results
- convergedbool
CCSD converged or not
- 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)
-
MP2
(*args, **kwargs)¶
-
TDA
(*args, **kwargs)¶
-
TDHF
(*args, **kwargs)¶
-
analyze
(verbose=None, with_meta_lowdin=True, **kwargs)¶ Analyze the given SCF object: print orbital energies, occupancies; print orbital coefficients; Mulliken population analysis; Diople moment.
-
canonicalize
(mo_coeff_kpts, mo_occ_kpts, fock=None)¶ Canonicalization diagonalizes the UHF Fock matrix within occupied, virtual subspaces separatedly (without change occupancy).
-
conv_tol
= 1e-07¶
-
conv_tol_grad
= None¶
-
convert_from_
(mf)¶ Convert given mean-field object to KUHF
-
dip_moment
(cell=None, dm=None, unit='Debye', verbose=3, **kwargs)¶ Dipole moment in the unit cell.
- Args:
cell : an instance of
Cell
dm_kpts (two lists of ndarrays) : KUHF density matrices of k-points
- Return:
A list: the dipole moment on x, y and z components
-
direct_scf
= False¶
-
dump_flags
(verbose=None)¶
-
eig
(h_kpts, s_kpts)¶ Solver for generalized eigenvalue problem
\[HC = SCE\]
-
energy_elec
(dm_kpts=None, h1e_kpts=None, vhf_kpts=None)¶ Following pyscf.scf.hf.energy_elec()
-
gen_response
(mo_coeff=None, mo_occ=None, with_j=True, hermi=0, max_memory=None)¶ Generate a function to compute the product of UHF response function and UHF density matrices.
-
get_bands
(kpts_band, cell=None, dm_kpts=None, kpts=None)¶ Get energy bands at the given (arbitrary) ‘band’ k-points.
- Returns:
- mo_energy(nmo,) ndarray or a list of (nmo,) ndarray
Bands energies E_n(k)
- mo_coeff(nao, nmo) ndarray or a list of (nao,nmo) ndarray
Band orbitals psi_n(k)
-
get_fermi
(mo_energy_kpts=None, mo_occ_kpts=None)¶ A pair of Fermi level for spin-up and spin-down orbitals
-
get_fock
(h1e=None, s1e=None, vhf=None, dm=None, cycle=- 1, diis=None, diis_start_cycle=None, level_shift_factor=None, damp_factor=None)¶ F = h^{core} + V^{HF}
Special treatment (damping, DIIS, or level shift) will be applied to the Fock matrix if diis and cycle is specified (The two parameters are passed to get_fock function during the SCF iteration)
- Kwargs:
- h1e2D ndarray
Core hamiltonian
- s1e2D ndarray
Overlap matrix, for DIIS
- vhf2D ndarray
HF potential matrix
- dm2D ndarray
Density matrix, for DIIS
- cycleint
Then present SCF iteration step, for DIIS
- diisan object of
SCF.DIIS
class DIIS object to hold intermediate Fock and error vectors
- diis_start_cycleint
The step to start DIIS. Default is 0.
- level_shift_factorfloat or int
Level shift (in AU) for virtual space. Default is 0.
-
get_grad
(mo_coeff_kpts, mo_occ_kpts, fock=None)¶ returns 1D array of gradients, like non K-pt version note that occ and virt indices of different k pts now occur in sequential patches of the 1D array
-
get_init_guess
(cell=None, key='minao')¶
-
get_occ
(mo_energy_kpts=None, mo_coeff_kpts=None)¶ Label the occupancies for each orbital for sampled k-points.
This is a k-point version of scf.hf.SCF.get_occ
-
get_rho
(dm=None, grids=None, kpts=None)¶ Compute density in real space
-
get_veff
(cell=None, dm_kpts=None, dm_last=0, vhf_last=0, hermi=1, kpts=None, kpts_band=None)¶ Hartree-Fock potential matrix for the given density matrix. See
scf.hf.get_veff()
andscf.hf.RHF.get_veff()
-
init_guess_by_1e
(cell=None)¶ Generate initial guess density matrix from core hamiltonian
- Returns:
Density matrix, 2D ndarray
-
init_guess_by_atom
(mol=None, breaksym=True)¶ Generate initial guess density matrix from superposition of atomic HF density matrix. The atomic HF is occupancy averaged RHF
- Returns:
Density matrix, 2D ndarray
-
init_guess_by_chkfile
(chk=None, project=True, kpts=None)¶ Read the HF results from checkpoint file, then project it to the basis defined by
mol
- Returns:
Density matrix, 2D ndarray
-
init_guess_by_huckel
(mol=None, breaksym=True)¶ Generate initial guess density matrix from a Huckel calculation based on occupancy averaged atomic RHF calculations, doi:10.1021/acs.jctc.8b01089
- Returns:
Density matrix, 2D ndarray
-
init_guess_by_minao
(mol=None, breaksym=True)¶ Initial guess in terms of the overlap to minimal basis.
-
make_rdm1
(mo_coeff_kpts=None, mo_occ_kpts=None, **kwargs)¶ One-particle density matrix
- Returns:
A list of 2D ndarrays for alpha and beta spins
-
mulliken_meta
(cell=None, dm=None, verbose=5, pre_orth_method='ANO', s=None)¶ A modified Mulliken population analysis, based on meta-Lowdin AOs.
Note this function only computes the Mulliken population for the gamma point density matrix.
-
mulliken_pop
()¶ Mulliken population analysis
\[M_{ij} = D_{ij} S_{ji}\]Mulliken charges
\[\delta_i = \sum_j M_{ij}\]- Returns:
A list : pop, charges
- popnparray
Mulliken population on each atomic orbitals
- chargesnparray
Mulliken charges
-
property
nelec
¶
-
nuc_grad_method
()¶ Hook to create object for analytical nuclear gradients.
-
spin_square
(mo_coeff=None, s=None)¶ Spin square and multiplicity of UHF determinant
\[S^2 = \frac{1}{2}(S_+ S_- + S_- S_+) + S_z^2\]where \(S_+ = \sum_i S_{i+}\) is effective for all beta occupied orbitals; \(S_- = \sum_i S_{i-}\) is effective for all alpha occupied orbitals.
- There are two possibilities for \(S_+ S_-\)
same electron \(S_+ S_- = \sum_i s_{i+} s_{i-}\),
\[\sum_i \langle UHF|s_{i+} s_{i-}|UHF\rangle = \sum_{pq}\langle p|s_+s_-|q\rangle \gamma_{qp} = n_\alpha\]2) different electrons \(S_+ S_- = \sum s_{i+} s_{j-}, (i\neq j)\). There are in total \(n(n-1)\) terms. As a two-particle operator,
\[\langle S_+ S_- \rangle = \langle ij|s_+ s_-|ij\rangle - \langle ij|s_+ s_-|ji\rangle = -\langle i^\alpha|j^\beta\rangle \langle j^\beta|i^\alpha\rangle\]
- Similarly, for \(S_- S_+\)
same electron
\[\sum_i \langle s_{i-} s_{i+}\rangle = n_\beta\]different electrons
\[\langle S_- S_+ \rangle = -\langle i^\beta|j^\alpha\rangle \langle j^\alpha|i^\beta\rangle\]
- For \(S_z^2\)
same electron
\[\langle s_z^2\rangle = \frac{1}{4}(n_\alpha + n_\beta)\]different electrons
\[\begin{split}&\frac{1}{2}\sum_{ij}(\langle ij|2s_{z1}s_{z2}|ij\rangle -\langle ij|2s_{z1}s_{z2}|ji\rangle) \\ &=\frac{1}{4}(\langle i^\alpha|i^\alpha\rangle \langle j^\alpha|j^\alpha\rangle - \langle i^\alpha|i^\alpha\rangle \langle j^\beta|j^\beta\rangle - \langle i^\beta|i^\beta\rangle \langle j^\alpha|j^\alpha\rangle + \langle i^\beta|i^\beta\rangle \langle j^\beta|j^\beta\rangle) \\ &-\frac{1}{4}(\langle i^\alpha|j^\alpha\rangle \langle j^\alpha|i^\alpha\rangle + \langle i^\beta|j^\beta\rangle\langle j^\beta|i^\beta\rangle) \\ &=\frac{1}{4}(n_\alpha^2 - n_\alpha n_\beta - n_\beta n_\alpha + n_\beta^2) -\frac{1}{4}(n_\alpha + n_\beta) \\ &=\frac{1}{4}((n_\alpha-n_\beta)^2 - (n_\alpha+n_\beta))\end{split}\]
In total
\[\begin{split}\langle S^2\rangle &= \frac{1}{2} (n_\alpha-\sum_{ij}\langle i^\alpha|j^\beta\rangle \langle j^\beta|i^\alpha\rangle +n_\beta -\sum_{ij}\langle i^\beta|j^\alpha\rangle\langle j^\alpha|i^\beta\rangle) + \frac{1}{4}(n_\alpha-n_\beta)^2 \\\end{split}\]- Args:
- moa list of 2 ndarrays
Occupied alpha and occupied beta orbitals
- Kwargs:
- sndarray
AO overlap
- Returns:
A list of two floats. The first is the expectation value of S^2. The second is the corresponding 2S+1
Examples:
>>> mol = gto.M(atom='O 0 0 0; H 0 0 1; H 0 1 0', basis='ccpvdz', charge=1, spin=1, verbose=0) >>> mf = scf.UHF(mol) >>> mf.kernel() -75.623975516256706 >>> mo = (mf.mo_coeff[0][:,mf.mo_occ[0]>0], mf.mo_coeff[1][:,mf.mo_occ[1]>0]) >>> print('S^2 = %.7f, 2S+1 = %.7f' % spin_square(mo, mol.intor('int1e_ovlp_sph'))) S^2 = 0.7570150, 2S+1 = 2.0070027
-
stability
(internal=True, external=False, verbose=None)¶ Stability analysis for RHF/RKS method.
See also pyscf.scf.stability.uhf_stability function.
- Args:
mf : UHF or UKS object
- Kwargs:
- internalbool
Internal stability, within the UHF space.
- externalbool
External stability. Including the UHF -> GHF and real -> complex stability analysis.
- Returns:
New orbitals that are more close to the stable condition. The return value includes two set of orbitals. The first corresponds to the internal stability and the second corresponds to the external stability.
-
-
pyscf.pbc.scf.kuhf.
canonicalize
(mf, mo_coeff_kpts, mo_occ_kpts, fock=None)¶ Canonicalization diagonalizes the UHF Fock matrix within occupied, virtual subspaces separatedly (without change occupancy).
-
pyscf.pbc.scf.kuhf.
dip_moment
(cell, dm_kpts, unit='Debye', verbose=3, grids=None, rho=None, kpts=array([[0., 0., 0.]]))¶ Dipole moment in the unit cell.
- Args:
cell : an instance of
Cell
dm_kpts (two lists of ndarrays) : KUHF density matrices of k-points
- Return:
A list: the dipole moment on x, y and z components
-
pyscf.pbc.scf.kuhf.
energy_elec
(mf, dm_kpts=None, h1e_kpts=None, vhf_kpts=None)¶ Following pyscf.scf.hf.energy_elec()
-
pyscf.pbc.scf.kuhf.
get_fermi
(mf, mo_energy_kpts=None, mo_occ_kpts=None)¶ A pair of Fermi level for spin-up and spin-down orbitals
-
pyscf.pbc.scf.kuhf.
get_fock
(mf, h1e=None, s1e=None, vhf=None, dm=None, cycle=- 1, diis=None, diis_start_cycle=None, level_shift_factor=None, damp_factor=None)¶
-
pyscf.pbc.scf.kuhf.
get_occ
(mf, mo_energy_kpts=None, mo_coeff_kpts=None)¶ Label the occupancies for each orbital for sampled k-points.
This is a k-point version of scf.hf.SCF.get_occ
-
pyscf.pbc.scf.kuhf.
init_guess_by_chkfile
(cell, chkfile_name, project=None, kpts=None)¶ Read the KHF results from checkpoint file, then project it to the basis defined by
cell
- Returns:
Density matrix, 3D ndarray
-
pyscf.pbc.scf.kuhf.
make_rdm1
(mo_coeff_kpts, mo_occ_kpts, **kwargs)¶ Alpha and beta spin one particle density matrices for all k-points.
- Returns:
dm_kpts : (2, nkpts, nao, nao) ndarray
-
pyscf.pbc.scf.kuhf.
mulliken_meta
(cell, dm_ao_kpts, verbose=5, pre_orth_method='ANO', s=None)¶ A modified Mulliken population analysis, based on meta-Lowdin AOs.
Note this function only computes the Mulliken population for the gamma point density matrix.
pyscf.pbc.scf.newton_ah module¶
Co-iterative augmented hessian second order SCF solver (CIAH-SOSCF)
-
pyscf.pbc.scf.newton_ah.
gen_g_hop_rhf
(mf, mo_coeff, mo_occ, fock_ao=None, h1e=None)¶
-
pyscf.pbc.scf.newton_ah.
gen_g_hop_rohf
(mf, mo_coeff, mo_occ, fock_ao=None, h1e=None)¶
-
pyscf.pbc.scf.newton_ah.
gen_g_hop_uhf
(mf, mo_coeff, mo_occ, fock_ao=None, h1e=None)¶
-
pyscf.pbc.scf.newton_ah.
newton
(mf)¶
pyscf.pbc.scf.rohf module¶
Restricted open-shell Hartree-Fock for periodic systems at a single k-point
- See Also:
pyscf/pbc/scf/khf.py : Hartree-Fock for periodic systems with k-point sampling
-
class
pyscf.pbc.scf.rohf.
ROHF
(cell, kpt=array([0., 0., 0.]), exxdiv='ewald')¶ Bases:
pyscf.pbc.scf.hf.RHF
,pyscf.scf.rohf.ROHF
ROHF class for PBCs.
-
CCSD
= None¶
-
CISD
= None¶
-
MP2
= None¶
-
TDA
= None¶
-
TDHF
= None¶
-
analyze
(verbose=None, with_meta_lowdin=True, **kwargs)¶ Analyze the given SCF object: print orbital energies, occupancies; print orbital coefficients; Mulliken population analysis
-
canonicalize
(mo_coeff, mo_occ, fock=None)¶ Canonicalization diagonalizes the Fock matrix within occupied, open, virtual subspaces separatedly (without change occupancy).
-
convert_from_
(mf)¶ Convert given mean-field object to RHF/ROHF
-
dip_moment
(cell=None, dm=None, unit='Debye', verbose=3, **kwargs)¶ Dipole moment in the unit cell.
- Args:
cell : an instance of
Cell
dm_kpts (a list of ndarrays) : density matrices of k-points
- Return:
A list: the dipole moment on x, y and z components
-
direct_scf
= False¶
-
dump_flags
(verbose=None)¶
-
eig
(fock, s)¶ Solver for generalized eigenvalue problem
\[HC = SCE\]
-
energy_elec
(dm=None, h1e=None, vhf=None)¶ Electronic part of Hartree-Fock energy, for given core hamiltonian and HF potential
… math:
E = \sum_{ij}h_{ij} \gamma_{ji} + \frac{1}{2}\sum_{ijkl} \gamma_{ji}\gamma_{lk} \langle ik||jl\rangle
Note this function has side effects which cause mf.scf_summary updated.
- Args:
mf : an instance of SCF class
- Kwargs:
- dm2D ndarray
one-partical density matrix
- h1e2D ndarray
Core hamiltonian
- vhf2D ndarray
HF potential
- Returns:
Hartree-Fock electronic energy and the Coulomb energy
Examples:
>>> from pyscf import gto, scf >>> mol = gto.M(atom='H 0 0 0; H 0 0 1.1') >>> mf = scf.RHF(mol) >>> mf.scf() >>> dm = mf.make_rdm1() >>> scf.hf.energy_elec(mf, dm) (-1.5176090667746334, 0.60917167853723675) >>> mf.energy_elec(dm) (-1.5176090667746334, 0.60917167853723675)
-
get_bands
(kpts_band, cell=None, dm=None, kpt=None)¶ Get energy bands at the given (arbitrary) ‘band’ k-points.
- Returns:
- mo_energy(nmo,) ndarray or a list of (nmo,) ndarray
Bands energies E_n(k)
- mo_coeff(nao, nmo) ndarray or a list of (nao,nmo) ndarray
Band orbitals psi_n(k)
-
get_fock
(h1e=None, s1e=None, vhf=None, dm=None, cycle=- 1, diis=None, diis_start_cycle=None, level_shift_factor=None, damp_factor=None)¶ Build fock matrix based on Roothaan’s effective fock. See also
get_roothaan_fock()
-
get_grad
(mo_coeff, mo_occ, fock=None)¶ ROHF gradients is the off-diagonal block [co + cv + ov], where [ cc co cv ] [ oc oo ov ] [ vc vo vv ]
-
get_init_guess
(cell=None, key='minao')¶
-
get_occ
(mo_energy=None, mo_coeff=None)¶ Label the occupancies for each orbital. NOTE the occupancies are not assigned based on the orbital energy ordering. The first N orbitals are assigned to be occupied orbitals.
Examples:
>>> mol = gto.M(atom='H 0 0 0; O 0 0 1.1', spin=1) >>> mf = scf.hf.SCF(mol) >>> energy = numpy.array([-10., -1., 1, -2., 0, -3]) >>> mf.get_occ(energy) array([2, 2, 2, 2, 1, 0])
-
get_rho
(dm=None, grids=None, kpt=None)¶ Compute density in real space
-
get_veff
(cell=None, dm=None, dm_last=0, vhf_last=0, hermi=1, kpt=None, kpts_band=None)¶ Hartree-Fock potential matrix for the given density matrix. See
scf.hf.get_veff()
andscf.hf.RHF.get_veff()
-
init_guess_by_1e
(cell=None)¶ Generate initial guess density matrix from core hamiltonian
- Returns:
Density matrix, 2D ndarray
-
init_guess_by_atom
(mol=None)¶ Generate initial guess density matrix from superposition of atomic HF density matrix. The atomic HF is occupancy averaged RHF
- Returns:
Density matrix, 2D ndarray
-
init_guess_by_chkfile
(chk=None, project=True, kpt=None)¶ Read the HF results from checkpoint file, then project it to the basis defined by
mol
- Returns:
Density matrix, 2D ndarray
-
init_guess_by_huckel
(mol=None)¶ Generate initial guess density matrix from a Huckel calculation based on occupancy averaged atomic RHF calculations, doi:10.1021/acs.jctc.8b01089
- Returns:
Density matrix, 2D ndarray
-
init_guess_by_minao
(mol=None)¶ Generate initial guess density matrix based on ANO basis, then project the density matrix to the basis set defined by
mol
- Returns:
Density matrix, 2D ndarray
Examples:
>>> from pyscf import gto, scf >>> mol = gto.M(atom='H 0 0 0; H 0 0 1.1') >>> scf.hf.init_guess_by_minao(mol) array([[ 0.94758917, 0.09227308], [ 0.09227308, 0.94758917]])
-
make_rdm1
(mo_coeff=None, mo_occ=None, **kwargs)¶ One-particle densit matrix. mo_occ is a 1D array, with occupancy 1 or 2.
-
property
nelec
¶
-
spin_square
(mo_coeff=None, s=None)¶ Spin square and multiplicity of RHF determinant
-
stability
(internal=True, external=False, verbose=None)¶ ROHF/ROKS stability analysis.
See also pyscf.scf.stability.rohf_stability function.
- Kwargs:
- internalbool
Internal stability, within the RHF optimization space.
- externalbool
External stability. It is not available in current version.
- Returns:
The return value includes two set of orbitals which are more close to the required stable condition.
-
pyscf.pbc.scf.scfint module¶
- SCF (Hartree-Fock and DFT) tools for periodic systems at a single k-point,
using analytical GTO integrals instead of PWs.
- See Also:
kscf.py : SCF tools for periodic systems with k-point sampling.
-
pyscf.pbc.scf.scfint.
get_hcore
(cell, kpt=None)¶ Get the core Hamiltonian AO matrix, following
dft.rks.get_veff_()
.
-
pyscf.pbc.scf.scfint.
get_int1e
(intor, cell, kpt=None)¶ Get the one-electron integral defined by intor using lattice sums.
-
pyscf.pbc.scf.scfint.
get_int1e_cross
(intor, cell1, cell2, kpt=None, comp=1)¶ 1-electron integrals from two molecules like
\[\langle \mu | intor | \nu \rangle, \mu \in cell1, \nu \in cell2\]
-
pyscf.pbc.scf.scfint.
get_ovlp
(cell, kpt=None)¶ Get the overlap AO matrix.
-
pyscf.pbc.scf.scfint.
get_t
(cell, kpt=None)¶ Get the kinetic energy AO matrix.
pyscf.pbc.scf.stability module¶
Wave Function Stability Analysis
Ref. JCP 66, 3045 (1977); DOI:10.1063/1.434318 JCP 104, 9047 (1996); DOI:10.1063/1.471637
See also tddft/rhf.py and scf/newton_ah.py
-
pyscf.pbc.scf.stability.
rhf_external
(mf, verbose=None)¶
-
pyscf.pbc.scf.stability.
rhf_internal
(mf, verbose=None)¶
-
pyscf.pbc.scf.stability.
rhf_stability
(mf, internal=True, external=False, verbose=None)¶
-
pyscf.pbc.scf.stability.
uhf_external
(mf, verbose=None)¶
-
pyscf.pbc.scf.stability.
uhf_internal
(mf, verbose=None)¶
-
pyscf.pbc.scf.stability.
uhf_stability
(mf, internal=True, external=False, verbose=None)¶
pyscf.pbc.scf.uhf module¶
Unrestricted Hartree-Fock for periodic systems at a single k-point
- See Also:
pyscf/pbc/scf/khf.py : Hartree-Fock for periodic systems with k-point sampling
-
class
pyscf.pbc.scf.uhf.
UHF
(cell, kpt=array([0., 0., 0.]), exxdiv='ewald')¶ Bases:
pyscf.pbc.scf.hf.SCF
,pyscf.scf.uhf.UHF
UHF class for PBCs.
-
CCSD
(*args, **kwargs)¶ 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 stablize 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() >>> # freeze 2 core orbitals and 3 high lying unoccupied orbitals >>> mycc.set(frozen = [0,1,16,17,18]).run()
Saved results
- convergedbool
CCSD converged or not
- 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)
-
CISD
(*args, **kwargs)¶
-
MP2
(*args, **kwargs)¶
-
TDA
(*args, **kwargs)¶
-
TDHF
(*args, **kwargs)¶
-
analyze
(verbose=None, with_meta_lowdin=True, **kwargs)¶ Analyze the given SCF object: print orbital energies, occupancies; print orbital coefficients; Mulliken population analysis; Diople moment.
-
canonicalize
(mo_coeff, mo_occ, fock=None)¶ Canonicalization diagonalizes the UHF Fock matrix within occupied, virtual subspaces separatedly (without change occupancy).
-
convert_from_
(mf)¶ Convert given mean-field object to RHF/ROHF
-
dip_moment
(cell=None, dm=None, unit='Debye', verbose=3, **kwargs)¶ Dipole moment in the unit cell.
- Args:
cell : an instance of
Cell
dm_kpts (a list of ndarrays) : density matrices of k-points
- Return:
A list: the dipole moment on x, y and z components
-
direct_scf
= False¶
-
dump_flags
(verbose=None)¶
-
eig
(fock, s)¶ Solver for generalized eigenvalue problem
\[HC = SCE\]
-
energy_elec
(dm=None, h1e=None, vhf=None)¶ Electronic energy of Unrestricted Hartree-Fock
Note this function has side effects which cause mf.scf_summary updated.
- Returns:
Hartree-Fock electronic energy and the 2-electron part contribution
-
get_bands
(kpts_band, cell=None, dm=None, kpt=None)¶ Get energy bands at the given (arbitrary) ‘band’ k-points.
- Returns:
- mo_energy(nmo,) ndarray or a list of (nmo,) ndarray
Bands energies E_n(k)
- mo_coeff(nao, nmo) ndarray or a list of (nao,nmo) ndarray
Band orbitals psi_n(k)
-
get_fock
(h1e=None, s1e=None, vhf=None, dm=None, cycle=- 1, diis=None, diis_start_cycle=None, level_shift_factor=None, damp_factor=None)¶ F = h^{core} + V^{HF}
Special treatment (damping, DIIS, or level shift) will be applied to the Fock matrix if diis and cycle is specified (The two parameters are passed to get_fock function during the SCF iteration)
- Kwargs:
- h1e2D ndarray
Core hamiltonian
- s1e2D ndarray
Overlap matrix, for DIIS
- vhf2D ndarray
HF potential matrix
- dm2D ndarray
Density matrix, for DIIS
- cycleint
Then present SCF iteration step, for DIIS
- diisan object of
SCF.DIIS
class DIIS object to hold intermediate Fock and error vectors
- diis_start_cycleint
The step to start DIIS. Default is 0.
- level_shift_factorfloat or int
Level shift (in AU) for virtual space. Default is 0.
-
get_grad
(mo_coeff, mo_occ, fock=None)¶ RHF orbital gradients
- Args:
- mo_coeff2D ndarray
Obital coefficients
- mo_occ1D ndarray
Orbital occupancy
- fock_ao2D ndarray
Fock matrix in AO representation
- Returns:
Gradients in MO representation. It’s a num_occ*num_vir vector.
-
get_init_guess
(cell=None, key='minao')¶
-
get_occ
(mo_energy=None, mo_coeff=None)¶ Label the occupancies for each orbital
- Kwargs:
- mo_energy1D ndarray
Obital energies
- mo_coeff2D ndarray
Obital coefficients
Examples:
>>> from pyscf import gto, scf >>> mol = gto.M(atom='H 0 0 0; F 0 0 1.1') >>> mf = scf.hf.SCF(mol) >>> energy = numpy.array([-10., -1., 1, -2., 0, -3]) >>> mf.get_occ(energy) array([2, 2, 0, 2, 2, 2])
-
get_rho
(dm=None, grids=None, kpt=None)¶ Compute density in real space
-
get_veff
(cell=None, dm=None, dm_last=0, vhf_last=0, hermi=1, kpt=None, kpts_band=None)¶ Hartree-Fock potential matrix for the given density matrix. See
scf.hf.get_veff()
andscf.hf.RHF.get_veff()
-
init_guess_by_1e
(cell=None)¶ Generate initial guess density matrix from core hamiltonian
- Returns:
Density matrix, 2D ndarray
-
init_guess_by_atom
(mol=None, breaksym=True)¶ Generate initial guess density matrix from superposition of atomic HF density matrix. The atomic HF is occupancy averaged RHF
- Returns:
Density matrix, 2D ndarray
-
init_guess_by_chkfile
(chk=None, project=True, kpt=None)¶ Read the HF results from checkpoint file, then project it to the basis defined by
mol
- Returns:
Density matrix, 2D ndarray
-
init_guess_by_huckel
(mol=None, breaksym=True)¶ Generate initial guess density matrix from a Huckel calculation based on occupancy averaged atomic RHF calculations, doi:10.1021/acs.jctc.8b01089
- Returns:
Density matrix, 2D ndarray
-
init_guess_by_minao
(mol=None, breaksym=True)¶ Initial guess in terms of the overlap to minimal basis.
-
make_rdm1
(mo_coeff=None, mo_occ=None, **kwargs)¶ One-particle density matrix
- Returns:
A list of 2D ndarrays for alpha and beta spins
-
mulliken_meta
(mol=None, dm=None, verbose=5, pre_orth_method='ANO', s=None)¶ Mulliken population analysis, based on meta-Lowdin AOs. In the meta-lowdin, the AOs are grouped in three sets: core, valence and Rydberg, the orthogonalization are carreid out within each subsets.
- Args:
mol : an instance of
Mole
- dmndarray or 2-item list of ndarray
Density matrix. ROHF dm is a 2-item list of 2D array
- Kwargs:
verbose : int or instance of
lib.logger.Logger
- pre_orth_methodstr
Pre-orthogonalization, which localized GTOs for each atom. To obtain the occupied and unoccupied atomic shells, there are three methods
‘ano’ : Project GTOs to ANO basis‘minao’ : Project GTOs to MINAO basis‘scf’ : Symmetry-averaged fractional occupation atomic RHF
- Returns:
A list : pop, charges
- popnparray
Mulliken population on each atomic orbitals
- chargesnparray
Mulliken charges
-
mulliken_pop
(mol=None, dm=None, s=None, verbose=5)¶ Mulliken population analysis
\[M_{ij} = D_{ij} S_{ji}\]Mulliken charges
\[\delta_i = \sum_j M_{ij}\]- Returns:
A list : pop, charges
- popnparray
Mulliken population on each atomic orbitals
- chargesnparray
Mulliken charges
-
property
nelec
¶
-
spin_square
(mo_coeff=None, s=None)¶ Spin square and multiplicity of UHF determinant
\[S^2 = \frac{1}{2}(S_+ S_- + S_- S_+) + S_z^2\]where \(S_+ = \sum_i S_{i+}\) is effective for all beta occupied orbitals; \(S_- = \sum_i S_{i-}\) is effective for all alpha occupied orbitals.
- There are two possibilities for \(S_+ S_-\)
same electron \(S_+ S_- = \sum_i s_{i+} s_{i-}\),
\[\sum_i \langle UHF|s_{i+} s_{i-}|UHF\rangle = \sum_{pq}\langle p|s_+s_-|q\rangle \gamma_{qp} = n_\alpha\]2) different electrons \(S_+ S_- = \sum s_{i+} s_{j-}, (i\neq j)\). There are in total \(n(n-1)\) terms. As a two-particle operator,
\[\langle S_+ S_- \rangle = \langle ij|s_+ s_-|ij\rangle - \langle ij|s_+ s_-|ji\rangle = -\langle i^\alpha|j^\beta\rangle \langle j^\beta|i^\alpha\rangle\]
- Similarly, for \(S_- S_+\)
same electron
\[\sum_i \langle s_{i-} s_{i+}\rangle = n_\beta\]different electrons
\[\langle S_- S_+ \rangle = -\langle i^\beta|j^\alpha\rangle \langle j^\alpha|i^\beta\rangle\]
- For \(S_z^2\)
same electron
\[\langle s_z^2\rangle = \frac{1}{4}(n_\alpha + n_\beta)\]different electrons
\[\begin{split}&\frac{1}{2}\sum_{ij}(\langle ij|2s_{z1}s_{z2}|ij\rangle -\langle ij|2s_{z1}s_{z2}|ji\rangle) \\ &=\frac{1}{4}(\langle i^\alpha|i^\alpha\rangle \langle j^\alpha|j^\alpha\rangle - \langle i^\alpha|i^\alpha\rangle \langle j^\beta|j^\beta\rangle - \langle i^\beta|i^\beta\rangle \langle j^\alpha|j^\alpha\rangle + \langle i^\beta|i^\beta\rangle \langle j^\beta|j^\beta\rangle) \\ &-\frac{1}{4}(\langle i^\alpha|j^\alpha\rangle \langle j^\alpha|i^\alpha\rangle + \langle i^\beta|j^\beta\rangle\langle j^\beta|i^\beta\rangle) \\ &=\frac{1}{4}(n_\alpha^2 - n_\alpha n_\beta - n_\beta n_\alpha + n_\beta^2) -\frac{1}{4}(n_\alpha + n_\beta) \\ &=\frac{1}{4}((n_\alpha-n_\beta)^2 - (n_\alpha+n_\beta))\end{split}\]
In total
\[\begin{split}\langle S^2\rangle &= \frac{1}{2} (n_\alpha-\sum_{ij}\langle i^\alpha|j^\beta\rangle \langle j^\beta|i^\alpha\rangle +n_\beta -\sum_{ij}\langle i^\beta|j^\alpha\rangle\langle j^\alpha|i^\beta\rangle) + \frac{1}{4}(n_\alpha-n_\beta)^2 \\\end{split}\]- Args:
- moa list of 2 ndarrays
Occupied alpha and occupied beta orbitals
- Kwargs:
- sndarray
AO overlap
- Returns:
A list of two floats. The first is the expectation value of S^2. The second is the corresponding 2S+1
Examples:
>>> mol = gto.M(atom='O 0 0 0; H 0 0 1; H 0 1 0', basis='ccpvdz', charge=1, spin=1, verbose=0) >>> mf = scf.UHF(mol) >>> mf.kernel() -75.623975516256706 >>> mo = (mf.mo_coeff[0][:,mf.mo_occ[0]>0], mf.mo_coeff[1][:,mf.mo_occ[1]>0]) >>> print('S^2 = %.7f, 2S+1 = %.7f' % spin_square(mo, mol.intor('int1e_ovlp_sph'))) S^2 = 0.7570150, 2S+1 = 2.0070027
-
stability
(internal=True, external=False, verbose=None)¶ Stability analysis for RHF/RKS method.
See also pyscf.scf.stability.uhf_stability function.
- Args:
mf : UHF or UKS object
- Kwargs:
- internalbool
Internal stability, within the UHF space.
- externalbool
External stability. Including the UHF -> GHF and real -> complex stability analysis.
- Returns:
New orbitals that are more close to the stable condition. The return value includes two set of orbitals. The first corresponds to the internal stability and the second corresponds to the external stability.
-
-
pyscf.pbc.scf.uhf.
dip_moment
(cell, dm, unit='Debye', verbose=3, grids=None, rho=None, kpt=array([0., 0., 0.]))¶ Dipole moment in the unit cell.
- Args:
cell : an instance of
Cell
dm_kpts (a list of ndarrays) : density matrices of k-points
- Return:
A list: the dipole moment on x, y and z components
-
pyscf.pbc.scf.uhf.
init_guess_by_chkfile
(cell, chkfile_name, project=None, kpt=None)¶ Read the HF results from checkpoint file and make the density matrix for UHF initial guess.
- Returns:
Density matrix, (nao,nao) ndarray
Module contents¶
Hartree-Fock for periodic systems
-
pyscf.pbc.scf.
HF
(cell, *args, **kwargs)¶
-
pyscf.pbc.scf.
KHF
(cell, *args, **kwargs)¶
-
pyscf.pbc.scf.
KKS
(cell, *args, **kwargs)¶
-
pyscf.pbc.scf.
KRKS
(cell, *args, **kwargs)¶
-
pyscf.pbc.scf.
KROKS
(cell, *args, **kwargs)¶
-
pyscf.pbc.scf.
KS
(cell, *args, **kwargs)¶
-
pyscf.pbc.scf.
KUKS
(cell, *args, **kwargs)¶
-
pyscf.pbc.scf.
RHF
(cell, *args, **kwargs)¶
-
pyscf.pbc.scf.
RKS
(cell, *args, **kwargs)¶
-
pyscf.pbc.scf.
ROKS
(cell, *args, **kwargs)¶
-
pyscf.pbc.scf.
UKS
(cell, *args, **kwargs)¶