Rivet  3.1.5
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 YODA::Scatter2D 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 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  double mean() const {
283  if (_sumW < 1e-10) return 0;
284  return _sumWX / _sumW;
285  }
286 
288  double sumW() const {
289  return _sumW;
290  }
291 
293  double sumW2() const {
294  return _sumW2;
295  }
296 
298  double sumWX() const {
299  return _sumWX;
300  }
301 
303  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  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  vector<CorSingleBin> getBins() const {
366  return bins;
367  }
368 
370  template<class T=CorBinBase>
371  vector<T*> getBinPtrs() {
372  vector<T*> ret(bins.size());
373  transform(bins.begin(), bins.end(), ret.begin(), [](CorSingleBin& b) {return &b;});
374  return ret;
375  }
376 
377  private:
378 
379  vector<CorSingleBin> bins;
380  size_t binIndex;
381  size_t nBins;
382 
383  };
384 
385 
386  public:
387 
392  class ECorrelator {
393  public:
394 
400  //ECorrelator(vector<int> h) : h1(h), h2({}),
401  // binX(0), binContent(0), reference() {
402  //};
403 
409  ECorrelator(vector<int> h, vector<double> binIn)
410  : h1(h), h2({}), binX(binIn), binContent(binIn.size() - 1), reference()
411  { }
412 
417  ECorrelator(vector<int> h1In, vector<int> h2In, vector<double> binIn)
418  : h1(h1In), h2(h2In), binX(binIn), binContent(binIn.size() - 1), reference()
419  { }
420 
422  void fill(const double& obs, const Correlators& c, double weight=1.0) {
423  int index = getBinIndex(obs);
424  if (index < 0) return;
425  binContent[index].fill(c.intCorrelator(h1), weight);
426  }
427 
431  void fill(const double& obs, const Correlators& c1, const Correlators& c2, double weight=1.0) {
432  if (!h2.size()) {
433  cout << "Trying to fill gapped correlator, but harmonics behind the gap (h2) are not given!" << endl;
434  return;
435  }
436  int index = getBinIndex(obs);
437  if (index < 0) return;
438  binContent[index].fill(c1.intCorrelatorGap(c2, h1, h2), weight);
439  }
440 
445  void fill(const Correlators& c, const double& weight=1.0) {
446  vector< pair<double, double> > diffCorr = c.pTBinnedCorrelators(h1);
447  // We always skip overflow when calculating the all-event average.
448  if (diffCorr.size() != binX.size() - 1)
449  cout << "Tried to fill event with wrong binning (ungapped)" << endl;
450  for (size_t i = 0; i < diffCorr.size(); ++i) {
451  int index = getBinIndex(binX[i]);
452  if (index < 0) return;
453  binContent[index].fill(diffCorr[i], weight);
454  }
455  reference.fill(c.intCorrelator(h1), weight);
456  }
457 
462  void fill(const Correlators& c1, const Correlators& c2, const double& weight = 1.0) {
463  if (!h2.size()) {
464  cout << "Trying to fill gapped correlator, but harmonics behind "
465  "the gap (h2) are not given!" << endl;
466  return;
467  }
468  vector< pair<double, double> > diffCorr = c1.pTBinnedCorrelatorsGap(c2, h1, h2);
469  // We always skip overflow when calculating the all event average.
470  if (diffCorr.size() != binX.size() - 1)
471  cout << "Tried to fill event with wrong binning (gapped)" << endl;
472  for (size_t i = 0; i < diffCorr.size(); ++i) {
473  int index = getBinIndex(binX[i]);
474  if (index < 0) return;
475  binContent[index].fill(diffCorr[i], weight);
476  }
477  reference.fill(c1.intCorrelatorGap(c2, h1, h2), weight);
478  }
479 
481  vector<CorBin> getBins() const {
482  return binContent;
483  }
484 
486  vector<CorBinBase*> getBinPtrs() {
487  vector<CorBinBase*> ret(binContent.size());
488  transform(binContent.begin(), binContent.end(), ret.begin(), [](CorBin& b) {return &b;});
489  return ret;
490  }
491 
493  vector<double> getBinX() const {
494  return binX;
495  }
496 
498  vector<int> getH1() const {
499  return h1;
500  }
501 
503  vector<int> getH2() const {
504  return h2;
505  }
506 
508  void setReference(CorBin refIn) {
509  reference = refIn;
510  }
511 
513  CorBin getReference() const {
514  if (reference.mean() < 1e-10)
515  cout << "Warning: ECorrelator, reference bin is zero." << endl;
516  return reference;
517  }
518 
521  void setProfs(vector<string> prIn) {
522  profs = prIn;
523  }
524 
526  bool fillFromProfile(YODA::AnalysisObjectPtr yao, string name) {
527  auto refs = reference.getBinPtrs<CorSingleBin>();
528  for (size_t i = 0; i < profs.size(); ++i) {
529  if (yao->path() == "/RAW/"+name+"/TMP/"+profs[i]) {
530  YODA::Profile1DPtr pPtr = dynamic_pointer_cast<YODA::Profile1D>(yao);
531  for (size_t j = 0; j < binX.size() - 1; ++j) {
532  const YODA::ProfileBin1D& pBin = pPtr->binAt(binX[j]);
533  auto tmp = binContent[j].getBinPtrs<CorSingleBin>();
534  tmp[i]->addContent(pBin.numEntries(), pBin.sumW(), pBin.sumW2(),
535  pBin.sumWY());
536  }
537  // Get the reference flow from the underflow bin of the histogram.
538  const YODA::Dbn2D& uBin = pPtr->underflow();
539  refs[i]->addContent(uBin.numEntries(), uBin.sumW(), uBin.sumW2(),
540  uBin.sumWY());
541  return true;
542  }
543  } // End loop of bootstrapped correlators.
544  return false;
545  }
546 
547  private:
548 
549  // Get correct bin index for a given @param obs value
550  int getBinIndex(const double& obs) const {
551  // Find the correct index of binContent.
552  // If we are in overflow, just skip.
553  if (obs >= binX.back()) return -1;
554  // If we are in underflow, ditto.
555  if (obs < binX[0]) return -1;
556  int index = 0;
557  for (int i = 0, N = binX.size() - 1; i < N; ++i, ++index)
558  if (obs >= binX[i] && obs < binX[i + 1]) break;
559  return index;
560  }
561 
562  // The harmonics vectors.
563  vector<int> h1;
564  vector<int> h2;
565 
566  // The bins.
567  vector<double> binX;
568  vector<CorBin> binContent;
569  // The reference flow.
570  CorBin reference;
571  public:
572  // The profile histograms associated with the CorBins for streaming.
573  vector <string> profs;
574 
575  };
576 
577 
579  const pair<int, int> getMaxValues() const {
580  vector< vector<int>> harmVecs;
581  for ( auto eItr = eCorrPtrs.begin(); eItr != eCorrPtrs.end(); ++eItr) {
582  vector<int> h1 = (*eItr)->getH1();
583  vector<int> h2 = (*eItr)->getH2();
584  if (h1.size() > 0) harmVecs.push_back(h1);
585  if (h2.size() > 0) harmVecs.push_back(h2);
586  }
587  if (harmVecs.size() == 0) {
588  cout << "Warning: You tried to extract max values from harmonic "
589  "vectors, but have not booked any." << endl;
590  return pair<int,int>();
591  }
592  return Correlators::getMaxValues(harmVecs);
593  }
594 
596  typedef shared_ptr<ECorrelator> ECorrPtr;
597 
600  ECorrPtr bookECorrelator(const string name, const vector<int>& h, const YODA::Scatter2D& hIn) {
601  vector<double> binIn;
602  for (auto b : hIn.points()) binIn.push_back(b.xMin());
603  binIn.push_back(hIn.points().back().xMax());
604  return bookECorrelator(name, h, binIn);
605  }
606 
609  ECorrPtr bookECorrelator(const string name, const vector<int>& h, vector<double>& binIn) {
610  ECorrPtr ecPtr = ECorrPtr(new ECorrelator(h, binIn));
611  vector<string> eCorrProfs;
612  for (int i = 0; i < BOOT_BINS; ++i) {
613  Profile1DPtr tmp;
614  book(tmp,"TMP/"+name+"-"+to_string(i),binIn);
615  eCorrProfs.push_back(name+"-"+to_string(i));
616  }
617  ecPtr->setProfs(eCorrProfs);
618  eCorrPtrs.push_back(ecPtr);
619  return ecPtr;
620  }
621 
624  ECorrPtr bookECorrelator(const string name, const vector<int>& h1,
625  const vector<int>& h2, vector<double>& binIn) {
626  ECorrPtr ecPtr = ECorrPtr(new ECorrelator(h1, h2, binIn));
627  vector<string> eCorrProfs;
628  for (int i = 0; i < BOOT_BINS; ++i) {
629  Profile1DPtr tmp;
630  book(tmp,"TMP/"+name+"-"+to_string(i),binIn);
631  eCorrProfs.push_back(name+"-"+to_string(i));
632  }
633  ecPtr->setProfs(eCorrProfs);
634  eCorrPtrs.push_back(ecPtr);
635  return ecPtr;
636  }
637 
640  ECorrPtr bookECorrelator(const string name, const vector<int>& h1,
641  const vector<int>& h2, const YODA::Scatter2D& hIn ) {
642  vector<double> binIn;
643  for (auto b : hIn.points()) binIn.push_back(b.xMin());
644  binIn.push_back(hIn.points().back().xMax());
645  return bookECorrelator(name, h1, h2, binIn);
646  }
647 
652  ECorrPtr bookECorrelatorGap(const string name, const vector<int>& h,
653  const YODA::Scatter2D& hIn) {
654  const vector<int> h1(h.begin(), h.begin() + h.size() / 2);
655  const vector<int> h2(h.begin() + h.size() / 2, h.end());
656  return bookECorrelator(name, h1, h2, hIn);
657  }
658 
663  template<unsigned int N, unsigned int M>
664  ECorrPtr bookECorrelator(const string name, vector<double> binIn) {
665  return bookECorrelator(name, Correlators::hVec(N, M), binIn);
666  }
667 
672  template<unsigned int N, unsigned int M>
673  ECorrPtr bookECorrelator(const string name, const YODA::Scatter2D& hIn) {
674  return bookECorrelator(name, Correlators::hVec(N, M), hIn);
675  }
676 
681  template<unsigned int N, unsigned int M>
682  ECorrPtr bookECorrelatorGap(const string name, const YODA::Scatter2D& hIn) {
683  const vector<int> h = Correlators::hVec(N,M);
684  const vector<int> h1(h.begin(), h.begin() + h.size() / 2);
685  const vector<int> h2(h.begin() + h.size() / 2, h.end());
686  return bookECorrelator(name, h1, h2, hIn);
687 
688  }
689 
690  protected:
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  // Test if we have proper bins from a booked histogram.
717  bool hasBins = (h->points().size() > 0);
718  for (int i = 0, N = binx.size() - 1; i < N; ++i) {
719  double xMid = (binx[i] + binx[i + 1]) / 2.0;
720  double xeMin = fabs(xMid - binx[i]);
721  double xeMax = fabs(xMid - binx[i + 1]);
722  if (hasBins) {
723  xMid = h->points()[i].x();
724  xeMin = h->points()[i].xErrMinus();
725  xeMax = h->points()[i].xErrPlus();
726  }
727  double yVal = func(i);
728  if (std::isnan(yVal)) yVal = 0.;
729  double yErr = 0;
730  points.push_back(YODA::Point2D(xMid, yVal, xeMin, xeMax, yErr, yErr));
731  }
732  h->reset();
733  h->points().clear();
734  for (int i = 0, N = points.size(); i < N; ++i) h->addPoint(points[i]);
735  }
736 
744  template<typename F>
745  void fillScatter(Scatter2DPtr h, vector<double>& binx, F func,
746  vector<pair<double, double> >& yErr) const {
747  vector<YODA::Point2D> points;
748  // Test if we have proper bins from a booked histogram.
749  bool hasBins = (h->points().size() > 0);
750  for (int i = 0, N = binx.size() - 1; i < N; ++i) {
751  double xMid = (binx[i] + binx[i + 1]) / 2.0;
752  double xeMin = fabs(xMid - binx[i]);
753  double xeMax = fabs(xMid - binx[i + 1]);
754  if (hasBins) {
755  xMid = h->points()[i].x();
756  xeMin = h->points()[i].xErrMinus();
757  xeMax = h->points()[i].xErrPlus();
758  }
759  double yVal = func(i);
760  if (std::isnan(yVal))
761  points.push_back(YODA::Point2D(xMid, 0., xeMin, xeMax,0., 0.));
762  else
763  points.push_back(YODA::Point2D(xMid, yVal, xeMin, xeMax,
764  yErr[i].first, yErr[i].second));
765  }
766  h->reset();
767  h->points().clear();
768 
769  for (int i = 0, N = points.size(); i < N; ++i)
770  h->addPoint(points[i]);
771  }
772 
773 
777  static void nthPow(Scatter2DPtr hOut, const Scatter2DPtr hIn, const double& n, const double& k = 1.0) {
778  if (n == 0 || n == 1) {
779  cout << "Error: Do not take the 0th or 1st power of a Scatter2D,"
780  " use scale instead." << endl;
781  return;
782  }
783  if (hIn->points().size() != hOut->points().size()) {
784  cout << "nthRoot: Scatterplots: " << hIn->name() << " and " <<
785  hOut->name() << " not initialized with same length." << endl;
786  return;
787  }
788  vector<YODA::Point2D> points;
789  // The error pre-factor is k^(1/n) / n by Taylors formula.
790  double eFac = pow(k,1./n)/n;
791  for (auto b : hIn->points()) {
792  double yVal = pow(k * b.y(),n);
793  if (std::isnan(yVal))
794  points.push_back(YODA::Point2D(b.x(), 0., b.xErrMinus(),
795  b.xErrPlus(), 0, 0 ));
796  else {
797  double yemin = abs(eFac * pow(yVal,1./(n - 1.))) * b.yErrMinus();
798  if (std::isnan(yemin)) yemin = b.yErrMinus();
799  double yemax = abs(eFac * pow(yVal,1./(n - 1.))) * b.yErrPlus();
800  if (std::isnan(yemax)) yemax = b.yErrPlus();
801  points.push_back(YODA::Point2D(b.x(), yVal, b.xErrMinus(),
802  b.xErrPlus(), yemin, yemax ));
803  }
804  }
805  hOut->reset();
806  hOut->points().clear();
807  for (int i = 0, N = points.size(); i < N; ++i)
808  hOut->addPoint(points[i]);
809  }
810 
811 
815  static void nthPow(Scatter2DPtr h, const double& n, const double& k = 1.0) {
816  if (n == 0 || n == 1) {
817  cout << "Error: Do not take the 0th or 1st power of a Scatter2D,"
818  " use scale instead." << endl;
819  return;
820  }
821  vector<YODA::Point2D> points;
822  vector<YODA::Point2D> pIn = h->points();
823  // The error pre-factor is k^(1/n) / n by Taylors formula.
824  double eFac = pow(k,1./n)/n;
825  for (auto b : pIn) {
826  double yVal = pow(k * b.y(),n);
827  if (std::isnan(yVal))
828  points.push_back(YODA::Point2D(b.x(), 0., b.xErrMinus(),
829  b.xErrPlus(), 0, 0 ));
830  else {
831  double yemin = abs(eFac * pow(yVal,1./(n - 1.))) * b.yErrMinus();
832  if (std::isnan(yemin)) yemin = b.yErrMinus();
833  double yemax = abs(eFac * pow(yVal,1./(n - 1.))) * b.yErrPlus();
834  if (std::isnan(yemax)) yemax = b.yErrPlus();
835  points.push_back(YODA::Point2D(b.x(), yVal, b.xErrMinus(),
836  b.xErrPlus(), yemin, yemax ));
837  }
838  }
839  h->reset();
840  h->points().clear();
841  for (int i = 0, N = points.size(); i < N; ++i) h->addPoint(points[i]);
842  }
843 
844 
849  template<typename T>
850  static pair<double, double> sampleVariance(T func) {
851  // First we calculate the mean (two pass calculation).
852  double avg = 0.;
853  for (int i = 0; i < BOOT_BINS; ++i) avg += func(i);
854  avg /= double(BOOT_BINS);
855  // Then we find the variance.
856  double var = 0.;
857  for (int i = 0; i < BOOT_BINS; ++i) var += pow(func(i) - avg, 2.);
858  var /= (double(BOOT_BINS) - 1);
859  return pair<double, double>(var, var);
860  }
861 
862 
867  template<typename T>
868  static pair<double, double> sampleEnvelope(T func) {
869  // First we calculate the mean.
870  double avg = 0.;
871  for (int i = 0; i < BOOT_BINS; ++i) avg += func(i);
872  avg /= double(BOOT_BINS);
873  double yMax = avg;
874  double yMin = avg;
875  // Then we find the envelope using the mean as initial value.
876  for (int i = 0; i < BOOT_BINS; ++i) {
877  double yVal = func(i);
878  if (yMin > yVal) yMin = yVal;
879  else if (yMax < yVal) yMax = yVal;
880  }
881  return pair<double, double>(fabs(avg - yMin), fabs(yMax - avg));
882  }
883 
884 
886  template<typename T>
887  const pair<double, double> sampleError(T func) const {
888  if (errorMethod == VARIANCE) return sampleVariance(func);
889  else if (errorMethod == ENVELOPE) return sampleEnvelope(func);
890  else
891  cout << "Error: Error method not found!" << endl;
892  return pair<double, double>(0.,0.);
893  }
894 
895 
897  void cnTwoInt(Scatter2DPtr h, ECorrPtr e2) const {
898  vector<CorBin> bins = e2->getBins();
899  vector<double> binx = e2->getBinX();
900  // Assert bin size.
901  if (binx.size() - 1 != bins.size()){
902  cout << "cnTwoInt: Bin size (x,y) differs!" << endl;
903  return;
904  }
905  vector<CorBinBase*> binPtrs;
906  // The mean value of the cumulant.
907  auto cn = [&] (int i) { return binPtrs[i]->mean(); };
908  // Error calculation.
909  vector<pair<double, double> > yErr;
910  for (int j = 0, N = bins.size(); j < N; ++j) {
911  binPtrs = bins[j].getBinPtrs();
912  yErr.push_back(sampleError(cn));
913  }
914  binPtrs = e2->getBinPtrs();
915  fillScatter(h, binx, cn, yErr);
916  }
917 
918 
920  void vnTwoInt(Scatter2DPtr h, ECorrPtr e2) const {
921  cnTwoInt(h, e2);
922  nthPow(h, 0.5);
923  }
924 
925 
929  void corrPlot(Scatter2DPtr h, ECorrPtr e) const {
930  cnTwoInt(h, e);
931  }
932 
933 
934 
935 
936  // TODO Use full path for lookup, change to single AU in output, rename.
937  void rawHookIn(YODA::AnalysisObjectPtr yao) final {
938  // Fill the corresponding ECorrelator.
939  for (auto ec : eCorrPtrs) if(ec->fillFromProfile(yao, name())) break;;
940  }
941 
946  void rawHookOut(vector<MultiweightAOPtr> raos, size_t iW) final {
947  // Loop over the correlators and extract the numbers.
948  for (auto ec : eCorrPtrs) {
949  vector<CorBin> corBins = ec->getBins();
950  vector<double> binx = ec->getBinX();
951  auto ref = ec->getReference();
952  auto refBins = ref.getBinPtrs<CorSingleBin>();
953  // Assert bin size.
954  if (binx.size() - 1 != corBins.size()){
955  cout << "corrPlot: Bin size (x,y) differs!" << endl;
956  return;
957  }
958  // Loop over the booked histograms using their names.
959  for (int i = 0, N = ec->profs.size(); i < N; ++i) {
960  for (auto rao : raos) {
961  if (rao->path() == "/"+name()+"/TMP/"+ec->profs[i]) {
962  // Get a pointer to the active profile.
963  rao.get()->setActiveWeightIdx(iW);
964  YODA::Profile1DPtr pPtr = dynamic_pointer_cast<YODA::Profile1D>(
965  rao.get()->activeYODAPtr());
966  // New bins.
967  vector<YODA::ProfileBin1D> profBins;
968  // Numbers for the summary distribution
969  double ne = 0., sow = 0., sow2 = 0.;
970  for (size_t j = 0, N = binx.size() - 1; j < N; ++j) {
971  vector<CorSingleBin*> binPtrs =
972  corBins[j].getBinPtrs<CorSingleBin>();
973  // Construct bin of the profiled quantities. We have no information
974  // (and no desire to add it) of sumWX of the profile, so really
975  // we should use a Dbn1D - but that does not work for Profile1D's.
976  profBins.push_back( YODA::ProfileBin1D(pPtr->bin(j).xEdges(),
977  YODA::Dbn2D( binPtrs[i]->numEntries(), binPtrs[i]->sumW(),
978  binPtrs[i]->sumW2(), 0., 0., binPtrs[i]->sumWX(), 0, 0)));
979  ne += binPtrs[i]->numEntries();
980  sow += binPtrs[i]->sumW();
981  sow2 += binPtrs[i]->sumW2();
982  }
983  // Put the ECorrelator into the raw histogram.
984  pPtr->reset();
985  pPtr->bins().clear();
986  // Add the bins.
987  pPtr->addBins(profBins);
988  // Set the total distribution.
989  pPtr->setTotalDbn(YODA::Dbn2D(ne,sow,sow2,0.,0.,0.,0.,0.));
990  // And reference flow in the underflow bin.
991  pPtr->setUnderflow(YODA::Dbn2D(refBins[i]->numEntries(),
992  refBins[i]->sumW(), refBins[i]->sumW2(), 0., 0.,
993  refBins[i]->sumWX(), 0., 0.));
994  }
995  }
996  }
997  }
998  }
999 
1000  // @brief Four particle integrated cn.
1001  void cnFourInt(Scatter2DPtr h, ECorrPtr e2, ECorrPtr e4) const {
1002  auto e2bins = e2->getBins();
1003  auto e4bins = e4->getBins();
1004  auto binx = e2->getBinX();
1005  if (binx.size() - 1 != e2bins.size()){
1006  cout << "cnFourInt: Bin size (x,y) differs!" << endl;
1007  return;
1008  }
1009  if (binx != e4->getBinX()) {
1010  cout << "Error in cnFourInt: Correlator x-binning differs!" << endl;
1011  return;
1012  }
1013  vector<CorBinBase*> e2binPtrs;
1014  vector<CorBinBase*> e4binPtrs;
1015  auto cn = [&] (int i) {
1016  double e22 = e2binPtrs[i]->mean() * e2binPtrs[i]->mean();
1017  return e4binPtrs[i]->mean() - 2. * e22;
1018  };
1019  // Error calculation.
1020  vector<pair<double, double> > yErr;
1021  for (int j = 0, N = e2bins.size(); j < N; ++j) {
1022  e2binPtrs = e2bins[j].getBinPtrs();
1023  e4binPtrs = e4bins[j].getBinPtrs();
1024  yErr.push_back(sampleError(cn));
1025  }
1026  // Put the bin ptrs back in place.
1027  e2binPtrs = e2->getBinPtrs();
1028  e4binPtrs = e4->getBinPtrs();
1029  fillScatter(h, binx, cn, yErr);
1030  }
1031 
1032 
1034  void vnFourInt(Scatter2DPtr h, ECorrPtr e2, ECorrPtr e4) const {
1035  cnFourInt(h, e2, e4);
1036  nthPow(h, 0.25, -1.0);
1037  }
1038 
1039 
1041  void cnSixInt(Scatter2DPtr h, ECorrPtr e2, ECorrPtr e4,
1042  ECorrPtr e6) const {
1043  auto e2bins = e2->getBins();
1044  auto e4bins = e4->getBins();
1045  auto e6bins = e6->getBins();
1046  auto binx = e2->getBinX();
1047  if (binx.size() - 1 != e2bins.size()){
1048  cout << "cnSixInt: Bin size (x,y) differs!" << endl;
1049  return;
1050  }
1051  if (binx != e4->getBinX() || binx != e6->getBinX()) {
1052  cout << "Error in cnSixInt: Correlator x-binning differs!" << endl;
1053  return;
1054  }
1055  vector<CorBinBase*> e2binPtrs;
1056  vector<CorBinBase*> e4binPtrs;
1057  vector<CorBinBase*> e6binPtrs;
1058  auto cn = [&] (int i) {
1059  double e23 = pow(e2binPtrs[i]->mean(), 3.0);
1060  return e6binPtrs[i]->mean() - 9.*e2binPtrs[i]->mean()*e4binPtrs[i]->mean() +
1061  12.*e23;
1062  };
1063  // Error calculation.
1064  vector<pair<double, double> > yErr;
1065  for (int j = 0, N = e2bins.size(); j < N; ++j) {
1066  e2binPtrs = e2bins[j].getBinPtrs();
1067  e4binPtrs = e4bins[j].getBinPtrs();
1068  e6binPtrs = e6bins[j].getBinPtrs();
1069  yErr.push_back(sampleError(cn));
1070  }
1071  // Put the bin ptrs back in place.
1072  e2binPtrs = e2->getBinPtrs();
1073  e4binPtrs = e4->getBinPtrs();
1074  e6binPtrs = e6->getBinPtrs();
1075  fillScatter(h, binx, cn, yErr);
1076  }
1077 
1078 
1080  void vnSixInt(Scatter2DPtr h, ECorrPtr e2, ECorrPtr e4,
1081  ECorrPtr e6) const {
1082  cnSixInt(h, e2, e4, e6);
1083  nthPow(h, 1./6., 0.25);
1084  }
1085 
1086 
1088  void cnEightInt(Scatter2DPtr h, ECorrPtr e2, ECorrPtr e4,
1089  ECorrPtr e6, ECorrPtr e8) const {
1090  auto e2bins = e2->getBins();
1091  auto e4bins = e4->getBins();
1092  auto e6bins = e6->getBins();
1093  auto e8bins = e8->getBins();
1094  auto binx = e2->getBinX();
1095  if (binx.size() - 1 != e2bins.size()){
1096  cout << "cnEightInt: Bin size (x,y) differs!" << endl;
1097  return;
1098  }
1099  if (binx != e4->getBinX() || binx != e6->getBinX() ||
1100  binx != e8->getBinX()) {
1101  cout << "Error in cnEightInt: Correlator x-binning differs!" << endl;
1102  return;
1103  }
1104  vector<CorBinBase*> e2binPtrs;
1105  vector<CorBinBase*> e4binPtrs;
1106  vector<CorBinBase*> e6binPtrs;
1107  vector<CorBinBase*> e8binPtrs;
1108  auto cn = [&] (int i ) {
1109  double e22 = e2binPtrs[i]->mean() * e2binPtrs[i]->mean();
1110  double e24 = e22 * e22;
1111  double e42 = e4binPtrs[i]->mean() * e4binPtrs[i]->mean();
1112  return e8binPtrs[i]->mean() - 16. * e6binPtrs[i]->mean() *
1113  e2binPtrs[i]->mean() - 18. * e42 + 144. * e4binPtrs[i]->mean()*e22 - 144. * e24;
1114  };
1115  // Error calculation.
1116  vector<pair<double, double> > yErr;
1117  for (int j = 0, N = e2bins.size(); j < N; ++j) {
1118  e2binPtrs = e2bins[j].getBinPtrs();
1119  e4binPtrs = e4bins[j].getBinPtrs();
1120  e6binPtrs = e6bins[j].getBinPtrs();
1121  e8binPtrs = e8bins[j].getBinPtrs();
1122  yErr.push_back(sampleError(cn));
1123  }
1124  // Put the bin ptrs back in place.
1125  e2binPtrs = e2->getBinPtrs();
1126  e4binPtrs = e4->getBinPtrs();
1127  e6binPtrs = e6->getBinPtrs();
1128  e8binPtrs = e8->getBinPtrs();
1129  fillScatter(h, binx, cn, yErr);
1130  }
1131 
1132 
1134  void vnEightInt(Scatter2DPtr h, ECorrPtr e2, ECorrPtr e4, ECorrPtr e6, ECorrPtr e8) const {
1135  cnEightInt(h, e2, e4, e6, e8);
1136  nthPow(h, 1./8., -1./33.);
1137  }
1138 
1139 
1141  void vnTwoDiff(Scatter2DPtr h, ECorrPtr e2Dif) const {
1142  auto e2bins = e2Dif->getBins();
1143  auto ref = e2Dif->getReference();
1144  auto binx = e2Dif->getBinX();
1145  if (binx.size() -1 != e2bins.size()) {
1146  cout << "vnTwoDif: Bin size (x,y) differs!" << endl;
1147  return;
1148  }
1149  vector<CorBinBase*> e2binPtrs;
1150  vector<CorBinBase*> refPtrs;
1151  auto vn = [&] (int i) {
1152  // Test reference flow.
1153  if (ref.mean() <= 0) return 0.;
1154  return e2binPtrs[i]->mean() / sqrt(ref.mean());
1155  };
1156  // We need here a separate error function, as we don't iterate over the reference flow.
1157  auto vnerr = [&] (int i) {
1158  // Test reference flow.
1159  if (refPtrs[i]->mean() <=0) return 0.;
1160  return e2binPtrs[i]->mean() / sqrt(refPtrs[i]->mean());
1161  };
1162  // Error calculation.
1163  vector<pair<double, double> > yErr;
1164  refPtrs = ref.getBinPtrs();
1165  for (int j = 0, N = e2bins.size(); j < N; ++j) {
1166  e2binPtrs = e2bins[j].getBinPtrs();
1167  yErr.push_back(sampleError(vnerr));
1168  }
1169  // Put the e2binPtrs back in place.
1170  e2binPtrs = e2Dif->getBinPtrs();
1171  fillScatter(h, binx, vn);
1172  }
1173 
1174 
1176  void vnFourDiff(Scatter2DPtr h, ECorrPtr e2Dif, ECorrPtr e4Dif) const {
1177  auto e2bins = e2Dif->getBins();
1178  auto e4bins = e4Dif->getBins();
1179  auto ref2 = e2Dif->getReference();
1180  auto ref4 = e4Dif->getReference();
1181  auto binx = e2Dif->getBinX();
1182  if (binx.size() - 1 != e2bins.size()){
1183  cout << "vnFourDif: Bin size (x,y) differs!" << endl;
1184  return;
1185  }
1186  if (binx != e4Dif->getBinX()) {
1187  cout << "Error in vnFourDif: Correlator x-binning differs!" << endl;
1188  return;
1189  }
1190  vector<CorBinBase*> e2binPtrs;
1191  vector<CorBinBase*> e4binPtrs;
1192  vector<CorBinBase*> ref2Ptrs;
1193  vector<CorBinBase*> ref4Ptrs;
1194  double denom = 2 * ref2.mean() * ref2.mean() - ref4.mean();
1195  auto vn = [&] (int i) {
1196  // Test denominator.
1197  if (denom <= 0 ) return 0.;
1198  return ((2 * ref2.mean() * e2bins[i].mean() - e4bins[i].mean()) / pow(denom, 0.75));
1199  };
1200  // We need here a separate error function, as we don't iterate over the reference flow.
1201  auto vnerr = [&] (int i) {
1202  double denom2 = 2 * ref2Ptrs[i]->mean() * ref2Ptrs[i]->mean() -
1203  ref4Ptrs[i]->mean();
1204  // Test denominator.
1205  if (denom2 <= 0) return 0.;
1206  return ((2 * ref2Ptrs[i]->mean() * e2binPtrs[i]->mean() - e4binPtrs[i]->mean()) / pow(denom2, 0.75));
1207  };
1208  // Error calculation.
1209  vector<pair<double, double> > yErr;
1210  ref2Ptrs = ref2.getBinPtrs();
1211  ref4Ptrs = ref4.getBinPtrs();
1212  for (int j = 0, N = e2bins.size(); j < N; ++j) {
1213  e2binPtrs = e2bins[j].getBinPtrs();
1214  e4binPtrs = e4bins[j].getBinPtrs();
1215  yErr.push_back(sampleError(vnerr));
1216  }
1217  // Put the binPtrs back in place.
1218  e2binPtrs = e2Dif->getBinPtrs();
1219  e4binPtrs = e4Dif->getBinPtrs();
1220  fillScatter(h, binx, vn, yErr);
1221  }
1222 
1223  };
1224 
1225 
1226 }
1227 
1228 #endif
ECorrPtr bookECorrelator(const string name, const YODA::Scatter2D &hIn)
Templated version of correlator booking which takes N desired harmonic and M number of particles...
Definition: Correlators.hh:673
Definition: MC_Cent_pPb.hh:10
static pair< double, double > sampleEnvelope(T func)
Calculate the bootstrapped sample envelope.
Definition: Correlators.hh:868
vector< CorBinBase * > getBinPtrs()
Return the bins as pointers to the base class.
Definition: Correlators.hh:486
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:416
ECorrPtr bookECorrelator(const string name, const vector< int > &h, const YODA::Scatter2D &hIn)
Book an ECorrelator in the same way as a histogram.
Definition: Correlators.hh:600
void rawHookOut(vector< MultiweightAOPtr > raos, size_t iW) final
Transform RAW ECorrelator Profiles to have content before writing them. Overloaded method from Analys...
Definition: Correlators.hh:946
CmpState compare(const Projection &p) const
Definition: Correlators.hh:145
ECorrPtr bookECorrelator(const string name, const vector< int > &h, vector< double > &binIn)
Book an ECorrelator in the same way as a histogram.
Definition: Correlators.hh:609
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 cnTwoInt(Scatter2DPtr h, ECorrPtr e2) const
Two-particle integrated cn.
Definition: Correlators.hh:897
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:462
vector< CorBin > getBins() const
Get a copy of the bin contents.
Definition: Correlators.hh:481
CorBin getReference() const
Extract the reference flow from a differential event averaged correlator.
Definition: Correlators.hh:513
vector< int > getH2() const
Get a copy of the h2 harmonic vector.
Definition: Correlators.hh:503
static pair< double, double > sampleVariance(T func)
Calculate the bootstrapped sample variance.
Definition: Correlators.hh:850
shared_ptr< ECorrelator > ECorrPtr
Typedef of shared pointer to ECorrelator.
Definition: Correlators.hh:596
ECorrPtr bookECorrelator(const string name, const vector< int > &h1, const vector< int > &h2, const YODA::Scatter2D &hIn)
Book a gapped ECorrelator with two harmonic vectors.
Definition: Correlators.hh:640
void vnSixInt(Scatter2DPtr h, ECorrPtr e2, ECorrPtr e4, ECorrPtr e6) const
Six particle integrated vn.
Definition: Correlators.hh:1080
ECorrelator(vector< int > h1In, vector< int > h2In, vector< double > binIn)
Constructor for gapped correlator.
Definition: Correlators.hh:417
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:392
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:431
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
ECorrPtr bookECorrelator(const string name, vector< double > binIn)
Templated version of correlator booking which takes N desired harmonic and M number of particles...
Definition: Correlators.hh:664
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()
ECorrPtr bookECorrelator(const string name, const vector< int > &h1, const vector< int > &h2, vector< double > &binIn)
Book a gapped ECorrelator with two harmonic vectors.
Definition: Correlators.hh:624
void fill(const Correlators &c, const double &weight=1.0)
Fill the bins with the appropriate correlator.
Definition: Correlators.hh:445
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:422
void setReference(CorBin refIn)
Replace reference flow bin with another, e.g. calculated in another phase space or with other pid...
Definition: Correlators.hh:508
void setProfs(vector< string > prIn)
Set the prIn list of profile histograms associated with the internal bins.
Definition: Correlators.hh:521
void vnFourDiff(Scatter2DPtr h, ECorrPtr e2Dif, ECorrPtr e4Dif) const
Four particle differential vn.
Definition: Correlators.hh:1176
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:745
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
void cnEightInt(Scatter2DPtr h, ECorrPtr e2, ECorrPtr e4, ECorrPtr e6, ECorrPtr e8) const
Eight particle integrated cn.
Definition: Correlators.hh:1088
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
vector< int > getH1() const
Get a copy of the h1 harmonic vector.
Definition: Correlators.hh:498
bool fillFromProfile(YODA::AnalysisObjectPtr yao, string name)
Fill bins with content from preloaded histograms.
Definition: Correlators.hh:526
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:815
void vnTwoDiff(Scatter2DPtr h, ECorrPtr e2Dif) const
Two particle differential vn.
Definition: Correlators.hh:1141
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:887
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()
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
vector< double > getBinX() const
Get a copy of the bin x-values.
Definition: Correlators.hh:493
const pair< int, int > getMaxValues() const
Get the correct max N and max P for the set of booked correlators.
Definition: Correlators.hh:579
CumulantAnalysis(const string &n)
Constructor.
Definition: Correlators.hh:702
void project(const Event &e)
void vnTwoInt(Scatter2DPtr h, ECorrPtr e2) const
Two particle integrated vn.
Definition: Correlators.hh:920
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:777
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:409
void vnEightInt(Scatter2DPtr h, ECorrPtr e2, ECorrPtr e4, ECorrPtr e6, ECorrPtr e8) const
Eight particle integrated vn.
Definition: Correlators.hh:1134
Base class for all Rivet projections.
Definition: Projection.hh:29
void corrPlot(Scatter2DPtr h, ECorrPtr e) const
Put an event-averaged correlator into a Scatter2D.
Definition: Correlators.hh:929
std::enable_if< std::is_arithmetic< NUM >::value, double >::type mean(const vector< NUM > &sample)
Definition: MathUtils.hh:458
ECorrPtr bookECorrelatorGap(const string name, const YODA::Scatter2D &hIn)
Templated version of gapped correlator booking which takes N desired harmonic and M number of particl...
Definition: Correlators.hh:682
void vnFourInt(Scatter2DPtr h, ECorrPtr e2, ECorrPtr e4) const
Four particle integrated vn.
Definition: Correlators.hh:1034
void cnSixInt(Scatter2DPtr h, ECorrPtr e2, ECorrPtr e4, ECorrPtr e6) const
Six particle integrated cn.
Definition: Correlators.hh:1041
ECorrPtr bookECorrelatorGap(const string name, const vector< int > &h, const YODA::Scatter2D &hIn)
Definition: Correlators.hh:652