momentGW.pbc.rpa
Construct RPA moments with periodic boundary conditions.
Module Contents
- class momentGW.pbc.rpa.dRPA(gw, nmom_max, integrals, mo_energy=None, mo_occ=None)
Bases:
momentGW.pbc.tda.dTDA
,momentGW.rpa.dRPA
Compute the self-energy moments using dRPA and numerical integration with periodic boundary conditions.
- Parameters:
gw (BaseKGW) – GW object.
nmom_max (int) – Maximum moment number to calculate.
integrals (KIntegrals) – Density-fitted integrals at each k-point.
mo_energy (dict, optional) – Molecular orbital energies. Keys are “g” and “w” for the Green’s function and screened Coulomb interaction, respectively. If None, use gw.mo_energy for both. Default value is None.
mo_occ (dict, optional) – Molecular orbital occupancies. Keys are “g” and “w” for the Green’s function and screened Coulomb interaction, respectively. If None, use gw.mo_occ for both. Default value is None.
Notes
See momentGW.tda.dTDA.__init__ for initialisation details and momentGW.tda.dTDA.kernel for calculation run details.
- property nov
Get the number of ov states in W.
- property kpts
Get the k-points.
- property nkpts
Get the number of k-points.
- property nmo
Get the number of MOs.
- property naux
Get the number of auxiliaries.
- integrate()
Optimise the quadrature and perform the integration for a given set of k-points for the zeroth moment.
- Returns:
integral – Integral array, including the offset part.
- Return type:
numpy.ndarray
- build_dd_moments(integral=None)
Build the moments of the density-density response.
- Parameters:
integral (numpy.ndarray, optional) – Integral array, including the offset part. If None, calculate from scratch. Default is None.
- Returns:
moments – Moments of the density-density response.
- Return type:
numpy.ndarray
- optimise_offset_quad(d, diag_eri, name='main')
Optimise the grid spacing of Gauss-Laguerre quadrature for the offset integral.
- Parameters:
d (numpy.ndarray) – Orbital energy differences at each k-point.
diag_eri (numpy.ndarray) – Diagonal of the ERIs at each k-point.
name (str, optional) – Name of the integral. Default value is “main”.
- Returns:
points (numpy.ndarray) – The quadrature points.
weights (numpy.ndarray) – The quadrature weights.
- eval_diag_offset_integral(quad, d, diag_eri)
Evaluate the diagonal of the offset integral.
- Parameters:
quad (tuple) – The quadrature points and weights.
d (numpy.ndarray) – Orbital energy differences at each k-point.
diag_eri (numpy.ndarray) – Diagonal of the ERIs at each k-point.
- Returns:
integral – Offset integral.
- Return type:
numpy.ndarray
- eval_offset_integral(quad, d, Lia=None)
Evaluate the offset integral.
- Parameters:
quad (tuple) – The quadrature points and weights.
d (numpy.ndarray) – Orbital energy differences at each k-point.
Lia (dict of numpy.ndarray) – Dict. with keys that are pairs of k-point indices (Nkpt, Nkpt) with an array of form (aux, W occ, W vir) at this k-point pair. The 1st Nkpt is defined by the difference between k-points and the second index’s kpoint. If None, use self.integrals.Lia.
Liad (dict of numpy.ndarray) – Product of Lia and the orbital energy differences at each k-point.
- Returns:
integral – Offset integral.
- Return type:
numpy.ndarray
- optimise_main_quad(d, diag_eri, name='main')
Optimise the grid spacing of Clenshaw-Curtis quadrature for the main integral.
- Parameters:
d (numpy.ndarray) – Orbital energy differences at each k-point.
diag_eri (numpy.ndarray) – Diagonal of the ERIs at each k-point.
- Returns:
points (numpy.ndarray) – The quadrature points.
weights (numpy.ndarray) – The quadrature weights.
- eval_diag_main_integral(quad, d, d_sq, d_eri, d_sq_eri)
Evaluate the diagonal of the main integral.
- Parameters:
quad (tuple) – The quadrature points and weights.
d (numpy.ndarray) – Orbital energy differences at each k-point.
d_sq (numpy.ndarray) – Orbital energy differences squared at each k-point. See “optimise_main_quad” for more details.
d_eri (numpy.ndarray) – Orbital energy differences times the diagonal of the ERIs at each k-point. See “optimise_main_quad” for more details.
d_sq_eri (numpy.ndarray) – Orbital energy differences times the diagonal of the ERIs plus the orbital energy differences at each k-point. See “optimise_main_quad” for more details.
- Returns:
integral – Main integral.
- Return type:
numpy.ndarray
- eval_main_integral(quad, d, Lia=None)
Evaluate the main integral.
- Parameters:
quad (tuple) – The quadrature points and weights.
Variables –
---------- –
d (numpy.ndarray) – Orbital energy differences at each k-point.
Lia (dict of numpy.ndarray) – Dict. with keys that are pairs of k-point indices (Nkpt, Nkpt) with an array of form (aux, W occ, W vir) at this k-point pair. The 1st Nkpt is defined by the difference between k-points and the second index’s kpoint. If None, use self.integrals.Lia.
Liad (dict of numpy.ndarray) – Product of Lia and the orbital energy differences at each k-point.
- Returns:
integral – Offset integral.
- Return type:
numpy.ndarray
- kernel(exact=False)
Run the polarizability calculation to compute moments of the self-energy.
- Parameters:
exact (bool, optional) – Has no effect and is only present for compatibility with dRPA. Default value is False.
- Returns:
moments_occ (numpy.ndarray) – Moments of the occupied self-energy at each k-point.
moments_vir (numpy.ndarray) – Moments of the virtual self-energy at each k-point.
- convolve(eta, mo_energy_g=None, mo_occ_g=None)
Handle the convolution of the moments of the Green’s function and screened Coulomb interaction.
- Parameters:
eta (numpy.ndarray) – Moments of the density-density response partly transformed into moments of the screened Coulomb interaction at each k-point.
mo_energy_g (numpy.ndarray, optional) – Energies of the Green’s function at each k-point. If None, use self.mo_energy_g. Default value is None.
mo_occ_g (numpy.ndarray, optional) – Occupancies of the Green’s function at each k-point. If None, use self.mo_occ_g. Default value is None.
- Returns:
moments_occ (numpy.ndarray) – Moments of the occupied self-energy at each k-point.
moments_vir (numpy.ndarray) – Moments of the virtual self-energy at each k-point.
- build_se_moments(moments_dd)
Build the moments of the self-energy via convolution.
- Parameters:
moments_dd (numpy.ndarray) – Moments of the density-density response at each k-point.
- Returns:
moments_occ (numpy.ndarray) – Moments of the occupied self-energy at each k-point.
moments_vir (numpy.ndarray) – Moments of the virtual self-energy at each k-point.
- build_dp_moments()
Build the moments of the dynamic polarizability for optical spectra calculations.
- Returns:
moments – Moments of the dynamic polarizability.
- Return type:
numpy.ndarray
- build_dd_moment_inv()
Build the first inverse (n=-1) moment of the density-density response.
- Returns:
moment – First inverse (n=-1) moment of the density-density response.
- Return type:
numpy.ndarray
Notes
This is not the full n=-1 moment, which is
\[\begin{split}D^{-1} - D^{-1} V^\dagger (I + V D^{-1} V^\dagger)^{-1} \\ V D^{-1}\end{split}\]but rather
\[(I + V D^{-1} V^\dagger)^{-1} V D^{-1}\]which ensures that the function scales properly. The final contractions are done when constructing the matrix-vector product.
- mpi_slice(n)
Return the start and end index for the current process for total size n.
- Parameters:
n (int) – Total size.
- Returns:
p0 (int) – Start index for current process.
p1 (int) – End index for current process.
- mpi_size(n)
Return the number of states in the current process for total size n.
- build_dd_moments_exact()
Build the exact moments of the density-density response.
- Returns:
moments – Moments of the density-density response.
- Return type:
numpy.ndarray
- static rescale_quad(bare_quad, a)
Rescale quadrature for grid space a.
- get_optimal_quad(bare_quad, integrand, exact, name=None)
Get the optimal quadrature.
- Parameters:
- Returns:
points (numpy.ndarray) – The quadrature points.
weights (numpy.ndarray) – The quadrature weights.
- gen_clencur_quad_inf(even=False)
Generate quadrature points and weights for Clenshaw-Curtis quadrature over an
(-inf, +inf)
.- Parameters:
even (bool, optional) – Whether to assume an even grid. Default is False.
- Returns:
points (numpy.ndarray) – Quadrature points.
weights (numpy.ndarray) – Quadrature weights.
- gen_gausslag_quad_semiinf()
Generate quadrature points and weights for Gauss-Laguerre quadrature over an
(0, +inf)
.- Returns:
points (numpy.ndarray) – Quadrature points.
weights (numpy.ndarray) – Quadrature weights.
- estimate_error_clencur(i4, i2, imag_tol=1e-10)
Estimate the quadrature error for Clenshaw-Curtis quadrature.
- Parameters:
i4 (numpy.ndarray) – Integral at one-quarter the number of points.
i2 (numpy.ndarray) – Integral at one-half the number of points.
imag_tol (float, optional) – Threshold to consider the imaginary part of a root to be zero. Default value is 1e-10.
- Returns:
error – Estimated error.
- Return type:
numpy.ndarray