pyerrors.correlators

   1import warnings
   2from itertools import permutations
   3import numpy as np
   4import autograd.numpy as anp
   5import matplotlib.pyplot as plt
   6import scipy.linalg
   7from .obs import Obs, reweight, correlate, CObs
   8from .misc import dump_object, _assert_equal_properties
   9from .fits import least_squares
  10from .roots import find_root
  11
  12
  13class Corr:
  14    """The class for a correlator (time dependent sequence of pe.Obs).
  15
  16    Everything, this class does, can be achieved using lists or arrays of Obs.
  17    But it is simply more convenient to have a dedicated object for correlators.
  18    One often wants to add or multiply correlators of the same length at every timeslice and it is inconvenient
  19    to iterate over all timeslices for every operation. This is especially true, when dealing with matrices.
  20
  21    The correlator can have two types of content: An Obs at every timeslice OR a GEVP
  22    matrix at every timeslice. Other dependency (eg. spatial) are not supported.
  23
  24    """
  25
  26    def __init__(self, data_input, padding=[0, 0], prange=None):
  27        """ Initialize a Corr object.
  28
  29        Parameters
  30        ----------
  31        data_input : list or array
  32            list of Obs or list of arrays of Obs or array of Corrs
  33        padding : list, optional
  34            List with two entries where the first labels the padding
  35            at the front of the correlator and the second the padding
  36            at the back.
  37        prange : list, optional
  38            List containing the first and last timeslice of the plateau
  39            region indentified for this correlator.
  40        """
  41
  42        if isinstance(data_input, np.ndarray):
  43
  44            # This only works, if the array fulfills the conditions below
  45            if not len(data_input.shape) == 2 and data_input.shape[0] == data_input.shape[1]:
  46                raise Exception("Incompatible array shape")
  47            if not all([isinstance(item, Corr) for item in data_input.flatten()]):
  48                raise Exception("If the input is an array, its elements must be of type pe.Corr")
  49            if not all([item.N == 1 for item in data_input.flatten()]):
  50                raise Exception("Can only construct matrix correlator from single valued correlators")
  51            if not len(set([item.T for item in data_input.flatten()])) == 1:
  52                raise Exception("All input Correlators must be defined over the same timeslices.")
  53
  54            T = data_input[0, 0].T
  55            N = data_input.shape[0]
  56            input_as_list = []
  57            for t in range(T):
  58                if any([(item.content[t] is None) for item in data_input.flatten()]):
  59                    if not all([(item.content[t] is None) for item in data_input.flatten()]):
  60                        warnings.warn("Input ill-defined at different timeslices. Conversion leads to data loss!", RuntimeWarning)
  61                    input_as_list.append(None)
  62                else:
  63                    array_at_timeslace = np.empty([N, N], dtype="object")
  64                    for i in range(N):
  65                        for j in range(N):
  66                            array_at_timeslace[i, j] = data_input[i, j][t]
  67                    input_as_list.append(array_at_timeslace)
  68            data_input = input_as_list
  69
  70        if isinstance(data_input, list):
  71
  72            if all([isinstance(item, (Obs, CObs)) or item is None for item in data_input]):
  73                _assert_equal_properties([o for o in data_input if o is not None])
  74                self.content = [np.asarray([item]) if item is not None else None for item in data_input]
  75                self.N = 1
  76
  77            elif all([isinstance(item, np.ndarray) or item is None for item in data_input]) and any([isinstance(item, np.ndarray) for item in data_input]):
  78                self.content = data_input
  79                noNull = [a for a in self.content if not (a is None)]  # To check if the matrices are correct for all undefined elements
  80                self.N = noNull[0].shape[0]
  81                if self.N > 1 and noNull[0].shape[0] != noNull[0].shape[1]:
  82                    raise Exception("Smearing matrices are not NxN")
  83                if (not all([item.shape == noNull[0].shape for item in noNull])):
  84                    raise Exception("Items in data_input are not of identical shape." + str(noNull))
  85            else:
  86                raise Exception("data_input contains item of wrong type")
  87        else:
  88            raise Exception("Data input was not given as list or correct array")
  89
  90        self.tag = None
  91
  92        # An undefined timeslice is represented by the None object
  93        self.content = [None] * padding[0] + self.content + [None] * padding[1]
  94        self.T = len(self.content)
  95        self.prange = prange
  96
  97    def __getitem__(self, idx):
  98        """Return the content of timeslice idx"""
  99        if self.content[idx] is None:
 100            return None
 101        elif len(self.content[idx]) == 1:
 102            return self.content[idx][0]
 103        else:
 104            return self.content[idx]
 105
 106    @property
 107    def reweighted(self):
 108        bool_array = np.array([list(map(lambda x: x.reweighted, o)) for o in [x for x in self.content if x is not None]])
 109        if np.all(bool_array == 1):
 110            return True
 111        elif np.all(bool_array == 0):
 112            return False
 113        else:
 114            raise Exception("Reweighting status of correlator corrupted.")
 115
 116    def gamma_method(self, **kwargs):
 117        """Apply the gamma method to the content of the Corr."""
 118        for item in self.content:
 119            if not (item is None):
 120                if self.N == 1:
 121                    item[0].gamma_method(**kwargs)
 122                else:
 123                    for i in range(self.N):
 124                        for j in range(self.N):
 125                            item[i, j].gamma_method(**kwargs)
 126
 127    def projected(self, vector_l=None, vector_r=None, normalize=False):
 128        """We need to project the Correlator with a Vector to get a single value at each timeslice.
 129
 130        The method can use one or two vectors.
 131        If two are specified it returns v1@G@v2 (the order might be very important.)
 132        By default it will return the lowest source, which usually means unsmeared-unsmeared (0,0), but it does not have to
 133        """
 134        if self.N == 1:
 135            raise Exception("Trying to project a Corr, that already has N=1.")
 136
 137        if vector_l is None:
 138            vector_l, vector_r = np.asarray([1.] + (self.N - 1) * [0.]), np.asarray([1.] + (self.N - 1) * [0.])
 139        elif (vector_r is None):
 140            vector_r = vector_l
 141        if isinstance(vector_l, list) and not isinstance(vector_r, list):
 142            if len(vector_l) != self.T:
 143                raise Exception("Length of vector list must be equal to T")
 144            vector_r = [vector_r] * self.T
 145        if isinstance(vector_r, list) and not isinstance(vector_l, list):
 146            if len(vector_r) != self.T:
 147                raise Exception("Length of vector list must be equal to T")
 148            vector_l = [vector_l] * self.T
 149
 150        if not isinstance(vector_l, list):
 151            if not vector_l.shape == vector_r.shape == (self.N,):
 152                raise Exception("Vectors are of wrong shape!")
 153            if normalize:
 154                vector_l, vector_r = vector_l / np.sqrt((vector_l @ vector_l)), vector_r / np.sqrt(vector_r @ vector_r)
 155            newcontent = [None if _check_for_none(self, item) else np.asarray([vector_l.T @ item @ vector_r]) for item in self.content]
 156
 157        else:
 158            # There are no checks here yet. There are so many possible scenarios, where this can go wrong.
 159            if normalize:
 160                for t in range(self.T):
 161                    vector_l[t], vector_r[t] = vector_l[t] / np.sqrt((vector_l[t] @ vector_l[t])), vector_r[t] / np.sqrt(vector_r[t] @ vector_r[t])
 162
 163            newcontent = [None if (_check_for_none(self, self.content[t]) or vector_l[t] is None or vector_r[t] is None) else np.asarray([vector_l[t].T @ self.content[t] @ vector_r[t]]) for t in range(self.T)]
 164        return Corr(newcontent)
 165
 166    def item(self, i, j):
 167        """Picks the element [i,j] from every matrix and returns a correlator containing one Obs per timeslice.
 168
 169        Parameters
 170        ----------
 171        i : int
 172            First index to be picked.
 173        j : int
 174            Second index to be picked.
 175        """
 176        if self.N == 1:
 177            raise Exception("Trying to pick item from projected Corr")
 178        newcontent = [None if (item is None) else item[i, j] for item in self.content]
 179        return Corr(newcontent)
 180
 181    def plottable(self):
 182        """Outputs the correlator in a plotable format.
 183
 184        Outputs three lists containing the timeslice index, the value on each
 185        timeslice and the error on each timeslice.
 186        """
 187        if self.N != 1:
 188            raise Exception("Can only make Corr[N=1] plottable")
 189        x_list = [x for x in range(self.T) if not self.content[x] is None]
 190        y_list = [y[0].value for y in self.content if y is not None]
 191        y_err_list = [y[0].dvalue for y in self.content if y is not None]
 192
 193        return x_list, y_list, y_err_list
 194
 195    def symmetric(self):
 196        """ Symmetrize the correlator around x0=0."""
 197        if self.N != 1:
 198            raise Exception('symmetric cannot be safely applied to multi-dimensional correlators.')
 199        if self.T % 2 != 0:
 200            raise Exception("Can not symmetrize odd T")
 201
 202        if np.argmax(np.abs(self.content)) != 0:
 203            warnings.warn("Correlator does not seem to be symmetric around x0=0.", RuntimeWarning)
 204
 205        newcontent = [self.content[0]]
 206        for t in range(1, self.T):
 207            if (self.content[t] is None) or (self.content[self.T - t] is None):
 208                newcontent.append(None)
 209            else:
 210                newcontent.append(0.5 * (self.content[t] + self.content[self.T - t]))
 211        if (all([x is None for x in newcontent])):
 212            raise Exception("Corr could not be symmetrized: No redundant values")
 213        return Corr(newcontent, prange=self.prange)
 214
 215    def anti_symmetric(self):
 216        """Anti-symmetrize the correlator around x0=0."""
 217        if self.N != 1:
 218            raise Exception('anti_symmetric cannot be safely applied to multi-dimensional correlators.')
 219        if self.T % 2 != 0:
 220            raise Exception("Can not symmetrize odd T")
 221
 222        test = 1 * self
 223        test.gamma_method()
 224        if not all([o.is_zero_within_error(3) for o in test.content[0]]):
 225            warnings.warn("Correlator does not seem to be anti-symmetric around x0=0.", RuntimeWarning)
 226
 227        newcontent = [self.content[0]]
 228        for t in range(1, self.T):
 229            if (self.content[t] is None) or (self.content[self.T - t] is None):
 230                newcontent.append(None)
 231            else:
 232                newcontent.append(0.5 * (self.content[t] - self.content[self.T - t]))
 233        if (all([x is None for x in newcontent])):
 234            raise Exception("Corr could not be symmetrized: No redundant values")
 235        return Corr(newcontent, prange=self.prange)
 236
 237    def is_matrix_symmetric(self):
 238        """Checks whether a correlator matrices is symmetric on every timeslice."""
 239        if self.N == 1:
 240            raise Exception("Only works for correlator matrices.")
 241        for t in range(self.T):
 242            if self[t] is None:
 243                continue
 244            for i in range(self.N):
 245                for j in range(i + 1, self.N):
 246                    if self[t][i, j] is self[t][j, i]:
 247                        continue
 248                    if hash(self[t][i, j]) != hash(self[t][j, i]):
 249                        return False
 250        return True
 251
 252    def matrix_symmetric(self):
 253        """Symmetrizes the correlator matrices on every timeslice."""
 254        if self.N == 1:
 255            raise Exception("Trying to symmetrize a correlator matrix, that already has N=1.")
 256        if self.is_matrix_symmetric():
 257            return 1.0 * self
 258        else:
 259            transposed = [None if _check_for_none(self, G) else G.T for G in self.content]
 260            return 0.5 * (Corr(transposed) + self)
 261
 262    def GEVP(self, t0, ts=None, sort="Eigenvalue", **kwargs):
 263        r'''Solve the generalized eigenvalue problem on the correlator matrix and returns the corresponding eigenvectors.
 264
 265        The eigenvectors are sorted according to the descending eigenvalues, the zeroth eigenvector(s) correspond to the
 266        largest eigenvalue(s). The eigenvector(s) for the individual states can be accessed via slicing
 267        ```python
 268        C.GEVP(t0=2)[0]  # Ground state vector(s)
 269        C.GEVP(t0=2)[:3]  # Vectors for the lowest three states
 270        ```
 271
 272        Parameters
 273        ----------
 274        t0 : int
 275            The time t0 for the right hand side of the GEVP according to $G(t)v_i=\lambda_i G(t_0)v_i$
 276        ts : int
 277            fixed time $G(t_s)v_i=\lambda_i G(t_0)v_i$ if sort=None.
 278            If sort="Eigenvector" it gives a reference point for the sorting method.
 279        sort : string
 280            If this argument is set, a list of self.T vectors per state is returned. If it is set to None, only one vector is returned.
 281            - "Eigenvalue": The eigenvector is chosen according to which eigenvalue it belongs individually on every timeslice.
 282            - "Eigenvector": Use the method described in arXiv:2004.10472 to find the set of v(t) belonging to the state.
 283              The reference state is identified by its eigenvalue at $t=t_s$.
 284
 285        Other Parameters
 286        ----------------
 287        state : int
 288           Returns only the vector(s) for a specified state. The lowest state is zero.
 289        '''
 290
 291        if self.N == 1:
 292            raise Exception("GEVP methods only works on correlator matrices and not single correlators.")
 293        if ts is not None:
 294            if (ts <= t0):
 295                raise Exception("ts has to be larger than t0.")
 296
 297        if "sorted_list" in kwargs:
 298            warnings.warn("Argument 'sorted_list' is deprecated, use 'sort' instead.", DeprecationWarning)
 299            sort = kwargs.get("sorted_list")
 300
 301        if self.is_matrix_symmetric():
 302            symmetric_corr = self
 303        else:
 304            symmetric_corr = self.matrix_symmetric()
 305
 306        G0 = np.vectorize(lambda x: x.value)(symmetric_corr[t0])
 307        np.linalg.cholesky(G0)  # Check if matrix G0 is positive-semidefinite.
 308
 309        if sort is None:
 310            if (ts is None):
 311                raise Exception("ts is required if sort=None.")
 312            if (self.content[t0] is None) or (self.content[ts] is None):
 313                raise Exception("Corr not defined at t0/ts.")
 314            Gt = np.vectorize(lambda x: x.value)(symmetric_corr[ts])
 315            reordered_vecs = _GEVP_solver(Gt, G0)
 316
 317        elif sort in ["Eigenvalue", "Eigenvector"]:
 318            if sort == "Eigenvalue" and ts is not None:
 319                warnings.warn("ts has no effect when sorting by eigenvalue is chosen.", RuntimeWarning)
 320            all_vecs = [None] * (t0 + 1)
 321            for t in range(t0 + 1, self.T):
 322                try:
 323                    Gt = np.vectorize(lambda x: x.value)(symmetric_corr[t])
 324                    all_vecs.append(_GEVP_solver(Gt, G0))
 325                except Exception:
 326                    all_vecs.append(None)
 327            if sort == "Eigenvector":
 328                if ts is None:
 329                    raise Exception("ts is required for the Eigenvector sorting method.")
 330                all_vecs = _sort_vectors(all_vecs, ts)
 331
 332            reordered_vecs = [[v[s] if v is not None else None for v in all_vecs] for s in range(self.N)]
 333        else:
 334            raise Exception("Unkown value for 'sort'.")
 335
 336        if "state" in kwargs:
 337            return reordered_vecs[kwargs.get("state")]
 338        else:
 339            return reordered_vecs
 340
 341    def Eigenvalue(self, t0, ts=None, state=0, sort="Eigenvalue"):
 342        """Determines the eigenvalue of the GEVP by solving and projecting the correlator
 343
 344        Parameters
 345        ----------
 346        state : int
 347            The state one is interested in ordered by energy. The lowest state is zero.
 348
 349        All other parameters are identical to the ones of Corr.GEVP.
 350        """
 351        vec = self.GEVP(t0, ts=ts, sort=sort)[state]
 352        return self.projected(vec)
 353
 354    def Hankel(self, N, periodic=False):
 355        """Constructs an NxN Hankel matrix
 356
 357        C(t) c(t+1) ... c(t+n-1)
 358        C(t+1) c(t+2) ... c(t+n)
 359        .................
 360        C(t+(n-1)) c(t+n) ... c(t+2(n-1))
 361
 362        Parameters
 363        ----------
 364        N : int
 365            Dimension of the Hankel matrix
 366        periodic : bool, optional
 367            determines whether the matrix is extended periodically
 368        """
 369
 370        if self.N != 1:
 371            raise Exception("Multi-operator Prony not implemented!")
 372
 373        array = np.empty([N, N], dtype="object")
 374        new_content = []
 375        for t in range(self.T):
 376            new_content.append(array.copy())
 377
 378        def wrap(i):
 379            while i >= self.T:
 380                i -= self.T
 381            return i
 382
 383        for t in range(self.T):
 384            for i in range(N):
 385                for j in range(N):
 386                    if periodic:
 387                        new_content[t][i, j] = self.content[wrap(t + i + j)][0]
 388                    elif (t + i + j) >= self.T:
 389                        new_content[t] = None
 390                    else:
 391                        new_content[t][i, j] = self.content[t + i + j][0]
 392
 393        return Corr(new_content)
 394
 395    def roll(self, dt):
 396        """Periodically shift the correlator by dt timeslices
 397
 398        Parameters
 399        ----------
 400        dt : int
 401            number of timeslices
 402        """
 403        return Corr(list(np.roll(np.array(self.content, dtype=object), dt)))
 404
 405    def reverse(self):
 406        """Reverse the time ordering of the Corr"""
 407        return Corr(self.content[:: -1])
 408
 409    def thin(self, spacing=2, offset=0):
 410        """Thin out a correlator to suppress correlations
 411
 412        Parameters
 413        ----------
 414        spacing : int
 415            Keep only every 'spacing'th entry of the correlator
 416        offset : int
 417            Offset the equal spacing
 418        """
 419        new_content = []
 420        for t in range(self.T):
 421            if (offset + t) % spacing != 0:
 422                new_content.append(None)
 423            else:
 424                new_content.append(self.content[t])
 425        return Corr(new_content)
 426
 427    def correlate(self, partner):
 428        """Correlate the correlator with another correlator or Obs
 429
 430        Parameters
 431        ----------
 432        partner : Obs or Corr
 433            partner to correlate the correlator with.
 434            Can either be an Obs which is correlated with all entries of the
 435            correlator or a Corr of same length.
 436        """
 437        if self.N != 1:
 438            raise Exception("Only one-dimensional correlators can be safely correlated.")
 439        new_content = []
 440        for x0, t_slice in enumerate(self.content):
 441            if _check_for_none(self, t_slice):
 442                new_content.append(None)
 443            else:
 444                if isinstance(partner, Corr):
 445                    if _check_for_none(partner, partner.content[x0]):
 446                        new_content.append(None)
 447                    else:
 448                        new_content.append(np.array([correlate(o, partner.content[x0][0]) for o in t_slice]))
 449                elif isinstance(partner, Obs):  # Should this include CObs?
 450                    new_content.append(np.array([correlate(o, partner) for o in t_slice]))
 451                else:
 452                    raise Exception("Can only correlate with an Obs or a Corr.")
 453
 454        return Corr(new_content)
 455
 456    def reweight(self, weight, **kwargs):
 457        """Reweight the correlator.
 458
 459        Parameters
 460        ----------
 461        weight : Obs
 462            Reweighting factor. An Observable that has to be defined on a superset of the
 463            configurations in obs[i].idl for all i.
 464        all_configs : bool
 465            if True, the reweighted observables are normalized by the average of
 466            the reweighting factor on all configurations in weight.idl and not
 467            on the configurations in obs[i].idl.
 468        """
 469        if self.N != 1:
 470            raise Exception("Reweighting only implemented for one-dimensional correlators.")
 471        new_content = []
 472        for t_slice in self.content:
 473            if _check_for_none(self, t_slice):
 474                new_content.append(None)
 475            else:
 476                new_content.append(np.array(reweight(weight, t_slice, **kwargs)))
 477        return Corr(new_content)
 478
 479    def T_symmetry(self, partner, parity=+1):
 480        """Return the time symmetry average of the correlator and its partner
 481
 482        Parameters
 483        ----------
 484        partner : Corr
 485            Time symmetry partner of the Corr
 486        partity : int
 487            Parity quantum number of the correlator, can be +1 or -1
 488        """
 489        if self.N != 1:
 490            raise Exception("T_symmetry only implemented for one-dimensional correlators.")
 491        if not isinstance(partner, Corr):
 492            raise Exception("T partner has to be a Corr object.")
 493        if parity not in [+1, -1]:
 494            raise Exception("Parity has to be +1 or -1.")
 495        T_partner = parity * partner.reverse()
 496
 497        t_slices = []
 498        test = (self - T_partner)
 499        test.gamma_method()
 500        for x0, t_slice in enumerate(test.content):
 501            if t_slice is not None:
 502                if not t_slice[0].is_zero_within_error(5):
 503                    t_slices.append(x0)
 504        if t_slices:
 505            warnings.warn("T symmetry partners do not agree within 5 sigma on time slices " + str(t_slices) + ".", RuntimeWarning)
 506
 507        return (self + T_partner) / 2
 508
 509    def deriv(self, variant="symmetric"):
 510        """Return the first derivative of the correlator with respect to x0.
 511
 512        Parameters
 513        ----------
 514        variant : str
 515            decides which definition of the finite differences derivative is used.
 516            Available choice: symmetric, forward, backward, improved, log, default: symmetric
 517        """
 518        if self.N != 1:
 519            raise Exception("deriv only implemented for one-dimensional correlators.")
 520        if variant == "symmetric":
 521            newcontent = []
 522            for t in range(1, self.T - 1):
 523                if (self.content[t - 1] is None) or (self.content[t + 1] is None):
 524                    newcontent.append(None)
 525                else:
 526                    newcontent.append(0.5 * (self.content[t + 1] - self.content[t - 1]))
 527            if (all([x is None for x in newcontent])):
 528                raise Exception('Derivative is undefined at all timeslices')
 529            return Corr(newcontent, padding=[1, 1])
 530        elif variant == "forward":
 531            newcontent = []
 532            for t in range(self.T - 1):
 533                if (self.content[t] is None) or (self.content[t + 1] is None):
 534                    newcontent.append(None)
 535                else:
 536                    newcontent.append(self.content[t + 1] - self.content[t])
 537            if (all([x is None for x in newcontent])):
 538                raise Exception("Derivative is undefined at all timeslices")
 539            return Corr(newcontent, padding=[0, 1])
 540        elif variant == "backward":
 541            newcontent = []
 542            for t in range(1, self.T):
 543                if (self.content[t - 1] is None) or (self.content[t] is None):
 544                    newcontent.append(None)
 545                else:
 546                    newcontent.append(self.content[t] - self.content[t - 1])
 547            if (all([x is None for x in newcontent])):
 548                raise Exception("Derivative is undefined at all timeslices")
 549            return Corr(newcontent, padding=[1, 0])
 550        elif variant == "improved":
 551            newcontent = []
 552            for t in range(2, self.T - 2):
 553                if (self.content[t - 2] is None) or (self.content[t - 1] is None) or (self.content[t + 1] is None) or (self.content[t + 2] is None):
 554                    newcontent.append(None)
 555                else:
 556                    newcontent.append((1 / 12) * (self.content[t - 2] - 8 * self.content[t - 1] + 8 * self.content[t + 1] - self.content[t + 2]))
 557            if (all([x is None for x in newcontent])):
 558                raise Exception('Derivative is undefined at all timeslices')
 559            return Corr(newcontent, padding=[2, 2])
 560        elif variant == 'log':
 561            newcontent = []
 562            for t in range(self.T):
 563                if (self.content[t] is None) or (self.content[t] <= 0):
 564                    newcontent.append(None)
 565                else:
 566                    newcontent.append(np.log(self.content[t]))
 567            if (all([x is None for x in newcontent])):
 568                raise Exception("Log is undefined at all timeslices")
 569            logcorr = Corr(newcontent)
 570            return self * logcorr.deriv('symmetric')
 571        else:
 572            raise Exception("Unknown variant.")
 573
 574    def second_deriv(self, variant="symmetric"):
 575        """Return the second derivative of the correlator with respect to x0.
 576
 577        Parameters
 578        ----------
 579        variant : str
 580            decides which definition of the finite differences derivative is used.
 581            Available choice: symmetric, improved, log, default: symmetric
 582        """
 583        if self.N != 1:
 584            raise Exception("second_deriv only implemented for one-dimensional correlators.")
 585        if variant == "symmetric":
 586            newcontent = []
 587            for t in range(1, self.T - 1):
 588                if (self.content[t - 1] is None) or (self.content[t + 1] is None):
 589                    newcontent.append(None)
 590                else:
 591                    newcontent.append((self.content[t + 1] - 2 * self.content[t] + self.content[t - 1]))
 592            if (all([x is None for x in newcontent])):
 593                raise Exception("Derivative is undefined at all timeslices")
 594            return Corr(newcontent, padding=[1, 1])
 595        elif variant == "improved":
 596            newcontent = []
 597            for t in range(2, self.T - 2):
 598                if (self.content[t - 2] is None) or (self.content[t - 1] is None) or (self.content[t] is None) or (self.content[t + 1] is None) or (self.content[t + 2] is None):
 599                    newcontent.append(None)
 600                else:
 601                    newcontent.append((1 / 12) * (-self.content[t + 2] + 16 * self.content[t + 1] - 30 * self.content[t] + 16 * self.content[t - 1] - self.content[t - 2]))
 602            if (all([x is None for x in newcontent])):
 603                raise Exception("Derivative is undefined at all timeslices")
 604            return Corr(newcontent, padding=[2, 2])
 605        elif variant == 'log':
 606            newcontent = []
 607            for t in range(self.T):
 608                if (self.content[t] is None) or (self.content[t] <= 0):
 609                    newcontent.append(None)
 610                else:
 611                    newcontent.append(np.log(self.content[t]))
 612            if (all([x is None for x in newcontent])):
 613                raise Exception("Log is undefined at all timeslices")
 614            logcorr = Corr(newcontent)
 615            return self * (logcorr.second_deriv('symmetric') + (logcorr.deriv('symmetric'))**2)
 616        else:
 617            raise Exception("Unknown variant.")
 618
 619    def m_eff(self, variant='log', guess=1.0):
 620        """Returns the effective mass of the correlator as correlator object
 621
 622        Parameters
 623        ----------
 624        variant : str
 625            log : uses the standard effective mass log(C(t) / C(t+1))
 626            cosh, periodic : Use periodicitiy of the correlator by solving C(t) / C(t+1) = cosh(m * (t - T/2)) / cosh(m * (t + 1 - T/2)) for m.
 627            sinh : Use anti-periodicitiy of the correlator by solving C(t) / C(t+1) = sinh(m * (t - T/2)) / sinh(m * (t + 1 - T/2)) for m.
 628            See, e.g., arXiv:1205.5380
 629            arccosh : Uses the explicit form of the symmetrized correlator (not recommended)
 630            logsym: uses the symmetric effective mass log(C(t-1) / C(t+1))/2
 631        guess : float
 632            guess for the root finder, only relevant for the root variant
 633        """
 634        if self.N != 1:
 635            raise Exception('Correlator must be projected before getting m_eff')
 636        if variant == 'log':
 637            newcontent = []
 638            for t in range(self.T - 1):
 639                if ((self.content[t] is None) or (self.content[t + 1] is None)) or (self.content[t + 1][0].value == 0):
 640                    newcontent.append(None)
 641                elif self.content[t][0].value / self.content[t + 1][0].value < 0:
 642                    newcontent.append(None)
 643                else:
 644                    newcontent.append(self.content[t] / self.content[t + 1])
 645            if (all([x is None for x in newcontent])):
 646                raise Exception('m_eff is undefined at all timeslices')
 647
 648            return np.log(Corr(newcontent, padding=[0, 1]))
 649
 650        elif variant == 'logsym':
 651            newcontent = []
 652            for t in range(1, self.T - 1):
 653                if ((self.content[t - 1] is None) or (self.content[t + 1] is None)) or (self.content[t + 1][0].value == 0):
 654                    newcontent.append(None)
 655                elif self.content[t - 1][0].value / self.content[t + 1][0].value < 0:
 656                    newcontent.append(None)
 657                else:
 658                    newcontent.append(self.content[t - 1] / self.content[t + 1])
 659            if (all([x is None for x in newcontent])):
 660                raise Exception('m_eff is undefined at all timeslices')
 661
 662            return np.log(Corr(newcontent, padding=[1, 1])) / 2
 663
 664        elif variant in ['periodic', 'cosh', 'sinh']:
 665            if variant in ['periodic', 'cosh']:
 666                func = anp.cosh
 667            else:
 668                func = anp.sinh
 669
 670            def root_function(x, d):
 671                return func(x * (t - self.T / 2)) / func(x * (t + 1 - self.T / 2)) - d
 672
 673            newcontent = []
 674            for t in range(self.T - 1):
 675                if (self.content[t] is None) or (self.content[t + 1] is None) or (self.content[t + 1][0].value == 0):
 676                    newcontent.append(None)
 677                # Fill the two timeslices in the middle of the lattice with their predecessors
 678                elif variant == 'sinh' and t in [self.T / 2, self.T / 2 - 1]:
 679                    newcontent.append(newcontent[-1])
 680                elif self.content[t][0].value / self.content[t + 1][0].value < 0:
 681                    newcontent.append(None)
 682                else:
 683                    newcontent.append(np.abs(find_root(self.content[t][0] / self.content[t + 1][0], root_function, guess=guess)))
 684            if (all([x is None for x in newcontent])):
 685                raise Exception('m_eff is undefined at all timeslices')
 686
 687            return Corr(newcontent, padding=[0, 1])
 688
 689        elif variant == 'arccosh':
 690            newcontent = []
 691            for t in range(1, self.T - 1):
 692                if (self.content[t] is None) or (self.content[t + 1] is None) or (self.content[t - 1] is None) or (self.content[t][0].value == 0):
 693                    newcontent.append(None)
 694                else:
 695                    newcontent.append((self.content[t + 1] + self.content[t - 1]) / (2 * self.content[t]))
 696            if (all([x is None for x in newcontent])):
 697                raise Exception("m_eff is undefined at all timeslices")
 698            return np.arccosh(Corr(newcontent, padding=[1, 1]))
 699
 700        else:
 701            raise Exception('Unknown variant.')
 702
 703    def fit(self, function, fitrange=None, silent=False, **kwargs):
 704        r'''Fits function to the data
 705
 706        Parameters
 707        ----------
 708        function : obj
 709            function to fit to the data. See fits.least_squares for details.
 710        fitrange : list
 711            Two element list containing the timeslices on which the fit is supposed to start and stop.
 712            Caution: This range is inclusive as opposed to standard python indexing.
 713            `fitrange=[4, 6]` corresponds to the three entries 4, 5 and 6.
 714            If not specified, self.prange or all timeslices are used.
 715        silent : bool
 716            Decides whether output is printed to the standard output.
 717        '''
 718        if self.N != 1:
 719            raise Exception("Correlator must be projected before fitting")
 720
 721        if fitrange is None:
 722            if self.prange:
 723                fitrange = self.prange
 724            else:
 725                fitrange = [0, self.T - 1]
 726        else:
 727            if not isinstance(fitrange, list):
 728                raise Exception("fitrange has to be a list with two elements")
 729            if len(fitrange) != 2:
 730                raise Exception("fitrange has to have exactly two elements [fit_start, fit_stop]")
 731
 732        xs = [x for x in range(fitrange[0], fitrange[1] + 1) if not self.content[x] is None]
 733        ys = [self.content[x][0] for x in range(fitrange[0], fitrange[1] + 1) if not self.content[x] is None]
 734        result = least_squares(xs, ys, function, silent=silent, **kwargs)
 735        return result
 736
 737    def plateau(self, plateau_range=None, method="fit", auto_gamma=False):
 738        """ Extract a plateau value from a Corr object
 739
 740        Parameters
 741        ----------
 742        plateau_range : list
 743            list with two entries, indicating the first and the last timeslice
 744            of the plateau region.
 745        method : str
 746            method to extract the plateau.
 747                'fit' fits a constant to the plateau region
 748                'avg', 'average' or 'mean' just average over the given timeslices.
 749        auto_gamma : bool
 750            apply gamma_method with default parameters to the Corr. Defaults to None
 751        """
 752        if not plateau_range:
 753            if self.prange:
 754                plateau_range = self.prange
 755            else:
 756                raise Exception("no plateau range provided")
 757        if self.N != 1:
 758            raise Exception("Correlator must be projected before getting a plateau.")
 759        if (all([self.content[t] is None for t in range(plateau_range[0], plateau_range[1] + 1)])):
 760            raise Exception("plateau is undefined at all timeslices in plateaurange.")
 761        if auto_gamma:
 762            self.gamma_method()
 763        if method == "fit":
 764            def const_func(a, t):
 765                return a[0]
 766            return self.fit(const_func, plateau_range)[0]
 767        elif method in ["avg", "average", "mean"]:
 768            returnvalue = np.mean([item[0] for item in self.content[plateau_range[0]:plateau_range[1] + 1] if item is not None])
 769            return returnvalue
 770
 771        else:
 772            raise Exception("Unsupported plateau method: " + method)
 773
 774    def set_prange(self, prange):
 775        """Sets the attribute prange of the Corr object."""
 776        if not len(prange) == 2:
 777            raise Exception("prange must be a list or array with two values")
 778        if not ((isinstance(prange[0], int)) and (isinstance(prange[1], int))):
 779            raise Exception("Start and end point must be integers")
 780        if not (0 <= prange[0] <= self.T and 0 <= prange[1] <= self.T and prange[0] < prange[1]):
 781            raise Exception("Start and end point must define a range in the interval 0,T")
 782
 783        self.prange = prange
 784        return
 785
 786    def show(self, x_range=None, comp=None, y_range=None, logscale=False, plateau=None, fit_res=None, ylabel=None, save=None, auto_gamma=False, hide_sigma=None, references=None, title=None):
 787        """Plots the correlator using the tag of the correlator as label if available.
 788
 789        Parameters
 790        ----------
 791        x_range : list
 792            list of two values, determining the range of the x-axis e.g. [4, 8].
 793        comp : Corr or list of Corr
 794            Correlator or list of correlators which are plotted for comparison.
 795            The tags of these correlators are used as labels if available.
 796        logscale : bool
 797            Sets y-axis to logscale.
 798        plateau : Obs
 799            Plateau value to be visualized in the figure.
 800        fit_res : Fit_result
 801            Fit_result object to be visualized.
 802        ylabel : str
 803            Label for the y-axis.
 804        save : str
 805            path to file in which the figure should be saved.
 806        auto_gamma : bool
 807            Apply the gamma method with standard parameters to all correlators and plateau values before plotting.
 808        hide_sigma : float
 809            Hides data points from the first value on which is consistent with zero within 'hide_sigma' standard errors.
 810        references : list
 811            List of floating point values that are displayed as horizontal lines for reference.
 812        title : string
 813            Optional title of the figure.
 814        """
 815        if self.N != 1:
 816            raise Exception("Correlator must be projected before plotting")
 817
 818        if auto_gamma:
 819            self.gamma_method()
 820
 821        if x_range is None:
 822            x_range = [0, self.T - 1]
 823
 824        fig = plt.figure()
 825        ax1 = fig.add_subplot(111)
 826
 827        x, y, y_err = self.plottable()
 828        if hide_sigma:
 829            hide_from = np.argmax((hide_sigma * np.array(y_err[1:])) > np.abs(y[1:])) - 1
 830        else:
 831            hide_from = None
 832        ax1.errorbar(x[:hide_from], y[:hide_from], y_err[:hide_from], label=self.tag)
 833        if logscale:
 834            ax1.set_yscale('log')
 835        else:
 836            if y_range is None:
 837                try:
 838                    y_min = min([(x[0].value - x[0].dvalue) for x in self.content[x_range[0]: x_range[1] + 1] if (x is not None) and x[0].dvalue < 2 * np.abs(x[0].value)])
 839                    y_max = max([(x[0].value + x[0].dvalue) for x in self.content[x_range[0]: x_range[1] + 1] if (x is not None) and x[0].dvalue < 2 * np.abs(x[0].value)])
 840                    ax1.set_ylim([y_min - 0.1 * (y_max - y_min), y_max + 0.1 * (y_max - y_min)])
 841                except Exception:
 842                    pass
 843            else:
 844                ax1.set_ylim(y_range)
 845        if comp:
 846            if isinstance(comp, (Corr, list)):
 847                for corr in comp if isinstance(comp, list) else [comp]:
 848                    if auto_gamma:
 849                        corr.gamma_method()
 850                    x, y, y_err = corr.plottable()
 851                    if hide_sigma:
 852                        hide_from = np.argmax((hide_sigma * np.array(y_err[1:])) > np.abs(y[1:])) - 1
 853                    else:
 854                        hide_from = None
 855                    ax1.errorbar(x[:hide_from], y[:hide_from], y_err[:hide_from], label=corr.tag, mfc=plt.rcParams['axes.facecolor'])
 856            else:
 857                raise Exception("'comp' must be a correlator or a list of correlators.")
 858
 859        if plateau:
 860            if isinstance(plateau, Obs):
 861                if auto_gamma:
 862                    plateau.gamma_method()
 863                ax1.axhline(y=plateau.value, linewidth=2, color=plt.rcParams['text.color'], alpha=0.6, marker=',', ls='--', label=str(plateau))
 864                ax1.axhspan(plateau.value - plateau.dvalue, plateau.value + plateau.dvalue, alpha=0.25, color=plt.rcParams['text.color'], ls='-')
 865            else:
 866                raise Exception("'plateau' must be an Obs")
 867
 868        if references:
 869            if isinstance(references, list):
 870                for ref in references:
 871                    ax1.axhline(y=ref, linewidth=1, color=plt.rcParams['text.color'], alpha=0.6, marker=',', ls='--')
 872            else:
 873                raise Exception("'references' must be a list of floating pint values.")
 874
 875        if self.prange:
 876            ax1.axvline(self.prange[0], 0, 1, ls='-', marker=',')
 877            ax1.axvline(self.prange[1], 0, 1, ls='-', marker=',')
 878
 879        if fit_res:
 880            x_samples = np.arange(x_range[0], x_range[1] + 1, 0.05)
 881            ax1.plot(x_samples,
 882                     fit_res.fit_function([o.value for o in fit_res.fit_parameters], x_samples),
 883                     ls='-', marker=',', lw=2)
 884
 885        ax1.set_xlabel(r'$x_0 / a$')
 886        if ylabel:
 887            ax1.set_ylabel(ylabel)
 888        ax1.set_xlim([x_range[0] - 0.5, x_range[1] + 0.5])
 889
 890        handles, labels = ax1.get_legend_handles_labels()
 891        if labels:
 892            ax1.legend()
 893
 894        if title:
 895            plt.title(title)
 896
 897        plt.draw()
 898
 899        if save:
 900            if isinstance(save, str):
 901                fig.savefig(save, bbox_inches='tight')
 902            else:
 903                raise Exception("'save' has to be a string.")
 904
 905    def spaghetti_plot(self, logscale=True):
 906        """Produces a spaghetti plot of the correlator suited to monitor exceptional configurations.
 907
 908        Parameters
 909        ----------
 910        logscale : bool
 911            Determines whether the scale of the y-axis is logarithmic or standard.
 912        """
 913        if self.N != 1:
 914            raise Exception("Correlator needs to be projected first.")
 915
 916        mc_names = list(set([item for sublist in [sum(map(o[0].e_content.get, o[0].mc_names), []) for o in self.content if o is not None] for item in sublist]))
 917        x0_vals = [n for (n, o) in zip(np.arange(self.T), self.content) if o is not None]
 918
 919        for name in mc_names:
 920            data = np.array([o[0].deltas[name] + o[0].r_values[name] for o in self.content if o is not None]).T
 921
 922            fig = plt.figure()
 923            ax = fig.add_subplot(111)
 924            for dat in data:
 925                ax.plot(x0_vals, dat, ls='-', marker='')
 926
 927            if logscale is True:
 928                ax.set_yscale('log')
 929
 930            ax.set_xlabel(r'$x_0 / a$')
 931            plt.title(name)
 932            plt.draw()
 933
 934    def dump(self, filename, datatype="json.gz", **kwargs):
 935        """Dumps the Corr into a file of chosen type
 936        Parameters
 937        ----------
 938        filename : str
 939            Name of the file to be saved.
 940        datatype : str
 941            Format of the exported file. Supported formats include
 942            "json.gz" and "pickle"
 943        path : str
 944            specifies a custom path for the file (default '.')
 945        """
 946        if datatype == "json.gz":
 947            from .input.json import dump_to_json
 948            if 'path' in kwargs:
 949                file_name = kwargs.get('path') + '/' + filename
 950            else:
 951                file_name = filename
 952            dump_to_json(self, file_name)
 953        elif datatype == "pickle":
 954            dump_object(self, filename, **kwargs)
 955        else:
 956            raise Exception("Unknown datatype " + str(datatype))
 957
 958    def print(self, print_range=None):
 959        print(self.__repr__(print_range))
 960
 961    def __repr__(self, print_range=None):
 962        if print_range is None:
 963            print_range = [0, None]
 964
 965        content_string = ""
 966        content_string += "Corr T=" + str(self.T) + " N=" + str(self.N) + "\n"  # +" filled with"+ str(type(self.content[0][0])) there should be a good solution here
 967
 968        if self.tag is not None:
 969            content_string += "Description: " + self.tag + "\n"
 970        if self.N != 1:
 971            return content_string
 972        if isinstance(self[0], CObs):
 973            return content_string
 974
 975        if print_range[1]:
 976            print_range[1] += 1
 977        content_string += 'x0/a\tCorr(x0/a)\n------------------\n'
 978        for i, sub_corr in enumerate(self.content[print_range[0]:print_range[1]]):
 979            if sub_corr is None:
 980                content_string += str(i + print_range[0]) + '\n'
 981            else:
 982                content_string += str(i + print_range[0])
 983                for element in sub_corr:
 984                    content_string += '\t' + ' ' * int(element >= 0) + str(element)
 985                content_string += '\n'
 986        return content_string
 987
 988    def __str__(self):
 989        return self.__repr__()
 990
 991    # We define the basic operations, that can be performed with correlators.
 992    # While */+- get defined here, they only work for Corr*Obs and not Obs*Corr.
 993    # This is because Obs*Corr checks Obs.__mul__ first and does not catch an exception.
 994    # One could try and tell Obs to check if the y in __mul__ is a Corr and
 995
 996    def __add__(self, y):
 997        if isinstance(y, Corr):
 998            if ((self.N != y.N) or (self.T != y.T)):
 999                raise Exception("Addition of Corrs with different shape")
1000            newcontent = []
1001            for t in range(self.T):
1002                if _check_for_none(self, self.content[t]) or _check_for_none(y, y.content[t]):
1003                    newcontent.append(None)
1004                else:
1005                    newcontent.append(self.content[t] + y.content[t])
1006            return Corr(newcontent)
1007
1008        elif isinstance(y, (Obs, int, float, CObs)):
1009            newcontent = []
1010            for t in range(self.T):
1011                if _check_for_none(self, self.content[t]):
1012                    newcontent.append(None)
1013                else:
1014                    newcontent.append(self.content[t] + y)
1015            return Corr(newcontent, prange=self.prange)
1016        elif isinstance(y, np.ndarray):
1017            if y.shape == (self.T,):
1018                return Corr(list((np.array(self.content).T + y).T))
1019            else:
1020                raise ValueError("operands could not be broadcast together")
1021        else:
1022            raise TypeError("Corr + wrong type")
1023
1024    def __mul__(self, y):
1025        if isinstance(y, Corr):
1026            if not ((self.N == 1 or y.N == 1 or self.N == y.N) and self.T == y.T):
1027                raise Exception("Multiplication of Corr object requires N=N or N=1 and T=T")
1028            newcontent = []
1029            for t in range(self.T):
1030                if _check_for_none(self, self.content[t]) or _check_for_none(y, y.content[t]):
1031                    newcontent.append(None)
1032                else:
1033                    newcontent.append(self.content[t] * y.content[t])
1034            return Corr(newcontent)
1035
1036        elif isinstance(y, (Obs, int, float, CObs)):
1037            newcontent = []
1038            for t in range(self.T):
1039                if _check_for_none(self, self.content[t]):
1040                    newcontent.append(None)
1041                else:
1042                    newcontent.append(self.content[t] * y)
1043            return Corr(newcontent, prange=self.prange)
1044        elif isinstance(y, np.ndarray):
1045            if y.shape == (self.T,):
1046                return Corr(list((np.array(self.content).T * y).T))
1047            else:
1048                raise ValueError("operands could not be broadcast together")
1049        else:
1050            raise TypeError("Corr * wrong type")
1051
1052    def __truediv__(self, y):
1053        if isinstance(y, Corr):
1054            if not ((self.N == 1 or y.N == 1 or self.N == y.N) and self.T == y.T):
1055                raise Exception("Multiplication of Corr object requires N=N or N=1 and T=T")
1056            newcontent = []
1057            for t in range(self.T):
1058                if _check_for_none(self, self.content[t]) or _check_for_none(y, y.content[t]):
1059                    newcontent.append(None)
1060                else:
1061                    newcontent.append(self.content[t] / y.content[t])
1062            for t in range(self.T):
1063                if _check_for_none(self, newcontent[t]):
1064                    continue
1065                if np.isnan(np.sum(newcontent[t]).value):
1066                    newcontent[t] = None
1067
1068            if all([item is None for item in newcontent]):
1069                raise Exception("Division returns completely undefined correlator")
1070            return Corr(newcontent)
1071
1072        elif isinstance(y, (Obs, CObs)):
1073            if isinstance(y, Obs):
1074                if y.value == 0:
1075                    raise Exception('Division by zero will return undefined correlator')
1076            if isinstance(y, CObs):
1077                if y.is_zero():
1078                    raise Exception('Division by zero will return undefined correlator')
1079
1080            newcontent = []
1081            for t in range(self.T):
1082                if _check_for_none(self, self.content[t]):
1083                    newcontent.append(None)
1084                else:
1085                    newcontent.append(self.content[t] / y)
1086            return Corr(newcontent, prange=self.prange)
1087
1088        elif isinstance(y, (int, float)):
1089            if y == 0:
1090                raise Exception('Division by zero will return undefined correlator')
1091            newcontent = []
1092            for t in range(self.T):
1093                if _check_for_none(self, self.content[t]):
1094                    newcontent.append(None)
1095                else:
1096                    newcontent.append(self.content[t] / y)
1097            return Corr(newcontent, prange=self.prange)
1098        elif isinstance(y, np.ndarray):
1099            if y.shape == (self.T,):
1100                return Corr(list((np.array(self.content).T / y).T))
1101            else:
1102                raise ValueError("operands could not be broadcast together")
1103        else:
1104            raise TypeError('Corr / wrong type')
1105
1106    def __neg__(self):
1107        newcontent = [None if _check_for_none(self, item) else -1. * item for item in self.content]
1108        return Corr(newcontent, prange=self.prange)
1109
1110    def __sub__(self, y):
1111        return self + (-y)
1112
1113    def __pow__(self, y):
1114        if isinstance(y, (Obs, int, float, CObs)):
1115            newcontent = [None if _check_for_none(self, item) else item**y for item in self.content]
1116            return Corr(newcontent, prange=self.prange)
1117        else:
1118            raise TypeError('Type of exponent not supported')
1119
1120    def __abs__(self):
1121        newcontent = [None if _check_for_none(self, item) else np.abs(item) for item in self.content]
1122        return Corr(newcontent, prange=self.prange)
1123
1124    # The numpy functions:
1125    def sqrt(self):
1126        return self ** 0.5
1127
1128    def log(self):
1129        newcontent = [None if _check_for_none(self, item) else np.log(item) for item in self.content]
1130        return Corr(newcontent, prange=self.prange)
1131
1132    def exp(self):
1133        newcontent = [None if _check_for_none(self, item) else np.exp(item) for item in self.content]
1134        return Corr(newcontent, prange=self.prange)
1135
1136    def _apply_func_to_corr(self, func):
1137        newcontent = [None if _check_for_none(self, item) else func(item) for item in self.content]
1138        for t in range(self.T):
1139            if _check_for_none(self, newcontent[t]):
1140                continue
1141            tmp_sum = np.sum(newcontent[t])
1142            if hasattr(tmp_sum, "value"):
1143                if np.isnan(tmp_sum.value):
1144                    newcontent[t] = None
1145        if all([item is None for item in newcontent]):
1146            raise Exception('Operation returns undefined correlator')
1147        return Corr(newcontent)
1148
1149    def sin(self):
1150        return self._apply_func_to_corr(np.sin)
1151
1152    def cos(self):
1153        return self._apply_func_to_corr(np.cos)
1154
1155    def tan(self):
1156        return self._apply_func_to_corr(np.tan)
1157
1158    def sinh(self):
1159        return self._apply_func_to_corr(np.sinh)
1160
1161    def cosh(self):
1162        return self._apply_func_to_corr(np.cosh)
1163
1164    def tanh(self):
1165        return self._apply_func_to_corr(np.tanh)
1166
1167    def arcsin(self):
1168        return self._apply_func_to_corr(np.arcsin)
1169
1170    def arccos(self):
1171        return self._apply_func_to_corr(np.arccos)
1172
1173    def arctan(self):
1174        return self._apply_func_to_corr(np.arctan)
1175
1176    def arcsinh(self):
1177        return self._apply_func_to_corr(np.arcsinh)
1178
1179    def arccosh(self):
1180        return self._apply_func_to_corr(np.arccosh)
1181
1182    def arctanh(self):
1183        return self._apply_func_to_corr(np.arctanh)
1184
1185    # Right hand side operations (require tweak in main module to work)
1186    def __radd__(self, y):
1187        return self + y
1188
1189    def __rsub__(self, y):
1190        return -self + y
1191
1192    def __rmul__(self, y):
1193        return self * y
1194
1195    def __rtruediv__(self, y):
1196        return (self / y) ** (-1)
1197
1198    @property
1199    def real(self):
1200        def return_real(obs_OR_cobs):
1201            if isinstance(obs_OR_cobs.flatten()[0], CObs):
1202                return np.vectorize(lambda x: x.real)(obs_OR_cobs)
1203            else:
1204                return obs_OR_cobs
1205
1206        return self._apply_func_to_corr(return_real)
1207
1208    @property
1209    def imag(self):
1210        def return_imag(obs_OR_cobs):
1211            if isinstance(obs_OR_cobs.flatten()[0], CObs):
1212                return np.vectorize(lambda x: x.imag)(obs_OR_cobs)
1213            else:
1214                return obs_OR_cobs * 0  # So it stays the right type
1215
1216        return self._apply_func_to_corr(return_imag)
1217
1218    def prune(self, Ntrunc, tproj=3, t0proj=2, basematrix=None):
1219        r''' Project large correlation matrix to lowest states
1220
1221        This method can be used to reduce the size of an (N x N) correlation matrix
1222        to (Ntrunc x Ntrunc) by solving a GEVP at very early times where the noise
1223        is still small.
1224
1225        Parameters
1226        ----------
1227        Ntrunc: int
1228            Rank of the target matrix.
1229        tproj: int
1230            Time where the eigenvectors are evaluated, corresponds to ts in the GEVP method.
1231            The default value is 3.
1232        t0proj: int
1233            Time where the correlation matrix is inverted. Choosing t0proj=1 is strongly
1234            discouraged for O(a) improved theories, since the correctness of the procedure
1235            cannot be granted in this case. The default value is 2.
1236        basematrix : Corr
1237            Correlation matrix that is used to determine the eigenvectors of the
1238            lowest states based on a GEVP. basematrix is taken to be the Corr itself if
1239            is is not specified.
1240
1241        Notes
1242        -----
1243        We have the basematrix $C(t)$ and the target matrix $G(t)$. We start by solving
1244        the GEVP $$C(t) v_n(t, t_0) = \lambda_n(t, t_0) C(t_0) v_n(t, t_0)$$ where $t \equiv t_\mathrm{proj}$
1245        and $t_0 \equiv t_{0, \mathrm{proj}}$. The target matrix is projected onto the subspace of the
1246        resulting eigenvectors $v_n, n=1,\dots,N_\mathrm{trunc}$ via
1247        $$G^\prime_{i, j}(t) = (v_i, G(t) v_j)$$. This allows to reduce the size of a large
1248        correlation matrix and to remove some noise that is added by irrelevant operators.
1249        This may allow to use the GEVP on $G(t)$ at late times such that the theoretically motivated
1250        bound $t_0 \leq t/2$ holds, since the condition number of $G(t)$ is decreased, compared to $C(t)$.
1251        '''
1252
1253        if self.N == 1:
1254            raise Exception('Method cannot be applied to one-dimensional correlators.')
1255        if basematrix is None:
1256            basematrix = self
1257        if Ntrunc >= basematrix.N:
1258            raise Exception('Cannot truncate using Ntrunc <= %d' % (basematrix.N))
1259        if basematrix.N != self.N:
1260            raise Exception('basematrix and targetmatrix have to be of the same size.')
1261
1262        evecs = basematrix.GEVP(t0proj, tproj, sort=None)[:Ntrunc]
1263
1264        tmpmat = np.empty((Ntrunc, Ntrunc), dtype=object)
1265        rmat = []
1266        for t in range(basematrix.T):
1267            for i in range(Ntrunc):
1268                for j in range(Ntrunc):
1269                    tmpmat[i][j] = evecs[i].T @ self[t] @ evecs[j]
1270            rmat.append(np.copy(tmpmat))
1271
1272        newcontent = [None if (self.content[t] is None) else rmat[t] for t in range(self.T)]
1273        return Corr(newcontent)
1274
1275
1276def _sort_vectors(vec_set, ts):
1277    """Helper function used to find a set of Eigenvectors consistent over all timeslices"""
1278    reference_sorting = np.array(vec_set[ts])
1279    N = reference_sorting.shape[0]
1280    sorted_vec_set = []
1281    for t in range(len(vec_set)):
1282        if vec_set[t] is None:
1283            sorted_vec_set.append(None)
1284        elif not t == ts:
1285            perms = [list(o) for o in permutations([i for i in range(N)], N)]
1286            best_score = 0
1287            for perm in perms:
1288                current_score = 1
1289                for k in range(N):
1290                    new_sorting = reference_sorting.copy()
1291                    new_sorting[perm[k], :] = vec_set[t][k]
1292                    current_score *= abs(np.linalg.det(new_sorting))
1293                if current_score > best_score:
1294                    best_score = current_score
1295                    best_perm = perm
1296            sorted_vec_set.append([vec_set[t][k] for k in best_perm])
1297        else:
1298            sorted_vec_set.append(vec_set[t])
1299
1300    return sorted_vec_set
1301
1302
1303def _check_for_none(corr, entry):
1304    """Checks if entry for correlator corr is None"""
1305    return len(list(filter(None, np.asarray(entry).flatten()))) < corr.N ** 2
1306
1307
1308def _GEVP_solver(Gt, G0):
1309    """Helper function for solving the GEVP and sorting the eigenvectors.
1310
1311    The helper function assumes that both provided matrices are symmetric and
1312    only processes the lower triangular part of both matrices. In case the matrices
1313    are not symmetric the upper triangular parts are effectively discarded."""
1314    return scipy.linalg.eigh(Gt, G0, lower=True)[1].T[::-1]
class Corr:
  14class Corr:
  15    """The class for a correlator (time dependent sequence of pe.Obs).
  16
  17    Everything, this class does, can be achieved using lists or arrays of Obs.
  18    But it is simply more convenient to have a dedicated object for correlators.
  19    One often wants to add or multiply correlators of the same length at every timeslice and it is inconvenient
  20    to iterate over all timeslices for every operation. This is especially true, when dealing with matrices.
  21
  22    The correlator can have two types of content: An Obs at every timeslice OR a GEVP
  23    matrix at every timeslice. Other dependency (eg. spatial) are not supported.
  24
  25    """
  26
  27    def __init__(self, data_input, padding=[0, 0], prange=None):
  28        """ Initialize a Corr object.
  29
  30        Parameters
  31        ----------
  32        data_input : list or array
  33            list of Obs or list of arrays of Obs or array of Corrs
  34        padding : list, optional
  35            List with two entries where the first labels the padding
  36            at the front of the correlator and the second the padding
  37            at the back.
  38        prange : list, optional
  39            List containing the first and last timeslice of the plateau
  40            region indentified for this correlator.
  41        """
  42
  43        if isinstance(data_input, np.ndarray):
  44
  45            # This only works, if the array fulfills the conditions below
  46            if not len(data_input.shape) == 2 and data_input.shape[0] == data_input.shape[1]:
  47                raise Exception("Incompatible array shape")
  48            if not all([isinstance(item, Corr) for item in data_input.flatten()]):
  49                raise Exception("If the input is an array, its elements must be of type pe.Corr")
  50            if not all([item.N == 1 for item in data_input.flatten()]):
  51                raise Exception("Can only construct matrix correlator from single valued correlators")
  52            if not len(set([item.T for item in data_input.flatten()])) == 1:
  53                raise Exception("All input Correlators must be defined over the same timeslices.")
  54
  55            T = data_input[0, 0].T
  56            N = data_input.shape[0]
  57            input_as_list = []
  58            for t in range(T):
  59                if any([(item.content[t] is None) for item in data_input.flatten()]):
  60                    if not all([(item.content[t] is None) for item in data_input.flatten()]):
  61                        warnings.warn("Input ill-defined at different timeslices. Conversion leads to data loss!", RuntimeWarning)
  62                    input_as_list.append(None)
  63                else:
  64                    array_at_timeslace = np.empty([N, N], dtype="object")
  65                    for i in range(N):
  66                        for j in range(N):
  67                            array_at_timeslace[i, j] = data_input[i, j][t]
  68                    input_as_list.append(array_at_timeslace)
  69            data_input = input_as_list
  70
  71        if isinstance(data_input, list):
  72
  73            if all([isinstance(item, (Obs, CObs)) or item is None for item in data_input]):
  74                _assert_equal_properties([o for o in data_input if o is not None])
  75                self.content = [np.asarray([item]) if item is not None else None for item in data_input]
  76                self.N = 1
  77
  78            elif all([isinstance(item, np.ndarray) or item is None for item in data_input]) and any([isinstance(item, np.ndarray) for item in data_input]):
  79                self.content = data_input
  80                noNull = [a for a in self.content if not (a is None)]  # To check if the matrices are correct for all undefined elements
  81                self.N = noNull[0].shape[0]
  82                if self.N > 1 and noNull[0].shape[0] != noNull[0].shape[1]:
  83                    raise Exception("Smearing matrices are not NxN")
  84                if (not all([item.shape == noNull[0].shape for item in noNull])):
  85                    raise Exception("Items in data_input are not of identical shape." + str(noNull))
  86            else:
  87                raise Exception("data_input contains item of wrong type")
  88        else:
  89            raise Exception("Data input was not given as list or correct array")
  90
  91        self.tag = None
  92
  93        # An undefined timeslice is represented by the None object
  94        self.content = [None] * padding[0] + self.content + [None] * padding[1]
  95        self.T = len(self.content)
  96        self.prange = prange
  97
  98    def __getitem__(self, idx):
  99        """Return the content of timeslice idx"""
 100        if self.content[idx] is None:
 101            return None
 102        elif len(self.content[idx]) == 1:
 103            return self.content[idx][0]
 104        else:
 105            return self.content[idx]
 106
 107    @property
 108    def reweighted(self):
 109        bool_array = np.array([list(map(lambda x: x.reweighted, o)) for o in [x for x in self.content if x is not None]])
 110        if np.all(bool_array == 1):
 111            return True
 112        elif np.all(bool_array == 0):
 113            return False
 114        else:
 115            raise Exception("Reweighting status of correlator corrupted.")
 116
 117    def gamma_method(self, **kwargs):
 118        """Apply the gamma method to the content of the Corr."""
 119        for item in self.content:
 120            if not (item is None):
 121                if self.N == 1:
 122                    item[0].gamma_method(**kwargs)
 123                else:
 124                    for i in range(self.N):
 125                        for j in range(self.N):
 126                            item[i, j].gamma_method(**kwargs)
 127
 128    def projected(self, vector_l=None, vector_r=None, normalize=False):
 129        """We need to project the Correlator with a Vector to get a single value at each timeslice.
 130
 131        The method can use one or two vectors.
 132        If two are specified it returns v1@G@v2 (the order might be very important.)
 133        By default it will return the lowest source, which usually means unsmeared-unsmeared (0,0), but it does not have to
 134        """
 135        if self.N == 1:
 136            raise Exception("Trying to project a Corr, that already has N=1.")
 137
 138        if vector_l is None:
 139            vector_l, vector_r = np.asarray([1.] + (self.N - 1) * [0.]), np.asarray([1.] + (self.N - 1) * [0.])
 140        elif (vector_r is None):
 141            vector_r = vector_l
 142        if isinstance(vector_l, list) and not isinstance(vector_r, list):
 143            if len(vector_l) != self.T:
 144                raise Exception("Length of vector list must be equal to T")
 145            vector_r = [vector_r] * self.T
 146        if isinstance(vector_r, list) and not isinstance(vector_l, list):
 147            if len(vector_r) != self.T:
 148                raise Exception("Length of vector list must be equal to T")
 149            vector_l = [vector_l] * self.T
 150
 151        if not isinstance(vector_l, list):
 152            if not vector_l.shape == vector_r.shape == (self.N,):
 153                raise Exception("Vectors are of wrong shape!")
 154            if normalize:
 155                vector_l, vector_r = vector_l / np.sqrt((vector_l @ vector_l)), vector_r / np.sqrt(vector_r @ vector_r)
 156            newcontent = [None if _check_for_none(self, item) else np.asarray([vector_l.T @ item @ vector_r]) for item in self.content]
 157
 158        else:
 159            # There are no checks here yet. There are so many possible scenarios, where this can go wrong.
 160            if normalize:
 161                for t in range(self.T):
 162                    vector_l[t], vector_r[t] = vector_l[t] / np.sqrt((vector_l[t] @ vector_l[t])), vector_r[t] / np.sqrt(vector_r[t] @ vector_r[t])
 163
 164            newcontent = [None if (_check_for_none(self, self.content[t]) or vector_l[t] is None or vector_r[t] is None) else np.asarray([vector_l[t].T @ self.content[t] @ vector_r[t]]) for t in range(self.T)]
 165        return Corr(newcontent)
 166
 167    def item(self, i, j):
 168        """Picks the element [i,j] from every matrix and returns a correlator containing one Obs per timeslice.
 169
 170        Parameters
 171        ----------
 172        i : int
 173            First index to be picked.
 174        j : int
 175            Second index to be picked.
 176        """
 177        if self.N == 1:
 178            raise Exception("Trying to pick item from projected Corr")
 179        newcontent = [None if (item is None) else item[i, j] for item in self.content]
 180        return Corr(newcontent)
 181
 182    def plottable(self):
 183        """Outputs the correlator in a plotable format.
 184
 185        Outputs three lists containing the timeslice index, the value on each
 186        timeslice and the error on each timeslice.
 187        """
 188        if self.N != 1:
 189            raise Exception("Can only make Corr[N=1] plottable")
 190        x_list = [x for x in range(self.T) if not self.content[x] is None]
 191        y_list = [y[0].value for y in self.content if y is not None]
 192        y_err_list = [y[0].dvalue for y in self.content if y is not None]
 193
 194        return x_list, y_list, y_err_list
 195
 196    def symmetric(self):
 197        """ Symmetrize the correlator around x0=0."""
 198        if self.N != 1:
 199            raise Exception('symmetric cannot be safely applied to multi-dimensional correlators.')
 200        if self.T % 2 != 0:
 201            raise Exception("Can not symmetrize odd T")
 202
 203        if np.argmax(np.abs(self.content)) != 0:
 204            warnings.warn("Correlator does not seem to be symmetric around x0=0.", RuntimeWarning)
 205
 206        newcontent = [self.content[0]]
 207        for t in range(1, self.T):
 208            if (self.content[t] is None) or (self.content[self.T - t] is None):
 209                newcontent.append(None)
 210            else:
 211                newcontent.append(0.5 * (self.content[t] + self.content[self.T - t]))
 212        if (all([x is None for x in newcontent])):
 213            raise Exception("Corr could not be symmetrized: No redundant values")
 214        return Corr(newcontent, prange=self.prange)
 215
 216    def anti_symmetric(self):
 217        """Anti-symmetrize the correlator around x0=0."""
 218        if self.N != 1:
 219            raise Exception('anti_symmetric cannot be safely applied to multi-dimensional correlators.')
 220        if self.T % 2 != 0:
 221            raise Exception("Can not symmetrize odd T")
 222
 223        test = 1 * self
 224        test.gamma_method()
 225        if not all([o.is_zero_within_error(3) for o in test.content[0]]):
 226            warnings.warn("Correlator does not seem to be anti-symmetric around x0=0.", RuntimeWarning)
 227
 228        newcontent = [self.content[0]]
 229        for t in range(1, self.T):
 230            if (self.content[t] is None) or (self.content[self.T - t] is None):
 231                newcontent.append(None)
 232            else:
 233                newcontent.append(0.5 * (self.content[t] - self.content[self.T - t]))
 234        if (all([x is None for x in newcontent])):
 235            raise Exception("Corr could not be symmetrized: No redundant values")
 236        return Corr(newcontent, prange=self.prange)
 237
 238    def is_matrix_symmetric(self):
 239        """Checks whether a correlator matrices is symmetric on every timeslice."""
 240        if self.N == 1:
 241            raise Exception("Only works for correlator matrices.")
 242        for t in range(self.T):
 243            if self[t] is None:
 244                continue
 245            for i in range(self.N):
 246                for j in range(i + 1, self.N):
 247                    if self[t][i, j] is self[t][j, i]:
 248                        continue
 249                    if hash(self[t][i, j]) != hash(self[t][j, i]):
 250                        return False
 251        return True
 252
 253    def matrix_symmetric(self):
 254        """Symmetrizes the correlator matrices on every timeslice."""
 255        if self.N == 1:
 256            raise Exception("Trying to symmetrize a correlator matrix, that already has N=1.")
 257        if self.is_matrix_symmetric():
 258            return 1.0 * self
 259        else:
 260            transposed = [None if _check_for_none(self, G) else G.T for G in self.content]
 261            return 0.5 * (Corr(transposed) + self)
 262
 263    def GEVP(self, t0, ts=None, sort="Eigenvalue", **kwargs):
 264        r'''Solve the generalized eigenvalue problem on the correlator matrix and returns the corresponding eigenvectors.
 265
 266        The eigenvectors are sorted according to the descending eigenvalues, the zeroth eigenvector(s) correspond to the
 267        largest eigenvalue(s). The eigenvector(s) for the individual states can be accessed via slicing
 268        ```python
 269        C.GEVP(t0=2)[0]  # Ground state vector(s)
 270        C.GEVP(t0=2)[:3]  # Vectors for the lowest three states
 271        ```
 272
 273        Parameters
 274        ----------
 275        t0 : int
 276            The time t0 for the right hand side of the GEVP according to $G(t)v_i=\lambda_i G(t_0)v_i$
 277        ts : int
 278            fixed time $G(t_s)v_i=\lambda_i G(t_0)v_i$ if sort=None.
 279            If sort="Eigenvector" it gives a reference point for the sorting method.
 280        sort : string
 281            If this argument is set, a list of self.T vectors per state is returned. If it is set to None, only one vector is returned.
 282            - "Eigenvalue": The eigenvector is chosen according to which eigenvalue it belongs individually on every timeslice.
 283            - "Eigenvector": Use the method described in arXiv:2004.10472 to find the set of v(t) belonging to the state.
 284              The reference state is identified by its eigenvalue at $t=t_s$.
 285
 286        Other Parameters
 287        ----------------
 288        state : int
 289           Returns only the vector(s) for a specified state. The lowest state is zero.
 290        '''
 291
 292        if self.N == 1:
 293            raise Exception("GEVP methods only works on correlator matrices and not single correlators.")
 294        if ts is not None:
 295            if (ts <= t0):
 296                raise Exception("ts has to be larger than t0.")
 297
 298        if "sorted_list" in kwargs:
 299            warnings.warn("Argument 'sorted_list' is deprecated, use 'sort' instead.", DeprecationWarning)
 300            sort = kwargs.get("sorted_list")
 301
 302        if self.is_matrix_symmetric():
 303            symmetric_corr = self
 304        else:
 305            symmetric_corr = self.matrix_symmetric()
 306
 307        G0 = np.vectorize(lambda x: x.value)(symmetric_corr[t0])
 308        np.linalg.cholesky(G0)  # Check if matrix G0 is positive-semidefinite.
 309
 310        if sort is None:
 311            if (ts is None):
 312                raise Exception("ts is required if sort=None.")
 313            if (self.content[t0] is None) or (self.content[ts] is None):
 314                raise Exception("Corr not defined at t0/ts.")
 315            Gt = np.vectorize(lambda x: x.value)(symmetric_corr[ts])
 316            reordered_vecs = _GEVP_solver(Gt, G0)
 317
 318        elif sort in ["Eigenvalue", "Eigenvector"]:
 319            if sort == "Eigenvalue" and ts is not None:
 320                warnings.warn("ts has no effect when sorting by eigenvalue is chosen.", RuntimeWarning)
 321            all_vecs = [None] * (t0 + 1)
 322            for t in range(t0 + 1, self.T):
 323                try:
 324                    Gt = np.vectorize(lambda x: x.value)(symmetric_corr[t])
 325                    all_vecs.append(_GEVP_solver(Gt, G0))
 326                except Exception:
 327                    all_vecs.append(None)
 328            if sort == "Eigenvector":
 329                if ts is None:
 330                    raise Exception("ts is required for the Eigenvector sorting method.")
 331                all_vecs = _sort_vectors(all_vecs, ts)
 332
 333            reordered_vecs = [[v[s] if v is not None else None for v in all_vecs] for s in range(self.N)]
 334        else:
 335            raise Exception("Unkown value for 'sort'.")
 336
 337        if "state" in kwargs:
 338            return reordered_vecs[kwargs.get("state")]
 339        else:
 340            return reordered_vecs
 341
 342    def Eigenvalue(self, t0, ts=None, state=0, sort="Eigenvalue"):
 343        """Determines the eigenvalue of the GEVP by solving and projecting the correlator
 344
 345        Parameters
 346        ----------
 347        state : int
 348            The state one is interested in ordered by energy. The lowest state is zero.
 349
 350        All other parameters are identical to the ones of Corr.GEVP.
 351        """
 352        vec = self.GEVP(t0, ts=ts, sort=sort)[state]
 353        return self.projected(vec)
 354
 355    def Hankel(self, N, periodic=False):
 356        """Constructs an NxN Hankel matrix
 357
 358        C(t) c(t+1) ... c(t+n-1)
 359        C(t+1) c(t+2) ... c(t+n)
 360        .................
 361        C(t+(n-1)) c(t+n) ... c(t+2(n-1))
 362
 363        Parameters
 364        ----------
 365        N : int
 366            Dimension of the Hankel matrix
 367        periodic : bool, optional
 368            determines whether the matrix is extended periodically
 369        """
 370
 371        if self.N != 1:
 372            raise Exception("Multi-operator Prony not implemented!")
 373
 374        array = np.empty([N, N], dtype="object")
 375        new_content = []
 376        for t in range(self.T):
 377            new_content.append(array.copy())
 378
 379        def wrap(i):
 380            while i >= self.T:
 381                i -= self.T
 382            return i
 383
 384        for t in range(self.T):
 385            for i in range(N):
 386                for j in range(N):
 387                    if periodic:
 388                        new_content[t][i, j] = self.content[wrap(t + i + j)][0]
 389                    elif (t + i + j) >= self.T:
 390                        new_content[t] = None
 391                    else:
 392                        new_content[t][i, j] = self.content[t + i + j][0]
 393
 394        return Corr(new_content)
 395
 396    def roll(self, dt):
 397        """Periodically shift the correlator by dt timeslices
 398
 399        Parameters
 400        ----------
 401        dt : int
 402            number of timeslices
 403        """
 404        return Corr(list(np.roll(np.array(self.content, dtype=object), dt)))
 405
 406    def reverse(self):
 407        """Reverse the time ordering of the Corr"""
 408        return Corr(self.content[:: -1])
 409
 410    def thin(self, spacing=2, offset=0):
 411        """Thin out a correlator to suppress correlations
 412
 413        Parameters
 414        ----------
 415        spacing : int
 416            Keep only every 'spacing'th entry of the correlator
 417        offset : int
 418            Offset the equal spacing
 419        """
 420        new_content = []
 421        for t in range(self.T):
 422            if (offset + t) % spacing != 0:
 423                new_content.append(None)
 424            else:
 425                new_content.append(self.content[t])
 426        return Corr(new_content)
 427
 428    def correlate(self, partner):
 429        """Correlate the correlator with another correlator or Obs
 430
 431        Parameters
 432        ----------
 433        partner : Obs or Corr
 434            partner to correlate the correlator with.
 435            Can either be an Obs which is correlated with all entries of the
 436            correlator or a Corr of same length.
 437        """
 438        if self.N != 1:
 439            raise Exception("Only one-dimensional correlators can be safely correlated.")
 440        new_content = []
 441        for x0, t_slice in enumerate(self.content):
 442            if _check_for_none(self, t_slice):
 443                new_content.append(None)
 444            else:
 445                if isinstance(partner, Corr):
 446                    if _check_for_none(partner, partner.content[x0]):
 447                        new_content.append(None)
 448                    else:
 449                        new_content.append(np.array([correlate(o, partner.content[x0][0]) for o in t_slice]))
 450                elif isinstance(partner, Obs):  # Should this include CObs?
 451                    new_content.append(np.array([correlate(o, partner) for o in t_slice]))
 452                else:
 453                    raise Exception("Can only correlate with an Obs or a Corr.")
 454
 455        return Corr(new_content)
 456
 457    def reweight(self, weight, **kwargs):
 458        """Reweight the correlator.
 459
 460        Parameters
 461        ----------
 462        weight : Obs
 463            Reweighting factor. An Observable that has to be defined on a superset of the
 464            configurations in obs[i].idl for all i.
 465        all_configs : bool
 466            if True, the reweighted observables are normalized by the average of
 467            the reweighting factor on all configurations in weight.idl and not
 468            on the configurations in obs[i].idl.
 469        """
 470        if self.N != 1:
 471            raise Exception("Reweighting only implemented for one-dimensional correlators.")
 472        new_content = []
 473        for t_slice in self.content:
 474            if _check_for_none(self, t_slice):
 475                new_content.append(None)
 476            else:
 477                new_content.append(np.array(reweight(weight, t_slice, **kwargs)))
 478        return Corr(new_content)
 479
 480    def T_symmetry(self, partner, parity=+1):
 481        """Return the time symmetry average of the correlator and its partner
 482
 483        Parameters
 484        ----------
 485        partner : Corr
 486            Time symmetry partner of the Corr
 487        partity : int
 488            Parity quantum number of the correlator, can be +1 or -1
 489        """
 490        if self.N != 1:
 491            raise Exception("T_symmetry only implemented for one-dimensional correlators.")
 492        if not isinstance(partner, Corr):
 493            raise Exception("T partner has to be a Corr object.")
 494        if parity not in [+1, -1]:
 495            raise Exception("Parity has to be +1 or -1.")
 496        T_partner = parity * partner.reverse()
 497
 498        t_slices = []
 499        test = (self - T_partner)
 500        test.gamma_method()
 501        for x0, t_slice in enumerate(test.content):
 502            if t_slice is not None:
 503                if not t_slice[0].is_zero_within_error(5):
 504                    t_slices.append(x0)
 505        if t_slices:
 506            warnings.warn("T symmetry partners do not agree within 5 sigma on time slices " + str(t_slices) + ".", RuntimeWarning)
 507
 508        return (self + T_partner) / 2
 509
 510    def deriv(self, variant="symmetric"):
 511        """Return the first derivative of the correlator with respect to x0.
 512
 513        Parameters
 514        ----------
 515        variant : str
 516            decides which definition of the finite differences derivative is used.
 517            Available choice: symmetric, forward, backward, improved, log, default: symmetric
 518        """
 519        if self.N != 1:
 520            raise Exception("deriv only implemented for one-dimensional correlators.")
 521        if variant == "symmetric":
 522            newcontent = []
 523            for t in range(1, self.T - 1):
 524                if (self.content[t - 1] is None) or (self.content[t + 1] is None):
 525                    newcontent.append(None)
 526                else:
 527                    newcontent.append(0.5 * (self.content[t + 1] - self.content[t - 1]))
 528            if (all([x is None for x in newcontent])):
 529                raise Exception('Derivative is undefined at all timeslices')
 530            return Corr(newcontent, padding=[1, 1])
 531        elif variant == "forward":
 532            newcontent = []
 533            for t in range(self.T - 1):
 534                if (self.content[t] is None) or (self.content[t + 1] is None):
 535                    newcontent.append(None)
 536                else:
 537                    newcontent.append(self.content[t + 1] - self.content[t])
 538            if (all([x is None for x in newcontent])):
 539                raise Exception("Derivative is undefined at all timeslices")
 540            return Corr(newcontent, padding=[0, 1])
 541        elif variant == "backward":
 542            newcontent = []
 543            for t in range(1, self.T):
 544                if (self.content[t - 1] is None) or (self.content[t] is None):
 545                    newcontent.append(None)
 546                else:
 547                    newcontent.append(self.content[t] - self.content[t - 1])
 548            if (all([x is None for x in newcontent])):
 549                raise Exception("Derivative is undefined at all timeslices")
 550            return Corr(newcontent, padding=[1, 0])
 551        elif variant == "improved":
 552            newcontent = []
 553            for t in range(2, self.T - 2):
 554                if (self.content[t - 2] is None) or (self.content[t - 1] is None) or (self.content[t + 1] is None) or (self.content[t + 2] is None):
 555                    newcontent.append(None)
 556                else:
 557                    newcontent.append((1 / 12) * (self.content[t - 2] - 8 * self.content[t - 1] + 8 * self.content[t + 1] - self.content[t + 2]))
 558            if (all([x is None for x in newcontent])):
 559                raise Exception('Derivative is undefined at all timeslices')
 560            return Corr(newcontent, padding=[2, 2])
 561        elif variant == 'log':
 562            newcontent = []
 563            for t in range(self.T):
 564                if (self.content[t] is None) or (self.content[t] <= 0):
 565                    newcontent.append(None)
 566                else:
 567                    newcontent.append(np.log(self.content[t]))
 568            if (all([x is None for x in newcontent])):
 569                raise Exception("Log is undefined at all timeslices")
 570            logcorr = Corr(newcontent)
 571            return self * logcorr.deriv('symmetric')
 572        else:
 573            raise Exception("Unknown variant.")
 574
 575    def second_deriv(self, variant="symmetric"):
 576        """Return the second derivative of the correlator with respect to x0.
 577
 578        Parameters
 579        ----------
 580        variant : str
 581            decides which definition of the finite differences derivative is used.
 582            Available choice: symmetric, improved, log, default: symmetric
 583        """
 584        if self.N != 1:
 585            raise Exception("second_deriv only implemented for one-dimensional correlators.")
 586        if variant == "symmetric":
 587            newcontent = []
 588            for t in range(1, self.T - 1):
 589                if (self.content[t - 1] is None) or (self.content[t + 1] is None):
 590                    newcontent.append(None)
 591                else:
 592                    newcontent.append((self.content[t + 1] - 2 * self.content[t] + self.content[t - 1]))
 593            if (all([x is None for x in newcontent])):
 594                raise Exception("Derivative is undefined at all timeslices")
 595            return Corr(newcontent, padding=[1, 1])
 596        elif variant == "improved":
 597            newcontent = []
 598            for t in range(2, self.T - 2):
 599                if (self.content[t - 2] is None) or (self.content[t - 1] is None) or (self.content[t] is None) or (self.content[t + 1] is None) or (self.content[t + 2] is None):
 600                    newcontent.append(None)
 601                else:
 602                    newcontent.append((1 / 12) * (-self.content[t + 2] + 16 * self.content[t + 1] - 30 * self.content[t] + 16 * self.content[t - 1] - self.content[t - 2]))
 603            if (all([x is None for x in newcontent])):
 604                raise Exception("Derivative is undefined at all timeslices")
 605            return Corr(newcontent, padding=[2, 2])
 606        elif variant == 'log':
 607            newcontent = []
 608            for t in range(self.T):
 609                if (self.content[t] is None) or (self.content[t] <= 0):
 610                    newcontent.append(None)
 611                else:
 612                    newcontent.append(np.log(self.content[t]))
 613            if (all([x is None for x in newcontent])):
 614                raise Exception("Log is undefined at all timeslices")
 615            logcorr = Corr(newcontent)
 616            return self * (logcorr.second_deriv('symmetric') + (logcorr.deriv('symmetric'))**2)
 617        else:
 618            raise Exception("Unknown variant.")
 619
 620    def m_eff(self, variant='log', guess=1.0):
 621        """Returns the effective mass of the correlator as correlator object
 622
 623        Parameters
 624        ----------
 625        variant : str
 626            log : uses the standard effective mass log(C(t) / C(t+1))
 627            cosh, periodic : Use periodicitiy of the correlator by solving C(t) / C(t+1) = cosh(m * (t - T/2)) / cosh(m * (t + 1 - T/2)) for m.
 628            sinh : Use anti-periodicitiy of the correlator by solving C(t) / C(t+1) = sinh(m * (t - T/2)) / sinh(m * (t + 1 - T/2)) for m.
 629            See, e.g., arXiv:1205.5380
 630            arccosh : Uses the explicit form of the symmetrized correlator (not recommended)
 631            logsym: uses the symmetric effective mass log(C(t-1) / C(t+1))/2
 632        guess : float
 633            guess for the root finder, only relevant for the root variant
 634        """
 635        if self.N != 1:
 636            raise Exception('Correlator must be projected before getting m_eff')
 637        if variant == 'log':
 638            newcontent = []
 639            for t in range(self.T - 1):
 640                if ((self.content[t] is None) or (self.content[t + 1] is None)) or (self.content[t + 1][0].value == 0):
 641                    newcontent.append(None)
 642                elif self.content[t][0].value / self.content[t + 1][0].value < 0:
 643                    newcontent.append(None)
 644                else:
 645                    newcontent.append(self.content[t] / self.content[t + 1])
 646            if (all([x is None for x in newcontent])):
 647                raise Exception('m_eff is undefined at all timeslices')
 648
 649            return np.log(Corr(newcontent, padding=[0, 1]))
 650
 651        elif variant == 'logsym':
 652            newcontent = []
 653            for t in range(1, self.T - 1):
 654                if ((self.content[t - 1] is None) or (self.content[t + 1] is None)) or (self.content[t + 1][0].value == 0):
 655                    newcontent.append(None)
 656                elif self.content[t - 1][0].value / self.content[t + 1][0].value < 0:
 657                    newcontent.append(None)
 658                else:
 659                    newcontent.append(self.content[t - 1] / self.content[t + 1])
 660            if (all([x is None for x in newcontent])):
 661                raise Exception('m_eff is undefined at all timeslices')
 662
 663            return np.log(Corr(newcontent, padding=[1, 1])) / 2
 664
 665        elif variant in ['periodic', 'cosh', 'sinh']:
 666            if variant in ['periodic', 'cosh']:
 667                func = anp.cosh
 668            else:
 669                func = anp.sinh
 670
 671            def root_function(x, d):
 672                return func(x * (t - self.T / 2)) / func(x * (t + 1 - self.T / 2)) - d
 673
 674            newcontent = []
 675            for t in range(self.T - 1):
 676                if (self.content[t] is None) or (self.content[t + 1] is None) or (self.content[t + 1][0].value == 0):
 677                    newcontent.append(None)
 678                # Fill the two timeslices in the middle of the lattice with their predecessors
 679                elif variant == 'sinh' and t in [self.T / 2, self.T / 2 - 1]:
 680                    newcontent.append(newcontent[-1])
 681                elif self.content[t][0].value / self.content[t + 1][0].value < 0:
 682                    newcontent.append(None)
 683                else:
 684                    newcontent.append(np.abs(find_root(self.content[t][0] / self.content[t + 1][0], root_function, guess=guess)))
 685            if (all([x is None for x in newcontent])):
 686                raise Exception('m_eff is undefined at all timeslices')
 687
 688            return Corr(newcontent, padding=[0, 1])
 689
 690        elif variant == 'arccosh':
 691            newcontent = []
 692            for t in range(1, self.T - 1):
 693                if (self.content[t] is None) or (self.content[t + 1] is None) or (self.content[t - 1] is None) or (self.content[t][0].value == 0):
 694                    newcontent.append(None)
 695                else:
 696                    newcontent.append((self.content[t + 1] + self.content[t - 1]) / (2 * self.content[t]))
 697            if (all([x is None for x in newcontent])):
 698                raise Exception("m_eff is undefined at all timeslices")
 699            return np.arccosh(Corr(newcontent, padding=[1, 1]))
 700
 701        else:
 702            raise Exception('Unknown variant.')
 703
 704    def fit(self, function, fitrange=None, silent=False, **kwargs):
 705        r'''Fits function to the data
 706
 707        Parameters
 708        ----------
 709        function : obj
 710            function to fit to the data. See fits.least_squares for details.
 711        fitrange : list
 712            Two element list containing the timeslices on which the fit is supposed to start and stop.
 713            Caution: This range is inclusive as opposed to standard python indexing.
 714            `fitrange=[4, 6]` corresponds to the three entries 4, 5 and 6.
 715            If not specified, self.prange or all timeslices are used.
 716        silent : bool
 717            Decides whether output is printed to the standard output.
 718        '''
 719        if self.N != 1:
 720            raise Exception("Correlator must be projected before fitting")
 721
 722        if fitrange is None:
 723            if self.prange:
 724                fitrange = self.prange
 725            else:
 726                fitrange = [0, self.T - 1]
 727        else:
 728            if not isinstance(fitrange, list):
 729                raise Exception("fitrange has to be a list with two elements")
 730            if len(fitrange) != 2:
 731                raise Exception("fitrange has to have exactly two elements [fit_start, fit_stop]")
 732
 733        xs = [x for x in range(fitrange[0], fitrange[1] + 1) if not self.content[x] is None]
 734        ys = [self.content[x][0] for x in range(fitrange[0], fitrange[1] + 1) if not self.content[x] is None]
 735        result = least_squares(xs, ys, function, silent=silent, **kwargs)
 736        return result
 737
 738    def plateau(self, plateau_range=None, method="fit", auto_gamma=False):
 739        """ Extract a plateau value from a Corr object
 740
 741        Parameters
 742        ----------
 743        plateau_range : list
 744            list with two entries, indicating the first and the last timeslice
 745            of the plateau region.
 746        method : str
 747            method to extract the plateau.
 748                'fit' fits a constant to the plateau region
 749                'avg', 'average' or 'mean' just average over the given timeslices.
 750        auto_gamma : bool
 751            apply gamma_method with default parameters to the Corr. Defaults to None
 752        """
 753        if not plateau_range:
 754            if self.prange:
 755                plateau_range = self.prange
 756            else:
 757                raise Exception("no plateau range provided")
 758        if self.N != 1:
 759            raise Exception("Correlator must be projected before getting a plateau.")
 760        if (all([self.content[t] is None for t in range(plateau_range[0], plateau_range[1] + 1)])):
 761            raise Exception("plateau is undefined at all timeslices in plateaurange.")
 762        if auto_gamma:
 763            self.gamma_method()
 764        if method == "fit":
 765            def const_func(a, t):
 766                return a[0]
 767            return self.fit(const_func, plateau_range)[0]
 768        elif method in ["avg", "average", "mean"]:
 769            returnvalue = np.mean([item[0] for item in self.content[plateau_range[0]:plateau_range[1] + 1] if item is not None])
 770            return returnvalue
 771
 772        else:
 773            raise Exception("Unsupported plateau method: " + method)
 774
 775    def set_prange(self, prange):
 776        """Sets the attribute prange of the Corr object."""
 777        if not len(prange) == 2:
 778            raise Exception("prange must be a list or array with two values")
 779        if not ((isinstance(prange[0], int)) and (isinstance(prange[1], int))):
 780            raise Exception("Start and end point must be integers")
 781        if not (0 <= prange[0] <= self.T and 0 <= prange[1] <= self.T and prange[0] < prange[1]):
 782            raise Exception("Start and end point must define a range in the interval 0,T")
 783
 784        self.prange = prange
 785        return
 786
 787    def show(self, x_range=None, comp=None, y_range=None, logscale=False, plateau=None, fit_res=None, ylabel=None, save=None, auto_gamma=False, hide_sigma=None, references=None, title=None):
 788        """Plots the correlator using the tag of the correlator as label if available.
 789
 790        Parameters
 791        ----------
 792        x_range : list
 793            list of two values, determining the range of the x-axis e.g. [4, 8].
 794        comp : Corr or list of Corr
 795            Correlator or list of correlators which are plotted for comparison.
 796            The tags of these correlators are used as labels if available.
 797        logscale : bool
 798            Sets y-axis to logscale.
 799        plateau : Obs
 800            Plateau value to be visualized in the figure.
 801        fit_res : Fit_result
 802            Fit_result object to be visualized.
 803        ylabel : str
 804            Label for the y-axis.
 805        save : str
 806            path to file in which the figure should be saved.
 807        auto_gamma : bool
 808            Apply the gamma method with standard parameters to all correlators and plateau values before plotting.
 809        hide_sigma : float
 810            Hides data points from the first value on which is consistent with zero within 'hide_sigma' standard errors.
 811        references : list
 812            List of floating point values that are displayed as horizontal lines for reference.
 813        title : string
 814            Optional title of the figure.
 815        """
 816        if self.N != 1:
 817            raise Exception("Correlator must be projected before plotting")
 818
 819        if auto_gamma:
 820            self.gamma_method()
 821
 822        if x_range is None:
 823            x_range = [0, self.T - 1]
 824
 825        fig = plt.figure()
 826        ax1 = fig.add_subplot(111)
 827
 828        x, y, y_err = self.plottable()
 829        if hide_sigma:
 830            hide_from = np.argmax((hide_sigma * np.array(y_err[1:])) > np.abs(y[1:])) - 1
 831        else:
 832            hide_from = None
 833        ax1.errorbar(x[:hide_from], y[:hide_from], y_err[:hide_from], label=self.tag)
 834        if logscale:
 835            ax1.set_yscale('log')
 836        else:
 837            if y_range is None:
 838                try:
 839                    y_min = min([(x[0].value - x[0].dvalue) for x in self.content[x_range[0]: x_range[1] + 1] if (x is not None) and x[0].dvalue < 2 * np.abs(x[0].value)])
 840                    y_max = max([(x[0].value + x[0].dvalue) for x in self.content[x_range[0]: x_range[1] + 1] if (x is not None) and x[0].dvalue < 2 * np.abs(x[0].value)])
 841                    ax1.set_ylim([y_min - 0.1 * (y_max - y_min), y_max + 0.1 * (y_max - y_min)])
 842                except Exception:
 843                    pass
 844            else:
 845                ax1.set_ylim(y_range)
 846        if comp:
 847            if isinstance(comp, (Corr, list)):
 848                for corr in comp if isinstance(comp, list) else [comp]:
 849                    if auto_gamma:
 850                        corr.gamma_method()
 851                    x, y, y_err = corr.plottable()
 852                    if hide_sigma:
 853                        hide_from = np.argmax((hide_sigma * np.array(y_err[1:])) > np.abs(y[1:])) - 1
 854                    else:
 855                        hide_from = None
 856                    ax1.errorbar(x[:hide_from], y[:hide_from], y_err[:hide_from], label=corr.tag, mfc=plt.rcParams['axes.facecolor'])
 857            else:
 858                raise Exception("'comp' must be a correlator or a list of correlators.")
 859
 860        if plateau:
 861            if isinstance(plateau, Obs):
 862                if auto_gamma:
 863                    plateau.gamma_method()
 864                ax1.axhline(y=plateau.value, linewidth=2, color=plt.rcParams['text.color'], alpha=0.6, marker=',', ls='--', label=str(plateau))
 865                ax1.axhspan(plateau.value - plateau.dvalue, plateau.value + plateau.dvalue, alpha=0.25, color=plt.rcParams['text.color'], ls='-')
 866            else:
 867                raise Exception("'plateau' must be an Obs")
 868
 869        if references:
 870            if isinstance(references, list):
 871                for ref in references:
 872                    ax1.axhline(y=ref, linewidth=1, color=plt.rcParams['text.color'], alpha=0.6, marker=',', ls='--')
 873            else:
 874                raise Exception("'references' must be a list of floating pint values.")
 875
 876        if self.prange:
 877            ax1.axvline(self.prange[0], 0, 1, ls='-', marker=',')
 878            ax1.axvline(self.prange[1], 0, 1, ls='-', marker=',')
 879
 880        if fit_res:
 881            x_samples = np.arange(x_range[0], x_range[1] + 1, 0.05)
 882            ax1.plot(x_samples,
 883                     fit_res.fit_function([o.value for o in fit_res.fit_parameters], x_samples),
 884                     ls='-', marker=',', lw=2)
 885
 886        ax1.set_xlabel(r'$x_0 / a$')
 887        if ylabel:
 888            ax1.set_ylabel(ylabel)
 889        ax1.set_xlim([x_range[0] - 0.5, x_range[1] + 0.5])
 890
 891        handles, labels = ax1.get_legend_handles_labels()
 892        if labels:
 893            ax1.legend()
 894
 895        if title:
 896            plt.title(title)
 897
 898        plt.draw()
 899
 900        if save:
 901            if isinstance(save, str):
 902                fig.savefig(save, bbox_inches='tight')
 903            else:
 904                raise Exception("'save' has to be a string.")
 905
 906    def spaghetti_plot(self, logscale=True):
 907        """Produces a spaghetti plot of the correlator suited to monitor exceptional configurations.
 908
 909        Parameters
 910        ----------
 911        logscale : bool
 912            Determines whether the scale of the y-axis is logarithmic or standard.
 913        """
 914        if self.N != 1:
 915            raise Exception("Correlator needs to be projected first.")
 916
 917        mc_names = list(set([item for sublist in [sum(map(o[0].e_content.get, o[0].mc_names), []) for o in self.content if o is not None] for item in sublist]))
 918        x0_vals = [n for (n, o) in zip(np.arange(self.T), self.content) if o is not None]
 919
 920        for name in mc_names:
 921            data = np.array([o[0].deltas[name] + o[0].r_values[name] for o in self.content if o is not None]).T
 922
 923            fig = plt.figure()
 924            ax = fig.add_subplot(111)
 925            for dat in data:
 926                ax.plot(x0_vals, dat, ls='-', marker='')
 927
 928            if logscale is True:
 929                ax.set_yscale('log')
 930
 931            ax.set_xlabel(r'$x_0 / a$')
 932            plt.title(name)
 933            plt.draw()
 934
 935    def dump(self, filename, datatype="json.gz", **kwargs):
 936        """Dumps the Corr into a file of chosen type
 937        Parameters
 938        ----------
 939        filename : str
 940            Name of the file to be saved.
 941        datatype : str
 942            Format of the exported file. Supported formats include
 943            "json.gz" and "pickle"
 944        path : str
 945            specifies a custom path for the file (default '.')
 946        """
 947        if datatype == "json.gz":
 948            from .input.json import dump_to_json
 949            if 'path' in kwargs:
 950                file_name = kwargs.get('path') + '/' + filename
 951            else:
 952                file_name = filename
 953            dump_to_json(self, file_name)
 954        elif datatype == "pickle":
 955            dump_object(self, filename, **kwargs)
 956        else:
 957            raise Exception("Unknown datatype " + str(datatype))
 958
 959    def print(self, print_range=None):
 960        print(self.__repr__(print_range))
 961
 962    def __repr__(self, print_range=None):
 963        if print_range is None:
 964            print_range = [0, None]
 965
 966        content_string = ""
 967        content_string += "Corr T=" + str(self.T) + " N=" + str(self.N) + "\n"  # +" filled with"+ str(type(self.content[0][0])) there should be a good solution here
 968
 969        if self.tag is not None:
 970            content_string += "Description: " + self.tag + "\n"
 971        if self.N != 1:
 972            return content_string
 973        if isinstance(self[0], CObs):
 974            return content_string
 975
 976        if print_range[1]:
 977            print_range[1] += 1
 978        content_string += 'x0/a\tCorr(x0/a)\n------------------\n'
 979        for i, sub_corr in enumerate(self.content[print_range[0]:print_range[1]]):
 980            if sub_corr is None:
 981                content_string += str(i + print_range[0]) + '\n'
 982            else:
 983                content_string += str(i + print_range[0])
 984                for element in sub_corr:
 985                    content_string += '\t' + ' ' * int(element >= 0) + str(element)
 986                content_string += '\n'
 987        return content_string
 988
 989    def __str__(self):
 990        return self.__repr__()
 991
 992    # We define the basic operations, that can be performed with correlators.
 993    # While */+- get defined here, they only work for Corr*Obs and not Obs*Corr.
 994    # This is because Obs*Corr checks Obs.__mul__ first and does not catch an exception.
 995    # One could try and tell Obs to check if the y in __mul__ is a Corr and
 996
 997    def __add__(self, y):
 998        if isinstance(y, Corr):
 999            if ((self.N != y.N) or (self.T != y.T)):
1000                raise Exception("Addition of Corrs with different shape")
1001            newcontent = []
1002            for t in range(self.T):
1003                if _check_for_none(self, self.content[t]) or _check_for_none(y, y.content[t]):
1004                    newcontent.append(None)
1005                else:
1006                    newcontent.append(self.content[t] + y.content[t])
1007            return Corr(newcontent)
1008
1009        elif isinstance(y, (Obs, int, float, CObs)):
1010            newcontent = []
1011            for t in range(self.T):
1012                if _check_for_none(self, self.content[t]):
1013                    newcontent.append(None)
1014                else:
1015                    newcontent.append(self.content[t] + y)
1016            return Corr(newcontent, prange=self.prange)
1017        elif isinstance(y, np.ndarray):
1018            if y.shape == (self.T,):
1019                return Corr(list((np.array(self.content).T + y).T))
1020            else:
1021                raise ValueError("operands could not be broadcast together")
1022        else:
1023            raise TypeError("Corr + wrong type")
1024
1025    def __mul__(self, y):
1026        if isinstance(y, Corr):
1027            if not ((self.N == 1 or y.N == 1 or self.N == y.N) and self.T == y.T):
1028                raise Exception("Multiplication of Corr object requires N=N or N=1 and T=T")
1029            newcontent = []
1030            for t in range(self.T):
1031                if _check_for_none(self, self.content[t]) or _check_for_none(y, y.content[t]):
1032                    newcontent.append(None)
1033                else:
1034                    newcontent.append(self.content[t] * y.content[t])
1035            return Corr(newcontent)
1036
1037        elif isinstance(y, (Obs, int, float, CObs)):
1038            newcontent = []
1039            for t in range(self.T):
1040                if _check_for_none(self, self.content[t]):
1041                    newcontent.append(None)
1042                else:
1043                    newcontent.append(self.content[t] * y)
1044            return Corr(newcontent, prange=self.prange)
1045        elif isinstance(y, np.ndarray):
1046            if y.shape == (self.T,):
1047                return Corr(list((np.array(self.content).T * y).T))
1048            else:
1049                raise ValueError("operands could not be broadcast together")
1050        else:
1051            raise TypeError("Corr * wrong type")
1052
1053    def __truediv__(self, y):
1054        if isinstance(y, Corr):
1055            if not ((self.N == 1 or y.N == 1 or self.N == y.N) and self.T == y.T):
1056                raise Exception("Multiplication of Corr object requires N=N or N=1 and T=T")
1057            newcontent = []
1058            for t in range(self.T):
1059                if _check_for_none(self, self.content[t]) or _check_for_none(y, y.content[t]):
1060                    newcontent.append(None)
1061                else:
1062                    newcontent.append(self.content[t] / y.content[t])
1063            for t in range(self.T):
1064                if _check_for_none(self, newcontent[t]):
1065                    continue
1066                if np.isnan(np.sum(newcontent[t]).value):
1067                    newcontent[t] = None
1068
1069            if all([item is None for item in newcontent]):
1070                raise Exception("Division returns completely undefined correlator")
1071            return Corr(newcontent)
1072
1073        elif isinstance(y, (Obs, CObs)):
1074            if isinstance(y, Obs):
1075                if y.value == 0:
1076                    raise Exception('Division by zero will return undefined correlator')
1077            if isinstance(y, CObs):
1078                if y.is_zero():
1079                    raise Exception('Division by zero will return undefined correlator')
1080
1081            newcontent = []
1082            for t in range(self.T):
1083                if _check_for_none(self, self.content[t]):
1084                    newcontent.append(None)
1085                else:
1086                    newcontent.append(self.content[t] / y)
1087            return Corr(newcontent, prange=self.prange)
1088
1089        elif isinstance(y, (int, float)):
1090            if y == 0:
1091                raise Exception('Division by zero will return undefined correlator')
1092            newcontent = []
1093            for t in range(self.T):
1094                if _check_for_none(self, self.content[t]):
1095                    newcontent.append(None)
1096                else:
1097                    newcontent.append(self.content[t] / y)
1098            return Corr(newcontent, prange=self.prange)
1099        elif isinstance(y, np.ndarray):
1100            if y.shape == (self.T,):
1101                return Corr(list((np.array(self.content).T / y).T))
1102            else:
1103                raise ValueError("operands could not be broadcast together")
1104        else:
1105            raise TypeError('Corr / wrong type')
1106
1107    def __neg__(self):
1108        newcontent = [None if _check_for_none(self, item) else -1. * item for item in self.content]
1109        return Corr(newcontent, prange=self.prange)
1110
1111    def __sub__(self, y):
1112        return self + (-y)
1113
1114    def __pow__(self, y):
1115        if isinstance(y, (Obs, int, float, CObs)):
1116            newcontent = [None if _check_for_none(self, item) else item**y for item in self.content]
1117            return Corr(newcontent, prange=self.prange)
1118        else:
1119            raise TypeError('Type of exponent not supported')
1120
1121    def __abs__(self):
1122        newcontent = [None if _check_for_none(self, item) else np.abs(item) for item in self.content]
1123        return Corr(newcontent, prange=self.prange)
1124
1125    # The numpy functions:
1126    def sqrt(self):
1127        return self ** 0.5
1128
1129    def log(self):
1130        newcontent = [None if _check_for_none(self, item) else np.log(item) for item in self.content]
1131        return Corr(newcontent, prange=self.prange)
1132
1133    def exp(self):
1134        newcontent = [None if _check_for_none(self, item) else np.exp(item) for item in self.content]
1135        return Corr(newcontent, prange=self.prange)
1136
1137    def _apply_func_to_corr(self, func):
1138        newcontent = [None if _check_for_none(self, item) else func(item) for item in self.content]
1139        for t in range(self.T):
1140            if _check_for_none(self, newcontent[t]):
1141                continue
1142            tmp_sum = np.sum(newcontent[t])
1143            if hasattr(tmp_sum, "value"):
1144                if np.isnan(tmp_sum.value):
1145                    newcontent[t] = None
1146        if all([item is None for item in newcontent]):
1147            raise Exception('Operation returns undefined correlator')
1148        return Corr(newcontent)
1149
1150    def sin(self):
1151        return self._apply_func_to_corr(np.sin)
1152
1153    def cos(self):
1154        return self._apply_func_to_corr(np.cos)
1155
1156    def tan(self):
1157        return self._apply_func_to_corr(np.tan)
1158
1159    def sinh(self):
1160        return self._apply_func_to_corr(np.sinh)
1161
1162    def cosh(self):
1163        return self._apply_func_to_corr(np.cosh)
1164
1165    def tanh(self):
1166        return self._apply_func_to_corr(np.tanh)
1167
1168    def arcsin(self):
1169        return self._apply_func_to_corr(np.arcsin)
1170
1171    def arccos(self):
1172        return self._apply_func_to_corr(np.arccos)
1173
1174    def arctan(self):
1175        return self._apply_func_to_corr(np.arctan)
1176
1177    def arcsinh(self):
1178        return self._apply_func_to_corr(np.arcsinh)
1179
1180    def arccosh(self):
1181        return self._apply_func_to_corr(np.arccosh)
1182
1183    def arctanh(self):
1184        return self._apply_func_to_corr(np.arctanh)
1185
1186    # Right hand side operations (require tweak in main module to work)
1187    def __radd__(self, y):
1188        return self + y
1189
1190    def __rsub__(self, y):
1191        return -self + y
1192
1193    def __rmul__(self, y):
1194        return self * y
1195
1196    def __rtruediv__(self, y):
1197        return (self / y) ** (-1)
1198
1199    @property
1200    def real(self):
1201        def return_real(obs_OR_cobs):
1202            if isinstance(obs_OR_cobs.flatten()[0], CObs):
1203                return np.vectorize(lambda x: x.real)(obs_OR_cobs)
1204            else:
1205                return obs_OR_cobs
1206
1207        return self._apply_func_to_corr(return_real)
1208
1209    @property
1210    def imag(self):
1211        def return_imag(obs_OR_cobs):
1212            if isinstance(obs_OR_cobs.flatten()[0], CObs):
1213                return np.vectorize(lambda x: x.imag)(obs_OR_cobs)
1214            else:
1215                return obs_OR_cobs * 0  # So it stays the right type
1216
1217        return self._apply_func_to_corr(return_imag)
1218
1219    def prune(self, Ntrunc, tproj=3, t0proj=2, basematrix=None):
1220        r''' Project large correlation matrix to lowest states
1221
1222        This method can be used to reduce the size of an (N x N) correlation matrix
1223        to (Ntrunc x Ntrunc) by solving a GEVP at very early times where the noise
1224        is still small.
1225
1226        Parameters
1227        ----------
1228        Ntrunc: int
1229            Rank of the target matrix.
1230        tproj: int
1231            Time where the eigenvectors are evaluated, corresponds to ts in the GEVP method.
1232            The default value is 3.
1233        t0proj: int
1234            Time where the correlation matrix is inverted. Choosing t0proj=1 is strongly
1235            discouraged for O(a) improved theories, since the correctness of the procedure
1236            cannot be granted in this case. The default value is 2.
1237        basematrix : Corr
1238            Correlation matrix that is used to determine the eigenvectors of the
1239            lowest states based on a GEVP. basematrix is taken to be the Corr itself if
1240            is is not specified.
1241
1242        Notes
1243        -----
1244        We have the basematrix $C(t)$ and the target matrix $G(t)$. We start by solving
1245        the GEVP $$C(t) v_n(t, t_0) = \lambda_n(t, t_0) C(t_0) v_n(t, t_0)$$ where $t \equiv t_\mathrm{proj}$
1246        and $t_0 \equiv t_{0, \mathrm{proj}}$. The target matrix is projected onto the subspace of the
1247        resulting eigenvectors $v_n, n=1,\dots,N_\mathrm{trunc}$ via
1248        $$G^\prime_{i, j}(t) = (v_i, G(t) v_j)$$. This allows to reduce the size of a large
1249        correlation matrix and to remove some noise that is added by irrelevant operators.
1250        This may allow to use the GEVP on $G(t)$ at late times such that the theoretically motivated
1251        bound $t_0 \leq t/2$ holds, since the condition number of $G(t)$ is decreased, compared to $C(t)$.
1252        '''
1253
1254        if self.N == 1:
1255            raise Exception('Method cannot be applied to one-dimensional correlators.')
1256        if basematrix is None:
1257            basematrix = self
1258        if Ntrunc >= basematrix.N:
1259            raise Exception('Cannot truncate using Ntrunc <= %d' % (basematrix.N))
1260        if basematrix.N != self.N:
1261            raise Exception('basematrix and targetmatrix have to be of the same size.')
1262
1263        evecs = basematrix.GEVP(t0proj, tproj, sort=None)[:Ntrunc]
1264
1265        tmpmat = np.empty((Ntrunc, Ntrunc), dtype=object)
1266        rmat = []
1267        for t in range(basematrix.T):
1268            for i in range(Ntrunc):
1269                for j in range(Ntrunc):
1270                    tmpmat[i][j] = evecs[i].T @ self[t] @ evecs[j]
1271            rmat.append(np.copy(tmpmat))
1272
1273        newcontent = [None if (self.content[t] is None) else rmat[t] for t in range(self.T)]
1274        return Corr(newcontent)

The class for a correlator (time dependent sequence of pe.Obs).

Everything, this class does, can be achieved using lists or arrays of Obs. But it is simply more convenient to have a dedicated object for correlators. One often wants to add or multiply correlators of the same length at every timeslice and it is inconvenient to iterate over all timeslices for every operation. This is especially true, when dealing with matrices.

The correlator can have two types of content: An Obs at every timeslice OR a GEVP matrix at every timeslice. Other dependency (eg. spatial) are not supported.

Corr(data_input, padding=[0, 0], prange=None)
27    def __init__(self, data_input, padding=[0, 0], prange=None):
28        """ Initialize a Corr object.
29
30        Parameters
31        ----------
32        data_input : list or array
33            list of Obs or list of arrays of Obs or array of Corrs
34        padding : list, optional
35            List with two entries where the first labels the padding
36            at the front of the correlator and the second the padding
37            at the back.
38        prange : list, optional
39            List containing the first and last timeslice of the plateau
40            region indentified for this correlator.
41        """
42
43        if isinstance(data_input, np.ndarray):
44
45            # This only works, if the array fulfills the conditions below
46            if not len(data_input.shape) == 2 and data_input.shape[0] == data_input.shape[1]:
47                raise Exception("Incompatible array shape")
48            if not all([isinstance(item, Corr) for item in data_input.flatten()]):
49                raise Exception("If the input is an array, its elements must be of type pe.Corr")
50            if not all([item.N == 1 for item in data_input.flatten()]):
51                raise Exception("Can only construct matrix correlator from single valued correlators")
52            if not len(set([item.T for item in data_input.flatten()])) == 1:
53                raise Exception("All input Correlators must be defined over the same timeslices.")
54
55            T = data_input[0, 0].T
56            N = data_input.shape[0]
57            input_as_list = []
58            for t in range(T):
59                if any([(item.content[t] is None) for item in data_input.flatten()]):
60                    if not all([(item.content[t] is None) for item in data_input.flatten()]):
61                        warnings.warn("Input ill-defined at different timeslices. Conversion leads to data loss!", RuntimeWarning)
62                    input_as_list.append(None)
63                else:
64                    array_at_timeslace = np.empty([N, N], dtype="object")
65                    for i in range(N):
66                        for j in range(N):
67                            array_at_timeslace[i, j] = data_input[i, j][t]
68                    input_as_list.append(array_at_timeslace)
69            data_input = input_as_list
70
71        if isinstance(data_input, list):
72
73            if all([isinstance(item, (Obs, CObs)) or item is None for item in data_input]):
74                _assert_equal_properties([o for o in data_input if o is not None])
75                self.content = [np.asarray([item]) if item is not None else None for item in data_input]
76                self.N = 1
77
78            elif all([isinstance(item, np.ndarray) or item is None for item in data_input]) and any([isinstance(item, np.ndarray) for item in data_input]):
79                self.content = data_input
80                noNull = [a for a in self.content if not (a is None)]  # To check if the matrices are correct for all undefined elements
81                self.N = noNull[0].shape[0]
82                if self.N > 1 and noNull[0].shape[0] != noNull[0].shape[1]:
83                    raise Exception("Smearing matrices are not NxN")
84                if (not all([item.shape == noNull[0].shape for item in noNull])):
85                    raise Exception("Items in data_input are not of identical shape." + str(noNull))
86            else:
87                raise Exception("data_input contains item of wrong type")
88        else:
89            raise Exception("Data input was not given as list or correct array")
90
91        self.tag = None
92
93        # An undefined timeslice is represented by the None object
94        self.content = [None] * padding[0] + self.content + [None] * padding[1]
95        self.T = len(self.content)
96        self.prange = prange

Initialize a Corr object.

Parameters
  • data_input (list or array): list of Obs or list of arrays of Obs or array of Corrs
  • padding (list, optional): List with two entries where the first labels the padding at the front of the correlator and the second the padding at the back.
  • prange (list, optional): List containing the first and last timeslice of the plateau region indentified for this correlator.
def gamma_method(self, **kwargs):
117    def gamma_method(self, **kwargs):
118        """Apply the gamma method to the content of the Corr."""
119        for item in self.content:
120            if not (item is None):
121                if self.N == 1:
122                    item[0].gamma_method(**kwargs)
123                else:
124                    for i in range(self.N):
125                        for j in range(self.N):
126                            item[i, j].gamma_method(**kwargs)

Apply the gamma method to the content of the Corr.

def projected(self, vector_l=None, vector_r=None, normalize=False):
128    def projected(self, vector_l=None, vector_r=None, normalize=False):
129        """We need to project the Correlator with a Vector to get a single value at each timeslice.
130
131        The method can use one or two vectors.
132        If two are specified it returns v1@G@v2 (the order might be very important.)
133        By default it will return the lowest source, which usually means unsmeared-unsmeared (0,0), but it does not have to
134        """
135        if self.N == 1:
136            raise Exception("Trying to project a Corr, that already has N=1.")
137
138        if vector_l is None:
139            vector_l, vector_r = np.asarray([1.] + (self.N - 1) * [0.]), np.asarray([1.] + (self.N - 1) * [0.])
140        elif (vector_r is None):
141            vector_r = vector_l
142        if isinstance(vector_l, list) and not isinstance(vector_r, list):
143            if len(vector_l) != self.T:
144                raise Exception("Length of vector list must be equal to T")
145            vector_r = [vector_r] * self.T
146        if isinstance(vector_r, list) and not isinstance(vector_l, list):
147            if len(vector_r) != self.T:
148                raise Exception("Length of vector list must be equal to T")
149            vector_l = [vector_l] * self.T
150
151        if not isinstance(vector_l, list):
152            if not vector_l.shape == vector_r.shape == (self.N,):
153                raise Exception("Vectors are of wrong shape!")
154            if normalize:
155                vector_l, vector_r = vector_l / np.sqrt((vector_l @ vector_l)), vector_r / np.sqrt(vector_r @ vector_r)
156            newcontent = [None if _check_for_none(self, item) else np.asarray([vector_l.T @ item @ vector_r]) for item in self.content]
157
158        else:
159            # There are no checks here yet. There are so many possible scenarios, where this can go wrong.
160            if normalize:
161                for t in range(self.T):
162                    vector_l[t], vector_r[t] = vector_l[t] / np.sqrt((vector_l[t] @ vector_l[t])), vector_r[t] / np.sqrt(vector_r[t] @ vector_r[t])
163
164            newcontent = [None if (_check_for_none(self, self.content[t]) or vector_l[t] is None or vector_r[t] is None) else np.asarray([vector_l[t].T @ self.content[t] @ vector_r[t]]) for t in range(self.T)]
165        return Corr(newcontent)

We need to project the Correlator with a Vector to get a single value at each timeslice.

The method can use one or two vectors. If two are specified it returns v1@G@v2 (the order might be very important.) By default it will return the lowest source, which usually means unsmeared-unsmeared (0,0), but it does not have to

def item(self, i, j):
167    def item(self, i, j):
168        """Picks the element [i,j] from every matrix and returns a correlator containing one Obs per timeslice.
169
170        Parameters
171        ----------
172        i : int
173            First index to be picked.
174        j : int
175            Second index to be picked.
176        """
177        if self.N == 1:
178            raise Exception("Trying to pick item from projected Corr")
179        newcontent = [None if (item is None) else item[i, j] for item in self.content]
180        return Corr(newcontent)

Picks the element [i,j] from every matrix and returns a correlator containing one Obs per timeslice.

Parameters
  • i (int): First index to be picked.
  • j (int): Second index to be picked.
def plottable(self):
182    def plottable(self):
183        """Outputs the correlator in a plotable format.
184
185        Outputs three lists containing the timeslice index, the value on each
186        timeslice and the error on each timeslice.
187        """
188        if self.N != 1:
189            raise Exception("Can only make Corr[N=1] plottable")
190        x_list = [x for x in range(self.T) if not self.content[x] is None]
191        y_list = [y[0].value for y in self.content if y is not None]
192        y_err_list = [y[0].dvalue for y in self.content if y is not None]
193
194        return x_list, y_list, y_err_list

Outputs the correlator in a plotable format.

Outputs three lists containing the timeslice index, the value on each timeslice and the error on each timeslice.

def symmetric(self):
196    def symmetric(self):
197        """ Symmetrize the correlator around x0=0."""
198        if self.N != 1:
199            raise Exception('symmetric cannot be safely applied to multi-dimensional correlators.')
200        if self.T % 2 != 0:
201            raise Exception("Can not symmetrize odd T")
202
203        if np.argmax(np.abs(self.content)) != 0:
204            warnings.warn("Correlator does not seem to be symmetric around x0=0.", RuntimeWarning)
205
206        newcontent = [self.content[0]]
207        for t in range(1, self.T):
208            if (self.content[t] is None) or (self.content[self.T - t] is None):
209                newcontent.append(None)
210            else:
211                newcontent.append(0.5 * (self.content[t] + self.content[self.T - t]))
212        if (all([x is None for x in newcontent])):
213            raise Exception("Corr could not be symmetrized: No redundant values")
214        return Corr(newcontent, prange=self.prange)

Symmetrize the correlator around x0=0.

def anti_symmetric(self):
216    def anti_symmetric(self):
217        """Anti-symmetrize the correlator around x0=0."""
218        if self.N != 1:
219            raise Exception('anti_symmetric cannot be safely applied to multi-dimensional correlators.')
220        if self.T % 2 != 0:
221            raise Exception("Can not symmetrize odd T")
222
223        test = 1 * self
224        test.gamma_method()
225        if not all([o.is_zero_within_error(3) for o in test.content[0]]):
226            warnings.warn("Correlator does not seem to be anti-symmetric around x0=0.", RuntimeWarning)
227
228        newcontent = [self.content[0]]
229        for t in range(1, self.T):
230            if (self.content[t] is None) or (self.content[self.T - t] is None):
231                newcontent.append(None)
232            else:
233                newcontent.append(0.5 * (self.content[t] - self.content[self.T - t]))
234        if (all([x is None for x in newcontent])):
235            raise Exception("Corr could not be symmetrized: No redundant values")
236        return Corr(newcontent, prange=self.prange)

Anti-symmetrize the correlator around x0=0.

def is_matrix_symmetric(self):
238    def is_matrix_symmetric(self):
239        """Checks whether a correlator matrices is symmetric on every timeslice."""
240        if self.N == 1:
241            raise Exception("Only works for correlator matrices.")
242        for t in range(self.T):
243            if self[t] is None:
244                continue
245            for i in range(self.N):
246                for j in range(i + 1, self.N):
247                    if self[t][i, j] is self[t][j, i]:
248                        continue
249                    if hash(self[t][i, j]) != hash(self[t][j, i]):
250                        return False
251        return True

Checks whether a correlator matrices is symmetric on every timeslice.

def matrix_symmetric(self):
253    def matrix_symmetric(self):
254        """Symmetrizes the correlator matrices on every timeslice."""
255        if self.N == 1:
256            raise Exception("Trying to symmetrize a correlator matrix, that already has N=1.")
257        if self.is_matrix_symmetric():
258            return 1.0 * self
259        else:
260            transposed = [None if _check_for_none(self, G) else G.T for G in self.content]
261            return 0.5 * (Corr(transposed) + self)

Symmetrizes the correlator matrices on every timeslice.

def GEVP(self, t0, ts=None, sort='Eigenvalue', **kwargs):
263    def GEVP(self, t0, ts=None, sort="Eigenvalue", **kwargs):
264        r'''Solve the generalized eigenvalue problem on the correlator matrix and returns the corresponding eigenvectors.
265
266        The eigenvectors are sorted according to the descending eigenvalues, the zeroth eigenvector(s) correspond to the
267        largest eigenvalue(s). The eigenvector(s) for the individual states can be accessed via slicing
268        ```python
269        C.GEVP(t0=2)[0]  # Ground state vector(s)
270        C.GEVP(t0=2)[:3]  # Vectors for the lowest three states
271        ```
272
273        Parameters
274        ----------
275        t0 : int
276            The time t0 for the right hand side of the GEVP according to $G(t)v_i=\lambda_i G(t_0)v_i$
277        ts : int
278            fixed time $G(t_s)v_i=\lambda_i G(t_0)v_i$ if sort=None.
279            If sort="Eigenvector" it gives a reference point for the sorting method.
280        sort : string
281            If this argument is set, a list of self.T vectors per state is returned. If it is set to None, only one vector is returned.
282            - "Eigenvalue": The eigenvector is chosen according to which eigenvalue it belongs individually on every timeslice.
283            - "Eigenvector": Use the method described in arXiv:2004.10472 to find the set of v(t) belonging to the state.
284              The reference state is identified by its eigenvalue at $t=t_s$.
285
286        Other Parameters
287        ----------------
288        state : int
289           Returns only the vector(s) for a specified state. The lowest state is zero.
290        '''
291
292        if self.N == 1:
293            raise Exception("GEVP methods only works on correlator matrices and not single correlators.")
294        if ts is not None:
295            if (ts <= t0):
296                raise Exception("ts has to be larger than t0.")
297
298        if "sorted_list" in kwargs:
299            warnings.warn("Argument 'sorted_list' is deprecated, use 'sort' instead.", DeprecationWarning)
300            sort = kwargs.get("sorted_list")
301
302        if self.is_matrix_symmetric():
303            symmetric_corr = self
304        else:
305            symmetric_corr = self.matrix_symmetric()
306
307        G0 = np.vectorize(lambda x: x.value)(symmetric_corr[t0])
308        np.linalg.cholesky(G0)  # Check if matrix G0 is positive-semidefinite.
309
310        if sort is None:
311            if (ts is None):
312                raise Exception("ts is required if sort=None.")
313            if (self.content[t0] is None) or (self.content[ts] is None):
314                raise Exception("Corr not defined at t0/ts.")
315            Gt = np.vectorize(lambda x: x.value)(symmetric_corr[ts])
316            reordered_vecs = _GEVP_solver(Gt, G0)
317
318        elif sort in ["Eigenvalue", "Eigenvector"]:
319            if sort == "Eigenvalue" and ts is not None:
320                warnings.warn("ts has no effect when sorting by eigenvalue is chosen.", RuntimeWarning)
321            all_vecs = [None] * (t0 + 1)
322            for t in range(t0 + 1, self.T):
323                try:
324                    Gt = np.vectorize(lambda x: x.value)(symmetric_corr[t])
325                    all_vecs.append(_GEVP_solver(Gt, G0))
326                except Exception:
327                    all_vecs.append(None)
328            if sort == "Eigenvector":
329                if ts is None:
330                    raise Exception("ts is required for the Eigenvector sorting method.")
331                all_vecs = _sort_vectors(all_vecs, ts)
332
333            reordered_vecs = [[v[s] if v is not None else None for v in all_vecs] for s in range(self.N)]
334        else:
335            raise Exception("Unkown value for 'sort'.")
336
337        if "state" in kwargs:
338            return reordered_vecs[kwargs.get("state")]
339        else:
340            return reordered_vecs

Solve the generalized eigenvalue problem on the correlator matrix and returns the corresponding eigenvectors.

The eigenvectors are sorted according to the descending eigenvalues, the zeroth eigenvector(s) correspond to the largest eigenvalue(s). The eigenvector(s) for the individual states can be accessed via slicing

C.GEVP(t0=2)[0]  # Ground state vector(s)
C.GEVP(t0=2)[:3]  # Vectors for the lowest three states
Parameters
  • t0 (int): The time t0 for the right hand side of the GEVP according to $G(t)v_i=\lambda_i G(t_0)v_i$
  • ts (int): fixed time $G(t_s)v_i=\lambda_i G(t_0)v_i$ if sort=None. If sort="Eigenvector" it gives a reference point for the sorting method.
  • sort (string): If this argument is set, a list of self.T vectors per state is returned. If it is set to None, only one vector is returned.
    • "Eigenvalue": The eigenvector is chosen according to which eigenvalue it belongs individually on every timeslice.
    • "Eigenvector": Use the method described in arXiv:2004.10472 to find the set of v(t) belonging to the state. The reference state is identified by its eigenvalue at $t=t_s$.
Other Parameters
  • state (int): Returns only the vector(s) for a specified state. The lowest state is zero.
def Eigenvalue(self, t0, ts=None, state=0, sort='Eigenvalue'):
342    def Eigenvalue(self, t0, ts=None, state=0, sort="Eigenvalue"):
343        """Determines the eigenvalue of the GEVP by solving and projecting the correlator
344
345        Parameters
346        ----------
347        state : int
348            The state one is interested in ordered by energy. The lowest state is zero.
349
350        All other parameters are identical to the ones of Corr.GEVP.
351        """
352        vec = self.GEVP(t0, ts=ts, sort=sort)[state]
353        return self.projected(vec)

Determines the eigenvalue of the GEVP by solving and projecting the correlator

Parameters
  • state (int): The state one is interested in ordered by energy. The lowest state is zero.
  • All other parameters are identical to the ones of Corr.GEVP.
def Hankel(self, N, periodic=False):
355    def Hankel(self, N, periodic=False):
356        """Constructs an NxN Hankel matrix
357
358        C(t) c(t+1) ... c(t+n-1)
359        C(t+1) c(t+2) ... c(t+n)
360        .................
361        C(t+(n-1)) c(t+n) ... c(t+2(n-1))
362
363        Parameters
364        ----------
365        N : int
366            Dimension of the Hankel matrix
367        periodic : bool, optional
368            determines whether the matrix is extended periodically
369        """
370
371        if self.N != 1:
372            raise Exception("Multi-operator Prony not implemented!")
373
374        array = np.empty([N, N], dtype="object")
375        new_content = []
376        for t in range(self.T):
377            new_content.append(array.copy())
378
379        def wrap(i):
380            while i >= self.T:
381                i -= self.T
382            return i
383
384        for t in range(self.T):
385            for i in range(N):
386                for j in range(N):
387                    if periodic:
388                        new_content[t][i, j] = self.content[wrap(t + i + j)][0]
389                    elif (t + i + j) >= self.T:
390                        new_content[t] = None
391                    else:
392                        new_content[t][i, j] = self.content[t + i + j][0]
393
394        return Corr(new_content)

Constructs an NxN Hankel matrix

C(t) c(t+1) ... c(t+n-1) C(t+1) c(t+2) ... c(t+n) ................. C(t+(n-1)) c(t+n) ... c(t+2(n-1))

Parameters
  • N (int): Dimension of the Hankel matrix
  • periodic (bool, optional): determines whether the matrix is extended periodically
def roll(self, dt):
396    def roll(self, dt):
397        """Periodically shift the correlator by dt timeslices
398
399        Parameters
400        ----------
401        dt : int
402            number of timeslices
403        """
404        return Corr(list(np.roll(np.array(self.content, dtype=object), dt)))

Periodically shift the correlator by dt timeslices

Parameters
  • dt (int): number of timeslices
def reverse(self):
406    def reverse(self):
407        """Reverse the time ordering of the Corr"""
408        return Corr(self.content[:: -1])

Reverse the time ordering of the Corr

def thin(self, spacing=2, offset=0):
410    def thin(self, spacing=2, offset=0):
411        """Thin out a correlator to suppress correlations
412
413        Parameters
414        ----------
415        spacing : int
416            Keep only every 'spacing'th entry of the correlator
417        offset : int
418            Offset the equal spacing
419        """
420        new_content = []
421        for t in range(self.T):
422            if (offset + t) % spacing != 0:
423                new_content.append(None)
424            else:
425                new_content.append(self.content[t])
426        return Corr(new_content)

Thin out a correlator to suppress correlations

Parameters
  • spacing (int): Keep only every 'spacing'th entry of the correlator
  • offset (int): Offset the equal spacing
def correlate(self, partner):
428    def correlate(self, partner):
429        """Correlate the correlator with another correlator or Obs
430
431        Parameters
432        ----------
433        partner : Obs or Corr
434            partner to correlate the correlator with.
435            Can either be an Obs which is correlated with all entries of the
436            correlator or a Corr of same length.
437        """
438        if self.N != 1:
439            raise Exception("Only one-dimensional correlators can be safely correlated.")
440        new_content = []
441        for x0, t_slice in enumerate(self.content):
442            if _check_for_none(self, t_slice):
443                new_content.append(None)
444            else:
445                if isinstance(partner, Corr):
446                    if _check_for_none(partner, partner.content[x0]):
447                        new_content.append(None)
448                    else:
449                        new_content.append(np.array([correlate(o, partner.content[x0][0]) for o in t_slice]))
450                elif isinstance(partner, Obs):  # Should this include CObs?
451                    new_content.append(np.array([correlate(o, partner) for o in t_slice]))
452                else:
453                    raise Exception("Can only correlate with an Obs or a Corr.")
454
455        return Corr(new_content)

Correlate the correlator with another correlator or Obs

Parameters
  • partner (Obs or Corr): partner to correlate the correlator with. Can either be an Obs which is correlated with all entries of the correlator or a Corr of same length.
def reweight(self, weight, **kwargs):
457    def reweight(self, weight, **kwargs):
458        """Reweight the correlator.
459
460        Parameters
461        ----------
462        weight : Obs
463            Reweighting factor. An Observable that has to be defined on a superset of the
464            configurations in obs[i].idl for all i.
465        all_configs : bool
466            if True, the reweighted observables are normalized by the average of
467            the reweighting factor on all configurations in weight.idl and not
468            on the configurations in obs[i].idl.
469        """
470        if self.N != 1:
471            raise Exception("Reweighting only implemented for one-dimensional correlators.")
472        new_content = []
473        for t_slice in self.content:
474            if _check_for_none(self, t_slice):
475                new_content.append(None)
476            else:
477                new_content.append(np.array(reweight(weight, t_slice, **kwargs)))
478        return Corr(new_content)

Reweight the correlator.

Parameters
  • weight (Obs): Reweighting factor. An Observable that has to be defined on a superset of the configurations in obs[i].idl for all i.
  • all_configs (bool): if True, the reweighted observables are normalized by the average of the reweighting factor on all configurations in weight.idl and not on the configurations in obs[i].idl.
def T_symmetry(self, partner, parity=1):
480    def T_symmetry(self, partner, parity=+1):
481        """Return the time symmetry average of the correlator and its partner
482
483        Parameters
484        ----------
485        partner : Corr
486            Time symmetry partner of the Corr
487        partity : int
488            Parity quantum number of the correlator, can be +1 or -1
489        """
490        if self.N != 1:
491            raise Exception("T_symmetry only implemented for one-dimensional correlators.")
492        if not isinstance(partner, Corr):
493            raise Exception("T partner has to be a Corr object.")
494        if parity not in [+1, -1]:
495            raise Exception("Parity has to be +1 or -1.")
496        T_partner = parity * partner.reverse()
497
498        t_slices = []
499        test = (self - T_partner)
500        test.gamma_method()
501        for x0, t_slice in enumerate(test.content):
502            if t_slice is not None:
503                if not t_slice[0].is_zero_within_error(5):
504                    t_slices.append(x0)
505        if t_slices:
506            warnings.warn("T symmetry partners do not agree within 5 sigma on time slices " + str(t_slices) + ".", RuntimeWarning)
507
508        return (self + T_partner) / 2

Return the time symmetry average of the correlator and its partner

Parameters
  • partner (Corr): Time symmetry partner of the Corr
  • partity (int): Parity quantum number of the correlator, can be +1 or -1
def deriv(self, variant='symmetric'):
510    def deriv(self, variant="symmetric"):
511        """Return the first derivative of the correlator with respect to x0.
512
513        Parameters
514        ----------
515        variant : str
516            decides which definition of the finite differences derivative is used.
517            Available choice: symmetric, forward, backward, improved, log, default: symmetric
518        """
519        if self.N != 1:
520            raise Exception("deriv only implemented for one-dimensional correlators.")
521        if variant == "symmetric":
522            newcontent = []
523            for t in range(1, self.T - 1):
524                if (self.content[t - 1] is None) or (self.content[t + 1] is None):
525                    newcontent.append(None)
526                else:
527                    newcontent.append(0.5 * (self.content[t + 1] - self.content[t - 1]))
528            if (all([x is None for x in newcontent])):
529                raise Exception('Derivative is undefined at all timeslices')
530            return Corr(newcontent, padding=[1, 1])
531        elif variant == "forward":
532            newcontent = []
533            for t in range(self.T - 1):
534                if (self.content[t] is None) or (self.content[t + 1] is None):
535                    newcontent.append(None)
536                else:
537                    newcontent.append(self.content[t + 1] - self.content[t])
538            if (all([x is None for x in newcontent])):
539                raise Exception("Derivative is undefined at all timeslices")
540            return Corr(newcontent, padding=[0, 1])
541        elif variant == "backward":
542            newcontent = []
543            for t in range(1, self.T):
544                if (self.content[t - 1] is None) or (self.content[t] is None):
545                    newcontent.append(None)
546                else:
547                    newcontent.append(self.content[t] - self.content[t - 1])
548            if (all([x is None for x in newcontent])):
549                raise Exception("Derivative is undefined at all timeslices")
550            return Corr(newcontent, padding=[1, 0])
551        elif variant == "improved":
552            newcontent = []
553            for t in range(2, self.T - 2):
554                if (self.content[t - 2] is None) or (self.content[t - 1] is None) or (self.content[t + 1] is None) or (self.content[t + 2] is None):
555                    newcontent.append(None)
556                else:
557                    newcontent.append((1 / 12) * (self.content[t - 2] - 8 * self.content[t - 1] + 8 * self.content[t + 1] - self.content[t + 2]))
558            if (all([x is None for x in newcontent])):
559                raise Exception('Derivative is undefined at all timeslices')
560            return Corr(newcontent, padding=[2, 2])
561        elif variant == 'log':
562            newcontent = []
563            for t in range(self.T):
564                if (self.content[t] is None) or (self.content[t] <= 0):
565                    newcontent.append(None)
566                else:
567                    newcontent.append(np.log(self.content[t]))
568            if (all([x is None for x in newcontent])):
569                raise Exception("Log is undefined at all timeslices")
570            logcorr = Corr(newcontent)
571            return self * logcorr.deriv('symmetric')
572        else:
573            raise Exception("Unknown variant.")

Return the first derivative of the correlator with respect to x0.

Parameters
  • variant (str): decides which definition of the finite differences derivative is used. Available choice: symmetric, forward, backward, improved, log, default: symmetric
def second_deriv(self, variant='symmetric'):
575    def second_deriv(self, variant="symmetric"):
576        """Return the second derivative of the correlator with respect to x0.
577
578        Parameters
579        ----------
580        variant : str
581            decides which definition of the finite differences derivative is used.
582            Available choice: symmetric, improved, log, default: symmetric
583        """
584        if self.N != 1:
585            raise Exception("second_deriv only implemented for one-dimensional correlators.")
586        if variant == "symmetric":
587            newcontent = []
588            for t in range(1, self.T - 1):
589                if (self.content[t - 1] is None) or (self.content[t + 1] is None):
590                    newcontent.append(None)
591                else:
592                    newcontent.append((self.content[t + 1] - 2 * self.content[t] + self.content[t - 1]))
593            if (all([x is None for x in newcontent])):
594                raise Exception("Derivative is undefined at all timeslices")
595            return Corr(newcontent, padding=[1, 1])
596        elif variant == "improved":
597            newcontent = []
598            for t in range(2, self.T - 2):
599                if (self.content[t - 2] is None) or (self.content[t - 1] is None) or (self.content[t] is None) or (self.content[t + 1] is None) or (self.content[t + 2] is None):
600                    newcontent.append(None)
601                else:
602                    newcontent.append((1 / 12) * (-self.content[t + 2] + 16 * self.content[t + 1] - 30 * self.content[t] + 16 * self.content[t - 1] - self.content[t - 2]))
603            if (all([x is None for x in newcontent])):
604                raise Exception("Derivative is undefined at all timeslices")
605            return Corr(newcontent, padding=[2, 2])
606        elif variant == 'log':
607            newcontent = []
608            for t in range(self.T):
609                if (self.content[t] is None) or (self.content[t] <= 0):
610                    newcontent.append(None)
611                else:
612                    newcontent.append(np.log(self.content[t]))
613            if (all([x is None for x in newcontent])):
614                raise Exception("Log is undefined at all timeslices")
615            logcorr = Corr(newcontent)
616            return self * (logcorr.second_deriv('symmetric') + (logcorr.deriv('symmetric'))**2)
617        else:
618            raise Exception("Unknown variant.")

Return the second derivative of the correlator with respect to x0.

Parameters
  • variant (str): decides which definition of the finite differences derivative is used. Available choice: symmetric, improved, log, default: symmetric
def m_eff(self, variant='log', guess=1.0):
620    def m_eff(self, variant='log', guess=1.0):
621        """Returns the effective mass of the correlator as correlator object
622
623        Parameters
624        ----------
625        variant : str
626            log : uses the standard effective mass log(C(t) / C(t+1))
627            cosh, periodic : Use periodicitiy of the correlator by solving C(t) / C(t+1) = cosh(m * (t - T/2)) / cosh(m * (t + 1 - T/2)) for m.
628            sinh : Use anti-periodicitiy of the correlator by solving C(t) / C(t+1) = sinh(m * (t - T/2)) / sinh(m * (t + 1 - T/2)) for m.
629            See, e.g., arXiv:1205.5380
630            arccosh : Uses the explicit form of the symmetrized correlator (not recommended)
631            logsym: uses the symmetric effective mass log(C(t-1) / C(t+1))/2
632        guess : float
633            guess for the root finder, only relevant for the root variant
634        """
635        if self.N != 1:
636            raise Exception('Correlator must be projected before getting m_eff')
637        if variant == 'log':
638            newcontent = []
639            for t in range(self.T - 1):
640                if ((self.content[t] is None) or (self.content[t + 1] is None)) or (self.content[t + 1][0].value == 0):
641                    newcontent.append(None)
642                elif self.content[t][0].value / self.content[t + 1][0].value < 0:
643                    newcontent.append(None)
644                else:
645                    newcontent.append(self.content[t] / self.content[t + 1])
646            if (all([x is None for x in newcontent])):
647                raise Exception('m_eff is undefined at all timeslices')
648
649            return np.log(Corr(newcontent, padding=[0, 1]))
650
651        elif variant == 'logsym':
652            newcontent = []
653            for t in range(1, self.T - 1):
654                if ((self.content[t - 1] is None) or (self.content[t + 1] is None)) or (self.content[t + 1][0].value == 0):
655                    newcontent.append(None)
656                elif self.content[t - 1][0].value / self.content[t + 1][0].value < 0:
657                    newcontent.append(None)
658                else:
659                    newcontent.append(self.content[t - 1] / self.content[t + 1])
660            if (all([x is None for x in newcontent])):
661                raise Exception('m_eff is undefined at all timeslices')
662
663            return np.log(Corr(newcontent, padding=[1, 1])) / 2
664
665        elif variant in ['periodic', 'cosh', 'sinh']:
666            if variant in ['periodic', 'cosh']:
667                func = anp.cosh
668            else:
669                func = anp.sinh
670
671            def root_function(x, d):
672                return func(x * (t - self.T / 2)) / func(x * (t + 1 - self.T / 2)) - d
673
674            newcontent = []
675            for t in range(self.T - 1):
676                if (self.content[t] is None) or (self.content[t + 1] is None) or (self.content[t + 1][0].value == 0):
677                    newcontent.append(None)
678                # Fill the two timeslices in the middle of the lattice with their predecessors
679                elif variant == 'sinh' and t in [self.T / 2, self.T / 2 - 1]:
680                    newcontent.append(newcontent[-1])
681                elif self.content[t][0].value / self.content[t + 1][0].value < 0:
682                    newcontent.append(None)
683                else:
684                    newcontent.append(np.abs(find_root(self.content[t][0] / self.content[t + 1][0], root_function, guess=guess)))
685            if (all([x is None for x in newcontent])):
686                raise Exception('m_eff is undefined at all timeslices')
687
688            return Corr(newcontent, padding=[0, 1])
689
690        elif variant == 'arccosh':
691            newcontent = []
692            for t in range(1, self.T - 1):
693                if (self.content[t] is None) or (self.content[t + 1] is None) or (self.content[t - 1] is None) or (self.content[t][0].value == 0):
694                    newcontent.append(None)
695                else:
696                    newcontent.append((self.content[t + 1] + self.content[t - 1]) / (2 * self.content[t]))
697            if (all([x is None for x in newcontent])):
698                raise Exception("m_eff is undefined at all timeslices")
699            return np.arccosh(Corr(newcontent, padding=[1, 1]))
700
701        else:
702            raise Exception('Unknown variant.')

Returns the effective mass of the correlator as correlator object

Parameters
  • variant (str): log : uses the standard effective mass log(C(t) / C(t+1)) cosh, periodic : Use periodicitiy of the correlator by solving C(t) / C(t+1) = cosh(m * (t - T/2)) / cosh(m * (t + 1 - T/2)) for m. sinh : Use anti-periodicitiy of the correlator by solving C(t) / C(t+1) = sinh(m * (t - T/2)) / sinh(m * (t + 1 - T/2)) for m. See, e.g., arXiv:1205.5380 arccosh : Uses the explicit form of the symmetrized correlator (not recommended) logsym: uses the symmetric effective mass log(C(t-1) / C(t+1))/2
  • guess (float): guess for the root finder, only relevant for the root variant
def fit(self, function, fitrange=None, silent=False, **kwargs):
704    def fit(self, function, fitrange=None, silent=False, **kwargs):
705        r'''Fits function to the data
706
707        Parameters
708        ----------
709        function : obj
710            function to fit to the data. See fits.least_squares for details.
711        fitrange : list
712            Two element list containing the timeslices on which the fit is supposed to start and stop.
713            Caution: This range is inclusive as opposed to standard python indexing.
714            `fitrange=[4, 6]` corresponds to the three entries 4, 5 and 6.
715            If not specified, self.prange or all timeslices are used.
716        silent : bool
717            Decides whether output is printed to the standard output.
718        '''
719        if self.N != 1:
720            raise Exception("Correlator must be projected before fitting")
721
722        if fitrange is None:
723            if self.prange:
724                fitrange = self.prange
725            else:
726                fitrange = [0, self.T - 1]
727        else:
728            if not isinstance(fitrange, list):
729                raise Exception("fitrange has to be a list with two elements")
730            if len(fitrange) != 2:
731                raise Exception("fitrange has to have exactly two elements [fit_start, fit_stop]")
732
733        xs = [x for x in range(fitrange[0], fitrange[1] + 1) if not self.content[x] is None]
734        ys = [self.content[x][0] for x in range(fitrange[0], fitrange[1] + 1) if not self.content[x] is None]
735        result = least_squares(xs, ys, function, silent=silent, **kwargs)
736        return result

Fits function to the data

Parameters
  • function (obj): function to fit to the data. See fits.least_squares for details.
  • fitrange (list): Two element list containing the timeslices on which the fit is supposed to start and stop. Caution: This range is inclusive as opposed to standard python indexing. fitrange=[4, 6] corresponds to the three entries 4, 5 and 6. If not specified, self.prange or all timeslices are used.
  • silent (bool): Decides whether output is printed to the standard output.
def plateau(self, plateau_range=None, method='fit', auto_gamma=False):
738    def plateau(self, plateau_range=None, method="fit", auto_gamma=False):
739        """ Extract a plateau value from a Corr object
740
741        Parameters
742        ----------
743        plateau_range : list
744            list with two entries, indicating the first and the last timeslice
745            of the plateau region.
746        method : str
747            method to extract the plateau.
748                'fit' fits a constant to the plateau region
749                'avg', 'average' or 'mean' just average over the given timeslices.
750        auto_gamma : bool
751            apply gamma_method with default parameters to the Corr. Defaults to None
752        """
753        if not plateau_range:
754            if self.prange:
755                plateau_range = self.prange
756            else:
757                raise Exception("no plateau range provided")
758        if self.N != 1:
759            raise Exception("Correlator must be projected before getting a plateau.")
760        if (all([self.content[t] is None for t in range(plateau_range[0], plateau_range[1] + 1)])):
761            raise Exception("plateau is undefined at all timeslices in plateaurange.")
762        if auto_gamma:
763            self.gamma_method()
764        if method == "fit":
765            def const_func(a, t):
766                return a[0]
767            return self.fit(const_func, plateau_range)[0]
768        elif method in ["avg", "average", "mean"]:
769            returnvalue = np.mean([item[0] for item in self.content[plateau_range[0]:plateau_range[1] + 1] if item is not None])
770            return returnvalue
771
772        else:
773            raise Exception("Unsupported plateau method: " + method)

Extract a plateau value from a Corr object

Parameters
  • plateau_range (list): list with two entries, indicating the first and the last timeslice of the plateau region.
  • method (str): method to extract the plateau. 'fit' fits a constant to the plateau region 'avg', 'average' or 'mean' just average over the given timeslices.
  • auto_gamma (bool): apply gamma_method with default parameters to the Corr. Defaults to None
def set_prange(self, prange):
775    def set_prange(self, prange):
776        """Sets the attribute prange of the Corr object."""
777        if not len(prange) == 2:
778            raise Exception("prange must be a list or array with two values")
779        if not ((isinstance(prange[0], int)) and (isinstance(prange[1], int))):
780            raise Exception("Start and end point must be integers")
781        if not (0 <= prange[0] <= self.T and 0 <= prange[1] <= self.T and prange[0] < prange[1]):
782            raise Exception("Start and end point must define a range in the interval 0,T")
783
784        self.prange = prange
785        return

Sets the attribute prange of the Corr object.

def show( self, x_range=None, comp=None, y_range=None, logscale=False, plateau=None, fit_res=None, ylabel=None, save=None, auto_gamma=False, hide_sigma=None, references=None, title=None):
787    def show(self, x_range=None, comp=None, y_range=None, logscale=False, plateau=None, fit_res=None, ylabel=None, save=None, auto_gamma=False, hide_sigma=None, references=None, title=None):
788        """Plots the correlator using the tag of the correlator as label if available.
789
790        Parameters
791        ----------
792        x_range : list
793            list of two values, determining the range of the x-axis e.g. [4, 8].
794        comp : Corr or list of Corr
795            Correlator or list of correlators which are plotted for comparison.
796            The tags of these correlators are used as labels if available.
797        logscale : bool
798            Sets y-axis to logscale.
799        plateau : Obs
800            Plateau value to be visualized in the figure.
801        fit_res : Fit_result
802            Fit_result object to be visualized.
803        ylabel : str
804            Label for the y-axis.
805        save : str
806            path to file in which the figure should be saved.
807        auto_gamma : bool
808            Apply the gamma method with standard parameters to all correlators and plateau values before plotting.
809        hide_sigma : float
810            Hides data points from the first value on which is consistent with zero within 'hide_sigma' standard errors.
811        references : list
812            List of floating point values that are displayed as horizontal lines for reference.
813        title : string
814            Optional title of the figure.
815        """
816        if self.N != 1:
817            raise Exception("Correlator must be projected before plotting")
818
819        if auto_gamma:
820            self.gamma_method()
821
822        if x_range is None:
823            x_range = [0, self.T - 1]
824
825        fig = plt.figure()
826        ax1 = fig.add_subplot(111)
827
828        x, y, y_err = self.plottable()
829        if hide_sigma:
830            hide_from = np.argmax((hide_sigma * np.array(y_err[1:])) > np.abs(y[1:])) - 1
831        else:
832            hide_from = None
833        ax1.errorbar(x[:hide_from], y[:hide_from], y_err[:hide_from], label=self.tag)
834        if logscale:
835            ax1.set_yscale('log')
836        else:
837            if y_range is None:
838                try:
839                    y_min = min([(x[0].value - x[0].dvalue) for x in self.content[x_range[0]: x_range[1] + 1] if (x is not None) and x[0].dvalue < 2 * np.abs(x[0].value)])
840                    y_max = max([(x[0].value + x[0].dvalue) for x in self.content[x_range[0]: x_range[1] + 1] if (x is not None) and x[0].dvalue < 2 * np.abs(x[0].value)])
841                    ax1.set_ylim([y_min - 0.1 * (y_max - y_min), y_max + 0.1 * (y_max - y_min)])
842                except Exception:
843                    pass
844            else:
845                ax1.set_ylim(y_range)
846        if comp:
847            if isinstance(comp, (Corr, list)):
848                for corr in comp if isinstance(comp, list) else [comp]:
849                    if auto_gamma:
850                        corr.gamma_method()
851                    x, y, y_err = corr.plottable()
852                    if hide_sigma:
853                        hide_from = np.argmax((hide_sigma * np.array(y_err[1:])) > np.abs(y[1:])) - 1
854                    else:
855                        hide_from = None
856                    ax1.errorbar(x[:hide_from], y[:hide_from], y_err[:hide_from], label=corr.tag, mfc=plt.rcParams['axes.facecolor'])
857            else:
858                raise Exception("'comp' must be a correlator or a list of correlators.")
859
860        if plateau:
861            if isinstance(plateau, Obs):
862                if auto_gamma:
863                    plateau.gamma_method()
864                ax1.axhline(y=plateau.value, linewidth=2, color=plt.rcParams['text.color'], alpha=0.6, marker=',', ls='--', label=str(plateau))
865                ax1.axhspan(plateau.value - plateau.dvalue, plateau.value + plateau.dvalue, alpha=0.25, color=plt.rcParams['text.color'], ls='-')
866            else:
867                raise Exception("'plateau' must be an Obs")
868
869        if references:
870            if isinstance(references, list):
871                for ref in references:
872                    ax1.axhline(y=ref, linewidth=1, color=plt.rcParams['text.color'], alpha=0.6, marker=',', ls='--')
873            else:
874                raise Exception("'references' must be a list of floating pint values.")
875
876        if self.prange:
877            ax1.axvline(self.prange[0], 0, 1, ls='-', marker=',')
878            ax1.axvline(self.prange[1], 0, 1, ls='-', marker=',')
879
880        if fit_res:
881            x_samples = np.arange(x_range[0], x_range[1] + 1, 0.05)
882            ax1.plot(x_samples,
883                     fit_res.fit_function([o.value for o in fit_res.fit_parameters], x_samples),
884                     ls='-', marker=',', lw=2)
885
886        ax1.set_xlabel(r'$x_0 / a$')
887        if ylabel:
888            ax1.set_ylabel(ylabel)
889        ax1.set_xlim([x_range[0] - 0.5, x_range[1] + 0.5])
890
891        handles, labels = ax1.get_legend_handles_labels()
892        if labels:
893            ax1.legend()
894
895        if title:
896            plt.title(title)
897
898        plt.draw()
899
900        if save:
901            if isinstance(save, str):
902                fig.savefig(save, bbox_inches='tight')
903            else:
904                raise Exception("'save' has to be a string.")

Plots the correlator using the tag of the correlator as label if available.

Parameters
  • x_range (list): list of two values, determining the range of the x-axis e.g. [4, 8].
  • comp (Corr or list of Corr): Correlator or list of correlators which are plotted for comparison. The tags of these correlators are used as labels if available.
  • logscale (bool): Sets y-axis to logscale.
  • plateau (Obs): Plateau value to be visualized in the figure.
  • fit_res (Fit_result): Fit_result object to be visualized.
  • ylabel (str): Label for the y-axis.
  • save (str): path to file in which the figure should be saved.
  • auto_gamma (bool): Apply the gamma method with standard parameters to all correlators and plateau values before plotting.
  • hide_sigma (float): Hides data points from the first value on which is consistent with zero within 'hide_sigma' standard errors.
  • references (list): List of floating point values that are displayed as horizontal lines for reference.
  • title (string): Optional title of the figure.
def spaghetti_plot(self, logscale=True):
906    def spaghetti_plot(self, logscale=True):
907        """Produces a spaghetti plot of the correlator suited to monitor exceptional configurations.
908
909        Parameters
910        ----------
911        logscale : bool
912            Determines whether the scale of the y-axis is logarithmic or standard.
913        """
914        if self.N != 1:
915            raise Exception("Correlator needs to be projected first.")
916
917        mc_names = list(set([item for sublist in [sum(map(o[0].e_content.get, o[0].mc_names), []) for o in self.content if o is not None] for item in sublist]))
918        x0_vals = [n for (n, o) in zip(np.arange(self.T), self.content) if o is not None]
919
920        for name in mc_names:
921            data = np.array([o[0].deltas[name] + o[0].r_values[name] for o in self.content if o is not None]).T
922
923            fig = plt.figure()
924            ax = fig.add_subplot(111)
925            for dat in data:
926                ax.plot(x0_vals, dat, ls='-', marker='')
927
928            if logscale is True:
929                ax.set_yscale('log')
930
931            ax.set_xlabel(r'$x_0 / a$')
932            plt.title(name)
933            plt.draw()

Produces a spaghetti plot of the correlator suited to monitor exceptional configurations.

Parameters
  • logscale (bool): Determines whether the scale of the y-axis is logarithmic or standard.
def dump(self, filename, datatype='json.gz', **kwargs):
935    def dump(self, filename, datatype="json.gz", **kwargs):
936        """Dumps the Corr into a file of chosen type
937        Parameters
938        ----------
939        filename : str
940            Name of the file to be saved.
941        datatype : str
942            Format of the exported file. Supported formats include
943            "json.gz" and "pickle"
944        path : str
945            specifies a custom path for the file (default '.')
946        """
947        if datatype == "json.gz":
948            from .input.json import dump_to_json
949            if 'path' in kwargs:
950                file_name = kwargs.get('path') + '/' + filename
951            else:
952                file_name = filename
953            dump_to_json(self, file_name)
954        elif datatype == "pickle":
955            dump_object(self, filename, **kwargs)
956        else:
957            raise Exception("Unknown datatype " + str(datatype))

Dumps the Corr into a file of chosen type

Parameters
  • filename (str): Name of the file to be saved.
  • datatype (str): Format of the exported file. Supported formats include "json.gz" and "pickle"
  • path (str): specifies a custom path for the file (default '.')
def print(self, print_range=None):
959    def print(self, print_range=None):
960        print(self.__repr__(print_range))
def sqrt(self):
1126    def sqrt(self):
1127        return self ** 0.5
def log(self):
1129    def log(self):
1130        newcontent = [None if _check_for_none(self, item) else np.log(item) for item in self.content]
1131        return Corr(newcontent, prange=self.prange)
def exp(self):
1133    def exp(self):
1134        newcontent = [None if _check_for_none(self, item) else np.exp(item) for item in self.content]
1135        return Corr(newcontent, prange=self.prange)
def sin(self):
1150    def sin(self):
1151        return self._apply_func_to_corr(np.sin)
def cos(self):
1153    def cos(self):
1154        return self._apply_func_to_corr(np.cos)
def tan(self):
1156    def tan(self):
1157        return self._apply_func_to_corr(np.tan)
def sinh(self):
1159    def sinh(self):
1160        return self._apply_func_to_corr(np.sinh)
def cosh(self):
1162    def cosh(self):
1163        return self._apply_func_to_corr(np.cosh)
def tanh(self):
1165    def tanh(self):
1166        return self._apply_func_to_corr(np.tanh)
def arcsin(self):
1168    def arcsin(self):
1169        return self._apply_func_to_corr(np.arcsin)
def arccos(self):
1171    def arccos(self):
1172        return self._apply_func_to_corr(np.arccos)
def arctan(self):
1174    def arctan(self):
1175        return self._apply_func_to_corr(np.arctan)
def arcsinh(self):
1177    def arcsinh(self):
1178        return self._apply_func_to_corr(np.arcsinh)
def arccosh(self):
1180    def arccosh(self):
1181        return self._apply_func_to_corr(np.arccosh)
def arctanh(self):
1183    def arctanh(self):
1184        return self._apply_func_to_corr(np.arctanh)
def prune(self, Ntrunc, tproj=3, t0proj=2, basematrix=None):
1219    def prune(self, Ntrunc, tproj=3, t0proj=2, basematrix=None):
1220        r''' Project large correlation matrix to lowest states
1221
1222        This method can be used to reduce the size of an (N x N) correlation matrix
1223        to (Ntrunc x Ntrunc) by solving a GEVP at very early times where the noise
1224        is still small.
1225
1226        Parameters
1227        ----------
1228        Ntrunc: int
1229            Rank of the target matrix.
1230        tproj: int
1231            Time where the eigenvectors are evaluated, corresponds to ts in the GEVP method.
1232            The default value is 3.
1233        t0proj: int
1234            Time where the correlation matrix is inverted. Choosing t0proj=1 is strongly
1235            discouraged for O(a) improved theories, since the correctness of the procedure
1236            cannot be granted in this case. The default value is 2.
1237        basematrix : Corr
1238            Correlation matrix that is used to determine the eigenvectors of the
1239            lowest states based on a GEVP. basematrix is taken to be the Corr itself if
1240            is is not specified.
1241
1242        Notes
1243        -----
1244        We have the basematrix $C(t)$ and the target matrix $G(t)$. We start by solving
1245        the GEVP $$C(t) v_n(t, t_0) = \lambda_n(t, t_0) C(t_0) v_n(t, t_0)$$ where $t \equiv t_\mathrm{proj}$
1246        and $t_0 \equiv t_{0, \mathrm{proj}}$. The target matrix is projected onto the subspace of the
1247        resulting eigenvectors $v_n, n=1,\dots,N_\mathrm{trunc}$ via
1248        $$G^\prime_{i, j}(t) = (v_i, G(t) v_j)$$. This allows to reduce the size of a large
1249        correlation matrix and to remove some noise that is added by irrelevant operators.
1250        This may allow to use the GEVP on $G(t)$ at late times such that the theoretically motivated
1251        bound $t_0 \leq t/2$ holds, since the condition number of $G(t)$ is decreased, compared to $C(t)$.
1252        '''
1253
1254        if self.N == 1:
1255            raise Exception('Method cannot be applied to one-dimensional correlators.')
1256        if basematrix is None:
1257            basematrix = self
1258        if Ntrunc >= basematrix.N:
1259            raise Exception('Cannot truncate using Ntrunc <= %d' % (basematrix.N))
1260        if basematrix.N != self.N:
1261            raise Exception('basematrix and targetmatrix have to be of the same size.')
1262
1263        evecs = basematrix.GEVP(t0proj, tproj, sort=None)[:Ntrunc]
1264
1265        tmpmat = np.empty((Ntrunc, Ntrunc), dtype=object)
1266        rmat = []
1267        for t in range(basematrix.T):
1268            for i in range(Ntrunc):
1269                for j in range(Ntrunc):
1270                    tmpmat[i][j] = evecs[i].T @ self[t] @ evecs[j]
1271            rmat.append(np.copy(tmpmat))
1272
1273        newcontent = [None if (self.content[t] is None) else rmat[t] for t in range(self.T)]
1274        return Corr(newcontent)

Project large correlation matrix to lowest states

This method can be used to reduce the size of an (N x N) correlation matrix to (Ntrunc x Ntrunc) by solving a GEVP at very early times where the noise is still small.

Parameters
  • Ntrunc (int): Rank of the target matrix.
  • tproj (int): Time where the eigenvectors are evaluated, corresponds to ts in the GEVP method. The default value is 3.
  • t0proj (int): Time where the correlation matrix is inverted. Choosing t0proj=1 is strongly discouraged for O(a) improved theories, since the correctness of the procedure cannot be granted in this case. The default value is 2.
  • basematrix (Corr): Correlation matrix that is used to determine the eigenvectors of the lowest states based on a GEVP. basematrix is taken to be the Corr itself if is is not specified.
Notes

We have the basematrix $C(t)$ and the target matrix $G(t)$. We start by solving the GEVP $$C(t) v_n(t, t_0) = \lambda_n(t, t_0) C(t_0) v_n(t, t_0)$$ where $t \equiv t_\mathrm{proj}$ and $t_0 \equiv t_{0, \mathrm{proj}}$. The target matrix is projected onto the subspace of the resulting eigenvectors $v_n, n=1,\dots,N_\mathrm{trunc}$ via $$G^\prime_{i, j}(t) = (v_i, G(t) v_j)$$. This allows to reduce the size of a large correlation matrix and to remove some noise that is added by irrelevant operators. This may allow to use the GEVP on $G(t)$ at late times such that the theoretically motivated bound $t_0 \leq t/2$ holds, since the condition number of $G(t)$ is decreased, compared to $C(t)$.