AliPhysics  88b7ad0 (88b7ad0)
FastCentEstimators.C
Go to the documentation of this file.
1 #ifndef FASTCENTESTIMATOR_C
2 #define FASTCENTESTIMATOR_C
3 #ifndef __CINT__
4 # include <TObject.h>
5 # include <TString.h>
6 # include <TH1.h>
7 # include <TMath.h>
8 # include <TCollection.h>
9 # include <TTree.h>
10 # include <TParticle.h>
11 # include <TError.h>
12 # include "FastShortHeader.C"
13 # include <TH2.h>
14 #else
15 class TH1;
16 // class TH2;
17 class TCollection;
18 class TTree;
19 class TParticle;
20 class FastShortHeader;
21 #endif
22 
23 enum {
24  kNeutron = 2112,
25  kProton = 2212
26 };
27 
28 //====================================================================
32 struct FastCentEstimator : public TObject
33 {
41  FastCentEstimator(const char* name="")
42  : TObject(), fName(name), fVerbose(false)
43  {}
47  virtual ~FastCentEstimator() {}
53  const char* GetName() const { return fName.Data(); }
54  void SetVerbose(Bool_t verb) { fVerbose = verb; }
65  virtual void Setup(TCollection* out, TTree* tree, UShort_t sNN,
66  Bool_t tgtA, Bool_t projA) = 0;
71  virtual void PreEvent() {}
72  virtual void ProcessHeader(FastShortHeader&) {}
80  virtual void Process(const TParticle* p) = 0;
84  virtual void PostEvent() {}
90  virtual void Terminate(TCollection* out) = 0;
91  virtual void Print(Option_t* option="") const
92  {
93  Printf("%s: %s", ClassName(), GetName());
94  }
102  static Double_t Theta(const TParticle* p)
103  {
104  Double_t pT = p->Pt();
105  Double_t pZ = p->Pz();
106  Double_t theta = TMath::ATan2(pT, pZ);
107  return theta;
108  }
116  static Double_t Eta(const TParticle* p)
117  {
118  Double_t theta = Theta(p);
119  Double_t tanth = TMath::Tan(theta/2);
120  if (tanth < 1e-100) return +1e10;
121  if (tanth > 1e100) return -1e10;
122  Double_t eta = -TMath::Log(tanth);
123  return eta;
124  }
132  static Double_t Phi(const TParticle* p)
133  {
134  Double_t px = p->Px();
135  Double_t py = p->Py();
136  Double_t phi = TMath::ATan2(px, py);
137  return phi;
138  }
146  static Bool_t IsPrimary(const TParticle* p)
147  {
148  return p->TestBit(BIT(14));
149  }
157  static Bool_t IsWeakDecay(const TParticle* p)
158  {
159  return p->TestBit(BIT(15));
160  }
168  static Bool_t IsCharged(const TParticle* p)
169  {
170  return p->TestBit(BIT(16));
171  }
172  ClassDef(FastCentEstimator,1);
173 };
174 
175 
176 
177 
178 //____________________________________________________________________
183 {
185  ULong64_t fCache;
197  Fast1DCentEstimator(const char* name="")
198  : FastCentEstimator(name), fHistogram(0), fFromTop(true)
199  {}
203  virtual ~Fast1DCentEstimator() {}
204  virtual const char* MultSpec() const { return "l"; }
212  virtual void Setup(TCollection* l, TTree* tree, UShort_t,
213  Bool_t, Bool_t)
214  {
215  if (fHistogram && l) {
216  Info("Setup", "Adding histogram %s to output",
217  fHistogram->GetName());
218  l->Add(fHistogram);
219  }
220  TString leaves; leaves.Form("value/%s", MultSpec());
221  if (tree) tree->Branch(GetName(), &fCache, leaves.Data());
222  }
226  virtual void PreEvent()
227  {
228  fCache = 0;
229  }
233  virtual void PostEvent()
234  {
235  if (fVerbose) Info("PostEvent", " Got %lld %s particles",
236  fCache, GetName());
237  fHistogram->Fill(fCache);
238  }
239  virtual TH1* GetHistogram(TCollection* l) = 0;
249  virtual void Terminate(TCollection* out)
250  {
251  TH1* h = GetHistogram(out);
252  TH1* cent = static_cast<TH1*>(h->Clone(GetName()));
253  cent->SetDirectory(0);
254  cent->SetYTitle("Centrality [%]");
255  cent->SetTitle(Form("%s mapping", GetName()));
256  cent->Reset();
257  out->Add(cent);
258 
259  Int_t nX = h->GetNbinsX();
260  Double_t total = h->Integral(1,nX);
261  Int_t dBin = (fFromTop ? -1 : 1);
262  Int_t end = (fFromTop ? 0 : nX+1);
263  Int_t start = (fFromTop ? nX : 1);
264  if (fVerbose)
265  Info("Teminate", "Integrating %s from bin %d to 1",
266  h->GetName(), nX);
267  for (Int_t i = start; i != end; i += dBin) {
268  Double_t curInt = (fFromTop ?
269  h->Integral(i, start) :
270  h->Integral(start,i));
271  if (curInt < 0) continue;
272  Double_t curCent = curInt / total * 100;
273  cent->SetBinContent(i, curCent);
274  if (fVerbose)
275  Info("Terminate", "Bin %3d -> %9f/%9f -> %5.1f%%",
276  i, curInt, total, curCent);
277  }
278  }
279  virtual void Print(Option_t* option="nh") const
280  {
281  TString opt(option); opt.ToLower();
282  if (opt.Contains("n"))
283  Printf("1D Estimator: %s (%screasing)",
284  GetName(),(fFromTop ? "de" : "in"));
285  if (opt.Contains("h") && fHistogram) {
286  Int_t nBin = fHistogram->GetNbinsX();
287  Printf(" %d bins between %f and %f",nBin,
288  fHistogram->GetXaxis()->GetXmin(),
289  fHistogram->GetXaxis()->GetXmax());
290  }
291  }
292  ClassDef(Fast1DCentEstimator,1);
293 };
294 
295 //____________________________________________________________________
297 {
298  // TH2* fBvsC;
301  : Fast1DCentEstimator("B"),
302  // fBvsC(0),
303  fkFactor(1000)
304  {}
314  virtual void Setup(TCollection* out, TTree* tree, UShort_t sNN,
315  Bool_t tgtA, Bool_t projA)
316  {
317  fHistogram = MakeHistogram(sNN, tgtA, projA);
318  Fast1DCentEstimator::Setup(out, tree, sNN, tgtA, projA);
319  // if (tree) tree->Branch(GetName(), &fB, "value/D");
320  }
327  {
328  return static_cast<TH1*>(l->FindObject(Form("raw%s",GetName())));
329  }
331  {
332  TArrayD cents;
333  TArrayD bins; // In B
334  if (tgtA && projA) { // Pb-Pb
335  if (sNN == 2760) {
336  // PbPb @ 2.76TeV only
337  // Updated 4th of November 2014 from
338  // cern.ch/twiki/bin/view/ALICE/CentStudies
339  // #Tables_with_centrality_bins_AN1
340  Double_t bs[] = { 0, 1.57, 2.22, 2.71, 3.13,
341  3.50, 4.94, 6.05, 6.98, 7.81,
342  8.55, 9.23, 9.88, 10.47, 11.04,
343  11.58, 12.09, 12.58, 13.05, 13.52,
344  13.97, 14.43, 14.96, 15.67, 20.00 };
345  Double_t cs[] = { 0.5, 1.5, 2.5, 3.5, 4.5,
346  7.5, 12.5, 17.5, 22.5, 27.5,
347  32.5, 37.5, 42.5, 47.5, 52.5,
348  57.5, 62.5, 67.5, 72.5, 77.5,
349  82.5, 87.5, 92.5, 97.5 };
350  cents.Set(24,cs);
351  bins.Set(25,bs);
352  }
353  else if (sNN == 5023) {
354  // PbPb @ 5.02TeV only
355  // https://twiki.cern.ch/twiki/bin/viewauth/ALICE/CentralityCodeSnippets
356  Double_t bs[] = { 0.00, 1.56, 2.22, 2.71, 3.13,
357  3.51, 3.84, 4.15, 4.43, 4.71,
358  4.96, 6.08, 7.01, 7.84, 8.59,
359  9.27, 9.92, 10.5, 11.1, 11.6,
360  12.1, 12.6, 13.1, 13.6, 14.0,
361  14.5, 15.0, 15.7, 19.6 }; // 29
362  Double_t cs[] = { 0.5, 1.5, 2.5, 3.5, 4.5,
363  5.5, 6.5, 7.5, 8.5, 9.5,
364  12.5, 17.5, 22.5, 27.5, 32.5,
365  37.5, 42.5, 47.5, 52.5, 57.5,
366  62.5, 67.5, 72.5, 77.5, 82.5,
367  87.5, 92.5, 97.5 };
368  cents.Set(28,cs);
369  bins.Set(29,bs);
370  }
371  else if (sNN == 5440) {
372  // Xe-Xe @ 5.44TeV
373  // https://twiki.cern.ch/twiki/bin/viewauth/ALICE/XeXeCentStudies
374  Double_t bs[] = {
375  0.00, 2.12, 3.00, 3.66, 4.23, 5.18, 5.98, 6.68, 7.32, 7.91, 8.46,
376  8.97, 9.45, 9.91, 10.4, 10.8, 11.2, 11.6, 12.0, 12.4, 12.9, 18.3
377  };
378  Double_t cs[] = {
379  1.25, 3.75, 6.25, 8.75, 12.50, 17.50, 22.50, 27.50, 32.50, 37.50,
380  42.50, 47.50, 52.50, 57.50, 62.50, 67.50, 72.50, 77.50, 82.50,
381  87.50, 95.50
382  };
383  cents.Set(21, cs);
384  bins.Set(22, bs);
385  }
386  }
387  else if (tgtA || projA) { // p-Pb or Pb-p
388  if (sNN == 5023) {
389  Double_t cs[] = { 2.5, 7.5, 15., 30.,
390  50., 70., 90. };
391  Double_t bs[] = { 0, 1.83675, 2.59375, 3.66875,
392  5.18625, 6.35475, 7.40225, 13.8577 };
393  cents.Set(7, cs);
394  bins.Set(8, bs);
395  }
396  }
397  if (bins.GetSize() <= 0 || cents.GetSize() <= 0 ) {
398  // Nothing defined
399  Warning("MakeHistogram", "No bins defined for sNN=%d (%c-%c)",
400  sNN, tgtA ? 'A' : 'p', projA ? 'A' : 'p');
401  return 0;
402  }
403  printf("b bins: ");
404  for (Int_t i = 0; i < bins.GetSize(); i++) {
405  bins[i] *= fkFactor; // Scale to 1/1000 fm
406  printf("%s%7.1f", i!=0 ? "-" : "", bins[i]);
407  }
408  Printf("");
409  TH1* h = new TH1D(Form("raw%s",GetName()), "B to Centrality",
410  bins.GetSize()-1, bins.GetArray());
411  h->SetDirectory(0);
412  h->SetXTitle("b\\hbox{ [10^{3}fm]}");
413  h->SetYTitle("c\\hbox{ [\\%]}");
414  h->SetBinContent(0,1);
415  for (Int_t i = 1; i <= cents.GetSize(); i++) {
416  h->SetBinContent(i, cents[i-1]);
417  }
418  return h;
419  }
424  void PreEvent() { fCache = 100*fkFactor; };
431  {
432  // In Ap (EPOS) we have many spec in the target, meaning they will
433  // be detected on the C side (p is the projectile, A is the target)
434  // if (!fSpectators) return;
435  fCache = h.fB * fkFactor;
436  // Info("", "Cache=%ld", fCache);
437  // if (fBvsC) fBvsC->Fill(fCache, h.fC);
438  }
442  virtual void Process(const TParticle*) {};
447  void PostEvent() {};
453  virtual void Terminate(TCollection* out)
454  {
455  TH1* h = GetHistogram(out);
456  if (!h) {
457  Warning("Terminate", "No histogram on input");
458  out->ls();
459  return;
460  }
461  TH1* cent = static_cast<TH1*>(h->Clone(GetName()));
462  cent->SetDirectory(0);
463  cent->SetYTitle("Centrality [%]");
464  cent->SetTitle(Form("%s mapping", GetName()));
465  out->Add(cent);
466  // Scale to number of workers
467  Double_t scale = cent->GetBinContent(0);
468  cent->Scale(1./scale);
469  }
477  {
478  Double_t ret = (fHistogram ?
479  fHistogram->GetBinContent(fHistogram->FindBin(b*fkFactor)) :
480  200);
481  // Info("", "Look-up of %f (%f) -> %5.1f%%", b*fkFactor, b, ret);
482  return ret;
483  }
484  ClassDef(BCentEstimator,2);
485 };
486 
487 //____________________________________________________________________
492 {
498  FastNchCentEstimator(const char* name="")
499  : Fast1DCentEstimator(name)
500  {}
510  virtual void Process(const TParticle* p)
511  {
512  if (!IsCharged(p)) return;
513  if (!Accept(p)) return;
514  fCache++;
515  }
523  virtual Bool_t Accept(const TParticle* p) = 0;
524  virtual void Print(Option_t* option="nh") const
525  {
526  TString opt(option); opt.ToLower();
527  if (opt.Contains("n"))
528  Printf("1D Nch Estimator: %s (%screasing)",
529  GetName(),(fFromTop ? "de" : "in"));
530  opt.ReplaceAll("n", "");
532  }
533  ClassDef(FastNchCentEstimator,1);
534 };
535 
536 //____________________________________________________________________
541 {
555  V0CentEstimator(Short_t mode=0, Bool_t onlyPrimary=false)
556  : FastNchCentEstimator(Form("%s%s",
557  (mode < 0 ? "V0C" :
558  mode > 0 ? "V0A" : "V0M"),
559  (onlyPrimary ? "P" : ""))),
560  fMode(mode),
561  fOnlyPrimary(onlyPrimary)
562  {
563  fAEtaMin = +2.8;
564  fAEtaMax = +5.1;
565  fCEtaMin = -3.7;
566  fCEtaMax = -1.7;
567  }
568  void Flip(Bool_t onlySign=true)
569  {
570  if (!onlySign) {
571  std::swap(fAEtaMin, fCEtaMax); // a=[-1.7,*] c=[*,+2.8]
572  std::swap(fAEtaMax, fCEtaMin); // a=[*,-3.7] c=[+5.1,*]
573  }
574  else {
575  std::swap(fAEtaMin, fAEtaMax); // a=[+5.1,+2.8]
576  std::swap(fCEtaMin, fCEtaMax); // c=[-1.7,-3.7]
577  }
578  // onlySign !onlySign
579  fAEtaMin *= -1; // a=[-5.1,*] a=[+1.7,*]
580  fAEtaMax *= -1; // a=[*,-2.8] a=[*,+3.7]
581  fCEtaMin *= -1; // c=[+1.7,*] c=[-5.1,*]
582  fCEtaMax *= -1; // c=[*,+3.7] c=[*,-2.8]
583  Printf("Flipped %s", (onlySign ? "sign" : "acceptance"));
584  }
595  void Setup(TCollection* l, TTree* tree, UShort_t sNN,
596  Bool_t tgtA, Bool_t projA)
597  {
598  Bool_t isAA = (tgtA && projA);
599  Bool_t isPA = (tgtA ^ projA); // XOR
600  UInt_t max = (isAA ? 13000 : isPA ? 800 : 300);
601  UInt_t dBin = (isAA ? 10 : isPA ? 1 : 1);
602  Color_t color = (fMode < 0 ? kRed : fMode > 0 ? kBlue : kGreen)+2;
603  fHistogram = new TH1D(Form("raw%s",GetName()),
604  Form("%s #it{N}_{ch} distribution", GetName()),
605  max/dBin, 0, (fMode == 0 ? 2 : 1)*max);
606  fHistogram->SetXTitle("#it{N}_{ch}");
607  fHistogram->SetYTitle("Raw #it{P}(#it{N}_{ch})");
608  fHistogram->SetDirectory(0);
609  fHistogram->SetLineColor(color);
610  fHistogram->SetFillColor(color);
611  fHistogram->SetMarkerColor(color);
612  fHistogram->SetMarkerStyle(20);
613  fHistogram->SetFillStyle(3002);
614 
615  Fast1DCentEstimator::Setup(l, tree, sNN, tgtA, projA);
616  }
625  Bool_t Accept(const TParticle* p)
626  {
627  if (fOnlyPrimary && !IsPrimary(p)) return false;
628  Double_t eta = Eta(p);
629  Bool_t v0A = ((eta >= fAEtaMin) && (eta <= fAEtaMax));
630  Bool_t v0C = ((eta >= fCEtaMin) && (eta <= fCEtaMax));
631  if (fMode < 0) return v0C;
632  if (fMode > 0) return v0A;
633  return v0A || v0C;
634  }
641  {
642  return static_cast<TH1*>(l->FindObject(Form("raw%s",GetName())));
643  }
644  virtual void Print(Option_t* option="nah") const
645  {
646  TString opt(option); opt.ToLower();
647  if (opt.Contains("n"))
648  Printf("V0 Estimator: %s (%screasing) mode=%d",
649  GetName(),(fFromTop ? "de" : "in"), fMode);
650  if (opt.Contains("a"))
651  Printf(" A: eta=[%f,%f], C: eta=[%f,%f]",
652  fAEtaMin, fAEtaMax, fCEtaMin, fCEtaMax);
653  opt.ReplaceAll("n", "");
655  }
656  ClassDef(V0CentEstimator,1);
657 };
658 //____________________________________________________________________
663 {
672  : FastNchCentEstimator(Form("RefMult%02dd%02d",
673  Int_t(etaCut), Int_t(100*etaCut)%100)),
674  fEtaCut(etaCut)
675  {}
686  void Setup(TCollection* l, TTree* tree, UShort_t sNN,
687  Bool_t tgtA, Bool_t projA)
688  {
689  Bool_t isAA = (tgtA && projA);
690  Bool_t isPA = (tgtA ^ projA); // XOR
691  UInt_t max = (isAA ? 15000 : isPA ? 900 : 200);
692  UInt_t dBin = (isAA ? 10 : isPA ? 1 : 1);
693  Color_t color = kMagenta;
694  fHistogram = new TH1D(Form("raw%s",GetName()),
695  Form("#it{N}_{ch} |#it{#eta}|<%5.2f distribution",
696  fEtaCut), max/dBin, 0, max);
697  fHistogram->SetXTitle("#it{N}_{ch}");
698  fHistogram->SetYTitle("Raw #it{P}(#it{N}_{ch})");
699  fHistogram->SetDirectory(0);
700  fHistogram->SetLineColor(color);
701  fHistogram->SetFillColor(color);
702  fHistogram->SetMarkerColor(color);
703  fHistogram->SetMarkerStyle(20);
704  fHistogram->SetFillStyle(3002);
705 
706  Fast1DCentEstimator::Setup(l, tree, sNN, tgtA, projA);
707  }
716  Bool_t Accept(const TParticle* p)
717  {
718  Double_t eta = Eta(p);
719  if (!IsPrimary(p)) return false;
720  if (TMath::Abs(eta) > fEtaCut) return false;
721  return true;
722  }
729  {
730  return static_cast<TH1*>(l->FindObject(Form("raw%s",GetName())));
731  }
732  ClassDef(RefMultEstimator,1);
733 };
734 //____________________________________________________________________
739 {
740  // What to put in
749  ULong64_t fNSpec;
759  Bool_t neutrons=true,
760  Bool_t spectators=false,
761  Bool_t primary=false)
762  : Fast1DCentEstimator(Form("Z%c%c%c%c",
763  (neutrons ? 'N' : 'P'),
764  (mode < 0 ? 'C' : (mode > 0 ? 'A' : 'M')),
765  (spectators ? 'S' : 'E'),
766  (primary ? 'P' : 'A'))),
767  fMode(mode),
768  fNeutrons(neutrons),
769  fSpectators(spectators),
770  fPrimary(primary),
771  fMinEta(0),
772  fMaxEta(0),
773  fMaxPhi(-1)
774  {
775  fFromTop = (fSpectators ? false : true);
776  if (fNeutrons) {
777  const Double_t zN = 1161.3; // distance to ZN cm
778  const Double_t dN = 7; // Area of ZN
779  const Double_t rN = TMath::Sqrt(2*dN*dN); // Radius of ZN
780  const Double_t tN = TMath::ATan2(rN,zN); // Largest angle ZN
781  fMinEta = -TMath::Log(TMath::Tan(tN/2));
782  fMaxEta = TMath::Infinity();
783  fMaxPhi = -1;
784  }
785  else {
786  const Double_t zP = 1156.3; // distance to ZC cm
787  const Double_t wP = 20.8; // Width of cal
788  const Double_t hP = 12; // Width of cal
789  const Double_t oP = 19;
790  const Double_t rO = TMath::Sqrt(TMath::Power(oP+wP/2,2)+hP*hP/4);
791  const Double_t rI = TMath::Sqrt(TMath::Power(oP-wP/2,2)+hP*hP/4);
792  const Double_t tO = TMath::ATan2(rO, zP);
793  const Double_t tI = TMath::ATan2(rI, zP);
794  fMinEta = -TMath::Log(TMath::Tan(tO/2));
795  fMaxEta = -TMath::Log(TMath::Tan(tI/2));
796  fMaxPhi = TMath::ATan2(hP/2,oP);
797  }
798  }
809  void Setup(TCollection* l, TTree* tree, UShort_t sNN,
810  Bool_t tgtA, Bool_t projA)
811  {
812  Bool_t isAA = (tgtA && projA);
813  Bool_t isPA = (tgtA ^ projA); // XOR
814  UInt_t max = (isAA ? 2000 : isPA ? 300 : 30);
815  UInt_t dBin = (isAA ? 10 : isPA ? 1 : 1);
816  Color_t color = (fMode < 0 ? kRed : fMode > 0 ? kBlue : kGreen)+2;
817 
818  TString nTxt;
819  nTxt.Form("#it{N}_{%s%c} (%s)",
820  fSpectators ? "spec," : "", fNeutrons ? 'n' : 'p',
821  fPrimary ? "primary" : "all");
822  fHistogram = new TH1D(Form("raw%s",GetName()),
823  Form("%s %s distribution", GetName(), nTxt.Data()),
824  max/dBin, 0, (fMode == 0 ? 2 : 1)*max);
825  fHistogram->SetXTitle(nTxt);
826  fHistogram->SetYTitle(Form("Raw #it{P}(%s)", nTxt.Data()));
827  fHistogram->SetDirectory(0);
828  fHistogram->SetLineColor(color);
829  fHistogram->SetFillColor(color);
830  fHistogram->SetMarkerColor(color);
831  fHistogram->SetMarkerStyle(20);
832  fHistogram->SetFillStyle(3002);
833 
834  Fast1DCentEstimator::Setup(l, tree, sNN, tgtA, projA);
835  }
841  virtual void Process(const TParticle* p)
842  {
843  // Check if we should count emitted
844  if (fSpectators) return;
845 
846  // Check if we should count only primaries
847  if (fPrimary && !IsPrimary(p)) return;
848  else if (!fPrimary && p->GetStatusCode()==4) return;
849  // if ((fFlags & kPrimary) != 0 && !IsPrimary(p)) return;
850 
851  // Check if this is a spectator
852  // Int_t status = p->GetStatusCode();
853  // if ((fFlags & kSpectators) == 0 && (status==13 || status==14)) return;
854 
855  // Check particle type
856  Int_t aPdg = TMath::Abs(p->GetPdgCode());
857  if (fNeutrons && aPdg != kNeutron)
858  // Looking for neutrons
859  return;
860  if (!fNeutrons && aPdg != kProton)
861  // Looking for protons
862  return;
863 
864  // Check side
865  Double_t eta = Eta(p);
866  if (fMode < 0 && eta > 0) return; // Wrong side
867  if (fMode > 0 && eta < 0) return; // Wrong side
868 
869  // Check acceptance
870  Double_t aeta = TMath::Abs(eta);
871  if (aeta > fMaxEta || aeta < fMinEta) return; // Not acceptance
872 
873  if (fMaxPhi > 0) {
874  Double_t phi = TMath::Abs(Phi(p) - (eta < 0 ? 0 : TMath::Pi()));
875  if (phi > fMaxPhi) return; // Not acceptance
876  }
877 
878  // Decrement
879  fCache++;
880  // Printf("%s: increment from %d (%f) -> %lld",
881  // GetName(), aPdg, eta, fCache);
882  }
884  {
885  // In Ap (EPOS) we have many spec in the target, meaning they will
886  // be detected on the C side (p is the projectile, A is the target)
887  // if (!fSpectators) return;
888  fCache = 0;
889  fNSpec = 0;
890  Int_t nC = 0;
891  Int_t nA = 0;
892  if (fNeutrons) {
893  nC = h.fNSpecNtgt;
894  nA = h.fNSpecNproj;
895  }
896  else {
897  nC = h.fNSpecPtgt;
898  nA = h.fNSpecPproj;
899  }
900  if (fMode < 0) fNSpec = nC;
901  else if (fMode > 0) fNSpec = nA;
902  else fNSpec = nC + nA;
903  // Printf("%s: Initial cache value: %d", GetName(), fNSpec);
904  }
905  virtual void PostEvent()
906  {
907  // Info(GetName(), "Nspec=%lld Nneu=%lld", fNSpec, fCache);
908  if (fSpectators) fCache = fNSpec;
909  else {
910  if (fNSpec > fCache) {
911  // Warning("PostEvent", "Nspec (%d) > Nneu (%d)", fNSpec, fCache);
912  fNSpec = fCache;
913  }
914  if (!fPrimary) fCache -= fNSpec;
915  }
917  }
924  {
925  return static_cast<TH1*>(l->FindObject(Form("raw%s",GetName())));
926  }
927  virtual void Print(Option_t* option="nah") const
928  {
929  TString opt(option); opt.ToLower();
930  if (opt.Contains("n"))
931  Printf("ZN Estimator: %s (%screasing) mode=%d - %s %s",
932  GetName(),(fFromTop ? "de" : "in"), fMode,
933  (fSpectators ? "spectator" : "emitted"),
934  (fPrimary ? "primary" : "all"),
935  (fNeutrons ? "neutrons" : "protons"));
936  if (opt.Contains("a"))
937  Printf(" eta=[%f,%f], phi=%f", fMinEta, fMaxEta, fMaxPhi);
938  opt.ReplaceAll("n", "");
940  }
941  ClassDef(ZNCentEstimator,1);
942 };
943 #endif
944 //
945 // EOF
946 //
Int_t color[]
print message on plot with ok/not ok
double Double_t
Definition: External.C:58
virtual void Terminate(TCollection *out)
virtual void PostEvent()
virtual void Setup(TCollection *out, TTree *tree, UShort_t sNN, Bool_t tgtA, Bool_t projA)
static Bool_t IsPrimary(const TParticle *p)
virtual void PostEvent()
static Double_t Eta(const TParticle *p)
virtual void Terminate(TCollection *out)=0
virtual void Process(const TParticle *p)
const char * GetName() const
virtual void ProcessHeader(FastShortHeader &)
virtual void Print(Option_t *option="nah") const
const Double_t fkFactor
virtual void Print(Option_t *option="") const
TH1 * MakeHistogram(UShort_t sNN, Bool_t tgtA, Bool_t projA)
static Bool_t IsWeakDecay(const TParticle *p)
void ProcessHeader(FastShortHeader &h)
int Int_t
Definition: External.C:63
unsigned int UInt_t
Definition: External.C:33
void Flip(Bool_t onlySign=true)
static Bool_t IsCharged(const TParticle *p)
virtual void Print(Option_t *option="nh") const
Double_t GetCentrality(Double_t b) const
Definition: External.C:212
virtual void Print(Option_t *option="nah") const
Bool_t Accept(const TParticle *p)
virtual TH1 * GetHistogram(TCollection *l)
V0CentEstimator(Short_t mode=0, Bool_t onlyPrimary=false)
Int_t mode
Definition: anaM.C:41
virtual TH1 * GetHistogram(TCollection *l)
short Short_t
Definition: External.C:23
virtual void PostEvent()
void Setup(TCollection *l, TTree *tree, UShort_t sNN, Bool_t tgtA, Bool_t projA)
void SetVerbose(Bool_t verb)
FastCentEstimator(const char *name="")
static Double_t Phi(const TParticle *p)
ZNCentEstimator(Short_t mode=0, Bool_t neutrons=true, Bool_t spectators=false, Bool_t primary=false)
virtual void Process(const TParticle *p)=0
virtual TH1 * GetHistogram(TCollection *l)
virtual const char * MultSpec() const
RefMultEstimator(Double_t etaCut=0.8)
virtual void Terminate(TCollection *out)
static Double_t Theta(const TParticle *p)
void Setup(TCollection *l, TTree *tree, UShort_t sNN, Bool_t tgtA, Bool_t projA)
virtual ~FastCentEstimator()
unsigned short UShort_t
Definition: External.C:28
void Setup(TCollection *l, TTree *tree, UShort_t sNN, Bool_t tgtA, Bool_t projA)
const char Option_t
Definition: External.C:48
void ProcessHeader(FastShortHeader &h)
bool Bool_t
Definition: External.C:53
FastNchCentEstimator(const char *name="")
virtual void Setup(TCollection *l, TTree *tree, UShort_t, Bool_t, Bool_t)
virtual void PreEvent()
Bool_t Accept(const TParticle *p)
virtual void Process(const TParticle *p)
virtual void Print(Option_t *option="nh") const
void swap(PWGJE::EMCALJetTasks::AliAnalysisTaskEmcalJetHPerformance &first, PWGJE::EMCALJetTasks::AliAnalysisTaskEmcalJetHPerformance &second)
virtual void Setup(TCollection *out, TTree *tree, UShort_t sNN, Bool_t tgtA, Bool_t projA)=0
virtual void Process(const TParticle *)
Definition: External.C:196
virtual TH1 * GetHistogram(TCollection *l)
Fast1DCentEstimator(const char *name="")