Lattice Green’s functions
- triqs_tprf.lattice.lattice_dyson_g0_wk()
Signature : (float mu, triqs_tprf::e_k_cvt e_k, gf_mesh<triqs::mesh::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
\[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})\)
- 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})\)
- 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
\[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})\)
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.
\[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_{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.
\[\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_contiguous_view<std::complex<double>,4> U) -> triqs_tprf::chi_wk_t 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}, i\omega_n)\)
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}, i\omega_n)\)
- 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})\).
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
\[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)
\[W_{abcd}(i\omega_n, \mathbf{k}) = W^{(full)}_{abcd}(i\omega_n, \mathbf{k}) - V_{abcd}(\mathbf{k})\]- Parameters:
PI_wk :
polarization bubble \(\Pi_{abcd}(i\omega_n, \mathbf{k})\)
V_k :
static interaction \(V_{abcd}(\mathbf{k})\)
- Returns:
out :
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
\[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)
\[W_{abcd}(i\omega_n, \mathbf{k}) = W^{(full)}_{abcd}(i\omega_n, \mathbf{k}) - V_{abcd}(\mathbf{k})\]- Parameters:
chi_wk :
polarization bubble \(\Pi_{abcd}(i\omega_n, \mathbf{k})\)
V_k :
static interaction \(V_{abcd}(\mathbf{k})\)
- Returns:
out :
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
\[\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.
\[G_{ab}(\tau, \mathbf{r}) = \mathcal{F}^{-1} \left\{ G_{ab}(i\omega_n, \mathbf{k}) \right\}\]\[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
\[\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
\[\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(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
- 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.
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
- 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
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 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
the real diagonal entries
the real part of the upper triangular entries
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.
- 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.
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
- 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 ofParameterCollection
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