Private names index

These are not considered API and therefore not covered by any semver promises.

Core.IntMethod

Get prefactor n for the Matsubara frequency ω = n*π/β

source
Core.UnionMethod
(polyFT::PiecewiseLegendreFT)(ω)

Obtain Fourier transform of polynomial for given MatsubaraFreq ω.

source
SparseIR.AbstractAugmentationType
AbstractAugmentation

Scalar function in imaginary time/frequency.

This represents a single function in imaginary time and frequency, together with some auxiliary methods that make it suitable for augmenting a basis.

See also: AugmentedBasis

source
SparseIR.AbstractBasisType
AbstractBasis

Abstract base class for bases on the imaginary-time axis.

Let basis be an abstract basis. Then we can expand a two-point propagator G(τ), where τ is imaginary time, into a set of basis functions:

G(τ) == sum(basis.u[l](τ) * g[l] for l in 1:length(basis)) + ϵ(τ),

where basis.u[l] is the l-th basis function, g[l] is the associated expansion coefficient and ϵ(τ) is an error term. Similarly, the Fourier transform Ĝ(n), where n is now a Matsubara frequency, can be expanded as follows:

Ĝ(n) == sum(basis.uhat[l](n) * g[l] for l in 1:length(basis)) + ϵ(n),

where basis.uhat[l] is now the Fourier transform of the basis function.

source
SparseIR.AbstractKernelType
(kernel::AbstractKernel)(x, y[, x₊, x₋])

Evaluate kernel at point (x, y).

The parameters x₊ and x₋, if given, shall contain the values of x - xₘᵢₙ and xₘₐₓ - x, respectively. This is useful if either difference is to be formed and cancellation expected.

source
SparseIR.AbstractKernelType
AbstractKernel

Integral kernel K(x, y).

Abstract base type for an integral kernel, i.e. a AbstractFloat binary function $K(x, y)$ used in a Fredhold integral equation of the first kind:

\[ u(x) = ∫ K(x, y) v(y) dy\]

where $x ∈ [x_\mathrm{min}, x_\mathrm{max}]$ and $y ∈ [y_\mathrm{min}, y_\mathrm{max}]$. For its SVE to exist, the kernel must be square-integrable, for its singular values to decay exponentially, it must be smooth.

In general, the kernel is applied to a scaled spectral function $ρ'(y)$ as:

\[ ∫ K(x, y) ρ'(y) dy,\]

where $ρ'(y) = w(y) ρ(y)$.

source
SparseIR.AbstractSamplingType
AbstractSampling

Abstract type for sparse sampling.

Encodes the "basis transformation" of a propagator from the truncated IR basis coefficients G_ir[l] to time/frequency sampled on sparse points G(x[i]) together with its inverse, a least squares fit:

     ________________                   ___________________
    |                |    evaluate     |                   |
    |     Basis      |---------------->|     Value on      |
    |  coefficients  |<----------------|  sampling points  |
    |________________|      fit        |___________________|
source
SparseIR.CentrosymmSVEType
CentrosymmSVE <: AbstractSVE

SVE of centrosymmetric kernel in block-diagonal (even/odd) basis.

For a centrosymmetric kernel K, i.e., a kernel satisfying: K(x, y) == K(-x, -y), one can make the following ansatz for the singular functions:

u[l](x) = ured[l](x) + sign[l] * ured[l](-x)
v[l](y) = vred[l](y) + sign[l] * ured[l](-y)

where sign[l] is either +1 or -1. This means that the singular value expansion can be block-diagonalized into an even and an odd part by (anti-)symmetrizing the kernel:

K_even = K(x, y) + K(x, -y)
K_odd  = K(x, y) - K(x, -y)

The lth basis function, restricted to the positive interval, is then the singular function of one of these kernels. If the kernel generates a Chebyshev system [1], then even and odd basis functions alternate.

[1]: A. Karlin, Total Positivity (1968).

source
SparseIR.LogisticKernelOddType
LogisticKernelOdd <: AbstractReducedKernel

Fermionic analytical continuation kernel, odd.

In dimensionless variables $x = 2τ/β - 1$, $y = βω/Λ$, the fermionic integral kernel is a function on $[-1, 1] × [-1, 1]$:

\[ K(x, y) = -\frac{\sinh(Λ x y / 2)}{\cosh(Λ y / 2)}\]

source
SparseIR.PiecewiseLegendreFTType
PiecewiseLegendreFT <: Function

Fourier transform of a piecewise Legendre polynomial.

For a given frequency index n, the Fourier transform of the Legendre function is defined as:

    p̂(n) == ∫ dx exp(im * π * n * x / (xmax - xmin)) p(x)

The polynomial is continued either periodically (freq=:even), in which case n must be even, or antiperiodically (freq=:odd), in which case n must be odd.

source
SparseIR.PiecewiseLegendrePolyType
PiecewiseLegendrePoly <: Function

Piecewise Legendre polynomial.

Models a function on the interval $[xmin, xmax]$ as a set of segments on the intervals $S[i] = [a[i], a[i+1]]$, where on each interval the function is expanded in scaled Legendre polynomials.

source
SparseIR.PowerModelType
PowerModel

Model from a high-frequency series expansion::

A(iω) == sum(A[n] / (iω)^(n+1) for n in 1:N)

where $iω == i * π/2 * wn$ is a reduced imaginary frequency, i.e., $wn$ is an odd/even number for fermionic/bosonic frequencies.

source
SparseIR.ReducedKernelType
ReducedKernel

Restriction of centrosymmetric kernel to positive interval.

For a kernel $K$ on $[-1, 1] × [-1, 1]$ that is centrosymmetric, i.e. $K(x, y) = K(-x, -y)$, it is straight-forward to show that the left/right singular vectors can be chosen as either odd or even functions.

Consequentially, they are singular functions of a reduced kernel $K_\mathrm{red}$ on $[0, 1] × [0, 1]$ that is given as either:

\[ K_\mathrm{red}(x, y) = K(x, y) \pm K(x, -y)\]

This kernel is what this type represents. The full singular functions can be reconstructed by (anti-)symmetrically continuing them to the negative axis.

source
SparseIR.RegularizedBoseKernelOddType
RegularizedBoseKernelOdd <: AbstractReducedKernel

Bosonic analytical continuation kernel, odd.

In dimensionless variables $x = 2 τ / β - 1$, $y = β ω / Λ$, the fermionic integral kernel is a function on $[-1, 1] × [-1, 1]$:

\[ K(x, y) = -y \frac{\sinh(Λ x y / 2)}{\sinh(Λ y / 2)}\]

source
SparseIR.RuleType
Rule{T<:AbstractFloat}

Quadrature rule.

Approximation of an integral over [a, b] by a sum over discrete points x with weights w:

\[ ∫ f(x) ω(x) dx ≈ ∑_i f(x_i) w_i\]

where we generally have superexponential convergence for smooth $f(x)$ in the number of quadrature points.

source
SparseIR.SVEResultMethod
SVEResult(kernel::AbstractKernel;
    Twork=nothing, ε=nothing, lmax=typemax(Int),
    n_gauss=nothing, svd_strat=:auto,
    sve_strat=iscentrosymmetric(kernel) ? CentrosymmSVE : SamplingSVE
)

Perform truncated singular value expansion of a kernel.

Perform a truncated singular value expansion (SVE) of an integral kernel kernel : [xmin, xmax] x [ymin, ymax] -> ℝ:

kernel(x, y) == sum(s[l] * u[l](x) * v[l](y) for l in (1, 2, 3, ...)),

where s[l] are the singular values, which are ordered in non-increasing fashion, u[l](x) are the left singular functions, which form an orthonormal system on [xmin, xmax], and v[l](y) are the right singular functions, which form an orthonormal system on [ymin, ymax].

The SVE is mapped onto the singular value decomposition (SVD) of a matrix by expanding the kernel in piecewise Legendre polynomials (by default by using a collocation).

Arguments

  • K::AbstractKernel: Integral kernel to take SVE from.

  • ε::Real: Accuracy target for the basis: attempt to have singular values down to a relative magnitude of ε, and have each singular value and singular vector be accurate to ε. A Twork with a machine epsilon of ε^2 or lower is required to satisfy this. Defaults to 2.2e-16 if xprec is available, and 1.5e-8 otherwise.

  • cutoff::Real: Relative cutoff for the singular values. A Twork with machine epsilon of cutoff is required to satisfy this. Defaults to a small multiple of the machine epsilon.

    Note that cutoff and ε serve distinct purposes. cutoff reprsents the accuracy to which the kernel is reproduced, whereas ε is the accuracy to which the singular values and vectors are guaranteed.

  • lmax::Integer: Maximum basis size. If given, only at most the lmax most significant singular values and associated singular functions are returned.

  • `n_gauss (int): Order of Legendre polynomials. Defaults to kernel hinted value.

  • Twork: Working data type. Defaults to a data type with machine epsilon of at mostε^2and at mostcutoff`, or otherwise most accurate data type available.

  • sve_strat::AbstractSVE: SVE to SVD translation strategy. Defaults to SamplingSVE, optionally wrapped inside of a CentrosymmSVE if the kernel is centrosymmetric.

  • svd_strat ('fast' or 'default' or 'accurate'): SVD solver. Defaults to fast (ID/RRQR) based solution when accuracy goals are moderate, and more accurate Jacobi-based algorithm otherwise.

Returns: An SVEResult containing the truncated singular value expansion.

source
SparseIR.SamplingSVEType
SamplingSVE <: AbstractSVE

SVE to SVD translation by sampling technique [1].

Maps the singular value expansion (SVE) of a kernel kernel onto the singular value decomposition of a matrix A. This is achieved by choosing two sets of Gauss quadrature rules: (x, wx) and (y, wy) and approximating the integrals in the SVE equations by finite sums. This implies that the singular values of the SVE are well-approximated by the singular values of the following matrix:

A[i, j] = √(wx[i]) * K(x[i], y[j]) * √(wy[j])

and the values of the singular functions at the Gauss sampling points can be reconstructed from the singular vectors u and v as follows:

u[l,i] ≈ √(wx[i]) u[l](x[i])
v[l,j] ≈ √(wy[j]) u[l](y[j])

[1] P. Hansen, Discrete Inverse Problems, Ch. 3.1

source
SparseIR.accuracyFunction
accuracy(basis::AbstractBasis)

Accuracy of the basis.

Upper bound to the relative error of reprensenting a propagator with the given number of basis functions (number between 0 and 1).

source
SparseIR.canonicalize!Method
canonicalize!(u, v)

Canonicalize basis.

Each SVD (u[l], v[l]) pair is unique only up to a global phase, which may differ from implementation to implementation and also platform. We fix that gauge by demanding u[l](1) > 0. This ensures a diffeomorphic connection to the Legendre polynomials as Λ → 0.

source
SparseIR.conv_radiusFunction
conv_radius(kernel)

Convergence radius of the Matsubara basis asymptotic model.

For improved relative numerical accuracy, the IR basis functions on the Matsubara axis uhat(basis, n) can be evaluated from an asymptotic expression for abs(n) > conv_radius. If isinf(conv_radius), then the asymptotics are unused (the default).

source
SparseIR.default_matsubara_sampling_pointsFunction
default_matsubara_sampling_points(basis::AbstractBasis; positive_only=false)

Default sampling points on the imaginary frequency axis.

Arguments

  • positive_only::Bool: Only return non-negative frequencies. This is useful if the object to be fitted is symmetric in Matsubura frequency, ĝ(ω) == conj(ĝ(-ω)), or, equivalently, real in imaginary time.
source
SparseIR.eval_matrixFunction
eval_matrix(T, basis, x)

Return evaluation matrix from coefficients to sampling points. T <: AbstractSampling.

source
SparseIR.find_extremaMethod
find_extrema(polyFT::PiecewiseLegendreFT; part=nothing, grid=DEFAULT_GRID)

Obtain extrema of Fourier-transformed polynomial.

source
SparseIR.finite_temp_basesFunction
finite_temp_bases(β::Real, ωmax::Real, ε=nothing;
                  kernel=LogisticKernel(β * ωmax), sve_result=SVEResult(kernel; ε))

Construct FiniteTempBasis objects for fermion and bosons using the same LogisticKernel instance.

source
SparseIR.from_IRFunction
from_IR(dlr::DiscreteLehmannRepresentation, gl::AbstractArray, dims=1)

From IR to DLR. gl`: Expansion coefficients in IR.

source
SparseIR.get_symmetrizedMethod
get_symmetrized(kernel, sign)

Construct a symmetrized version of kernel, i.e. kernel(x, y) + sign * kernel(x, -y).

Beware!

By default, this returns a simple wrapper over the current instance which naively performs the sum. You may want to override this to avoid cancellation.

source
SparseIR.get_tnlMethod
get_tnl(l, w)

Fourier integral of the l-th Legendre polynomial::

Tₗ(ω) == ∫ dx exp(iωx) Pₗ(x)
source
SparseIR.giwMethod
giw(polyFT, wn)

Return model Green's function for reduced frequencies

source
SparseIR.iscentrosymmetricFunction
iscentrosymmetric(kernel)

Return true if kernel(x, y) == kernel(-x, -y) for all values of x and y in range. This allows the kernel to be block-diagonalized, speeding up the singular value expansion by a factor of 4. Defaults to false.

source
SparseIR.matop!Method
matop!(buffer, mat, arr::AbstractArray, op, dim)

Apply the operator op to the matrix mat and to the array arr along the first dimension (dim=1) or the last dimension (dim=N).

source
SparseIR.matop_along_dim!Method
matop_along_dim!(buffer, mat, arr::AbstractArray, dim::Integer, op)

Apply the operator op to the matrix mat and to the array arr along the dimension dim, writing the result to buffer.

source
SparseIR.movedimMethod
movedim(arr::AbstractArray, src => dst)

Move arr's dimension at src to dst while keeping the order of the remaining dimensions unchanged.

source
SparseIR.nsvalsMethod
nsvals(hints)

Upper bound for number of singular values.

Upper bound on the number of singular values above the given threshold, i.e. where s[l] ≥ ε * first(s).

source
SparseIR.phase_stableMethod
phase_stable(poly, wn)

Phase factor for the piecewise Legendre to Matsubara transform.

Compute the following phase factor in a stable way:

exp.(iπ/2 * wn * cumsum(poly.Δx))
source
SparseIR.piecewiseMethod
piecewise(rule, edges)

Piecewise quadrature with the same quadrature rule, but scaled.

source
SparseIR.rescaleMethod
rescale(basis::FiniteTempBasis, new_β)

Return a basis for different temperature.

Uses the same kernel with the same $ε$, but a different temperature. Note that this implies a different UV cutoff $ωmax$, since $Λ == β * ωmax$ stays constant.

source
SparseIR.segments_xMethod
segments_x(sve_hints::AbstractSVEHints[, T])

Segments for piecewise polynomials on the $x$ axis.

List of segments on the $x$ axis for the associated piecewise polynomial. Should reflect the approximate position of roots of a high-order singular function in $x$.

source
SparseIR.segments_yMethod
segments_y(sve_hints::AbstractSVEHints[, T])

Segments for piecewise polynomials on the $y$ axis.

List of segments on the $y$ axis for the associated piecewise polynomial. Should reflect the approximate position of roots of a high-order singular function in $y$.

source
SparseIR.shift_xmidMethod
shift_xmid(knots, Δx)

Return midpoint relative to the nearest integer plus a shift.

Return the midpoints xmid of the segments, as pair (diff, shift), where shift is in (0, 1, -1) and diff is a float such that xmid == shift + diff to floating point accuracy.

source
SparseIR.significanceFunction
significance(basis::AbstractBasis)

Return vector σ, where 0 ≤ σ[i] ≤ 1 is the significance level of the i-th basis function. If ϵ is the desired accuracy to which to represent a propagator, then any basis function where σ[i] < ϵ can be neglected.

For the IR basis, we simply have that σ[i] = s[i] / first(s).

source
SparseIR.splitMethod
split(poly, x)

Split segment.

Find segment of poly's domain that covers x.

source
SparseIR.statisticsMethod
statistics(basis::AbstractBasis)

Quantum statistic (Statistics instance, Fermionic() or Bosonic()).

source
SparseIR.sve_hintsFunction
sve_hints(kernel, ε)

Provide discretisation hints for the SVE routines.

Advises the SVE routines of discretisation parameters suitable in tranforming the (infinite) SVE into an (finite) SVD problem.

See also AbstractSVEHints.

source
SparseIR.to_IRFunction
to_IR(dlr::DiscreteLehmannRepresentation, g_dlr::AbstractArray, dims=1)

From DLR to IR. g_dlr`: Expansion coefficients in DLR.

source
SparseIR.truncateMethod
truncate(u, s, v; rtol=0.0, lmax=typemax(Int))

Truncate singular value expansion.

Arguments

- `u`, `s`, `v`: Thin singular value expansion
- `rtol`: Only singular values satisfying `s[l]/s[1] > rtol` are retained.
- `lmax`: At most the `lmax` most significant singular values are retained.
source
SparseIR.weight_funcFunction
weight_func(kernel, statistics::Statistics)

Return the weight function for the given statistics.

  • Fermion: w(x) == 1
  • Boson: w(y) == 1/tanh(Λ*y/2)
source
SparseIR.xrangeFunction
xrange(kernel)

Return a tuple $(x_\mathrm{min}, x_\mathrm{max})$ delimiting the range of allowed x values.

source
SparseIR.yrangeFunction
yrange(kernel)

Return a tuple $(y_\mathrm{min}, y_\mathrm{max})$ delimiting the range of allowed y values.

source
SparseIR.ΛFunction
Λ(basis::AbstractBasis)
lambda(basis::AbstractBasis)

Basis cutoff parameter, Λ = β * ωmax, or None if not present

source
SparseIR.βMethod
β(basis::AbstractBasis)
beta(basis::AbstractBasis)

Inverse temperature or nothing if unscaled basis.

source
SparseIR.ωmaxFunction
ωmax(basis::AbstractBasis)
wmax(basis::AbstractBasis)

Real frequency cutoff or nothing if unscaled basis.

source
SparseIR._LinAlg.rrqr!Method

Truncated rank-revealing QR decomposition with full column pivoting.

Decomposes a (m, n) matrix A into the product:

A[:,piv] == Q * R

where Q is an (m, k) isometric matrix, R is a (k, n) upper triangular matrix, piv is a permutation vector, and k is chosen such that the relative tolerance tol is met in the equality above.

source
SparseIR._LinAlg.svd2x2Method

Perform the SVD of an arbitrary two-by-two matrix:

  [ a11  a12 ]  =  [  cu  -su ] [ smax     0 ] [  cv   sv ]
  [ a21  a22 ]     [  su   cu ] [    0  smin ] [ -sv   cv ]

Note that smax and smin can be negative.

source
SparseIR._LinAlg.svd2x2Method

Perform the SVD of upper triangular two-by-two matrix:

  [ f    g   ]  =  [  cu  -su ] [ smax     0 ] [  cv   sv ]
  [ 0    h   ]     [  su   cu ] [    0  smin ] [ -sv   cv ]

Note that smax and smin can be negative.

source
SparseIR._LinAlg.tsvd!Method

Truncated singular value decomposition.

Decomposes an (m, n) matrix A into the product:

A == U * (s .* VT)

where U is a (m, k) matrix with orthogonal columns, VT is a (k, n) matrix with orthogonal rows and s are the singular values, a set of k nonnegative numbers in non-ascending order. The SVD is truncated in the sense that singular values below tol are discarded.

source