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