Private names index
These are not considered API and therefore not covered by any semver promises.
Core.Int — MethodGet prefactor n for the Matsubara frequency ω = n*π/β
Core.Integer — MethodGet prefactor n for the Matsubara frequency ω = n*π/β
Core.Union — Method(polyFT::PiecewiseLegendreFT)(ω)Obtain Fourier transform of polynomial for given MatsubaraFreq ω.
SparseIR.AbstractAugmentation — TypeAbstractAugmentationScalar 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
SparseIR.AbstractBasis — TypeAbstractBasisAbstract 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.
SparseIR.AbstractKernel — Type(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.
SparseIR.AbstractKernel — TypeAbstractKernelIntegral 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)$.
SparseIR.AbstractSVEHints — TypeAbstractSVEHintsDiscretization hints for singular value expansion of a given kernel.
SparseIR.AbstractSampling — TypeAbstractSamplingAbstract 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 |___________________|SparseIR.CentrosymmSVE — TypeCentrosymmSVE <: AbstractSVESVE 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).
SparseIR.LogisticKernelOdd — TypeLogisticKernelOdd <: AbstractReducedKernelFermionic 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)}\]
SparseIR.PiecewiseLegendreFT — TypePiecewiseLegendreFT <: FunctionFourier 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.
SparseIR.PiecewiseLegendrePoly — TypePiecewiseLegendrePoly <: FunctionPiecewise 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.
SparseIR.PiecewiseLegendrePolyVector — TypePiecewiseLegendrePolyVectorAlias for Vector{PiecewiseLegendrePoly}.
SparseIR.PowerModel — TypePowerModelModel 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.
SparseIR.ReducedKernel — TypeReducedKernelRestriction 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.
SparseIR.RegularizedBoseKernelOdd — TypeRegularizedBoseKernelOdd <: AbstractReducedKernelBosonic 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)}\]
SparseIR.Rule — TypeRule{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.
SparseIR.SVEResult — MethodSVEResult(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ε. ATworkwith a machine epsilon ofε^2or lower is required to satisfy this. Defaults to2.2e-16if xprec is available, and1.5e-8otherwise.cutoff::Real: Relative cutoff for the singular values. ATworkwith machine epsilon ofcutoffis required to satisfy this. Defaults to a small multiple of the machine epsilon.Note that
cutoffandεserve distinct purposes.cutoffreprsents 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 thelmaxmost 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 toSamplingSVE, optionally wrapped inside of aCentrosymmSVEif 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.
SparseIR.SamplingSVE — TypeSamplingSVE <: AbstractSVESVE 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
SparseIR.Statistics — TypeStatistics(zeta)Abstract type for quantum statistics (fermionic/bosonic/etc.)
SparseIR.accuracy — Functionaccuracy(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).
SparseIR.canonicalize! — Methodcanonicalize!(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.
SparseIR.choose_accuracy — Methodchoose_accuracy(ε, Twork[, svd_strat])Choose work type and accuracy based on specs and defaults
SparseIR.compute_unl_inner — Methodcompute_unl_inner(poly, wn)Compute piecewise Legendre to Matsubara transform.
SparseIR.conv_radius — Functionconv_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).
SparseIR.default_matsubara_sampling_points — Functiondefault_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.
SparseIR.default_tau_sampling_points — Functiondefault_tau_sampling_points(basis::AbstractBasis)Default sampling points on the imaginary time/x axis.
SparseIR.deriv — Methodderiv(poly[, ::Val{n}=Val(1)])Get polynomial for the nth derivative.
SparseIR.eval_matrix — Functioneval_matrix(T, basis, x)Return evaluation matrix from coefficients to sampling points. T <: AbstractSampling.
SparseIR.find_extrema — Methodfind_extrema(polyFT::PiecewiseLegendreFT; part=nothing, grid=DEFAULT_GRID)Obtain extrema of Fourier-transformed polynomial.
SparseIR.finite_temp_bases — Functionfinite_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.
SparseIR.from_IR — Functionfrom_IR(dlr::DiscreteLehmannRepresentation, gl::AbstractArray, dims=1)From IR to DLR. gl`: Expansion coefficients in IR.
SparseIR.get_symmetrized — Methodget_symmetrized(kernel, sign)Construct a symmetrized version of kernel, i.e. kernel(x, y) + sign * kernel(x, -y).
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.
SparseIR.get_tnl — Methodget_tnl(l, w)Fourier integral of the l-th Legendre polynomial::
Tₗ(ω) == ∫ dx exp(iωx) Pₗ(x)SparseIR.giw — Methodgiw(polyFT, wn)Return model Green's function for reduced frequencies
SparseIR.iscentrosymmetric — Functioniscentrosymmetric(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.
SparseIR.iswellconditioned — Methodiswellconditioned(basis::AbstractBasis)Returns true if the sampling is expected to be well-conditioned.
SparseIR.joinrules — Methodjoinrules(rules)Join multiple Gauss quadratures together.
SparseIR.legder — MethodlegderSparseIR.legendre — Methodlegendre(n[, T])Gauss-Legendre quadrature with n points on [-1, 1].
SparseIR.legendre_collocation — Functionlegendre_collocation(rule, n=length(rule.x))Generate collocation matrix from Gauss-Legendre rule.
SparseIR.legvander — Methodlegvander(x, deg)
Pseudo-Vandermonde matrix of degree deg.
SparseIR.matop! — Methodmatop!(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).
SparseIR.matop_along_dim! — Methodmatop_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.
SparseIR.matrices — Methodmatrices(sve::AbstractSVE)SVD problems underlying the SVE.
SparseIR.matrix_from_gauss — Methodmatrix_from_gauss(kernel, gauss_x, gauss_y)Compute matrix for kernel from Gauss rules.
SparseIR.movedim — Methodmovedim(arr::AbstractArray, src => dst)Move arr's dimension at src to dst while keeping the order of the remaining dimensions unchanged.
SparseIR.ngauss — Functionngauss(hints)Gauss-Legendre order to use to guarantee accuracy.
SparseIR.nsvals — Methodnsvals(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).
SparseIR.phase_stable — Methodphase_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))SparseIR.piecewise — Methodpiecewise(rule, edges)Piecewise quadrature with the same quadrature rule, but scaled.
SparseIR.postprocess — Methodpostprocess(sve::AbstractSVE, u, s, v)Construct the SVE result from the SVD.
SparseIR.rescale — Methodrescale(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.
SparseIR.reseat — Methodreseat(rule, a, b)Reseat quadrature rule to new domain.
SparseIR.roots — Methodroots(poly)Find all roots of the piecewise polynomial poly.
SparseIR.scale — Methodscale(rule, factor)Scale weights by factor.
SparseIR.segments_x — Methodsegments_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$.
SparseIR.segments_y — Methodsegments_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$.
SparseIR.shift_xmid — Methodshift_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.
SparseIR.significance — Functionsignificance(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).
SparseIR.split — Methodsplit(poly, x)Split segment.
Find segment of poly's domain that covers x.
SparseIR.statistics — Methodstatistics(basis::AbstractBasis)Quantum statistic (Statistics instance, Fermionic() or Bosonic()).
SparseIR.sve_hints — Functionsve_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.
SparseIR.to_IR — Functionto_IR(dlr::DiscreteLehmannRepresentation, g_dlr::AbstractArray, dims=1)From DLR to IR. g_dlr`: Expansion coefficients in DLR.
SparseIR.truncate — Methodtruncate(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.SparseIR.value — MethodGet value of the Matsubara frequency ω = n*π/β
SparseIR.valueim — MethodGet complex value of the Matsubara frequency iω = iπ/β * n
SparseIR.weight_func — Functionweight_func(kernel, statistics::Statistics)Return the weight function for the given statistics.
- Fermion:
w(x) == 1 - Boson:
w(y) == 1/tanh(Λ*y/2)
SparseIR.workarrlength — Methodworkarrlength(smpl::AbstractSampling, al; dim=1)Return length of workarr for fit!.
SparseIR.xrange — Functionxrange(kernel)Return a tuple $(x_\mathrm{min}, x_\mathrm{max})$ delimiting the range of allowed x values.
SparseIR.ypower — Functionypower(kernel)Power with which the $y$ coordinate scales.
SparseIR.yrange — Functionyrange(kernel)Return a tuple $(y_\mathrm{min}, y_\mathrm{max})$ delimiting the range of allowed y values.
SparseIR.zeta — MethodGet statistics ζ for Matsubara frequency ω = (2*m+ζ)*π/β
SparseIR.Λ — FunctionΛ(basis::AbstractBasis)
lambda(basis::AbstractBasis)Basis cutoff parameter, Λ = β * ωmax, or None if not present
SparseIR.β — Methodβ(basis::AbstractBasis)
beta(basis::AbstractBasis)Inverse temperature or nothing if unscaled basis.
SparseIR.ωmax — Functionωmax(basis::AbstractBasis)
wmax(basis::AbstractBasis)Real frequency cutoff or nothing if unscaled basis.
SparseIR._LinAlg.givens_lmul — MethodApply Givens rotation to vector:
[ a ] = [ c s ] [ x ]
[ b ] [ -s c ] [ y ]SparseIR._LinAlg.givens_params — MethodCompute Givens rotation R matrix that satisfies:
[ c s ] [ f ] [ r ]
[ -s c ] [ g ] = [ 0 ]SparseIR._LinAlg.rrqr! — MethodTruncated rank-revealing QR decomposition with full column pivoting.
Decomposes a (m, n) matrix A into the product:
A[:,piv] == Q * Rwhere 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.
SparseIR._LinAlg.rrqr — MethodTruncated rank-revealing QR decomposition with full column pivoting.
SparseIR._LinAlg.svd2x2 — MethodPerform 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.
SparseIR._LinAlg.svd2x2 — MethodPerform 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.
SparseIR._LinAlg.svd_jacobi! — MethodSingular value decomposition using Jacobi rotations.
SparseIR._LinAlg.svd_jacobi — MethodSingular value decomposition using Jacobi rotations.
SparseIR._LinAlg.truncate_qr_result — MethodTruncate RRQR result low-rank
SparseIR._LinAlg.tsvd! — MethodTruncated 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.
SparseIR._LinAlg.tsvd — MethodTruncated singular value decomposition.