momentGW.uhf.qsgw

Spin-unrestricted quasiparticle self-consistent GW via self-energy moment constraints for molecular systems.

Module Contents

class momentGW.uhf.qsgw.qsUGW(mf, **kwargs)

Bases: momentGW.uhf.UGW, momentGW.qsgw.qsGW

Spin-unrestricted quasiparticle self-consistent GW via self-energy moment constraints for molecules.

Parameters:
  • mf (pyscf.scf.SCF) – PySCF mean-field class.

  • diagonal_se (bool, optional) – If True, use a diagonal approximation in the self-energy. Default value is False.

  • polarizability (str, optional) – Type of polarizability to use, can be one of (“drpa”, “drpa-exact”, “dtda”, “thc-dtda”). Default value is `”drpa”.

  • npoints (int, optional) – Number of numerical integration points. Default value is 48.

  • optimise_chempot (bool, optional) – If True, optimise the chemical potential by shifting the position of the poles in the self-energy relative to those in the Green’s function. Default value is False.

  • fock_loop (bool, optional) – If True, self-consistently renormalise the density matrix according to the updated Green’s function. Default value is False.

  • fock_opts (dict, optional) – Dictionary of options passed to the Fock loop. For more details see momentGW.fock.

  • compression (str, optional) – Blocks of the ERIs to use as a metric for compression. Can be one or more of (“oo”, “ov”, “vv”, “ia”) which can be passed as a comma-separated string. “oo”, “ov” and “vv” refer to compression on the initial ERIs, whereas “ia” refers to compression on the ERIs entering RPA, which may change under a self-consistent scheme. Default value is “ia”.

  • compression_tol (float, optional) – Tolerance for the compression. Default value is 1e-10.

  • thc_opts (dict, optional) – Dictionary of options to be used for THC calculations. Current implementation requires a filepath to import the THC integrals.

  • max_cycle (int, optional) – Maximum number of iterations. Default value is 50.

  • max_cycle_qp (int, optional) – Maximum number of iterations in the quasiparticle equation loop. Default value is 50.

  • conv_tol (float, optional) – Convergence threshold in the change in the HOMO and LUMO. Default value is 1e-8.

  • conv_tol_moms (float, optional) – Convergence threshold in the change in the moments. Default value is 1e-8.

  • conv_tol_qp (float, optional) – Convergence threshold in the change in the density matrix in the quasiparticle equation loop. Default value is 1e-8.

  • conv_logical (callable, optional) – Function that takes an iterable of booleans as input indicating whether the individual conv_tol, conv_tol_moms, conv_tol_qp have been satisfied, respectively, and returns a boolean indicating overall convergence. For example, the function all requires both metrics to be met, and any requires just one. Default value is all.

  • diis_space (int, optional) – Size of the DIIS extrapolation space. Default value is 8.

  • diis_space_qp (int, optional) – Size of the DIIS extrapolation space in the quasiparticle loop. Default value is 8.

  • damping (float, optional) – Damping parameter. Default value is 0.0.

  • eta (float, optional) – Small value to regularise the self-energy. Default value is 1e-1.

  • srg (float, optional) – If non-zero, use the similarity renormalisation group approach of Marie and Loos in place of the eta regularisation. For value recommendations refer to their paper (arXiv:2303.05984). Default value is 0.0.

  • solver (BaseGW, optional) – Solver to use to obtain the self-energy. Compatible with any BaseGW-like class. Default value is momentGW.gw.GW.

  • solver_options (dict, optional) – Keyword arguments to pass to the solver. Default value is an empty dict.

property name

Get the method name.

property qp_energy

Get the quasiparticle energies.

Notes

For most GW methods, this simply consists of the poles of the self.gf that best overlap with the MOs, in order. In some methods such as qsGW, these two quantities are not the same.

property has_fock_loop

Get a boolean indicating whether the solver requires a Fock loop.

Notes

For most GW methods, this is simply self.fock_loop. In some methods such as qsGW, a Fock loop is required with or without self.fock_loop for the quasiparticle self-consistency, with this property acting as a hook to indicate this.

property mol

Get the molecule object.

property with_df

Get the density fitting object.

property nao

Get the number of atomic orbitals.

property nmo

Get the number of molecular orbitals.

property nocc

Get the number of occupied molecular orbitals.

property active

Get the mask to remove frozen orbitals.

property mo_energy

Get the molecular orbital energies.

property mo_energy_with_frozen

Get the molecular orbital energies with frozen orbitals.

property mo_coeff

Get the molecular orbital coefficients.

property mo_coeff_with_frozen

Get the molecular orbital coefficients with frozen orbitals.

property mo_occ

Get the molecular orbital occupation numbers.

property mo_occ_with_frozen

Get the molecular orbital occupation numbers with frozen orbitals.

static project_basis(matrix, ovlp, mo1, mo2)

Project a matrix from one basis to another.

Parameters:
  • matrix (numpy.ndarray or tuple of dyson.Lehmann) – Matrix to project for each spin channel. Can also be a tuple of dyson.Lehmann objects, in which case the couplings attributes are projected.

  • ovlp (numpy.ndarray) – Overlap matrix in the shared (AO) basis.

  • mo1 (numpy.ndarray) – First basis, rotates from the shared (AO) basis into the basis of matrix for each spin channel.

  • mo2 (numpy.ndarray) – Second basis, rotates from the shared (AO) basis into the desired basis of the output for each spin channel.

Returns:

proj – Matrix projected into the desired basis for each spin channel.

Return type:

numpy.ndarray or tuple of dyson.Lehmann

static self_energy_to_moments(se, nmom_max)

Return the hole and particle moments for a self-energy.

Parameters:
  • se (tuple of dyson.Lehmann) – Self-energy to compute the moments of for each spin channel.

  • nmom_max (int) – Maximum moment number to calculate.

Returns:

  • th (numpy.ndarray) – Hole moments for each spin channel.

  • tp (numpy.ndarray) – Particle moments for each spin channel.

build_static_potential(mo_energy, se)

Build the static potential approximation to the self-energy.

Parameters:
  • mo_energy (numpy.ndarray) – Molecular orbital energies for each spin channel.

  • se (tuple of dyson.Lehmann) – Self-energy to approximate for each spin channel.

Returns:

se_qp – Static potential approximation to the self-energy for each spin channel.

Return type:

numpy.ndarray

build_se_static(integrals)

Build the static part of the self-energy, including the Fock matrix.

Parameters:

integrals (UIntegrals) – Integrals object.

Returns:

se_static – Static part of the self-energy for each spin channel. If self.diagonal_se, non-diagonal elements are set to zero.

Return type:

numpy.ndarray

build_se_moments(nmom_max, integrals, **kwargs)

Build the moments of the self-energy.

Parameters:
  • nmom_max (int) – Maximum moment number to calculate.

  • integrals (UIntegrals) – Integrals object.

  • **kwargs (dict, optional) – Additional keyword arguments passed to polarizability class.

Returns:

  • se_moments_hole (tuple of numpy.ndarray) – Moments of the hole self-energy for each spin channel. If self.diagonal_se, non-diagonal elements are set to zero.

  • se_moments_part (tuple of numpy.ndarray) – Moments of the particle self-energy for each spin channel. If self.diagonal_se, non-diagonal elements are set to zero.

ao2mo(transform=True)

Get the integrals object.

Parameters:

transform (bool, optional) – Whether to transform the integrals object.

Returns:

integrals – Integrals object.

Return type:

UIntegrals

solve_dyson(se_moments_hole, se_moments_part, se_static, integrals=None)

Solve the Dyson equation due to a self-energy resulting from a list of hole and particle moments, along with a static contribution for each spin channel.

Also finds a chemical potential best satisfying the physical number of electrons. If self.optimise_chempot, this will shift the self-energy poles relative to the Green’s function, which is a partial self-consistency that better conserves the particle number.

If self.fock_loop, this function will also require that the outputted Green’s function is self-consistent with respect to the corresponding density and Fock matrix.

Parameters:
  • se_moments_hole (tuple of numpy.ndarray) – Moments of the hole self-energy for each spin channel.

  • se_moments_part (tuple of numpy.ndarray) – Moments of the particle self-energy for each spin channel.

  • se_static (tuple of numpy.ndarray) – Static part of the self-energy for each spin channel.

  • integrals (UIntegrals) – Integrals object. Required if self.fock_loop is True. Default value is None.

Returns:

  • gf (tuple of dyson.Lehmann) – Green’s function for each spin channel.

  • se (tuple of dyson.Lehmann) – Self-energy for each spin channel.

kernel(nmom_max, moments=None, integrals=None)

Driver for the method.

Parameters:
  • nmom_max (int) – Maximum moment number to calculate.

  • moments (tuple of numpy.ndarray, optional) – Tuple of (hole, particle) moments for each spin channel, if passed then they will be used instead of calculating them. Default value is None.

  • integrals (UIntegrals, optional) – Integrals object. If None, generate from scratch. Default value is None.

Returns:

  • converged (bool) – Whether the solver converged. For single-shot calculations, this is always True.

  • gf (tuple of dyson.Lehmann) – Green’s function object for each spin channel.

  • se (tuple of dyson.Lehmann) – Self-energy object for each spin channel.

  • qp_energy (NoneType) – Quasiparticle energies. For most GW methods, this is None.

make_rdm1(gf=None)

Get the first-order reduced density matrix.

Parameters:

gf (tuple of dyson.Lehmann, optional) – Green’s function for each spin channel. If None, use either self.gf, or the mean-field Green’s function. Default value is None.

Returns:

rdm1 – First-order reduced density matrix for each spin channel.

Return type:

numpy.ndarray

energy_hf(gf=None, integrals=None)

Calculate the one-body (Hartree–Fock) energy.

Parameters:
  • gf (tuple of dyson.Lehmann, optional) – Green’s function for each spin channel. If None, use either self.gf, or the mean-field Green’s function. Default value is None.

  • integrals (UIntegrals, optional) – Integrals object. If None, generate from scratch. Default value is None.

Returns:

e_1b – One-body energy.

Return type:

float

energy_gm(gf=None, se=None, g0=True)

Calculate the two-body (Galitskii–Migdal) energy.

Parameters:
  • gf (tuple of dyson.Lehmann, optional) – Green’s function for each spin channel. If None, use self.gf. Default value is None.

  • se (tuple of dyson.Lehmann, optional) – Self-energy for each spin channel. If None, use self.se. Default value is None.

  • g0 (bool, optional) – If True, use the mean-field Green’s function. Default value is True.

Returns:

e_2b – Two-body energy.

Return type:

float

init_gf(mo_energy=None)

Initialise the mean-field Green’s function.

Parameters:

mo_energy (tuple of numpy.ndarray, optional) – Molecular orbital energies for each spin channel. Default value is self.mo_energy.

Returns:

gf – Mean-field Green’s function for each spin channel.

Return type:

tuple of dyson.Lehmann

run(*args, **kwargs)

Alias for kernel, instead returning self.

Parameters:
  • *args (tuple) – Positional arguments to pass to kernel.

  • **kwargs (dict) – Keyword arguments to pass to kernel.

Returns:

self – The solver object.

Return type:

BaseGW

moment_error(se_moments_hole, se_moments_part, se)

Return the error in the moments.

Parameters:
  • se_moments_hole (numpy.ndarray) – Moments of the hole self-energy.

  • se_moments_part (numpy.ndarray) – Moments of the particle self-energy.

  • se (dyson.Lehmann) – Self-energy object.

Returns:

  • eh (float) – Error in the hole moments.

  • ep (float) – Error in the particle moments.

energy_nuc()

Calculate the nuclear repulsion energy.

Returns:

e_nuc – Nuclear repulsion energy.

Return type:

float