Rivet 3.1.9
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
27namespace 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
62 using Projection::operator =;
63
64
72 const pair<double,double> intCorrelator(vector<int> n) const;
73
78 const vector<pair<double,double>> pTBinnedCorrelators(vector<int> n,
79 bool overflow = false) const;
80
89 const pair<double,double> intCorrelatorGap(const Correlators& other,
90 vector<int> n1, vector<int> n2) const;
91
103 const vector<pair<double,double>>
104 pTBinnedCorrelatorsGap(const Correlators& other, vector<int> n1, vector<int> n2, bool overflow=false) const;
105
110 static vector<int> hVec(int n, int m) {
111 if (m % 2 != 0) {
112 cout << "Harmonic Vector: Number of particles must be an even number." << endl;
113 return {};
114 }
115 vector<int> ret = {};
116 for (int i = 0; i < m; ++i) {
117 if (i < m/2) ret.push_back(n);
118 else ret.push_back(-n);
119 }
120 return ret;
121 }
122
125 static pair<int,int> getMaxValues(vector< vector<int> >& hList){
126 int nMax = 0, pMax = 0;
127 for (vector<int> h : hList) {
128 int tmpN = 0, tmpP = 0;
129 for ( int i = 0; i < int(h.size()); ++i) {
130 tmpN += abs(h[i]);
131 ++tmpP;
132 }
133 if (tmpN > nMax) nMax = tmpN;
134 if (tmpP > pMax) pMax = tmpP;
135 }
136 return make_pair(nMax,pMax);
137 }
138
139 // Clone on the heap.
140 DEFAULT_RIVET_PROJ_CLONE(Correlators);
141
142
143 protected:
144
145 // @brief Loop over array and calculates Q and P vectors if needed
146 void project(const Event& e);
147
148 // @brief Compare to other projection, testing harmonics, pT bins and underlying final state similarity
149 CmpState compare(const Projection& p) const {
150 const Correlators* other = dynamic_cast<const Correlators*>(&p);
151 if (nMax != other->nMax) return CmpState::NEQ;
152 if (pMax != other->pMax) return CmpState::NEQ;
153 if (pTbinEdges != other->pTbinEdges) return CmpState::NEQ;
154 return mkPCmp(p, "FS");
155 }
156
157 // @brief Calculate correlators from one particle
158 void fillCorrelators(const Particle& p, const double& weight);
159
160 // @brief Return a Q-vector.
161 const complex<double> getQ(int n, int p) const {
162 bool isNeg = (n < 0);
163 if (isNeg) return conj( qVec[abs(n)][p] );
164 else return qVec[n][p];
165 }
166
167 // @brief Return a P-vector
168 const complex<double> getP(int n, int p, double pT = 0.) const {
169 bool isNeg = (n < 0);
170 map<double, Vec2D>::const_iterator pTitr = pVec.lower_bound(pT);
171 if (pTitr == pVec.end()) return DBL_NAN;
172 if (isNeg) return conj( pTitr->second[abs(n)][p] );
173 else return pTitr->second[n][p];
174 }
175
176
177 private:
178
179 // @brief Find correlators by recursion
180 //
181 // Order = M (# of particles), n's are harmonics, p's are the powers of the weights
182 const complex<double> recCorr(int order, vector<int> n,
183 vector<int> p, bool useP, double pT = 0.) const;
184
185 // @brief Two-particle correlator
186 //
187 // Cf. eq. (19) p6. Flag if p-vectors or q-vectors should be used to
188 // calculate the correlator.
189 const complex<double> twoPartCorr(int n1, int n2, int p1 = 1,
190 int p2 = 1, double pT = 0., bool useP = false) const;
191
192 // Set elements in vectors to zero
193 void setToZero();
194
195 // Shorthands for setting and comparing to zero
196 const complex<double> _ZERO = {0., 0.};
197 const double _TINY = 1e-10;
198
199 // Shorthand typedefs for vec<vec<complex>>.
200 typedef vector< vector<complex<double>> > Vec2D;
201
202 // Define Q-vectors and P-vectors
203 Vec2D qVec; // Q[n][p]
204 map<double, Vec2D> pVec; // p[pT][n][p]
205
206 // The max values of n and p to be calculated.
207 int nMax, pMax;
208
209 // pT bin-edges
210 vector<double> pTbinEdges;
211
212 bool isPtDiff;
213
214 };
215
216
217
225 class CumulantAnalysis : public Analysis {
226 private:
227
228 // Number of bins used for bootstrap calculation of statistical
229 // uncertainties. It is hard coded, and shout NOT be changed unless there
230 // are good reasons to do so.
231 static const int BOOT_BINS = 9;
232
233 // Enum for choosing the method of error calculation.
234 enum Error {
235 VARIANCE,
236 ENVELOPE
237 };
238
239 // The desired error method. Can be changed in the analysis constructor
240 // by setting it appropriately.
241 Error errorMethod;
242
243
245 class CorBinBase {
246 public:
247 CorBinBase() {}
248 virtual ~CorBinBase() {};
249 // Derived class should have fill and mean defined.
250 virtual void fill(const pair<double, double>& cor, const double& weight = 1.0) = 0;
251 virtual double mean() const = 0;
252 };
253
254
260 class CorSingleBin : public CorBinBase {
261 public:
262
264 CorSingleBin()
265 : _sumWX(0.), _sumW(0.), _sumW2(0.), _numEntries(0.)
266 { }
267
268 // Destructor does nothing but must be implemented (?)
269 ~CorSingleBin() {}
270
274 void fill(const pair<double, double>& cor, const double& weight = 1.0) {
275 // Test if denominator for the single event average is zero.
276 if (cor.second < 1e-10) return;
277 // The single event average <M> is then cor.first / cor.second.
278 // With weights this becomes just:
279 _sumWX += cor.first * weight;
280 _sumW += weight * cor.second;
281 _sumW2 += weight * weight * cor.second * cor.second;
282 _numEntries += 1.;
283 }
284
286 double mean() const {
287 if (_sumW < 1e-10) return 0;
288 return _sumWX / _sumW;
289 }
290
292 double sumW() const {
293 return _sumW;
294 }
295
297 double sumW2() const {
298 return _sumW2;
299 }
300
302 double sumWX() const {
303 return _sumWX;
304 }
305
307 double numEntries() const {
308 return _numEntries;
309 }
310
312 void addContent(double ne, double sw, double sw2, double swx) {
313 _numEntries += ne;
314 _sumW += sw;
315 _sumW2 += sw2;
316 _sumWX += swx;
317 }
318
319
320 private:
321
322 double _sumWX, _sumW, _sumW2, _numEntries;
323
324 };
325
326
331 class CorBin : public CorBinBase {
332 public:
333
339 CorBin() : binIndex(0), nBins(BOOT_BINS) {
340 for(size_t i = 0; i < nBins; ++i) bins.push_back(CorSingleBin());
341 }
342
343 // Destructor does nothing but must be implemented (?)
344 ~CorBin() {}
345
347 void fill(const pair<double, double>& cor, const double& weight = 1.0) {
348 // Test if denominator for the single event average is zero.
349 if (cor.second < 1e-10) return;
350 // Fill the correct bin.
351 bins[binIndex].fill(cor, weight);
352 if (binIndex == nBins - 1) binIndex = 0;
353 else ++binIndex;
354 }
355
357 double mean() const {
358 double sow = 0;
359 double sowx = 0;
360 for (auto b : bins) {
361 if (b.sumW() < 1e-10) continue;
362 sow += b.sumW();
363 sowx += b.sumWX();
364 }
365 return sowx / sow;
366 }
367
369 vector<CorSingleBin> getBins() const {
370 return bins;
371 }
372
374 template<class T=CorBinBase>
375 vector<T*> getBinPtrs() {
376 vector<T*> ret(bins.size());
377 transform(bins.begin(), bins.end(), ret.begin(), [](CorSingleBin& b) {return &b;});
378 return ret;
379 }
380
381 private:
382
383 vector<CorSingleBin> bins;
384 size_t binIndex;
385 size_t nBins;
386
387 };
388
389
390 public:
391
397 public:
398
404 //ECorrelator(vector<int> h) : h1(h), h2({}),
405 // binX(0), binContent(0), reference() {
406 //};
407
413 ECorrelator(vector<int> h, vector<double> binIn)
414 : h1(h), h2({}), binX(binIn), binContent(binIn.size() - 1), reference()
415 { }
416
421 ECorrelator(vector<int> h1In, vector<int> h2In, vector<double> binIn)
422 : h1(h1In), h2(h2In), binX(binIn), binContent(binIn.size() - 1), reference()
423 { }
424
426 void fill(const double& obs, const Correlators& c, double weight=1.0) {
427 int index = getBinIndex(obs);
428 if (index < 0) return;
429 binContent[index].fill(c.intCorrelator(h1), weight);
430 }
431
435 void fill(const double& obs, const Correlators& c1, const Correlators& c2, double weight=1.0) {
436 if (!h2.size()) {
437 cout << "Trying to fill gapped correlator, but harmonics behind the gap (h2) are not given!" << endl;
438 return;
439 }
440 int index = getBinIndex(obs);
441 if (index < 0) return;
442 binContent[index].fill(c1.intCorrelatorGap(c2, h1, h2), weight);
443 }
444
449 void fill(const Correlators& c, const double& weight=1.0) {
450 vector< pair<double, double> > diffCorr = c.pTBinnedCorrelators(h1);
451 // We always skip overflow when calculating the all-event average.
452 if (diffCorr.size() != binX.size() - 1)
453 cout << "Tried to fill event with wrong binning (ungapped)" << endl;
454 for (size_t i = 0; i < diffCorr.size(); ++i) {
455 int index = getBinIndex(binX[i]);
456 if (index < 0) return;
457 binContent[index].fill(diffCorr[i], weight);
458 }
459 reference.fill(c.intCorrelator(h1), weight);
460 }
461
466 void fill(const Correlators& c1, const Correlators& c2, const double& weight = 1.0) {
467 if (!h2.size()) {
468 cout << "Trying to fill gapped correlator, but harmonics behind "
469 "the gap (h2) are not given!" << endl;
470 return;
471 }
472 vector< pair<double, double> > diffCorr = c1.pTBinnedCorrelatorsGap(c2, h1, h2);
473 // We always skip overflow when calculating the all event average.
474 if (diffCorr.size() != binX.size() - 1)
475 cout << "Tried to fill event with wrong binning (gapped)" << endl;
476 for (size_t i = 0; i < diffCorr.size(); ++i) {
477 int index = getBinIndex(binX[i]);
478 if (index < 0) return;
479 binContent[index].fill(diffCorr[i], weight);
480 }
481 reference.fill(c1.intCorrelatorGap(c2, h1, h2), weight);
482 }
483
485 vector<CorBin> getBins() const {
486 return binContent;
487 }
488
490 vector<CorBinBase*> getBinPtrs() {
491 vector<CorBinBase*> ret(binContent.size());
492 transform(binContent.begin(), binContent.end(), ret.begin(), [](CorBin& b) {return &b;});
493 return ret;
494 }
495
497 vector<double> getBinX() const {
498 return binX;
499 }
500
502 vector<int> getH1() const {
503 return h1;
504 }
505
507 vector<int> getH2() const {
508 return h2;
509 }
510
512 void setReference(CorBin refIn) {
513 reference = refIn;
514 }
515
517 CorBin getReference() const {
518 if (reference.mean() < 1e-10)
519 cout << "Warning: ECorrelator, reference bin is zero." << endl;
520 return reference;
521 }
522
525 void setProfs(vector<string> prIn) {
526 profs = prIn;
527 }
528
530 bool fillFromProfile(YODA::AnalysisObjectPtr yao, string name) {
531 auto refs = reference.getBinPtrs<CorSingleBin>();
532 for (size_t i = 0; i < profs.size(); ++i) {
533 if (yao->path() == "/RAW/"+name+"/TMP/"+profs[i]) {
534 YODA::Profile1DPtr pPtr = dynamic_pointer_cast<YODA::Profile1D>(yao);
535 for (size_t j = 0; j < binX.size() - 1; ++j) {
536 const YODA::ProfileBin1D& pBin = pPtr->binAt(binX[j]);
537 auto tmp = binContent[j].getBinPtrs<CorSingleBin>();
538 tmp[i]->addContent(pBin.numEntries(), pBin.sumW(), pBin.sumW2(),
539 pBin.sumWY());
540 }
541 // Get the reference flow from the underflow bin of the histogram.
542 const YODA::Dbn2D& uBin = pPtr->underflow();
543 refs[i]->addContent(uBin.numEntries(), uBin.sumW(), uBin.sumW2(),
544 uBin.sumWY());
545 return true;
546 }
547 } // End loop of bootstrapped correlators.
548 return false;
549 }
550
551 private:
552
553 // Get correct bin index for a given @param obs value
554 int getBinIndex(const double& obs) const {
555 // Find the correct index of binContent.
556 // If we are in overflow, just skip.
557 if (obs >= binX.back()) return -1;
558 // If we are in underflow, ditto.
559 if (obs < binX[0]) return -1;
560 int index = 0;
561 for (int i = 0, N = binX.size() - 1; i < N; ++i, ++index)
562 if (obs >= binX[i] && obs < binX[i + 1]) break;
563 return index;
564 }
565
566 // The harmonics vectors.
567 vector<int> h1;
568 vector<int> h2;
569
570 // The bins.
571 vector<double> binX;
572 vector<CorBin> binContent;
573 // The reference flow.
574 CorBin reference;
575 public:
576 // The profile histograms associated with the CorBins for streaming.
577 vector <string> profs;
578
579 };
580
581
583 const pair<int, int> getMaxValues() const {
584 vector< vector<int>> harmVecs;
585 for ( auto eItr = eCorrPtrs.begin(); eItr != eCorrPtrs.end(); ++eItr) {
586 vector<int> h1 = (*eItr)->getH1();
587 vector<int> h2 = (*eItr)->getH2();
588 if (h1.size() > 0) harmVecs.push_back(h1);
589 if (h2.size() > 0) harmVecs.push_back(h2);
590 }
591 if (harmVecs.size() == 0) {
592 cout << "Warning: You tried to extract max values from harmonic "
593 "vectors, but have not booked any." << endl;
594 return pair<int,int>();
595 }
596 return Correlators::getMaxValues(harmVecs);
597 }
598
600 typedef shared_ptr<ECorrelator> ECorrPtr;
601
604 ECorrPtr bookECorrelator(const string name, const vector<int>& h, const YODA::Scatter2D& hIn) {
605 vector<double> binIn;
606 for (auto b : hIn.points()) binIn.push_back(b.xMin());
607 binIn.push_back(hIn.points().back().xMax());
608 return bookECorrelator(name, h, binIn);
609 }
610
613 ECorrPtr bookECorrelator(const string name, const vector<int>& h, vector<double>& binIn) {
614 ECorrPtr ecPtr = ECorrPtr(new ECorrelator(h, binIn));
615 vector<string> eCorrProfs;
616 for (int i = 0; i < BOOT_BINS; ++i) {
617 Profile1DPtr tmp;
618 book(tmp,"TMP/"+name+"-"+to_string(i),binIn);
619 eCorrProfs.push_back(name+"-"+to_string(i));
620 }
621 ecPtr->setProfs(eCorrProfs);
622 eCorrPtrs.push_back(ecPtr);
623 return ecPtr;
624 }
625
628 ECorrPtr bookECorrelator(const string name, const vector<int>& h1,
629 const vector<int>& h2, vector<double>& binIn) {
630 ECorrPtr ecPtr = ECorrPtr(new ECorrelator(h1, h2, binIn));
631 vector<string> eCorrProfs;
632 for (int i = 0; i < BOOT_BINS; ++i) {
633 Profile1DPtr tmp;
634 book(tmp,"TMP/"+name+"-"+to_string(i),binIn);
635 eCorrProfs.push_back(name+"-"+to_string(i));
636 }
637 ecPtr->setProfs(eCorrProfs);
638 eCorrPtrs.push_back(ecPtr);
639 return ecPtr;
640 }
641
644 ECorrPtr bookECorrelator(const string name, const vector<int>& h1,
645 const vector<int>& h2, const YODA::Scatter2D& hIn ) {
646 vector<double> binIn;
647 for (auto b : hIn.points()) binIn.push_back(b.xMin());
648 binIn.push_back(hIn.points().back().xMax());
649 return bookECorrelator(name, h1, h2, binIn);
650 }
651
656 ECorrPtr bookECorrelatorGap(const string name, const vector<int>& h,
657 const YODA::Scatter2D& hIn) {
658 const vector<int> h1(h.begin(), h.begin() + h.size() / 2);
659 const vector<int> h2(h.begin() + h.size() / 2, h.end());
660 return bookECorrelator(name, h1, h2, hIn);
661 }
662
667 template<unsigned int N, unsigned int M>
668 ECorrPtr bookECorrelator(const string name, vector<double> binIn) {
669 return bookECorrelator(name, Correlators::hVec(N, M), binIn);
670 }
671
676 template<unsigned int N, unsigned int M>
677 ECorrPtr bookECorrelator(const string name, const YODA::Scatter2D& hIn) {
678 return bookECorrelator(name, Correlators::hVec(N, M), hIn);
679 }
680
685 template<unsigned int N, unsigned int M>
686 ECorrPtr bookECorrelatorGap(const string name, const YODA::Scatter2D& hIn) {
687 const vector<int> h = Correlators::hVec(N,M);
688 const vector<int> h1(h.begin(), h.begin() + h.size() / 2);
689 const vector<int> h2(h.begin() + h.size() / 2, h.end());
690 return bookECorrelator(name, h1, h2, hIn);
691
692 }
693
694 protected:
695
696 // Bookkeeping of the event averaged correlators.
697 list<ECorrPtr> eCorrPtrs;
698
699
700 public:
701
706 CumulantAnalysis(const string& n)
707 : Analysis(n), errorMethod(VARIANCE)
708 { }
709
717 template<typename T>
718 static void fillScatter(Scatter2DPtr h, vector<double>& binx, T func) {
719 vector<YODA::Point2D> points;
720 // Test if we have proper bins from a booked histogram.
721 bool hasBins = (h->points().size() > 0);
722 for (int i = 0, N = binx.size() - 1; i < N; ++i) {
723 double xMid = (binx[i] + binx[i + 1]) / 2.0;
724 double xeMin = fabs(xMid - binx[i]);
725 double xeMax = fabs(xMid - binx[i + 1]);
726 if (hasBins) {
727 xMid = h->points()[i].x();
728 xeMin = h->points()[i].xErrMinus();
729 xeMax = h->points()[i].xErrPlus();
730 }
731 double yVal = func(i);
732 if (std::isnan(yVal)) yVal = 0.;
733 double yErr = 0;
734 points.push_back(YODA::Point2D(xMid, yVal, xeMin, xeMax, yErr, yErr));
735 }
736 h->reset();
737 h->points().clear();
738 for (int i = 0, N = points.size(); i < N; ++i) h->addPoint(points[i]);
739 }
740
748 template<typename F>
749 void fillScatter(Scatter2DPtr h, vector<double>& binx, F func,
750 vector<pair<double, double> >& yErr) const {
751 vector<YODA::Point2D> points;
752 // Test if we have proper bins from a booked histogram.
753 bool hasBins = (h->points().size() > 0);
754 for (int i = 0, N = binx.size() - 1; i < N; ++i) {
755 double xMid = (binx[i] + binx[i + 1]) / 2.0;
756 double xeMin = fabs(xMid - binx[i]);
757 double xeMax = fabs(xMid - binx[i + 1]);
758 if (hasBins) {
759 xMid = h->points()[i].x();
760 xeMin = h->points()[i].xErrMinus();
761 xeMax = h->points()[i].xErrPlus();
762 }
763 double yVal = func(i);
764 if (std::isnan(yVal))
765 points.push_back(YODA::Point2D(xMid, 0., xeMin, xeMax,0., 0.));
766 else
767 points.push_back(YODA::Point2D(xMid, yVal, xeMin, xeMax,
768 yErr[i].first, yErr[i].second));
769 }
770 h->reset();
771 h->points().clear();
772
773 for (int i = 0, N = points.size(); i < N; ++i)
774 h->addPoint(points[i]);
775 }
776
777
781 static void nthPow(Scatter2DPtr hOut, const Scatter2DPtr hIn, const double& n, const double& k = 1.0) {
782 if (n == 0 || n == 1) {
783 cout << "Error: Do not take the 0th or 1st power of a Scatter2D,"
784 " use scale instead." << endl;
785 return;
786 }
787 if (hIn->points().size() != hOut->points().size()) {
788 cout << "nthRoot: Scatterplots: " << hIn->name() << " and " <<
789 hOut->name() << " not initialized with same length." << endl;
790 return;
791 }
792 vector<YODA::Point2D> points;
793 // The error pre-factor is k^(1/n) / n by Taylors formula.
794 double eFac = pow(k,1./n)/n;
795 for (auto b : hIn->points()) {
796 double yVal = pow(k * b.y(),n);
797 if (std::isnan(yVal))
798 points.push_back(YODA::Point2D(b.x(), 0., b.xErrMinus(),
799 b.xErrPlus(), 0, 0 ));
800 else {
801 double yemin = abs(eFac * pow(yVal,1./(n - 1.))) * b.yErrMinus();
802 if (std::isnan(yemin)) yemin = b.yErrMinus();
803 double yemax = abs(eFac * pow(yVal,1./(n - 1.))) * b.yErrPlus();
804 if (std::isnan(yemax)) yemax = b.yErrPlus();
805 points.push_back(YODA::Point2D(b.x(), yVal, b.xErrMinus(),
806 b.xErrPlus(), yemin, yemax ));
807 }
808 }
809 hOut->reset();
810 hOut->points().clear();
811 for (int i = 0, N = points.size(); i < N; ++i)
812 hOut->addPoint(points[i]);
813 }
814
815
819 static void nthPow(Scatter2DPtr h, const double& n, const double& k = 1.0) {
820 if (n == 0 || n == 1) {
821 cout << "Error: Do not take the 0th or 1st power of a Scatter2D,"
822 " use scale instead." << endl;
823 return;
824 }
825 vector<YODA::Point2D> points;
826 vector<YODA::Point2D> pIn = h->points();
827 // The error pre-factor is k^(1/n) / n by Taylors formula.
828 double eFac = pow(k,1./n)/n;
829 for (auto b : pIn) {
830 double yVal = pow(k * b.y(),n);
831 if (std::isnan(yVal))
832 points.push_back(YODA::Point2D(b.x(), 0., b.xErrMinus(),
833 b.xErrPlus(), 0, 0 ));
834 else {
835 double yemin = abs(eFac * pow(yVal,1./(n - 1.))) * b.yErrMinus();
836 if (std::isnan(yemin)) yemin = b.yErrMinus();
837 double yemax = abs(eFac * pow(yVal,1./(n - 1.))) * b.yErrPlus();
838 if (std::isnan(yemax)) yemax = b.yErrPlus();
839 points.push_back(YODA::Point2D(b.x(), yVal, b.xErrMinus(),
840 b.xErrPlus(), yemin, yemax ));
841 }
842 }
843 h->reset();
844 h->points().clear();
845 for (int i = 0, N = points.size(); i < N; ++i) h->addPoint(points[i]);
846 }
847
848
853 template<typename T>
854 static pair<double, double> sampleVariance(T func) {
855 // First we calculate the mean (two pass calculation).
856 double avg = 0.;
857 for (int i = 0; i < BOOT_BINS; ++i) avg += func(i);
858 avg /= double(BOOT_BINS);
859 // Then we find the variance.
860 double var = 0.;
861 for (int i = 0; i < BOOT_BINS; ++i) var += pow(func(i) - avg, 2.);
862 var /= (double(BOOT_BINS) - 1);
863 return pair<double, double>(var, var);
864 }
865
866
871 template<typename T>
872 static pair<double, double> sampleEnvelope(T func) {
873 // First we calculate the mean.
874 double avg = 0.;
875 for (int i = 0; i < BOOT_BINS; ++i) avg += func(i);
876 avg /= double(BOOT_BINS);
877 double yMax = avg;
878 double yMin = avg;
879 // Then we find the envelope using the mean as initial value.
880 for (int i = 0; i < BOOT_BINS; ++i) {
881 double yVal = func(i);
882 if (yMin > yVal) yMin = yVal;
883 else if (yMax < yVal) yMax = yVal;
884 }
885 return pair<double, double>(fabs(avg - yMin), fabs(yMax - avg));
886 }
887
888
890 template<typename T>
891 const pair<double, double> sampleError(T func) const {
892 if (errorMethod == VARIANCE) return sampleVariance(func);
893 else if (errorMethod == ENVELOPE) return sampleEnvelope(func);
894 else
895 cout << "Error: Error method not found!" << endl;
896 return pair<double, double>(0.,0.);
897 }
898
899
901 void cnTwoInt(Scatter2DPtr h, ECorrPtr e2) const {
902 vector<CorBin> bins = e2->getBins();
903 vector<double> binx = e2->getBinX();
904 // Assert bin size.
905 if (binx.size() - 1 != bins.size()){
906 cout << "cnTwoInt: Bin size (x,y) differs!" << endl;
907 return;
908 }
909 vector<CorBinBase*> binPtrs;
910 // The mean value of the cumulant.
911 auto cn = [&] (int i) { return binPtrs[i]->mean(); };
912 // Error calculation.
913 vector<pair<double, double> > yErr;
914 for (int j = 0, N = bins.size(); j < N; ++j) {
915 binPtrs = bins[j].getBinPtrs();
916 yErr.push_back(sampleError(cn));
917 }
918 binPtrs = e2->getBinPtrs();
919 fillScatter(h, binx, cn, yErr);
920 }
921
922
924 void vnTwoInt(Scatter2DPtr h, ECorrPtr e2) const {
925 cnTwoInt(h, e2);
926 nthPow(h, 0.5);
927 }
928
929
933 void corrPlot(Scatter2DPtr h, ECorrPtr e) const {
934 cnTwoInt(h, e);
935 }
936
937
938
939
940 // TODO Use full path for lookup, change to single AU in output, rename.
941 void rawHookIn(YODA::AnalysisObjectPtr yao) final {
942 // Fill the corresponding ECorrelator.
943 for (auto ec : eCorrPtrs) if(ec->fillFromProfile(yao, name())) break;;
944 }
945
950 void rawHookOut(vector<MultiweightAOPtr> raos, size_t iW) final {
951 // Loop over the correlators and extract the numbers.
952 for (auto ec : eCorrPtrs) {
953 vector<CorBin> corBins = ec->getBins();
954 vector<double> binx = ec->getBinX();
955 auto ref = ec->getReference();
956 auto refBins = ref.getBinPtrs<CorSingleBin>();
957 // Assert bin size.
958 if (binx.size() - 1 != corBins.size()){
959 cout << "corrPlot: Bin size (x,y) differs!" << endl;
960 return;
961 }
962 // Loop over the booked histograms using their names.
963 for (int i = 0, N = ec->profs.size(); i < N; ++i) {
964 for (auto rao : raos) {
965 if (rao->path() == "/"+name()+"/TMP/"+ec->profs[i]) {
966 // Get a pointer to the active profile.
967 rao.get()->setActiveWeightIdx(iW);
968 YODA::Profile1DPtr pPtr = dynamic_pointer_cast<YODA::Profile1D>(
969 rao.get()->activeYODAPtr());
970 // New bins.
971 vector<YODA::ProfileBin1D> profBins;
972 // Numbers for the summary distribution
973 double ne = 0., sow = 0., sow2 = 0.;
974 for (size_t j = 0, N = binx.size() - 1; j < N; ++j) {
975 vector<CorSingleBin*> binPtrs =
976 corBins[j].getBinPtrs<CorSingleBin>();
977 // Construct bin of the profiled quantities. We have no information
978 // (and no desire to add it) of sumWX of the profile, so really
979 // we should use a Dbn1D - but that does not work for Profile1D's.
980 profBins.push_back( YODA::ProfileBin1D(pPtr->bin(j).xEdges(),
981 YODA::Dbn2D( binPtrs[i]->numEntries(), binPtrs[i]->sumW(),
982 binPtrs[i]->sumW2(), 0., 0., binPtrs[i]->sumWX(), 0, 0)));
983 ne += binPtrs[i]->numEntries();
984 sow += binPtrs[i]->sumW();
985 sow2 += binPtrs[i]->sumW2();
986 }
987 // Put the ECorrelator into the raw histogram.
988 pPtr->reset();
989 pPtr->bins().clear();
990 // Add the bins.
991 pPtr->addBins(profBins);
992 // Set the total distribution.
993 pPtr->setTotalDbn(YODA::Dbn2D(ne,sow,sow2,0.,0.,0.,0.,0.));
994 // And reference flow in the underflow bin.
995 pPtr->setUnderflow(YODA::Dbn2D(refBins[i]->numEntries(),
996 refBins[i]->sumW(), refBins[i]->sumW2(), 0., 0.,
997 refBins[i]->sumWX(), 0., 0.));
998 }
999 }
1000 }
1001 }
1002 }
1003
1004 // @brief Four particle integrated cn.
1005 void cnFourInt(Scatter2DPtr h, ECorrPtr e2, ECorrPtr e4) const {
1006 auto e2bins = e2->getBins();
1007 auto e4bins = e4->getBins();
1008 auto binx = e2->getBinX();
1009 if (binx.size() - 1 != e2bins.size()){
1010 cout << "cnFourInt: Bin size (x,y) differs!" << endl;
1011 return;
1012 }
1013 if (binx != e4->getBinX()) {
1014 cout << "Error in cnFourInt: Correlator x-binning differs!" << endl;
1015 return;
1016 }
1017 vector<CorBinBase*> e2binPtrs;
1018 vector<CorBinBase*> e4binPtrs;
1019 auto cn = [&] (int i) {
1020 double e22 = e2binPtrs[i]->mean() * e2binPtrs[i]->mean();
1021 return e4binPtrs[i]->mean() - 2. * e22;
1022 };
1023 // Error calculation.
1024 vector<pair<double, double> > yErr;
1025 for (int j = 0, N = e2bins.size(); j < N; ++j) {
1026 e2binPtrs = e2bins[j].getBinPtrs();
1027 e4binPtrs = e4bins[j].getBinPtrs();
1028 yErr.push_back(sampleError(cn));
1029 }
1030 // Put the bin ptrs back in place.
1031 e2binPtrs = e2->getBinPtrs();
1032 e4binPtrs = e4->getBinPtrs();
1033 fillScatter(h, binx, cn, yErr);
1034 }
1035
1036
1038 void vnFourInt(Scatter2DPtr h, ECorrPtr e2, ECorrPtr e4) const {
1039 cnFourInt(h, e2, e4);
1040 nthPow(h, 0.25, -1.0);
1041 }
1042
1043
1045 void cnSixInt(Scatter2DPtr h, ECorrPtr e2, ECorrPtr e4,
1046 ECorrPtr e6) const {
1047 auto e2bins = e2->getBins();
1048 auto e4bins = e4->getBins();
1049 auto e6bins = e6->getBins();
1050 auto binx = e2->getBinX();
1051 if (binx.size() - 1 != e2bins.size()){
1052 cout << "cnSixInt: Bin size (x,y) differs!" << endl;
1053 return;
1054 }
1055 if (binx != e4->getBinX() || binx != e6->getBinX()) {
1056 cout << "Error in cnSixInt: Correlator x-binning differs!" << endl;
1057 return;
1058 }
1059 vector<CorBinBase*> e2binPtrs;
1060 vector<CorBinBase*> e4binPtrs;
1061 vector<CorBinBase*> e6binPtrs;
1062 auto cn = [&] (int i) {
1063 double e23 = pow(e2binPtrs[i]->mean(), 3.0);
1064 return e6binPtrs[i]->mean() - 9.*e2binPtrs[i]->mean()*e4binPtrs[i]->mean() +
1065 12.*e23;
1066 };
1067 // Error calculation.
1068 vector<pair<double, double> > yErr;
1069 for (int j = 0, N = e2bins.size(); j < N; ++j) {
1070 e2binPtrs = e2bins[j].getBinPtrs();
1071 e4binPtrs = e4bins[j].getBinPtrs();
1072 e6binPtrs = e6bins[j].getBinPtrs();
1073 yErr.push_back(sampleError(cn));
1074 }
1075 // Put the bin ptrs back in place.
1076 e2binPtrs = e2->getBinPtrs();
1077 e4binPtrs = e4->getBinPtrs();
1078 e6binPtrs = e6->getBinPtrs();
1079 fillScatter(h, binx, cn, yErr);
1080 }
1081
1082
1084 void vnSixInt(Scatter2DPtr h, ECorrPtr e2, ECorrPtr e4,
1085 ECorrPtr e6) const {
1086 cnSixInt(h, e2, e4, e6);
1087 nthPow(h, 1./6., 0.25);
1088 }
1089
1090
1092 void cnEightInt(Scatter2DPtr h, ECorrPtr e2, ECorrPtr e4,
1093 ECorrPtr e6, ECorrPtr e8) const {
1094 auto e2bins = e2->getBins();
1095 auto e4bins = e4->getBins();
1096 auto e6bins = e6->getBins();
1097 auto e8bins = e8->getBins();
1098 auto binx = e2->getBinX();
1099 if (binx.size() - 1 != e2bins.size()){
1100 cout << "cnEightInt: Bin size (x,y) differs!" << endl;
1101 return;
1102 }
1103 if (binx != e4->getBinX() || binx != e6->getBinX() ||
1104 binx != e8->getBinX()) {
1105 cout << "Error in cnEightInt: Correlator x-binning differs!" << endl;
1106 return;
1107 }
1108 vector<CorBinBase*> e2binPtrs;
1109 vector<CorBinBase*> e4binPtrs;
1110 vector<CorBinBase*> e6binPtrs;
1111 vector<CorBinBase*> e8binPtrs;
1112 auto cn = [&] (int i ) {
1113 double e22 = e2binPtrs[i]->mean() * e2binPtrs[i]->mean();
1114 double e24 = e22 * e22;
1115 double e42 = e4binPtrs[i]->mean() * e4binPtrs[i]->mean();
1116 return e8binPtrs[i]->mean() - 16. * e6binPtrs[i]->mean() *
1117 e2binPtrs[i]->mean() - 18. * e42 + 144. * e4binPtrs[i]->mean()*e22 - 144. * e24;
1118 };
1119 // Error calculation.
1120 vector<pair<double, double> > yErr;
1121 for (int j = 0, N = e2bins.size(); j < N; ++j) {
1122 e2binPtrs = e2bins[j].getBinPtrs();
1123 e4binPtrs = e4bins[j].getBinPtrs();
1124 e6binPtrs = e6bins[j].getBinPtrs();
1125 e8binPtrs = e8bins[j].getBinPtrs();
1126 yErr.push_back(sampleError(cn));
1127 }
1128 // Put the bin ptrs back in place.
1129 e2binPtrs = e2->getBinPtrs();
1130 e4binPtrs = e4->getBinPtrs();
1131 e6binPtrs = e6->getBinPtrs();
1132 e8binPtrs = e8->getBinPtrs();
1133 fillScatter(h, binx, cn, yErr);
1134 }
1135
1136
1138 void vnEightInt(Scatter2DPtr h, ECorrPtr e2, ECorrPtr e4, ECorrPtr e6, ECorrPtr e8) const {
1139 cnEightInt(h, e2, e4, e6, e8);
1140 nthPow(h, 1./8., -1./33.);
1141 }
1142
1143
1145 void vnTwoDiff(Scatter2DPtr h, ECorrPtr e2Dif) const {
1146 auto e2bins = e2Dif->getBins();
1147 auto ref = e2Dif->getReference();
1148 auto binx = e2Dif->getBinX();
1149 if (binx.size() -1 != e2bins.size()) {
1150 cout << "vnTwoDif: Bin size (x,y) differs!" << endl;
1151 return;
1152 }
1153 vector<CorBinBase*> e2binPtrs;
1154 vector<CorBinBase*> refPtrs;
1155 auto vn = [&] (int i) {
1156 // Test reference flow.
1157 if (ref.mean() <= 0) return 0.;
1158 return e2binPtrs[i]->mean() / sqrt(ref.mean());
1159 };
1160 // We need here a separate error function, as we don't iterate over the reference flow.
1161 auto vnerr = [&] (int i) {
1162 // Test reference flow.
1163 if (refPtrs[i]->mean() <=0) return 0.;
1164 return e2binPtrs[i]->mean() / sqrt(refPtrs[i]->mean());
1165 };
1166 // Error calculation.
1167 vector<pair<double, double> > yErr;
1168 refPtrs = ref.getBinPtrs();
1169 for (int j = 0, N = e2bins.size(); j < N; ++j) {
1170 e2binPtrs = e2bins[j].getBinPtrs();
1171 yErr.push_back(sampleError(vnerr));
1172 }
1173 // Put the e2binPtrs back in place.
1174 e2binPtrs = e2Dif->getBinPtrs();
1175 fillScatter(h, binx, vn);
1176 }
1177
1178
1180 void vnFourDiff(Scatter2DPtr h, ECorrPtr e2Dif, ECorrPtr e4Dif) const {
1181 auto e2bins = e2Dif->getBins();
1182 auto e4bins = e4Dif->getBins();
1183 auto ref2 = e2Dif->getReference();
1184 auto ref4 = e4Dif->getReference();
1185 auto binx = e2Dif->getBinX();
1186 if (binx.size() - 1 != e2bins.size()){
1187 cout << "vnFourDif: Bin size (x,y) differs!" << endl;
1188 return;
1189 }
1190 if (binx != e4Dif->getBinX()) {
1191 cout << "Error in vnFourDif: Correlator x-binning differs!" << endl;
1192 return;
1193 }
1194 vector<CorBinBase*> e2binPtrs;
1195 vector<CorBinBase*> e4binPtrs;
1196 vector<CorBinBase*> ref2Ptrs;
1197 vector<CorBinBase*> ref4Ptrs;
1198 double denom = 2 * ref2.mean() * ref2.mean() - ref4.mean();
1199 auto vn = [&] (int i) {
1200 // Test denominator.
1201 if (denom <= 0 ) return 0.;
1202 return ((2 * ref2.mean() * e2bins[i].mean() - e4bins[i].mean()) / pow(denom, 0.75));
1203 };
1204 // We need here a separate error function, as we don't iterate over the reference flow.
1205 auto vnerr = [&] (int i) {
1206 double denom2 = 2 * ref2Ptrs[i]->mean() * ref2Ptrs[i]->mean() -
1207 ref4Ptrs[i]->mean();
1208 // Test denominator.
1209 if (denom2 <= 0) return 0.;
1210 return ((2 * ref2Ptrs[i]->mean() * e2binPtrs[i]->mean() - e4binPtrs[i]->mean()) / pow(denom2, 0.75));
1211 };
1212 // Error calculation.
1213 vector<pair<double, double> > yErr;
1214 ref2Ptrs = ref2.getBinPtrs();
1215 ref4Ptrs = ref4.getBinPtrs();
1216 for (int j = 0, N = e2bins.size(); j < N; ++j) {
1217 e2binPtrs = e2bins[j].getBinPtrs();
1218 e4binPtrs = e4bins[j].getBinPtrs();
1219 yErr.push_back(sampleError(vnerr));
1220 }
1221 // Put the binPtrs back in place.
1222 e2binPtrs = e2Dif->getBinPtrs();
1223 e4binPtrs = e4Dif->getBinPtrs();
1224 fillScatter(h, binx, vn, yErr);
1225 }
1226
1227 };
1228
1229
1230}
1231
1232#endif
This is the base class of all analysis classes in Rivet.
Definition Analysis.hh:65
double sumW2() const
Get the sum of squared event weights seen (via the analysis handler).
double sumW() const
Get the sum of event weights seen (via the analysis handler).
Projection for calculating correlators for flow measurements.
Definition Correlators.hh:38
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()
CmpState compare(const Projection &p) const
Definition Correlators.hh:149
void project(const Event &e)
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 pair< double, double > intCorrelator(vector< int > n) const
Integrated correlator of n harmonic, with the number of powers being the size of n....
Correlators(const ParticleFinder &fsp, int nMaxIn=2, int pMaxIn=0, vector< double > pTbinEdgesIn={})
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()
static vector< int > hVec(int n, int m)
Construct a harmonic vectors from n harmonics and m number of particles.
Definition Correlators.hh:110
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,...
Definition Correlators.hh:125
A helper class to calculate all event averages of correlators.
Definition Correlators.hh:396
vector< double > getBinX() const
Get a copy of the bin x-values.
Definition Correlators.hh:497
void setProfs(vector< string > prIn)
Set the prIn list of profile histograms associated with the internal bins.
Definition Correlators.hh:525
void setReference(CorBin refIn)
Replace reference flow bin with another, e.g. calculated in another phase space or with other pid.
Definition Correlators.hh:512
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:413
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:466
vector< CorBin > getBins() const
Get a copy of the bin contents.
Definition Correlators.hh:485
vector< int > getH2() const
Get a copy of the h2 harmonic vector.
Definition Correlators.hh:507
ECorrelator(vector< int > h1In, vector< int > h2In, vector< double > binIn)
Constructor for gapped correlator.
Definition Correlators.hh:421
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:426
bool fillFromProfile(YODA::AnalysisObjectPtr yao, string name)
Fill bins with content from preloaded histograms.
Definition Correlators.hh:530
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:435
CorBin getReference() const
Extract the reference flow from a differential event averaged correlator.
Definition Correlators.hh:517
void fill(const Correlators &c, const double &weight=1.0)
Fill the bins with the appropriate correlator.
Definition Correlators.hh:449
vector< int > getH1() const
Get a copy of the h1 harmonic vector.
Definition Correlators.hh:502
vector< CorBinBase * > getBinPtrs()
Return the bins as pointers to the base class.
Definition Correlators.hh:490
Tools for flow analyses.
Definition Correlators.hh:225
void vnTwoDiff(Scatter2DPtr h, ECorrPtr e2Dif) const
Two particle differential vn.
Definition Correlators.hh:1145
void vnSixInt(Scatter2DPtr h, ECorrPtr e2, ECorrPtr e4, ECorrPtr e6) const
Six particle integrated vn.
Definition Correlators.hh:1084
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:628
void vnFourDiff(Scatter2DPtr h, ECorrPtr e2Dif, ECorrPtr e4Dif) const
Four particle differential vn.
Definition Correlators.hh:1180
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:668
const pair< double, double > sampleError(T func) const
Selection method for which sample error to use, given in the constructor.
Definition Correlators.hh:891
const pair< int, int > getMaxValues() const
Get the correct max N and max P for the set of booked correlators.
Definition Correlators.hh:583
void cnSixInt(Scatter2DPtr h, ECorrPtr e2, ECorrPtr e4, ECorrPtr e6) const
Six particle integrated cn.
Definition Correlators.hh:1045
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:686
CumulantAnalysis(const string &n)
Constructor.
Definition Correlators.hh:706
void vnFourInt(Scatter2DPtr h, ECorrPtr e2, ECorrPtr e4) const
Four particle integrated vn.
Definition Correlators.hh:1038
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:749
static void fillScatter(Scatter2DPtr h, vector< double > &binx, T func)
Helper method for turning correlators into Scatter2Ds.
Definition Correlators.hh:718
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:644
void corrPlot(Scatter2DPtr h, ECorrPtr e) const
Put an event-averaged correlator into a Scatter2D.
Definition Correlators.hh:933
void cnTwoInt(Scatter2DPtr h, ECorrPtr e2) const
Two-particle integrated cn.
Definition Correlators.hh:901
static pair< double, double > sampleVariance(T func)
Calculate the bootstrapped sample variance.
Definition Correlators.hh:854
static pair< double, double > sampleEnvelope(T func)
Calculate the bootstrapped sample envelope.
Definition Correlators.hh:872
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:781
ECorrPtr bookECorrelatorGap(const string name, const vector< int > &h, const YODA::Scatter2D &hIn)
Definition Correlators.hh:656
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:677
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:613
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:819
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:604
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:950
void vnTwoInt(Scatter2DPtr h, ECorrPtr e2) const
Two particle integrated vn.
Definition Correlators.hh:924
void vnEightInt(Scatter2DPtr h, ECorrPtr e2, ECorrPtr e4, ECorrPtr e6, ECorrPtr e8) const
Eight particle integrated vn.
Definition Correlators.hh:1138
shared_ptr< ECorrelator > ECorrPtr
Typedef of shared pointer to ECorrelator.
Definition Correlators.hh:600
void cnEightInt(Scatter2DPtr h, ECorrPtr e2, ECorrPtr e4, ECorrPtr e6, ECorrPtr e8) const
Eight particle integrated cn.
Definition Correlators.hh:1092
Representation of a HepMC event, and enabler of Projection caching.
Definition Event.hh:22
Base class for projections which return subsets of an event's particles.
Definition ParticleFinder.hh:11
Particle representation, either from a HepMC::GenEvent or reconstructed.
Definition Particle.hh:53
Base class for all Rivet projections.
Definition Projection.hh:29
Cmp< Projection > mkPCmp(const Projection &otherparent, const std::string &pname) const
CounterPtr & book(CounterPtr &, const std::string &name)
Book a counter.
virtual std::string name() const
Get the name of the analysis.
Definition Analysis.hh:131
double pT(const ParticleBase &p)
Unbound function access to pT.
Definition ParticleBaseUtils.hh:687
double p(const ParticleBase &p)
Unbound function access to p.
Definition ParticleBaseUtils.hh:684
Definition MC_Cent_pPb.hh:10
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:456
static constexpr double DBL_NAN
Convenient const for getting the double NaN value.
Definition Utils.hh:37
std::enable_if< std::is_arithmetic< NUM >::value, double >::type mean(const vector< NUM > &sample)
Definition MathUtils.hh:498