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