SumK DFT

class triqs_dft_tools.sumk_dft.SumkDFT(hdf_file, h_field=0.0, use_dft_blocks=False, dft_data='dft_input', symmcorr_data='dft_symmcorr_input', parproj_data='dft_parproj_input', symmpar_data='dft_symmpar_input', bands_data='dft_bands_input', transp_data='dft_transp_input', misc_data='dft_misc_input')[source]

Bases: object

This class provides a general SumK method for combining ab-initio code and pytriqs.

__init__(hdf_file, h_field=0.0, use_dft_blocks=False, dft_data='dft_input', symmcorr_data='dft_symmcorr_input', parproj_data='dft_parproj_input', symmpar_data='dft_symmpar_input', bands_data='dft_bands_input', transp_data='dft_transp_input', misc_data='dft_misc_input')[source]

Initialises the class from data previously stored into an hdf5 archive.

Parameters:

hdf_file : string

Name of hdf5 containing the data.

h_field : scalar, optional

The value of magnetic field to add to the DFT Hamiltonian. The contribution -h_field*sigma is added to diagonal elements of the Hamiltonian. It cannot be used with the spin-orbit coupling on; namely h_field is set to 0 if self.SO=True.

use_dft_blocks : boolean, optional

If True, the local Green’s function matrix for each spin is divided into smaller blocks

with the block structure determined from the DFT density matrix of the corresponding correlated shell.

Alternatively and additionally, the block structure can be analysed using analyse_block_structure and manipulated using the SumkDFT.block_structre attribute (see BlockStructure).

dft_data : string, optional

Name of hdf5 subgroup in which DFT data for projector and lattice Green’s function construction are stored.

symmcorr_data : string, optional

Name of hdf5 subgroup in which DFT data on symmetries of correlated shells (symmetry operations, permutaion matrices etc.) are stored.

parproj_data : string, optional

Name of hdf5 subgroup in which DFT data on non-normalized projectors for non-correlated states (used in the partial density of states calculations) are stored.

symmpar_data : string, optional

Name of hdf5 subgroup in which DFT data on symmetries of the non-normalized projectors are stored.

bands_data : string, optional

Name of hdf5 subgroup in which DFT data necessary for band-structure/k-resolved spectral function calculations (projectors, DFT Hamiltonian for a chosen path in the Brillouin zone etc.) are stored.

transp_data : string, optional

Name of hdf5 subgroup in which DFT data necessary for transport calculations are stored.

misc_data : string, optional

Name of hdf5 subgroup in which miscellaneous DFT data are stored.

__weakref__

list of weak references to the object (if defined)

add_dc(iw_or_w='iw')[source]

Subtracts the double counting term from the impurity self energy.

Parameters:

iw_or_w : string, optional

  • iw_or_w = ‘iw’ for a imaginary-frequency self-energy
  • iw_or_w = ‘w’ for a real-frequency self-energy
Returns:

sigma_minus_dc : gf_struct_sumk like

Self-energy with a subtracted double-counting term.

analyse_block_structure(threshold=1e-05, include_shells=None, dm=None, hloc=None)[source]

Determines the block structure of local Green’s functions by analysing the structure of the corresponding density matrices and the local Hamiltonian. The resulting block structures for correlated shells are stored in the SumkDFT.block_structure attribute.

Parameters:

threshold : real, optional

If the difference between density matrix / hloc elements is below threshold, they are considered to be equal.

include_shells : list of integers, optional

List of correlated shells to be analysed. If include_shells is not provided all correlated shells will be analysed.

dm : list of dict, optional

List of density matrices from which block stuctures are to be analysed. Each density matrix is a dict {block names: 2d numpy arrays}. If not provided, dm will be calculated from the DFT Hamiltonian by a simple-point BZ integration.

hloc : list of dict, optional

List of local Hamiltonian matrices from which block stuctures are to be analysed Each Hamiltonian is a dict {block names: 2d numpy arrays}. If not provided, it will be calculated using eff_atomic_levels.

analyse_block_structure_from_gf(G, threshold=1e-05, include_shells=None, analyse_deg_shells=True)[source]

Determines the block structure of local Green’s functions by analysing the structure of the corresponding non-interacting Green’s function. The resulting block structures for correlated shells are stored in the SumkDFT.block_structure attribute.

This is a safer alternative to analyse_block_structure, because the full non-interacting Green’s function is taken into account and not just the density matrix and Hloc.

Parameters:

G : list of BlockGf of GfImFreq, GfImTime, GfReFreq or GfReTime

the non-interacting Green’s function for each inequivalent correlated shell

threshold : real, optional

If the difference between matrix elements is below threshold, they are considered to be equal.

include_shells : list of integers, optional

List of correlated shells to be analysed. If include_shells is not provided all correlated shells will be analysed.

analyse_deg_shells : bool

Whether to call the analyse_deg_shells function after having finished the block structure analysis

Returns:

G : list of BlockGf of GfImFreq or GfImTime

the Green’s function transformed into the new block structure

analyse_deg_shells(G, threshold=1e-05, include_shells=None)[source]

Determines the degenerate shells of local Green’s functions by analysing the structure of the corresponding non-interacting Green’s function. The results are stored in the SumkDFT.block_structure attribute.

Due to the implementation and numerics, the maximum difference between two matrix elements that are detected as equal can be a bit higher (e.g. a factor of two) than the actual threshold.

Parameters:

G : list of BlockGf of GfImFreq or GfImTime

the non-interacting Green’s function for each inequivalent correlated shell

threshold : real, optional

If the difference between matrix elements is below threshold, they are considered to be equal.

include_shells : list of integers, optional

List of correlated shells to be analysed. If include_shells is not provided all correlated shells will be analysed.

calc_dc(dens_mat, orb=0, U_interact=None, J_hund=None, use_dc_formula=0, use_dc_value=None)[source]

Calculates and sets the double counting corrections.

If ‘use_dc_value’ is provided the double-counting term is uniformly initialized with this constant and ‘U_interact’ and ‘J_hund’ are ignored.

If ‘use_dc_value’ is None the correction is evaluated according to one of the following formulae:

  • use_dc_formula = 0: fully-localised limit (FLL)
  • use_dc_formula = 1: Held’s formula, i.e. mean-field formula for the Kanamori
    type of the interaction Hamiltonian
  • use_dc_formula = 2: around mean-field (AMF)

Note that FLL and AMF formulae were derived assuming a full Slater-type interaction term and should be thus used accordingly. For the Kanamori-type interaction one should use formula 1.

The double-counting self-energy term is stored in self.dc_imp and the energy correction in self.dc_energ.

Parameters:

dens_mat : gf_struct_solver like

Density matrix for the specified correlated shell.

orb : int, optional

Index of an inequivalent shell.

U_interact : float, optional

Value of interaction parameter U.

J_hund : float, optional

Value of interaction parameter J.

use_dc_formula : int, optional

Type of double-counting correction (see description).

use_dc_value : float, optional

Value of the double-counting correction. If specified U_interact, J_hund and use_dc_formula are ignored.

calc_dc_for_density(orb, dc_init, dens_mat, density=None, precision=0.01)[source]

Searches for DC in order to fulfill charge neutrality. If density is given, then DC is set such that the LOCAL charge of orbital orb coincides with the given density.

calc_density_correction(filename=None, dm_type='wien2k')[source]

Calculates the charge density correction and stores it into a file.

The charge density correction is needed for charge-self-consistent DFT+DMFT calculations. It represents a density matrix of the interacting system defined in Bloch basis and it is calculated from the sum over Matsubara frequecies of the full GF,

..math:: N_{nunu’}(k) = sum_{iomega_{n}} G_{nunu’}(k, iomega_{n})

The density matrix for every k-point is stored into a file.

Parameters:

filename : string

Name of the file to store the charge density correction.

Returns:

(deltaN, dens) : tuple

Returns a tuple containing the density matrix deltaN and the corresponing total charge dens.

calc_mu(precision=0.01, iw_or_w='iw', broadening=None, delta=0.5)[source]

Searches for the chemical potential that gives the DFT total charge. A simple bisection method is used.

Parameters:

precision : float, optional

A desired precision of the resulting total charge.

iw_or_w : string, optional

  • iw_or_w = ‘iw’ for a imaginary-frequency self-energy
  • iw_or_w = ‘w’ for a real-frequency self-energy

broadening : float, optional

Imaginary shift for the axis along which the real-axis GF is calculated. If not provided, broadening will be set to double of the distance between mesh points in ‘mesh’. Only relevant for real-frequency GF.

Returns:

mu : float

Value of the chemical potential giving the DFT total charge within specified precision.

check_projectors()[source]

Calculated the density matrix from projectors (DM = P Pdagger) to check that it is correct and specifically that it matches DFT.

density_matrix(method='using_gf', beta=40.0)[source]

Calculate density matrices in one of two ways.

Parameters:

method : string, optional

  • if ‘using_gf’: First get lattice gf (g_loc is not set up), then density matrix.
    It is useful for Hubbard I, and very quick. No assumption on the hopping structure is made (ie diagonal or not).
  • if ‘using_point_integration’: Only works for diagonal hopping matrix (true in wien2k).

beta : float, optional

Inverse temperature.

Returns:

dens_mat : list of dicts

Density matrix for each spin in each correlated shell.

downfold(ik, ish, bname, gf_to_downfold, gf_inp, shells='corr', ir=None)[source]

Downfolds a block of the Green’s function for a given shell and k-point using the corresponding projector matrices.

Parameters:

ik : integer

k-point index for which the downfolding is to be done.

ish : integer

Shell index of GF to be downfolded.

  • if shells=’corr’: ish labels all correlated shells (equivalent or not)
  • if shells=’all’: ish labels only representative (inequivalent) non-correlated shells

bname : string

Block name of the target block of the lattice Green’s function.

gf_to_downfold : Gf

Block of the Green’s function that is to be downfolded.

gf_inp : Gf

FIXME

shells : string, optional

  • if shells=’corr’: orthonormalized projectors for correlated shells are used for the downfolding.
  • if shells=’all’: non-normalized projectors for all included shells are used for the downfolding.

ir : integer, optional

Index of equivalent site in the non-correlated shell ‘ish’, only used if shells=’all’.

Returns:

gf_downfolded : Gf

Downfolded block of the lattice Green’s function.

eff_atomic_levels()[source]

Calculates the effective local Hamiltonian required as an input for the Hubbard I Solver. The local Hamiltonian (effective atomic levels) is calculated by projecting the on-site Bloch Hamiltonian:

\[H^{loc}_{m m'} = \sum_{k} P_{m \nu}(k) H_{\nu\nu'}(k) P^{*}_{\nu' m'}(k),\]

where

\[H_{\nu\nu'}(k) = [\epsilon_{\nu k} - h_{z} \sigma_{z}] \delta_{\nu\nu'}.\]
Parameters:

None :

Returns:

eff_atlevels : gf_struct_sumk like

Effective local Hamiltonian \(H^{loc}_{m m'}\) for each inequivalent correlated shell.

extract_G_loc(mu=None, iw_or_w='iw', with_Sigma=True, with_dc=True, broadening=None)[source]

Extracts the local downfolded Green function by the Brillouin-zone integration of the lattice Green’s function.

Parameters:

mu : real, optional

Input chemical potential. If not provided the value of self.chemical_potential is used as mu.

with_Sigma : boolean, optional

If True then the local GF is calculated with the self-energy self.Sigma_imp.

with_dc : boolean, optional

If True then the double-counting correction is subtracted from the self-energy in calculating the GF.

broadening : float, optional

Imaginary shift for the axis along which the real-axis GF is calculated. If not provided, broadening will be set to double of the distance between mesh points in ‘mesh’. Only relevant for real-frequency GF.

Returns:

G_loc_inequiv : list of BlockGf (Green’s function) objects

List of the local Green’s functions for all inequivalent correlated shells, rotated into the corresponding local frames.

init_dc()[source]

Initializes the double counting terms.

Parameters:None :
lattice_gf(ik, mu=None, iw_or_w='iw', beta=40, broadening=None, mesh=None, with_Sigma=True, with_dc=True)[source]

Calculates the lattice Green function for a given k-point from the DFT Hamiltonian and the self energy.

Parameters:

ik : integer

k-point index.

mu : real, optional

Chemical potential for which the Green’s function is to be calculated. If not provided, self.chemical_potential is used for mu.

iw_or_w : string, optional

  • iw_or_w = ‘iw’ for a imaginary-frequency self-energy
  • iw_or_w = ‘w’ for a real-frequency self-energy

beta : real, optional

Inverse temperature.

broadening : real, optional

Imaginary shift for the axis along which the real-axis GF is calculated. If not provided, broadening will be set to double of the distance between mesh points in ‘mesh’.

mesh : list, optional

Data defining mesh on which the real-axis GF will be calculated, given in the form (om_min,om_max,n_points), where om_min is the minimum omega, om_max is the maximum omega and n_points is the number of points.

with_Sigma : boolean, optional

If True the GF will be calculated with the self-energy stored in self.Sigmaimp_(w/iw), for real/Matsubara GF, respectively. In this case the mesh is taken from the self.Sigma_imp object. If with_Sigma=True but self.Sigmaimp_(w/iw) is not present, with_Sigma is reset to False.

with_dc : boolean, optional

if True and with_Sigma=True, the dc correction is substracted from the self-energy before it is included into GF.

Returns:

G_latt : BlockGf

Lattice Green’s function.

load(things_to_load, subgrp='user_data')[source]

Loads user data from the HDF file. Raises an exeption if a requested dataset is not found.

Parameters:

things_to_read : list of strings

List of datasets to be read from the hdf5 file.

subgrp : string, optional

Name of hdf5 file subgroup from which the data are to be read.

Returns:

list_to_return : list

A list containing data read from hdf5.

number_of_atoms(shells)[source]

Determine the number of inequivalent atoms.

put_Sigma(Sigma_imp)[source]

Inserts the impurity self-energies into the sumk_dft class.

Parameters:

Sigma_imp : list of BlockGf (Green’s function) objects

List containing impurity self-energy for all inequivalent correlated shells. Self-energies for equivalent shells are then automatically set by this function. The self-energies can be of the real or imaginary-frequency type.

read_input_from_hdf(subgrp, things_to_read)[source]

Reads data from the HDF file. Prints a warning if a requested dataset is not found.

Parameters:

subgrp : string

Name of hdf5 file subgroup from which the data are to be read.

things_to_read : list of strings

List of datasets to be read from the hdf5 file.

Returns:

subgroup_present : boolean

Is the subgrp is present in hdf5 file?

value_read : boolean

Did the reading of requested datasets succeed?

rotloc(ish, gf_to_rotate, direction, shells='corr')[source]

Rotates a block of the local Green’s function from the local frame to the global frame and vice versa.

Parameters:

ish : integer

Shell index of GF to be rotated.

  • if shells=’corr’: ish labels all correlated shells (equivalent or not)
  • if shells=’all’: ish labels only representative (inequivalent) non-correlated shells

gf_to_rotate : Gf

Block of the Green’s function that is to be rotated.

direction : string

The direction of rotation can be either

  • ‘toLocal’ : global -> local transformation,
  • ‘toGlobal’ : local -> global transformation.

shells : string, optional

  • if shells=’corr’: the rotation matrix for the correlated shell ‘ish’ is used,
  • if shells=’all’: the rotation matrix for the generic (non-correlated) shell ‘ish’ is used.
Returns:

gf_rotated : Gf

Rotated block of the local Green’s function.

save(things_to_save, subgrp='user_data')[source]

Saves data from a list into the HDF file. Prints a warning if a requested data is not found in SumkDFT object.

Parameters:

things_to_save : list of strings

List of datasets to be saved into the hdf5 file.

subgrp : string, optional

Name of hdf5 file subgroup in which the data are to be stored.

set_dc(dc_imp, dc_energ)[source]

Sets double counting corrections to given values.

Parameters:

dc_imp : gf_struct_sumk like

Double-counting self-energy term.

dc_energ : list of floats

Double-counting energy corrections for each correlated shell.

set_mu(mu)[source]

Sets a new chemical potential.

Parameters:

mu : float

New value of the chemical potential.

sorts_of_atoms(shells)[source]

Determine the number of inequivalent sorts.

symm_deg_gf(gf_to_symm, orb)[source]

Averages a GF over degenerate shells.

Degenerate shells of an inequivalent correlated shell are defined by self.deg_shells. This function enforces corresponding degeneracies in the input GF.

Parameters:

gf_to_symm : gf_struct_solver like

Input and output GF (i.e., it gets overwritten)

orb : int

Index of an inequivalent shell.

total_density(mu=None, iw_or_w='iw', with_Sigma=True, with_dc=True, broadening=None)[source]

Calculates the total charge within the energy window for a given chemical potential. The chemical potential is either given by parameter mu or, if it is not specified, taken from self.chemical_potential.

The total charge is calculated from the trace of the GF in the Bloch basis. By default, a full interacting GF is used. To use the non-interacting GF, set parameter with_Sigma = False.

The number of bands within the energy windows generally depends on k. The trace is therefore calculated separately for each k-point.

Since in general n_orbitals depends on k, the calculation is done in the following order:

\[n_{tot} = \sum_{k} n(k),\]

with

\[n(k) = Tr G_{\nu\nu'}(k, i\omega_{n}). \]

The calculation is done in the global coordinate system, if distinction is made between local/global.

Parameters:

mu : float, optional

Input chemical potential. If not specified, self.chemical_potential is used instead.

iw_or_w : string, optional

  • iw_or_w = ‘iw’ for a imaginary-frequency self-energy
  • iw_or_w = ‘w’ for a real-frequency self-energy

with_Sigma : boolean, optional

If True the full interacing GF is evaluated, otherwise the self-energy is not included and the charge would correspond to a non-interacting system.

with_dc : boolean, optional

Whether or not to subtract the double-counting term from the self-energy.

broadening : float, optional

Imaginary shift for the axis along which the real-axis GF is calculated. If not provided, broadening will be set to double of the distance between mesh points in ‘mesh’. Only relevant for real-frequency GF.

Returns:

dens : float

Total charge \(n_{tot}\).

upfold(ik, ish, bname, gf_to_upfold, gf_inp, shells='corr', ir=None)[source]

Upfolds a block of the Green’s function for a given shell and k-point using the corresponding projector matrices.

Parameters:

ik : integer

k-point index for which the upfolding is to be done.

ish : integer

Shell index of GF to be upfolded.

  • if shells=’corr’: ish labels all correlated shells (equivalent or not)
  • if shells=’all’: ish labels only representative (inequivalent) non-correlated shells

bname : string

Block name of the target block of the lattice Green’s function.

gf_to_upfold : Gf

Block of the Green’s function that is to be upfolded.

gf_inp : Gf

FIXME

shells : string, optional

  • if shells=’corr’: orthonormalized projectors for correlated shells are used for the upfolding.
  • if shells=’all’: non-normalized projectors for all included shells are used for the upfolding.

ir : integer, optional

Index of equivalent site in the non-correlated shell ‘ish’, only used if shells=’all’.

Returns:

gf_upfolded : Gf

Upfolded block of the lattice Green’s function.