Sparse sampling

Sparse sampling is a way to transition between the truncated IR representation of a propagator (useful for convergence analysis) and a sparse set of sampling points in imaginary time or Matsubara frequency. This is mediated by two classes:

  • sparse_ir.TauSampling: sparse sampling in imaginary time, useful for, e.g., constructing Feynman diagrams with a spontaneous interation term.

  • sparse_ir.MatsubaraSampling: sparse sampling in Matsubara frequencies, useful for, e.g., solving diagrammatic equations such as the Dyson equation.

All sampling classes contain sampling_points, which are the corresponding sampling points in time or frequency, and a method evaluate(), which allows one to go from coefficients to sampling points, and a method fit() to go back:

 ________________                   ___________________
|                |    evaluate()   |                   |
|     Basis      |---------------->|     Value on      |
|  coefficients  |<----------------|  sampling_points  |
|________________|      fit()      |___________________|

Warning

When storing data in sparse time/frequency, always store the sampling points together with the data. The exact location of the sampling points may be different from between platforms and/or between releases.

Sparse sampling transformers

class sparse_ir.TauSampling(basis, sampling_points=None)

Sparse sampling in imaginary time.

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

evaluate(al, axis=None, *, points=None)

Evaluate the basis coefficients at sampling points.

Parameters:
  • al (array) – Array where the l-th item along axis corresponds to the l-th basis coefficient

  • axis (integer) – Axis or dimension of al along which to evaluate the function. Defaults to the last, i.e., rightmost axis.

  • points (vector) –

    Points on which the results should be evaluated. Defaults to the sampling points for which the sampling objects was created.

    Added in version 1.1.

Returns:

Array where the n-th item along axis corresponds to the value on the n-th sampling point (or value on point[n], if given.)

Note

If points is given, a new sampling is created at each invocation, which can result in a performance hit. Consider caching sampling objects or simply using the .u() and .uhat() methods of the underlying basis.

fit(ax, axis=None, *, points=None)

Fit the basis coefficients from the sampling points.

Parameters:
  • ax (array) – Array where the n-th item along axis corresponds to the value on the n-th sampling point (or value on point[n], if given.)

  • axis (integer) – Axis or dimension of ax along which to fit the function. Defaults to the last, i.e., rightmost axis.

  • points (vector) –

    Points on which the ax is given. Defaults to the sampling points for which the sampling objects was created.

    Added in version 1.1.

Returns:

Array where the l-th item along axis corresponds to the l-th basis coefficient

Note

If points is given, a new sampling is created at each invocation, which can result in a performance hit. Consider caching sampling objects.

property tau

Sampling points in (reduced) imaginary time

class sparse_ir.MatsubaraSampling(basis, sampling_points=None, *, positive_only=False)

Sparse sampling in Matsubara frequencies.

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

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

Ghat(iv) == Ghat(-iv).conj()

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.

evaluate(al, axis=None, *, points=None)

Evaluate the basis coefficients at sampling points.

Parameters:
  • al (array) – Array where the l-th item along axis corresponds to the l-th basis coefficient

  • axis (integer) – Axis or dimension of al along which to evaluate the function. Defaults to the last, i.e., rightmost axis.

  • points (vector) –

    Points on which the results should be evaluated. Defaults to the sampling points for which the sampling objects was created.

    Added in version 1.1.

Returns:

Array where the n-th item along axis corresponds to the value on the n-th sampling point (or value on point[n], if given.)

Note

If points is given, a new sampling is created at each invocation, which can result in a performance hit. Consider caching sampling objects or simply using the .u() and .uhat() methods of the underlying basis.

fit(ax, axis=None, *, points=None)

Fit the basis coefficients from the sampling points.

Parameters:
  • ax (array) – Array where the n-th item along axis corresponds to the value on the n-th sampling point (or value on point[n], if given.)

  • axis (integer) – Axis or dimension of ax along which to fit the function. Defaults to the last, i.e., rightmost axis.

  • points (vector) –

    Points on which the ax is given. Defaults to the sampling points for which the sampling objects was created.

    Added in version 1.1.

Returns:

Array where the l-th item along axis corresponds to the l-th basis coefficient

Note

If points is given, a new sampling is created at each invocation, which can result in a performance hit. Consider caching sampling objects.

property wn

Sampling points as (reduced) Matsubara frequencies

Base classes

class sparse_ir.sampling.AbstractSampling

Base class 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        |___________________|
property basis

Basis instance

property cond

Condition number of the fitting problem

evaluate(al, axis=None, *, points=None)

Evaluate the basis coefficients at sampling points.

Parameters:
  • al (array) – Array where the l-th item along axis corresponds to the l-th basis coefficient

  • axis (integer) – Axis or dimension of al along which to evaluate the function. Defaults to the last, i.e., rightmost axis.

  • points (vector) –

    Points on which the results should be evaluated. Defaults to the sampling points for which the sampling objects was created.

    Added in version 1.1.

Returns:

Array where the n-th item along axis corresponds to the value on the n-th sampling point (or value on point[n], if given.)

Note

If points is given, a new sampling is created at each invocation, which can result in a performance hit. Consider caching sampling objects or simply using the .u() and .uhat() methods of the underlying basis.

fit(ax, axis=None, *, points=None)

Fit the basis coefficients from the sampling points.

Parameters:
  • ax (array) – Array where the n-th item along axis corresponds to the value on the n-th sampling point (or value on point[n], if given.)

  • axis (integer) – Axis or dimension of ax along which to fit the function. Defaults to the last, i.e., rightmost axis.

  • points (vector) –

    Points on which the ax is given. Defaults to the sampling points for which the sampling objects was created.

    Added in version 1.1.

Returns:

Array where the l-th item along axis corresponds to the l-th basis coefficient

Note

If points is given, a new sampling is created at each invocation, which can result in a performance hit. Consider caching sampling objects.

property matrix

Evaluation matrix is decomposed form

property sampling_points

Set of sampling points

class sparse_ir.sampling.DecomposedMatrix(a, svd_result=None)

Matrix in SVD decomposed form for fast and accurate fitting.

Stores a matrix A together with its thin SVD form:

A == (u * s) @ vH.

This allows for fast and accurate least squares fits using A.lstsq(x).

property a

Full matrix

property cond

Condition number of matrix

lstsq(x, axis=None)

Return y such that np.linalg.norm(A @ y - x) is minimal

matmul(x, axis=None)

Compute A @ x (optionally along specified axis of x)

property s

Most significant, nonzero singular values

property u

Left singular vectors, aranged column-wise

property vH

Right singular vectors, transposed

class sparse_ir.sampling.SplitDecomposedMatrix(a, ssvd_result)

Matrix in “split” SVD decomposed form for fast and accurate fitting.

Stores a matrix A together with its “split SVD” form:

A == u * s @ vT

where vT is a real matrix and u is a complex matrix. The “split” SVD form differs from the SVD in that the least squares fit has to be constructed as follows:

fit(A, x) == vT.T / s @ (u.conj().T * x).real

This again allows for fast and accurate least squares fits using A.lstsq(x). This is useful in the case where.

property a

Full matrix

property cond

Condition number of matrix

lstsq(x, axis=None)

Return y such that np.linalg.norm(A @ y - x) is minimal

matmul(x, axis=None)

Compute A @ x (optionally along specified axis of x)

property s

Most significant, nonzero singular values

property u

Split left singular vectors, aranged column-wise

property vH

Right singular vectors, transposed

class sparse_ir.sampling.ConditioningWarning

Warns about a poorly conditioned problem.

This warning is issued if the library detects a poorly conditioned fitting problem. This essentially means there is a high degree of ambiguity in how to choose the solution. One must therefore expect to lose significant precision in the parameter values.