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