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