Lattice Green’s functions

triqs_tprf.lattice.lattice_dyson_g0_wk()

Construct a non-interacting Matsubara frequency lattice Green’s function \(G^{(0)}_{a\bar{b}}(i\omega_n, \mathbf{k})\)

Computes

\[G^{(0)}_{a\bar{b}}(i\omega_n, \mathbf{k}) = \left[ (i\omega_n + \mu ) \cdot \mathbf{1} - \epsilon(\mathbf{k}) \right]^{-1}_{a\bar{b}},\]

using a discretized dispersion \(\epsilon_{\bar{a}b}(\mathbf{k})\), chemical potential \(\mu\), and a Matsubara frequency Green’s function mesh.

Parameters:

mu :

chemical potential \(\mu\)

e_k

discretized lattice dispersion \(\epsilon_{\bar{a}b}(\mathbf{k})\)

mesh

imaginary frequency mesh

Returns:

out :

Matsubara frequency lattice Green’s function \(G^{(0)}_{a\bar{b}}(i\omega_n, \mathbf{k})\)

Signature : (float mu, triqs_tprf::e_k_cvt e_k, triqs::mesh::dlr_imfreq mesh) -> triqs_tprf::g_Dwk_t

triqs_tprf.lattice.lattice_dyson_g0_fk()

Signature : (float mu, triqs_tprf::e_k_cvt e_k, triqs::mesh::refreq mesh, float delta) -> triqs_tprf::g_fk_t Construct a non-interacting real frequency lattice Green’s function \(G^{(0)}_{a\bar{b}}(\omega, \mathbf{k})\)

Computes

\[G^{(0)}_{a\bar{b}}(\omega, \mathbf{k}) = \left[ (\omega + i\delta + \mu ) \cdot \mathbf{1} - \epsilon(\mathbf{k}) \right]^{-1}_{a\bar{b}},\]

using a discretized dispersion \(\epsilon_{\bar{a}b}(\mathbf{k})\), chemical potential \(\mu\), broadening \(\delta\), and a real frequency Green’s function mesh.

Parameters:

mu :

chemical potential \(\mu\)

e_k :

discretized lattice dispersion \(\epsilon_{\bar{a}b}(\mathbf{k})\)

mesh :

real frequency mesh

delta :

broadening \(\delta\)

Returns:

out :

Matsubara frequency lattice Green’s function \(G^{(0)}_{a\bar{b}}(\omega, \mathbf{k})\)

triqs_tprf.lattice.lattice_dyson_g_wk()

Construct an interacting Matsubara frequency lattice Green’s function \(G_{a\bar{b}}(i\omega_n, \mathbf{k})\)

Computes

\[G_{a\bar{b}}(i\omega_n, \mathbf{k}) = \left[ (i\omega_n + \mu ) \cdot \mathbf{1} - \epsilon(\mathbf{k}) - \Sigma(i\omega_n) \right]^{-1}_{a\bar{b}},\]

using a discretized dispersion \(\epsilon_{\bar{a}b}(\mathbf{k})\), chemical potential \(\mu\), and a momentum independent Matsubara frequency self energy \(\Sigma_{\bar{a}b}(i\omega_n)\).

Parameters:

mu :

chemical potential \(\mu\)

e_k

discretized lattice dispersion \(\epsilon_{\bar{a}b}(\mathbf{k})\)

sigma_wk

imaginary frequency self-energy \(\Sigma_{\bar{a}b}(i\omega_n, \mathbf{k})\)

Returns:

out :

Matsubara frequency lattice Green’s function \(G_{a\bar{b}}(i\omega_n, \mathbf{k})\)

Signature : (float mu, triqs_tprf::e_k_cvt e_k, triqs_tprf::g_Dwk_cvt sigma_wk) -> triqs_tprf::g_Dwk_t

triqs_tprf.lattice.lattice_dyson_g_fk()

Construct an interacting real frequency lattice Green’s function \(G_{a\bar{b}}(\omega, \mathbf{k})\)

Computes

\[G_{a\bar{b}}(\omega, \mathbf{k}) = \left[ (\omega + i\delta + \mu ) \cdot \mathbf{1} - \epsilon(\mathbf{k}) - \Sigma(\omega, \mathbf{k}) \right]^{-1}_{a\bar{b}},\]

using a discretized dispersion \(\epsilon_{\bar{a}b}(\mathbf{k})\), chemical potential \(\mu\), broadening \(\delta\), and a real frequency self energy \(\Sigma_{\bar{a}b}(\omega, \mathbf{k})\).

Parameters:

mu :

chemical potential \(\mu\)

e_k

discretized lattice dispersion \(\epsilon_{\bar{a}b}(\mathbf{k})\)

sigma_fk

real frequency self-energy \(\Sigma_{\bar{a}b}(\omega)\)

delta

broadening \(\delta\)

Returns:

out :

real frequency lattice Green’s function \(G_{a\bar{b}}(\omega, \mathbf{k})\)

triqs_tprf.lattice.lattice_dyson_g_f()

Signature : (float mu, triqs_tprf::e_k_cvt e_k, triqs_tprf::g_f_cvt sigma_f, float delta) -> triqs_tprf::g_f_t Construct an interacting real frequency local (\(\mathbf{r}=\mathbf{0}\)) lattice Green’s function \(G_{a\bar{b}}(\omega)\)

Computes

\[G_{a\bar{b}}(\omega) = \frac{1}{N_k} \sum_\mathbf{k} \left[ (i\omega_n + i\delta + \mu ) \cdot \mathbf{1} - \epsilon(\mathbf{k}) - \Sigma(\omega) \right]^{-1}_{a\bar{b}},\]

using a discretized dispersion \(\epsilon_{\bar{a}b}(\mathbf{k})\), chemical potential \(\mu\), and a momentum independent real frequency self energy \(\Sigma_{\bar{a}b}(\omega)\).

Parameters:

mu :

chemical potential \(\mu\)

e_k :

discretized lattice dispersion \(\epsilon_{\bar{a}b}(\mathbf{k})\)

sigma_f :

real frequency self-energy \(\Sigma_{\bar{a}b}(\omega)\)

delta :

broadening \(\delta\)

Returns:

out :

Real frequency lattice Green’s function \(G_{a\bar{b}}(\omega, \mathbf{k})\)

triqs_tprf.lattice.lattice_dyson_g_w()

Construct an interacting Matsubara frequency local (\(\mathbf{r}=\mathbf{0}\)) lattice Green’s function \(G_{a\bar{b}}(i\omega_n)\)

Computes

\[G_{a\bar{b}}(i\omega_n) = \frac{1}{N_k} \sum_\mathbf{k} \left[ (i\omega_n + \mu ) \cdot \mathbf{1} - \epsilon(\mathbf{k}) - \Sigma(i\omega_n) \right]^{-1}_{a\bar{b}},\]

using a discretized dispersion \(\epsilon_{\bar{a}b}(\mathbf{k})\), chemical potential \(\mu\), and a momentum independent Matsubara frequency self energy \(\Sigma_{\bar{a}b}(i\omega_n)\).

Parameters:

mu :

chemical potential \(\mu\)

e_k

discretized lattice dispersion \(\epsilon_{\bar{a}b}(\mathbf{k})\)

sigma_w

imaginary frequency self-energy \(\Sigma_{\bar{a}b}(i\omega_n)\)

Returns:

out :

Matsubara frequency lattice Green’s function \(G_{a\bar{b}}(i\omega_n, \mathbf{k})\)

Signature : (float mu, triqs_tprf::e_k_cvt e_k, triqs_tprf::g_Dw_cvt sigma_w) -> triqs_tprf::g_Dw_t

Non-interacting generalized susceptibility

triqs_tprf.lattice_utils.imtime_bubble_chi0_wk(g_wk, nw=1, save_memory=False, verbose=True)[source]
triqs_tprf.lattice.lindhard_chi00()

Generalized Lindhard susceptibility in the particle-hole channel \(\chi^{(00)}_{\bar{a}b\bar{c}d}(i\omega_n, \mathbf{q})\).

Analytic calculation of the generalized (non-interacting) Lindhard susceptibility in the particle-hole channel. The analytic expression is obtained using residue calculus to explicitly evaluate the matsubara sum of the fourier transformed imaginary time bubble product of two non-interacting single-particle Green’s functions.

\[G^{(0)}_{a\bar{b}}(\mathbf{k}, i\omega_n) = \left[ i\omega_n \cdot \mathbf{1} - \epsilon(\mathbf{k}) \right]^{-1} .\]

The analytic evaluation of the bubble diagram gives

\[\begin{split}\chi^{(00)}_{\bar{a}b\bar{c}d}(i\omega_n, \mathbf{q}) \equiv \mathcal{F} \left\{ - G^{(0)}_{d\bar{a}}(\tau, \mathbf{r}) G^{(0)}_{b\bar{c}}(-\tau, -\mathbf{r}) \right\} = - \frac{1}{N_k} \sum_{\nu} \sum_{\mathbf{k}} G^{(0)}_{d\bar{a}}(\nu, \mathbf{k}) G^{(0)}_{b\bar{c}}(\nu + \omega, \mathbf{k} + \mathbf{q}) \\ = - \frac{1}{N_k} \sum_{\nu} \sum_{\mathbf{k}} \left( \sum_{i} U^\dagger_{di}(\mathbf{k}) \frac{1}{i\nu - \epsilon_{\mathbf{k}, i}} U_{i\bar{a}}(\mathbf{k}) \right) \left( \sum_j U^\dagger_{bj}(\mathbf{k} + \mathbf{q}) \frac{1}{i\nu + i\omega - \epsilon_{\mathbf{k} + \mathbf{q}, j}} U_{j\bar{c}}(\mathbf{k} + \mathbf{q}) \right) \\ = \frac{1}{N_k} \sum_{\mathbf{k}} \sum_{ij} \left( [1 - \delta_{0, \omega_n} \delta_{\epsilon_{\mathbf{k},i},\epsilon_{\mathbf{k}+\mathbf{q}, j}})] \frac{ f(\epsilon_{\mathbf{k}, i}) - f(\epsilon_{\mathbf{k}+\mathbf{q}, j}) } {i\omega_n + \epsilon_{\mathbf{k} + \mathbf{q}, j} - \epsilon_{\mathbf{k}, i}} + \delta_{0, \omega_n} \delta_{\epsilon_{\mathbf{k},i},\epsilon_{\mathbf{k}+\mathbf{q}, j}} \frac{\beta}{4 \cosh^2 (\beta \epsilon_{\mathbf{k}, i} / 2) } \right) \\ \times U_{\bar{a}i}(\mathbf{k}) U^\dagger_{id}(\mathbf{k}) U_{\bar{c}j}(\mathbf{k} + \mathbf{q}) U^\dagger_{jb}(\mathbf{k} + \mathbf{q})\end{split}\]

where the \(U(\mathbf{k})\) matrices are the diagonalizing unitary transform of the matrix valued dispersion relation \(\epsilon_{\bar{a}b}(\mathbf{k})\), i.e.

Parameters:

e_k :

discretized lattice dispersion \(\epsilon_{\bar{a}b}(\mathbf{k})\)

mesh

real frequency mesh

beta

inverse temperature

mu

chemical potential \(\mu\)

delta

broadening \(\delta\)

Returns:

out :

real frequency generalized Lindhard susceptibility in the particle-hole channel \(\chi^{(00)}_{\bar{a}b\bar{c}d}(\omega, \mathbf{q})\)

Random Phase Approximation

triqs_tprf.lattice.solve_rpa_PH()

Random Phase Approximation (RPA) in the particle-hole channel

Computes the equation

\[\chi(\bar{a}b\bar{c}d) = \big( \mathbb{1} - \chi^{(0)}(\bar{a}b\bar{B}A) U(A\bar{B}D\bar{C}) \big)^{-1} \chi^{(0)}(\bar{C}D\bar{c}d)\,.\]
Parameters:

chi0 :

bare particle-hole bubble \(\chi^{(0)}_{\bar{a}b\bar{c}d}(\mathbf{k}, \omega)\)

U

RPA static vertex as obtained from triqs_tprf.rpa_tensor.get_rpa_tensor \(U_{a\bar{b}c\bar{d}}\)

Returns:

out :

RPA suceptibility \(\chi_{\bar{a}b\bar{c}d}(\mathbf{k}, \omega)\)

triqs_tprf.rpa_tensor.kanamori_quartic_tensor(norb, U, Up, J, Jp)[source]

Return Kanamori interaction as a quartic tensor

\[\begin{split}\hat{U}_{\text { Kanamori }} = U \sum_{i} \hat{n}_{i, \uparrow} \hat{n}_{i, \downarrow}+ \sum_{i>j, s, s^{\prime}}\left(U^{\prime}-J \delta_{\sigma, \sigma^{\prime}}\right) \hat{n}_{i, \sigma} \hat{n}_{j, \sigma^{\prime}} - \\ J \sum_{i \neq j}\left(\hat{c}_{i, \downarrow}^{\dagger} \hat{c}_{j, \uparrow}^{\dagger} \hat{c}_{j, \downarrow} \hat{c}_{i, \uparrow}+\hat{c}_{j, \uparrow}^{\dagger} \hat{c}_{j, \downarrow}^{\dagger} \hat{c}_{i, \uparrow} \hat{c}_{i, \downarrow}+ \mathrm{h.c.}\right)\end{split}\]
Parameters:

norb : int,

Number of orbitals excluding spin.

U : complex,

Strength of intra-orbital interaction.

Up : complex,

Strength of inter-orbital interaction.

J : complex,

Strength of Hound’s coupling.

Jp : complex,

Strength pair hopping and spin-flip.

Returns:

U : np.ndarray,

shape = (2*norb, 2*norb, 2*norb, 2*norb)

Impurity susceptibility and Bethe-Salpeter Equation

triqs_tprf.chi_from_gg2.chi0_from_gg2_PH()

Signature : (triqs_tprf::g_iw_vt g, triqs_tprf::g2_iw_vt g2) -> triqs_tprf::g2_iw_t Bubble susceptibility \(\chi^{(0)} = GG\) in the Particle-Hole channel

Computes

\[\chi^{(0)}_{\bar{a}b\bar{c}d}(\omega, \nu, \nu') = - \beta \delta_{\nu, \nu'} G_{da}(\nu) \cdot G_{bc}(\omega + \nu)\]
Parameters:

g :

single particle Green’s function \(G_{ab}(\nu)\)

g2 :

two-particle Green’s function with the mesh to use for \(\chi^{(0)}\)

Returns:

out :

chi0 particle-hole bubble \(\chi^{(0)}_{\bar{a}b\bar{c}d}(\omega, \nu,\nu')\)

triqs_tprf.chi_from_gg2.chi0_from_gg2_PP()

Signature : (triqs_tprf::g_iw_vt g, triqs_tprf::g2_iw_vt g2) -> triqs_tprf::g2_iw_t Bubble susceptibility \(\chi^{(0)} = GG\) in the Particle-Particle

channel

Computes

\[\chi^{(0)}_{\bar{a}b\bar{c}d}(\omega, \nu, \nu') = - \beta \delta_{\nu, \nu'} G_{da}(\nu) \cdot G_{bc}(\omega - \nu)\]
Parameters:

g :

single particle Green’s function \(G_{ab}(\nu)\)

g2 :

two-particle Green’s function with the mesh to use for \(\chi^{(0)}\)

Returns:

out :

chi0 particle-particle bubble \(\chi^{(0)}_{\bar{a}b\bar{c}d}(\omega, \nu,\nu')\)

triqs_tprf.chi_from_gg2.chi_from_gg2_PH()

Signature : (triqs_tprf::g_iw_vt g, triqs_tprf::g2_iw_vt g2) -> triqs_tprf::g2_iw_t Generalized susceptibility \(\chi^{(0)} = G^{(2)} - GG\) in the

Particle-Hole channel

Computes

\[\chi_{\bar{a}b\bar{c}d}(\omega, \nu, \nu') = G^{(2)}_{\bar{a}b\bar{c}d}(\omega, \nu, \nu') - \beta \delta_{\omega} G_{ba}(\nu) \cdot G_{dc}(\nu')\]
Parameters:

g :

single particle Green’s function \(G_{ab}(\nu)\)

g2 :

two-particle Green’s function \(G^{(2)}_{\bar{a}b\bar{c}d}(\omega, \nu, \nu')\)

Returns:

out :

chi generalized particle-hole susceptibility \(\chi_{\bar{a}b\bar{c}d}(\omega, \nu,\nu')\)

triqs_tprf.chi_from_gg2.chi_from_gg2_PP()

Signature : (triqs_tprf::g_iw_vt g, triqs_tprf::g2_iw_vt g2) -> triqs_tprf::g2_iw_t Generalized susceptibility \(\chi^{(0)} = G^{(2)} - GG\) in the

Particle-Particle channel

Computes

\[\chi_{\bar{a}b\bar{c}d}(\omega, \nu, \nu') = G^{(2)}_{\bar{a}b\bar{c}d}(\omega, \nu, \nu') - \beta \delta_{\nu + \nu' - \omega} G_{ba}(\nu) \cdot G_{dc}(\nu')\]
Parameters:

g :

single particle Green’s function \(G_{ab}(\nu)\)

g2 :

two-particle Green’s function \(G^{(2)}_{\bar{a}b\bar{c}d}(\omega, \nu, \nu')\)

Returns:

out :

chi generalized particle-hole susceptibility \(\chi_{\bar{a}b\bar{c}d}(\omega, \nu,\nu')\)

triqs_tprf.bse.solve_local_bse(chi0_wnn, chi_wnn)[source]

Solve the Bethe-Salpeter equation for the local vertex function \(\Gamma_{abcd}(\omega, \nu, \nu')\).

Computes:

\[\Gamma_{abcd}(\omega, \nu, \nu') = [\chi^{(0)}]^{-1} - \chi^{-1}\]

where the inverses are taken in the particle-hole channel pairing of fermionic frequencies \(\nu\) and \(\nu'\) and orbital indices.

Parameters:

chi0_wnn : Gerealized local bubble susceptibility

\(\chi^{(0)}_{abcd}(\omega, \nu, \nu')\)

chi_wnn : Generalized local susceptibility

\(\chi_{abcd}(\omega, \nu, \nu')\)

Returns:

gamma_wnn : Particle-hole vertex function

\(\Gamma_{abcd}(\omega, \nu, \nu')\)

Lattice Bethe-Salpeter Equation

triqs_tprf.bse.solve_lattice_bse(g_wk, gamma_wnn)[source]

Compute the generalized lattice susceptibility \(\chi_{\bar{a}b\bar{c}d}(\mathbf{k}, \omega_n)\) using the Bethe-Salpeter equation (BSE).

Parameters:

g_wk : Gf,

Single-particle Green’s function \(G_{a\bar{b}}(i\nu_n, \mathbf{k})\).

gamma_wnn : Gf,

Local particle-hole vertex function \(\Gamma_{a\bar{b}c\bar{d}}(i\omega_n, i\nu_n, i\nu_n')\).

Returns:

chi_kw : Gf,

Generalized lattice susceptibility \(\chi_{\bar{a}b\bar{c}d}(\mathbf{k}, i\omega_n)\).

chi0_kw : Gf,

Generalized bare lattice susceptibility \(\chi^0_{\bar{a}b\bar{c}d}(\mathbf{k}, i\omega_n)\).

triqs_tprf.bse.solve_lattice_bse_at_specific_w(g_wk, gamma_wnn, nw_index)[source]

Compute the generalized lattice susceptibility \(\chi_{\bar{a}b\bar{c}d}(i\omega_{n=\mathrm{nw\_index}}, \mathbf{k})\) using the Bethe-Salpeter equation (BSE) for a specific \(i\omega_{n=\mathrm{nw\_index}}\).

Parameters:

g_wk : Gf,

Single-particle Green’s function \(G_{a\bar{b}}(i\nu_n, \mathbf{k})\).

gamma_wnn : Gf,

Local particle-hole vertex function \(\Gamma_{a\bar{b}c\bar{d}}(i\omega_n, i\nu_n, i\nu_n')\).

nw_index : int,

The bosonic Matsubara frequency index \(i\omega_{n=\mathrm{nw\_index}}\) at which the BSE is solved.

Returns:

chi_k : Gf,

Generalized lattice susceptibility \(\chi_{\bar{a}b\bar{c}d}(i\omega_{n=\mathrm{nw\_index}}, \mathbf{k})\).

chi0_k : Gf,

Generalized bare lattice susceptibility \(\chi^0_{\bar{a}b\bar{c}d}(i\omega_{n=\mathrm{nw\_index}}, \mathbf{k})\).

triqs_tprf.bse.get_chi0_wnk(g_wk, nw=1, nwf=None)[source]

Compute the generalized bare lattice susceptibility \(\chi^{0}_{\bar{a}b\bar{c}d}(i\omega_n, i\nu_n, \mathbf{k})\) from the single-particle Green’s function \(G_{a\bar{b}}(i\nu_n, \mathbf{k})\).

Parameters:

g_wk : Gf,

Single-particle Green’s function \(G_{a\bar{b}}(i\nu_n, \mathbf{k})\).

nw : int,

Number of bosonic frequencies in \(\chi^0\).

nwf : int,

Number of fermionic frequencies in \(\chi^0\).

Returns:

chi0_wnk : Gf,

Generalized bare lattice susceptibility \(\chi^{0}_{\bar{a}b\bar{c}d}(i\omega_n, i\nu_n, \mathbf{k})\).

triqs_tprf.bse.get_chi0_nk_at_specific_w(g_wk, nw_index=1, nwf=None)[source]

Compute the generalized bare lattice susceptibility \(\chi^{0}_{\bar{a}b\bar{c}d}(i\omega_{n=\mathrm{nw\_index}}, i\nu_n, \mathbf{k})\) from the single-particle Green’s function \(G_{a\bar{b}}(i\nu_n, \mathbf{k})\) for a specific \(i\omega_{n=\mathrm{nw\_index}}\).

Parameters:

g_wk : Gf,

Single-particle Green’s function \(G_{a\bar{b}}(i\nu_n, \mathbf{k})\).

nw_index : int,

The bosonic Matsubara frequency index \(i\omega_{n=\mathrm{nw\_index}}\) at which \(\chi^0\) is calculated.

nwf : int,

Number of fermionic frequencies in \(\chi^0\).

Returns:

chi0_nk : Gf,

Generalized bare lattice susceptibility \(\chi^{0}_{\bar{a}b\bar{c}d}(i\omega_{n=\mathrm{nw\_index}}, i\nu_n, \mathbf{k})\).

Dual Bethe-Salpeter Equation

triqs_tprf.dbse.solve_lattice_dbse(g_wk, F_wnn, L_wn, chi_imp_w)[source]

Compute the generalized lattice susceptibility \(\chi_{\bar{a}b\bar{c}d}(\mathbf{k}, \omega_n)\) using the dual Bethe-Salpeter equation (DBSE).

Parameters:

g_wk : Gf,

Single-particle Green’s function \(G_{a\bar{b}}(i\nu_n, \mathbf{k})\).

F_wnn : Gf,

Local particle-hole reducible vertex function \(F_{a\bar{b}c\bar{d}}(i\omega_n, i\nu_n, i\nu_n')\).

L_wn : Gf,

Local particle-hole reducible triangle vertex function \(L_{a\bar{b}c\bar{d}}(i\omega_n, i\nu_n)\).

chi_imp_w : Gf,

Generalized DMFT impurity susceptibility \(\chi_{a\bar{b}c\bar{d}}(i\omega_n)\).

Returns:

chi_kw : Gf,

Generalized lattice susceptibility \(\chi_{\bar{a}b\bar{c}d}(\mathbf{k}, i\omega_n)\).

triqs_tprf.dbse.impurity_reducible_vertex_F(g_w, g2_wnn)[source]

Compute the impurity reducible vertex function \(F_{abcd}(\omega, \nu, \nu')\).

Computes:

\[F_{abcd}(\omega, \nu, \nu') = [\chi^{(0)}]^{-1} (\chi - \chi^{(0)} ) [\chi^{(0)}]^{-1} \]

where the inverses are taken in the particle-hole channel pairing of fermionic frequencies \(\nu\) and \(\nu'\) and orbital indices.

Parameters:

g_w : Single particle Green’s function

\(G_{ab}(\nu)\)

g2_wnn : Two-particle Green’s function

\(G^{(2)}_{abcd}(\omega, \nu, \nu')\)

Returns:

F_wnn : Particle-hole reducible vertex function

\(F_{abcd}(\omega, \nu, \nu')\)

GW approximation

class triqs_tprf.gw_solver.GWSolver(e_k, V_k, wmesh, mu=None, g_wk=None, N_fix=False, N_tol=1e-05, mu_bracket=None)[source]
triqs_tprf.lattice.dynamical_screened_interaction_W()

Dynamical screened interaction \(W(i\omega_n, \mathbf{k})\) calculator for static momentum-dependent bare interactions \(V(\mathbf{k})\).

The full screened interaction \(W(i\omega_n, \mathbf{k})\) is given by

\[W^{(full)}_{abcd}(i\omega_n, \mathbf{k}) = V_{abcd}(\mathbf{k}) + \sum_{efgh} V_{abef}(\mathbf{k}) \cdot \Pi_{fegh}(i\omega_n, \mathbf{k}) \cdot W^{(full)}_{hgcd}(i\omega_n, \mathbf{k})\]
Parameters:

PI_fk :

polarization bubble \(\Pi_{abcd}(\omega, \mathbf{k})\)

V_fk

dynamic bare interaction \(V_{abcd}(\omega, \mathbf{k})\)

Returns:

out :

dynamical screened interaction \(W_{abcd}(\omega, \mathbf{k})\)

triqs_tprf.lattice.dynamical_screened_interaction_W_from_generalized_susceptibility()

Dynamical screened interaction \(W(i\omega_n, \mathbf{k})\) calculator for static momentum-dependent bare interactions \(V(\mathbf{k})\) and known generalized susceptibility \(\chi(i\omega_n, \mathbf{k})\)

The full screened interaction \(W(i\omega_n, \mathbf{k})\) is given by

\[W^{(full)}_{abcd}(i\omega_n, \mathbf{k}) = V_{abcd}(\mathbf{k}) + \sum_{efgh} V_{abef}(\mathbf{k}) \cdot \chi_{fegh}(i\omega_n, \mathbf{k}) \cdot V_{hgcd}(\mathbf{k})\]
Parameters:

chi_fk :

generalized susceptibility \(\chi_{abcd}(\omega, \mathbf{k})\)

V_fk

dynamic bare interaction \(V_{abcd}(\omega, \mathbf{k})\)

Returns:

out :

dynamical screened interaction \(W_{abcd}(\omega, \mathbf{k})\)

triqs_tprf.gw.bubble_PI_wk(g_wk)[source]

Compute the particle-hole bubble from the single particle lattice Green’s function

\[\Pi_{abcd}(i\omega_n, k) = - \mathcal{F}_{\tau, \mathbf{r} \rightarrow i\omega_n, \mathbf{k}} \left\{ G_{d\bar{a}}(\tau, \mathbf{r}) G_{b\bar{c}}(-\tau, -\mathbf{r}) \right\}\]
Parameters:

g_wk : TRIQS Green’s function (rank 2) on Matsubara and Brillouinzone meshes

Single particle lattice Green’s function.

Returns:

PI_wk : TRIQS Green’s function (rank 4) on Matsubara and Brillouinzone meshes

Particle hole bubble.

triqs_tprf.gw.gw_sigma()

GW self energy \(\Sigma(i\omega_n, \mathbf{k})\) calculator for dynamic interactions

Splits the interaction into a dynamic and a static part

\[W_{abcd}(i\omega_n, \mathbf{k}) = W^{(dyn)}_{abcd}(i\omega_n, \mathbf{k}) + V_{abcd}(\mathbf{k})\]

by fitting the high-frequency tail.

Fourier transforms the dynamic part of the interaction and the single-particle Green’s function to imaginary time and real space.

\[G_{ab}(\tau, \mathbf{r}) = \mathcal{F}^{-1} \left\{ G_{ab}(i\omega_n, \mathbf{k}) \right\}\]
\[W^{(dyn)}_{abcd}(\tau, \mathbf{r}) = \mathcal{F}^{-1} \left\{ W^{(dyn)}_{abcd}(i\omega_n, \mathbf{k}) \right\}\]

computes the GW self-energy as the product

\[\Sigma^{(dyn)}_{ab}(\tau, \mathbf{r}) = - \sum_{cd} W^{(dyn)}_{abcd}(\tau, \mathbf{r}) G_{cd}(\tau, \mathbf{r})\]

and transforms back to frequency and momentum

\[\Sigma^{(dyn)}_{ab}(i\omega_n, \mathbf{k}) = \mathcal{F} \left\{ \Sigma^{(dyn)}_{ab}(\tau, \mathbf{r}) \right\}\]

The self-energy of the static part of the interaction is calculated as the sum

\[\Sigma^{(stat)}_{ab}(\mathbf{k}) = -\frac{1}{N_k} \sum_{\mathbf{q}} V_{abab}(\mathbf{k}) \rho(G(i\omega_n, \mathbf{k+q}))_{ab}\]

where \(\rho(G(i\omega_n, \mathbf{k+q}))\) is the density matrix of the single particle Green’s function.

The total GW self-energy is given by

\[\Sigma_{ab}(i\omega_n, \mathbf{k}) = \Sigma^{(dyn)}_{ab}(i\omega_n, \mathbf{k}) + \Sigma^{(stat)}_{ab}(\mathbf{k})\]
Parameters:

V_k :

static interaction \(V_{abcd}(\mathbf{k})\)

g_wk

single particle Green’s function \(G_{ab}(i\omega_n, \mathbf{k})\)

Returns:

out :

Static GW self-energy (Fock) \(\Sigma_{ab}(\mathbf{k})\)

triqs_tprf.gw.g0w_sigma()

Add some docs

Signature(float mu, float beta, triqs_tprf::e_k_cvt e_k, triqs_tprf::chi_k_cvt v_k, mesh::brzone::value_t kpoint) -> array<std::complex<double>, 2>

Add some docs

Signature(float mu, float beta, triqs_tprf::e_k_cvt e_k, triqs_tprf::chi_k_cvt v_k, mesh::brzone kmesh) -> triqs_tprf::e_k_t

Add some docs

Signature(float mu, float beta, triqs_tprf::e_k_cvt e_k, triqs_tprf::chi_k_cvt v_k) -> triqs_tprf::e_k_t

GW self energy \(\Sigma(\mathbf{k})\) calculator for static interactions

Computes the GW self-energy of a static interaction as the product

\[\Sigma_{ab}(\mathbf{k}) = \frac{-1}{N_k} \sum_{\mathbf{q}} \sum_{l} U_{la}(\mathbf{k}+\mathbf{q}) U^\dagger_{bl}(\mathbf{k}+\mathbf{q}) V_{abab}(\mathbf{q}) f(\epsilon_{\mathbf{k}+\mathbf{q}, l})\]

where the \(U(\mathbf{k})\) matrices are the diagonalizing unitary transform of the matrix valued dispersion relation \(\epsilon_{\bar{a}b}(\mathbf{k})\), i.e.

\[\sum_{\bar{a}b} U_{i\bar{a}}(\mathbf{k}) \epsilon_{\bar{a}b}(\mathbf{k}) U^\dagger_{bj} (\mathbf{k}) = \delta_{ij} \epsilon_{\mathbf{k}, i}\]
Parameters:

mu :

chemical potential \(\mu\)

beta

inverse temperature

e_k

discretized lattice dispersion \(\epsilon_{\bar{a}b}(\mathbf{k})\)

W_fk

fully screened interaction \(W_{abcd}(\omega, \mathbf{k})\)

V_k

bare interaction \(V_{abcd}(\mathbf{k})\)

delta

broadening \(\delta\)

Returns:

out :

real frequency GW self-energy \(\Sigma_{ab}(\omega, \mathbf{k})\)

Linearized Eliashberg equation

triqs_tprf.eliashberg.solve_eliashberg(Gamma_pp_wk, g_wk, initial_delta=None, Gamma_pp_const_k=None, tol=1e-10, product='FFT', solver='IRAM', symmetrize_fct=<function <lambda>>, k=6)[source]

Solve the linearized Eliashberg equation

Returns the biggest eigenvalues and corresponding eigenvectors of the linearized Eliashberg equation, for a particle-particle vertex in the random phase approximation, as described here Random phase approximation for the irreducible particle-particle vertex. The Eliashberg equation implementation is using fourier transformations for computational efficiency. The eigenvalues are found using an iterative algorithm from scipy.

Parameters:

Gamma_pp_wk : Gf,

Pairing vertex \(\Gamma(i\omega_n, \mathbf{k})\). The mesh attribute of the Gf must be a MeshProduct with the components (MeshImFreq, MeshBrZone).

g_wk : Gf,

Green’s function \(G(i\nu_n, \mathbf{k})\). The mesh attribute of the Gf must be a MeshProduct with the components (MeshImFreq, MeshBrZone).

initial_delta : Gf, optional

An initial anomalous self-energy \(\Delta(i\nu_n, \mathbf{k})\) to start an iterative solver, given as a Gf with MeshProduct with the components (MeshImFreq, MeshBrZone). If not given semi_random_initial_delta() will be called.

Gamma_pp_const_k : float or np.ndarray or Gf, optional

Part of the pairing vertex that is constant in Matsubara frequency space \(\Gamma(\mathbf{k})\). If given as a Gf its mesh attribute needs to be a MeshBrZone. If not given, the constant part will be fitted.

tol : float, optional

Relative accuracy for eigenvalues (stopping criterion).

product : str, [‘FFT’, ‘SUM’], optional

Which function of the Eliashberg product shall be used:

‘FFT’triqs_tprf.lattice.eliashberg_product_fft,

which uses Fourier transformation for optimal computational efficiency.

‘SUM’triqs_tprf.lattice.eliashberg_product, uses the explicit sum.

Restrictions : wmesh of Gamma_pp_wk must be atleast twice the size of the one of g_wk.

solver : str, [‘IRAM’, ‘PM’], optional

Which eigenvalue solver shall be used:

‘IRAM’ : Use the Implicitly Restarted Arnoldi Method implemented in implicitly_restarted_arnoldi_method().

‘PM’ : Use the Power Method implemented in power_method_LR().

symmetrize_fct : function, optional

A function that takes one parameter: A Green’s function \(G(i\nu_n, \mathbf{k})\). The mesh attribute of the Gf must be a MeshProduct with the components (MeshImFreq, MeshBrZone). This function is applied after every iteration of the eigenvalue solver and can be used to enforce a specific symmetry. If no symmetries are enforced, caution is need, because unphysical symmetries can occur.

k : int, optional

The number of leading superconducting gaps that shall be calculated. Does only have an effect, if ‘IRAM’ is used as a solver.

Returns:

Es : list of float,

Biggest eigenvalues of the linearized Eliashberg equation \(\lambda\).

eigen_modes : list of Gf,

Corresponding eigenvectors (anomalous self-energies) \(\Delta(i\nu_n, \mathbf{k})\) as Gf with MeshProduct with the components (MeshImFreq, MeshBrZone).

See also

Linearized Eliashberg Equation

Theory of the linearized Eliashberg equation.

triqs_tprf.eliashberg.preprocess_gamma_for_fft(Gamma_pp_wk, Gamma_pp_const_k=None)[source]

Prepare Gamma to be used with the FFT implementation

Parameters:

Gamma_pp_wk : Gf,

Pairing vertex \(\Gamma(i\omega_n, \mathbf{k})\). The mesh attribute of the Gf must be a MeshProduct with the components (MeshImFreq, MeshBrZone).

Gamma_pp_const_k : float or np.ndarray or Gf

Part of the pairing vertex that is constant in Matsubara frequency space \(\Gamma(\mathbf{k})\). If given as a Gf its mesh attribute needs to be a MeshBrZone. If not given, the constant part will be fitted.

Returns:

Gamma_pp_dyn_tr : Gf,

The dynamic part of Gamma, which converges to zero for \(\omega_n \rightarrow \infty\), but now in \(\tau\)-space. Its mesh attribute is MeshProduct with the components (MeshImTime, MeshCycLat).

Gamma_pp_const_r : Gf,

The constant part of Gamma with mesh attribute MeshCycLat.

triqs_tprf.eliashberg.semi_random_initial_delta(g_wk, nr_factor=0.5, seed=None)[source]

Create a delta based on the GF with random elements

Returns an anomalous self-energy that can be used as an inital input for the iterative solvers. The momentum space is random, while the Matsubara space is only partialy randomized to ensure working tail fits for the Fourier transformations.

Parameters:

g_wk : Gf,

Green’s function \(G(i\nu_n, \mathbf{k})\). The mesh attribute of the Gf must be a MeshProduct with the components (MeshImFreq, MeshBrZone).

nr_factor : float, optional

Percentage of \(\omega\) points which shall not be randomized. This is needed to assure a working tail fit for the Fourier transformations. The default is 0.5, meaning that 50% of the \(\omega\) points will not be randomized.

seed : int, optional

Set a np.random.seed to enforce predictable results.

Returns:

delta : Gf,

An initial anomalous self-energy \(\Delta(i\nu_n, \mathbf{k})\) to start an iterative solver, given as a Gf with MeshProduct with the components (MeshImFreq, MeshBrZone).

triqs_tprf.eliashberg.power_method_LR(matvec, init, tol=1e-10, max_it=100000.0)[source]

Find the eigenvalue with the largest real value via the power method

Parameters:

matvec : callable f(v),

Returns A*v.

init : np.ndarray,

The array representation of the anomalous self-energy to start the iterative method with. Restriction: len(init.shape) == 1.

tol : float, optional

The tolerance at which the iterative scheme is considered to be converged.

max_it : float, optional

The maximum number of iterations that shall be done before a error is raised.

Returns:

norm : float,

The eigenvalue with the largest positive real part.

v_k : np.ndarray,

The corresponding eigenvector.

triqs_tprf.eliashberg.implicitly_restarted_arnoldi_method(matvec, init, tol=1e-10, k=6)[source]

Find the eigenvalue with the largest real value via the Implicitly Restarted Arnoldi Method

Parameters:

matvec : callable f(v),

Returns A*v.

init : np.ndarray,

The array representation of the anomalous self-energy to start the iterative method with. Restriction: len(init.shape) == 1.

tol : float, optional

The tolerance at which the iterative scheme is considered to be converged.

k : int, optional

The number of eigenvalues and eigenvectors desired.

Returns:

Es : list of float,

The eigenvalues with the largest positive real part.

U : list of np.ndarray,

The corresponding eigenvectors.

Notes

scipy.sparse.linalg.eigs

scipy.sparse.linalg.LinearOperator

triqs_tprf.eliashberg.construct_gamma_singlet_rpa(U_d, U_m, phi_d_wk, phi_m_wk)[source]

Construct the irreducible singlet vertex in the RPA limit

The irreducible singlet vertex in the random phase approximation limit for a symmetrized calculations of the Eliashberg equation is given by

\[\Gamma^{\text{s}}_{a\overline{b}c\overline{d}}(Q=0, K, K') \equiv \frac{1}{2}U_{a\overline{b}c\overline{d}}^{\mathrm{d}} + \frac{3}{2}U_{a\overline{b}c\overline{d}}^{\mathrm{m}} + \text{Complex Conjugate} \left[ 3 \Phi^{\text{m}}_{c\overline{b}a\overline{d}}(K-K') - \Phi^{\text{d}}_{c\overline{b}a\overline{d}}(K-K') \right] \,.\]

For more details see Linearized Eliashberg Equation.

Parameters:

U_d : np.ndarray,

The local static interaction in the density channel.

U_m : np.ndarray,

The local static interaction in the magnetic channel.

phi_d_wk : Gf,

The reducible ladder vertex in the density channel \(\Phi^{\mathrm{d}}(i\omega_n, \mathbf{q})\). The mesh attribute of the Gf must be a MeshProduct with the components (MeshImFreq, MeshBrZone).

phi_m_wk : Gf,

The reducible ladder vertex in the magnetic channel \(\Phi^{\mathrm{m}}(i\omega_n, \mathbf{q})\). The mesh attribute of the Gf must be a MeshProduct with the components (MeshImFreq, MeshBrZone).

Returns:

gamma_singlet : Gf,

The irreducible singlet vertex in the RPA limit for a symmetrized calculation of the Eliashberg equation \(\Gamma^{\mathrm{s}}(i\omega_n,\mathbf{q})\).

triqs_tprf.eliashberg.construct_gamma_triplet_rpa(U_d, U_m, phi_d_wk, phi_m_wk)[source]

Construct the irreducible triplet vertex in the RPA limit

The irreducible triplet vertex in the random phase approximation limit for a symmetrized calculations of the Eliashberg equation is given by

\[\Gamma^{\text{t}}_{a\overline{b}c\overline{d}}(Q=0, K, K') \equiv - \frac{1}{2}U_{a\overline{b}c\overline{d}}^{\mathrm{d}} + \frac{1}{2}U_{a\overline{b}c\overline{d}}^{\mathrm{m}} + \text{Complex Conjuate} \left[ - \Phi^{\text{m}}_{c\overline{b}a\overline{d}}(K-K') - \Phi^{\text{d}}_{c\overline{b}a\overline{d}}(K-K') \right] \,.\]

For more details see Linearized Eliashberg Equation.

Parameters:

U_d : np.ndarray,

The local static interaction in the density channel.

U_m : np.ndarray,

The local static interaction in the magnetic channel.

phi_d_wk : Gf,

The reducible ladder vertex in the density channel \(\Phi^{\mathrm{d}}(i\omega_n, \mathbf{q})\). The mesh attribute of the Gf must be a MeshProduct with the components (MeshImFreq, MeshBrZone).

phi_m_wk : Gf,

The reducible ladder vertex in the magnetic channel \(\Phi^{\mathrm{m}}(i\omega_n, \mathbf{q})\). The mesh attribute of the Gf must be a MeshProduct with the components (MeshImFreq, MeshBrZone).

Returns:

gamma_triplet : Gf,

The irreducible triplet vertex in the RPA limit for a symmetrized calculation of the Eliashberg equation \(\Gamma^{\mathrm{t}}(i\omega_n,\mathbf{q})\).

Hubbard atom analytic response functions

triqs_tprf.analytic_hubbard_atom.analytic_hubbard_atom(beta, U, nw, nwf, nwf_gf)[source]

Compute dynamical response functions for the Hubbard atom at half filling.

This function returns an object that contains the single-particle Greens function \(G(\omega)\), the magnetic two-particle generalized susceptibility \(\chi_m(\omega, \nu, \nu')\), and the corresponding bare bubble \(\chi^{(0)}_m(\omega, \nu, \nu')\), and the magnetic vertex function \(\Gamma_m(\omega, \nu, \nu')\).

This is implemented using analytical formulas from Thunstrom et al. [PRB 98, 235107 (2018)] please cite the paper if you use this function!

In particular this is one exact solution to the Bethe-Salpeter equation, that is the infinite matrix inverse problem:

\[\Gamma_m = [\chi^{(0)}_m]^{-1} - \chi_m^{-1}\]
Parameters:

beta : float

Inverse temperature.

U : float

Hubbard U interaction parameter.

nw : int

Number of bosonic Matsubara frequencies in the computed two-particle response functions.

nwf : int

Number of fermionic Matsubara frequencies in the computed two-particle response functions.

nwf_gf : int

Number of fermionic Matsubara frequencies in the computed single-particle Greens function.

Returns:

p : ParameterCollection

Object containing all the response functions and some other observables, p.G_iw, p.chi_m, p.chi0_m, p.gamma_m, p.Z, p.m2, p.chi_m_static.

Two-particle response function linear-algebra operations

triqs_tprf.linalg.inverse_PH()

Two-particle response-function inversion \([g]^{-1}\) in the particle-hole channel (PH).

The two-particle response function \(g_{abcd}(\omega, \nu, \nu')\) is cast to matrix form and inverted

\[[g]^{-1} = [ g_{\{\nu\alpha\beta\}, \{\nu'\gamma\delta\}}(\omega) ]^{-1}\]

where the mapping of target-space indices \(\{a, b, c, d \}\) to \(\{\alpha, \beta\}, \{\gamma, \delta\}\) is channel dependent.

Storage is allocated and the inverse is returned by value.

Parameters:

g :

two-particle response function to invert, \(g \equiv g_{abcd}(\omega, \nu, \nu')\)

Returns:

out :

\([g]^{-1}\) in the given channel

Signature : (triqs_tprf::g2_nn_vt g) -> triqs_tprf::g2_nn_t

triqs_tprf.linalg.inverse_PP()

Two-particle response-function inversion \([g]^{-1}\) in the particle-particle channel (PP).

The two-particle response function \(g_{abcd}(\omega, \nu, \nu')\) is cast to matrix form and inverted

\[[g]^{-1} = [ g_{\{\nu\alpha\beta\}, \{\nu'\gamma\delta\}}(\omega) ]^{-1}\]

where the mapping of target-space indices \(\{a, b, c, d \}\) to \(\{\alpha, \beta\}, \{\gamma, \delta\}\) is channel dependent.

Storage is allocated and the inverse is returned by value.

Parameters:

g :

two-particle response function to invert, \(g \equiv g_{abcd}(\omega, \nu, \nu')\)

Returns:

out :

\([g]^{-1}\) in the given channel

Signature : (triqs_tprf::g2_nn_vt g) -> triqs_tprf::g2_nn_t

triqs_tprf.linalg.inverse_PH_bar()

Two-particle response-function inversion \([g]^{-1}\) in the particle-hole-bar channel (PH-bar).

The two-particle response function \(g_{abcd}(\omega, \nu, \nu')\) is cast to matrix form and inverted

\[[g]^{-1} = [ g_{\{\nu\alpha\beta\}, \{\nu'\gamma\delta\}}(\omega) ]^{-1}\]

where the mapping of target-space indices \(\{a, b, c, d \}\) to \(\{\alpha, \beta\}, \{\gamma, \delta\}\) is channel dependent.

Storage is allocated and the inverse is returned by value.

Parameters:

g :

two-particle response function to invert, \(g \equiv g_{abcd}(\omega, \nu, \nu')\)

Returns:

out :

\([g]^{-1}\) in the given channel

Signature : (triqs_tprf::g2_nn_vt g) -> triqs_tprf::g2_nn_t

triqs_tprf.linalg.product_PH()

Two-particle response-function product \(A * B\) in the particle-hole channel (PH).

The two-particle response functions \(A \equiv A_{abcd}(\omega, \nu, \nu')\) and \(B \equiv B_{abcd}(\omega, \nu, \nu')\) are cast to matrix form and their product is computed

\[(A * B)_{\{\nu\alpha\beta\}, \{\nu'\gamma\delta\}}(\omega) = \sum_{\bar{\nu}ab} A_{\{\nu\alpha\beta\}, \{\bar{\nu}ab\}}(\omega) B_{\{\bar{\nu}ab\}, \{\nu'\gamma\delta\}}(\omega)\]

where the mapping of target-space indices \(\{a, b, c, d \}\) to \(\{\alpha, \beta\}, \{\gamma, \delta\}\) is channel dependent.

Storage is allocated and the product is returned by value.

Parameters:

A :

two-particle response function \(A \equiv A_{abcd}(\omega, \nu, \nu')\)

B

two-particle response function \(B \equiv A_{abcd}(\omega, \nu, \nu')\)

Returns:

out :

\((A * B)\) in the given channel

Signature : (triqs_tprf::g2_nn_vt A, triqs_tprf::g2_nn_vt B) -> triqs_tprf::g2_nn_t

triqs_tprf.linalg.product_PP()

Two-particle response-function product \(A * B\) in the particle-particle channel (PP).

The two-particle response functions \(A \equiv A_{abcd}(\omega, \nu, \nu')\) and \(B \equiv B_{abcd}(\omega, \nu, \nu')\) are cast to matrix form and their product is computed

\[(A * B)_{\{\nu\alpha\beta\}, \{\nu'\gamma\delta\}}(\omega) = \sum_{\bar{\nu}ab} A_{\{\nu\alpha\beta\}, \{\bar{\nu}ab\}}(\omega) B_{\{\bar{\nu}ab\}, \{\nu'\gamma\delta\}}(\omega)\]

where the mapping of target-space indices \(\{a, b, c, d \}\) to \(\{\alpha, \beta\}, \{\gamma, \delta\}\) is channel dependent.

Storage is allocated and the product is returned by value.

Parameters:

A :

two-particle response function \(A \equiv A_{abcd}(\omega, \nu, \nu')\)

B

two-particle response function \(B \equiv A_{abcd}(\omega, \nu, \nu')\)

Returns:

out :

\((A * B)\) in the given channel

Signature : (triqs_tprf::g2_nn_vt A, triqs_tprf::g2_nn_vt B) -> triqs_tprf::g2_nn_t

triqs_tprf.linalg.product_PH_bar()

Two-particle response-function product \(A * B\) in the particle-hole-bar channel (PH-bar).

The two-particle response functions \(A \equiv A_{abcd}(\omega, \nu, \nu')\) and \(B \equiv B_{abcd}(\omega, \nu, \nu')\) are cast to matrix form and their product is computed

\[(A * B)_{\{\nu\alpha\beta\}, \{\nu'\gamma\delta\}}(\omega) = \sum_{\bar{\nu}ab} A_{\{\nu\alpha\beta\}, \{\bar{\nu}ab\}}(\omega) B_{\{\bar{\nu}ab\}, \{\nu'\gamma\delta\}}(\omega)\]

where the mapping of target-space indices \(\{a, b, c, d \}\) to \(\{\alpha, \beta\}, \{\gamma, \delta\}\) is channel dependent.

Storage is allocated and the product is returned by value.

Parameters:

A :

two-particle response function \(A \equiv A_{abcd}(\omega, \nu, \nu')\)

B

two-particle response function \(B \equiv A_{abcd}(\omega, \nu, \nu')\)

Returns:

out :

\((A * B)\) in the given channel

Signature : (triqs_tprf::g2_nn_vt A, triqs_tprf::g2_nn_vt B) -> triqs_tprf::g2_nn_t

triqs_tprf.linalg.identity_PH()

Two-particle response-function identity operator \(\mathbf{1}\) in the particle-hole channel (PH).

Constructs the unity-operator in the given channel

\[\mathbf{1}_{abcd}(\omega,\nu,\nu') = \mathbf{1}_{\{\nu\alpha\beta\}, \{\nu'\gamma\delta\}}(\omega) \equiv \delta_{\nu\nu'} \delta_{\alpha\gamma} \delta_{\beta\delta}\]

where the mapping of target-space indices \(\{a, b, c, d \}\) to \(\{\alpha, \beta\}, \{\gamma, \delta\}\) is channel dependent.

Storage is allocated and the result is returned by value.

Parameters:

A :

two-particle response function \(A \equiv A_{abcd}(\omega, \nu, \nu')\) determinig the shape and size of the unity operator

Returns:

out :

the unity operator \(\mathbf{1}\), in the given channel

Signature : (triqs_tprf::g2_nn_vt g) -> triqs_tprf::g2_nn_t

triqs_tprf.linalg.identity_PP()

Two-particle response-function identity operator \(\mathbf{1}\) in the particle-particle channel (PP).

Constructs the unity-operator in the given channel

\[\mathbf{1}_{abcd}(\omega,\nu,\nu') = \mathbf{1}_{\{\nu\alpha\beta\}, \{\nu'\gamma\delta\}}(\omega) \equiv \delta_{\nu\nu'} \delta_{\alpha\gamma} \delta_{\beta\delta}\]

where the mapping of target-space indices \(\{a, b, c, d \}\) to \(\{\alpha, \beta\}, \{\gamma, \delta\}\) is channel dependent.

Storage is allocated and the result is returned by value.

Parameters:

A :

two-particle response function \(A \equiv A_{abcd}(\omega, \nu, \nu')\) determinig the shape and size of the unity operator

Returns:

out :

the unity operator \(\mathbf{1}\), in the given channel

Signature : (triqs_tprf::g2_nn_vt g) -> triqs_tprf::g2_nn_t

triqs_tprf.linalg.identity_PH_bar()

Two-particle response-function identity operator \(\mathbf{1}\) in the particle-hole-bar channel (PH-bar).

Constructs the unity-operator in the given channel

\[\mathbf{1}_{abcd}(\omega,\nu,\nu') = \mathbf{1}_{\{\nu\alpha\beta\}, \{\nu'\gamma\delta\}}(\omega) \equiv \delta_{\nu\nu'} \delta_{\alpha\gamma} \delta_{\beta\delta}\]

where the mapping of target-space indices \(\{a, b, c, d \}\) to \(\{\alpha, \beta\}, \{\gamma, \delta\}\) is channel dependent.

Storage is allocated and the result is returned by value.

Parameters:

A :

two-particle response function \(A \equiv A_{abcd}(\omega, \nu, \nu')\) determinig the shape and size of the unity operator

Returns:

out :

the unity operator \(\mathbf{1}\), in the given channel

Signature : (triqs_tprf::g2_nn_vt g) -> triqs_tprf::g2_nn_t

Wannier90 tight binding parsers

triqs_tprf.wannier90.parse_hopping_from_wannier90_hr_dat(filename)[source]

Wannier90 real space hopping parser of *_hr.dat files.

Returns a dictionary where the keys are the real-space hopping vectors, in terms of multiples of the lattice vectors, and the values are num_wann * num_wann numpy ndarrays with the hopping integrals.

Parameters:

filename : str

Wannier90 *_hr.dat file to parse.

Returns:

hopp_dict : dict

Dictionary of real space hoppings.

num_wann : int

Total number of Wannier functions per unit-cell.

triqs_tprf.wannier90.parse_lattice_vectors_from_wannier90_wout(filename)[source]

Wannier90 real space lattice vectors parser of *.wout files.

Parameters:

filename : str

Wannier90 *.wout file to parse.

Returns:

vectors : list of three three-tuples of floats

Lattice vectors.

triqs_tprf.wannier90.parse_reciprocal_lattice_vectors_from_wannier90_wout(filename)[source]

Wannier90 reciprocal space lattice vectors parser of *.wout files.

Parameters:

filename : str

Wannier90 *.wout file to parse.

Returns:

vectors : list of three three-tuples of floats

Reciprocal lattice vectors.

triqs_tprf.wannier90.parse_band_structure_from_wannier90_band_dat(filename)[source]

Wannier90 band structure parser of *_band.dat` files.

Parameters:

filename : str

Wannier90 *_band.dat file to parse.

Returns:

E : ndarray (2D)

Band energies.

w : ndarray (1D)

k-space path points.

Tight binding lattice model

class triqs_tprf.tight_binding.TBLattice(units, hoppings={}, orbital_positions=[(0, 0, 0)], orbital_names=None, hopping=None)[source]

Class describing a tight binding hamiltonian on top of a bravais lattice. Has objects of type BravaisLattice, BrillouinZone and TightBinding as attributes, and exposes part of their interfaces.

Attributes:

bl : BravaisLattice

The associated Bravais Lattice

bz : BrillouinZone

The associated Brillouin Zone

tb : TightBinding

The tight-binding Hamiltonian

dispersion(arg)[source]

Evaluate the dispersion relation for a momentum vector k in units of the reciprocal lattice vectors

Signature(k_cvt K) -> nda::array<double, 1>

Evaluate the dispersion relation for a momentum vector k in units of the reciprocal lattice vectors

Signature(nda::array_const_view<double,2> K) -> nda::array<double, 2>

Evaluate the dispersion relation for an array of momentum vectors k in units of the reciprocal lattice vectors

Signature(mesh::brzone k_mesh) -> gf<mesh::brzone, tensor_real_valued<1>>

Evaluate the dispersion relation on the k_mesh and return the associated Green-function object

fourier(arg)[source]

Evaluate the fourier transform for a momentum vector k in units of the reciprocal lattice vectors

Signature(k_cvt K) -> matrix<dcomplex>

Evaluate the fourier transform for a momentum vector k in units of the reciprocal lattice vectors

Signature(nda::array_const_view<double,2> K) -> nda::array<dcomplex, 3>

Evaluate the fourier transform for an array of momentum vectors k in units of the reciprocal lattice vectors

Signature(mesh::brzone k_mesh) -> gf<mesh::brzone, matrix_valued>

Evaluate the fourier transform on the k_mesh and return the associated Green-function object

get_kmesh(n_k)[source]

Return a mesh on the Brillouin zone with a given discretization

Parameters:

n_k : int or three-tuple of int

The linear dimension(s)

Returns:

MeshBrZone :

The mesh on the Brillouin zone

get_rmesh(n_r)[source]

Return a mesh on the Bravais lattice with a given periodicity

Parameters:

n_r : int or three-tuple of int

The periodicity in each dimension

Returns:

MeshCycLat :

The cyclic lattice mesh

lattice_to_real_coordinates(x)[source]

Signature : (r_t x) -> r_t Transform into real coordinates.

property n_orbitals

Number of orbitals in the unit cell

property ndim

Number of dimensions of the lattice

property orbital_names

The list of orbital names

property orbital_positions

The list of orbital positions

property units

Number of dimensions of the lattice

triqs_tprf.tight_binding.create_square_lattice(norb, t, tp=0.0, zeeman=0.0, spin=False, **kwargs)[source]

Retuns TBLattice that represents a model on a square lattice

The model is described by the Hamiltonian

\[H=-t \sum_{\langle j, l\rangle \sigma}\left(c_{j \sigma}^{\dagger} c_{l \sigma}+c_{l \sigma}^{\dagger} c_{j \sigma}\right) - t' \sum_{\langle\langle j, l\rangle\rangle \sigma}\left(c_{j \sigma}^{\dagger} c_{l \sigma}+c_{l \sigma}^{\dagger} c_{j \sigma}\right) + \xi \sum_j \left(n_{j \uparrow} - n_{j \downarrow} \right)\,,\]

where the angular bracket describes a sum over nearest neighbors and the double angular bracket over next-nearest neighbors.

Parameters:

norb : int,

Number of orbitals excluding spin

t : complex,

Kinetic energy of nearest neighbor hopping. Corresponds to \(t\).

tp : complex, optional

Kinetic energy of next-nearest neighbor hopping. Corresponds to \(t'\).

zeeman : complex, optional

Strength of Zeeman term. Corresponds to \(\xi\).

spin : bool,

True if spin index should be used explicitly, False otherwise. The Zeeman term can only be applied if spin=True.

Returns:

square_lattice : TBLattice

class triqs_tprf.super_lattice.TBSuperLattice(tb_lattice, super_lattice_units, cluster_sites=None, remove_internal_hoppings=False)[source]

Builds a superlattice on top of a base TBLattice.

Parameters:

tb_lattice : TBLattice instance

The base tight binding lattice.

super_lattice_units : ndarray (2D)

The unit vectors of the superlattice in the tb_lattice (integer) coordinates.

cluster_sites : :

Coordinates of the cluster in tb_lattice coordinates. If None, an automatic computation of cluster positions is made as follows: it takes all points whose coordinates in the basis of the superlattice are in [0, 1[^dimension.

remove_internal_hoppings : bool

If true, the hopping terms are removed inside the cluster. Useful to add Hartree Fock terms at the boundary of a cluster, e.g.

change_coordinates_L_to_SL(x)[source]

Given a point on the tb_lattice in lattice coordinates, returns its coordinates (R, alpha) in the Superlattice

change_coordinates_SL_to_L(R, alpha)[source]

Given a point in the supercell R, site (number) alpha, it computes its position on the tb_lattice in lattice coordinates

cluster_sites()[source]

Generate the position of the cluster site in the tb_lattice coordinates.

fold(D1, remove_internal=False, create_zero=None)[source]

Input: a function r-> f(r) on the tb_lattice given as a dictionnary Output: the function R-> F(R) folded on the superlattice. Only requirement is that f(r)[orbital1, orbital2] is properly defined. Hence f(r) can be a numpy, a GFBloc, etc…

pack_index_site_orbital(n_site, n_orbital)[source]

nsite and n_orbital must start at 0

unpack_index_site_orbital(index)[source]

Inverse of pack_index_site_orbital

Hartree-Fock and Hartree solvers

class triqs_tprf.hf_solver.HartreeFockSolver(e_k, beta, H_int=None, gf_struct=None, mu0=0.0, mu_max=None, mu_min=None)[source]

Hartree-Fock solver for local interactions.

Parameters:

e_k : TRIQS Greens function on a Brillouin zone mesh

Single-particle dispersion.

beta : float

inverse temperature

H_int : TRIQS Operator instance

Local interaction Hamiltonian

gf_struct : list of pairs of block index and its size

gf_struct fixing orbital order between e_k and H_int

mu0 : float

chemical potential

mu_min, mu_max : float, optional

range for chemical potential search.

mat2vec(mat)[source]

Converts a unitary matrix to a vector representation with the order

  1. the real diagonal entries

  2. the real part of the upper triangular entries

  3. the imaginary part of the upper triangular entries

solve_iter(N_target, M0=None, mu0=None, nitermax=100, mixing=0.5, tol=1e-09)[source]

Solve the HF-equations using forward recursion at fixed density.

Parameters:

N_target : float

Total density.

M0 : ndarray (2D), optional

Initial mean-field (0 if None).

mu0 : float, optional

Initial chemical potential.

nitermax : int, optional

Maximal number of self consistent iterations.

mixing : float, optional

Linear mixing parameter.

tol : float, optional

Convergence in relative change of the density matrix.

Returns:

rho : ndarray (2D)

Local density matrix.

solve_newton(N_target, M0=None, mu0=None)[source]

Solve the HF-equations using a quasi Newton method at fixed density.

Parameters:

N_tot : float

Total density.

M0 : ndarray (2D), optional

Initial mean-field (0 if None).

mu0 : float, optional

Initial chemical potential.

solve_newton_mu(mu, M0=None)[source]

Solve the HF-equations using a quasi Newton method at fixed chemical potential.

Parameters:

mu0 : float

Initial chemical potential.

M0 : ndarray (2D), optional

Initial mean-field (0 if None).

vec2mat(vec)[source]

Converts from vector representation to a unitary matrix see mat2vec(…) for details.

class triqs_tprf.hf_response.HartreeFockResponse(hartree_fock_solver, eps=1e-09)[source]

Hartree-Fock linear response calculator.

Parameters:

hartree_fock_solver : HartreeFockSolver instance

Converged Hartree-Fock solver.

eps : float

Step size in finite difference linear response calculation.

class triqs_tprf.hf_solver.HartreeSolver(e_k, beta, H_int=None, gf_struct=None, mu0=0.0, mu_max=None, mu_min=None)[source]

Hartree solver for local interactions.

Parameters:

e_k : TRIQS Greens function on a Brillouin zone mesh

Single-particle dispersion.

beta : float

inverse temperature

H_int : TRIQS Operator instance

Local interaction Hamiltonian

gf_struct : list of lists of block index and list of internal orbital indices

gf_struct fixing orbital order between e_k and H_int

mu0 : float

chemical potential

mu_min, mu_max : float, optional

range for chemical potential search.

mat2vec(mat)[source]

Converts a unitary matrix to a vector representation with the order

  1. the real diagonal entries

  2. the real part of the upper triangular entries

  3. the imaginary part of the upper triangular entries

vec2mat(vec)[source]

Converts from vector representation to a unitary matrix see mat2vec(…) for details.

class triqs_tprf.hf_response.HartreeResponse(hartree_solver, eps=1e-09)[source]

Hartree linear response calculator.

Parameters:

hartree_solver : HartreeSolver instance

Converged Hartree solver.

eps : float

Step size in finite difference linear response calculation.

Utility functions

triqs_tprf.lattice_utils.k_space_path(paths, num=100, bz=None, relative_coordinates=False)[source]

Generate an array of k-vectors along a path defined by a list of pairs of k-vectors

Parameters:

paths : list of pairs of three-vectors of floats

List of pairs of k-vectors in reciprocal units to create a path in-between.

num : int, default=100

Number of k-vectors along each segment of the overall path

bz : brillouin_zone, optional

When a Brillouin Zone is passed, calculate distance in absolute units

relative_coordinates : bool, optional

Return k-vectors in reciprocal units. (Default False)

Returns:

kvecs: numpy.ndarray [shape=(len(paths)*num,3)] :

Two-dimensional numpy array containing the path vectors (in reciprocal units) as rows

dist: numpy.ndarray [shape=(kvecs.shape[0])] :

One-dimensional numpy array containing, for each element in kvecs, the distance travelled along the path. Useful for plotting. If bz is provided, calculate the distance in absolute units. The distances for the relevant k-points in paths can be obtained with dist[num::num].

ticks : numpy.ndarray [shape=(len(paths)+1)]

Array with tick points, i.e. distances at the interfaces of path segments. Only returned if return_ticks is True.

triqs_tprf.lattice_utils.chi_contraction(chi, op1, op2)[source]

Contract a susceptibility with two operators

Parameters:

chi : Gf,

Susceptibility \(\chi(i\omega_n, \mathbf{k})\). The mesh attribute of the Gf must be a MeshProduct with the components (MeshImFreq, MeshBrZone) and its target_rank 4.

op1, op2 : np.ndarray,

Operators in matrix representation.

Returns:

Gf, :

Susceptibility :math:`chi(iomega_n, mathbf{k})`. With a target_rank of 0. :

Parameter collections

class triqs_tprf.ParameterCollection.ParameterCollection(**kwargs)[source]

Helper class for handing collections of parameters.

Parameters:

kwargs : dict

Key-word argument list of parameters.

Examples

A ParameterCollection has any number of attributes, accessible with the dot operator.

>>> p = ParameterCollection(beta=10., U=1.0, t=1.0)
>>> print p
U = 1.0
beta = 10.0
t = 1.0
>>> print p.beta
10.0
>>> p.W = 1.2
>>> print p
U = 1.0
W = 1.2
beta = 10.0
t = 1.0

and can be stored and loaded to/from TRIQS hdf archives.

>>> from h5 import HDFArchive
>>> with HDFArchive('data.h5', 'w') as arch: arch['p'] = p
>>> with HDFArchive('data.h5', 'r') as arch: p_ref = arch['p']
>>> print p_ref
U = 1.0
beta = 10.0
t = 1.0
alter(**kwargs)[source]

Change or add attributes

Returns:

p : ParameterCollection

convert_keys_from_string_to_python(dict_key)[source]

h5.HDFArchive incorrectly mangles tuple keys to string running this on the affected dict tries to revert this by running eval on the string representation. UGLY FIX…

copy()[source]

Shallow copy

class triqs_tprf.ParameterCollection.ParameterCollections(objects=[])[source]

Helper class for handing a series of collections of parameters.

Parameters:

objects : list

List of ParameterCollection instances.

Examples

The ParameterCollections class makes it easy to get vectors of parameters from a list of ParameterCollection objects.

>>> p1 = ParameterCollection(beta=10., U=1.0, t=1.0)
>>> p2 = ParameterCollection(beta=5., U=2.0, t=1.337)
>>> ps = ParameterCollections(objects=[p1, p2])
>>> print ps.beta
[10.  5.]
>>> print ps.U
[1. 2.]
triqs_tprf.ParameterCollection.parameter_scan(p, **kwargs)[source]

Return ParameterCollections with copies of ParameterCollection for different parameters

Uses a given ParameterCollection as a template to create copies of it with one or more parameters changing. Stores all of these copies in a ParameterCollections for easy access.

Parameters:

p : ParameterCollection,

The ParameterCollection that shall be used as a template for all the others

**kwargs : Sequence,

The keyword gives the parameter name and the Sequence the values that shall be scanned through.

Returns:

ParameterCollections :

Examples

>>> p = ParameterCollection(beta=10., U=1.0, t=1.0)
>>> ps = parameter_scan(p, U=[1.0, 1.5, 2.0])
>>> print ps[0]
U = 1.0
beta = 10.0
t = 1.0
>>> print ps[1]
U = 1.5
beta = 10.0
t = 1.0
>>> print ps[2]
U = 2.0
beta = 10.0
t = 1.0