Tools for GFs

triqs.gf.tools.conjugate(x)[source]

Return the conjugate of a Green’s function

triqs.gf.tools.delta(g)[source]

Compute Delta_iw from G0_iw. CAUTION: This function assumes the following properties of g

  • The diagonal components of g should decay as 1/iOmega

  • g should fullfill the property g[iw][i,j] = conj(g[-iw][j,i])

Parameters:

g (BlockGf (of GfImFreq) or GfImFreq) – Non-interacting Green’s function.

Returns:

delta_iw – Hybridization function.

Return type:

BlockGf (of GfImFreq) or GfImFreq

triqs.gf.tools.discretize_bath(delta_in, Nb, eps0=3, V0=None, tol=1e-15, maxiter=10000, cmplx=False, method='BFGS')[source]

Discretize a given hybridization function using Nb bath sites.

The discretized hybridization is constructed as .. math:: Delta_{kl}^{disc} (i omega_n) = sum_{j=1}^{Nb} V_{kj} S V_{jl}^* . where S is either: .. math:: [i omega_n - eps_j]^{-1} for MeshImFreq or .. math:: - exp(-tau * eps_j) / (1 + exp(-eta * eps_j) ) for MeshImTime.

The hoppings V and energies eps are chosen to minimize the norm .. math:: left[

rac{1}{sqrt(N)} sum_{i omega_n}^{N} | Delta^{disc} (i omega_n) - Delta (i omega_n) |^2 ight]^{ rac{1}{2}}

and for MeshImTime .. math:: left[

rac{1}{sqrt(N)} sum_{ au}^{N} | Delta^{disc} ( au) - Delta ( au) |^2 ight]^{ rac{1}{2}}

This minimization is performed with the given tolerance using scipy.optimize.minimize or the scipy.optimize.basinhopping frontend.

delta_inGf or BlockGf

Matsubara or imaginary-time hybridization function to discretize

Nbint

Number of bath sites per Gf block

eps0: float or list(float), default=3.0

Approximate bandwith or initial guesses for bath energies.

V0float or np.ndarray (shape norb X Nb), optional

If float: initial guess used for all hopping values. If np.ndarray: initial guess for V. Otherwise use the cholesky decomposition of .. math:: lim_{omega->infty} iomega*Delta(iomega) or .. math:: -Delta( au=0^+) - Delta( au=eta^-) to obtain an initial guess for V.

tolfloat, default=1e-15

Tolerance for scipy minimize on data to optimize (xatol / ftol)

maxiterint, default=10000

Maximum number of optimization steps

complxbool, default=False

Allow the hoppings V to be complex

methodstring, default=BFGS

Method for minimizing the function. Should be one of ‘BFGS’, ‘Nelder-Mead’ and ‘basinhopping’.

V_optnp.array (shape norb x Nb) or list thereof

Optimized bath hoppings

eps_optlist(float) or list(list(float)

Optimized bath energies (sorted)

delta_discGf or BlockGf

Discretized hybridization function

triqs.gf.tools.dyson(**kwargs)[source]

Solve Dyson’s equation for given two of G0_iw, G_iw and Sigma_iw to yield the third.

Parameters:
  • G0_iw (Gf, optional) – Non-interacting Green’s function.

  • G_iw (Gf, optional) – Interacting Green’s function.

  • Sigma_iw (Gf, optional) – Self-energy.

triqs.gf.tools.fit_legendre(g_t, order=10)[source]

General fit of a noisy imaginary time Green’s function to a low order Legendre expansion in imaginary time.

Only Hermiticity is imposed on the fit, so discontinuities has to be fixed separately (see the method enforce_discontinuity)

Author: Hugo U.R. Strand

Parameters:
  • g_t (TRIQS imaginary time Green's function (matrix valued)) – Imaginary time Green’s function to fit (possibly noisy binned data)

  • order (int) – Maximal order of the fitted Legendre expansion

Returns:

g_l – Fitted Legendre Green’s function with order order

Return type:

TRIQS Legendre polynomial Green’s function (matrix valued)

triqs.gf.tools.inverse(x)[source]

Return the inverse of a Green’s function

triqs.gf.tools.make_delta(V, eps, mesh, block_names=None)[source]

Create a hybridization function from given hoppings and bath energies as .. math:: Delta_{kl}^{disc} (i omega_n) = sum_{j=1}^{Nb} V_{kj} S V_{jl}^* . where S is either .. math:: [i omega_n - eps_j]^{-1} for MeshImFreq or .. math:: - exp(-tau * eps_j) / (1 + exp(-eta * eps_j) ) for MeshImTime

Parameters:
  • V (np.array (shape Norb x NB) or list thereof) – Bath hopping matrix or matrix for each Gf block

  • eps (list(float) or list(list(float))) – Bath energies or energies for each Gf block

  • mesh (MeshImFreq or MeshImTime) – Mesh of the hybridization function

  • block_names (list(str)) – List of block names, used if V, eps are lists

Returns:

delta – Hybridization function on given mesh

Return type:

Gf or BlockGf

triqs.gf.tools.make_zero_tail(g, n_moments=10)[source]

Return a container for the high-frequency coefficients of a given Green function initialized to zero.

Parameters:
  • g (GfImFreq or GfReFreq or GfImTime or GfReTime) – The real/imaginary frequency/time Green’s function that we create the tail-array for.

  • [default=10] (n_moments) –

triqs.gf.tools.read_gf_from_txt(block_txtfiles, block_name)[source]

Read a GfReFreq from text files with the format (w, Re(G), Im(G)) for a single block.

Notes

A BlockGf must be constructed from multiple GfReFreq objects if desired. The mesh must be the same for all files read in. Non-uniform meshes are not supported.

Parameters:
  • block_txtfiles (Rank 2 square np.array(str) or list[list[str]]) –

    The text files containing the GF data that need to read for the block. e.g. [[‘up_eg1.dat’]] for a one-dimensional block and

    [[‘up_eg1_1.dat’,’up_eg2_1.dat’],

    [‘up_eg1_2.dat’,’up_eg2_2.dat’]] for a 2x2 block.

  • block_name (str) – Name of the block.

Returns:

g – The real frequency Green’s function read in.

Return type:

GfReFreq

triqs.gf.tools.transpose(x)[source]

Return the transpose of a Green’s function

triqs.gf.tools.write_gf_to_txt(g)[source]

Write a GfReFreq or GfImFreq to in the format (w/iw, Re(G), Im(G)) for a single block.

Parameters:

g (GfReFreq or GfImFreq) – The real/imaginary frequency Green’s function to be written out.