SparseIR

Documentation for SparseIR.

SparseIR.SparseIRModule

Intermediate representation (IR) for many-body propagators

SparseIR.DimensionlessBasisType
DimensionlessBasis <: AbstractBasis

Intermediate representation (IR) basis in reduced variables.

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(x, y) ≈ sum(u[l](x) * s[l] * v[l](y) for l in range(L))

The functions are given in reduced variables, x = 2τ/β - 1 and y = ω/ωmax, which scales both sides to the interval [-1, 1]. The kernel then only depends on a cutoff parameter Λ = β * ωmax.

Examples

The following example code assumes the spectral function is a single pole at x = 0.2. We first compute an IR basis suitable for fermions and β*W ≤ 42. Then we get G(iw) on the first few Matsubara frequencies:

julia> using SparseIR

julia> basis = DimensionlessBasis(fermion, 42);

julia> gl = basis.s .* basis.v(0.2);

julia> giw = transpose(basis.uhat([1, 3, 5, 7])) * gl

Fields

  • u::PiecewiseLegendrePolyVector: Set of IR basis functions on the reduced imaginary time (x) 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::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 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::PiecewiseLegendrePolyVector: Set of IR basis functions on the reduced real frequency (y) 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 y, you can call the function v(y). To obtain a single basis function, a slice or a subset l, you can use v[l].

See also FiniteTempBasis for a basis directly in time/frequency.

SparseIR.DimensionlessBasisType
DimensionlessBasis(statistics, Λ, ε=nothing; kernel=LogisticKernel(Λ), sve_result=compute_sve(kernel; ε))

Construct an IR basis suitable for the given statistics and cutoff Λ.

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

Examples

The following example code assumes the spectral function is a single pole at ω = 2.5. We first compute an IR basis suitable for fermions and β = 10, W ≤ 4.2. Then we get G(iw) on the first few Matsubara frequencies:

julia> using SparseIR

julia> basis = FiniteTempBasis(fermion, 42, 4.2);

julia> gl = basis.s .* basis.v(2.5);

julia> giw = transpose(basis.uhat([1, 3, 5, 7])) * gl

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

SparseIR.FiniteTempBasisType
FiniteTempBasis(statistics, β, wmax, ε=nothing; kernel=LogisticKernel(β * wmax), sve_result=compute_sve(kernel; ε))

Construct a finite temperature basis suitable for the given statistics and cutoffs β and wmax.

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
  • beta::Float64: Inverse temperature
  • wmax::Float64: Cut-off frequency
  • 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
SparseIR.FiniteTempBasisSetMethod
FiniteTempBasisSet(β, wmax, ε; sve_result=compute_sve(LogisticKernel(β * wmax); ε))

Create basis sets for fermion and boson and associated sampling objects. Fermion and bosonic bases are constructed by SVE of the logistic kernel.

SparseIR.LegendreBasisType

Legendre basis

In the original paper [L. Boehnke et al., PRB 84, 075145 (2011)], they used:

G(\tau) = \sum_{l=0} \sqrt{2l+1} P_l[x(\tau)] G_l/beta,

where P_l[x] is the $l$-th Legendre polynomial.

In this type, the basis functions are defined by

U_l(\tau) \equiv c_l (\sqrt{2l+1}/beta) * P_l[x(\tau)],

where cl are additional l-dependent constant factors. By default, we take cl = 1, which reduces to the original definition.

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

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.

SparseIR.MatsubaraSamplingType
MatsubaraSampling(basis[, sampling_points])

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

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

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.

SparseIR.TauSamplingType
TauSampling(basis[, sampling_points])

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

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.

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

Evaluate the basis coefficients al at the sparse sampling points.

SparseIR.fit!Method
fit!(buffer, sampling, al::Array{T,N}; dim=1)

Like fit, but write the result to buffer. Please use dim = 1 or N to avoid allocating large temporary arrays internally. The length of workarry cannot be smaller than the returned value of workarrlengthfit.

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.

SparseIR.overlapMethod
overlap(poly::PiecewiseLegendrePoly, f)

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.

SparseIR.to_IRFunction

From SPR to IR

g_spr: Expansion coefficients in SPR