<li>automatic differentiation for exact linear error propagation as suggested in <ahref="https://arxiv.org/abs/1809.01289">arXiv:1809.01289</a> (partly based on the <ahref="https://github.com/HIPS/autograd">autograd</a> package).</li>
<li>treatment of slow modes in the simulation as suggested in <ahref="https://arxiv.org/abs/1009.5228">arXiv:1009.5228</a>.</li>
<li>coherent error propagation for data from different Markov chains.</li>
<li>non-linear fits with x- and y-errors and exact linear error propagation based on automatic differentiation as introduced in <ahref="https://arxiv.org/abs/1809.01289">arXiv:1809.01289</a>.</li>
<li>real and complex matrix operations and their error propagation based on automatic differentiation (Matrix inverse, Cholesky decomposition, calculation of eigenvalues and eigenvectors, singular value decomposition...).</li>
<p>More detailed examples can found in the <ahref="https://github.com/fjosw/pyerrors/tree/develop/examples">GitHub repository</a><ahref="https://mybinder.org/v2/gh/fjosw/pyerrors/HEAD?labpath=examples"><imgsrc="https://img.shields.io/badge/-try%20it%20out-579ACA.svg?logo=data:image/png;base64,iVBORw0KGgoAAAANSUhEUgAAAFkAAABZCAMAAABi1XidAAAB8lBMVEX///9XmsrmZYH1olJXmsr1olJXmsrmZYH1olJXmsr1olJXmsrmZYH1olL1olJXmsr1olJXmsrmZYH1olL1olJXmsrmZYH1olJXmsr1olL1olJXmsrmZYH1olL1olJXmsrmZYH1olL1olL0nFf1olJXmsrmZYH1olJXmsq8dZb1olJXmsrmZYH1olJXmspXmspXmsr1olL1olJXmsrmZYH1olJXmsr1olL1olJXmsrmZYH1olL1olLeaIVXmsrmZYH1olL1olL1olJXmsrmZYH1olLna31Xmsr1olJXmsr1olJXmsrmZYH1olLqoVr1olJXmsr1olJXmsrmZYH1olL1olKkfaPobXvviGabgadXmsqThKuofKHmZ4Dobnr1olJXmsr1olJXmspXmsr1olJXmsrfZ4TuhWn1olL1olJXmsqBi7X1olJXmspZmslbmMhbmsdemsVfl8ZgmsNim8Jpk8F0m7R4m7F5nLB6jbh7jbiDirOEibOGnKaMhq+PnaCVg6qWg6qegKaff6WhnpKofKGtnomxeZy3noG6dZi+n3vCcpPDcpPGn3bLb4/Mb47UbIrVa4rYoGjdaIbeaIXhoWHmZYHobXvpcHjqdHXreHLroVrsfG/uhGnuh2bwj2Hxk17yl1vzmljzm1j0nlX1olL3AJXWAAAAbXRSTlMAEBAQHx8gICAuLjAwMDw9PUBAQEpQUFBXV1hgYGBkcHBwcXl8gICAgoiIkJCQlJicnJ2goKCmqK+wsLC4usDAwMjP0NDQ1NbW3Nzg4ODi5+3v8PDw8/T09PX29vb39/f5+fr7+/z8/Pz9/v7+zczCxgAABC5JREFUeAHN1ul3k0UUBvCb1CTVpmpaitAGSLSpSuKCLWpbTKNJFGlcSMAFF63iUmRccNG6gLbuxkXU66JAUef/9LSpmXnyLr3T5AO/rzl5zj137p136BISy44fKJXuGN/d19PUfYeO67Znqtf2KH33Id1psXoFdW30sPZ1sMvs2D060AHqws4FHeJojLZqnw53cmfvg+XR8mC0OEjuxrXEkX5ydeVJLVIlV0e10PXk5k7dYeHu7Cj1j+49uKg7uLU61tGLw1lq27ugQYlclHC4bgv7VQ+TAyj5Zc/UjsPvs1sd5cWryWObtvWT2EPa4rtnWW3JkpjggEpbOsPr7F7EyNewtpBIslA7p43HCsnwooXTEc3UmPmCNn5lrqTJxy6nRmcavGZVt/3Da2pD5NHvsOHJCrdc1G2r3DITpU7yic7w/7Rxnjc0kt5GC4djiv2Sz3Fb2iEZg41/ddsFDoyuYrIkmFehz0HR2thPgQqMyQYb2OtB0WxsZ3BeG3+wpRb1vzl2UYBog8FfGhttFKjtAclnZYrRo9ryG9uG/FZQU4AEg8ZE9LjGMzTmqKXPLnlWVnIlQQTvxJf8ip7VgjZjyVPrjw1te5otM7RmP7xm+sK2Gv9I8Gi++BRbEkR9EBw8zRUcKxwp73xkaLiqQb+kGduJTNHG72zcW9LoJgqQxpP3/Tj//c3yB0tqzaml05/+orHLksVO+95kX7/7qgJvnjlrfr2Ggsyx0eoy9uPzN5SPd86aXggOsEKW2Prz7du3VID3/tzs/sSRs2w7ovVHKtjrX2pd7ZMlTxAYfBAL9jiDwfLkq55Tm7ifhMlTGPyCAs7RFRhn47JnlcB9RM5T97ASuZXIcVNuUDIndpDbdsfrqsOppeXl5Y+XVKdjFCTh+zGaVuj0d9zy05PPK3QzBamxdwtTCrzyg/2Rvf2EstUjordGwa/kx9mSJLr8mLLtCW8HHGJc2R5hS219IiF6PnTusOqcMl57gm0Z8kanKMAQg0qSyuZfn7zItsbGyO9QlnxY0eCuD1XL2ys/MsrQhltE7Ug0uFOzufJFE2PxBo/YAx8XPPdDwWN0MrDRYIZF0mSMKCNHgaIVFoBbNoLJ7tEQDKxGF0kcLQimojCZopv0OkNOyWCCg9XMVAi7ARJzQdM2QUh0gmBozjc3Skg6dSBRqDGYSUOu66Zg+I2fNZs/M3/f/Grl/XnyF1Gw3VKCez0PN5IUfFLqvgUN4C0qNqYs5YhPL+aVZYDE4IpUk57oSFnJm4FyCqqOE0jhY2SMyLFoo56zyo6becOS5UVDdj7Vih0zp+tcMhwRpBeLyqtIjlJKAIZSbI8SGSF3k0pA3mR5tHuwPFoa7N7reoq2bqCsAk1HqCu5uvI1n6JuRXI+S1Mco54YmYTwcn6Aeic+kssXi8XpXC4V3t7/ADuTNKaQJdScAAAAAElFTkSuQmCC"alt="badge"/></a>.</p>
<p>If you use <code><ahref="">pyerrors</a></code> for research that leads to a publication please consider citing:
Fabian Joswig, Simon Kuberski, Justus T. Kuhlmann, Jan Neuendorf, <em>pyerrors: a python framework for error analysis of Monte Carlo data</em>. [arXiv:2209.14371 [hep-lat]].</p>
<li>Ulli Wolff, <em>Monte Carlo errors with less errors</em>. Comput.Phys.Commun. 156 (2004) 143-153, Comput.Phys.Commun. 176 (2007) 383 (erratum).</li>
<li>Alberto Ramos, <em>Automatic differentiation for error analysis of Monte Carlo data</em>. Comput.Phys.Commun. 238 (2019) 19-35.</li>
</ul>
<p>and</p>
<ul>
<li>Stefan Schaefer, Rainer Sommer, Francesco Virotta, <em>Critical slowing down and error analysis in lattice QCD simulations</em>. Nucl.Phys.B 845 (2011) 93-119.</li>
<p>There exist similar publicly available implementations of gamma method error analysis suites in <ahref="https://gitlab.ift.uam-csic.es/alberto/aderrors">Fortran</a>, <ahref="https://gitlab.ift.uam-csic.es/alberto/aderrors.jl">Julia</a> and <ahref="https://github.com/mbruno46/pyobs">Python</a>.</p>
<spanclass="n">my_obs</span><spanclass="o">=</span><spanclass="n">pe</span><spanclass="o">.</span><spanclass="n">Obs</span><spanclass="p">([</span><spanclass="n">samples</span><spanclass="p">],</span><spanclass="p">[</span><spanclass="s1">'ensemble_name'</span><spanclass="p">])</span><spanclass="c1"># Initialize an Obs object</span>
<spanclass="n">my_new_obs</span><spanclass="o">.</span><spanclass="n">gamma_method</span><spanclass="p">()</span><spanclass="c1"># Estimate the statistical error</span>
<spanclass="nb">print</span><spanclass="p">(</span><spanclass="n">my_new_obs</span><spanclass="p">)</span><spanclass="c1"># Print the result to stdout</span>
<p><code><ahref="">pyerrors</a></code> introduces a new datatype, <code>Obs</code>, which simplifies error propagation and estimation for auto- and cross-correlated data.
An <code>Obs</code> object can be initialized with two arguments, the first is a list containing the samples for an observable from a Monte Carlo chain.
The samples can either be provided as python list or as numpy array.
The second argument is a list containing the names of the respective Monte Carlo chains as strings. These strings uniquely identify a Monte Carlo chain/ensemble.</p>
<p>When performing mathematical operations on <code>Obs</code> objects the correct error propagation is intrinsically taken care of using a first order Taylor expansion
as introduced in <ahref="https://arxiv.org/abs/hep-lat/0306017">arXiv:hep-lat/0306017</a>.
The required derivatives $\bar{f}_\alpha$ are evaluated up to machine precision via automatic differentiation as suggested in <ahref="https://arxiv.org/abs/1809.01289">arXiv:1809.01289</a>.</p>
<p>The error estimation within <code><ahref="">pyerrors</a></code> is based on the gamma method introduced in <ahref="https://arxiv.org/abs/hep-lat/0306017">arXiv:hep-lat/0306017</a>.
<p>We use the following definition of the integrated autocorrelation time established in <ahref="https://link.springer.com/article/10.1007/BF01022990">Madras & Sokal 1988</a>
The window $W$ is determined via the automatic windowing procedure described in <ahref="https://arxiv.org/abs/hep-lat/0306017">arXiv:hep-lat/0306017</a>.
The standard value for the parameter $S$ of this automatic windowing procedure is $S=2$. Other values for $S$ can be passed to the <code>gamma_method</code> as parameter.</p>
<p>The integrated autocorrelation time $\tau_\mathrm{int}$ and the autocorrelation function $\rho(W)$ can be monitored via the methods <code><ahref="pyerrors/obs.html#Obs.plot_tauint">pyerrors.obs.Obs.plot_tauint</a></code> and <code><ahref="pyerrors/obs.html#Obs.plot_rho">pyerrors.obs.Obs.plot_rho</a></code>.</p>
<p>Slow modes in the Monte Carlo history can be accounted for by attaching an exponential tail to the autocorrelation function $\rho$ as suggested in <ahref="https://arxiv.org/abs/1009.5228">arXiv:1009.5228</a>. The longest autocorrelation time in the history, $\tau_\mathrm{exp}$, can be passed to the <code>gamma_method</code> as parameter. In this case the automatic windowing procedure is vacated and the parameter $S$ does not affect the error estimate.</p>
<p>Error propagation for multiple ensembles (Markov chains with different simulation parameters) is handled automatically. Ensembles are uniquely identified by their <code>name</code>.</p>
<p><code><ahref="">pyerrors</a></code> identifies multiple replica (independent Markov chains with identical simulation parameters) by the vertical bar <code>|</code> in the name of the data set.</p>
<p>In order to keep track of different error analysis parameters for different ensembles one can make use of global dictionaries as detailed in the following example.</p>
<p>In case the <code>gamma_method</code> is called without any parameters it will use the values specified in the dictionaries for the respective ensembles.
Passing arguments to the <code>gamma_method</code> still dominates over the dictionaries.</p>
<p><code>Obs</code> objects defined on regular and irregular histories of the same ensemble can be combined with each other and the correct error propagation and estimation is automatically taken care of.</p>
Make sure to check the autocorrelation time with e.g. <code><ahref="pyerrors/obs.html#Obs.plot_rho">pyerrors.obs.Obs.plot_rho</a></code> or <code><ahref="pyerrors/obs.html#Obs.plot_tauint">pyerrors.obs.Obs.plot_tauint</a></code>.</p>
<p>When one is not interested in single observables but correlation functions, <code><ahref="">pyerrors</a></code> offers the <code>Corr</code> class which simplifies the corresponding error propagation and provides the user with a set of standard methods. In order to initialize a <code>Corr</code> objects one needs to arrange the data as a list of <code>Obs</code></p>
<p>In case the correlation functions are not defined on the outermost timeslices, for example because of fixed boundary conditions, a padding can be introduced.</p>
<p>Error propagation with the <code>Corr</code> class works very similar to <code>Obs</code> objects. Mathematical operations are overloaded and <code>Corr</code> objects can be computed together with other <code>Corr</code> objects, <code>Obs</code> objects or real numbers and integers.</p>
<li><code>Corr.gamma_method</code> applies the gamma method to all entries of the correlator.</li>
<li><code>Corr.m_eff</code> to construct effective masses. Various variants for periodic and fixed temporal boundary conditions are available.</li>
<li><code>Corr.deriv</code> returns the first derivative of the correlator as <code>Corr</code>. Different discretizations of the numerical derivative are available.</li>
<li><code>Corr.second_deriv</code> returns the second derivative of the correlator as <code>Corr</code>. Different discretizations of the numerical derivative are available.</li>
<li><code>Corr.symmetric</code> symmetrizes parity even correlations functions, assuming periodic boundary conditions.</li>
<li><code>Corr.T_symmetry</code> averages a correlator with its time symmetry partner, assuming fixed boundary conditions.</li>
<li><code>Corr.plateau</code> extracts a plateau value from the correlator in a given range.</li>
<li><code>Corr.roll</code> periodically shifts the correlator.</li>
<li><code>Corr.reverse</code> reverses the time ordering of the correlator.</li>
<li><code>Corr.correlate</code> constructs a disconnected correlation function from the correlator and another <code>Corr</code> or <code>Obs</code> object.</li>
<li><code>Corr.reweight</code> reweights the correlator.</li>
<p><code><ahref="">pyerrors</a></code> can also handle matrices of correlation functions and extract energy states from these matrices via a generalized eigenvalue problem (see <code><ahref="pyerrors/correlators.html#Corr.GEVP">pyerrors.correlators.Corr.GEVP</a></code>).</p>
<p><code><ahref="">pyerrors</a></code> can handle complex valued observables via the class <code><ahref="pyerrors/obs.html#CObs">pyerrors.obs.CObs</a></code>.
<p>In many projects, auxiliary data that is not based on Monte Carlo chains enters. Examples are experimentally determined mesons masses which are used to set the scale or renormalization constants. These numbers come with an error that has to be propagated through the analysis. The <code>Covobs</code> class allows to define such quantities in <code><ahref="">pyerrors</a></code>. Furthermore, external input might consist of correlated quantities. An example are the parameters of an interpolation formula, which are defined via mean values and a covariance matrix between all parameters. The contribution of the interpolation formula to the error of a derived quantity therefore might depend on the complete covariance matrix.</p>
<p>This concept is built into the definition of <code>Covobs</code>. In <code><ahref="">pyerrors</a></code>, external input is defined by $M$ mean values, a $M\times M$ covariance matrix, where $M=1$ is permissible, and a name that uniquely identifies the covariance matrix. Below, we define the pion mass, based on its mean value and error, 134.9768(5). Note, that the square of the error enters <code>cov_Obs</code>, since the second argument of this function is the covariance matrix of the <code>Covobs</code>.</p>
<p>The resulting object <code>mpi</code> is an <code>Obs</code> that contains a <code>Covobs</code>. In the following, it may be handled as any other <code>Obs</code>. The contribution of the covariance matrix to the error of an <code>Obs</code> is determined from the $M \times M$ covariance matrix $\Sigma$ and the gradient of the <code>Obs</code> with respect to the external quantities, which is the $1\times M$ Jacobian matrix $J$, via
<p>Since the gradient of a derived observable with respect to an external covariance matrix is propagated through the entire analysis, the <code>Covobs</code> class allows to quote the derivative of a result with respect to the external quantities. If these derivatives are published together with the result, small shifts in the definition of external quantities, e.g., the definition of the physical point, can be performed a posteriori based on the published information. This may help to compare results of different groups. The gradient of an <code>Obs</code><code>o</code> with respect to a covariance matrix with the identifying string <code>k</code> may be accessed via</p>
<p><code><ahref="">pyerrors</a></code> supports exact linear error propagation for iterative algorithms like various variants of non-linear least squares fits or root finding. The derivatives required for the error propagation are calculated as described in <ahref="https://arxiv.org/abs/1809.01289">arXiv:1809.01289</a>.</p>
<p>Standard non-linear least square fits with errors on the dependent but not the independent variables can be performed with <code><ahref="pyerrors/fits.html#least_squares">pyerrors.fits.least_squares</a></code>. As default solver the Levenberg-Marquardt algorithm implemented in <ahref="https://docs.scipy.org/doc/scipy/reference/generated/scipy.optimize.least_squares.html">scipy</a> is used.</p>
<p>Fit functions have to be of the following form</p>
<p><strong>It is important that numerical functions refer to <code>autograd.numpy</code> instead of <code>numpy</code> for the automatic differentiation in iterative algorithms to work properly.</strong></p>
<p>where x is a <code>list</code> or <code>numpy.array</code> of <code>floats</code> and y is a <code>list</code> or <code>numpy.array</code> of <code>Obs</code>.</p>
<p>Data stored in <code>Corr</code> objects can be fitted directly using the <code>Corr.fit</code> method.</p>
Details about how the required covariance matrix is estimated can be found in <code><ahref="pyerrors/obs.html#covariance">pyerrors.obs.covariance</a></code>.</p>
<p>Direct visualizations of the performed fits can be triggered via <code>resplot=True</code> or <code>qqplot=True</code>. For all available options see <code><ahref="pyerrors/fits.html#least_squares">pyerrors.fits.least_squares</a></code>.</p>
<p><code><ahref="">pyerrors</a></code> can also fit data with errors on both the dependent and independent variables using the total least squares method also referred to orthogonal distance regression as implemented in <ahref="https://docs.scipy.org/doc/scipy/reference/odr.html">scipy</a>, see <code><ahref="pyerrors/fits.html#least_squares">pyerrors.fits.least_squares</a></code>. The syntax is identical to the standard least squares case, the only difference being that <code>x</code> also has to be a <code>list</code> or <code>numpy.array</code> of <code>Obs</code>.</p>
<p>For the full API see <code><ahref="pyerrors/fits.html">pyerrors.fits</a></code> for fits and <code><ahref="pyerrors/roots.html">pyerrors.roots</a></code> for finding roots of functions.</p>
<p><code><ahref="">pyerrors</a></code> provides wrappers for <code>Obs</code>- and <code>CObs</code>-valued matrix operations based on <code>numpy.linalg</code>. The supported functions include:</p>
<p>The preferred exported file format within <code><ahref="">pyerrors</a></code> is json.gz. Files written to this format are valid JSON files that have been compressed using gzip. The structure of the content is inspired by the dobs format of the ALPHA collaboration. The aim of the format is to facilitate the storage of data in a self-contained way such that, even years after the creation of the file, it is possible to extract all necessary information:</p>
<p>This can be achieved by storing all information in one single file. The export routines of <code><ahref="">pyerrors</a></code> are written such that as much information as possible is written automatically as described in the following example</p>
<spanclass="n">pe</span><spanclass="o">.</span><spanclass="n">input</span><spanclass="o">.</span><spanclass="n">json</span><spanclass="o">.</span><spanclass="n">dump_to_json</span><spanclass="p">(</span><spanclass="n">my_obs</span><spanclass="p">,</span><spanclass="s2">"test_output_file"</span><spanclass="p">,</span><spanclass="n">description</span><spanclass="o">=</span><spanclass="s2">"This file contains a test observable"</span><spanclass="p">)</span>
<spanclass="c1"># For a single observable one can equivalently use the class method dump</span>
<spanclass="n">my_obs</span><spanclass="o">.</span><spanclass="n">dump</span><spanclass="p">(</span><spanclass="s2">"test_output_file"</span><spanclass="p">,</span><spanclass="n">description</span><spanclass="o">=</span><spanclass="s2">"This file contains a test observable"</span><spanclass="p">)</span>
<p>The format also allows to directly write out the content of <code>Corr</code> objects or lists and arrays of <code>Obs</code> objects by passing the desired data to <code><ahref="pyerrors/input/json.html#dump_to_json">pyerrors.input.json.dump_to_json</a></code>.</p>
<li><code>description</code> contains information on the content of the file. This field is not filled automatically in <code><ahref="">pyerrors</a></code>. The user is advised to provide as detailed information as possible in this field. Examples are: Input files of measurements or simulations, LaTeX formulae or references to publications to specify how the observables have been computed, details on the analysis strategy, ... This field may be any valid JSON type. Strings, arrays or objects (equivalent to dicts in python) are well suited to provide information.</li>
<p>Each entry of the array belongs to a single structure of observables. Currently, these structures can be either of <code>Obs</code>, <code>list</code>, <code>numpy.ndarray</code>, <code>Corr</code>. All <code>Obs</code> inside a structure (with dimension > 0) have to be defined on the same set of configurations. Different structures, that are represented by entries of the array <code>obsdata</code>, are treated independently. Each entry of the array <code>obsdata</code> has the following required entries:</p>
<li><code>type</code> is a string that specifies the type of the structure. This allows to parse the content to the correct form after reading the file. It is always possible to interpret the content as list of Obs.</li>
<li><code>value</code> is an array that contains the mean values of the Obs inside the structure.
The following entries are optional:</li>
<li><code>layout</code> is a string that specifies the layout of multi-dimensional structures. Examples are "2, 2" for a 2x2 dimensional matrix or "64, 4, 4" for a Corr with $T=64$ and 4x4 matrices on each time slices. "1" denotes a single Obs. Multi-dimensional structures are stored in row-major format (see below).</li>
<li><code>tag</code> is any JSON type. It contains additional information concerning the structure. The <code>tag</code> of an <code>Obs</code> in <code><ahref="">pyerrors</a></code> is written here.</li>
<li><code>cdata</code> is an array that contains the data from external quantities with an error (<code>Covobs</code> in <code><ahref="">pyerrors</a></code>). We will define it below.</li>
<p>Each entry in <code>deltas</code> corresponds to one configuration of the replica and has $1+N$ many entries. The first entry is an integer that specifies the configuration number that, together with ensemble and replica name, may be used to uniquely identify the configuration on which the data has been obtained. The following N entries specify the deltas, i.e., the deviation of the observable from the mean value on this configuration, of each <code>Obs</code> inside the structure. Multi-dimensional structures are stored in a row-major format. For primary observables, such as correlation functions, $value + delta_i$ matches the primary data obtained on the configuration.</p>
<p>The array <code>cdata</code> contains information about the contribution of auxiliary observables, represented by <code>Covobs</code> in <code><ahref="">pyerrors</a></code>, to the total error of the observables. Each entry of the array belongs to one auxiliary covariance matrix and contains:</p>
<li><code>id</code>, a string that identifies the covariance matrix</li>
<li><code>layout</code>, a string that defines the dimensions of the $M\times M$ covariance matrix (has to be "M, M" or "1").</li>
<li><code>cov</code>, an array that contains the $M\times M$ many entries of the covariance matrix, stored in row-major format.</li>
<li><code>grad</code>, an array that contains N entries, one for each <code>Obs</code> inside the structure. Each entry itself is an array, that contains the M gradients of the Nth observable with respect to the quantity that corresponds to the Mth diagonal entry of the covariance matrix.</li>
</ul>
<p>A JSON schema that may be used to verify the correctness of a file with respect to the format definition is stored in ./examples/json_schema.json. The schema is a self-descriptive format definition and contains an exemplary file.</p>
<p>Julia I/O routines for the json.gz format, compatible with <ahref="https://gitlab.ift.uam-csic.es/alberto/aderrors.jl">ADerrors.jl</a>, can be found <ahref="https://github.com/fjosw/ADjson.jl">here</a>.</p>
</span><spanid="L-2"><ahref="#L-2"><spanclass="linenos"> 2</span></a><spanclass="sd"># What is pyerrors?</span>
</span><spanid="L-3"><ahref="#L-3"><spanclass="linenos"> 3</span></a><spanclass="sd">`pyerrors` is a python package for error computation and propagation of Markov chain Monte Carlo data.</span>
</span><spanid="L-4"><ahref="#L-4"><spanclass="linenos"> 4</span></a><spanclass="sd">It is based on the gamma method [arXiv:hep-lat/0306017](https://arxiv.org/abs/hep-lat/0306017). Some of its features are:</span>
</span><spanid="L-5"><ahref="#L-5"><spanclass="linenos"> 5</span></a><spanclass="sd">- automatic differentiation for exact linear error propagation as suggested in [arXiv:1809.01289](https://arxiv.org/abs/1809.01289) (partly based on the [autograd](https://github.com/HIPS/autograd) package).</span>
</span><spanid="L-6"><ahref="#L-6"><spanclass="linenos"> 6</span></a><spanclass="sd">- treatment of slow modes in the simulation as suggested in [arXiv:1009.5228](https://arxiv.org/abs/1009.5228).</span>
</span><spanid="L-7"><ahref="#L-7"><spanclass="linenos"> 7</span></a><spanclass="sd">- coherent error propagation for data from different Markov chains.</span>
</span><spanid="L-8"><ahref="#L-8"><spanclass="linenos"> 8</span></a><spanclass="sd">- non-linear fits with x- and y-errors and exact linear error propagation based on automatic differentiation as introduced in [arXiv:1809.01289](https://arxiv.org/abs/1809.01289).</span>
</span><spanid="L-9"><ahref="#L-9"><spanclass="linenos"> 9</span></a><spanclass="sd">- real and complex matrix operations and their error propagation based on automatic differentiation (Matrix inverse, Cholesky decomposition, calculation of eigenvalues and eigenvectors, singular value decomposition...).</span>
</span><spanid="L-11"><ahref="#L-11"><spanclass="linenos"> 11</span></a><spanclass="sd">More detailed examples can found in the [GitHub repository](https://github.com/fjosw/pyerrors/tree/develop/examples) [](https://mybinder.org/v2/gh/fjosw/pyerrors/HEAD?labpath=examples).</span>
</span><spanid="L-13"><ahref="#L-13"><spanclass="linenos"> 13</span></a><spanclass="sd">If you use `pyerrors` for research that leads to a publication please consider citing:</span>
</span><spanid="L-14"><ahref="#L-14"><spanclass="linenos"> 14</span></a><spanclass="sd">Fabian Joswig, Simon Kuberski, Justus T. Kuhlmann, Jan Neuendorf, *pyerrors: a python framework for error analysis of Monte Carlo data*. [arXiv:2209.14371 [hep-lat]].</span>
</span><spanid="L-15"><ahref="#L-15"><spanclass="linenos"> 15</span></a><spanclass="sd">- Ulli Wolff, *Monte Carlo errors with less errors*. Comput.Phys.Commun. 156 (2004) 143-153, Comput.Phys.Commun. 176 (2007) 383 (erratum).</span>
</span><spanid="L-16"><ahref="#L-16"><spanclass="linenos"> 16</span></a><spanclass="sd">- Alberto Ramos, *Automatic differentiation for error analysis of Monte Carlo data*. Comput.Phys.Commun. 238 (2019) 19-35.</span>
</span><spanid="L-20"><ahref="#L-20"><spanclass="linenos"> 20</span></a><spanclass="sd">- Stefan Schaefer, Rainer Sommer, Francesco Virotta, *Critical slowing down and error analysis in lattice QCD simulations*. Nucl.Phys.B 845 (2011) 93-119.</span>
</span><spanid="L-24"><ahref="#L-24"><spanclass="linenos"> 24</span></a><spanclass="sd">There exist similar publicly available implementations of gamma method error analysis suites in [Fortran](https://gitlab.ift.uam-csic.es/alberto/aderrors), [Julia](https://gitlab.ift.uam-csic.es/alberto/aderrors.jl) and [Python](https://github.com/mbruno46/pyobs).</span>
</span><spanid="L-41"><ahref="#L-41"><spanclass="linenos"> 41</span></a><spanclass="sd">`pyerrors` introduces a new datatype, `Obs`, which simplifies error propagation and estimation for auto- and cross-correlated data.</span>
</span><spanid="L-42"><ahref="#L-42"><spanclass="linenos"> 42</span></a><spanclass="sd">An `Obs` object can be initialized with two arguments, the first is a list containing the samples for an observable from a Monte Carlo chain.</span>
</span><spanid="L-43"><ahref="#L-43"><spanclass="linenos"> 43</span></a><spanclass="sd">The samples can either be provided as python list or as numpy array.</span>
</span><spanid="L-44"><ahref="#L-44"><spanclass="linenos"> 44</span></a><spanclass="sd">The second argument is a list containing the names of the respective Monte Carlo chains as strings. These strings uniquely identify a Monte Carlo chain/ensemble.</span>
</span><spanid="L-54"><ahref="#L-54"><spanclass="linenos"> 54</span></a><spanclass="sd">When performing mathematical operations on `Obs` objects the correct error propagation is intrinsically taken care of using a first order Taylor expansion</span>
</span><spanid="L-56"><ahref="#L-56"><spanclass="linenos"> 56</span></a><spanclass="sd">as introduced in [arXiv:hep-lat/0306017](https://arxiv.org/abs/hep-lat/0306017).</span>
</span><spanid="L-57"><ahref="#L-57"><spanclass="linenos"> 57</span></a><spanclass="sd">The required derivatives $\bar{f}_\alpha$ are evaluated up to machine precision via automatic differentiation as suggested in [arXiv:1809.01289](https://arxiv.org/abs/1809.01289).</span>
</span><spanid="L-59"><ahref="#L-59"><spanclass="linenos"> 59</span></a><spanclass="sd">The `Obs` class is designed such that mathematical numpy functions can be used on `Obs` just as for regular floats.</span>
</span><spanid="L-73"><ahref="#L-73"><spanclass="linenos"> 73</span></a><spanclass="sd"># Check that value and fluctuations are zero within machine precision</span>
</span><spanid="L-80"><ahref="#L-80"><spanclass="linenos"> 80</span></a><spanclass="sd">The error estimation within `pyerrors` is based on the gamma method introduced in [arXiv:hep-lat/0306017](https://arxiv.org/abs/hep-lat/0306017).</span>
</span><spanid="L-81"><ahref="#L-81"><spanclass="linenos"> 81</span></a><spanclass="sd">After having arrived at the derived quantity of interest the `gamma_method` can be called as detailed in the following example.</span>
</span><spanid="L-95"><ahref="#L-95"><spanclass="linenos"> 95</span></a><spanclass="sd">We use the following definition of the integrated autocorrelation time established in [Madras & Sokal 1988](https://link.springer.com/article/10.1007/BF01022990)</span>
</span><spanid="L-97"><ahref="#L-97"><spanclass="linenos"> 97</span></a><spanclass="sd">The window $W$ is determined via the automatic windowing procedure described in [arXiv:hep-lat/0306017](https://arxiv.org/abs/hep-lat/0306017).</span>
</span><spanid="L-98"><ahref="#L-98"><spanclass="linenos"> 98</span></a><spanclass="sd">The standard value for the parameter $S$ of this automatic windowing procedure is $S=2$. Other values for $S$ can be passed to the `gamma_method` as parameter.</span>
</span><spanid="L-110"><ahref="#L-110"><spanclass="linenos">110</span></a><spanclass="sd">The integrated autocorrelation time $\tau_\mathrm{int}$ and the autocorrelation function $\rho(W)$ can be monitored via the methods `pyerrors.obs.Obs.plot_tauint` and `pyerrors.obs.Obs.plot_rho`.</span>
</span><spanid="L-112"><ahref="#L-112"><spanclass="linenos">112</span></a><spanclass="sd">If the parameter $S$ is set to zero it is assumed that the dataset does not exhibit any autocorrelation and the window size is chosen to be zero.</span>
</span><spanid="L-113"><ahref="#L-113"><spanclass="linenos">113</span></a><spanclass="sd">In this case the error estimate is identical to the sample standard error.</span>
</span><spanid="L-117"><ahref="#L-117"><spanclass="linenos">117</span></a><spanclass="sd">Slow modes in the Monte Carlo history can be accounted for by attaching an exponential tail to the autocorrelation function $\rho$ as suggested in [arXiv:1009.5228](https://arxiv.org/abs/1009.5228). The longest autocorrelation time in the history, $\tau_\mathrm{exp}$, can be passed to the `gamma_method` as parameter. In this case the automatic windowing procedure is vacated and the parameter $S$ does not affect the error estimate.</span>
</span><spanid="L-132"><ahref="#L-132"><spanclass="linenos">132</span></a><spanclass="sd">Error propagation for multiple ensembles (Markov chains with different simulation parameters) is handled automatically. Ensembles are uniquely identified by their `name`.</span>
</span><spanid="L-146"><ahref="#L-146"><spanclass="linenos">146</span></a><spanclass="sd">`pyerrors` identifies multiple replica (independent Markov chains with identical simulation parameters) by the vertical bar `|` in the name of the data set.</span>
</span><spanid="L-163"><ahref="#L-163"><spanclass="linenos">163</span></a><spanclass="sd">In order to keep track of different error analysis parameters for different ensembles one can make use of global dictionaries as detailed in the following example.</span>
</span><spanid="L-171"><ahref="#L-171"><spanclass="linenos">171</span></a><spanclass="sd">In case the `gamma_method` is called without any parameters it will use the values specified in the dictionaries for the respective ensembles.</span>
</span><spanid="L-172"><ahref="#L-172"><spanclass="linenos">172</span></a><spanclass="sd">Passing arguments to the `gamma_method` still dominates over the dictionaries.</span>
</span><spanid="L-177"><ahref="#L-177"><spanclass="linenos">177</span></a><spanclass="sd">`Obs` objects defined on irregular Monte Carlo chains can be initialized with the parameter `idl`.</span>
</span><spanid="L-187"><ahref="#L-187"><spanclass="linenos">187</span></a><spanclass="sd"># Observable defined on every second configuration between 5 and 1003</span>
</span><spanid="L-190"><ahref="#L-190"><spanclass="linenos">190</span></a><spanclass="sd">> Result 9.99100712e-01</span>
</span><spanid="L-191"><ahref="#L-191"><spanclass="linenos">191</span></a><spanclass="sd">> 500 samples in 1 ensemble:</span>
</span><spanid="L-192"><ahref="#L-192"><spanclass="linenos">192</span></a><spanclass="sd">> · Ensemble 'ensemble1' : 500 configurations (from 5 to 1003 in steps of 2)</span>
</span><spanid="L-194"><ahref="#L-194"><spanclass="linenos">194</span></a><spanclass="sd"># Observable defined on configurations 2, 9, 28, 29 and 501</span>
</span><spanid="L-203"><ahref="#L-203"><spanclass="linenos">203</span></a><spanclass="sd">`Obs` objects defined on regular and irregular histories of the same ensemble can be combined with each other and the correct error propagation and estimation is automatically taken care of.</span>
</span><spanid="L-205"><ahref="#L-205"><spanclass="linenos">205</span></a><spanclass="sd">**Warning:** Irregular Monte Carlo chains can result in odd patterns in the autocorrelation functions.</span>
</span><spanid="L-206"><ahref="#L-206"><spanclass="linenos">206</span></a><spanclass="sd">Make sure to check the autocorrelation time with e.g. `pyerrors.obs.Obs.plot_rho` or `pyerrors.obs.Obs.plot_tauint`.</span>
</span><spanid="L-211"><ahref="#L-211"><spanclass="linenos">211</span></a><spanclass="sd">When one is not interested in single observables but correlation functions, `pyerrors` offers the `Corr` class which simplifies the corresponding error propagation and provides the user with a set of standard methods. In order to initialize a `Corr` objects one needs to arrange the data as a list of `Obs`</span>
</span><spanid="L-222"><ahref="#L-222"><spanclass="linenos">222</span></a><spanclass="sd">In case the correlation functions are not defined on the outermost timeslices, for example because of fixed boundary conditions, a padding can be introduced.</span>
</span><spanid="L-235"><ahref="#L-235"><spanclass="linenos">235</span></a><spanclass="sd">The individual entries of a correlator can be accessed via slicing</span>
</span><spanid="L-240"><ahref="#L-240"><spanclass="linenos">240</span></a><spanclass="sd">Error propagation with the `Corr` class works very similar to `Obs` objects. Mathematical operations are overloaded and `Corr` objects can be computed together with other `Corr` objects, `Obs` objects or real numbers and integers.</span>
</span><spanid="L-245"><ahref="#L-245"><spanclass="linenos">245</span></a><spanclass="sd">`pyerrors` provides the user with a set of regularly used methods for the manipulation of correlator objects:</span>
</span><spanid="L-246"><ahref="#L-246"><spanclass="linenos">246</span></a><spanclass="sd">- `Corr.gamma_method` applies the gamma method to all entries of the correlator.</span>
</span><spanid="L-247"><ahref="#L-247"><spanclass="linenos">247</span></a><spanclass="sd">- `Corr.m_eff` to construct effective masses. Various variants for periodic and fixed temporal boundary conditions are available.</span>
</span><spanid="L-248"><ahref="#L-248"><spanclass="linenos">248</span></a><spanclass="sd">- `Corr.deriv` returns the first derivative of the correlator as `Corr`. Different discretizations of the numerical derivative are available.</span>
</span><spanid="L-249"><ahref="#L-249"><spanclass="linenos">249</span></a><spanclass="sd">- `Corr.second_deriv` returns the second derivative of the correlator as `Corr`. Different discretizations of the numerical derivative are available.</span>
</span><spanid="L-252"><ahref="#L-252"><spanclass="linenos">252</span></a><spanclass="sd">- `Corr.T_symmetry` averages a correlator with its time symmetry partner, assuming fixed boundary conditions.</span>
</span><spanid="L-253"><ahref="#L-253"><spanclass="linenos">253</span></a><spanclass="sd">- `Corr.plateau` extracts a plateau value from the correlator in a given range.</span>
</span><spanid="L-254"><ahref="#L-254"><spanclass="linenos">254</span></a><spanclass="sd">- `Corr.roll` periodically shifts the correlator.</span>
</span><spanid="L-255"><ahref="#L-255"><spanclass="linenos">255</span></a><spanclass="sd">- `Corr.reverse` reverses the time ordering of the correlator.</span>
</span><spanid="L-256"><ahref="#L-256"><spanclass="linenos">256</span></a><spanclass="sd">- `Corr.correlate` constructs a disconnected correlation function from the correlator and another `Corr` or `Obs` object.</span>
</span><spanid="L-257"><ahref="#L-257"><spanclass="linenos">257</span></a><spanclass="sd">- `Corr.reweight` reweights the correlator.</span>
</span><spanid="L-259"><ahref="#L-259"><spanclass="linenos">259</span></a><spanclass="sd">`pyerrors` can also handle matrices of correlation functions and extract energy states from these matrices via a generalized eigenvalue problem (see `pyerrors.correlators.Corr.GEVP`).</span>
</span><spanid="L-265"><ahref="#L-265"><spanclass="linenos">265</span></a><spanclass="sd">`pyerrors` can handle complex valued observables via the class `pyerrors.obs.CObs`.</span>
</span><spanid="L-266"><ahref="#L-266"><spanclass="linenos">266</span></a><spanclass="sd">`CObs` are initialized with a real and an imaginary part which both can be `Obs` valued.</span>
</span><spanid="L-278"><ahref="#L-278"><spanclass="linenos">278</span></a><spanclass="sd">Elementary mathematical operations are overloaded and samples are properly propagated as for the `Obs` class.</span>
</span><spanid="L-286"><ahref="#L-286"><spanclass="linenos">286</span></a><spanclass="sd"># The `Covobs` class</span>
</span><spanid="L-287"><ahref="#L-287"><spanclass="linenos">287</span></a><spanclass="sd">In many projects, auxiliary data that is not based on Monte Carlo chains enters. Examples are experimentally determined mesons masses which are used to set the scale or renormalization constants. These numbers come with an error that has to be propagated through the analysis. The `Covobs` class allows to define such quantities in `pyerrors`. Furthermore, external input might consist of correlated quantities. An example are the parameters of an interpolation formula, which are defined via mean values and a covariance matrix between all parameters. The contribution of the interpolation formula to the error of a derived quantity therefore might depend on the complete covariance matrix.</span>
</span><spanid="L-289"><ahref="#L-289"><spanclass="linenos">289</span></a><spanclass="sd">This concept is built into the definition of `Covobs`. In `pyerrors`, external input is defined by $M$ mean values, a $M\times M$ covariance matrix, where $M=1$ is permissible, and a name that uniquely identifies the covariance matrix. Below, we define the pion mass, based on its mean value and error, 134.9768(5). Note, that the square of the error enters `cov_Obs`, since the second argument of this function is the covariance matrix of the `Covobs`.</span>
</span><spanid="L-302"><ahref="#L-302"><spanclass="linenos">302</span></a><spanclass="sd">The resulting object `mpi` is an `Obs` that contains a `Covobs`. In the following, it may be handled as any other `Obs`. The contribution of the covariance matrix to the error of an `Obs` is determined from the $M \times M$ covariance matrix $\Sigma$ and the gradient of the `Obs` with respect to the external quantities, which is the $1\times M$ Jacobian matrix $J$, via</span>
</span><spanid="L-304"><ahref="#L-304"><spanclass="linenos">304</span></a><spanclass="sd">where the Jacobian is computed for each derived quantity via automatic differentiation.</span>
</span><spanid="L-306"><ahref="#L-306"><spanclass="linenos">306</span></a><spanclass="sd">Correlated auxiliary data is defined similarly to above, e.g., via</span>
</span><spanid="L-312"><ahref="#L-312"><spanclass="linenos">312</span></a><spanclass="sd">where `RAP` now is a list of two `Obs` that contains the two correlated parameters.</span>
</span><spanid="L-314"><ahref="#L-314"><spanclass="linenos">314</span></a><spanclass="sd">Since the gradient of a derived observable with respect to an external covariance matrix is propagated through the entire analysis, the `Covobs` class allows to quote the derivative of a result with respect to the external quantities. If these derivatives are published together with the result, small shifts in the definition of external quantities, e.g., the definition of the physical point, can be performed a posteriori based on the published information. This may help to compare results of different groups. The gradient of an `Obs` `o` with respect to a covariance matrix with the identifying string `k` may be accessed via</span>
</span><spanid="L-321"><ahref="#L-321"><spanclass="linenos">321</span></a><spanclass="sd">`pyerrors` supports exact linear error propagation for iterative algorithms like various variants of non-linear least squares fits or root finding. The derivatives required for the error propagation are calculated as described in [arXiv:1809.01289](https://arxiv.org/abs/1809.01289).</span>
</span><spanid="L-325"><ahref="#L-325"><spanclass="linenos">325</span></a><spanclass="sd">Standard non-linear least square fits with errors on the dependent but not the independent variables can be performed with `pyerrors.fits.least_squares`. As default solver the Levenberg-Marquardt algorithm implemented in [scipy](https://docs.scipy.org/doc/scipy/reference/generated/scipy.optimize.least_squares.html) is used.</span>
</span><spanid="L-334"><ahref="#L-334"><spanclass="linenos">334</span></a><spanclass="sd">**It is important that numerical functions refer to `autograd.numpy` instead of `numpy` for the automatic differentiation in iterative algorithms to work properly.**</span>
</span><spanid="L-352"><ahref="#L-352"><spanclass="linenos">352</span></a><spanclass="sd">where x is a `list` or `numpy.array` of `floats` and y is a `list` or `numpy.array` of `Obs`.</span>
</span><spanid="L-354"><ahref="#L-354"><spanclass="linenos">354</span></a><spanclass="sd">Data stored in `Corr` objects can be fitted directly using the `Corr.fit` method.</span>
</span><spanid="L-359"><ahref="#L-359"><spanclass="linenos">359</span></a><spanclass="sd">this can simplify working with absolute fit ranges and takes care of gaps in the data automatically.</span>
</span><spanid="L-361"><ahref="#L-361"><spanclass="linenos">361</span></a><spanclass="sd">For fit functions with multiple independent variables the fit function can be of the form</span>
</span><spanid="L-369"><ahref="#L-369"><spanclass="linenos">369</span></a><spanclass="sd">`pyerrors` also supports correlated fits which can be triggered via the parameter `correlated_fit=True`.</span>
</span><spanid="L-370"><ahref="#L-370"><spanclass="linenos">370</span></a><spanclass="sd">Details about how the required covariance matrix is estimated can be found in `pyerrors.obs.covariance`.</span>
</span><spanid="L-372"><ahref="#L-372"><spanclass="linenos">372</span></a><spanclass="sd">Direct visualizations of the performed fits can be triggered via `resplot=True` or `qqplot=True`. For all available options see `pyerrors.fits.least_squares`.</span>
</span><spanid="L-374"><ahref="#L-374"><spanclass="linenos">374</span></a><spanclass="sd">## Total least squares fits</span>
</span><spanid="L-375"><ahref="#L-375"><spanclass="linenos">375</span></a><spanclass="sd">`pyerrors` can also fit data with errors on both the dependent and independent variables using the total least squares method also referred to orthogonal distance regression as implemented in [scipy](https://docs.scipy.org/doc/scipy/reference/odr.html), see `pyerrors.fits.least_squares`. The syntax is identical to the standard least squares case, the only difference being that `x` also has to be a `list` or `numpy.array` of `Obs`.</span>
</span><spanid="L-377"><ahref="#L-377"><spanclass="linenos">377</span></a><spanclass="sd">For the full API see `pyerrors.fits` for fits and `pyerrors.roots` for finding roots of functions.</span>
</span><spanid="L-380"><ahref="#L-380"><spanclass="linenos">380</span></a><spanclass="sd">`pyerrors` provides wrappers for `Obs`- and `CObs`-valued matrix operations based on `numpy.linalg`. The supported functions include:</span>
</span><spanid="L-381"><ahref="#L-381"><spanclass="linenos">381</span></a><spanclass="sd">- `inv` for the matrix inverse.</span>
</span><spanid="L-382"><ahref="#L-382"><spanclass="linenos">382</span></a><spanclass="sd">- `cholseky` for the Cholesky decomposition.</span>
</span><spanid="L-383"><ahref="#L-383"><spanclass="linenos">383</span></a><spanclass="sd">- `det` for the matrix determinant.</span>
</span><spanid="L-384"><ahref="#L-384"><spanclass="linenos">384</span></a><spanclass="sd">- `eigh` for eigenvalues and eigenvectors of hermitean matrices.</span>
</span><spanid="L-385"><ahref="#L-385"><spanclass="linenos">385</span></a><spanclass="sd">- `eig` for eigenvalues of general matrices.</span>
</span><spanid="L-386"><ahref="#L-386"><spanclass="linenos">386</span></a><spanclass="sd">- `pinv` for the Moore-Penrose pseudoinverse.</span>
</span><spanid="L-387"><ahref="#L-387"><spanclass="linenos">387</span></a><spanclass="sd">- `svd` for the singular-value-decomposition.</span>
</span><spanid="L-393"><ahref="#L-393"><spanclass="linenos">393</span></a><spanclass="sd">The preferred exported file format within `pyerrors` is json.gz. Files written to this format are valid JSON files that have been compressed using gzip. The structure of the content is inspired by the dobs format of the ALPHA collaboration. The aim of the format is to facilitate the storage of data in a self-contained way such that, even years after the creation of the file, it is possible to extract all necessary information:</span>
</span><spanid="L-394"><ahref="#L-394"><spanclass="linenos">394</span></a><spanclass="sd">- What observables are stored? Possibly: How exactly are they defined.</span>
</span><spanid="L-395"><ahref="#L-395"><spanclass="linenos">395</span></a><spanclass="sd">- How does each single ensemble or external quantity contribute to the error of the observable?</span>
</span><spanid="L-396"><ahref="#L-396"><spanclass="linenos">396</span></a><spanclass="sd">- Who did write the file when and on which machine?</span>
</span><spanid="L-398"><ahref="#L-398"><spanclass="linenos">398</span></a><spanclass="sd">This can be achieved by storing all information in one single file. The export routines of `pyerrors` are written such that as much information as possible is written automatically as described in the following example</span>
</span><spanid="L-403"><ahref="#L-403"><spanclass="linenos">403</span></a><spanclass="sd">pe.input.json.dump_to_json(my_obs, "test_output_file", description="This file contains a test observable")</span>
</span><spanid="L-404"><ahref="#L-404"><spanclass="linenos">404</span></a><spanclass="sd"># For a single observable one can equivalently use the class method dump</span>
</span><spanid="L-405"><ahref="#L-405"><spanclass="linenos">405</span></a><spanclass="sd">my_obs.dump("test_output_file", description="This file contains a test observable")</span>
</span><spanid="L-412"><ahref="#L-412"><spanclass="linenos">412</span></a><spanclass="sd">The format also allows to directly write out the content of `Corr` objects or lists and arrays of `Obs` objects by passing the desired data to `pyerrors.input.json.dump_to_json`.</span>
</span><spanid="L-414"><ahref="#L-414"><spanclass="linenos">414</span></a><spanclass="sd">## json.gz format specification</span>
</span><spanid="L-415"><ahref="#L-415"><spanclass="linenos">415</span></a><spanclass="sd">The first entries of the file provide optional auxiliary information:</span>
</span><spanid="L-416"><ahref="#L-416"><spanclass="linenos">416</span></a><spanclass="sd">- `program` is a string that indicates which program was used to write the file.</span>
</span><spanid="L-417"><ahref="#L-417"><spanclass="linenos">417</span></a><spanclass="sd">- `version` is a string that specifies the version of the format.</span>
</span><spanid="L-418"><ahref="#L-418"><spanclass="linenos">418</span></a><spanclass="sd">- `who` is a string that specifies the user name of the creator of the file.</span>
</span><spanid="L-419"><ahref="#L-419"><spanclass="linenos">419</span></a><spanclass="sd">- `date` is a string and contains the creation date of the file.</span>
</span><spanid="L-420"><ahref="#L-420"><spanclass="linenos">420</span></a><spanclass="sd">- `host` is a string and contains the hostname of the machine where the file has been written.</span>
</span><spanid="L-421"><ahref="#L-421"><spanclass="linenos">421</span></a><spanclass="sd">- `description` contains information on the content of the file. This field is not filled automatically in `pyerrors`. The user is advised to provide as detailed information as possible in this field. Examples are: Input files of measurements or simulations, LaTeX formulae or references to publications to specify how the observables have been computed, details on the analysis strategy, ... This field may be any valid JSON type. Strings, arrays or objects (equivalent to dicts in python) are well suited to provide information.</span>
</span><spanid="L-426"><ahref="#L-426"><spanclass="linenos">426</span></a><spanclass="sd">Each entry of the array belongs to a single structure of observables. Currently, these structures can be either of `Obs`, `list`, `numpy.ndarray`, `Corr`. All `Obs` inside a structure (with dimension > 0) have to be defined on the same set of configurations. Different structures, that are represented by entries of the array `obsdata`, are treated independently. Each entry of the array `obsdata` has the following required entries:</span>
</span><spanid="L-427"><ahref="#L-427"><spanclass="linenos">427</span></a><spanclass="sd">- `type` is a string that specifies the type of the structure. This allows to parse the content to the correct form after reading the file. It is always possible to interpret the content as list of Obs.</span>
</span><spanid="L-428"><ahref="#L-428"><spanclass="linenos">428</span></a><spanclass="sd">- `value` is an array that contains the mean values of the Obs inside the structure.</span>
</span><spanid="L-429"><ahref="#L-429"><spanclass="linenos">429</span></a><spanclass="sd">The following entries are optional:</span>
</span><spanid="L-430"><ahref="#L-430"><spanclass="linenos">430</span></a><spanclass="sd">- `layout` is a string that specifies the layout of multi-dimensional structures. Examples are "2, 2" for a 2x2 dimensional matrix or "64, 4, 4" for a Corr with $T=64$ and 4x4 matrices on each time slices. "1" denotes a single Obs. Multi-dimensional structures are stored in row-major format (see below).</span>
</span><spanid="L-431"><ahref="#L-431"><spanclass="linenos">431</span></a><spanclass="sd">- `tag` is any JSON type. It contains additional information concerning the structure. The `tag` of an `Obs` in `pyerrors` is written here.</span>
</span><spanid="L-432"><ahref="#L-432"><spanclass="linenos">432</span></a><spanclass="sd">- `reweighted` is a Bool that may be used to specify, whether the `Obs` in the structure have been reweighted.</span>
</span><spanid="L-433"><ahref="#L-433"><spanclass="linenos">433</span></a><spanclass="sd">- `data` is an array that contains the data from MC chains. We will define it below.</span>
</span><spanid="L-434"><ahref="#L-434"><spanclass="linenos">434</span></a><spanclass="sd">- `cdata` is an array that contains the data from external quantities with an error (`Covobs` in `pyerrors`). We will define it below.</span>
</span><spanid="L-436"><ahref="#L-436"><spanclass="linenos">436</span></a><spanclass="sd">The array `data` contains the data from MC chains. Each entry of the array corresponds to one ensemble and contains:</span>
</span><spanid="L-437"><ahref="#L-437"><spanclass="linenos">437</span></a><spanclass="sd">- `id`, a string that contains the name of the ensemble</span>
</span><spanid="L-438"><ahref="#L-438"><spanclass="linenos">438</span></a><spanclass="sd">- `replica`, an array that contains an entry per replica of the ensemble.</span>
</span><spanid="L-440"><ahref="#L-440"><spanclass="linenos">440</span></a><spanclass="sd">Each entry of `replica` contains</span>
</span><spanid="L-441"><ahref="#L-441"><spanclass="linenos">441</span></a><spanclass="sd">`name`, a string that contains the name of the replica</span>
</span><spanid="L-442"><ahref="#L-442"><spanclass="linenos">442</span></a><spanclass="sd">`deltas`, an array that contains the actual data.</span>
</span><spanid="L-444"><ahref="#L-444"><spanclass="linenos">444</span></a><spanclass="sd">Each entry in `deltas` corresponds to one configuration of the replica and has $1+N$ many entries. The first entry is an integer that specifies the configuration number that, together with ensemble and replica name, may be used to uniquely identify the configuration on which the data has been obtained. The following N entries specify the deltas, i.e., the deviation of the observable from the mean value on this configuration, of each `Obs` inside the structure. Multi-dimensional structures are stored in a row-major format. For primary observables, such as correlation functions, $value + delta_i$ matches the primary data obtained on the configuration.</span>
</span><spanid="L-446"><ahref="#L-446"><spanclass="linenos">446</span></a><spanclass="sd">The array `cdata` contains information about the contribution of auxiliary observables, represented by `Covobs` in `pyerrors`, to the total error of the observables. Each entry of the array belongs to one auxiliary covariance matrix and contains:</span>
</span><spanid="L-447"><ahref="#L-447"><spanclass="linenos">447</span></a><spanclass="sd">- `id`, a string that identifies the covariance matrix</span>
</span><spanid="L-448"><ahref="#L-448"><spanclass="linenos">448</span></a><spanclass="sd">- `layout`, a string that defines the dimensions of the $M\times M$ covariance matrix (has to be "M, M" or "1").</span>
</span><spanid="L-449"><ahref="#L-449"><spanclass="linenos">449</span></a><spanclass="sd">- `cov`, an array that contains the $M\times M$ many entries of the covariance matrix, stored in row-major format.</span>
</span><spanid="L-450"><ahref="#L-450"><spanclass="linenos">450</span></a><spanclass="sd">- `grad`, an array that contains N entries, one for each `Obs` inside the structure. Each entry itself is an array, that contains the M gradients of the Nth observable with respect to the quantity that corresponds to the Mth diagonal entry of the covariance matrix.</span>
</span><spanid="L-452"><ahref="#L-452"><spanclass="linenos">452</span></a><spanclass="sd">A JSON schema that may be used to verify the correctness of a file with respect to the format definition is stored in ./examples/json_schema.json. The schema is a self-descriptive format definition and contains an exemplary file.</span>
</span><spanid="L-454"><ahref="#L-454"><spanclass="linenos">454</span></a><spanclass="sd">Julia I/O routines for the json.gz format, compatible with [ADerrors.jl](https://gitlab.ift.uam-csic.es/alberto/aderrors.jl), can be found [here](https://github.com/fjosw/ADjson.jl).</span>