11. API documentation of ana_cont

11.1. Continuation module

class ana_cont.continuation.AnalyticContinuationProblem(im_axis=None, re_axis=None, im_data=None, kernel_mode=None, beta=None)

Specification of an analytic continuation problem.

This class is designed to hold the information that specifies the analytic continuation problem:

  • Imaginary axis of the data

  • Real axis on which the result is anticipated

  • Imaginary-axis data (real-valued array)

  • Type of kernel (fermionic, bosonic, time, frequency)

  • Inverse temperature beta (only necessary for time kernels)

Furthermore the class provides an interface to call the solver.

Create instance of AnalyticContinuationProblem

Parameters
  • im_axis (numpy.ndarray) – One-dimensional numpy array of type float Matsubara frequencies or imaginary time points of the input data.

  • re_axis (numpy.ndarray) – One-dimensional numpy array of type float Real-frequency grid to use for the analytic continuation.

  • im_data (numpy.ndarray) – One-dimensional numpy array of type float or complex Imaginary-axis data, e.g. Matsubara Green’s function

  • kernel_mode ({‘freq_fermionic’, ‘freq_bosonic’, ‘time_fermionic’, ‘time_bosonic’}) –

    • ‘freq_fermionic’ fermionic Matsubara Greens function

    • ’freq_bosonic’ bosonic Matsubara Greens function

    • ’time_fermionic’ fermionic Green’s function in imaginary time

    • ’time_bosonic’ bosonic Green’s function (susceptibility) in imaginary time

  • beta (float) – Inverse temperature (only necessary for time kernels)

partial_solution(method='', **kwargs)

Maxent optimization at a specific value of alpha.

solve(method='', **kwargs)

Interface function for solving the analytic continuation problem.

Parameters
  • method – Analytic continuation method, possible choices are ‘maxent_svd’, ‘pade’

  • kwargs – Further keyword arguments, specific to the solution method

Returns

OptimizationResult object

This function first creates an instance of an AnalyticContinuationSolver object, then the respective solve function is called.

class ana_cont.continuation.GreensFunction(spectrum=None, wgrid=None, kind='')

This class defines a GreensFunction object. The main use of this is to calculate a full Green’s function with real- and imaginary part from a spectrum.

Initialize with spectral function and real-frequency grid.

Parameters
  • spectrum (np.ndarray) – numpy array containing the spectral function

  • wgrid (np.ndarray) – one-dimensional numpy array of type float real-frequency values corresponding to the spectrum

  • kind ({‘fermionic’, ‘bosonic’}) – specify if the Green’s function is fermionic or bosonic.

kkt()

Kramers Kronig transformation

Obtain full complex Green’s function from its imaginary part.

11.2. Solvers module

class ana_cont.solvers.AnalyticContinuationSolver

Abstract base class for solver classes.

Each Solver class has to inherit from this class. The purpose of this is to ensure that each child class has a method ‘solve’ implemented.

class ana_cont.solvers.MaxentSolverSVD(im_axis, re_axis, im_data, kernel_mode='', model=None, stdev=None, cov=None, offdiag=False, preblur=False, blur_width=0.0, optimizer='newton', verbose=True, **kwargs)

Maxent with singular-value decomposition.

The singular-value decomposition of the kernel leads to a basis in which the dimensionality of the analytic continuation problem is much reduced. This makes computations faster and enforces some constraints.

This class is never instantiated directly by the user, but instead by the solve method of continuation.AnalyticContinuationProblem.

Create a Maxent solver object.

Here we pass the data to the solver and do some precomputations. This makes the process of solving much faster.

Parameters
  • im_axis (numpy.ndarray) – One-dimensional numpy array of type float Matsubara frequencies or imaginary time points of the input data.

  • re_axis (numpy.ndarray) – One-dimensional numpy array of type float Real-frequency grid to use for the analytic continuation.

  • im_data (numpy.ndarray) – One-dimensional numpy array of type float or complex Imaginary-axis data, e.g. Matsubara Green’s function

  • kernel_mode ({‘freq_fermionic’, ‘freq_bosonic’, ‘time_fermionic’, ‘time_bosonic’}) –

    • ‘freq_fermionic’ fermionic Matsubara Greens function

    • ’freq_bosonic’ bosonic Matsubara Greens function

    • ’time_fermionic’ fermionic Green’s function in imaginary time

    • ’time_bosonic’ bosonic Green’s function (susceptibility) in imaginary time

    Additionally there are less established special-purpose options [‘freq_bosonic_xyz’, ‘freq_fermionic_phsym’, ‘time_fermionic_phsym’].

  • model (numpy.ndarray) – One-dimensional numpy array of type float. (Values must be greater or equal to zero.) Default model of the Maxent calculation, i.e. the entropy of the spectral function is computed with respect to this model. Thus the shape has to match re_axis.

  • stdev (numpy.ndarray) – One-dimensional numpy array of positive float values. stdev are the standard deviation values of the measured data points, thus they have to be larger than zero. The shape has to match im_data The keywords stdev and cov are mutually exclusive, i.e. only one of them can be used. For complex data (e.g. Matsubara Green’s function) the same standard deviation is used for both the real and imaginary part, such that the array passed to this keyword is taken as the standard deviation of the real part (not the absolute).

  • cov (numpy.ndarray) – Two-dimensional numpy array. cov is the covariance matrix of the input data. If it is diagonal, the diagonal elements are the variance of the data, i.e. the square of the standard deviation. The keywords stdev and cov are mutually exclusive, i.e. only one of them can be used. For complex data, the covariance matrix can be complex (see Kappl et al., PRB 102 (8), 085124)

  • offdiag (bool, default: False) – Whether the input data are offdiagonal elements of a matrix. For diagonal elements, the spectral function is positive, for offdiagonal ones it has positive and negative values, and the integral is 0.

  • preblur (bool, default: False) – Whether to use preblur in the calculation

  • blur_width (float, default: 0.) – If preblur=True, specify a width of the preblur that is larger than 0.

  • optimizer ({'scipy_lm', 'newton'}, default: 'newton') – Which optimizer to use at each value of alpha to solve the minimization / root finding problem. Usually the custom implementation of the newton algorithm is much faster.

  • verbose (bool, default: True) – If set to False, no output is printed. This can be useful when doing a very large number of computations, but generally I recommend leaving it on True, to see what is happening.

backtransform(A)

Backtransformation from real to imaginary axis.

\(G(i\omega_n) = \int d\nu \; K(i\omega_n, \nu) A(\nu)\)

Also at this place we return from the ‘diagonal-covariance space’ Note: this function is not a bottleneck.

Parameters

A (numpy.ndarray) – Spectral function

Returns

Back-transformed Green’s function on imaginary axis

Return type

np.ndarray

bayes_conv(A, entr, alpha)

Bayesian convergence criterion for classic maxent (maximum of probablility distribution)

Parameters
  • A (numpy.ndarray) – spectral function

  • entr (float) – entropy

  • alpha (float) – weight factor of the entropy

Returns

  • ng (float) – “number of good data points”

  • tr (float) – trace of the gamma matrix

  • conv (float) – convergence criterion (1 -> converged)

bayes_conv_offdiag(A, entr, alpha)

Bayesian convergence criterion for classic maxent, offdiagonal version

Parameters
  • A (numpy.ndarray) – spectral function

  • entr (float) – entropy

  • alpha (float) – weight factor of the entropy

Returns

  • ng (float) – “number of good data points”

  • tr (float) – trace of the gamma matrix

  • conv (float) – convergence criterion (1 -> converged)

chi2(A)

compute chi-squared-deviation

Compute the log-likelihood function or chi-squared-deviation of the spectral function: \(\sum_n \frac{|G(i\omega_n) - \int K(i\omega_n, \nu) A(\nu)|^2}{\sigma_n^2}\)

Parameters

A (numpy.ndarray) – Spectral function

Returns

chi-squared deviation

Return type

float

compute_f_J_diag(u, alpha)

This function evaluates the function whose root we want to find.

The function \(f_m(u)\) is defined as the singular value decomposition of the derivative \(dQ[A]/dA_m\). Since we want to minimize \(Q[A]\), we have to find the root of the vector-valued function \(f\), i.e. \(f_m(u) = SVD(dQ/dA)_m = 0\). For more efficient root finding, we also need the Jacobian \(J\). It is directly computed in singular space, \(J_{mi}=df_m/du_i\).

Parameters
  • u (numpy.ndarray) – singular-space vector that parametrizes the spectral function

  • alpha (float) – (positive) weight factor of the entropy

Returns

  • f (numpy.ndarray) – value of the function whose zero we want to find

  • J (numpy.ndarray) – Jacobian at the current position

compute_f_J_offdiag(u, alpha)

The analogue to compute_f_J_diag for offdiagonal elements.

Parameters
  • u (numpy.ndarray) – singular-space vector that parametrizes the spectral function

  • alpha (float) – (positive) weight factor of the entropy

Returns

  • f (numpy.ndarray) – value of the function whose zero we want to find

  • J (numpy.ndarray) – Jacobian at the current position

entropy_pos(A, u)

Compute entropy for positive definite spectral function.

Parameters
  • A (numpy.ndarray) – Spectral function

  • u (numpy.ndarray) – Singular-space vector representing the same spectral function

Returns

entropy

Return type

float

entropy_posneg(A, u)

Compute “positive-negative entropy” for spectral function with norm 0.

Parameters
  • A (numpy.ndarray) – Spectral function

  • u (numpy.ndarray) – Singular-space vector representing the same spectral function

Returns

entropy

Return type

float

error_propagation(obs, args)

Calculate the deviation for a callable function obs. This feature is still experimental

log(msg)

Print log messages if verbose=True was passed to the solver.

Parameters

msg (str) – Message to print

maxent_optimization(alpha, ustart, iterfac=10000000, use_bayes=False, **kwargs)

optimization of maxent functional for a given value of alpha

Since a priori the best value of \(\alpha\) is unknown, this function has to be called several times in order to find a good value.

Parameters
  • alpha (float) – weight factor of the entropy

  • ustart (numpy.ndarray) – singular-space vector used as a starting value for the optimization. For the first optimization, done at large alpha, we use zeros, which corresponds to the default model. Then we use the result of the previous optimization as a starting value.

  • iterfac (int, default: 10000000) – Control parameter for maximum number of iterations in scipy_lm, which is <number of singular values> * iterfac. It has no effect when using the newton optimizer.

  • use_bayes (bool, default=False) – Whether to use the Bayesian inference parameters for alpha.

Returns

Object that holds the results of the optimization, e.g. spectral function, chi-squared deviation.

Return type

OptimizationResult

posterior_probability(A, alpha, entr, chisq)

Bayesian a-posteriori probability for alpha after optimization of A

Parameters
  • A (numpy.ndarray) – spectral function

  • entr (float) – entropy

  • alpha (float) – weight factor of the entropy

  • chisq (float) – chi-squared deviation

Returns

Probability

Return type

float

singular_to_realspace_diag(u)

Go from singular to real space.

Transform the singular space vector u into real-frequency space (spectral function) by \(A(\omega) = D(\omega) e^{V u}\), where \(D\) is the default model and \(V\) is the matrix from the SVD.

Parameters

u (numpy.ndarray) – singular-space vector that parametrizes the spectral function

Returns

Spectral function \(A\) obtained from u

Return type

numpy.ndarray

singular_to_realspace_offdiag(u)

Go from singular to real space.

Transform the singular space vector u into real-frequency space in the case of an offdiagonal element.

Parameters

u (numpy.ndarray) – singular-space vector that parametrizes the spectral function

Returns

Spectral function \(A\) obtained from u

Return type

numpy.ndarray

solve(**kwargs)

Wrapper function for solve, which calls the chosen method of alpha_determination.

solve_bryan(alphastart=500, alphadiv=1.1, interactive=False)

Bryan’s method of determining the optimal spectrum.

Bryan’s maxent calculates an average of spectral functions, weighted by their Bayesian probability

Parameters
  • alphastart (float, default=500) – Starting value of alpha. This is the largest value of alpha, it is decreased during the calculation.

  • alphadiv (float, default=1.1) – After each optimization, the current alpha is divided by this number. Hence, the number has to be larger than 1.

  • interactive (bool, default=False) – Whether to show a plot of the probability. (Needs matplotlib)

Returns

  • OptimizationResult – Contains the weighted average of all results

  • list – List of OptimizationResult objects for all used values of alpha

solve_chi2kink(alpha_start=1000000000.0, alpha_end=0.001, alpha_div=10.0, fit_position=2.5, interactive=False, **kwargs)

Determine optimal alpha by searching for the kink in log(chi2) (log(alpha))

We start with an optimization at a large value of alpha (alpha_start), where we should get only the default model. Then, alpha is decreased step-by-step, where \(\alpha_{n+1} = \alpha_n / \alpha_{div}\), until the minimal value of alpha_end is reached. Then, we fit a function \(\phi(x; a, b, c, d) = a + b / [1 + exp(-d*(x-c))]\), from which the optimal alpha is determined by x_opt = c - fit_position / d; alpha_opt = 10^x_opt.

Parameters
  • alpha_start (float, default=1e9) – Value of alpha where to start the procedure. This is the largest value of alpha, it is subsequently decreased in the algorithm.

  • alpha_end (float, default=1e-3) – Last (smallest) value of alpha that is considered in the algorithm

  • alpha_div (float, default=10.) – After each optimization, alpha is divided by alpha_div. Thus it has to be larger than 1. The default value of 10 is a good compromise of function and speed. You can take smaller values if you are unsure or want to make fancy plots.

  • fit_position (float, default=2.5) – Control parameter for under/overfitting. In my experience, good values are usually between 2 and 2.5. Smaller values lead to underfitting, which is sometimes desirable. Larger values lead to overfitting, which should be avoided.

  • interactive (bool, default=False) – Decide whether to show a plot of chi2 vs alpha.

Returns

  • OptimizationResult – Result of the optimization at “best” alpha value

  • list – List of OptimizationResult objects for all used values of alpha

solve_classic()

Classic maxent alpha determination.

Classic maxent uses Bayes statistics to approximately determine the most probable value of alpha We start at a large value of alpha, where the optimization yields basically the default model, therefore u_opt is only a few steps away from ustart=0 (=default model) Then we gradually decrease alpha, step by step moving away from the default model towards data fitting. Using u_opt as ustart for the next (smaller) alpha brings a great speedup into this procedure.

Returns

  • OptimizationResult – Result of the optimization at “best” alpha value

  • list – List of OptimizationResult objects for all used values of alpha

solve_historic()

Historic Maxent: choose alpha in a way that chi^2 pprox N

Returns

  • OptimizationResult – Result of the optimization at “best” alpha value

  • list – List of OptimizationResult objects for all used values of alpha

class ana_cont.solvers.NewtonOptimizer(opt_size, max_hist=1, max_iter=20000, initial_guess=None)

Newton root finding.

Parameters
  • opt_size – number of variables (integer)

  • max_hist – maximal history for mixing (integer)

  • max_iter – maximum number of iterations for root finding

  • initial_guess – initial guess for the root finding.

__call__(function_and_jacobian, alpha)

Main function of Newton optimization.

This function implements the self-consistent iteration of the root finding.

Returns

sol.x is the result sol.nfev is the number of function evaluations

Return type

collections.namedtuple

get_proposal(mixing=0.5)

Propose a new solution by linear mixing.

iteration_function(proposal, function_vector, jacobian_matrix)

The function, whose fixed point we are searching.

This function generates the iteration procedure in the Newton root finding. Basically, it computes the “result” from the “proposal”. It has the form result = proposal + increment, but since the increment may be large, we apply a reduction of step width in such cases.

class ana_cont.solvers.OptimizationResult(u_opt=None, A_opt=None, chi2=None, backtransform=None, entropy=None, n_good=None, probability=None, alpha=None, convergence=None, trace=None, Q=None, norm=None, blur_width=None, numerator=None, denominator=None, numerator_function=None, denominator_function=None, check=None, ivcheck=None, g_ret=None)

Object for holding the result of an optimization.

This class has no methods except the constructor, it is thus essentially a collection of output numbers.

All member variables have None as default value, different solvers override different variables. A_opt is always set, since it is the main result of analytic continuation.

class ana_cont.solvers.PadeSolver(im_axis, re_axis, im_data)

Pade solver

Parameters
  • im_axis (numpy.ndarray) – Matsubara frequencies which are used for the continuation

  • re_axis (numpy.ndarray) – Real-frequency points at which the Pade interpolant is evaluated

  • im_data (Green’s function values at the given points im_axis) –

check(im_axis_fine=None)

Sanity check for Pade approximant

Evaluate the Pade approximant on the imaginary axis, however not only at Matsubara frequencies, but on a dense grid. If the approximant is good, then this should yield a smooth interpolating curve on the Matsubara axis. On the other hand, little ‘spikes’ are a hint for pole-zero pairs close to the imaginary axis. This is usually a result of noise and a different choice of Matsubara frequencies may help.

Parameters

im_axis_fine (numpy.ndarray, default=None) – Imaginary-axis points where the approximant is evaluated for checking. If not specified, an array of length 500 is generated, which goes up to twice the maximum frequency that was used for constructing the approximant.

Returns

Values of the Pade approximant at the points of im_axis_fine.

Return type

numpy.ndarray

solve(show_plot=False)

Compute the Pade approximation on the real axis.

The main numerically heavy computation is done in the Cython module pade.pyx. Here we just call the functions. In the Pade method, the numerator and denominator approximants are generated separately, and then the division is done. As an additional feature we add the callable numerator_function and denominator_function to the OptimizationResult object.

11.3. Kernels module

class ana_cont.kernels.Kernel(kind=None, re_axis=None, im_axis=None)

This class handles the kernel of the analytic continuation.

The kernel matrix is calculated from the real-frequency vector re_axis and the imaginary frequency or time vector im_axis. It is important to note that the Kernel does not know the inverse temperature beta. In case of the frequency kernels (‘freq_fermionic’ or ‘freq_bosonic’) the temperature is implicitly contained through the Matsubara frequencies. In case of the time kernels (‘time_fermionic’ or ‘time_bosonic’) the im_axis list has to be rescaled such that the temperature would be 1. Usually this is done by imag_time = imag_time / beta.

For regular users of the library this is no problem, since the Kernel class is only instantiated inside the AnalyticContinuationProblem class, where the rescaling is done automatically. (Also the results are scaled back automatically there.)

Parameters
  • kind (str, default=None) – Which kind of kernel to use. Possible options: ‘freq_bosonic’, ‘time_bosonic’, ‘freq_bosonic_xyz’, ‘freq_fermionic’, ‘time_fermionic’, ‘freq_fermionic_phsym’, ‘time_fermionic_phsym’

  • re_axis (numpy.ndarray, default=None) – real-frequency axis to generate the kernel matrix

  • im_axis (numpy.ndarray, default=None) – imaginary axis (time, Matsubara frequency) to generate the kernel matrix

blur(hidden_spectrum)

Convert hidden spectral function to spectral function.

convolve_kernel()

Convolve bosonic or fermionic kernel with a Gaussian.

The Gaussian is \(g(x) = \frac{1}{b \sqrt{2 \pi}} exp(-x^2/(2b^2))\). In the fermionic case, the convolution can be written as \(K_{preblur}(i\nu_n, \omega) = \int_{-5b}^{5b} dx\; \frac{g(x)}{i\nu_n - x - \omega}\)

In the bosonic case, the convolution can be written as \(K_{preblur}(i\omega_n, \nu) = \frac{1}{2} \int_{-5b}^{5b} dx\; g(x) [\frac{(x+\nu)^2 }{ (x+\nu)^2 + \omega_n^2} + \frac{(x-\nu)^2 }{ (x-\nu)^2 + \omega_n^2}]\)

Integration over the Gaussian from \(-5b\) to \(5b\) is certainly sufficient. Thus the Gaussian has to be computed only once and integration by scipy.integrate.simps gives very accurate results even for tiny values of b.

Note: More direct calculation by, e.g., scipy.integrate.quad is possible but unstable and extremely slow.

kernel_matrix()

Compute the kernel matrix. If you want to implement another kernel, you just have to add another ‘elif’ here.

preblur(blur_width)

Preblur for the kernel, if applicable.

real_matrix()

Return real and imaginary part one after another.

rotate_to_cov_eb(ucov)

Rotate first axis of kernel to eigenbasis of covariance matrix.

11.4. Graphical user interface

class gui.gui_backend.InputData(fname=None, iter_type=None, iter_num=None, data_type=None, atom=None, orbital=None, spin=None, num_mats=None, ignore_real_part=False)

Input data for the analytic continuation of a w2dynamics result.

Initialize the input data object.

fname – Name (str) of a valid w2dynamics DMFT output file.

A w2dynamics version from 2020 or newer should be used.

iter_type – iteration type: ‘dmft’ or ‘stat’ iter_num – iteration number: integer or ‘last’ (-1 also points to last) data_type – “Green’s function” or “Self-energy” atom – which inequivalent atom (one-based integer), e.g. 1 leads to ‘ineq-001’ orbital – which orbital component to load (one-based integer) spin – which spin projection to load: ‘up’, ‘down’, ‘average’ num_mats – number of Matsubara frequencies for continuation (integer) ignore_real_part – if True, the real part is set to zero.

generate_mats_freq()

Generate the Matsubara frequency grid.

get_iteration()

Compose the whole iteration string, e.g. ‘dmft-003’.

load_data()

Load Matsubara-frequency data.

Load w2dynamics data on Matsubara frequency axis. This is either a self-energy, or a Green’s function. For the self-energy, there are two different possible formats.

load_giw()

Load Matsubara Green’s function data.

Here we have only one possible format, i.e. (orb, spin, orb, spin, freq). If no QMC error is present, we set it to a constant value err_const.

load_siw_1()

Load self-energy with jackknife error, type 1.

In the first implementation, the self-energy with jackknife error was written in the group ‘siw-full-jkerr’, where the frequency axis is the first, followed by orbital and spin: (freq, orb, spin, orb, spin) The Hartree term is subtracted automatically.

load_siw_2()

Load self-energy with jackknife error, type 2.

In newer versions, the self-energy with jackknife error is stored in siw-full, in standard format (orb, spin, orb, spin, freq). If the QMC error is not present, it is set to a constant value err_const. The Hartree term is subtracted automatically.

class gui.gui_backend.OutputData

Output data of the analytic continuation.

interpolate(original_grid, original_function, n_reg)

Spline interpolation of real-frequency data.

save(interpolate=False, n_reg=None)

Save output data to text file.

Here we have to distinguish several different cases. * For a self-energy, we first do a Kramers-Kronig transformation

to get also the real part. Then the Hartree term is added. Remember: when reading input data from a text file, the Hartree term has to be handeled manually, the program treats it as zero. The whole function, i.e. real and imaginary part of the self-energy, is saved.

  • In case of a Green’s function, we save the spectral function,

    i.e. -Im[G(omega)] / pi. The real part is not computed.

  • If desired, an interpolation to a regular-spaced grid is done.

Output format: text file. For spectral function: two colums: frequency, spectrum For self-energy: three columns: frequency, real part, imaginary part.

update(w, spec, input_data)

Update the output data with the latest input and results.

class gui.gui_backend.RealFrequencyGrid(wmax=None, nw=None, type=None)

Class for real-frequency grids.

Our real-frequency grids are always symmetric around 0. Thus they cover an interval [-wmax, wmax], where the endpoint is included. If the number of grid-points nw is an odd number, zero is included. We can create two types of grids: equispaced or non-equispaced (centered). The latter is more dense in the low-frequency region, and more sparse in the high-frequency region.

Initialize the real-frequency grid.

wmax – border of the real-frequency interval [-wmax, wmax] nw – total number of real frequency grid points type – grid type: ‘equispaced grid’ or ‘centered grid’

create_grid()

Create the real-frequency grid.

There are four possible types of real-frequency grids. In each case, the endpoint is included in the grid.

  • An ‘equispaced symmetric’ grid is simply a linspace, containing also the endpoint.

  • The ‘centered symmetric’ grid is created by a tangent function, such that the grid points are denser around zero.

  • ‘equispaced positive’ is a simple linspace, starting from zero and containing only positive values.

  • ‘centered positive’ also starts from zero and contains only positive values, but close to zero they are lying denser.

class gui.gui_backend.TextInputData(fname=None, data_type=None, num_mats=None, n_skip=None)

Input data for analytic continuation, from a generic text file.

Initialize the text file input.

fname – file name of the text file data_type – either “Green’s function” or “Self-energy” num_mats – number of Matsubara frequencies n_skip – number of lines to skip at the beginning of the text file.

Note the following points: * The Hartree term has to be subtracted beforehand,

i.e. both the real and imaginary part of the data in the input file have to approach zero in the high-frequency limit,

  • The input file must contain only data on the positive half of the Matsubara axis,

  • data_type will not have any impact on the analytic continuation itself,

    but if it is a self-energy, the real part is calculated by Kramers-Kronig after the continuation, and the full self-energy is stored.

  • The Hartree term has to be handeled completely manually in pre- and postprocessing.

read_data()

Read the text file by np.loadtxt.

gui.gui_backend.input_data_plot(mats, value, error, datatype, mom1=None, hartree=None)

Generate a very simple plot of the Matsubara data.