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