Documentation

Solver class

class edipack2triqs.solver.EDIpackSolver(hamiltonian: Operator, fops_imp_up: list[tuple[int | str, int | str]], fops_imp_dn: list[tuple[int | str, int | str]], fops_bath_up: list[tuple[int | str, int | str]], fops_bath_dn: list[tuple[int | str, int | str]], **kwargs)[source]

This class represents the EDIpack exact diagonalization library and incapsulates its internal state. Its methods and attributes allow to perform ED calculations and to access their results.

Note

At most one instance of this type can exists at any time. An attempt to create more such objects will raise an AssertionError.

__init__(hamiltonian: Operator, fops_imp_up: list[tuple[int | str, int | str]], fops_imp_dn: list[tuple[int | str, int | str]], fops_bath_up: list[tuple[int | str, int | str]], fops_bath_dn: list[tuple[int | str, int | str]], **kwargs)[source]

Initialize internal state of the underlying EDIpack solver.

The fundamental operator sets (fops for short) define sizes of impurity and bath Hilbert spaces. The expected fops objects are lists of tuples (pairs). Each element (b, i) of such a list corresponds to a single fermionic degree of freedom created by the many-body operator c_dag(b, i). b and i are a (string or integer) block index and an index within a block respectively, which are used to construct output Green's function containers.

Parameters:
  • hamiltonian – Many-body electronic Hamiltonian to diagonalize. Symmetries of this Hamiltonian are automatically analyzed and dictate the choice of EDIpacks’s ED mode and bath geometry. This choice remains unchangeable throughout object’s lifetime.

  • fops_imp_up (list[tuple[int | str, int | str]]) – Fundamental operator set for spin-up impurity degrees of freedom.

  • fops_imp_dn (list[tuple[int | str, int | str]]) – Fundamental operator set for spin-down impurity degrees of freedom.

  • fops_bath_up (list[tuple[int | str, int | str]]) – Fundamental operator set for spin-up bath degrees of freedom.

  • fops_bath_dn (list[tuple[int | str, int | str]]) – Fundamental operator set for spin-down bath degrees of freedom.

  • input_file (str, optional) – Path to a custom input file compatible with EDIpack’s f/ed_input_vars/ed_read_input(). When specified, it is used to initialize EDIpack. This option is mutually exclusive with all other keyword arguments.

  • verbose (int, default=3) – Verbosity level of the solver: 0 for almost no output, 5 for the full output from EDIpack.

  • print_input_vars (bool, default=False) – Flag to toggle the printing of all the input variables and their values on the console.

  • cutoff (float, default=1e-9) – Spectrum cutoff used to determine the number of states to be retained.

  • gs_threshold (float, default=1e-9) – Energy threshold for ground state degeneracy loop.

  • ed_sparse_h (bool, default=True) – Switch between storing the sparse matrix \(\hat H\) (True) and direct on-the-fly evaluation of the product \(\hat H |v\rangle\) (False).

  • ed_total_ud (bool, default=False) – Force use of the total spin-up and spin-down occupancies as quantum numbers instead of the orbital-resolved occupancies.

  • lanc_method (str, default="arpack") – Select the method to be used in the determination of the spectrum, one of “arpack”, “lanczos” (simple Lanczos method, only works at zero temperature, can be useful in some rare cases, when ARPACK fails) and “dvdson” (no-MPI mode only).

  • lanc_nstates_sector (int, default=2) – Initial number of states per sector to be determined.

  • lanc_nstates_total (int, default=2) – Initial total number of states to be determined.

  • lanc_nstates_step (int, default=2) – Number of states added to the spectrum at each step.

  • lanc_ncv_factor (int, default=10) – Set the size of the block used in Lanczos-ARPACK by multiplying the required Neigen (NCV = lanc_ncv_factor * Neigen + lanc_ncv_add).

  • lanc_ncv_add (int, default=0) – Adds up to the size of the block to prevent it from becoming too small (NCV = lanc_ncv_factor * Neigen + lanc_ncv_add).

  • lanc_niter (int, default=512) – Number of Lanczos iterations in spectrum determination.

  • lanc_ngfiter (int, default=200) – Number of Lanczos iterations in GF determination (number of momenta).

  • lanc_tolerance (float, default=1e-18) – Tolerance for the Lanczos iterations as used in ARPACK.

  • lanc_dim_threshold (int, default=1024) – Minimal dimension threshold to use Lanczos determination of the spectrum rather than LAPACK-based exact diagonalization.

  • bath_fitting_params (BathFittingParams, optional) – Parameters used to perform bath fitting.

property hloc: ndarray

Read-only access to the matrix of the local impurity Hamiltonian \(\hat H_\text{loc}\).

property Uloc: ndarray

Local intra-orbital interaction \(U_\text{loc}\), one value per orbital, up to 5 values.

property Ust: float

Local inter-orbital interaction \(U_\text{st}\).

property Jh: float

Hund’s coupling \(J_h\).

property Jx: float

Spin-exchange coupling \(J_x\).

property Jp: float

Pair-hopping coupling \(J_p\).

property bath: Bath

Access to the current bath object stored in the solver. It is possible to assign a new bath object as long as it has the matching type and describes a bath of the same geometry.

solve(beta: float, *, n_iw: int = 4096, energy_window: tuple[float, float] = (-5.0, 5.0), n_w: int = 5000, broadening: float = 0.01)[source]

Solve the impurity problem and calculate the observables, Green’s function and self-energy.

Parameters:
  • beta (float) – Inverse temperature \(\beta = 1 / (k_B T)\) for observable and Green’s function calculations.

  • n_iw (int, default=4096) – Number of Matsubara frequencies for impurity GF calculations.

  • energy_window (tuple[float, float], default=(-5.0, 5.0)) – Energy window for real-frequency impurity GF calculations.

  • n_w (int, default=5000) – Number of real-frequency points for impurity GF calculations.

  • broadening (float, default=0.01) – Broadening (imaginary shift away from the real-frequency axis) for real-frequency impurity GF calculations.

property e_pot: float

Potential energy from interaction.

property e_kin: float

Kinetic energy.

property densities: ndarray

Impurity occupations, one element per orbital.

property double_occ: ndarray

Impurity double occupancy, one element per orbital.

property superconductive_phi: ndarray

Impurity superconductive \(\phi\)-matrix in orbital space.

property magnetization: ndarray

Cartesian components of impurity magnetization vectors, one row per orbital.

property g_iw: BlockGf

Matsubara impurity Green’s function.

property g_an_iw: BlockGf

Anomalous Matsubara impurity Green’s function.

property Sigma_iw: BlockGf

Matsubara impurity self-energy.

property Sigma_an_iw: BlockGf

Anomalous Matsubara impurity self-energy.

property g_w: BlockGf

Real-frequency impurity Green’s function.

property g_an_w: BlockGf

Anomalous real-frequency impurity Green’s function.

property Sigma_w: BlockGf

Real-frequency impurity self-energy.

property Sigma_an_w: BlockGf

Anomalous real-frequency impurity self-energy.

Objects representing sets of bath parameters

class edipack2triqs.bath.Bath[source]

Base class for all bath classes.

The following classes derived from Bath represent sets of bath parameters used by EDIpack. Each of the classes corresponds to one bath geometry (normal, hybrid or general). They support addition, subtraction and multiplication by a scalar, which translate into the respective operations with the stored parameters.

class edipack2triqs.bath.BathNormal(ed_mode: str, nspin: int, norb: int, nbath: int, data: ndarray | None = None)[source]

Bases: Bath

Parameters of a bath with normal topology. The normal topology means that each of norb impurity orbitals is connected to nbath independent bath levels. There are, therefore, norb * nbath bath levels in total.

Instances of this class are compatible with TRIQS’ HDF5 interface.

nbath: int

Number of bath levels per impurity orbital.

eps: ndarray

Bath energy levels as an array of the shape (nspin, norb, nbath).

U: ndarray | None

Spin off-diagonal hopping amplitudes between the impurity and the bath as an array of the shape (nspin, norb, nbath). Only present when ed_mode == "nonsu2".

Delta: ndarray | None

Local super-conducting order parameters of the bath as an array of the shape (nspin, norb, nbath). Only present when ed_mode == "superc".

V: ndarray

Spin-diagonal hopping amplitudes between the impurity and the bath as an array of the shape (nspin, norb, nbath).

property ed_mode

ED mode this bath object is usable with, one of “normal”, “nonsu2” and “superc”.

class edipack2triqs.bath.BathHybrid(ed_mode: str, nspin: int, norb: int, nbath: int, data: ndarray | None = None)[source]

Bases: Bath

Parameters of a bath with hybrid topology. In the hybrid topology there are nbath independent bath levels. Each of these levels is connected to each impurity orbital via hopping amplitudes V.

Instances of this class are compatible with TRIQS’ HDF5 interface.

nbath: int

Total number of bath levels.

eps: ndarray

Bath energy levels as an array of the shape (nspin, nbath).

U: ndarray | None

Spin off-diagonal hopping amplitudes between the impurity and the bath as an array of the shape (nspin, norb, nbath). Only present when ed_mode == "nonsu2".

Delta: ndarray | None

Local super-conducting order parameters of the bath as an array of the shape (nspin, nbath). Only present when ed_mode == "superc".

V: ndarray

Spin-diagonal hopping amplitudes between the impurity and the bath as an array of the shape (nspin, norb, nbath).

property ed_mode

ED mode this bath object is usable with, one of “normal”, “nonsu2” and “superc”.

class edipack2triqs.bath.BathGeneral(nspin: int, norb: int, nbath: int, hvec: ndarray, data: ndarray | None = None)[source]

Bases: Bath

Parameters of a bath with general topology. General bath is a set of nbath independent replicas of the impurity. Hamiltonian of each replica is constructed as a linear combination of nsym basis matrices \(\hat O_i\) with coefficients \(\lambda^\nu_i\). Each impurity orbital is coupled to the corresponding orbital of each replica via hopping amplitudes V.

Instances of this class are compatible with TRIQS’ HDF5 interface.

nbath: int

Number of replicas.

hvec: ndarray

Basis matrices \(\hat O_i\) as an array of the shape (nspin, nspin, norb, norb, nsym).

nsym: int

Number of bath basis matrices \(\hat O_i\).

V: list[ndarray]

Hopping amplitudes \(V^\nu_{\sigma,\alpha}\). Each of the nbath elements is an array of the shape (nspin, norb) corresponding to the \(\nu\)-th replica.

l: list[ndarray]

Coefficients of linear combinations \(\lambda^\nu_i\). Each of the nbath elements is an array of length nsym corresponding to the \(\nu\)-th replica.

Bath fitting

BathFittingParams is a collection of parameters used in EDIpack’s bath fitting procedure. The procedure can be invoked by calling a method of edipack2triqs.solver.EDIpackSolver.

class edipack2triqs.fit.BathFittingParams(*, scheme: str = 'weiss', method: str = 'minimize', grad: str = 'numeric', tol: float = 1e-05, stop: str = 'both', niter: int = 500, n_iw: int = 1000, weight: str = '1', norm: str = 'elemental', pow: int = 2, minimize_ver: bool = False, minimize_hh: float = 0.0001)[source]

Parameters of bath fitting.

scheme: str = 'weiss'

Fitting scheme: “delta” to fit the hybridization function, “weiss” to fit the Weiss field.

method: str = 'minimize'

Minimization routine to use: “CGnr” for an algorithm from Numerical Recipes, “minimize” for an older FORTRAN 77 minimization procedure.

grad: str = 'numeric'

Gradient evaluation method, either “analytic” or “numeric”.

tol: float = 1e-05

Tolerance level for the conjugate gradient method.

stop: str = 'both'

Stopping condition for the conjugate gradient method.

  • “target”: \(|F_{n-1} - F_n| < tol (1+F_n)\).

  • “vars”: \(||x_{n-1} - x_n|| < tol (1+||x_n||)\).

  • “both”: Both conditions are fulfilled.

niter: int = 500

Maximal number of iterations.

n_iw: int = 1000

The number of Matsubara frequencies used in the fit.

weight: str = '1'

Weight function for the conjugate gradient minimization.

  • “1”: \(1\)

  • “1/n”: \(1/n\), where \(n\) is a Matsubara frequency number.

  • “1/w_n”: \(1/\omega_n\), where \(\omega_n\) is a Matsubara frequency.

norm: str = 'elemental'

Matrix norm to use in optimization, either “elemental” (sum of \(\chi^2\)-norms for each matrix element) or “frobenius”.

pow: int = 2

Fit power for the calculation of the generalized distance as \(|G_0 - G_{0,\text{and}}| ^ \text{pow}\).

minimize_ver: bool = False

Use the old/Krauth (False) or the new/Lichtenstein (True) version of the minimization conjugate gradient procedure. Only relevant for method=”minimize”.

minimize_hh: float = 0.0001

An unknown parameter used in the conjugate gradient minimization procedure.

EDIpackSolver.chi2_fit_bath(g: BlockGf, f: BlockGf | None = None)

Perform bath parameter fit of a given Green’s function.

Parameters:
  • g (triqs.gf.block_gf.BlockGf) – Normal component of the function to fit (either the hybridization function or the Weiss field).

  • f – Anomalous component of the function to fit (either the hybridization function or the Weiss field). Required iff the bath is superconducting.

Returns:

  • A bath object that contains resulting parameters of the fit.

  • The normal component of the fitted function.

  • (optional) The anomalous component of the fitted function.

Return type:

tuple[Bath, triqs.gf.block_gf.BlockGf] or tuple[Bath, triqs.gf.block_gf.BlockGf, triqs.gf.block_gf.BlockGf]