pyscf.gw package#
Submodules#
pyscf.gw.gw_ac module#
Spin-restricted G0W0 approximation with analytic continuation This implementation has N^4 scaling, and is faster than GW-CD (N^4) and analytic GW (N^6) methods. GW-AC is recommended for valence states only, and is inaccurate for core states.
- Method:
See T. Zhu and G.K.-L. Chan, arxiv:2007.03148 (2020) for details Compute Sigma on imaginary frequency with density fitting, then analytically continued to real frequency
- Useful References:
J. Chem. Theory Comput. 12, 3623-3635 (2016) New J. Phys. 14, 053020 (2012)
- pyscf.gw.gw_ac.AC_pade_thiele_diag(sigma, omega)[source]#
Analytic continuation to real axis using a Pade approximation from Thiele’s reciprocal difference method Reference: J. Low Temp. Phys. 29, 179 (1977) Returns:
coeff: 2D array (ncoeff, norbs) omega: 2D array (norbs, npade)
- pyscf.gw.gw_ac.AC_twopole_diag(sigma, omega, orbs, nocc)[source]#
Analytic continuation to real axis using a two-pole model Returns:
coeff: 2D array (ncoeff, norbs)
- class pyscf.gw.gw_ac.GWAC(mf, frozen=None)[source]#
Bases:
StreamObject
- ac = 'pade'#
- get_frozen_mask()#
Get boolean mask for the restricted reference orbitals.
In the returned boolean (mask) array of frozen orbital indices, the element is False if it corresponds to the frozen orbital.
- get_nmo()#
- get_nocc()#
- kernel(mo_energy=None, mo_coeff=None, Lpq=None, orbs=None, nw=100, vhf_df=False)[source]#
- Input:
orbs: self-energy orbs nw: grid number vhf_df: whether using density fitting for HF exchange
- Output:
mo_energy: GW quasiparticle energy
- linearized = False#
- property nmo#
- property nocc#
- pyscf.gw.gw_ac.get_rho_response(omega, mo_energy, Lpq)[source]#
Compute density response function in auxiliary basis at freq iw
- pyscf.gw.gw_ac.get_sigma_diag(gw, orbs, Lpq, freqs, wts, iw_cutoff=None)[source]#
Compute GW correlation self-energy (diagonal elements) in MO basis on imaginary axis
pyscf.gw.gw_cd module#
Spin-restricted G0W0 approximation with contour deformation
This implementation has the same scaling (N^4) as GW-AC, more robust but slower. GW-CD is particularly recommended for accurate core and high-energy states.
- Method:
See T. Zhu and G.K.-L. Chan, arxiv:2007.03148 (2020) for details Compute Sigma directly on real axis with density fitting through a contour deformation method
- Useful References:
Chem. Theory Comput. 14, 4856-4869 (2018)
- class pyscf.gw.gw_cd.GWCD(mf, frozen=None)[source]#
Bases:
StreamObject
- eta = 0.001#
- get_frozen_mask()#
Get boolean mask for the restricted reference orbitals.
In the returned boolean (mask) array of frozen orbital indices, the element is False if it corresponds to the frozen orbital.
- get_nmo()#
- get_nocc()#
- kernel(mo_energy=None, mo_coeff=None, Lpq=None, orbs=None, nw=100, vhf_df=False)[source]#
- Input:
orbs: self-energy orbs nw: grid number
- Output:
mo_energy: GW quasiparticle energy
- linearized = False#
- property nmo#
- property nocc#
- pyscf.gw.gw_cd.get_WmnI_diag(gw, orbs, Lpq, freqs)[source]#
Compute W_mn(iw) on imarginary axis grids Return:
Wmn: (Nmo, Norbs, Nw)
- pyscf.gw.gw_cd.get_rho_response(omega, mo_energy, Lpq)[source]#
Compute density response function in auxiliary basis at freq iw
- pyscf.gw.gw_cd.get_rho_response_R(gw, omega, Lpq)[source]#
Compute density response function in auxiliary basis at poles
- pyscf.gw.gw_cd.get_sigmaI_diag(gw, omega, Wmn, sign, freqs, wts)[source]#
Compute self-energy by integrating on imaginary axis
- pyscf.gw.gw_cd.get_sigmaR_diag(gw, omega, orbp, ef, Lpq)[source]#
Compute self-energy for poles inside contour (more and more expensive away from Fermi surface)
pyscf.gw.gw_exact module#
Spin-restricted G0W0 approximation with exact frequency integration
- class pyscf.gw.gw_exact.GWExact(mf, frozen=None, tdmf=None)[source]#
Bases:
StreamObject
non-relativistic restricted GW
Saved results
- mo_energy :
Orbital energies
- mo_coeff
Orbital coefficients
- eta = 1e-08#
- get_frozen_mask()#
Get boolean mask for the restricted reference orbitals.
In the returned boolean (mask) array of frozen orbital indices, the element is False if it corresponds to the frozen orbital.
- get_nmo()#
- get_nocc()#
- kernel(mo_energy=None, mo_coeff=None, td_e=None, td_xy=None, eris=None, orbs=None)[source]#
Kernel function is the main driver of a method. Every method should define the kernel function as the entry of the calculation. Note the return value of kernel function is not strictly defined. It can be anything related to the method (such as the energy, the wave-function, the DFT mesh grids etc.).
- linearized = False#
- property nmo#
- property nocc#
pyscf.gw.gw_slow module#
This module implements the G0W0 approximation on top of pyscf.tdscf.rhf_slow and pyscf.tdscf.proxy TD implementations. Unlike gw.py, all integrals are stored in memory. Several variants of GW are available:
(this module) pyscf.gw_slow: the molecular implementation;
pyscf.pbc.gw.gw_slow: single-kpoint PBC (periodic boundary condition) implementation;
pyscf.pbc.gw.kgw_slow_supercell: a supercell approach to PBC implementation with multiple k-points. Runs the molecular code for a model with several k-points for the cost of discarding momentum conservation and using dense instead of sparse matrixes;
pyscf.pbc.gw.kgw_slow: a PBC implementation with multiple k-points;
- class pyscf.gw.gw_slow.AbstractIMDS(td, eri=None)[source]#
Bases:
object
- get_rhs(p)[source]#
The right-hand side of the quasiparticle equation. Args:
p (int, tuple): the orbital;
- Returns:
Right-hand sides of the quasiparticle equation
- get_sigma_element(omega, p, **kwargs)[source]#
The diagonal matrix element of the self-energy matrix. Args:
omega (float): the energy value; p (int, tuple): the orbital;
- Returns:
The diagonal matrix element.
- initial_guess(p)[source]#
Retrieves the initial guess for the quasiparticle energy for orbital p. Args:
p (int, tuple): the orbital;
- Returns:
The value of initial guess (float).
- orb_dims = 1#
- class pyscf.gw.gw_slow.IMDS(td, eri=None)[source]#
Bases:
AbstractIMDS
- property entire_space#
The entire orbital space. Returns:
An iterable of the entire orbital space.
- get_rhs(p)[source]#
The right-hand side of the quasiparticle equation. Args:
p (int, tuple): the orbital;
- Returns:
Right-hand sides of the quasiparticle equation
- pyscf.gw.gw_slow.corrected_moe(eri, p)[source]#
Calculates the corrected orbital energy. Args:
eri (PhysERI): a container with electron repulsion integrals; p (int): orbital;
- Returns:
The corrected orbital energy.
- pyscf.gw.gw_slow.kernel(imds, orbs=None, linearized=False, eta=0.001, tol=1e-09, method='fallback')[source]#
Calculates GW energies. Args:
imds (AbstractIMDS): GW intermediates; orbs (Iterable): indexes of MO orbitals to correct; linearized (bool): whether to apply a single-step linearized correction to energies instead of iterative procedure; eta (float): imaginary energy for the Green’s function; tol (float): tolerance for the search of zero; method (str): ‘bisect’ finds roots no matter what but, potentially, wrong ones, ‘newton’ finding roots close to the correct one but, potentially, failing during iterations, or ‘fallback’ using ‘newton’ and proceeding to ‘bisect’ in case of failure;
- Returns:
Corrected orbital energies.
pyscf.gw.rpa module#
Spin-restricted random phase approximation (direct RPA/dRPA in chemistry) with N^4 scaling
- Method:
Main routines are based on GW-AC method described in: T. Zhu and G.K.-L. Chan, J. Chem. Theory. Comput. 17, 727-741 (2021) X. Ren et al., New J. Phys. 14, 053020 (2012)
- class pyscf.gw.rpa.DirectRPA(mf, frozen=None, auxbasis=None)[source]#
Bases:
StreamObject
- get_frozen_mask()#
Get boolean mask for the restricted reference orbitals.
In the returned boolean (mask) array of frozen orbital indices, the element is False if it corresponds to the frozen orbital.
- get_nmo()#
- get_nocc()#
- kernel(mo_energy=None, mo_coeff=None, cderi_ov=None, nw=40, x0=0.5)[source]#
The kernel function for direct RPA
- make_dielectric_matrix(omega, e_ov=None, cderi_ov=None, blksize=None)[source]#
- Args:
omega : float, frequency e_ov : 1D array (nocc * nvir), orbital energy differences mo_coeff : (nao, nmo), mean-field mo coefficient cderi_ov : (naux, nocc, nvir), Cholesky decomposed ERI in OV subspace.
- Returns:
diel : 2D array (naux, naux), dielectric matrix
- property nmo#
- property nocc#
- pyscf.gw.rpa.kernel(rpa, mo_energy, mo_coeff, cderi_ov=None, nw=40, x0=0.5, verbose=3)[source]#
RPA correlation and total energy
- Args:
- cderi_ov:
Array-like object, Cholesky decomposed ERI in OV subspace.
- nw:
number of frequency point on imaginary axis.
- x0:
scaling factor for frequency grid.
- Returns:
- e_tot:
RPA total energy
- e_hf:
EXX energy
- e_corr:
RPA correlation energy
- pyscf.gw.rpa.make_dielectric_matrix(omega, e_ov, cderi_ov, blksize=None)[source]#
Compute dielectric matrix at a given frequency omega
- Args:
omega : float, frequency e_ov : 1D array (nocc * nvir), orbital energy differences cderi_ov : 2D array (naux, nocc * nvir), Cholesky decomposed ERI
in OV subspace.
- Returns:
diel : 2D array (naux, naux), dielectric matrix
pyscf.gw.ugw_ac module#
Spin-unrestricted G0W0 approximation with analytic continuation
This implementation has N^4 scaling, and is faster than GW-CD (N^4) and analytic GW (N^6) methods. GW-AC is recommended for valence states only, and is inaccurate for core states.
- Method:
See T. Zhu and G.K.-L. Chan, arxiv:2007.03148 (2020) for details Compute Sigma on imaginary frequency with density fitting, then analytically continued to real frequency
- Useful References:
J. Chem. Theory Comput. 12, 3623-3635 (2016) New J. Phys. 14, 053020 (2012)
- pyscf.gw.ugw_ac.AC_pade_thiele_diag(sigma, omega)[source]#
Analytic continuation to real axis using a Pade approximation from Thiele’s reciprocal difference method Reference: J. Low Temp. Phys. 29, 179 (1977) Returns:
coeff: 2D array (ncoeff, norbs) omega: 2D array (norbs, npade)
- pyscf.gw.ugw_ac.AC_twopole_diag(sigma, omega, orbs, nocc)[source]#
Analytic continuation to real axis using a two-pole model Returns:
coeff: 2D array (ncoeff, norbs)
- class pyscf.gw.ugw_ac.UGWAC(mf, frozen=None)[source]#
Bases:
StreamObject
- ac = 'pade'#
- get_frozen_mask()#
Get boolean mask for the unrestricted reference orbitals.
In the returned boolean (mask) array of frozen orbital indices, the element is False if it corresponds to the frozen orbital.
- get_nmo()#
- get_nocc()#
- kernel(mo_energy=None, mo_coeff=None, Lpq=None, orbs=None, nw=100, vhf_df=False)[source]#
- Input:
orbs: self-energy orbs nw: grid number vhf_df: whether using density fitting for HF exchange
- Output:
mo_energy: GW quasiparticle energy
- linearized = False#
- property nmo#
- property nocc#
- pyscf.gw.ugw_ac.get_rho_response(omega, mo_energy, Lpqa, Lpqb)[source]#
Compute density response function in auxiliary basis at freq iw
- pyscf.gw.ugw_ac.get_sigma_diag(gw, orbs, Lpq, freqs, wts, iw_cutoff=None)[source]#
Compute GW correlation self-energy (diagonal elements) in MO basis on imaginary axis
pyscf.gw.urpa module#
Spin-unrestricted random phase approximation (direct RPA/dRPA in chemistry) with N^4 scaling
- Method:
Main routines are based on GW-AC method described in: T. Zhu and G.K.-L. Chan, J. Chem. Theory. Comput. 17, 727-741 (2021) X. Ren et al., New J. Phys. 14, 053020 (2012)
- class pyscf.gw.urpa.URPA(mf, frozen=None, auxbasis=None)[source]#
Bases:
DirectRPA
- get_frozen_mask()#
Get boolean mask for the unrestricted reference orbitals.
In the returned boolean (mask) array of frozen orbital indices, the element is False if it corresponds to the frozen orbital.
- get_nmo()#
- get_nocc()#
- make_dielectric_matrix(omega, e_ov=None, cderi_ov=None, blksize=None)[source]#
- Args:
omega : float, frequency mo_energy : (2, nmo), mean-field mo energy mo_coeff : (2, nao, nmo), mean-field mo coefficient cderi_ov : (2, naux, nocc, nvir), Cholesky decomposed ERI in OV subspace.
- Returns:
diel : 2D array (naux, naux), dielectric matrix
Module contents#
G0W0 approximation