Solver

class triqs_cthyb.Solver(beta, gf_struct, n_iw=1025, n_tau=10001, n_l=30)[source]

Bases: triqs_cthyb.solver_core.SolverCore

__init__(beta, gf_struct, n_iw=1025, n_tau=10001, n_l=30)[source]

Initialise the solver.

Parameters:

beta : scalar

Inverse temperature.

gf_struct : list of pairs [ (str,[int,…]), …]

Structure of the Green’s functions. It must be a list of pairs, each containing the name of the Green’s function block as a string and a list of integer indices. For example: [ ('up', [0, 1, 2]), ('down', [0, 1, 2]) ].

n_iw : integer, optional

Number of Matsubara frequencies used for the Green’s functions.

n_tau : integer, optional

Number of imaginary time points used for the Green’s functions.

n_l : integer, optional

Number of legendre polynomials to use in accumulations of the Green’s functions.

Delta_infty

\(G_0^{-1}(i\omega_n = \infty)\) in Matsubara Frequency.

Delta_tau

\(\Delta(\tau)\) in imaginary time.

G0_iw

\(G_0(i\omega)\) in imaginary frequencies.

G2_iw

Two-particle Green’s function \(G^{(2)}(i\nu,i\nu',i\nu'')\) (three Fermionic frequencies)

G2_iw_nfft

Two-particle Green’s function \(G^{(2)}(i\nu,i\nu',i\nu'')\) (three Fermionic frequencies)

G2_iw_ph

Two-particle Green’s function \(G^{(2)}(i\omega,i\nu,i\nu')\) in the ph-channel (one bosonic matsubara and two fermionic)

G2_iw_ph_nfft

Two-particle Green’s function \(G^{(2)}(i\omega,i\nu,i\nu')\) in the ph-channel (one bosonic matsubara and two fermionic)

G2_iw_pp

Two-particle Green’s function \(G^{(2)}(i\omega,i\nu,i\nu')\) in the pp-channel (one bosonic matsubara and two fermionic)

G2_iw_pp_nfft

Two-particle Green’s function \(G^{(2)}(i\omega,i\nu,i\nu')\) in the pp-channel (one bosonic matsubara and two fermionic)

G2_iwll_ph

Two-particle Green’s function \(G^{(2)}(i\omega,l,l')\) in the ph-channel (one bosonic matsubara and two legendre)

G2_iwll_pp

Two-particle Green’s function \(G^{(2)}(i\omega,l,l')\) in the pp-channel (one bosonic matsubara and two legendre)

G2_tau

Two-particle Green’s function \(G^{(2)}(\tau_1,\tau_2,\tau_3)\) (three Fermionic imaginary times)

G_l

Single-particle Green’s function \(G_l\) in Legendre polynomial representation.

G_tau

Single-particle Green’s function \(G(\tau)\) in imaginary time.

G_tau_accum

Intermediate Green’s function to accumulate g(tau), either real or complex

O_tau

General operator Green’s function \(O(\tau)\) in imaginary time.

average_sign

Monte Carlo average sign.

constr_parameters
density_matrix

Accumulated density matrix.

h_loc

The local Hamiltonian of the problem: \(H_{loc}\) used in the last call to solve().

h_loc_diagonalization

Diagonalization of \(H_{loc}\).

last_constr_parameters

Set of parameters used in the construction of the solver_core class.

last_solve_parameters

Set of parameters used in the last call to solve().

performance_analysis

Histograms related to the performance analysis.

perturbation_order

Histograms of the perturbation order for each block.

perturbation_order_total

Histogram of the total perturbation order.

solve(**params_kw)[source]

Solve the impurity problem. If measure_g_tau (default = True), G_iw and Sigma_iw will be calculated and their tails fitted. In addition to the solver parameters, parameters to control the tail fitting can be provided.

Parameters:

params_kw : dict {‘param’:value} that is passed to the core solver.

Two required parameters are
  • h_int (Operator object): the local Hamiltonian of the impurity problem to be solved,
  • n_cycles (int): number of measurements to be made.

perform_post_proc : boolean, optional, default = True

Should G_iw and Sigma_iw be calculated?

perform_tail_fit : boolean, optional, default = False

Should the tails of Sigma_iw and G_iw be fitted?

fit_max_moment : integer, optional, default = 3

Highest moment to fit in the tail of Sigma_iw.

fit_known_moments : ndarray.shape[order, Sigma_iw[0].target_shape], optional, default = None

Known moments of Sigma_iw, given as an numpy ndarray

fit_min_n : integer, optional, default = int(0.8 * self.n_iw)

Index of iw from which to start fitting.

fit_max_n : integer, optional, default = n_iw

Index of iw to fit until.

solve_parameters
solve_status

Status of the solve() on exit.