fcdmft.gw.pbc.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 inaccurate for core states.

Method:

T. Zhu and G.K.-L. Chan, J. Chem. Theory. Comput. 17, 727-741 (2021) 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).

Note: MPI only works for GW code (DFT should run in serial)

class fcdmft.gw.pbc.krgw_ac.KRGWAC(mf, frozen=None)[source]#

Bases: StreamObject

Attributes:
nmo
nocc

Methods

__call__(*args, **kwargs)

Update the attributes of the current object.

apply(fn, *args, **kwargs)

Apply the fn to rest arguments: return fn(*args, **kwargs). The return value of method set is the object itself. This allows a series of functions/methods to be executed in pipe.

check_sanity()

Check input of class/object attributes, check whether a class method is overwritten.

copy()

Returns a shallow copy

kernel([orbs, kptlist])

Run a G0W0 calculation.

make_gf(omega, eta)

Get dynamical Green's function and self-energy.

make_rdm1([ao_repr])

Get GW density matrix from Green's function G(it=0).

post_kernel(envs)

A hook to be run after the main body of the kernel function.

pre_kernel(envs)

A hook to be run before the main body of kernel function is executed.

run(*args, **kwargs)

Call the kernel function of current object.

set(*args, **kwargs)

Update the attributes of the current object.

set_frozen_orbs()

Set .frozen attribute from frozen mask.

view(cls)

New view of object with the same attributes.

dump_flags

dump_flags(verbose=None)[source]#
kernel(orbs=None, kptlist=None)[source]#

Run a G0W0 calculation.

Parameters:
orbslist, optional

orbital list to calculate self-energy, by default None

kptlistlist, optional

k-point list to calculate self-energy, by default None

make_gf(omega, eta)#

Get dynamical Green’s function and self-energy.

Parameters:
gwKRGWAC

GW object, provides attributes: orbs, kptlist, ef, ac_coeff, omega_fit, vk, vxc, _scf.mo_energy

omegadouble or complex array

frequency grids

etadouble

broadening parameter

Returns:
gfcomplex ndarray

GW Green’s function

gf0complex ndarray

mean-field Green’s function

sigmacomplex ndarray

GW correlation self-energy

make_rdm1(ao_repr=False)#

Get GW density matrix from Green’s function G(it=0). G is from linear Dyson equation, which conserves particle number G = G0 + G0 Sigma G0 See equation 16 in 10.1021/acs.jctc.0c01264

Parameters:
gwKRGWAC

GW object, provides attributes: sigmaI, mol, _scf, freqs, wts, frozen, orbs, fc

ao_reprbool, optional

return density matrix in AO, by default False

Returns:
rdm1double ndarray

density matrix

property nmo#
property nocc#
set_frozen_orbs()#

Set .frozen attribute from frozen mask.

Parameters:
gwKRGWAC

unrestricted GW object

fcdmft.gw.pbc.krgw_ac.get_ef(kmf, mo_energy)[source]#

Get Fermi level. For gapped systems, Fermi level is computed as the average between HOMO and LUMO. For metallic systems, Fermi level is optmized according to mo_energy.

Parameters:
kmfpyscf.pbc.scf.rhf.RHF/pyscf.pbc.dft.rks.RKS

mean-field object, provides attributes: kpts, sigma, smearing_method

mo_energydouble array

orbital energy

Returns:
efdouble

Fermi level

fcdmft.gw.pbc.krgw_ac.get_g0_k(omega, mo_energy, eta)[source]#

Get non-interacting Green’s function.

Parameters:
omegadouble or complex ndarray

frequency grids

mo_energydouble ndarray

orbital energy

etadouble

broadening parameter

Returns:
gf0complex ndarray

non-interacting Green’s function

fcdmft.gw.pbc.krgw_ac.get_qij(gw, q, mo_energy, mo_coeff, uniform_grids=False)[source]#

Compute pair density matrix in the long-wavelength limit through kp perturbation theory qij = 1/Omega * |< psi_{ik} | e^{iqr} | psi_{ak-q} >|^2 equation 51 in 10.1021/acs.jctc.0c00704 Ref: Phys. Rev. B 83, 245122 (2011)

Parameters:
gwKRGWAC

gw object, provides attributes: nocc, nmo, kpts, mol

qdouble

q grid

mo_energydouble ndarray

orbital energy

mo_coeffcomplex ndarray

coefficient from AO to MO

uniform_gridsbool, optional

use uniform grids, by default False

Returns:
qijcomplex ndarray

pair density matrix in the long-wavelength limit

fcdmft.gw.pbc.krgw_ac.get_rho_response(omega, mo_energy, Lia, kidx)[source]#

Get Pi=PV. P is density-density response function. V is two-electron integral. See equation 24 in 10.1021/acs.jctc.0c00704

Parameters:
omegadouble

real position of imaginary frequency

mo_energydouble 2d array

orbital energy

Liacomplex 4d ndarray

occupied-virtual block of three-center density-fitting matrix in MO

kidxlist

momentum-conserved k-point list kj=kidx[ki]

Returns:
Picomplex ndarray

Pi in auxiliary basis at freq iw

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

Compute head (G=0, G’=0) density response function in auxiliary basis at freq iw. equation 48 in 10.1021/acs.jctc.0c00704

Parameters:
omegadouble

frequency point

mo_energydouble ndarray

orbital energy

qijcomplex ndarray

pair density matrix defined as equation 51 in 10.1021/acs.jctc.0c00704

Returns:
Pi_00complex

head response function

fcdmft.gw.pbc.krgw_ac.get_rho_response_metal(omega, mo_energy, mo_occ, Lpq, kidx)[source]#

Get Pi=PV for metallic systems. P is density-density response function. V is two-electron integral. See equation 24 in 10.1021/acs.jctc.0c00704

Parameters:
omegadouble

real position of imaginary frequency

mo_energydouble ndarray

orbital energy

mo_occdouble ndarray

occupation number

Lpqcomplex ndarray

three-center density-fitting matrix in MO

kidxlist

momentum-conserved k-point list kj=kidx[ki]

Returns:
Picomplex ndarray

Pi in auxiliary basis at freq iw

fcdmft.gw.pbc.krgw_ac.get_rho_response_wing(omega, mo_energy, Lia, qij)[source]#

Compute wing (G=P, G’=0) density response function in auxiliary basis at freq iw. equation 48 in 10.1021/acs.jctc.0c00704

Parameters:
omegadouble

frequency point

mo_energydouble 2d array

orbital energy

Liacomplex 4d array

occupied-virtual block of three-center density fitting matrix in MO

qijcomplex ndarray

pair density matrix defined as equation 51 in 10.1021/acs.jctc.0c00704

Returns:
Picomplex ndarray

wing response function

fcdmft.gw.pbc.krgw_ac.get_sigma(gw, freqs, wts, ef, mo_energy, orbs=None, kptlist=None, mo_coeff=None, mo_occ=None, iw_cutoff=None, fullsigma=False)[source]#

Get GW self-energy. See equation 27 in 10.1021/acs.jctc.0c00704

Parameters:
gwKRGWAC

GW objects, provides attributes: _scf, mol, frozen, nmo, nocc, kpts, nkpts, mo_coeff, mo_occ, fc, fc_grid, with_df

freqsdouble array

position of imaginary frequency

wtsdouble array

weight of frequency points

efdouble

Fermi level

mo_energydouble ndarray

non-frozen orbital energy

orbslist, optional

orbital index in non-frozen nmo to calculate self-energy, by default None

kptlistlist, optional

k-point index to calculate self-energy, by default None

mo_coeffcomplex ndarray, optional

coefficient from AO to non-frozen MO, by default None

mo_occdouble ndarray, optional

non-frozen occupation number, by default None

iw_cutoffcomplex, optional

imaginary grid cutoff for fitting, by default None

fullsigmabool, optional

calculate off-diagonal elements, by default False

Returns:
sigma: complex ndarray

self-energy on the imaginary axis

omega: complex ndarray

imaginary frequency grids of self-energy

fcdmft.gw.pbc.krgw_ac.kernel(gw)[source]#
fcdmft.gw.pbc.krgw_ac.make_gf(gw, omega, eta)[source]#

Get dynamical Green’s function and self-energy.

Parameters:
gwKRGWAC

GW object, provides attributes: orbs, kptlist, ef, ac_coeff, omega_fit, vk, vxc, _scf.mo_energy

omegadouble or complex array

frequency grids

etadouble

broadening parameter

Returns:
gfcomplex ndarray

GW Green’s function

gf0complex ndarray

mean-field Green’s function

sigmacomplex ndarray

GW correlation self-energy

fcdmft.gw.pbc.krgw_ac.make_rdm1_linear(gw, ao_repr=False)[source]#

Get GW density matrix from Green’s function G(it=0). G is from linear Dyson equation, which conserves particle number G = G0 + G0 Sigma G0 See equation 16 in 10.1021/acs.jctc.0c01264

Parameters:
gwKRGWAC

GW object, provides attributes: sigmaI, mol, _scf, freqs, wts, frozen, orbs, fc

ao_reprbool, optional

return density matrix in AO, by default False

Returns:
rdm1double ndarray

density matrix

fcdmft.gw.pbc.krgw_ac.set_frozen_orbs(gw)[source]#

Set .frozen attribute from frozen mask.

Parameters:
gwKRGWAC

unrestricted GW object