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

Apply the gamma method to the content of the Corr.

def gm(self, **kwargs):
119    def gamma_method(self, **kwargs):
120        """Apply the gamma method to the content of the Corr."""
121        for item in self.content:
122            if not (item is None):
123                if self.N == 1:
124                    item[0].gamma_method(**kwargs)
125                else:
126                    for i in range(self.N):
127                        for j in range(self.N):
128                            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):
132    def projected(self, vector_l=None, vector_r=None, normalize=False):
133        """We need to project the Correlator with a Vector to get a single value at each timeslice.
134
135        The method can use one or two vectors.
136        If two are specified it returns v1@G@v2 (the order might be very important.)
137        By default it will return the lowest source, which usually means unsmeared-unsmeared (0,0), but it does not have to
138        """
139        if self.N == 1:
140            raise Exception("Trying to project a Corr, that already has N=1.")
141
142        if vector_l is None:
143            vector_l, vector_r = np.asarray([1.] + (self.N - 1) * [0.]), np.asarray([1.] + (self.N - 1) * [0.])
144        elif (vector_r is None):
145            vector_r = vector_l
146        if isinstance(vector_l, list) and not isinstance(vector_r, list):
147            if len(vector_l) != self.T:
148                raise Exception("Length of vector list must be equal to T")
149            vector_r = [vector_r] * self.T
150        if isinstance(vector_r, list) and not isinstance(vector_l, list):
151            if len(vector_r) != self.T:
152                raise Exception("Length of vector list must be equal to T")
153            vector_l = [vector_l] * self.T
154
155        if not isinstance(vector_l, list):
156            if not vector_l.shape == vector_r.shape == (self.N,):
157                raise Exception("Vectors are of wrong shape!")
158            if normalize:
159                vector_l, vector_r = vector_l / np.sqrt((vector_l @ vector_l)), vector_r / np.sqrt(vector_r @ vector_r)
160            newcontent = [None if _check_for_none(self, item) else np.asarray([vector_l.T @ item @ vector_r]) for item in self.content]
161
162        else:
163            # There are no checks here yet. There are so many possible scenarios, where this can go wrong.
164            if normalize:
165                for t in range(self.T):
166                    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])
167
168            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)]
169        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):
171    def item(self, i, j):
172        """Picks the element [i,j] from every matrix and returns a correlator containing one Obs per timeslice.
173
174        Parameters
175        ----------
176        i : int
177            First index to be picked.
178        j : int
179            Second index to be picked.
180        """
181        if self.N == 1:
182            raise Exception("Trying to pick item from projected Corr")
183        newcontent = [None if (item is None) else item[i, j] for item in self.content]
184        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):
186    def plottable(self):
187        """Outputs the correlator in a plotable format.
188
189        Outputs three lists containing the timeslice index, the value on each
190        timeslice and the error on each timeslice.
191        """
192        if self.N != 1:
193            raise Exception("Can only make Corr[N=1] plottable")
194        x_list = [x for x in range(self.T) if not self.content[x] is None]
195        y_list = [y[0].value for y in self.content if y is not None]
196        y_err_list = [y[0].dvalue for y in self.content if y is not None]
197
198        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):
200    def symmetric(self):
201        """ Symmetrize the correlator around x0=0."""
202        if self.N != 1:
203            raise Exception('symmetric cannot be safely applied to multi-dimensional correlators.')
204        if self.T % 2 != 0:
205            raise Exception("Can not symmetrize odd T")
206
207        if self.content[0] is not None:
208            if np.argmax(np.abs([o[0].value if o is not None else 0 for o in self.content])) != 0:
209                warnings.warn("Correlator does not seem to be symmetric around x0=0.", RuntimeWarning)
210
211        newcontent = [self.content[0]]
212        for t in range(1, self.T):
213            if (self.content[t] is None) or (self.content[self.T - t] is None):
214                newcontent.append(None)
215            else:
216                newcontent.append(0.5 * (self.content[t] + self.content[self.T - t]))
217        if (all([x is None for x in newcontent])):
218            raise Exception("Corr could not be symmetrized: No redundant values")
219        return Corr(newcontent, prange=self.prange)

Symmetrize the correlator around x0=0.

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

Anti-symmetrize the correlator around x0=0.

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

Checks whether a correlator matrices is symmetric on every timeslice.

def trace(self):
258    def trace(self):
259        """Calculates the per-timeslice trace of a correlator matrix."""
260        if self.N == 1:
261            raise ValueError("Only works for correlator matrices.")
262        newcontent = []
263        for t in range(self.T):
264            if _check_for_none(self, self.content[t]):
265                newcontent.append(None)
266            else:
267                newcontent.append(np.trace(self.content[t]))
268        return Corr(newcontent)

Calculates the per-timeslice trace of a correlator matrix.

def matrix_symmetric(self):
270    def matrix_symmetric(self):
271        """Symmetrizes the correlator matrices on every timeslice."""
272        if self.N == 1:
273            raise Exception("Trying to symmetrize a correlator matrix, that already has N=1.")
274        if self.is_matrix_symmetric():
275            return 1.0 * self
276        else:
277            transposed = [None if _check_for_none(self, G) else G.T for G in self.content]
278            return 0.5 * (Corr(transposed) + self)

Symmetrizes the correlator matrices on every timeslice.

def GEVP(self, t0, ts=None, sort='Eigenvalue', **kwargs):
280    def GEVP(self, t0, ts=None, sort="Eigenvalue", **kwargs):
281        r'''Solve the generalized eigenvalue problem on the correlator matrix and returns the corresponding eigenvectors.
282
283        The eigenvectors are sorted according to the descending eigenvalues, the zeroth eigenvector(s) correspond to the
284        largest eigenvalue(s). The eigenvector(s) for the individual states can be accessed via slicing
285        ```python
286        C.GEVP(t0=2)[0]  # Ground state vector(s)
287        C.GEVP(t0=2)[:3]  # Vectors for the lowest three states
288        ```
289
290        Parameters
291        ----------
292        t0 : int
293            The time t0 for the right hand side of the GEVP according to $G(t)v_i=\lambda_i G(t_0)v_i$
294        ts : int
295            fixed time $G(t_s)v_i=\lambda_i G(t_0)v_i$ if sort=None.
296            If sort="Eigenvector" it gives a reference point for the sorting method.
297        sort : string
298            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.
299            - "Eigenvalue": The eigenvector is chosen according to which eigenvalue it belongs individually on every timeslice.
300            - "Eigenvector": Use the method described in arXiv:2004.10472 to find the set of v(t) belonging to the state.
301              The reference state is identified by its eigenvalue at $t=t_s$.
302
303        Other Parameters
304        ----------------
305        state : int
306           Returns only the vector(s) for a specified state. The lowest state is zero.
307        '''
308
309        if self.N == 1:
310            raise Exception("GEVP methods only works on correlator matrices and not single correlators.")
311        if ts is not None:
312            if (ts <= t0):
313                raise Exception("ts has to be larger than t0.")
314
315        if "sorted_list" in kwargs:
316            warnings.warn("Argument 'sorted_list' is deprecated, use 'sort' instead.", DeprecationWarning)
317            sort = kwargs.get("sorted_list")
318
319        if self.is_matrix_symmetric():
320            symmetric_corr = self
321        else:
322            symmetric_corr = self.matrix_symmetric()
323
324        G0 = np.vectorize(lambda x: x.value)(symmetric_corr[t0])
325        np.linalg.cholesky(G0)  # Check if matrix G0 is positive-semidefinite.
326
327        if sort is None:
328            if (ts is None):
329                raise Exception("ts is required if sort=None.")
330            if (self.content[t0] is None) or (self.content[ts] is None):
331                raise Exception("Corr not defined at t0/ts.")
332            Gt = np.vectorize(lambda x: x.value)(symmetric_corr[ts])
333            reordered_vecs = _GEVP_solver(Gt, G0)
334
335        elif sort in ["Eigenvalue", "Eigenvector"]:
336            if sort == "Eigenvalue" and ts is not None:
337                warnings.warn("ts has no effect when sorting by eigenvalue is chosen.", RuntimeWarning)
338            all_vecs = [None] * (t0 + 1)
339            for t in range(t0 + 1, self.T):
340                try:
341                    Gt = np.vectorize(lambda x: x.value)(symmetric_corr[t])
342                    all_vecs.append(_GEVP_solver(Gt, G0))
343                except Exception:
344                    all_vecs.append(None)
345            if sort == "Eigenvector":
346                if ts is None:
347                    raise Exception("ts is required for the Eigenvector sorting method.")
348                all_vecs = _sort_vectors(all_vecs, ts)
349
350            reordered_vecs = [[v[s] if v is not None else None for v in all_vecs] for s in range(self.N)]
351        else:
352            raise Exception("Unkown value for 'sort'.")
353
354        if "state" in kwargs:
355            return reordered_vecs[kwargs.get("state")]
356        else:
357            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'):
359    def Eigenvalue(self, t0, ts=None, state=0, sort="Eigenvalue"):
360        """Determines the eigenvalue of the GEVP by solving and projecting the correlator
361
362        Parameters
363        ----------
364        state : int
365            The state one is interested in ordered by energy. The lowest state is zero.
366
367        All other parameters are identical to the ones of Corr.GEVP.
368        """
369        vec = self.GEVP(t0, ts=ts, sort=sort)[state]
370        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):
372    def Hankel(self, N, periodic=False):
373        """Constructs an NxN Hankel matrix
374
375        C(t) c(t+1) ... c(t+n-1)
376        C(t+1) c(t+2) ... c(t+n)
377        .................
378        C(t+(n-1)) c(t+n) ... c(t+2(n-1))
379
380        Parameters
381        ----------
382        N : int
383            Dimension of the Hankel matrix
384        periodic : bool, optional
385            determines whether the matrix is extended periodically
386        """
387
388        if self.N != 1:
389            raise Exception("Multi-operator Prony not implemented!")
390
391        array = np.empty([N, N], dtype="object")
392        new_content = []
393        for t in range(self.T):
394            new_content.append(array.copy())
395
396        def wrap(i):
397            while i >= self.T:
398                i -= self.T
399            return i
400
401        for t in range(self.T):
402            for i in range(N):
403                for j in range(N):
404                    if periodic:
405                        new_content[t][i, j] = self.content[wrap(t + i + j)][0]
406                    elif (t + i + j) >= self.T:
407                        new_content[t] = None
408                    else:
409                        new_content[t][i, j] = self.content[t + i + j][0]
410
411        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):
413    def roll(self, dt):
414        """Periodically shift the correlator by dt timeslices
415
416        Parameters
417        ----------
418        dt : int
419            number of timeslices
420        """
421        return Corr(list(np.roll(np.array(self.content, dtype=object), dt, axis=0)))

Periodically shift the correlator by dt timeslices

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

Reverse the time ordering of the Corr

def thin(self, spacing=2, offset=0):
427    def thin(self, spacing=2, offset=0):
428        """Thin out a correlator to suppress correlations
429
430        Parameters
431        ----------
432        spacing : int
433            Keep only every 'spacing'th entry of the correlator
434        offset : int
435            Offset the equal spacing
436        """
437        new_content = []
438        for t in range(self.T):
439            if (offset + t) % spacing != 0:
440                new_content.append(None)
441            else:
442                new_content.append(self.content[t])
443        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):
445    def correlate(self, partner):
446        """Correlate the correlator with another correlator or Obs
447
448        Parameters
449        ----------
450        partner : Obs or Corr
451            partner to correlate the correlator with.
452            Can either be an Obs which is correlated with all entries of the
453            correlator or a Corr of same length.
454        """
455        if self.N != 1:
456            raise Exception("Only one-dimensional correlators can be safely correlated.")
457        new_content = []
458        for x0, t_slice in enumerate(self.content):
459            if _check_for_none(self, t_slice):
460                new_content.append(None)
461            else:
462                if isinstance(partner, Corr):
463                    if _check_for_none(partner, partner.content[x0]):
464                        new_content.append(None)
465                    else:
466                        new_content.append(np.array([correlate(o, partner.content[x0][0]) for o in t_slice]))
467                elif isinstance(partner, Obs):  # Should this include CObs?
468                    new_content.append(np.array([correlate(o, partner) for o in t_slice]))
469                else:
470                    raise Exception("Can only correlate with an Obs or a Corr.")
471
472        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):
474    def reweight(self, weight, **kwargs):
475        """Reweight the correlator.
476
477        Parameters
478        ----------
479        weight : Obs
480            Reweighting factor. An Observable that has to be defined on a superset of the
481            configurations in obs[i].idl for all i.
482        all_configs : bool
483            if True, the reweighted observables are normalized by the average of
484            the reweighting factor on all configurations in weight.idl and not
485            on the configurations in obs[i].idl.
486        """
487        if self.N != 1:
488            raise Exception("Reweighting only implemented for one-dimensional correlators.")
489        new_content = []
490        for t_slice in self.content:
491            if _check_for_none(self, t_slice):
492                new_content.append(None)
493            else:
494                new_content.append(np.array(reweight(weight, t_slice, **kwargs)))
495        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):
497    def T_symmetry(self, partner, parity=+1):
498        """Return the time symmetry average of the correlator and its partner
499
500        Parameters
501        ----------
502        partner : Corr
503            Time symmetry partner of the Corr
504        partity : int
505            Parity quantum number of the correlator, can be +1 or -1
506        """
507        if self.N != 1:
508            raise Exception("T_symmetry only implemented for one-dimensional correlators.")
509        if not isinstance(partner, Corr):
510            raise Exception("T partner has to be a Corr object.")
511        if parity not in [+1, -1]:
512            raise Exception("Parity has to be +1 or -1.")
513        T_partner = parity * partner.reverse()
514
515        t_slices = []
516        test = (self - T_partner)
517        test.gamma_method()
518        for x0, t_slice in enumerate(test.content):
519            if t_slice is not None:
520                if not t_slice[0].is_zero_within_error(5):
521                    t_slices.append(x0)
522        if t_slices:
523            warnings.warn("T symmetry partners do not agree within 5 sigma on time slices " + str(t_slices) + ".", RuntimeWarning)
524
525        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'):
527    def deriv(self, variant="symmetric"):
528        """Return the first derivative of the correlator with respect to x0.
529
530        Parameters
531        ----------
532        variant : str
533            decides which definition of the finite differences derivative is used.
534            Available choice: symmetric, forward, backward, improved, log, default: symmetric
535        """
536        if self.N != 1:
537            raise Exception("deriv only implemented for one-dimensional correlators.")
538        if variant == "symmetric":
539            newcontent = []
540            for t in range(1, self.T - 1):
541                if (self.content[t - 1] is None) or (self.content[t + 1] is None):
542                    newcontent.append(None)
543                else:
544                    newcontent.append(0.5 * (self.content[t + 1] - self.content[t - 1]))
545            if (all([x is None for x in newcontent])):
546                raise Exception('Derivative is undefined at all timeslices')
547            return Corr(newcontent, padding=[1, 1])
548        elif variant == "forward":
549            newcontent = []
550            for t in range(self.T - 1):
551                if (self.content[t] is None) or (self.content[t + 1] is None):
552                    newcontent.append(None)
553                else:
554                    newcontent.append(self.content[t + 1] - self.content[t])
555            if (all([x is None for x in newcontent])):
556                raise Exception("Derivative is undefined at all timeslices")
557            return Corr(newcontent, padding=[0, 1])
558        elif variant == "backward":
559            newcontent = []
560            for t in range(1, self.T):
561                if (self.content[t - 1] is None) or (self.content[t] is None):
562                    newcontent.append(None)
563                else:
564                    newcontent.append(self.content[t] - self.content[t - 1])
565            if (all([x is None for x in newcontent])):
566                raise Exception("Derivative is undefined at all timeslices")
567            return Corr(newcontent, padding=[1, 0])
568        elif variant == "improved":
569            newcontent = []
570            for t in range(2, self.T - 2):
571                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):
572                    newcontent.append(None)
573                else:
574                    newcontent.append((1 / 12) * (self.content[t - 2] - 8 * self.content[t - 1] + 8 * self.content[t + 1] - self.content[t + 2]))
575            if (all([x is None for x in newcontent])):
576                raise Exception('Derivative is undefined at all timeslices')
577            return Corr(newcontent, padding=[2, 2])
578        elif variant == 'log':
579            newcontent = []
580            for t in range(self.T):
581                if (self.content[t] is None) or (self.content[t] <= 0):
582                    newcontent.append(None)
583                else:
584                    newcontent.append(np.log(self.content[t]))
585            if (all([x is None for x in newcontent])):
586                raise Exception("Log is undefined at all timeslices")
587            logcorr = Corr(newcontent)
588            return self * logcorr.deriv('symmetric')
589        else:
590            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'):
592    def second_deriv(self, variant="symmetric"):
593        r"""Return the second derivative of the correlator with respect to x0.
594
595        Parameters
596        ----------
597        variant : str
598            decides which definition of the finite differences derivative is used.
599            Available choice:
600                - symmetric (default)
601                    $$\tilde{\partial}^2_0 f(x_0) = f(x_0+1)-2f(x_0)+f(x_0-1)$$
602                - big_symmetric
603                    $$\partial^2_0 f(x_0) = \frac{f(x_0+2)-2f(x_0)+f(x_0-2)}{4}$$
604                - improved
605                    $$\partial^2_0 f(x_0) = \frac{-f(x_0+2) + 16 * f(x_0+1) - 30 * f(x_0) + 16 * f(x_0-1) - f(x_0-2)}{12}$$
606                - log
607                    $$f(x) = \tilde{\partial}^2_0 log(f(x_0))+(\tilde{\partial}_0 log(f(x_0)))^2$$
608        """
609        if self.N != 1:
610            raise Exception("second_deriv only implemented for one-dimensional correlators.")
611        if variant == "symmetric":
612            newcontent = []
613            for t in range(1, self.T - 1):
614                if (self.content[t - 1] is None) or (self.content[t + 1] is None):
615                    newcontent.append(None)
616                else:
617                    newcontent.append((self.content[t + 1] - 2 * self.content[t] + self.content[t - 1]))
618            if (all([x is None for x in newcontent])):
619                raise Exception("Derivative is undefined at all timeslices")
620            return Corr(newcontent, padding=[1, 1])
621        elif variant == "big_symmetric":
622            newcontent = []
623            for t in range(2, self.T - 2):
624                if (self.content[t - 2] is None) or (self.content[t + 2] is None):
625                    newcontent.append(None)
626                else:
627                    newcontent.append((self.content[t + 2] - 2 * self.content[t] + self.content[t - 2]) / 4)
628            if (all([x is None for x in newcontent])):
629                raise Exception("Derivative is undefined at all timeslices")
630            return Corr(newcontent, padding=[2, 2])
631        elif variant == "improved":
632            newcontent = []
633            for t in range(2, self.T - 2):
634                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):
635                    newcontent.append(None)
636                else:
637                    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]))
638            if (all([x is None for x in newcontent])):
639                raise Exception("Derivative is undefined at all timeslices")
640            return Corr(newcontent, padding=[2, 2])
641        elif variant == 'log':
642            newcontent = []
643            for t in range(self.T):
644                if (self.content[t] is None) or (self.content[t] <= 0):
645                    newcontent.append(None)
646                else:
647                    newcontent.append(np.log(self.content[t]))
648            if (all([x is None for x in newcontent])):
649                raise Exception("Log is undefined at all timeslices")
650            logcorr = Corr(newcontent)
651            return self * (logcorr.second_deriv('symmetric') + (logcorr.deriv('symmetric'))**2)
652        else:
653            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 (default) $$\tilde{\partial}^2_0 f(x_0) = f(x_0+1)-2f(x_0)+f(x_0-1)$$ - big_symmetric $$\partial^2_0 f(x_0) = \frac{f(x_0+2)-2f(x_0)+f(x_0-2)}{4}$$ - improved $$\partial^2_0 f(x_0) = \frac{-f(x_0+2) + 16 * f(x_0+1) - 30 * f(x_0) + 16 * f(x_0-1) - f(x_0-2)}{12}$$ - log $$f(x) = \tilde{\partial}^2_0 log(f(x_0))+(\tilde{\partial}_0 log(f(x_0)))^2$$
def m_eff(self, variant='log', guess=1.0):
655    def m_eff(self, variant='log', guess=1.0):
656        """Returns the effective mass of the correlator as correlator object
657
658        Parameters
659        ----------
660        variant : str
661            log : uses the standard effective mass log(C(t) / C(t+1))
662            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.
663            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.
664            See, e.g., arXiv:1205.5380
665            arccosh : Uses the explicit form of the symmetrized correlator (not recommended)
666            logsym: uses the symmetric effective mass log(C(t-1) / C(t+1))/2
667        guess : float
668            guess for the root finder, only relevant for the root variant
669        """
670        if self.N != 1:
671            raise Exception('Correlator must be projected before getting m_eff')
672        if variant == 'log':
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                elif self.content[t][0].value / self.content[t + 1][0].value < 0:
678                    newcontent.append(None)
679                else:
680                    newcontent.append(self.content[t] / self.content[t + 1])
681            if (all([x is None for x in newcontent])):
682                raise Exception('m_eff is undefined at all timeslices')
683
684            return np.log(Corr(newcontent, padding=[0, 1]))
685
686        elif variant == 'logsym':
687            newcontent = []
688            for t in range(1, self.T - 1):
689                if ((self.content[t - 1] is None) or (self.content[t + 1] is None)) or (self.content[t + 1][0].value == 0):
690                    newcontent.append(None)
691                elif self.content[t - 1][0].value / self.content[t + 1][0].value < 0:
692                    newcontent.append(None)
693                else:
694                    newcontent.append(self.content[t - 1] / self.content[t + 1])
695            if (all([x is None for x in newcontent])):
696                raise Exception('m_eff is undefined at all timeslices')
697
698            return np.log(Corr(newcontent, padding=[1, 1])) / 2
699
700        elif variant in ['periodic', 'cosh', 'sinh']:
701            if variant in ['periodic', 'cosh']:
702                func = anp.cosh
703            else:
704                func = anp.sinh
705
706            def root_function(x, d):
707                return func(x * (t - self.T / 2)) / func(x * (t + 1 - self.T / 2)) - d
708
709            newcontent = []
710            for t in range(self.T - 1):
711                if (self.content[t] is None) or (self.content[t + 1] is None) or (self.content[t + 1][0].value == 0):
712                    newcontent.append(None)
713                # Fill the two timeslices in the middle of the lattice with their predecessors
714                elif variant == 'sinh' and t in [self.T / 2, self.T / 2 - 1]:
715                    newcontent.append(newcontent[-1])
716                elif self.content[t][0].value / self.content[t + 1][0].value < 0:
717                    newcontent.append(None)
718                else:
719                    newcontent.append(np.abs(find_root(self.content[t][0] / self.content[t + 1][0], root_function, guess=guess)))
720            if (all([x is None for x in newcontent])):
721                raise Exception('m_eff is undefined at all timeslices')
722
723            return Corr(newcontent, padding=[0, 1])
724
725        elif variant == 'arccosh':
726            newcontent = []
727            for t in range(1, self.T - 1):
728                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):
729                    newcontent.append(None)
730                else:
731                    newcontent.append((self.content[t + 1] + self.content[t - 1]) / (2 * self.content[t]))
732            if (all([x is None for x in newcontent])):
733                raise Exception("m_eff is undefined at all timeslices")
734            return np.arccosh(Corr(newcontent, padding=[1, 1]))
735
736        else:
737            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):
739    def fit(self, function, fitrange=None, silent=False, **kwargs):
740        r'''Fits function to the data
741
742        Parameters
743        ----------
744        function : obj
745            function to fit to the data. See fits.least_squares for details.
746        fitrange : list
747            Two element list containing the timeslices on which the fit is supposed to start and stop.
748            Caution: This range is inclusive as opposed to standard python indexing.
749            `fitrange=[4, 6]` corresponds to the three entries 4, 5 and 6.
750            If not specified, self.prange or all timeslices are used.
751        silent : bool
752            Decides whether output is printed to the standard output.
753        '''
754        if self.N != 1:
755            raise Exception("Correlator must be projected before fitting")
756
757        if fitrange is None:
758            if self.prange:
759                fitrange = self.prange
760            else:
761                fitrange = [0, self.T - 1]
762        else:
763            if not isinstance(fitrange, list):
764                raise Exception("fitrange has to be a list with two elements")
765            if len(fitrange) != 2:
766                raise Exception("fitrange has to have exactly two elements [fit_start, fit_stop]")
767
768        xs = np.array([x for x in range(fitrange[0], fitrange[1] + 1) if not self.content[x] is None])
769        ys = np.array([self.content[x][0] for x in range(fitrange[0], fitrange[1] + 1) if not self.content[x] is None])
770        result = least_squares(xs, ys, function, silent=silent, **kwargs)
771        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):
773    def plateau(self, plateau_range=None, method="fit", auto_gamma=False):
774        """ Extract a plateau value from a Corr object
775
776        Parameters
777        ----------
778        plateau_range : list
779            list with two entries, indicating the first and the last timeslice
780            of the plateau region.
781        method : str
782            method to extract the plateau.
783                'fit' fits a constant to the plateau region
784                'avg', 'average' or 'mean' just average over the given timeslices.
785        auto_gamma : bool
786            apply gamma_method with default parameters to the Corr. Defaults to None
787        """
788        if not plateau_range:
789            if self.prange:
790                plateau_range = self.prange
791            else:
792                raise Exception("no plateau range provided")
793        if self.N != 1:
794            raise Exception("Correlator must be projected before getting a plateau.")
795        if (all([self.content[t] is None for t in range(plateau_range[0], plateau_range[1] + 1)])):
796            raise Exception("plateau is undefined at all timeslices in plateaurange.")
797        if auto_gamma:
798            self.gamma_method()
799        if method == "fit":
800            def const_func(a, t):
801                return a[0]
802            return self.fit(const_func, plateau_range)[0]
803        elif method in ["avg", "average", "mean"]:
804            returnvalue = np.mean([item[0] for item in self.content[plateau_range[0]:plateau_range[1] + 1] if item is not None])
805            return returnvalue
806
807        else:
808            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):
810    def set_prange(self, prange):
811        """Sets the attribute prange of the Corr object."""
812        if not len(prange) == 2:
813            raise Exception("prange must be a list or array with two values")
814        if not ((isinstance(prange[0], int)) and (isinstance(prange[1], int))):
815            raise Exception("Start and end point must be integers")
816        if not (0 <= prange[0] <= self.T and 0 <= prange[1] <= self.T and prange[0] < prange[1]):
817            raise Exception("Start and end point must define a range in the interval 0,T")
818
819        self.prange = prange
820        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, fit_key=None, ylabel=None, save=None, auto_gamma=False, hide_sigma=None, references=None, title=None):
822    def show(self, x_range=None, comp=None, y_range=None, logscale=False, plateau=None, fit_res=None, fit_key=None, ylabel=None, save=None, auto_gamma=False, hide_sigma=None, references=None, title=None):
823        """Plots the correlator using the tag of the correlator as label if available.
824
825        Parameters
826        ----------
827        x_range : list
828            list of two values, determining the range of the x-axis e.g. [4, 8].
829        comp : Corr or list of Corr
830            Correlator or list of correlators which are plotted for comparison.
831            The tags of these correlators are used as labels if available.
832        logscale : bool
833            Sets y-axis to logscale.
834        plateau : Obs
835            Plateau value to be visualized in the figure.
836        fit_res : Fit_result
837            Fit_result object to be visualized.
838        fit_key : str
839            Key for the fit function in Fit_result.fit_function (for combined fits).
840        ylabel : str
841            Label for the y-axis.
842        save : str
843            path to file in which the figure should be saved.
844        auto_gamma : bool
845            Apply the gamma method with standard parameters to all correlators and plateau values before plotting.
846        hide_sigma : float
847            Hides data points from the first value on which is consistent with zero within 'hide_sigma' standard errors.
848        references : list
849            List of floating point values that are displayed as horizontal lines for reference.
850        title : string
851            Optional title of the figure.
852        """
853        if self.N != 1:
854            raise Exception("Correlator must be projected before plotting")
855
856        if auto_gamma:
857            self.gamma_method()
858
859        if x_range is None:
860            x_range = [0, self.T - 1]
861
862        fig = plt.figure()
863        ax1 = fig.add_subplot(111)
864
865        x, y, y_err = self.plottable()
866        if hide_sigma:
867            hide_from = np.argmax((hide_sigma * np.array(y_err[1:])) > np.abs(y[1:])) - 1
868        else:
869            hide_from = None
870        ax1.errorbar(x[:hide_from], y[:hide_from], y_err[:hide_from], label=self.tag)
871        if logscale:
872            ax1.set_yscale('log')
873        else:
874            if y_range is None:
875                try:
876                    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)])
877                    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)])
878                    ax1.set_ylim([y_min - 0.1 * (y_max - y_min), y_max + 0.1 * (y_max - y_min)])
879                except Exception:
880                    pass
881            else:
882                ax1.set_ylim(y_range)
883        if comp:
884            if isinstance(comp, (Corr, list)):
885                for corr in comp if isinstance(comp, list) else [comp]:
886                    if auto_gamma:
887                        corr.gamma_method()
888                    x, y, y_err = corr.plottable()
889                    if hide_sigma:
890                        hide_from = np.argmax((hide_sigma * np.array(y_err[1:])) > np.abs(y[1:])) - 1
891                    else:
892                        hide_from = None
893                    ax1.errorbar(x[:hide_from], y[:hide_from], y_err[:hide_from], label=corr.tag, mfc=plt.rcParams['axes.facecolor'])
894            else:
895                raise Exception("'comp' must be a correlator or a list of correlators.")
896
897        if plateau:
898            if isinstance(plateau, Obs):
899                if auto_gamma:
900                    plateau.gamma_method()
901                ax1.axhline(y=plateau.value, linewidth=2, color=plt.rcParams['text.color'], alpha=0.6, marker=',', ls='--', label=str(plateau))
902                ax1.axhspan(plateau.value - plateau.dvalue, plateau.value + plateau.dvalue, alpha=0.25, color=plt.rcParams['text.color'], ls='-')
903            else:
904                raise Exception("'plateau' must be an Obs")
905
906        if references:
907            if isinstance(references, list):
908                for ref in references:
909                    ax1.axhline(y=ref, linewidth=1, color=plt.rcParams['text.color'], alpha=0.6, marker=',', ls='--')
910            else:
911                raise Exception("'references' must be a list of floating pint values.")
912
913        if self.prange:
914            ax1.axvline(self.prange[0], 0, 1, ls='-', marker=',', color="black", zorder=0)
915            ax1.axvline(self.prange[1], 0, 1, ls='-', marker=',', color="black", zorder=0)
916
917        if fit_res:
918            x_samples = np.arange(x_range[0], x_range[1] + 1, 0.05)
919            if isinstance(fit_res.fit_function, dict):
920                if fit_key:
921                    ax1.plot(x_samples, fit_res.fit_function[fit_key]([o.value for o in fit_res.fit_parameters], x_samples), ls='-', marker=',', lw=2)
922                else:
923                    raise ValueError("Please provide a 'fit_key' for visualizing combined fits.")
924            else:
925                ax1.plot(x_samples, fit_res.fit_function([o.value for o in fit_res.fit_parameters], x_samples), ls='-', marker=',', lw=2)
926
927        ax1.set_xlabel(r'$x_0 / a$')
928        if ylabel:
929            ax1.set_ylabel(ylabel)
930        ax1.set_xlim([x_range[0] - 0.5, x_range[1] + 0.5])
931
932        handles, labels = ax1.get_legend_handles_labels()
933        if labels:
934            ax1.legend()
935
936        if title:
937            plt.title(title)
938
939        plt.draw()
940
941        if save:
942            if isinstance(save, str):
943                fig.savefig(save, bbox_inches='tight')
944            else:
945                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.
  • fit_key (str): Key for the fit function in Fit_result.fit_function (for combined fits).
  • 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):
947    def spaghetti_plot(self, logscale=True):
948        """Produces a spaghetti plot of the correlator suited to monitor exceptional configurations.
949
950        Parameters
951        ----------
952        logscale : bool
953            Determines whether the scale of the y-axis is logarithmic or standard.
954        """
955        if self.N != 1:
956            raise Exception("Correlator needs to be projected first.")
957
958        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]))
959        x0_vals = [n for (n, o) in zip(np.arange(self.T), self.content) if o is not None]
960
961        for name in mc_names:
962            data = np.array([o[0].deltas[name] + o[0].r_values[name] for o in self.content if o is not None]).T
963
964            fig = plt.figure()
965            ax = fig.add_subplot(111)
966            for dat in data:
967                ax.plot(x0_vals, dat, ls='-', marker='')
968
969            if logscale is True:
970                ax.set_yscale('log')
971
972            ax.set_xlabel(r'$x_0 / a$')
973            plt.title(name)
974            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):
976    def dump(self, filename, datatype="json.gz", **kwargs):
977        """Dumps the Corr into a file of chosen type
978        Parameters
979        ----------
980        filename : str
981            Name of the file to be saved.
982        datatype : str
983            Format of the exported file. Supported formats include
984            "json.gz" and "pickle"
985        path : str
986            specifies a custom path for the file (default '.')
987        """
988        if datatype == "json.gz":
989            from .input.json import dump_to_json
990            if 'path' in kwargs:
991                file_name = kwargs.get('path') + '/' + filename
992            else:
993                file_name = filename
994            dump_to_json(self, file_name)
995        elif datatype == "pickle":
996            dump_object(self, filename, **kwargs)
997        else:
998            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):
1000    def print(self, print_range=None):
1001        print(self.__repr__(print_range))
def sqrt(self):
1210    def sqrt(self):
1211        return self ** 0.5
def log(self):
1213    def log(self):
1214        newcontent = [None if _check_for_none(self, item) else np.log(item) for item in self.content]
1215        return Corr(newcontent, prange=self.prange)
def exp(self):
1217    def exp(self):
1218        newcontent = [None if _check_for_none(self, item) else np.exp(item) for item in self.content]
1219        return Corr(newcontent, prange=self.prange)
def sin(self):
1234    def sin(self):
1235        return self._apply_func_to_corr(np.sin)
def cos(self):
1237    def cos(self):
1238        return self._apply_func_to_corr(np.cos)
def tan(self):
1240    def tan(self):
1241        return self._apply_func_to_corr(np.tan)
def sinh(self):
1243    def sinh(self):
1244        return self._apply_func_to_corr(np.sinh)
def cosh(self):
1246    def cosh(self):
1247        return self._apply_func_to_corr(np.cosh)
def tanh(self):
1249    def tanh(self):
1250        return self._apply_func_to_corr(np.tanh)
def arcsin(self):
1252    def arcsin(self):
1253        return self._apply_func_to_corr(np.arcsin)
def arccos(self):
1255    def arccos(self):
1256        return self._apply_func_to_corr(np.arccos)
def arctan(self):
1258    def arctan(self):
1259        return self._apply_func_to_corr(np.arctan)
def arcsinh(self):
1261    def arcsinh(self):
1262        return self._apply_func_to_corr(np.arcsinh)
def arccosh(self):
1264    def arccosh(self):
1265        return self._apply_func_to_corr(np.arccosh)
def arctanh(self):
1267    def arctanh(self):
1268        return self._apply_func_to_corr(np.arctanh)
real
imag
def prune(self, Ntrunc, tproj=3, t0proj=2, basematrix=None):
1303    def prune(self, Ntrunc, tproj=3, t0proj=2, basematrix=None):
1304        r''' Project large correlation matrix to lowest states
1305
1306        This method can be used to reduce the size of an (N x N) correlation matrix
1307        to (Ntrunc x Ntrunc) by solving a GEVP at very early times where the noise
1308        is still small.
1309
1310        Parameters
1311        ----------
1312        Ntrunc: int
1313            Rank of the target matrix.
1314        tproj: int
1315            Time where the eigenvectors are evaluated, corresponds to ts in the GEVP method.
1316            The default value is 3.
1317        t0proj: int
1318            Time where the correlation matrix is inverted. Choosing t0proj=1 is strongly
1319            discouraged for O(a) improved theories, since the correctness of the procedure
1320            cannot be granted in this case. The default value is 2.
1321        basematrix : Corr
1322            Correlation matrix that is used to determine the eigenvectors of the
1323            lowest states based on a GEVP. basematrix is taken to be the Corr itself if
1324            is is not specified.
1325
1326        Notes
1327        -----
1328        We have the basematrix $C(t)$ and the target matrix $G(t)$. We start by solving
1329        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}$
1330        and $t_0 \equiv t_{0, \mathrm{proj}}$. The target matrix is projected onto the subspace of the
1331        resulting eigenvectors $v_n, n=1,\dots,N_\mathrm{trunc}$ via
1332        $$G^\prime_{i, j}(t) = (v_i, G(t) v_j)$$. This allows to reduce the size of a large
1333        correlation matrix and to remove some noise that is added by irrelevant operators.
1334        This may allow to use the GEVP on $G(t)$ at late times such that the theoretically motivated
1335        bound $t_0 \leq t/2$ holds, since the condition number of $G(t)$ is decreased, compared to $C(t)$.
1336        '''
1337
1338        if self.N == 1:
1339            raise Exception('Method cannot be applied to one-dimensional correlators.')
1340        if basematrix is None:
1341            basematrix = self
1342        if Ntrunc >= basematrix.N:
1343            raise Exception('Cannot truncate using Ntrunc <= %d' % (basematrix.N))
1344        if basematrix.N != self.N:
1345            raise Exception('basematrix and targetmatrix have to be of the same size.')
1346
1347        evecs = basematrix.GEVP(t0proj, tproj, sort=None)[:Ntrunc]
1348
1349        tmpmat = np.empty((Ntrunc, Ntrunc), dtype=object)
1350        rmat = []
1351        for t in range(basematrix.T):
1352            for i in range(Ntrunc):
1353                for j in range(Ntrunc):
1354                    tmpmat[i][j] = evecs[i].T @ self[t] @ evecs[j]
1355            rmat.append(np.copy(tmpmat))
1356
1357        newcontent = [None if (self.content[t] is None) else rmat[t] for t in range(self.T)]
1358        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)$.

N