Lattice Green’s functions

triqs_tprf.lattice.lattice_dyson_g0_wk()

Signature : (float mu, triqs_tprf::e_k_cvt e_k, gf_mesh<triqs::gfs::imfreq> mesh) -> triqs_tprf::g_wk_t Construct a non-interacting Matsubara frequency lattice Green’s function \(G^{(0)}_{a\bar{b}}(i\omega_n, \mathbf{k})\)

Computes

(1)\[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.

param mu:chemical potential \(\mu\)
param e_k:discretized lattice dispersion \(\epsilon_{\bar{a}b}(\mathbf{k})\)
param mesh:imaginary frequency mesh

@return Matsubara frequency lattice Green’s function \(G^{(0)}_{a\bar{b}}(i\omega_n, \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

(2)\[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)\).

param mu:chemical potential \(\mu\)
param e_k:discretized lattice dispersion \(\epsilon_{\bar{a}b}(\mathbf{k})\)
param sigma_w:imaginary frequency self-energy \(\Sigma_{\bar{a}b}(i\omega_n)\)

@return 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_w_cvt sigma_w) -> triqs_tprf::g_wk_t

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

Computes

(3)\[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)\).

param mu:chemical potential \(\mu\)
param e_k:discretized lattice dispersion \(\epsilon_{\bar{a}b}(\mathbf{k})\)
param sigma_w:imaginary frequency self-energy \(\Sigma_{\bar{a}b}(i\omega_n)\)

@return 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_wk_cvt sigma_wk) -> triqs_tprf::g_wk_t

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

Computes

(4)\[G_{a\bar{b}}(i\omega_n, \mathbf{k}) = \left[ (i\omega_n + \mu ) \cdot \mathbf{1} - \epsilon(\mathbf{k}) - \Sigma(i\omega_n, \mathbf{k}) \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, \mathbf{k})\).

param mu:chemical potential \(\mu\)
param e_k:discretized lattice dispersion \(\epsilon_{\bar{a}b}(\mathbf{k})\)
param sigma_wk:imaginary frequency self-energy \(\Sigma_{\bar{a}b}(i\omega_n, \mathbf{k})\)

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

triqs_tprf.lattice.lattice_dyson_g_w()

Signature : (float mu, triqs_tprf::e_k_cvt e_k, triqs_tprf::g_w_cvt sigma_w) -> triqs_tprf::g_w_t Construct an interacting Matsubara frequency local (\(\mathbf{r}=\mathbf{0}\)) lattice Green’s function \(G_{a\bar{b}}(i\omega_n)\)

Computes

(5)\[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)\).

param mu:chemical potential \(\mu\)
param e_k:discretized lattice dispersion \(\epsilon_{\bar{a}b}(\mathbf{k})\)
param sigma_w:imaginary frequency self-energy \(\Sigma_{\bar{a}b}(i\omega_n)\)

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

Lindhard non-interacting generalized susceptibility

triqs_tprf.lattice.lindhard_chi00_wk()

Signature : (triqs_tprf::e_k_cvt e_k, int nw, float beta, float mu) -> triqs_tprf::chi_wk_t 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.

(6)\[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

(7)\[\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_{i\bar{a}}(\mathbf{k}) U^\dagger_{di}(\mathbf{k}) U_{j\bar{c}}(\mathbf{k} + \mathbf{q}) U^\dagger_{bj}(\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.

(8)\[\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}\]

Note

The analytic formula is sub-optimal in terms of performance for higher temperatures. The evaluation scales as \(\mathcal{O}(N_k^2)\) which is worse than computing the bubble explicitly in imaginary time, with scaling \(\mathcal{O}(N_k N_\tau \log(N_k N_\tau)\) for \(N_k \gg N_\tau\).

Note

Care must be taken when evaluating the fermionic Matsubara frequency sum of the product of two simple poles. By extending the sum to an integral over the complex plane the standard expression for the Lindhard response is obtained when the poles are non-degenerate. The degenerate case produces an additional frequency independent contribution (the last term on the last row).

Random Phase Approximation

triqs_tprf.lattice.solve_rpa_PH()

Signature : (triqs_tprf::chi_wk_vt chi0, array_view<std::complex<double>,4> U) -> triqs_tprf::chi_wk_t Random Phase Approximation (RPA) in the particle-hole channel

Computes the equation

(9)\[\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)\,.\]
param chi0:bare particle-hole bubble \(\chi^{(0)}_{\bar{a}b\bar{c}d}(\mathbf{k}, i\omega_n)\)
param U:RPA static vertex as obtained from triqs_tprf.rpa_tensor.get_rpa_tensor \(U_{a\bar{b}c\bar{d}}\)

@return RPA suceptibility \(\chi_{\bar{a}b\bar{c}d}(\mathbf{k}, i\omega_n)\)

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

(10)\[\chi^{(0)}_{\bar{a}b\bar{c}d}(\omega, \nu, \nu') = - \beta \delta_{\nu, \nu'} G_{da}(\nu) \cdot G_{bc}(\omega + \nu)\]
param g:single particle Green’s function \(G_{ab}(\nu)\)
param g2:two-particle Green’s function with the mesh to use for
\(\chi^{(0)}\)
@return 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

(11)\[\chi^{(0)}_{\bar{a}b\bar{c}d}(\omega, \nu, \nu') = - \beta \delta_{\nu, \nu'} G_{da}(\nu) \cdot G_{bc}(\omega - \nu)\]
param g:single particle Green’s function \(G_{ab}(\nu)\)
param g2:two-particle Green’s function with the mesh to use for
\(\chi^{(0)}\)
@return 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

(12)\[\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')\]
param g:single particle Green’s function \(G_{ab}(\nu)\)
param g2:two-particle Green’s function
\(G^{(2)}_{\bar{a}b\bar{c}d}(\omega, \nu, \nu')\)
@return 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

(13)\[\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')\]
param g:single particle Green’s function \(G_{ab}(\nu)\)
param g2:two-particle Green’s function
\(G^{(2)}_{\bar{a}b\bar{c}d}(\omega, \nu, \nu')\)
@return 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:

(14)\[\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, tail_corr_nwf=None)[source]

Compute the generalized lattice susceptibility \(\chi_{abcd}(\omega, \mathbf{k})\) using the Bethe-Salpeter equation (BSE).

Parameters:

g_wk : Single-particle Green’s function \(G_{ab}(\omega, \mathbf{k})\).

gamma_wnn : Local particle-hole vertex function

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

tail_corr_nwf : Number of fermionic freqiencies to use in the

tail correction of the sum over fermionic frequencies.

Returns:

chi0_wk : Generalized lattice susceptibility

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

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

Compute the generalized lattice bubble susceptibility \(\chi^{(0)}_{abcd}(\omega, \nu, \mathbf{k})\) from the single-particle Green’s function \(G_{ab}(\omega, \mathbf{k})\).

Parameters:

g_wk : Single-particle Green’s function \(G_{ab}(\omega, \mathbf{k})\).

nw : Number of bosonic freqiencies in \(\chi\).

nwf : Number of fermionic freqiencies in \(\chi\).

Returns:

chi0_wnk : Generalized lattice bubble susceptibility

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

GW approximation

triqs_tprf.lattice.dynamical_screened_interaction_W_wk()

Signature : (triqs_tprf::chi_wk_cvt PI_wk, triqs_tprf::chi_k_cvt V_k) -> triqs_tprf::chi_wk_t Dynamical screened interaction \(W(i\omega_n, \mathbf{k})\) calculator

for static momentum-dependent interactions \(V(\mathbf{k})\).

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

(15)\[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})\]

Instead of returning \(W^{(full)}\) we return the dynamical/retarded part \(W^{(r)}\) (with zero high-frequency offset)

(16)\[W_{abcd}(i\omega_n, \mathbf{k}) = W^{(full)}_{abcd}(i\omega_n, \mathbf{k}) - V_{abcd}(\mathbf{k})\]
param PI_wk:polarization bubble \(\Pi_{abcd}(i\omega_n, \mathbf{k})\)
param V_k:static interaction \(V_{abcd}(\mathbf{k})\)

@return dynamical screened interaction \(W_{abcd}(i\omega_n, \mathbf{k})\)

triqs_tprf.lattice.dynamical_screened_interaction_W_wk_from_generalized_susceptibility()

Signature : (triqs_tprf::chi_wk_cvt chi_wk, triqs_tprf::chi_k_cvt V_k) -> triqs_tprf::chi_wk_t Dynamical screened interaction \(W(i\omega_n, \mathbf{k})\) calculator

for static momentum-dependent 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

(17)\[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})\]

Instead of returning \(W^{(full)}\) we return the dynamical/retarded part \(W^{(r)}\) (with zero high-frequency offset)

(18)\[W_{abcd}(i\omega_n, \mathbf{k}) = W^{(full)}_{abcd}(i\omega_n, \mathbf{k}) - V_{abcd}(\mathbf{k})\]
param chi_wk:polarization bubble \(\Pi_{abcd}(i\omega_n, \mathbf{k})\)
param V_k:static interaction \(V_{abcd}(\mathbf{k})\)

@return dynamical screened interaction \(W_{abcd}(i\omega_n, \mathbf{k})\)

triqs_tprf.gw.bubble_PI_wk(g_wk)[source]

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

(19)\[\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_wk(Wr_wk, g_wk, fft_flag=False)[source]

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

Fourier transforms the screened interaction and the single-particle Green’s function to imagiary time and real space.

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

computes the GW self-energy as the product

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

and transforms back to frequency and momentum

(23)\[\Sigma_{ab}(i\omega_n, \mathbf{k}) = \mathcal{F} \left\{ \Sigma_{ab}(\tau, \mathbf{r}) \right\}\]
Parameters:

V_k : TRIQS Green’s function (rank 4) on a Brillouinzone mesh

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

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

retarded screened interaction \(W^{(r)}_{abcd}(i\omega_n, \mathbf{k})\)

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

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

Returns:

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

GW self-energy \(\Sigma_{ab}(i\omega_n, \mathbf{k})\)

Linearized Eliashberg equation

triqs_tprf.eliashberg.solve_eliashberg_fft(Gamma_pp_wk, g_wk, Gamma_pp_const_k=None, tol=1e-10)[source]

Solve the linearized Eliashberg equation

Returns the biggest eigenvalues and corresponding eigenvectors of the linearized Eliashberg equation. 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, MeshBrillouinZone).

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, MeshBrillouinZone).

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 MeshBrillouinZone.

tol : float, optional

Relative accuracy for eigenvalues (stopping criterion).

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, MeshBrillouinZone).

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:

(24)\[\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

(25)\[[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.

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

@return \([g]^{-1}\) in the given channel @include tprf/linalg.hpp @note Assign to gf (g2_iw_t) yields move operation while assigning to gf_view (g2_iw_vt) causes extra copy operation

Signature : (triqs_tprf::g2_iw_vt g) -> triqs_tprf::g2_iw_t

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

(26)\[[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.

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

@return \([g]^{-1}\) in the given channel @include tprf/linalg.hpp @note Assign to gf (g2_iw_t) yields move operation while assigning to gf_view (g2_iw_vt) causes extra copy operation

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

(27)\[[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.

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

@return \([g]^{-1}\) in the given channel @include tprf/linalg.hpp @note Assign to gf (g2_iw_t) yields move operation while assigning to gf_view (g2_iw_vt) causes extra copy operation

Signature : (triqs_tprf::g2_iw_vt g) -> triqs_tprf::g2_iw_t

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

(28)\[[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.

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

@return \([g]^{-1}\) in the given channel @include tprf/linalg.hpp @note Assign to gf (g2_iw_t) yields move operation while assigning to gf_view (g2_iw_vt) causes extra copy operation

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

(29)\[[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.

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

@return \([g]^{-1}\) in the given channel @include tprf/linalg.hpp @note Assign to gf (g2_iw_t) yields move operation while assigning to gf_view (g2_iw_vt) causes extra copy operation

Signature : (triqs_tprf::g2_iw_vt g) -> triqs_tprf::g2_iw_t

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

(30)\[[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.

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

@return \([g]^{-1}\) in the given channel @include tprf/linalg.hpp @note Assign to gf (g2_iw_t) yields move operation while assigning to gf_view (g2_iw_vt) causes extra copy operation

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

(31)\[(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.

@tparam CH selects the two-particle channel :param A: two-particle response function \(A \equiv A_{abcd}(\omega, \nu, \nu')\) :param B: two-particle response function \(B \equiv A_{abcd}(\omega, \nu, \nu')\) @return \((A * B)\) in the given channel @include tprf/linalg.hpp @note Assign to gf (g2_iw_t) yields move operation while assigning to gf_view (g2_iw_vt) causes extra copy operation

Signature : (triqs_tprf::g2_iw_vt A, triqs_tprf::g2_iw_vt B) -> triqs_tprf::g2_iw_t

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

(32)\[(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.

@tparam CH selects the two-particle channel :param A: two-particle response function \(A \equiv A_{abcd}(\omega, \nu, \nu')\) :param B: two-particle response function \(B \equiv A_{abcd}(\omega, \nu, \nu')\) @return \((A * B)\) in the given channel @include tprf/linalg.hpp @note Assign to gf (g2_iw_t) yields move operation while assigning to gf_view (g2_iw_vt) causes extra copy operation

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

(33)\[(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.

@tparam CH selects the two-particle channel :param A: two-particle response function \(A \equiv A_{abcd}(\omega, \nu, \nu')\) :param B: two-particle response function \(B \equiv A_{abcd}(\omega, \nu, \nu')\) @return \((A * B)\) in the given channel @include tprf/linalg.hpp @note Assign to gf (g2_iw_t) yields move operation while assigning to gf_view (g2_iw_vt) causes extra copy operation

Signature : (triqs_tprf::g2_iw_vt A, triqs_tprf::g2_iw_vt B) -> triqs_tprf::g2_iw_t

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

(34)\[(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.

@tparam CH selects the two-particle channel :param A: two-particle response function \(A \equiv A_{abcd}(\omega, \nu, \nu')\) :param B: two-particle response function \(B \equiv A_{abcd}(\omega, \nu, \nu')\) @return \((A * B)\) in the given channel @include tprf/linalg.hpp @note Assign to gf (g2_iw_t) yields move operation while assigning to gf_view (g2_iw_vt) causes extra copy operation

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

(35)\[(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.

@tparam CH selects the two-particle channel :param A: two-particle response function \(A \equiv A_{abcd}(\omega, \nu, \nu')\) :param B: two-particle response function \(B \equiv A_{abcd}(\omega, \nu, \nu')\) @return \((A * B)\) in the given channel @include tprf/linalg.hpp @note Assign to gf (g2_iw_t) yields move operation while assigning to gf_view (g2_iw_vt) causes extra copy operation

Signature : (triqs_tprf::g2_iw_vt A, triqs_tprf::g2_iw_vt B) -> triqs_tprf::g2_iw_t

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

(36)\[(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.

@tparam CH selects the two-particle channel :param A: two-particle response function \(A \equiv A_{abcd}(\omega, \nu, \nu')\) :param B: two-particle response function \(B \equiv A_{abcd}(\omega, \nu, \nu')\) @return \((A * B)\) in the given channel @include tprf/linalg.hpp @note Assign to gf (g2_iw_t) yields move operation while assigning to gf_view (g2_iw_vt) causes extra copy operation

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

(37)\[\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.

@tparam CH selects the two-particle channel :param A: two-particle response function \(A \equiv A_{abcd}(\omega, \nu, \nu')\) determinig the shape and size of the unity operator @return the unity operator \(\mathbf{1}\), in the given channel @include tprf/linalg.hpp @note Assign to gf (g2_iw_t) yields move operation while assigning to gf_view (g2_iw_vt) causes extra copy operation

Signature : (triqs_tprf::g2_iw_vt g) -> triqs_tprf::g2_iw_t

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

Constructs the unity-operator in the given channel

(38)\[\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.

@tparam CH selects the two-particle channel :param A: two-particle response function \(A \equiv A_{abcd}(\omega, \nu, \nu')\) determinig the shape and size of the unity operator @return the unity operator \(\mathbf{1}\), in the given channel @include tprf/linalg.hpp @note Assign to gf (g2_iw_t) yields move operation while assigning to gf_view (g2_iw_vt) causes extra copy operation

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

(39)\[\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.

@tparam CH selects the two-particle channel :param A: two-particle response function \(A \equiv A_{abcd}(\omega, \nu, \nu')\) determinig the shape and size of the unity operator @return the unity operator \(\mathbf{1}\), in the given channel @include tprf/linalg.hpp @note Assign to gf (g2_iw_t) yields move operation while assigning to gf_view (g2_iw_vt) causes extra copy operation

Signature : (triqs_tprf::g2_iw_vt g) -> triqs_tprf::g2_iw_t

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

Constructs the unity-operator in the given channel

(40)\[\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.

@tparam CH selects the two-particle channel :param A: two-particle response function \(A \equiv A_{abcd}(\omega, \nu, \nu')\) determinig the shape and size of the unity operator @return the unity operator \(\mathbf{1}\), in the given channel @include tprf/linalg.hpp @note Assign to gf (g2_iw_t) yields move operation while assigning to gf_view (g2_iw_vt) causes extra copy operation

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

(41)\[\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.

@tparam CH selects the two-particle channel :param A: two-particle response function \(A \equiv A_{abcd}(\omega, \nu, \nu')\) determinig the shape and size of the unity operator @return the unity operator \(\mathbf{1}\), in the given channel @include tprf/linalg.hpp @note Assign to gf (g2_iw_t) yields move operation while assigning to gf_view (g2_iw_vt) causes extra copy operation

Signature : (triqs_tprf::g2_iw_vt g) -> triqs_tprf::g2_iw_t

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

Constructs the unity-operator in the given channel

(42)\[\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.

@tparam CH selects the two-particle channel :param A: two-particle response function \(A \equiv A_{abcd}(\omega, \nu, \nu')\) determinig the shape and size of the unity operator @return the unity operator \(\mathbf{1}\), in the given channel @include tprf/linalg.hpp @note Assign to gf (g2_iw_t) yields move operation while assigning to gf_view (g2_iw_vt) causes extra copy operation

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, hopping, orbital_positions=[(0, 0, 0)], orbital_names=[''])[source]

Class describing a tight binding lattice model.

This class is based in the TRIQS tight binding class and has been extended with some extra convienience methods and documentation.

Parameters:

units : list of three-tuples of floats

Basis vectors for the real space lattice.

hopping : dict

Dictionary with three tuple of integeers as keys, describing real space hoppings in multiples of the real space lattice basis vectors, and values being numpy ndarray hopping matrices in the orbital indices.

orbital_positions : list of three three-tuples of floats.

Internal orbital positions in the unit-cell.

orbital_names : list of strings

Names for each orbital.

on_mesh_brillouin_zone(n_k)[source]

Construct a discretization of the tight binding model in reciprocal space with given number of k-points.

Parameters:

n_k : three-tuple of ints

Number of k-points in every dimension.

Returns:

e_k : TRIQS Greens function on a Brillioun zone mesh

Reciprocal space tight binding dispersion.

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

Class building a superlattice on top of a base lattice (TBLattice).

The class inherits from TBLattice and has all the basic TBLattice methods, especially the on_mesh_brillouin_zone-method.

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=10, mu_min=-10.0)[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 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
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=10, mu_min=-10.0)[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.

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.

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 pytriqs.archive 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
convert_keys_from_string_to_python(dict_key)[source]

pytriqs.archive.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…

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.]