Public names index

SparseIR.AugmentedBasisType
AugmentedBasis <: AbstractBasis

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

Warning

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]

source
SparseIR.DiscreteLehmannRepresentationType
DiscreteLehmannRepresentation <: AbstractBasis

Discrete Lehmann representation (DLR) with poles selected according to extrema of IR.

This class implements a variant of the discrete Lehmann representation (DLR) 1. Instead of a truncated singular value expansion of the analytic continuation kernel $K$ like the IR, the discrete Lehmann representation is based on a "sketching" of $K$. The resulting basis is a linear combination of discrete set of poles on the real-frequency axis, continued to the imaginary-frequency axis:

 G(iv) == sum(a[i] / (iv - w[i]) for i in range(L))

Warning The poles on the real-frequency axis selected for the DLR are based on a rank-revealing decomposition, which offers accuracy guarantees. Here, we instead select the pole locations based on the zeros of the IR basis functions on the real axis, which is a heuristic. We do not expect that difference to matter, but please don't blame the DLR authors if we were wrong :-)

source
SparseIR.FiniteTempBasisType
FiniteTempBasis <: AbstractBasis

Intermediate 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 function u(x). To obtain a single basis function, a slice or a subset l, you can use u[l].

  • uhat::PiecewiseLegendreFT: 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 function uhat(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 subset l, you can use uhat[l].

  • s: Vector of singular values of the continuation kernel

  • v::PiecewiseLegendrePoly: 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 function v(w). To obtain a single basis function, a slice or a subset l, you can use v[l].

source
SparseIR.FiniteTempBasisMethod
FiniteTempBasis{S}(β, ωmax, ε=nothing; max_size=nothing, args...)

Construct a finite temperature basis suitable for the given S (Fermionic or Bosonic) and cutoffs β and ωmax.

source
SparseIR.FiniteTempBasisSetType
FiniteTempBasisSet

Type 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
source
SparseIR.LogisticKernelType
LogisticKernel <: AbstractKernel

Fermionic/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)}.\]

source
SparseIR.MatsubaraFreqType
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 π/β,
        0

where β 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): n even (periodic in β)
  • Fermionic frequency (S == Bosonic): n odd (anti-periodic in β)
source
SparseIR.MatsubaraSamplingType
MatsubaraSampling <: AbstractSampling

Sparse sampling in Matsubara frequencies.

Allows the transformation between the IR basis and a set of sampling points in (scaled/unscaled) imaginary frequencies.

source
SparseIR.MatsubaraSamplingMethod
MatsubaraSampling(basis; positive_only=false,
                  sampling_points=default_matsubara_sampling_points(basis; positive_only),
                  factorize=true)

Construct a MatsubaraSampling object. If not given, the sampling_points are chosen as the (discrete) extrema of the highest-order basis function in Matsubara. This turns out to be close to optimal with respect to conditioning for this size (within a few percent).

By setting positive_only=true, one assumes that functions to be fitted are symmetric in Matsubara frequency, i.e.:

\[ Ĝ(iν) = conj(Ĝ(-iν))\]

or equivalently, that they are purely real in imaginary time. In this case, sparse sampling is performed over non-negative frequencies only, cutting away half of the necessary sampling space. factorize controls whether the SVD decomposition is computed.

source
SparseIR.RegularizedBoseKernelType
RegularizedBoseKernel <: AbstractKernel

Regularized 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$.

source
SparseIR.TauConstType
TauConst <: AbstractAugmentation

Constant in imaginary time/discrete delta in frequency.

source
SparseIR.TauLinearType
TauLinear <: AbstractAugmentation

Linear function in imaginary time, antisymmetric around β/2.

source
SparseIR.TauSamplingType
TauSampling <: AbstractSampling

Sparse sampling in imaginary time.

Allows the transformation between the IR basis and a set of sampling points in (scaled/unscaled) imaginary time.

source
SparseIR.TauSamplingMethod
TauSampling(basis; sampling_points=default_tau_sampling_points(basis), factorize=true)

Construct a TauSampling object. If not given, the sampling_points are chosen as the extrema of the highest-order basis function in imaginary time. This turns out to be close to optimal with respect to conditioning for this size (within a few percent). factorize controls whether the SVD decomposition is computed.

source
SparseIR.evaluate!Method
evaluate!(buffer::AbstractArray{T,N}, sampling, al; dim=1) where {T,N}

Like evaluate, but write the result to buffer. Please use dim = 1 or N to avoid allocating large temporary arrays internally.

source
SparseIR.evaluateMethod
evaluate(sampling, al; dim=1)

Evaluate the basis coefficients al at the sparse sampling points.

source
SparseIR.fit!Method
fit!(buffer::Array{S,N}, smpl::AbstractSampling, al::Array{T,N}; 
    dim=1, workarr::Vector{S}) where {S,T,N}

Like fit, but write the result to buffer. Use dim = 1 or dim = N to avoid allocating large temporary arrays internally. The length of workarr cannot be smaller than SparseIR.workarrlength(smpl, al).

source
SparseIR.fitMethod
fit(sampling, al::AbstractArray{T,N}; dim=1)

Fit basis coefficients from the sparse sampling points Please use dim = 1 or N to avoid allocating large temporary arrays internally.

source
SparseIR.overlapMethod
overlap(poly::PiecewiseLegendrePoly, f; 
    rtol=eps(T), return_error=false, maxevals=10^4, points=T[])

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

source