fcdmft.utils.interpolate module#

fcdmft.utils.interpolate.get_bands(mf, kpts_band, cell=None, dm_kpts=None, kpts=None)[source]#

Get energy bands at the given (arbitrary) ‘band’ k-points.

Returns:
mo_energy(nmo,) ndarray or a list of (nmo,) ndarray

Bands energies E_n(k)

mo_coeff(nao, nmo) ndarray or a list of (nao,nmo) ndarray

Band orbitals psi_n(k)

fcdmft.utils.interpolate.get_phase_wigner_seitz(cell, kpts, R_vec_rel)[source]#
fcdmft.utils.interpolate.get_wigner_seitz_supercell(w90, ws_search_size=[2, 2, 2], ws_distance_tol=1e-06)[source]#

Adpated from pyWannier90 (hungpham2017/pyWannier90)

Return a grid that contains all the lattice within the Wigner-Seitz supercell Ref: the hamiltonian_wigner_seitz(count_pts) in wannier90/src/hamittonian.F90

fcdmft.utils.interpolate.interpolate_hf_diff(mf, kpts_band, fock_dft_band, mo_interpolate=False, mo_energy_hf=None, return_fock=True, C_ao_lo=None, w90=None, wigner_seitz=True, veff=None, cell=None, dm=None, kpts=None)[source]#

Get HF bands by interpolating the difference between HF and DFT Fock matrices.

Args:

mf : mean-field object kpts_band : (nkpts_band,) ndarray fock_dft_band : (nkpts_band, nbands, nbands), DFT Fock matrix in LO basis at kpts_band mo_interpolate : whether use HF mo_energy or veff for interpolation mo_energy_hf : HF/GW mo_energy

Returns:
mo_energy(nbands,) ndarray or a list of (nbands,) ndarray

Bands energies E_n(k)

mo_coeff(nao, nbands) ndarray or a list of (nao,nbands) ndarray

Band orbitals psi_n(k)

fock_band_lo : (nkpts_band, nbands, nbands)

fcdmft.utils.interpolate.interpolate_mf(mf, kpts_band, return_fock=True, C_ao_lo=None, w90=None, wigner_seitz=True, veff=None, hcore=None, cell=None, dm=None, kpts=None)[source]#

Get energy bands at the given (arbitrary) ‘band’ k-points through interpolating Fock matrix.

Args:

mf : mean-field object kpts_band : (nkpts_band,) ndarray, k-points for interpolation C_ao_lo : (nkpts, nao, nbands) ndarray, localized orbitals at kpts w90 : Wannier90 object wigner_seitz : use Wigner-Seitz supercell for interpolation

Returns:
mo_energy(nbands,) ndarray or a list of (nbands,) ndarray

Bands energies E_n(k)

mo_coeff(nao, nbands) ndarray or a list of (nao,nbands) ndarray

Band orbitals psi_n(k)

fock_band_lo : (nkpts_band, nbands, nbands)

fcdmft.utils.interpolate.interpolate_selfenergy(mf, kpts_band, sigma, C_ao_lo=None, w90=None, wigner_seitz=True, cell=None, kpts=None)[source]#

Get self-energy at kpts_band by interpolation.

Args:

mf : mean-field object kpts_band : (nkpts_band,) ndarray sigma : (nkpts, nmo, nmo, nw), self-energy in MO basis at kpts

Returns:

sigma_band_lo : (nkpts_band, nlo, nlo, nw), self-energy in LO basis at kpts_band