pyscf.pbc.gw package

Submodules

pyscf.pbc.gw.gw_slow module

This module implements the G0W0 approximation on top of pyscf.tdscf.rhf_slow TDHF implementation. Unlike gw.py, all integrals are stored in memory. Several variants of GW are available:

  • pyscf.gw_slow: the molecular implementation;

  • (this module) 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;

pyscf.pbc.gw.kgw_slow module

This module implements the G0W0 approximation on top of pyscf.tdscf.rhf_slow TDHF implementation. Unlike gw.py, all integrals are stored in memory. Several variants of GW are available:

  • 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;

  • (this module) pyscf.pbc.gw.kgw_slow: a PBC implementation with multiple k-points;

class pyscf.pbc.gw.kgw_slow.GW(td, eri=None)[source]

Bases: GW

base_imds

alias of IMDS

class pyscf.pbc.gw.kgw_slow.IMDS(td, eri=None)[source]

Bases: IMDS

construct_tdm()[source]
eri_ov(item)[source]
get_sigma_element(omega, p, eta, vir_sgn=1)[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.

pyscf.pbc.gw.kgw_slow_supercell module

This module implements the G0W0 approximation on top of pyscf.tdscf.rhf_slow TDHF implementation. Unlike gw.py, all integrals are stored in memory. Several variants of GW are available:

  • pyscf.gw_slow: the molecular implementation;

  • pyscf.pbc.gw.gw_slow: single-kpoint PBC (periodic boundary condition) implementation;

  • (this module) 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.pbc.gw.kgw_slow_supercell.GW(td, eri=None)[source]

Bases: GW

base_imds

alias of IMDS

class pyscf.pbc.gw.kgw_slow_supercell.IMDS(td, eri=None)[source]

Bases: IMDS

property entire_space

The entire orbital space. Returns:

An iterable of the entire orbital space.

eri_ov(item)[source]
get_rhs(p, components=False)[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, eta, vir_sgn=1)[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 = 2
pyscf.pbc.gw.kgw_slow_supercell.corrected_moe(eri, k, p)[source]

Calculates the corrected orbital energy. Args:

eri (PhysERI): a container with electron repulsion integrals; k (int): the k-point index; p (int): orbital;

Returns:

The corrected orbital energy.

pyscf.pbc.gw.krgw_ac module

PBC spin-restricted G0W0-AC QP eigenvalues with k-point sampling 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 inaccuarate 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. Gaussian density fitting must be used (FFTDF and MDF are not supported).

pyscf.pbc.gw.krgw_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.pbc.gw.krgw_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.pbc.gw.krgw_ac.KRGWAC(mf, frozen=0)[source]

Bases: StreamObject

ac = 'pade'
dump_flags()[source]
fc = True
get_frozen_mask()

Boolean mask for orbitals in k-point post-HF method.

Creates a boolean mask to remove frozen orbitals and keep other orbitals for post-HF calculations.

Args:

mp (MP2): An instantiation of an SCF or post-Hartree-Fock object.

Returns:

moidx (list of ndarray of bool): Boolean mask of orbitals to include.

get_nmo(per_kpoint=False)

Number of orbitals for k-point calculations.

Number of orbitals for use in a calculation with k-points, taking into account frozen orbitals.

Note:

If per_kpoint is False, then the number of orbitals here is equal to max(nocc) + max(nvir), where each max is done over all k-points. Otherwise the number of orbitals is returned as a list of number of orbitals at each k-point.

Args:

mp (MP2): An instantiation of an SCF or post-Hartree-Fock object. per_kpoint (bool, optional): True returns the number of orbitals at each k-point.

For a description of False, see Note.

Returns:
nmo (int, list of int): Number of orbitals. For return type, see description of arg

per_kpoint.

get_nocc(per_kpoint=False)

Number of occupied orbitals for k-point calculations.

Number of occupied orbitals for use in a calculation with k-points, taking into account frozen orbitals.

Args:

mp (MP2): An instantiation of an SCF or post-Hartree-Fock object. per_kpoint (bool, optional): True returns the number of occupied

orbitals at each k-point. False gives the max of this list.

Returns:
nocc (int, list of int): Number of occupied orbitals. For return type, see description of arg

per_kpoint.

kernel(mo_energy=None, mo_coeff=None, orbs=None, kptlist=None, nw=100)[source]
Input:

kptlist: self-energy k-points orbs: self-energy orbs nw: grid number

Output:

mo_energy: GW quasiparticle energy

linearized = False
property nmo
property nocc
pyscf.pbc.gw.krgw_ac.get_qij(gw, q, mo_coeff, uniform_grids=False)[source]

Compute qij = 1/Omega * |< psi_{ik} | e^{iqr} | psi_{ak-q} >|^2 at q: (nkpts, nocc, nvir) through kp perturbtation theory Ref: Phys. Rev. B 83, 245122 (2011)

pyscf.pbc.gw.krgw_ac.get_rho_response(gw, omega, mo_energy, Lpq, kL, kidx)[source]

Compute density response function in auxiliary basis at freq iw

pyscf.pbc.gw.krgw_ac.get_rho_response_head(gw, omega, mo_energy, qij)[source]

Compute head (G=0, G’=0) density response function in auxiliary basis at freq iw

pyscf.pbc.gw.krgw_ac.get_rho_response_wing(gw, omega, mo_energy, Lpq, qij)[source]

Compute wing (G=P, G’=0) density response function in auxiliary basis at freq iw

pyscf.pbc.gw.krgw_ac.get_sigma_diag(gw, orbs, kptlist, freqs, wts, iw_cutoff=None, max_memory=8000)[source]

Compute GW correlation self-energy (diagonal elements) in MO basis on imaginary axis

pyscf.pbc.gw.krgw_ac.kernel(gw, mo_energy, mo_coeff, orbs=None, kptlist=None, nw=None, verbose=3)[source]

GW-corrected quasiparticle orbital energies

Returns:

A list : converged, mo_energy, mo_coeff

pyscf.pbc.gw.krgw_ac.pade_thiele(freqs, zn, coeff)[source]
pyscf.pbc.gw.krgw_ac.thiele(fn, zn)[source]
pyscf.pbc.gw.krgw_ac.two_pole(freqs, coeff)[source]
pyscf.pbc.gw.krgw_ac.two_pole_fit(coeff, omega, sigma)[source]

pyscf.pbc.gw.krgw_cd module

PBC spin-restricted G0W0-CD QP eigenvalues with k-point sampling 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

class pyscf.pbc.gw.krgw_cd.KRGWCD(mf, frozen=0)[source]

Bases: StreamObject

dump_flags()[source]
eta = 0.001
fc = True
get_frozen_mask()

Boolean mask for orbitals in k-point post-HF method.

Creates a boolean mask to remove frozen orbitals and keep other orbitals for post-HF calculations.

Args:

mp (MP2): An instantiation of an SCF or post-Hartree-Fock object.

Returns:

moidx (list of ndarray of bool): Boolean mask of orbitals to include.

get_nmo(per_kpoint=False)

Number of orbitals for k-point calculations.

Number of orbitals for use in a calculation with k-points, taking into account frozen orbitals.

Note:

If per_kpoint is False, then the number of orbitals here is equal to max(nocc) + max(nvir), where each max is done over all k-points. Otherwise the number of orbitals is returned as a list of number of orbitals at each k-point.

Args:

mp (MP2): An instantiation of an SCF or post-Hartree-Fock object. per_kpoint (bool, optional): True returns the number of orbitals at each k-point.

For a description of False, see Note.

Returns:
nmo (int, list of int): Number of orbitals. For return type, see description of arg

per_kpoint.

get_nocc(per_kpoint=False)

Number of occupied orbitals for k-point calculations.

Number of occupied orbitals for use in a calculation with k-points, taking into account frozen orbitals.

Args:

mp (MP2): An instantiation of an SCF or post-Hartree-Fock object. per_kpoint (bool, optional): True returns the number of occupied

orbitals at each k-point. False gives the max of this list.

Returns:
nocc (int, list of int): Number of occupied orbitals. For return type, see description of arg

per_kpoint.

kernel(mo_energy=None, mo_coeff=None, orbs=None, kptlist=None, nw=100)[source]
Input:

kptlist: self-energy k-points orbs: self-energy orbs nw: grid number

Output:

mo_energy: GW quasiparticle energy

linearized = False
property nmo
property nocc
pyscf.pbc.gw.krgw_cd.get_WmnI_diag(gw, orbs, kptlist, freqs, max_memory=8000)[source]

Compute GW correlation self-energy (diagonal elements) in MO basis on imaginary axis

pyscf.pbc.gw.krgw_cd.get_qij(gw, q, mo_coeff, uniform_grids=False)[source]

Compute qij = 1/Omega * |< psi_{ik} | e^{iqr} | psi_{ak-q} >|^2 at q: (nkpts, nocc, nvir)

pyscf.pbc.gw.krgw_cd.get_rho_response(gw, omega, mo_energy, Lpq, kL, kidx)[source]

Compute density response function in auxiliary basis at freq iw

pyscf.pbc.gw.krgw_cd.get_rho_response_R(gw, omega, mo_energy, Lpq, kL, kidx)[source]

Compute density response function in auxiliary basis at freq iw

pyscf.pbc.gw.krgw_cd.get_rho_response_head(gw, omega, mo_energy, qij)[source]

Compute head (G=0, G’=0) density response function in auxiliary basis at freq iw

pyscf.pbc.gw.krgw_cd.get_rho_response_head_R(gw, omega, mo_energy, qij)[source]

Compute head (G=0, G’=0) density response function in auxiliary basis at freq w

pyscf.pbc.gw.krgw_cd.get_rho_response_wing(gw, omega, mo_energy, Lpq, qij)[source]

Compute wing (G=P, G’=0) density response function in auxiliary basis at freq iw

pyscf.pbc.gw.krgw_cd.get_rho_response_wing_R(gw, omega, mo_energy, Lpq, qij)[source]

Compute density response function in auxiliary basis at freq iw

pyscf.pbc.gw.krgw_cd.get_sigmaI_diag(gw, omega, kp, p, Wmn, Del_00, Del_P0, sign, freqs, wts)[source]

Compute self-energy by integrating on imaginary axis

pyscf.pbc.gw.krgw_cd.get_sigmaR_diag(gw, omega, kn, orbp, ef, freqs, qij, q_abs)[source]

Compute self-energy for poles inside coutour (more and more expensive away from Fermi surface)

pyscf.pbc.gw.krgw_cd.get_sigma_diag(gw, ep, kp, p, Wmn, Del_00, Del_P0, freqs, wts, qij, q_abs)[source]

Compute self-energy on real axis using contour deformation

pyscf.pbc.gw.krgw_cd.kernel(gw, mo_energy, mo_coeff, orbs=None, kptlist=None, nw=None, verbose=3)[source]

GW-corrected quasiparticle orbital energies

Returns:

A list : converged, mo_energy, mo_coeff

pyscf.pbc.gw.kugw_ac module

PBC spin-unrestricted G0W0-AC QP eigenvalues with k-point sampling GW-AC is recommended for valence states only, and is inaccuarate 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. Gaussian density fitting must be used (FFTDF and MDF are not supported).

pyscf.pbc.gw.kugw_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.pbc.gw.kugw_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.pbc.gw.kugw_ac.KUGWAC(mf, frozen=0)[source]

Bases: StreamObject

ac = 'pade'
dump_flags()[source]
fc = True
get_frozen_mask()

Boolean mask for orbitals in k-point post-HF method.

Creates a boolean mask to remove frozen orbitals and keep other orbitals for post-HF calculations.

Args:

mp (MP2): An instantiation of an SCF or post-Hartree-Fock object.

Returns:

moidx (list of ndarray of bool): Boolean mask of orbitals to include.

get_nmo(per_kpoint=False)

Number of orbitals for k-point calculations.

Number of orbitals for use in a calculation with k-points, taking into account frozen orbitals.

Note:

If per_kpoint is False, then the number of orbitals here is equal to max(nocc) + max(nvir), where each max is done over all k-points. Otherwise the number of orbitals is returned as a list of number of orbitals at each k-point.

Args:

mp (MP2): An instantiation of an SCF or post-Hartree-Fock object. per_kpoint (bool, optional): True returns the number of orbitals at each k-point.

For a description of False, see Note.

Returns:
nmo (int, list of int): Number of orbitals. For return type, see description of arg

per_kpoint.

get_nocc(per_kpoint=False)

Number of occupied orbitals for k-point calculations.

Number of occupied orbitals for use in a calculation with k-points, taking into account frozen orbitals.

Args:

mp (MP2): An instantiation of an SCF or post-Hartree-Fock object. per_kpoint (bool, optional): True returns the number of occupied

orbitals at each k-point. False gives the max of this list.

Returns:
nocc (int, list of int): Number of occupied orbitals. For return type, see description of arg

per_kpoint.

Notes:

For specifying frozen orbitals inside mp, the following options are accepted:

+=========================+========================================+===============================+ | Argument (Example) | Argument Meaning | Example Meaning | +=========================+========================================+===============================+ | int (1) | Freeze the same number of orbitals | Freeze one (lowest) orbital | | | regardless of spin and/or kpt | for all kpts and spin cases | +————————-+—————————————-+——————————-+ | 2-tuple of list of int | inner list: List of orbitals indices | Freeze the orbitals [0,4] for | | ([0, 4], [0, 5, 6]) | to freeze at all kpts | spin0, and orbitals [0,5,6] | | | outer list: Spin index | for spin1 at all kpts | +————————-+—————————————-+——————————-+ | list(2) of list of list | inner list: list of orbital indices to | Freeze orbital 0 for spin0 at | | ([[0,],[]], | freeze at each kpt for given spin | kpt0, and freeze orbital 0,1 | | [[0,1],[4]]) | outer list: spin index | for spin1 at kpt0 and orbital | | | | 4 at kpt1 | +————————-+—————————————-+——————————-+

kernel(mo_energy=None, mo_coeff=None, orbs=None, kptlist=None, nw=100)[source]
Input:

kptlist: self-energy k-points orbs: self-energy orbs nw: grid number

Output:

mo_energy: GW quasiparticle energy

linearized = False
property nmo
property nocc
pyscf.pbc.gw.kugw_ac.get_qij(gw, q, mo_coeff, uniform_grids=False)[source]

Compute qij = 1/Omega * |< psi_{ik} | e^{iqr} | psi_{ak-q} >|^2 at q: (nkpts, nocc, nvir) through kp perturbtation theory Ref: Phys. Rev. B 83, 245122 (2011)

pyscf.pbc.gw.kugw_ac.get_rho_response(gw, omega, mo_energy, Lpq, kL, kidx)[source]

Compute density response function in auxiliary basis at freq iw

pyscf.pbc.gw.kugw_ac.get_rho_response_head(gw, omega, mo_energy, qij)[source]

Compute head (G=0, G’=0) density response function in auxiliary basis at freq iw

pyscf.pbc.gw.kugw_ac.get_rho_response_wing(gw, omega, mo_energy, Lpq, qij)[source]

Compute wing (G=P, G’=0) density response function in auxiliary basis at freq iw

pyscf.pbc.gw.kugw_ac.get_sigma_diag(gw, orbs, kptlist, freqs, wts, iw_cutoff=None, max_memory=8000)[source]

Compute GW correlation self-energy (diagonal elements) in MO basis on imaginary axis

pyscf.pbc.gw.kugw_ac.kernel(gw, mo_energy, mo_coeff, orbs=None, kptlist=None, nw=None, verbose=3)[source]

GW-corrected quasiparticle orbital energies Returns:

A list : converged, mo_energy, mo_coeff

pyscf.pbc.gw.kugw_ac.pade_thiele(freqs, zn, coeff)[source]
pyscf.pbc.gw.kugw_ac.thiele(fn, zn)[source]
pyscf.pbc.gw.kugw_ac.two_pole(freqs, coeff)[source]
pyscf.pbc.gw.kugw_ac.two_pole_fit(coeff, omega, sigma)[source]

Module contents