fcdmft.gw.mol.gw_space_time module#

class fcdmft.gw.mol.gw_space_time.GWSpaceTime(mf, auxbasis=None)[source]#

Bases: GWAC

Attributes:
nmo
nocc

Methods

__call__(*args, **kwargs)

Update the attributes of the current object.

ao2mo([mo_coeff])

Transform density-fitting integral from AO to MO.

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

energy_tot()

Get GW correlation and total energy using Galitskii-Migdal formula on the imaginary time.

get_ef([mo_energy])

Get Fermi level.

get_frozen_mask()

Get boolean mask for the restricted reference orbitals.

get_sigma_exchange(mo_coeff)

Get exchange self-energy (EXX).

initialize_df([auxbasis])

Initialize density fitting.

kernel([mo_energy, mo_coeff, Lpq])

Run one-shot space-time GW calculation.

loop_ao2mo([mo_coeff, ijslice])

Transform density-fitting integral from AO to MO by block.

make_gf(omega[, eta, mode])

Get G0W0 Green's function by AC fitting.

make_rdm1([mode, ao_repr])

Get GW density matrix from 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.

setup_evaluation_grid([fallback_freqs, ...])

Set up self-energy grid, aka freqs2, aka gw.freqs.

shift_sigma(w[, E_min_scale, E_max_scale])

NOTE: experimental function.

view(cls)

New view of object with the same attributes.

dump_flags

get_nmo

get_nocc

dump_flags(verbose=None)[source]#
energy_tot()[source]#

Get GW correlation and total energy using Galitskii-Migdal formula on the imaginary time. See equation.13 in doi.org/10.1103/PhysRevB.101.205145

Returns:
e_totdouble

GW total energy

e_hfdouble

HF total energy

e_cdouble

GW correlation energy

kernel(mo_energy=None, mo_coeff=None, Lpq=None)[source]#

Run one-shot space-time GW calculation. This function is of O(N^4) scaling. No sparsity is used. Minimax grids were taken from CP2K. Reference [1] 10.1021/acs.jctc.0c01282 [2] 10.1021/acs.jpclett.7b02740 [3] 10.1103/PhysRevB.94.165109 [4] 10.1103/PhysRevB.101.205145 [5] 10.1021/ct5001268

Parameters:
mo_energydouble 1d array, optional

orbital energy, by default None

mo_coeffdouble 2d array, optional

coefficient from AO to MO, by default None

Lpqdouble 3d array, optional

3-center DF matrix in AO space, by default None

make_rdm1(mode='linear', ao_repr=False)[source]#

Get GW density matrix from G(it=0). The integration for G is separated into two parts: the Hartree-Fock part density matrix is simply diagonal (2, 2 … 2, 0, 0 …), the correlation part is Gc(it=0) = int Gc(iw) dw, where the correlation GW function is Gc(iw) = G(iw) - GHF(iw) and GHF(iw) = 1.0 / (iw - ef + H_HF). See equation 2-8 in 10.1103/PhysRevB.98.155143.

Parameters:
modestr, optional

mode for Dyson equation, by default ‘linear’

ao_reprbool, optional

return density matrix in AO, by default False

Returns:
rdm1double 2d array

one-particle density matrix

shift_sigma(w, E_min_scale=1.0, E_max_scale=1.0)[source]#

NOTE: experimental function. Get sigma(iw) on new grids by sigma(iw new) = FT[sigma(it old)]. Accuracy can be tuned by E_min_scale and E_max_scale.

Parameters:
wdouble 1d array

_description_

E_min_scaledouble, optional

scale factor of E_min, by default 1.0

E_max_scaledouble, optional

scale factor of E_max, by default 1.0

Returns:
sigma_wcomplex 3d array

sigma(iw) on new grids

fcdmft.gw.mol.gw_space_time.get_ft_weight(E_min, E_max, t, w, num=200, return_error=False)[source]#

Get Fourier transform weights between time and frequency domain.

Parameters:
E_mindouble

minimum frequency range

E_maxdouble

maximum frequency range

tdouble 1d array

scaled time grid

wdouble 1d array

scaled frequency grid

numint, optional

number of points per magnitude in least square fitting, by default 200

return_errorbool, optional

return least square error for FT weights, by default False

Returns:
wt_t_to_w_cosdouble 2d array

cosine t to w weight, [nt,nw]

wt_t_to_w_sindouble 2d array

sine t to w weight, [nt,nw]

wt_w_to_t_cosdouble 2d array

cosine w to t weight, [nw,nt]

wt_w_to_t_sindouble 2d array

sine w to t weight, [nw,nt]

error_t_to_w_cosdouble

cosine t to w error

error_t_to_w_sindouble

sine t to w error

error_w_to_t_cosdouble

cosine w to t error

error_w_to_t_sindouble

sine w to t error

fcdmft.gw.mol.gw_space_time.get_g0_t(nocc, mo_energy, t, ef, mo_coeff=None)[source]#

Get non-interacting Green’s function in the imaginary time domain. [G_mo]_{pq} (t) = delta_{pq} Theta(F-p) exp((mo_energy[p]-mu)t) t>0

Theta(p-F)-exp((mo_energy[p]-mu)t) t<0

the prefactor i is not included. if mo_coeff is not None, G_mo is rotated to G_ao.

Parameters:
noccint

number of occupied orbitals

mo_energydouble 1d array

orbital energy

tdouble 1d array

time grid

efdouble

Fermi level

mo_coeffdouble 2d array, optional

coefficient from AO to MO, by default None

Returns:
g0_t_posdouble 3d array

G0 on the positive imaginary time in MO/AO space

g0_t_negdouble 3d array

G0 on the negative imaginary time in MO/AO space

fcdmft.gw.mol.gw_space_time.get_imag_sigma_t(g_t_pos, g_t_neg, wc_t, Lpq, mo_coeff=None, fullsigma=False, fast_diag=True)[source]#

Get the correlation self-energy in the imaginary time domain. this function is similar to equation (28) in 10.1021/acs.jctc.0c01282 but we get sigma_c in AO first, then transform to the MO space

Parameters:
g_t_posdouble 3d array

Green’s function on the positive imaginary time

g_t_negdouble 3d array

Green’s function on the negative imaginary time

wc_tdouble 3d array

Wc matrix in the auxiliary space

Lpqdouble 3d array

3-center density-fitting matrix in the same space of G

mo_coeffdouble 2d array, optional

coefficient from initial space to final space, by default None

fullsigmabool, optional

calculate off-diagonal self-energy, by default False

fast_diagbool, optional

double the memory cost to speed up diagonal-only calculation, by default True

Returns:
sigma_t_posdouble 2d or 3d array

self-energy on the positive imaginary time in MO

sigma_t_negdouble 2d or 3d array

self-energy on the negative imaginary time in MO

fcdmft.gw.mol.gw_space_time.get_minimax_grid(k, E_min, E_max)[source]#

Get minimax grid for time and frequency.

Parameters:
kint

number of time/frequency grids

E_mindouble

minimum frequency range

E_maxdouble

maximum frequency range

Returns:
tdouble 1d array

scaled time grid

t_weightdouble 1d array

scaled time grid weight

twdouble 1d array

scaled frequency grid

w_weightdouble 1d array

scaled frequency grid weight

fcdmft.gw.mol.gw_space_time.get_pi_t(g_t_pos, g_t_neg, Lpq)[source]#

Get Pi(t)=P(t)V, where the response function P(t)=G(t)G(-t). G and Lpq must be in the same space (AO, MO or LO).

Parameters:
g_t_posdouble 3d array

Green’s function on positive imaginary time

g_t_negdouble 3d array

Green’s function on negative imaginary time

Lpqdouble 3d array

three-center density-fitting matrix

Returns:
Pi_tdouble 3d array

Pi on positive imaginary time domain in auxiliary space

fcdmft.gw.mol.gw_space_time.get_wc_w(Pi_w)[source]#

Get the correlation part screened interaction W in auxiliary orbital space and imaginary frequency space. Wc(w)=[I-Pi(w)]^-1 - I

Parameters:
Pi_wdouble 3d array

Pi on positive imaginary time domain in auxiliary space

Returns:
wc_wdouble 3d array

Wc on positive imaginary time domain in auxiliary space

fcdmft.gw.mol.gw_space_time.kernel(gw, mo_energy, mo_coeff, Lpq)[source]#

Run one-shot space-time GW calculation.

Parameters:
gwGWSpaceTime

gw object

mo_energydouble 1d array

orbital energy

mo_coeffdouble 2d array

coefficient from AO to MO

Lpqdouble 3d array

3-center density-fitting matrix in AO space