Public names index
SparseIR.AugmentedBasis — Type
AugmentedBasis <: AbstractBasisAugmented basis on the imaginary-time/frequency axis.
Groups a set of additional functions, augmentations, with a given basis. The augmented functions then form the first basis functions, while the rest is provided by the regular basis, i.e.:
u[l](x) == l < naug ? augmentations[l](x) : basis.u[l-naug](x),where naug = length(augmentations) is the number of added basis functions through augmentation. Similar expressions hold for Matsubara frequencies.
Augmentation is useful in constructing bases for vertex-like quantities such as self-energies [wallerberger2021] and when constructing a two-point kernel that serves as a base for multi-point functions [shinaoka2018].
Bases augmented with TauConst and TauLinear tend to be poorly conditioned. Care must be taken while fitting and compactness should be enforced if possible to regularize the problem.
While vertex bases, i.e. bases augmented with MatsubaraConst, stay reasonably well-conditioned, it is still good practice to treat the Hartree–Fock term separately rather than including it in the basis, if possible.
See also: MatsubaraConst for vertex basis [wallerberger2021], TauConst, TauLinear for multi-point [shinaoka2018]
SparseIR.Bosonic — Type
Bosonic statistics.
SparseIR.DiscreteLehmannRepresentation — Type
DiscreteLehmannRepresentation(basis::AbstractBasis, poles=nothing)Construct a DLR basis from an IR basis.
If poles is not provided, uses the default omega sampling points from the IR basis.
SparseIR.DiscreteLehmannRepresentation — Type
DiscreteLehmannRepresentation{S,B} <: AbstractBasis{S}Discrete Lehmann representation (DLR) with poles selected according to extrema of IR.
This type wraps the C API DLR functionality. The DLR basis is a variant of the IR basis that uses a "sketching" approach - representing functions as a linear combination of poles on the real-frequency axis:
G(iv) == sum(a[i] / (iv - w[i]) for i in 1:npoles)Fields
ptr::Ptr{spir_basis}: Pointer to the C DLR objectbasis::B: The underlying IR basispoles::Vector{Float64}: Pole locations on the real-frequency axis
SparseIR.Fermionic — Type
Fermionic statistics.
SparseIR.FiniteTempBasis — Type
FiniteTempBasis <: AbstractBasisIntermediate representation (IR) basis for given temperature.
For a continuation kernel K from real frequencies, ω ∈ [-ωmax, ωmax], to imaginary time, τ ∈ [0, β], this type stores the truncated singular value expansion or IR basis:
K(τ, ω) ≈ sum(u[l](τ) * s[l] * v[l](ω) for l in 1:L)This basis is inferred from a reduced form by appropriate scaling of the variables.
Fields
u::PiecewiseLegendrePolyVector: Set of IR basis functions on the imaginary time (tau) axis. These functions are stored as piecewise Legendre polynomials.To obtain the value of all basis functions at a point or a array of points
x, you can call the functionu(x). To obtain a single basis function, a slice or a subsetl, you can useu[l].uhat::PiecewiseLegendreFTVector: Set of IR basis functions on the Matsubara frequency (wn) axis. These objects are stored as a set of Bessel functions.To obtain the value of all basis functions at a Matsubara frequency or a array of points
wn, you can call the functionuhat(wn). Note that we expect reduced frequencies, which are simply even/odd numbers for bosonic/fermionic objects. To obtain a single basis function, a slice or a subsetl, you can useuhat[l].s: Vector of singular values of the continuation kernelv::PiecewiseLegendrePolyVector: Set of IR basis functions on the real frequency (w) axis. These functions are stored as piecewise Legendre polynomials.To obtain the value of all basis functions at a point or a array of points
w, you can call the functionv(w). To obtain a single basis function, a slice or a subsetl, you can usev[l].
SparseIR.FiniteTempBasis — Method
FiniteTempBasis{S}(β, ωmax, ε; kernel=LogisticKernel(β * ωmax), sve_result=SVEResult(kernel, ε), max_size=-1)Construct a finite temperature basis suitable for the given S (Fermionic or Bosonic) and cutoffs β and ωmax.
Arguments
β: Inverse temperature (must be positive)ωmax: Frequency cutoff (must be non-negative)ε: This parameter controls the number of basis functions. Only singular values ≥ ε * s[1] are kept. Typical values are 1e-6 to 1e-12 depending on the desired accuracy for your calculations. If ε is smaller than the square root of double precision machine epsilon (≈ 1.49e-8), the library will automatically use higher precision for the singular value decomposition, resulting in longer computation time for basis generation.
The number of basis functions grows logarithmically as log(1/ε) log (β * ωmax).
SparseIR.FiniteTempBasis — Method
FiniteTempBasis(stat::Statistics, β, ωmax, ε; kernel=LogisticKernel(β * ωmax), sve_result=SVEResult(kernel, ε), max_size=-1)Convenience constructor that matches SparseIR.jl signature.
Construct a finite temperature basis for the given statistics type and cutoffs.
Arguments
stat: Statistics type (Fermionic()orBosonic())β: Inverse temperature (must be positive)ωmax: Frequency cutoff (must be non-negative)ε: Accuracy target for the basis. This parameter controls the number of basis functions. Only singular values ≥ ε * s[1] are kept. Typical values are 1e-6 to 1e-12 depending on the desired accuracy for your calculations. If ε is smaller than the square root of double precision machine epsilon (≈ 1.49e-8), the library will automatically use higher precision for the singular value decomposition, resulting in longer computation time for basis generation.
The number of basis functions grows logarithmically as log(1/ε) log (β * ωmax).
SparseIR.FiniteTempBasisSet — Type
FiniteTempBasisSetType for holding IR bases and sparse-sampling objects.
An object of this type holds IR bases for fermions and bosons and associated sparse-sampling objects.
Fields
- basis_f::FiniteTempBasis: Fermion basis
- basis_b::FiniteTempBasis: Boson basis
- tau::Vector{Float64}: Sampling points in the imaginary-time domain
- wn_f::Vector{Int}: Sampling fermionic frequencies
- wn_b::Vector{Int}: Sampling bosonic frequencies
- smpltauf::TauSampling: Sparse sampling for tau & fermion
- smpltaub::TauSampling: Sparse sampling for tau & boson
- smplwnf::MatsubaraSampling: Sparse sampling for Matsubara frequency & fermion
- smplwnb::MatsubaraSampling: Sparse sampling for Matsubara frequency & boson
- sve_result::Tuple{PiecewiseLegendrePoly,Vector{Float64},PiecewiseLegendrePoly}: Results of SVE
Getters
- beta::Float64: Inverse temperature
- ωmax::Float64: Cut-off frequency
SparseIR.LogisticKernel — Type
LogisticKernel <: AbstractKernelFermionic/bosonic analytical continuation kernel.
In dimensionless variables $x = 2 τ/β - 1$, $y = β ω/Λ$, the integral kernel is a function on $[-1, 1] × [-1, 1]$:
\[ K(x, y) = \frac{e^{-Λ y (x + 1) / 2}}{1 + e^{-Λ y}}\]
LogisticKernel is a fermionic analytic continuation kernel. Nevertheless, one can model the $τ$ dependence of a bosonic correlation function as follows:
\[ ∫ \frac{e^{-Λ y (x + 1) / 2}}{1 - e^{-Λ y}} ρ(y) dy = ∫ K(x, y) ρ'(y) dy,\]
with
\[ ρ'(y) = w(y) ρ(y),\]
where the weight function is given by
\[ w(y) = \frac{1}{\tanh(Λ y/2)}.\]
SparseIR.MatsubaraConst — Type
MatsubaraConst{S} <: AbstractAugmentation{S}Constant in Matsubara, undefined in imaginary time.
Type Parameters
S: Statistics type (Fermionic or Bosonic). This is required for type consistency, though MatsubaraConst works identically for both statistics.
SparseIR.MatsubaraFreq — Type
MatsubaraFreq(n)Prefactor n of the Matsubara frequency ω = n*π/β
Struct representing the Matsubara frequency ω entering the Fourier transform of a propagator G(τ) on imaginary time τ to its Matsubara equivalent Ĝ(iω) on the imaginary-frequency axis:
β
Ĝ(iω) = ∫ dτ exp(iωτ) G(τ) with ω = n π/β,
0where β is inverse temperature and by convention we include the imaginary unit in the frequency argument, i.e, Ĝ(iω). The frequencies depend on the statistics of the propagator, i.e., we have that:
G(τ - β) = ± G(τ)where + is for bosons and - is for fermions. The frequencies are restricted accordingly.
- Bosonic frequency (
S == Fermionic):neven (periodic in β) - Fermionic frequency (
S == Bosonic):nodd (anti-periodic in β)
SparseIR.MatsubaraSampling — Type
MatsubaraSampling{T,B} <: AbstractSampling
Sparse sampling in Matsubara frequencies using the C API.
Allows transformation between IR basis coefficients and sampling points in Matsubara frequencies.
SparseIR.MatsubaraSampling — Method
MatsubaraSampling(basis::AbstractBasis; positive_only=false, sampling_points=nothing, factorize=true)Construct a MatsubaraSampling object from a basis. If sampling_points is not provided, the default Matsubara sampling points from the basis are used.
If positive_only=true, assumes functions are symmetric in Matsubara frequency.
SparseIR.RegularizedBoseKernel — Type
RegularizedBoseKernel <: AbstractKernelRegularized bosonic analytical continuation kernel.
In dimensionless variables $x = 2 τ/β - 1$, $y = β ω/Λ$, the fermionic integral kernel is a function on $[-1, 1] × [-1, 1]$:
\[ K(x, y) = y \frac{e^{-Λ y (x + 1) / 2}}{e^{-Λ y} - 1}\]
Care has to be taken in evaluating this expression around $y = 0$.
SparseIR.TauConst — Type
TauConst{S} <: AbstractAugmentation{S}Constant function in imaginary time with statistics-dependent periodicity.
Type Parameters
S: Statistics type (Fermionic or Bosonic)
SparseIR.TauLinear — Type
TauLinear{S} <: AbstractAugmentation{S}Linear function in imaginary time, antisymmetric around β/2, with statistics-dependent periodicity.
Type Parameters
S: Statistics type (Fermionic or Bosonic)
SparseIR.TauSampling — Type
TauSampling{T,B} <: AbstractSampling
Sparse sampling in imaginary time using the C API.
Allows transformation between IR basis coefficients and sampling points in imaginary time.
SparseIR.TauSampling — Method
TauSampling(basis::AbstractBasis; sampling_points=nothing, use_positive_taus=true)Construct a TauSampling object from a basis. If sampling_points is not provided, the default tau sampling points from the basis are used.
If use_positive_taus=true, the sampling points are folded to the positive tau domain [0, β) [default].
If use_positive_taus=false, the sampling points are in the range [-β/2, β/2].
SparseIR.default_omega_sampling_points — Method
default_omega_sampling_points(basis::AbstractBasis)Get the default real-frequency sampling points for a basis.
These are the extrema of the highest-order basis function on the real-frequency axis, which provide near-optimal conditioning for the DLR.
SparseIR.evaluate! — Method
evaluate!(output::Array, sampling::AbstractSampling, al::Array; dim=1)In-place version of evaluate. Write results to the pre-allocated output array.
SparseIR.evaluate — Method
evaluate(sampling::AbstractSampling, al::Array; dim=1)Evaluate basis coefficients at the sampling points using the C API.
For multidimensional arrays, dim specifies which dimension corresponds to the basis coefficients.
SparseIR.fit! — Method
fit!(output::Array, sampling::AbstractSampling, al::Array; dim=1)In-place version of fit. Write results to the pre-allocated output array.
SparseIR.fit — Method
fit(sampling::AbstractSampling, al::Array; dim=1)Fit basis coefficients from values at sampling points using the C API.
For multidimensional arrays, dim specifies which dimension corresponds to the sampling points.
SparseIR.from_IR — Method
from_IR(dlr::DiscreteLehmannRepresentation, gl::Array, dims=1)Transform from IR basis coefficients to DLR coefficients.
Arguments
dlr: The DLR basisgl: IR basis coefficientsdims: Dimension along which the basis coefficients are stored
Returns
DLR coefficients with the same shape as input, but with size length(dlr) along dimension dims.
SparseIR.get_poles — Method
get_poles(dlr::DiscreteLehmannRepresentation)Get the pole locations for the DLR basis.
Returns a vector of pole locations on the real-frequency axis.
SparseIR.iscentrosymmetric — Method
iscentrosymmetric(kernel::AbstractKernel)Return whether the kernel satisfies K(x, y) == K(-x, -y) for all values of x and y. Defaults to false.
A centrosymmetric kernel can be block-diagonalized, speeding up the singular value expansion by a factor of 4.
SparseIR.npoints — Method
npoints(sampling::AbstractSampling)Get the number of sampling points.
SparseIR.npoles — Method
npoles(dlr::DiscreteLehmannRepresentation)Get the number of poles in the DLR basis.
SparseIR.overlap — Method
overlap(poly::PiecewiseLegendrePoly, f, xmin::Float64, xmax::Float64;
rtol=eps(), return_error=false, maxevals=10^4, points=Float64[])Evaluate overlap integral of poly with arbitrary function f.
Given the function f, evaluate the integral
∫ dx f(x) poly(x)using adaptive Gauss-Legendre quadrature.
points is a sequence of break points in the integration interval where local difficulties of the integrand may occur (e.g. singularities, discontinuities).
SparseIR.overlap — Method
overlap(poly::PiecewiseLegendrePoly, f;
rtol=eps(), return_error=false, maxevals=10^4, points=Float64[])Evaluate overlap integral of poly with arbitrary function f using default range.
Given the function f, evaluate the integral
∫ dx f(x) poly(x)using adaptive Gauss-Legendre quadrature with the default integration range.
points is a sequence of break points in the integration interval where local difficulties of the integrand may occur (e.g. singularities, discontinuities).
SparseIR.overlap — Method
overlap(polys::PiecewiseLegendrePolyVector, f;
rtol=eps(), return_error=false, maxevals=10^4, points=Float64[])Evaluate overlap integral of polys with arbitrary function f using default range.
Given the function f, evaluate the integral
∫ dx f(x) polys[i](x)for each polynomial in the vector using adaptive Gauss-Legendre quadrature with the default integration range.
SparseIR.sampling_points — Method
sampling_points(sampling::AbstractSampling)Return sampling points.
SparseIR.to_IR — Method
to_IR(dlr::DiscreteLehmannRepresentation, g_dlr::Array, dims=1)Transform from DLR coefficients to IR basis coefficients.
Arguments
dlr: The DLR basisg_dlr: DLR coefficientsdims: Dimension along which the DLR coefficients are stored
Returns
IR basis coefficients with the same shape as input, but with size length(dlr.basis) along dimension dims.
- wallerberger2021https://doi.org/10.1103/PhysRevResearch.3.033168
- shinaoka2018https://doi.org/10.1103/PhysRevB.97.205111