Welcome to GfTool’s documentation!
- Release:
0.11.0+59.g354f390
- Date:
2024-04-27
- DOI:
This reference manual details functions, modules, and objects included in GfTool, describing what they are and what they do.
GfTool is a collection of commonly used Green’s functions and utilities. The main purpose of this module is to have a tested and thus reliable basis to do numerics. It happened to me too often, that I just made a mistake copying the Green’s function and was then wondering what was wrong with my algorithm.

Selection of lattice Green’s functions or rather the corresponding DOSs available in GfTool.
The main use case of GfTool was DMFT and its real space generalization, in particular using CT-QMC algorithms.
Getting started
GfTool’s main content
gftool
Collection of non-interacting Green’s functions and spectral functions, see also the
lattice
submodule.Utility functions like Matsubara frequencies and Fermi functions.
Reliable calculation of particle numbers via Matsubara sums.
cpa
/beb
Single site approximation to disorder.
Diagonal disorder only (CPA) and diagonal and off-diagonal (BEB).
Average local Green’s function and component Green’s functions.
fourier
Fourier transforms from Matsubara frequencies to imaginary time and back, including the handling of high-frequencies moments (especially important for transforms from Matsubara to imaginary time).
Laplace transform from real times to complex frequencies.
matrix
Helper for Green’s functions in matrix form.
pade
Analytic continuation via the Padé algorithm.
Calculates a rational polynomial as interpolation of the data points.
Note: the module is out-dated, so it’s not well documented at a bit awkward to use. Consider using
polepade
instead.
polepade
Analytic continuation via a pole-based variant of the Padé algorithm.
Based on explicit calculation of the pole in a least-squares sense.
Allows including uncertainties as weights.
siam
Basic Green’s functions for the non-interacting Single Impurity Anderson model.
Installation
The package is available on PyPI:
$ pip install gftool
For conda users, GfTool is also available on conda-forge
$ conda install -c conda-forge gftool
Alternatively, it can be installed via GitHub. You can install it using
$ pip install https://github.com/DerWeh/gftools/archive/VERSION.zip
where VERSION can be a release (e.g. 0.5.1) or a branch (e.g. develop). (As always, it is not advised to install it into your system Python, consider using pipenv, venv, conda, pyenv, or similar tools.) Of course, you can also clone or fork the project.
If you clone the project, you can locally build the documentation:
$ pip install -r requirements-doc.txt
$ python setup.py build_sphinx
Note on documentation
We try to follow numpy
broadcasting rules. Many functions acting on an axis
act like generalized ufuncs. In this case, a function can be called for
stacked arguments instead of looping over the specific arguments.
We indicate this by argument shapes containing an ellipse e.g. (…) or (…, N).
It must be possible for all ellipses to be broadcasted against each other.
A good example is the fourier
module.
We calculate the Fourier transforms iw2tau
for Green’s
functions with different on-site energies without looping:
>>> e_onsite = np.array([-0.5, 0, 0.5])
>>> beta = 10
>>> iws = gt.matsubara_frequencies(range(1024), beta=beta)
>>> gf_iw = gt.bethe_gf_z(iws - e_onsite[..., np.newaxis], half_bandwidth=1.0)
>>> gf_iw.shape
(3, 1024)
>>> from gftool import fourier
>>> gf_tau = fourier.iw2tau(gf_iw, beta=beta, moments=np.ones([1]))
>>> gf_tau.shape
(3, 2049)
The moments are automatically broadcasted. We can also explicitly give the second moments:
>>> moments = np.stack([np.ones([3]), e_onsite], axis=-1)
>>> gf_tau = fourier.iw2tau(gf_iw, beta=beta, moments=moments)
>>> gf_tau.shape
(3, 2049)
Tutorial
This tutorial explains some of the basic functionality. Throughout the tutorial we assume you have imported GfTool as
>>> import gftool as gt
and the packages numpy
and matplotlib
are imported as usual
>>> import numpy as np
>>> import matplotlib.pyplot as plt
Lattice Green’s functions
The package contains non-interacting Green’s functions for some tight-binding
lattices. They can be found in gftool.lattice
.
E.g. the dos
of the Bethe lattice, bethe
:
>>> ww = np.linspace(-1.1, 1.1, num=1000)
>>> dos_ww = gt.lattice.bethe.dos(ww, half_bandwidth=1.)
>>> __ = plt.plot(ww, dos_ww)
>>> __ = plt.xlabel(r"$\epsilon$")
>>> plt.show()

Typically, a shorthand for these functions exist in the top-level module e.g.
gftool.bethe_dos
>>> gt.bethe_dos is gt.lattice.bethe.dos
True
The corresponding Green’s functions are also available.
The Green’s functions can be evaluated for any complex frequency,
excluding the real axis, where it might become singular.
E.g. the gf_z
of the face-centered cubic (fcc) lattice,
fcc
, on a contour parallel to the real axis:
>>> ww = np.linspace(-0.9, 1.7, num=1000) + 1e-6j
>>> gf_ww = gt.lattice.fcc.gf_z(ww, half_bandwidth=1.)
>>> __ = plt.axhline(0, color="dimgray", linewidth=0.8)
>>> __ = plt.axvline(0, color="dimgray", linewidth=0.8)
>>> __ = plt.plot(ww.real, gf_ww.real)
>>> __ = plt.plot(ww.real, gf_ww.imag)
>>> __ = plt.xlabel(r"$\omega$")
>>> __ = plt.ylim(-7.0, 2.5)
>>> plt.show()

or on the Matsubara axis:
>>> beta = 50
>>> iws = gt.matsubara_frequencies(range(128), beta=beta)
>>> gf_iw = gt.lattice.fcc.gf_z(iws, half_bandwidth=1.)
>>> __ = plt.axhline(0, color="dimgray", linewidth=0.8)
>>> __ = plt.plot(gf_iw.real, "x--")
>>> __ = plt.plot(gf_iw.imag, "+--")
>>> __ = plt.xlabel("$n$")
>>> plt.show()

Density
We can also calculate the density (occupation number) from the imaginary axis for local Green’s function. We have the relation
To calculate the density for a given temperature from 1024 (fermionic)
Matsubara frequencies we use density_iw
:
>>> temperature = 0.02
>>> iws = gt.matsubara_frequencies(range(1024), beta=1./temperature)
>>> gf_iw = gt.bethe_gf_z(iws, half_bandwidth=1.)
>>> occ = gt.density_iw(iws, gf_iw, beta=1./temperature)
>>> occ
0.5
We can also search the chemical potential \(μ\) for a given occupation
using chemical_potential
.
To get, e.g., the Bethe lattice at quarter filling, we write:
>>> occ_quarter = 0.25
>>> def bethe_occ_diff(mu):
... """Calculate the difference to the desired occupation, note the sign."""
... gf_iw = gt.bethe_gf_z(iws + mu, half_bandwidth=1.)
... return gt.density_iw(iws, gf_iw, beta=1./temperature) - occ_quarter
...
>>> mu = gt.chemical_potential(bethe_occ_diff)
>>> mu
-0.406018...
Validate the result:
>>> gf_quarter_iw = gt.bethe_gf_z(iws + mu, half_bandwidth=1.)
>>> gt.density_iw(iws, gf_quarter_iw, beta=1./temperature)
0.249999...
Fourier transform
GfTool offers also accurate Fourier transformations between Matsubara frequencies
and imaginary time for local Green’s functions, see gftool.fourier
.
As a major part of the package, these functions are gu-functions.
This is indicated in the docstrings via the shapes (…, N). The ellipsis
stands for arbitrary leading dimensions. Let’s consider a simple example with magnetic splitting.
>>> beta = 20 # inverse temperature
>>> h = 0.3 # magnetic splitting
>>> eps = np.array([-0.5*h, +0.5*h]) # on-site energy
>>> iws = gt.matsubara_frequencies(range(1024), beta=beta)
We can calculate the Fourier transform using broadcasting, no need for any loops.
>>> gf_iw = gt.bethe_gf_z(iws + eps[:, np.newaxis], half_bandwidth=1)
>>> gf_iw.shape
(2, 1024)
>>> gf_tau = gt.fourier.iw2tau(gf_iw, beta=beta)
The Fourier transform generates the imaginary time Green’s function on the interval \(τ ∈ [0^+, β^-]\)
>>> tau = np.linspace(0, beta, num=gf_tau.shape[-1])
>>> __ = plt.plot(tau, gf_tau[0], label=r'$\sigma=\uparrow$')
>>> __ = plt.plot(tau, gf_tau[1], label=r'$\sigma=\downarrow$')
>>> __ = plt.xlabel(r'$\tau$')
>>> __ = plt.legend()
>>> plt.show()

We see the asymmetry due to the magnetic field. Let’s check the back transformation.
>>> gf_ft = gt.fourier.tau2iw(gf_tau, beta=beta)
>>> np.allclose(gf_ft, gf_iw, atol=2e-6)
True
Up to a certain threshold the transforms agree, they are not exact inverse transformations here. Accuracy can be improved e.g. by providing (or fitting) high-frequency moments.
Single site approximation of disorder
We also offer the single site approximation for disordered Hamiltonians,
namely cpa
and its extension to off-diagonal disorder beb
.
These methods treat substitutional disorder.
A multi-component system is considered, where each lattice site is randomly
occupied by one of the components.
The concentration of the components is known.
Coherent potential approximation (CPA)
We first consider the cpa
, where only the on-site energies depend on the component.
As example, we consider a system of three components.
We choose the on-site energies and concentrations (which should add to 1),
as lattice we consider a Bethe lattice with half-bandwidth 1:
>>> from functools import partial
>>> e_onsite = np.array([-0.3, -0.1, 0.4])
>>> concentration = np.array([0.3, 0.5, 0.2])
>>> g0 = partial(gt.bethe_gf_z, half_bandwidth=1)
The average local Green’s function and the component Green’s functions
(conditional average for local site fixed to a specific component) are calculated
in CPA using an effective medium.
The self-consistent effective medium is obtained via a root search
using solve_root
:
>>> ww = np.linspace(-1.5, 1.5, num=501) + 1e-6j
>>> self_cpa_ww = gt.cpa.solve_root(ww, e_onsite, concentration, hilbert_trafo=g0)
The average Green’s function is
>>> gf_coher_ww = g0(ww - self_cpa_ww)
>>> __ = plt.plot(ww.real, -1/np.pi*gf_coher_ww.imag)
>>> __ = plt.xlabel(r"$\omega$")
>>> plt.show()

For frequencies close to the real axis, issues might arise, that the conjugate
solution (advanced instead of retarded) is obtained.
The default restricted=True uses some heuristic to avoid this.
In this example we see at the left band-edge,
that for small imaginary part this can still fail.
In this case, it is enough to just increase the accuracy of the root search.
Additional keyword arguments are passed to scipy.optimize.root
:
>>> self_cpa_ww = gt.cpa.solve_root(ww, e_onsite, concentration, hilbert_trafo=g0,
... options=dict(fatol=1e-10))
>>> gf_coher_ww = g0(ww - self_cpa_ww)
>>> __ = plt.plot(ww.real, -1/np.pi*gf_coher_ww.imag)
>>> __ = plt.xlabel(r"$\omega$")
>>> plt.show()

Now, everything looks fine.
The component Green’s functions are calculated by gftool.cpa.gf_cmpt_z
.
The law of total expectation relates the component Green’s functions to the
average Green’s function: np.sum(concentration*gf_cmpt_ww, axis=-1) == gf_coher_ww:
>>> gf_cmpt_ww = gt.cpa.gf_cmpt_z(ww, self_cpa_ww, e_onsite, hilbert_trafo=g0)
>>> np.allclose(np.sum(concentration*gf_cmpt_ww, axis=-1), gf_coher_ww)
True
>>> for cmpt in range(3):
... __ = plt.plot(ww.real, -1/np.pi*gf_cmpt_ww[..., cmpt].imag, label=f"cmpt {cmpt}")
>>> __ = plt.plot(ww.real, -1/np.pi*gf_coher_ww.imag, linestyle=':', label="avg")
>>> __ = plt.legend()
>>> __ = plt.xlabel(r"$\omega$")
>>> plt.show()

Of course, it can be calculated for any lattice Hilbert transform.
Furthermore, the function is vectorized. Let’s consider a fcc
lattice, where one component has different on-site energies for up and down spin.
The on-site energies can simply be stacked as 2-dimensional array.
We can also take the previous self-energy as a starting guess self_cpa_z0:
>>> e_onsite = np.array([[-0.3, +0.15, +0.4],
... [-0.3, -0.35, +0.4]])
>>> concentration = np.array([0.3, 0.5, 0.2])
>>> g0 = partial(gt.fcc_gf_z, half_bandwidth=1)
>>> self_cpa_ww = gt.cpa.solve_root(ww[:, np.newaxis], e_onsite, concentration,
... hilbert_trafo=g0, options=dict(fatol=1e-8),
... self_cpa_z0=self_cpa_ww[:, np.newaxis])
>>> gf_cmpt_ww = gt.cpa.gf_cmpt_z(ww[:, np.newaxis], self_cpa_ww, e_onsite, hilbert_trafo=g0)
>>> __, axes = plt.subplots(nrows=2, sharex=True)
>>> for spin, ax in enumerate(axes):
... for cmpt in range(3):
... __ = ax.plot(ww.real, -1/np.pi*gf_cmpt_ww[:, spin, cmpt].imag, label=f"cmpt {cmpt}")
>>> __ = plt.legend()
>>> __ = plt.xlabel(r"$\omega$")
>>> plt.show()

Blackman, Esterling, Berk (BEB)
The beb
formalism is an extension of cpa
to off-diagonal disorder.
It allows us to provide different hopping amplitudes.
We have the additional parameter hopping which gives the relative hopping amplitudes.
The cpa
corresponds to hopping=np.ones([N, N]), where N is the number
of components.
The beb
module works very similar to cpa
:
We use solve_root
to get the effective medium,
in BEB, however, the effective medium is a matrix.
Next the component Green’s function are calculated using gf_loc_z
.
These are, however, already multiplied by the concentration.
So the average Green’s function is gf_loc_z.sum(axis=-1).
Let’s compare cpa
and beb
:
>>> from functools import partial
>>> e_onsite = np.array([-0.3, -0.1, 0.4])
>>> concentration = np.array([0.3, 0.5, 0.2])
>>> hopping = np.ones([3, 3])
>>> g0 = partial(gt.bethe_gf_z, half_bandwidth=1)
>>> ww = np.linspace(-1.5, 1.5, num=501) + 1e-5j
>>> self_cpa_ww = gt.cpa.solve_root(ww, e_onsite, concentration, hilbert_trafo=g0)
>>> gf_coher_ww = g0(ww - self_cpa_ww)
>>> self_beb_ww = gt.beb.solve_root(ww, e_onsite, concentration=concentration,
... hopping=hopping, hilbert_trafo=g0)
>>> gf_loc_ww = gt.beb.gf_loc_z(ww, self_beb_ww, hopping=hopping, hilbert_trafo=g0)
>>> __ = plt.plot(ww.real, -1/np.pi*gf_coher_ww.imag, label="CPA avg")
>>> __ = plt.plot(ww.real, -1/np.pi*gf_loc_ww.sum(axis=-1).imag,
... linestyle='--', label="BEB avg")
>>> __ = plt.xlabel(r"$\omega$")
>>> plt.show()

Of course, also the components match:
>>> gf_cmpt_ww = gt.cpa.gf_cmpt_z(ww, self_cpa_ww, e_onsite, hilbert_trafo=g0)
>>> c_gf_cmpt_ww = gf_cmpt_ww * concentration # to compare with BEB
>>> for cmpt in range(3):
... __ = plt.plot(ww.real, -1/np.pi*c_gf_cmpt_ww[..., cmpt].imag, label=f"CPA {cmpt}")
... __ = plt.plot(ww.real, -1/np.pi*gf_loc_ww[..., cmpt].imag, '--', label=f"BEB {cmpt}")
>>> __ = plt.legend()
>>> __ = plt.xlabel(r"$\omega$")
>>> plt.show()

The relevant case is when hopping differs from the CPA case. Then the components can have different band-widths and also the hopping between different components can be different. Let’s say we have two components ‘A’ and ‘B’. The values hopping=np.array([[1.0, 0.5], [0.5, 1.2]]) mean that the hopping amplitude between ‘B’ sites is 1.2 times the hopping amplitude between ‘A’ sites; the hopping amplitude from ‘A’ to ‘B’ is 0.5 times the hopping amplitude between ‘A’ sites.
>>> from functools import partial
>>> e_onsite = np.array([0.2, -0.2])
>>> concentration = np.array([0.3, 0.7])
>>> hopping = np.array([[1.0, 0.5],
... [0.5, 1.2]])
>>> g0 = partial(gt.bethe_gf_z, half_bandwidth=1)
>>> ww = np.linspace(-1.5, 1.5, num=501) + 1e-5j
>>> self_beb_ww = gt.beb.solve_root(ww, e_onsite, concentration=concentration,
... hopping=hopping, hilbert_trafo=g0)
>>> gf_loc_ww = gt.beb.gf_loc_z(ww, self_beb_ww, hopping=hopping, hilbert_trafo=g0)
>>> __ = plt.plot(ww.real, -1/np.pi*gf_loc_ww[..., 0].imag, label="A")
>>> __ = plt.plot(ww.real, -1/np.pi*gf_loc_ww[..., 1].imag, label="B")
>>> __ = plt.plot(ww.real, -1/np.pi*gf_loc_ww.sum(axis=-1).imag,
... linestyle='--', label="BEB avg")
>>> __ = plt.legend()
>>> __ = plt.xlabel(r"$\omega$")
>>> plt.show()

Additional diagnostic output is logged, you can get information on the convergence by setting:
>>> import logging
>>> logging.basicConfig()
>>> logging.getLogger('gftool.beb').setLevel(logging.DEBUG)
Matrix Green’s functions via diagonalization
The module gftool.matrix
contains some helper functions for matrix diagonalization.
A main use case is to calculate the one-particle Green’s function from the resolvent.
Instead of inverting the matrix for every frequency point,
we can diagonalize the Hamiltonian once:
Let’s consider the simple example of a 2D square lattice with nearest-neighbor hopping. The Hamiltonian can be easily constructed:
>>> N = 21 # system size in one dimension
>>> t = tx = ty = 0.5 # hopping amplitude
>>> hamilton = np.zeros([N]*4)
>>> diag = np.arange(N)
>>> hamilton[diag[1:], :, diag[:-1], :] = hamilton[diag[:-1], :, diag[1:], :] = -tx
>>> hamilton[:, diag[1:], :, diag[:-1]] = hamilton[:, diag[:-1], :, diag[1:]] = -ty
>>> ham_mat = hamilton.reshape(N**2, N**2) # turn in into a matrix
Let’s diagonalize it using the helper in gftool.matrix
and calculated the Green’s function
>>> dec = gt.matrix.decompose_her(ham_mat)
>>> ww = np.linspace(-2.5, 2.5, num=201) + 1e-1j # frequency match
>>> gf_ww = dec.reconstruct(1.0/(ww[:, np.newaxis] - dec.eig))
>>> gf_ww = gf_ww.reshape(ww.size, *[N]*4) # reshape for easy access
Where we used reconstruct
to calculate the matrix product for the given diagonal matrix.
Let’s check the local spectral function of the central lattice site:
>>> __ = plt.plot(ww.real, -1.0/np.pi*gf_ww.imag[:, N//2, N//2, N//2, N//2])
>>> __ = plt.plot(ww.real, -1.0/np.pi*gt.square_gf_z(ww, half_bandwidth=4*t).imag,
... color='black', linestyle='--')
>>> __ = plt.xlabel(r"$\omega$")
>>> plt.show()

Oftentimes we are only interested in the local Green’s functions and can avoid a large part of the computation, only calculating the diagonal elements. This can be done using the kind argument:
>>> gf_diag = dec.reconstruct(1.0/(ww[:, np.newaxis] - dec.eig), kind='diag')
>>> gf_diag = gf_diag.reshape(ww.size, N, N)
gftool
Collection of commonly used Green’s functions and utilities.
Main purpose is to have a tested base.
Submodules
Different function bases. |
|
Blackman, Esterling, and Berk (BEB) approach to off-diagonal disorder. |
|
Coherent cluster approximation (CPA) to substitutional disorder. |
|
Fourier transformations of Green's functions. |
|
Hermite-Padé approximants from Taylor expansion. |
|
Collection of different lattices and their Green's functions. |
|
Collection of linear algebra algorithms not contained in |
|
Linear prediction to extrapolated retarded Green's function. |
|
Functions to work with Green's functions in matrix from. |
|
Padé analytic continuation for Green's functions and self-energies. |
|
Padé based on robust pole finding. |
|
Basic functions for the (non-interacting) single impurity Anderson model (SIAM). |
Glossary
Green’s functions and lattices
1D
|
DOS of non-interacting 1D lattice. |
|
Calculate the m th moment of the 1D DOS. |
|
Local Green's function of the 1D lattice. |
|
Hilbert transform of non-interacting DOS of the 1D lattice. |
2D
|
DOS of non-interacting 2D square lattice. |
|
Calculate the m th moment of the square DOS. |
|
Local Green's function of the 2D square lattice. |
|
Hilbert transform of non-interacting DOS of the square lattice. |
|
DOS of non-interacting 2D triangular lattice. |
|
Calculate the m th moment of the triangular DOS. |
|
Local Green's function of the 2D triangular lattice. |
|
Hilbert transform of non-interacting DOS of the triangular lattice. |
|
DOS of non-interacting 2D honeycomb lattice. |
|
Calculate the m th moment of the honeycomb DOS. |
|
Local Green's function of the 2D honeycomb lattice. |
|
Hilbert transform of non-interacting DOS of the honeycomb lattice. |
3D
|
Local Green's function of 3D simple cubic lattice. |
|
Calculate the m th moment of the simple cubic DOS. |
|
Local Green's function of 3D simple cubic lattice. |
|
Hilbert transform of non-interacting DOS of the simple cubic lattice. |
|
DOS of non-interacting 3D body-centered cubic lattice. |
|
Calculate the m th moment of the body-centered cubic DOS. |
|
Local Green's function of 3D body-centered cubic (bcc) lattice. |
|
Hilbert transform of non-interacting DOS of the body-centered cubic lattice. |
|
DOS of non-interacting 3D face-centered cubic lattice. |
|
Calculate the m th moment of the face-centered cubic DOS. |
|
Local Green's function of the 3D face-centered cubic (fcc) lattice. |
|
Hilbert transform of non-interacting DOS of the face-centered cubic lattice. |
Miscellaneous
|
DOS of non-interacting Bethe lattice for infinite coordination number. |
|
Calculate the m th moment of the Bethe DOS. |
|
Local Green's function of Bethe lattice for infinite coordination number. |
|
First derivative of local Green's function of Bethe lattice for infinite coordination number. |
|
Second derivative of local Green's function of Bethe lattice for infinite coordination number. |
|
Hilbert transform of non-interacting DOS of the Bethe lattice. |
|
Green's function given by a finite number of poles. |
|
First derivative of Green's function given by a finite number of poles. |
|
High-frequency moments of the pole Green's function. |
|
Retarded time Green's function given by a finite number of poles. |
|
Imaginary time Green's function given by a finite number of poles. |
|
Bosonic imaginary time Green's function given by a finite number of poles. |
|
Self-energy in Hubbard-I approximation (atomic solution). |
|
Green's function for the two site Hubbard model on a dimer. |
|
Surface Green's function for stacked layers. |
Statistics and particle numbers
|
Return the Fermi function 1/(exp(βϵ)+1). |
|
Return the 1st derivative of the Fermi function. |
|
Inverse of the Fermi function. |
|
Return the Bose function 1/(exp(βϵ)-1). |
|
Return fermionic Matsubara frequencies \(iω_n\) for the points n_points. |
|
Return bosonic Matsubara frequencies \(iν_n\) for the points n_points. |
|
Return num fermionic Padé frequencies \(iz_p\). |
|
Calculate the number density of the Green's function gf_iw at finite temperature beta. |
|
Search chemical potential for a given occupation. |
|
Calculate the number density of the Green's function gf_iw at finite temperature beta. |
|
Return an estimate for the upper bound of the error in the density. |
|
Return an estimate for the upper bound of the error in the density. |
|
Return data for visual inspection of the density error. |
Utilities
Functions
Get version information or return default if unable to do so. |
Classes
|
What’s New
unreleased
Bug fixes
fix vectorization of
linearprediction
functions
0.11.0 (2022-04-29)
New Features
Add pole-base Padé analytic continuation
polepade
(41d57537).Allows determining number of poles avoiding overfitting.
Least-squares based approach allowing to include uncertainties of input data.
Add
linearprediction
module (b1c8f636).Linear prediction can be used to extend retarded-time Green’s functions.
Add Padé-Fourier approach to Laplace transform (fe1ac173).
Padé-Fourier allows to significantly reduce the truncation error. This allows for contours closer to real-axis for a given maximal time.
Linear Padé approximant
tt2z_pade
based on simple polesQuadratic Hermite-Padé approximant
tt2z_herm2
including quadratic branch cuts but introducing ambiguity which branch to choose.Module
hermpade
implements the necessary approximants.
Internal improvements
Use
numpy.testing.assert_allclose
for tests, providing more verbose output (dbb8fd7c).
Documentation
Start to adhere more closely to
numpydoc
(40d57d45).
0.10.1 (2021-12-01)
New Features
Internal improvements
Switch from Travis to GitHub actions #20 (23ba0a34)
This adds test for Mac and Windows
Documentation
Bug fixes
Accurately calculate
gftool.lattice.sc.dos
aroundeps=0
(5693184e). Previously, DOS was incorrect for tiny values, e.g.eps=1e-16
.
0.10.0 (2021-09-19)
Breaking Changes
Drop support for Python 3.6, minimal version is now 3.7
Content of
gftool.matrix
was renamed more appropriately:xi of
Decomposition
is noweig
, as it contains the eigenvaluesNew functions
decompose_mat
for general matrices,decompose_sym
for complex symmetric matrices, anddecompose_her
for Hermitian matrices.
Depreciations
Deprecate the
matrix
functionsdecompose_gf
,decompose_hamiltonian
,from_gf
, andfrom_hamiltonian
.
Documentation
New index page independent of README, separated Getting started page.
Improve Tutorial and
gftool.matrix
Generate PDF documentation on ReadTheDocs (3122e1ba)
Internal improvements
Use eigendecomposition instead of SVD in
gftool.beb
(0475c110)Drop slow
asfortranarray
ingftool.matrix
(4865cc05)Use pre-commit (6f4028d3)
0.9.1 (2021-06-01)
Bug fixes
CPA:
return scalar mu in
gftool.cpa.solve_fxdocc_root
(10fae4d)find mu more reliably
Other New Features
SIAM: add greater and lesser Green’s functions
gf0_loc_gr_t
andgf0_loc_le_t
(ea541f3)
0.9.0 (2021-05-09)
New Features
0.8.1 (2021-04-25)
New Features
The 3D cubic lattices were added:
body-centered cubic
gftool.lattice.bcc
(406acef8)face-centered cubic
gftool.lattice.fcc
(ddd559cb)
0.8.0 (2021-04-17)
New Features
The gftool.lattice
module was extended, especially regarding two-dimensional lattice.
There were also some enhancements, given DOS moments are now up to order 20,
and they should be accurate to machine precision.
The following lattices where added with full interface:
Simple cubic:
gftool.lattice.sc
(4e3021) by Andreas ÖstlinHoneycomb:
gftool.lattice.honeycomb
(7aa3133)Triangular:
gftool.lattice.triangular
(c56f33e)
Local Green’s function and DOS is now also available for the following lattices:
Lieb:
gftool.lattice.lieb
(c76e948)Kagome:
gftool.lattice.kagome
(28a41c0)Bethe lattice with general coordination:
gftool.lattice.bethez
(2648cf4)Rectangular:
gftool.lattice.rectangular
Other New Features
add retarded time Green’s function give by its poles
gftool.pole_gf_ret_t
added
gftool.siam
module with some basics for the non-interacting siam
Depreciations
gftool.density
is deprecated and will likely be discontinued. Consider the more flexiblegftool.density_iw
instead.
Documentation
Button to toggle the prompt (>>>) was added (46b6f39)
Internal improvements
0.7.0 (2020-10-18)
Breaking Changes
The
gftool.pade
module had a minor rework. The behavior of filters changed. Future breaking changes are to be expected, the module is not well-structured.
New Features
add
gftool.lattice.onedim
for Green’s function of one-dimensional latticeadd fitting of high-frequency moment to
gftool.fourier.iw2tau
(e2c92e2)
Other New Features
add
gftool.density_iw
function as common interface to calculate occupation number from Matsubara or Padé frequenciesallow calculation of
gftool.lattice.bethe
for Bethe lattice at complex points (note, that this is probably not a physically meaningful quantity) (ccbac7b)add stress tensor transformation
gftool.lattice.square.stress_trafo
for 2D (528fb21)
Bug fixes
Fix constant in
gftool.fourier.tau2iw_ft_lin
(e2163e3). This error most likely didn’t significantly affect any results for a reasonable number of tau-points.gftool.density
should work now with gu-style matrices (4deffdf)
Documentation
Functions exposed at the top level (
gftool
) should now properly appear in the documentation.