triqs.gf.dlr_crm_dyson_solver.minimize_dyson

triqs.gf.dlr_crm_dyson_solver.minimize_dyson(G0_dlr, G_dlr, Sigma_moments, method='trust-constr', options={'disp': True, 'finite_diff_rel_step': 1e-20, 'gtol': 1e-32, 'maxiter': 5000, 'xtol': 1e-100}, **kwargs)[source]

Contrained Residual Minimization Dyson solver as described in https://arxiv.org/abs/2310.01266

Defines the Dysons equation as an optimization problem:

G - G0 - G0*Σ*G = 0

and solves it using scipy.optimize.minimize using the DLR representation for the Green’s functions.

The solver optimizes only the dynamic part of the self-energy Σ_dyn(iν)= Σ(iν) - Σ_0. Here, Σ_0 is the Hartree shift. If provided the second moment Σ_1 is used as a non-linear constraint in the solver.

The moments can be explicitly calculated in the impurity solver, see for example the cthyb high frequency moments tutorial .

Alternatively the moments can be approximated by fitting the tail of the self-energy calculated via normal Dyson equation first:

>>> S_iw = inverse(G0_iw) - inverse(G_iw)
>>> tail, err = S_iw.fit_hermitian_tail()

and then used as input for the Dyson solver:

>>> S_iw_dlr, Sigma_HF, residual = minimize_dyson(G0_dlr=G0_dlr, G_dlr=G_dlr, Sigma_moments=tail[0:1])

The input G_dlr can for example obtained via fit_gf_dlr from a noisy imaginary time Green’s function or by directly setting the DLR mesh points from a full MeshImFreq G_iw object:

>>> for iwn in G_dlr_iw.mesh:
>>>     G_dlr_iw[iwn] = G_full_iw(iwn)
Parameters:
  • G0_dlr (triqs.gf.Gf or triqs.gf.BlockGf) – non-interacting Green’s function defined on a DLR, DLRImTime, or DLRImFreq mesh

  • G_dlr (triqs.gf.Gf or triqs.gf.BlockGf) – interacting Green’s function defined on a DLR, DLRImTime, or DLRImFreq mesh

  • Sigma_moments (list of numpy.ndarray or dict of list of numpy.ndarray) – moments of Σ. The first moment is the Hartree shift, i.e. the constant part of Σ. If provdided, use the second moment as a non-linear constraint for the Dyson solver.

  • method (str, optional) – optimization method, defaults to ‘trust-constr’ Note: For non-linear constraints this is one of the few available methods

  • options (dict, optional) – optimization options, defaults to dict(maxiter=5000, disp=True, gtol=1e-32, xtol=1e-100, finite_diff_rel_step=1e-20)

Returns:

  • Sigma_DLR (triqs.gf.Gf or triqs.gf.BlockGf) – optimized self-energy defined on a DLRImFreq mesh

  • Sigma_0 (numpy.ndarray) – Hartree shift

  • residual (float) – L2 norm of residual (G-G₀-G₀ΣG)