Rivet  3.1.0
Correlators.hh
1 // -*- C++ -*-
2 #ifndef RIVET_Correlators_HH
3 #define RIVET_Correlators_HH
4 
5 // Tools for calculating flow coefficients using correlators.
6 // Classes:
7 // Correlators: Calculates single event correlators of a given harmonic.
8 // Cumulants: An additional base class for flow analyses
9 // Use as: class MyAnalysis : public Analysis, Cumulants {};
10 // Includes a framework for calculating cumulants and flow coefficients
11 // from single event correlators, including automatic handling of
12 // statistical errors. Contains multiple internal sub-classes:
13 // CorBinBase: Base class for correlators binned in event or particle observables.
14 // CorSingleBin: A simple bin for correlators.
15 // CorBin: Has the interface of a simple bin, but does automatic calculation
16 // of statistical errors by a bootstrap method.
17 // ECorrelator: Data type for event averaged correlators.
18 //
19 // Authors: Vytautas Vislavicius, Christine O. Rasmussen, Christian Bierlich.
20 
21 #include "Rivet/Analysis.hh"
22 #include "Rivet/Projection.hh"
23 #include "Rivet/Projections/ParticleFinder.hh"
24 #include "YODA/Scatter2D.h"
25 #include <complex>
26 
27 namespace Rivet {
28  using std::complex;
29 
30 
38  class Correlators : public Projection {
39  public:
40 
54  Correlators(const ParticleFinder& fsp, int nMaxIn = 2,
55  int pMaxIn = 0, vector<double> pTbinEdgesIn = {});
56 
57  // Constructor which takes a Scatter2D to estimate bin edges.
58  Correlators(const ParticleFinder& fsp, int nMaxIn,
59  int pMaxIn, const Scatter2DPtr hIn);
60 
68  const pair<double,double> intCorrelator(vector<int> n) const;
69 
74  const vector<pair<double,double>> pTBinnedCorrelators(vector<int> n,
75  bool overflow = false) const;
76 
85  const pair<double,double> intCorrelatorGap(const Correlators& other,
86  vector<int> n1, vector<int> n2) const;
87 
99  const vector<pair<double,double>>
100  pTBinnedCorrelatorsGap(const Correlators& other, vector<int> n1, vector<int> n2, bool overflow=false) const;
101 
106  static vector<int> hVec(int n, int m) {
107  if (m % 2 != 0) {
108  cout << "Harmonic Vector: Number of particles must be an even number." << endl;
109  return {};
110  }
111  vector<int> ret = {};
112  for (int i = 0; i < m; ++i) {
113  if (i < m/2) ret.push_back(n);
114  else ret.push_back(-n);
115  }
116  return ret;
117  }
118 
121  static pair<int,int> getMaxValues(vector< vector<int> >& hList){
122  int nMax = 0, pMax = 0;
123  for (vector<int> h : hList) {
124  int tmpN = 0, tmpP = 0;
125  for ( int i = 0; i < int(h.size()); ++i) {
126  tmpN += abs(h[i]);
127  ++tmpP;
128  }
129  if (tmpN > nMax) nMax = tmpN;
130  if (tmpP > pMax) pMax = tmpP;
131  }
132  return make_pair(nMax,pMax);
133  }
134 
135  // Clone on the heap.
136  DEFAULT_RIVET_PROJ_CLONE(Correlators);
137 
138 
139  protected:
140 
141  // @brief Loop over array and calculates Q and P vectors if needed
142  void project(const Event& e);
143 
144  // @brief Compare to other projection, testing harmonics, pT bins and underlying final state similarity
145  CmpState compare(const Projection& p) const {
146  const Correlators* other = dynamic_cast<const Correlators*>(&p);
147  if (nMax != other->nMax) return CmpState::NEQ;
148  if (pMax != other->pMax) return CmpState::NEQ;
149  if (pTbinEdges != other->pTbinEdges) return CmpState::NEQ;
150  return mkPCmp(p, "FS");
151  }
152 
153  // @brief Calculate correlators from one particle
154  void fillCorrelators(const Particle& p, const double& weight);
155 
156  // @brief Return a Q-vector.
157  const complex<double> getQ(int n, int p) const {
158  bool isNeg = (n < 0);
159  if (isNeg) return conj( qVec[abs(n)][p] );
160  else return qVec[n][p];
161  }
162 
163  // @brief Return a P-vector
164  const complex<double> getP(int n, int p, double pT = 0.) const {
165  bool isNeg = (n < 0);
166  map<double, Vec2D>::const_iterator pTitr = pVec.lower_bound(pT);
167  if (pTitr == pVec.end()) return DBL_NAN;
168  if (isNeg) return conj( pTitr->second[abs(n)][p] );
169  else return pTitr->second[n][p];
170  }
171 
172 
173  private:
174 
175  // @brief Find correlators by recursion
176  //
177  // Order = M (# of particles), n's are harmonics, p's are the powers of the weights
178  const complex<double> recCorr(int order, vector<int> n,
179  vector<int> p, bool useP, double pT = 0.) const;
180 
181  // @brief Two-particle correlator
182  //
183  // Cf. eq. (19) p6. Flag if p-vectors or q-vectors should be used to
184  // calculate the correlator.
185  const complex<double> twoPartCorr(int n1, int n2, int p1 = 1,
186  int p2 = 1, double pT = 0., bool useP = false) const;
187 
188  // Set elements in vectors to zero
189  void setToZero();
190 
191  // Shorthands for setting and comparing to zero
192  const complex<double> _ZERO = {0., 0.};
193  const double _TINY = 1e-10;
194 
195  // Shorthand typedefs for vec<vec<complex>>.
196  typedef vector< vector<complex<double>> > Vec2D;
197 
198  // Define Q-vectors and P-vectors
199  Vec2D qVec; // Q[n][p]
200  map<double, Vec2D> pVec; // p[pT][n][p]
201 
202  // The max values of n and p to be calculated.
203  int nMax, pMax;
204 
205  // pT bin-edges
206  vector<double> pTbinEdges;
207 
208  bool isPtDiff;
209 
210  };
211 
212 
213 
221  class CumulantAnalysis : public Analysis {
222  private:
223 
224  // Number of bins used for bootstrap calculation of statistical
225  // uncertainties. It is hard coded, and shout NOT be changed unless there
226  // are good reasons to do so.
227  static const int BOOT_BINS = 9;
228 
229  // Enum for choosing the method of error calculation.
230  enum Error {
231  VARIANCE,
232  ENVELOPE
233  };
234 
235  // The desired error method. Can be changed in the analysis constructor
236  // by setting it appropriately.
237  Error errorMethod;
238 
239 
241  class CorBinBase {
242  public:
243  CorBinBase() {}
244  virtual ~CorBinBase() {};
245  // Derived class should have fill and mean defined.
246  virtual void fill(const pair<double, double>& cor, const double& weight = 1.0) = 0;
247  virtual const double mean() const = 0;
248  };
249 
250 
256  class CorSingleBin : public CorBinBase {
257  public:
258 
260  CorSingleBin()
261  : _sumWX(0.), _sumW(0.), _sumW2(0.), _numEntries(0.)
262  { }
263 
264  // Destructor does nothing but must be implemented (?)
265  ~CorSingleBin() {}
266 
270  void fill(const pair<double, double>& cor, const double& weight = 1.0) {
271  // Test if denominator for the single event average is zero.
272  if (cor.second < 1e-10) return;
273  // The single event average <M> is then cor.first / cor.second.
274  // With weights this becomes just:
275  _sumWX += cor.first * weight;
276  _sumW += weight * cor.second;
277  _sumW2 += weight * weight * cor.second * cor.second;
278  _numEntries += 1.;
279  }
280 
282  const double mean() const {
283  if (_sumW < 1e-10) return 0;
284  return _sumWX / _sumW;
285  }
286 
288  const double sumW() const {
289  return _sumW;
290  }
291 
293  const double sumW2() const {
294  return _sumW2;
295  }
296 
298  const double sumWX() const {
299  return _sumWX;
300  }
301 
303  const double numEntries() const {
304  return _numEntries;
305  }
306 
308  void addContent(double ne, double sw, double sw2, double swx) {
309  _numEntries += ne;
310  _sumW += sw;
311  _sumW2 += sw2;
312  _sumWX += swx;
313  }
314 
315 
316  private:
317 
318  double _sumWX, _sumW, _sumW2, _numEntries;
319 
320  };
321 
322 
327  class CorBin : public CorBinBase {
328  public:
329 
335  CorBin() : binIndex(0), nBins(BOOT_BINS) {
336  for(size_t i = 0; i < nBins; ++i) bins.push_back(CorSingleBin());
337  }
338 
339  // Destructor does nothing but must be implemented (?)
340  ~CorBin() {}
341 
343  void fill(const pair<double, double>& cor, const double& weight = 1.0) {
344  // Test if denominator for the single event average is zero.
345  if (cor.second < 1e-10) return;
346  // Fill the correct bin.
347  bins[binIndex].fill(cor, weight);
348  if (binIndex == nBins - 1) binIndex = 0;
349  else ++binIndex;
350  }
351 
353  const double mean() const {
354  double sow = 0;
355  double sowx = 0;
356  for (auto b : bins) {
357  if (b.sumW() < 1e-10) continue;
358  sow += b.sumW();
359  sowx += b.sumWX();
360  }
361  return sowx / sow;
362  }
363 
365  const vector<CorSingleBin> getBins() const {
366  return bins;
367  }
368 
370  template<class T=CorBinBase>
371  const vector<T*> getBinPtrs() {
372  vector<T*> ret(bins.size());
373  transform(bins.begin(), bins.end(), ret.begin(),
374  [](CorSingleBin& b) {return &b;});
375  return ret;
376  }
377 
378  private:
379 
380  vector<CorSingleBin> bins;
381  size_t binIndex;
382  size_t nBins;
383 
384  };
385 
386 
387  public:
388 
393  class ECorrelator {
394  public:
395 
401  //ECorrelator(vector<int> h) : h1(h), h2({}),
402  // binX(0), binContent(0), reference() {
403  //};
404 
410  ECorrelator(vector<int> h, vector<double> binIn)
411  : h1(h), h2({}), binX(binIn), binContent(binIn.size() - 1), reference()
412  { }
413 
418  ECorrelator(vector<int> h1In, vector<int> h2In, vector<double> binIn)
419  : h1(h1In), h2(h2In), binX(binIn), binContent(binIn.size() - 1), reference()
420  { }
421 
423  void fill(const double& obs, const Correlators& c, double weight=1.0) {
424  int index = getBinIndex(obs);
425  if (index < 0) return;
426  binContent[index].fill(c.intCorrelator(h1), weight);
427  }
428 
432  void fill(const double& obs, const Correlators& c1, const Correlators& c2, double weight=1.0) {
433  if (!h2.size()) {
434  cout << "Trying to fill gapped correlator, but harmonics behind the gap (h2) are not given!" << endl;
435  return;
436  }
437  int index = getBinIndex(obs);
438  if (index < 0) return;
439  binContent[index].fill(c1.intCorrelatorGap(c2, h1, h2), weight);
440  }
441 
446  void fill(const Correlators& c, const double& weight=1.0) {
447  vector< pair<double, double> > diffCorr = c.pTBinnedCorrelators(h1);
448  // We always skip overflow when calculating the all-event average.
449  if (diffCorr.size() != binX.size() - 1)
450  cout << "Tried to fill event with wrong binning (ungapped)" << endl;
451  for (size_t i = 0; i < diffCorr.size(); ++i) {
452  int index = getBinIndex(binX[i]);
453  if (index < 0) return;
454  binContent[index].fill(diffCorr[i], weight);
455  }
456  reference.fill(c.intCorrelator(h1), weight);
457  }
458 
463  void fill(const Correlators& c1, const Correlators& c2, const double& weight = 1.0) {
464  if (!h2.size()) {
465  cout << "Trying to fill gapped correlator, but harmonics behind "
466  "the gap (h2) are not given!" << endl;
467  return;
468  }
469  vector< pair<double, double> > diffCorr = c1.pTBinnedCorrelatorsGap(c2, h1, h2);
470  // We always skip overflow when calculating the all event average.
471  if (diffCorr.size() != binX.size() - 1)
472  cout << "Tried to fill event with wrong binning (gapped)" << endl;
473  for (size_t i = 0; i < diffCorr.size(); ++i) {
474  int index = getBinIndex(binX[i]);
475  if (index < 0) return;
476  binContent[index].fill(diffCorr[i], weight);
477  }
478  reference.fill(c1.intCorrelatorGap(c2, h1, h2), weight);
479  }
480 
482  const vector<CorBin> getBins() const {
483  return binContent;
484  }
485 
487  const vector<CorBinBase*> getBinPtrs() {
488  vector<CorBinBase*> ret(binContent.size());
489  transform(binContent.begin(), binContent.end(), ret.begin(), [](CorBin& b) {return &b;});
490  return ret;
491  }
492 
494  const vector<double> getBinX() const {
495  return binX;
496  }
497 
499  const vector<int> getH1() const {
500  return h1;
501  }
502 
504  const vector<int> getH2() const {
505  return h2;
506  }
507 
509  void setReference(CorBin refIn) {
510  reference = refIn;
511  }
512 
514  const CorBin getReference() const {
515  if (reference.mean() < 1e-10)
516  cout << "Warning: ECorrelator, reference bin is zero." << endl;
517  return reference;
518  }
519 
523  void setProfs(list<Profile1DPtr> prIn) {
524  profs = prIn;
525  }
526 
528  void fillFromProfs() {
529  list<Profile1DPtr>::iterator hItr = profs.begin();
530  auto refs = reference.getBinPtrs<CorSingleBin>();
531  for (size_t i = 0; i < profs.size(); ++i, ++hItr) {
532  for (size_t j = 0; j < binX.size() - 1; ++j) {
533  const YODA::ProfileBin1D& pBin = (*hItr)->binAt(binX[j]);
534  auto tmp = binContent[j].getBinPtrs<CorSingleBin>();
535  tmp[i]->addContent(pBin.numEntries(), pBin.sumW(), pBin.sumW2(),
536  pBin.sumWY());
537  }
538  // Get the reference flow from the underflow bin of the histogram.
539  const YODA::Dbn2D& uBin = (*hItr)->underflow();
540  refs[i]->addContent(uBin.numEntries(), uBin.sumW(), uBin.sumW2(),
541  uBin.sumWY());
542  } // End loop of bootstrapped correlators.
543  }
544 
546  list<Profile1DPtr>::iterator profBegin() {
547  return profs.begin();
548  }
549 
551  list<Profile1DPtr>::iterator profEnd() {
552  return profs.end();
553  }
554 
555 
556  private:
557 
558  // Get correct bin index for a given @param obs value
559  const int getBinIndex(const double& obs) const {
560  // Find the correct index of binContent.
561  // If we are in overflow, just skip.
562  if (obs >= binX.back()) return -1;
563  // If we are in underflow, ditto.
564  if (obs < binX[0]) return -1;
565  int index = 0;
566  for (int i = 0, N = binX.size() - 1; i < N; ++i, ++index)
567  if (obs >= binX[i] && obs < binX[i + 1]) break;
568  return index;
569  }
570 
571  // The harmonics vectors.
572  vector<int> h1;
573  vector<int> h2;
574 
575  // The bins.
576  vector<double> binX;
577  vector<CorBin> binContent;
578  // The reference flow.
579  CorBin reference;
580  // The profile histograms associated with the CorBins for streaming.
581  list<Profile1DPtr> profs;
582 
583  };
584 
585 
587  const pair<int, int> getMaxValues() const {
588  vector< vector<int>> harmVecs;
589  for ( auto eItr = eCorrPtrs.begin(); eItr != eCorrPtrs.end(); ++eItr) {
590  vector<int> h1 = (*eItr)->getH1();
591  vector<int> h2 = (*eItr)->getH2();
592  if (h1.size() > 0) harmVecs.push_back(h1);
593  if (h2.size() > 0) harmVecs.push_back(h2);
594  }
595  if (harmVecs.size() == 0) {
596  cout << "Warning: You tried to extract max values from harmonic "
597  "vectors, but have not booked any." << endl;
598  return pair<int,int>();
599  }
600  return Correlators::getMaxValues(harmVecs);
601  }
602 
604  typedef shared_ptr<ECorrelator> ECorrPtr;
605 
606 
609  ECorrPtr bookECorrelator(const string name, const vector<int>& h, const Scatter2DPtr hIn) {
610  vector<double> binIn;
611  for (auto b : hIn->points()) binIn.push_back(b.xMin());
612  binIn.push_back(hIn->points().back().xMax());
613  ECorrPtr ecPtr = ECorrPtr(new ECorrelator(h, binIn));
614  list<Profile1DPtr> eCorrProfs;
615  for (int i = 0; i < BOOT_BINS; ++i) {
616  Profile1DPtr tmp;
617  book(tmp, name+"-"+to_string(i),*hIn);
618  tmp->setPath(this->name()+"/FINAL/" + name+"-"+to_string(i));
619  //tmp->setPath(tmp->path()+"FINAL");
620  eCorrProfs.push_back(tmp);
621  }
622  ecPtr->setProfs(eCorrProfs);
623  eCorrPtrs.push_back(ecPtr);
624  return ecPtr;
625  }
626 
629  ECorrPtr bookECorrelator(const string name, const vector<int>& h1,
630  const vector<int>& h2, const Scatter2DPtr hIn ) {
631  vector<double> binIn;
632  for (auto b : hIn->points()) binIn.push_back(b.xMin());
633  binIn.push_back(hIn->points().back().xMax());
634  ECorrPtr ecPtr = ECorrPtr(new ECorrelator(h1, h2, binIn));
635  list<Profile1DPtr> eCorrProfs;
636  for (int i = 0; i < BOOT_BINS; ++i) {
637  Profile1DPtr tmp;
638  book(tmp, name+"-"+to_string(i),*hIn);
639  tmp->setPath(this->name()+"/FINAL/" + name+"-"+to_string(i));
640  //tmp->setPath(tmp->path()+"FINAL");
641  eCorrProfs.push_back(tmp);
642  }
643  ecPtr->setProfs(eCorrProfs);
644  eCorrPtrs.push_back(ecPtr);
645  return ecPtr;
646  }
647 
651  ECorrPtr bookECorrelatorGap (const string name, const vector<int>& h,
652  const Scatter2DPtr hIn) {
653  const vector<int> h1(h.begin(), h.begin() + h.size() / 2);
654  const vector<int> h2(h.begin() + h.size() / 2, h.end());
655  return bookECorrelator(name, h1, h2, hIn);
656  }
657 
662  template<unsigned int N, unsigned int M>
663  ECorrPtr bookECorrelator(const string name, const Scatter2DPtr hIn) {
664  return bookECorrelator(name, Correlators::hVec(N, M), hIn);
665  }
666 
671  template<unsigned int N, unsigned int M>
672  ECorrPtr bookECorrelatorGap(const string name, const Scatter2DPtr hIn) {
673  return bookECorrelatorGap(name, Correlators::hVec(N, M), hIn);
674  }
675 
676 
682  void stream() {
683  for (auto ecItr = eCorrPtrs.begin(); ecItr != eCorrPtrs.end(); ++ecItr){
684  (*ecItr)->fillFromProfs();
685  corrPlot(list<Profile1DPtr>((*ecItr)->profBegin(), (*ecItr)->profEnd()), *ecItr);
686  }
687  }
688 
689 
690  private:
691 
692  // Bookkeeping of the event averaged correlators.
693  list<ECorrPtr> eCorrPtrs;
694 
695 
696  public:
697 
702  CumulantAnalysis(const string& n)
703  : Analysis(n), errorMethod(VARIANCE)
704  { }
705 
713  template<typename T>
714  static void fillScatter(Scatter2DPtr h, vector<double>& binx, T func) {
715  vector<YODA::Point2D> points;
716  for (int i = 0, N = binx.size() - 1; i < N; ++i) {
717  double xMid = (binx[i] + binx[i + 1]) / 2.0;
718  double xeMin = fabs(xMid - binx[i]);
719  double xeMax = fabs(xMid - binx[i + 1]);
720  double yVal = func(i);
721  if (std::isnan(yVal)) yVal = 0.;
722  double yErr = 0;
723  points.push_back(YODA::Point2D(xMid, yVal, xeMin, xeMax, yErr, yErr));
724  }
725  h->reset();
726  h->points().clear();
727  for (int i = 0, N = points.size(); i < N; ++i) h->addPoint(points[i]);
728  }
729 
737  template<typename F>
738  const void fillScatter(Scatter2DPtr h, vector<double>& binx, F func,
739  vector<pair<double, double> >& yErr) const {
740  vector<YODA::Point2D> points;
741  for (int i = 0, N = binx.size() - 1; i < N; ++i) {
742  double xMid = (binx[i] + binx[i + 1]) / 2.0;
743  double xeMin = fabs(xMid - binx[i]);
744  double xeMax = fabs(xMid - binx[i + 1]);
745  double yVal = func(i);
746  if (std::isnan(yVal))
747  points.push_back(YODA::Point2D(xMid, 0., xeMin, xeMax,0., 0.));
748  else
749  points.push_back(YODA::Point2D(xMid, yVal, xeMin, xeMax,
750  yErr[i].first, yErr[i].second));
751  }
752  h->reset();
753  h->points().clear();
754  for (int i = 0, N = points.size(); i < N; ++i)
755  h->addPoint(points[i]);
756  }
757 
758 
762  static void nthPow(Scatter2DPtr hOut, const Scatter2DPtr hIn, const double& n, const double& k = 1.0) {
763  if (n == 0 || n == 1) {
764  cout << "Error: Do not take the 0th or 1st power of a Scatter2D,"
765  " use scale instead." << endl;
766  return;
767  }
768  if (hIn->points().size() != hOut->points().size()) {
769  cout << "nthRoot: Scatterplots: " << hIn->name() << " and " <<
770  hOut->name() << " not initialized with same length." << endl;
771  return;
772  }
773  vector<YODA::Point2D> points;
774  // The error pre-factor is k^(1/n) / n by Taylors formula.
775  double eFac = pow(k,1./n)/n;
776  for (auto b : hIn->points()) {
777  double yVal = pow(k * b.y(),n);
778  if (std::isnan(yVal))
779  points.push_back(YODA::Point2D(b.x(), 0., b.xErrMinus(),
780  b.xErrPlus(), 0, 0 ));
781  else {
782  double yemin = abs(eFac * pow(yVal,1./(n - 1.))) * b.yErrMinus();
783  if (std::isnan(yemin)) yemin = b.yErrMinus();
784  double yemax = abs(eFac * pow(yVal,1./(n - 1.))) * b.yErrPlus();
785  if (std::isnan(yemax)) yemax = b.yErrPlus();
786  points.push_back(YODA::Point2D(b.x(), yVal, b.xErrMinus(),
787  b.xErrPlus(), yemin, yemax ));
788  }
789  }
790  hOut->reset();
791  hOut->points().clear();
792  for (int i = 0, N = points.size(); i < N; ++i)
793  hOut->addPoint(points[i]);
794  }
795 
796 
800  static void nthPow(Scatter2DPtr h, const double& n, const double& k = 1.0) {
801  if (n == 0 || n == 1) {
802  cout << "Error: Do not take the 0th or 1st power of a Scatter2D,"
803  " use scale instead." << endl;
804  return;
805  }
806  vector<YODA::Point2D> points;
807  vector<YODA::Point2D> pIn = h->points();
808  // The error pre-factor is k^(1/n) / n by Taylors formula.
809  double eFac = pow(k,1./n)/n;
810  for (auto b : pIn) {
811  double yVal = pow(k * b.y(),n);
812  if (std::isnan(yVal))
813  points.push_back(YODA::Point2D(b.x(), 0., b.xErrMinus(),
814  b.xErrPlus(), 0, 0 ));
815  else {
816  double yemin = abs(eFac * pow(yVal,1./(n - 1.))) * b.yErrMinus();
817  if (std::isnan(yemin)) yemin = b.yErrMinus();
818  double yemax = abs(eFac * pow(yVal,1./(n - 1.))) * b.yErrPlus();
819  if (std::isnan(yemax)) yemax = b.yErrPlus();
820  points.push_back(YODA::Point2D(b.x(), yVal, b.xErrMinus(),
821  b.xErrPlus(), yemin, yemax ));
822  }
823  }
824  h->reset();
825  h->points().clear();
826  for (int i = 0, N = points.size(); i < N; ++i) h->addPoint(points[i]);
827  }
828 
829 
834  template<typename T>
835  static pair<double, double> sampleVariance(T func) {
836  // First we calculate the mean (two pass calculation).
837  double avg = 0.;
838  for (int i = 0; i < BOOT_BINS; ++i) avg += func(i);
839  avg /= double(BOOT_BINS);
840  // Then we find the variance.
841  double var = 0.;
842  for (int i = 0; i < BOOT_BINS; ++i) var += pow(func(i) - avg, 2.);
843  var /= (double(BOOT_BINS) - 1);
844  return pair<double, double>(var, var);
845  }
846 
847 
852  template<typename T>
853  static pair<double, double> sampleEnvelope(T func) {
854  // First we calculate the mean.
855  double avg = 0.;
856  for (int i = 0; i < BOOT_BINS; ++i) avg += func(i);
857  avg /= double(BOOT_BINS);
858  double yMax = avg;
859  double yMin = avg;
860  // Then we find the envelope using the mean as initial value.
861  for (int i = 0; i < BOOT_BINS; ++i) {
862  double yVal = func(i);
863  if (yMin > yVal) yMin = yVal;
864  else if (yMax < yVal) yMax = yVal;
865  }
866  return pair<double, double>(fabs(avg - yMin), fabs(yMax - avg));
867  }
868 
869 
871  template<typename T>
872  const pair<double, double> sampleError(T func) const {
873  if (errorMethod == VARIANCE) return sampleVariance(func);
874  else if (errorMethod == ENVELOPE) return sampleEnvelope(func);
875  else
876  cout << "Error: Error method not found!" << endl;
877  return pair<double, double>(0.,0.);
878  }
879 
880 
882  const void cnTwoInt(Scatter2DPtr h, ECorrPtr e2) const {
883  vector<CorBin> bins = e2->getBins();
884  vector<double> binx = e2->getBinX();
885  // Assert bin size.
886  if (binx.size() - 1 != bins.size()){
887  cout << "cnTwoInt: Bin size (x,y) differs!" << endl;
888  return;
889  }
890  vector<CorBinBase*> binPtrs;
891  // The mean value of the cumulant.
892  auto cn = [&] (int i) { return binPtrs[i]->mean(); };
893  // Error calculation.
894  vector<pair<double, double> > yErr;
895  for (int j = 0, N = bins.size(); j < N; ++j) {
896  binPtrs = bins[j].getBinPtrs();
897  yErr.push_back(sampleError(cn));
898  }
899  binPtrs = e2->getBinPtrs();
900  fillScatter(h, binx, cn, yErr);
901  }
902 
903 
905  const void vnTwoInt(Scatter2DPtr h, ECorrPtr e2) const {
906  cnTwoInt(h, e2);
907  nthPow(h, 0.5);
908  }
909 
910 
914  const void corrPlot(Scatter2DPtr h, ECorrPtr e) const {
915  cnTwoInt(h, e);
916  }
917 
918 
920  const void corrPlot(list<Profile1DPtr> profs, ECorrPtr e) const {
921  vector<CorBin> corBins = e->getBins();
922  vector<double> binx = e->getBinX();
923  auto ref = e->getReference();
924  auto refBins = ref.getBinPtrs<CorSingleBin>();
925  // Assert bin size.
926  if (binx.size() - 1 != corBins.size()){
927  cout << "corrPlot: Bin size (x,y) differs!" << endl;
928  return;
929  }
930  list<Profile1DPtr>::iterator hItr = profs.begin();
931  // Loop over the boostrapped correlators.
932  for (size_t i = 0; i < profs.size(); ++i, ++hItr) {
933  vector<YODA::ProfileBin1D> profBins;
934  // Numbers for the summary distribution
935  double ne = 0., sow = 0., sow2 = 0.;
936  for (size_t j = 0, N = binx.size() - 1; j < N; ++j) {
937  vector<CorSingleBin*> binPtrs =
938  corBins[j].getBinPtrs<CorSingleBin>();
939  // Construct bin of the profiled quantities. We have no information
940  // (and no desire to add it) of sumWX of the profile, so really
941  // we should use a Dbn1D - but that does not work for Profile1D's.
942  profBins.push_back( YODA::ProfileBin1D((*hItr)->bin(j).xEdges(),
943  YODA::Dbn2D( binPtrs[i]->numEntries(), binPtrs[i]->sumW(),
944  binPtrs[i]->sumW2(), 0., 0., binPtrs[i]->sumWX(), 0, 0)));
945  ne += binPtrs[i]->numEntries();
946  sow += binPtrs[i]->sumW();
947  sow2 += binPtrs[i]->sumW2();
948  }
949  // Reset the bins of the profiles.
950  (*hItr)->reset();
951  (*hItr)->bins().clear();
952  // Add our new bins.
953  // The total distribution
954  (*hItr)->setTotalDbn(YODA::Dbn2D(ne,sow,sow2,0.,0.,0.,0.,0.));
955  // The bins.
956  for (int j = 0, N = profBins.size(); j < N; ++j)
957  (*hItr)->addBin(profBins[j]);
958  // The reference flow in the underflow bin.
959  (*hItr)->setUnderflow(YODA::Dbn2D(refBins[i]->numEntries(),
960  refBins[i]->sumW(), refBins[i]->sumW2(), 0., 0.,
961  refBins[i]->sumWX(), 0., 0.));
962  } // End loop of bootstrapped correlators.
963  }
964  // @brief Four particle integrated cn.
965  const void cnFourInt(Scatter2DPtr h, ECorrPtr e2, ECorrPtr e4) const {
966  auto e2bins = e2->getBins();
967  auto e4bins = e4->getBins();
968  auto binx = e2->getBinX();
969  if (binx.size() - 1 != e2bins.size()){
970  cout << "cnFourInt: Bin size (x,y) differs!" << endl;
971  return;
972  }
973  if (binx != e4->getBinX()) {
974  cout << "Error in cnFourInt: Correlator x-binning differs!" << endl;
975  return;
976  }
977  vector<CorBinBase*> e2binPtrs;
978  vector<CorBinBase*> e4binPtrs;
979  auto cn = [&] (int i) {
980  double e22 = e2binPtrs[i]->mean() * e2binPtrs[i]->mean();
981  return e4binPtrs[i]->mean() - 2. * e22;
982  };
983  // Error calculation.
984  vector<pair<double, double> > yErr;
985  for (int j = 0, N = e2bins.size(); j < N; ++j) {
986  e2binPtrs = e2bins[j].getBinPtrs();
987  e4binPtrs = e4bins[j].getBinPtrs();
988  yErr.push_back(sampleError(cn));
989  }
990  // Put the bin ptrs back in place.
991  e2binPtrs = e2->getBinPtrs();
992  e4binPtrs = e4->getBinPtrs();
993  fillScatter(h, binx, cn, yErr);
994  }
995 
996 
998  const void vnFourInt(Scatter2DPtr h, ECorrPtr e2, ECorrPtr e4) const {
999  cnFourInt(h, e2, e4);
1000  nthPow(h, 0.25, -1.0);
1001  }
1002 
1003 
1005  const void cnSixInt(Scatter2DPtr h, ECorrPtr e2, ECorrPtr e4,
1006  ECorrPtr e6) const {
1007  auto e2bins = e2->getBins();
1008  auto e4bins = e4->getBins();
1009  auto e6bins = e6->getBins();
1010  auto binx = e2->getBinX();
1011  if (binx.size() - 1 != e2bins.size()){
1012  cout << "cnSixInt: Bin size (x,y) differs!" << endl;
1013  return;
1014  }
1015  if (binx != e4->getBinX() || binx != e6->getBinX()) {
1016  cout << "Error in cnSixInt: Correlator x-binning differs!" << endl;
1017  return;
1018  }
1019  vector<CorBinBase*> e2binPtrs;
1020  vector<CorBinBase*> e4binPtrs;
1021  vector<CorBinBase*> e6binPtrs;
1022  auto cn = [&] (int i) {
1023  double e23 = pow(e2binPtrs[i]->mean(), 3.0);
1024  return e6binPtrs[i]->mean() - 9.*e2binPtrs[i]->mean()*e4binPtrs[i]->mean() +
1025  12.*e23;
1026  };
1027  // Error calculation.
1028  vector<pair<double, double> > yErr;
1029  for (int j = 0, N = e2bins.size(); j < N; ++j) {
1030  e2binPtrs = e2bins[j].getBinPtrs();
1031  e4binPtrs = e4bins[j].getBinPtrs();
1032  e6binPtrs = e6bins[j].getBinPtrs();
1033  yErr.push_back(sampleError(cn));
1034  }
1035  // Put the bin ptrs back in place.
1036  e2binPtrs = e2->getBinPtrs();
1037  e4binPtrs = e4->getBinPtrs();
1038  e6binPtrs = e6->getBinPtrs();
1039  fillScatter(h, binx, cn, yErr);
1040  }
1041 
1042 
1044  const void vnSixInt(Scatter2DPtr h, ECorrPtr e2, ECorrPtr e4,
1045  ECorrPtr e6) const {
1046  cnSixInt(h, e2, e4, e6);
1047  nthPow(h, 1./6., 0.25);
1048  }
1049 
1050 
1052  const void cnEightInt(Scatter2DPtr h, ECorrPtr e2, ECorrPtr e4,
1053  ECorrPtr e6, ECorrPtr e8) const {
1054  auto e2bins = e2->getBins();
1055  auto e4bins = e4->getBins();
1056  auto e6bins = e6->getBins();
1057  auto e8bins = e8->getBins();
1058  auto binx = e2->getBinX();
1059  if (binx.size() - 1 != e2bins.size()){
1060  cout << "cnEightInt: Bin size (x,y) differs!" << endl;
1061  return;
1062  }
1063  if (binx != e4->getBinX() || binx != e6->getBinX() ||
1064  binx != e8->getBinX()) {
1065  cout << "Error in cnEightInt: Correlator x-binning differs!" << endl;
1066  return;
1067  }
1068  vector<CorBinBase*> e2binPtrs;
1069  vector<CorBinBase*> e4binPtrs;
1070  vector<CorBinBase*> e6binPtrs;
1071  vector<CorBinBase*> e8binPtrs;
1072  auto cn = [&] (int i ) {
1073  double e22 = e2binPtrs[i]->mean() * e2binPtrs[i]->mean();
1074  double e24 = e22 * e22;
1075  double e42 = e4binPtrs[i]->mean() * e4binPtrs[i]->mean();
1076  return e8binPtrs[i]->mean() - 16. * e6binPtrs[i]->mean() *
1077  e2binPtrs[i]->mean() - 18. * e42 + 144. * e4binPtrs[i]->mean()*e22
1078  - 144. * e24;
1079  };
1080  // Error calculation.
1081  vector<pair<double, double> > yErr;
1082  for (int j = 0, N = e2bins.size(); j < N; ++j) {
1083  e2binPtrs = e2bins[j].getBinPtrs();
1084  e4binPtrs = e4bins[j].getBinPtrs();
1085  e6binPtrs = e6bins[j].getBinPtrs();
1086  e8binPtrs = e8bins[j].getBinPtrs();
1087  yErr.push_back(sampleError(cn));
1088  }
1089  // Put the bin ptrs back in place.
1090  e2binPtrs = e2->getBinPtrs();
1091  e4binPtrs = e4->getBinPtrs();
1092  e6binPtrs = e6->getBinPtrs();
1093  e8binPtrs = e8->getBinPtrs();
1094  fillScatter(h, binx, cn, yErr);
1095  }
1096 
1097 
1099  const void vnEightInt(Scatter2DPtr h, ECorrPtr e2, ECorrPtr e4, ECorrPtr e6, ECorrPtr e8) const {
1100  cnEightInt(h, e2, e4, e6, e8);
1101  nthPow(h, 1./8., -1./33.);
1102  }
1103 
1104 
1106  const void vnTwoDiff(Scatter2DPtr h, ECorrPtr e2Dif) const {
1107  auto e2bins = e2Dif->getBins();
1108  auto ref = e2Dif->getReference();
1109  auto binx = e2Dif->getBinX();
1110  if (binx.size() -1 != e2bins.size()) {
1111  cout << "vnTwoDif: Bin size (x,y) differs!" << endl;
1112  return;
1113  }
1114  vector<CorBinBase*> e2binPtrs;
1115  vector<CorBinBase*> refPtrs;
1116  auto vn = [&] (int i) {
1117  // Test reference flow.
1118  if (ref.mean() <= 0) return 0.;
1119  return e2binPtrs[i]->mean() / sqrt(ref.mean());
1120  };
1121  // We need here a separate error function, as we don't iterate over the reference flow.
1122  auto vnerr = [&] (int i) {
1123  // Test reference flow.
1124  if (refPtrs[i]->mean() <=0) return 0.;
1125  return e2binPtrs[i]->mean() / sqrt(refPtrs[i]->mean());
1126  };
1127  // Error calculation.
1128  vector<pair<double, double> > yErr;
1129  refPtrs = ref.getBinPtrs();
1130  for (int j = 0, N = e2bins.size(); j < N; ++j) {
1131  e2binPtrs = e2bins[j].getBinPtrs();
1132  yErr.push_back(sampleError(vnerr));
1133  }
1134  // Put the e2binPtrs back in place.
1135  e2binPtrs = e2Dif->getBinPtrs();
1136  fillScatter(h, binx, vn);
1137  }
1138 
1139 
1141  const void vnFourDiff(Scatter2DPtr h, ECorrPtr e2Dif, ECorrPtr e4Dif) const {
1142  auto e2bins = e2Dif->getBins();
1143  auto e4bins = e4Dif->getBins();
1144  auto ref2 = e2Dif->getReference();
1145  auto ref4 = e4Dif->getReference();
1146  auto binx = e2Dif->getBinX();
1147  if (binx.size() - 1 != e2bins.size()){
1148  cout << "vnFourDif: Bin size (x,y) differs!" << endl;
1149  return;
1150  }
1151  if (binx != e4Dif->getBinX()) {
1152  cout << "Error in vnFourDif: Correlator x-binning differs!" << endl;
1153  return;
1154  }
1155  vector<CorBinBase*> e2binPtrs;
1156  vector<CorBinBase*> e4binPtrs;
1157  vector<CorBinBase*> ref2Ptrs;
1158  vector<CorBinBase*> ref4Ptrs;
1159  double denom = 2 * ref2.mean() * ref2.mean() - ref4.mean();
1160  auto vn = [&] (int i) {
1161  // Test denominator.
1162  if (denom <= 0 ) return 0.;
1163  return ((2 * ref2.mean() * e2bins[i].mean() - e4bins[i].mean()) / pow(denom, 0.75));
1164  };
1165  // We need here a separate error function, as we don't iterate over the reference flow.
1166  auto vnerr = [&] (int i) {
1167  double denom2 = 2 * ref2Ptrs[i]->mean() * ref2Ptrs[i]->mean() -
1168  ref4Ptrs[i]->mean();
1169  // Test denominator.
1170  if (denom2 <= 0) return 0.;
1171  return ((2 * ref2Ptrs[i]->mean() * e2binPtrs[i]->mean() - e4binPtrs[i]->mean()) / pow(denom2, 0.75));
1172  };
1173  // Error calculation.
1174  vector<pair<double, double> > yErr;
1175  ref2Ptrs = ref2.getBinPtrs();
1176  ref4Ptrs = ref4.getBinPtrs();
1177  for (int j = 0, N = e2bins.size(); j < N; ++j) {
1178  e2binPtrs = e2bins[j].getBinPtrs();
1179  e4binPtrs = e4bins[j].getBinPtrs();
1180  yErr.push_back(sampleError(vnerr));
1181  }
1182  // Put the binPtrs back in place.
1183  e2binPtrs = e2Dif->getBinPtrs();
1184  e4binPtrs = e4Dif->getBinPtrs();
1185  fillScatter(h, binx, vn, yErr);
1186  }
1187 
1188  };
1189 
1190 
1191 }
1192 
1193 #endif
Definition: MC_Cent_pPb.hh:10
static pair< double, double > sampleEnvelope(T func)
Calculate the bootstrapped sample envelope.
Definition: Correlators.hh:853
const vector< int > getH1() const
Get a copy of the h1 harmonic vector.
Definition: Correlators.hh:499
std::enable_if< std::is_arithmetic< NUM1 >::value &&std::is_arithmetic< NUM2 >::value, int >::type binIndex(NUM1 val, std::initializer_list< NUM2 > binedges, bool allow_overflow=false)
Return the bin index of the given value, val, given a vector of bin edges.
Definition: MathUtils.hh:365
void fillFromProfs()
Fill bins with content from preloaded histograms.
Definition: Correlators.hh:528
ECorrPtr bookECorrelator(const string name, const Scatter2DPtr hIn)
Templated version of correlator booking which takes N desired harmonic and M number of particles...
Definition: Correlators.hh:663
const vector< double > getBinX() const
Get a copy of the bin x-values.
Definition: Correlators.hh:494
const void vnTwoInt(Scatter2DPtr h, ECorrPtr e2) const
Two particle integrated vn.
Definition: Correlators.hh:905
const void cnSixInt(Scatter2DPtr h, ECorrPtr e2, ECorrPtr e4, ECorrPtr e6) const
Six particle integrated cn.
Definition: Correlators.hh:1005
CmpState compare(const Projection &p) const
Definition: Correlators.hh:145
static vector< int > hVec(int n, int m)
Construct a harmonic vectors from n harmonics and m number of particles.
Definition: Correlators.hh:106
void fill(const Correlators &c1, const Correlators &c2, const double &weight=1.0)
Fill bins with the appropriate correlator, and a rapidity gap between two Correlators.
Definition: Correlators.hh:463
const vector< CorBinBase * > getBinPtrs()
Return the bins as pointers to the base class.
Definition: Correlators.hh:487
const vector< int > getH2() const
Get a copy of the h2 harmonic vector.
Definition: Correlators.hh:504
const void vnFourDiff(Scatter2DPtr h, ECorrPtr e2Dif, ECorrPtr e4Dif) const
Four particle differential vn.
Definition: Correlators.hh:1141
static pair< double, double > sampleVariance(T func)
Calculate the bootstrapped sample variance.
Definition: Correlators.hh:835
ECorrPtr bookECorrelator(const string name, const vector< int > &h, const Scatter2DPtr hIn)
Book an ECorrelator in the same way as a histogram.
Definition: Correlators.hh:609
shared_ptr< ECorrelator > ECorrPtr
Typedef of shared pointer to ECorrelator.
Definition: Correlators.hh:604
const void cnEightInt(Scatter2DPtr h, ECorrPtr e2, ECorrPtr e4, ECorrPtr e6, ECorrPtr e8) const
Eight particle integrated cn.
Definition: Correlators.hh:1052
ECorrelator(vector< int > h1In, vector< int > h2In, vector< double > binIn)
Constructor for gapped correlator.
Definition: Correlators.hh:418
Particle representation, either from a HepMC::GenEvent or reconstructed.
Definition: Particle.hh:18
A helper class to calculate all event averages of correlators.
Definition: Correlators.hh:393
void fill(const double &obs, const Correlators &c1, const Correlators &c2, double weight=1.0)
Fill the appropriate bin given an input (per event) observable, e.g. centrality, with a rapidity gap ...
Definition: Correlators.hh:432
static constexpr double DBL_NAN
Convenient const for getting the double NaN value.
Definition: Utils.hh:46
Base class for projections which return subsets of an event&#39;s particles.
Definition: ParticleFinder.hh:11
Tools for flow analyses.
Definition: Correlators.hh:221
const vector< pair< double, double > > pTBinnedCorrelatorsGap(const Correlators &other, vector< int > n1, vector< int > n2, bool overflow=false) const
pT differential correlators of n1 harmonic, for number n1.size()
const void vnEightInt(Scatter2DPtr h, ECorrPtr e2, ECorrPtr e4, ECorrPtr e6, ECorrPtr e8) const
Eight particle integrated vn.
Definition: Correlators.hh:1099
void fill(const Correlators &c, const double &weight=1.0)
Fill the bins with the appropriate correlator.
Definition: Correlators.hh:446
void fill(const double &obs, const Correlators &c, double weight=1.0)
Fill the appropriate bin given an input (per event) observable, e.g. centrality.
Definition: Correlators.hh:423
void setReference(CorBin refIn)
Replace reference flow bin with another, e.g. calculated in another phase space or with other pid...
Definition: Correlators.hh:509
void stream()
Stream correlators to the yoda file.
Definition: Correlators.hh:682
This is the base class of all analysis classes in Rivet.
Definition: Analysis.hh:64
Representation of a HepMC event, and enabler of Projection caching.
Definition: Event.hh:22
Cmp< Projection > mkPCmp(const Projection &otherparent, const std::string &pname) const
const void fillScatter(Scatter2DPtr h, vector< double > &binx, F func, vector< pair< double, double > > &yErr) const
Helper method for turning correlators into Scatter2Ds with error estimates.
Definition: Correlators.hh:738
static pair< int, int > getMaxValues(vector< vector< int > > &hList)
Return the maximal values for n, p to be used in the constructor of Correlators(xxx, nMax, pMax, xxxx)
Definition: Correlators.hh:121
list< Profile1DPtr >::iterator profEnd()
End iterator for the list of associated profile histograms.
Definition: Correlators.hh:551
static void nthPow(Scatter2DPtr h, const double &n, const double &k=1.0)
Take the n th power of all points in h, and put the result back in the same Scatter2D.
Definition: Correlators.hh:800
static void fillScatter(Scatter2DPtr h, vector< double > &binx, T func)
Helper method for turning correlators into Scatter2Ds.
Definition: Correlators.hh:714
const pair< double, double > sampleError(T func) const
Selection method for which sample error to use, given in the constructor.
Definition: Correlators.hh:872
const vector< pair< double, double > > pTBinnedCorrelators(vector< int > n, bool overflow=false) const
pT differential correlator of n harmonic, for number of powers n.size()
const vector< CorBin > getBins() const
Get a copy of the bin contents.
Definition: Correlators.hh:482
ECorrPtr bookECorrelatorGap(const string name, const vector< int > &h, const Scatter2DPtr hIn)
Definition: Correlators.hh:651
virtual std::string name() const
Get the name of the projection.
Definition: Projection.hh:56
const pair< double, double > intCorrelatorGap(const Correlators &other, vector< int > n1, vector< int > n2) const
Integrated correlator of n1 harmonic, for number of powers n1.size()
Projection for calculating correlators for flow measurements.
Definition: Correlators.hh:38
const pair< int, int > getMaxValues() const
Get the correct max N and max P for the set of booked correlators.
Definition: Correlators.hh:587
const void corrPlot(list< Profile1DPtr > profs, ECorrPtr e) const
Put an event-averaged correlator into Profile1Ds, one for each bootstrapping bin. ...
Definition: Correlators.hh:920
CumulantAnalysis(const string &n)
Constructor.
Definition: Correlators.hh:702
void project(const Event &e)
list< Profile1DPtr >::iterator profBegin()
Begin iterator for the list of associated profile histograms.
Definition: Correlators.hh:546
const void vnSixInt(Scatter2DPtr h, ECorrPtr e2, ECorrPtr e4, ECorrPtr e6) const
Six particle integrated vn.
Definition: Correlators.hh:1044
const pair< double, double > intCorrelator(vector< int > n) const
Integrated correlator of n harmonic, with the number of powers being the size of n. E.G. n should be: <<2>>_2 => n = {2, -2} <<4>>_2 => n = {2, 2, -2, -2} <<2>>_4 => n = {4, -4} <<4>>_4 => n = {4, 4, -4, 4}, and so on.
static void nthPow(Scatter2DPtr hOut, const Scatter2DPtr hIn, const double &n, const double &k=1.0)
Take the n th power of all points in hIn and put the result in hOut.
Definition: Correlators.hh:762
void setProfs(list< Profile1DPtr > prIn)
Set the prIn list of profile histograms associated with the internal bins.
Definition: Correlators.hh:523
const void vnFourInt(Scatter2DPtr h, ECorrPtr e2, ECorrPtr e4) const
Four particle integrated vn.
Definition: Correlators.hh:998
const void vnTwoDiff(Scatter2DPtr h, ECorrPtr e2Dif) const
Two particle differential vn.
Definition: Correlators.hh:1106
Correlators(const ParticleFinder &fsp, int nMaxIn=2, int pMaxIn=0, vector< double > pTbinEdgesIn={})
ECorrelator(vector< int > h, vector< double > binIn)
Constructor. Takes as argument the desired harmonic and number of correlated particles as a generic f...
Definition: Correlators.hh:410
ECorrPtr bookECorrelator(const string name, const vector< int > &h1, const vector< int > &h2, const Scatter2DPtr hIn)
Book a gapped ECorrelator with two harmonic vectors.
Definition: Correlators.hh:629
ECorrPtr bookECorrelatorGap(const string name, const Scatter2DPtr hIn)
Templated version of gapped correlator booking which takes N desired harmonic and M number of particl...
Definition: Correlators.hh:672
Base class for all Rivet projections.
Definition: Projection.hh:29
const void cnTwoInt(Scatter2DPtr h, ECorrPtr e2) const
Two-particle integrated cn.
Definition: Correlators.hh:882
std::enable_if< std::is_arithmetic< NUM >::value, double >::type mean(const vector< NUM > &sample)
Definition: MathUtils.hh:407
const CorBin getReference() const
Extract the reference flow from a differential event averaged correlator.
Definition: Correlators.hh:514
const void corrPlot(Scatter2DPtr h, ECorrPtr e) const
Put an event-averaged correlator into a Scatter2D.
Definition: Correlators.hh:914