AliPhysics  dcc03a0 (dcc03a0)
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
FastCentEstimators.C
Go to the documentation of this file.
1 #ifndef __CINT__
2 # include <TObject.h>
3 # include <TString.h>
4 # include <TH1.h>
5 # include <TMath.h>
6 # include <TCollection.h>
7 # include <TTree.h>
8 # include <TParticle.h>
9 # include <TError.h>
10 # include "FastShortHeader.C"
11 # include <TH2.h>
12 #else
13 class TH1;
14 // class TH2;
15 class TCollection;
16 class TTree;
17 class TParticle;
18 class FastShortHeader;
19 #endif
20 
21 enum {
22  kNeutron = 2112,
23  kProton = 2212
24 };
25 
26 //====================================================================
30 struct FastCentEstimator : public TObject
31 {
39  FastCentEstimator(const char* name="")
40  : TObject(), fName(name), fVerbose(false)
41  {}
45  virtual ~FastCentEstimator() {}
51  const char* GetName() const { return fName.Data(); }
52  void SetVerbose(Bool_t verb) { fVerbose = verb; }
63  virtual void Setup(TCollection* out, TTree* tree, UShort_t sNN,
64  Bool_t tgtA, Bool_t projA) = 0;
69  virtual void PreEvent() {}
70  virtual void ProcessHeader(FastShortHeader&) {}
78  virtual void Process(const TParticle* p) = 0;
82  virtual void PostEvent() {}
88  virtual void Terminate(TCollection* out) = 0;
89  virtual void Print(Option_t* option="") const
90  {
91  Printf("%s: %s", ClassName(), GetName());
92  }
100  static Double_t Theta(const TParticle* p)
101  {
102  Double_t pT = p->Pt();
103  Double_t pZ = p->Pz();
104  Double_t theta = TMath::ATan2(pT, pZ);
105  return theta;
106  }
114  static Double_t Eta(const TParticle* p)
115  {
116  Double_t theta = Theta(p);
117  Double_t tanth = TMath::Tan(theta/2);
118  if (tanth < 1e-100) return +1e10;
119  if (tanth > 1e100) return -1e10;
120  Double_t eta = -TMath::Log(tanth);
121  return eta;
122  }
130  static Double_t Phi(const TParticle* p)
131  {
132  Double_t px = p->Px();
133  Double_t py = p->Py();
134  Double_t phi = TMath::ATan2(px, py);
135  return phi;
136  }
144  static Bool_t IsPrimary(const TParticle* p)
145  {
146  return p->TestBit(BIT(14));
147  }
155  static Bool_t IsWeakDecay(const TParticle* p)
156  {
157  return p->TestBit(BIT(15));
158  }
166  static Bool_t IsCharged(const TParticle* p)
167  {
168  return p->TestBit(BIT(16));
169  }
171 };
172 
173 
174 
175 
176 //____________________________________________________________________
181 {
183  ULong64_t fCache;
195  Fast1DCentEstimator(const char* name="")
196  : FastCentEstimator(name), fHistogram(0), fFromTop(true)
197  {}
201  virtual ~Fast1DCentEstimator() {}
202  virtual const char* MultSpec() const { return "l"; }
210  void Setup(TCollection* l, TTree* tree, UShort_t,
211  Bool_t, Bool_t)
212  {
213  if (fHistogram && l) l->Add(fHistogram);
214  TString leaves; leaves.Form("value/%s", MultSpec());
215  if (tree) tree->Branch(GetName(), &fCache, leaves.Data());
216  }
220  virtual void PreEvent()
221  {
222  fCache = 0;
223  }
227  virtual void PostEvent()
228  {
229  if (fVerbose) Info("PostEvent", " Got %lld %s particles",
230  fCache, GetName());
231  fHistogram->Fill(fCache);
232  }
233  virtual TH1* GetHistogram(TCollection* l) = 0;
243  virtual void Terminate(TCollection* out)
244  {
245  TH1* h = GetHistogram(out);
246  TH1* cent = static_cast<TH1*>(h->Clone(GetName()));
247  cent->SetDirectory(0);
248  cent->SetYTitle("Centrality [%]");
249  cent->SetTitle(Form("%s mapping", GetName()));
250  cent->Reset();
251  out->Add(cent);
252 
253  Int_t nX = h->GetNbinsX();
254  Double_t total = h->Integral(1,nX);
255  Int_t dBin = (fFromTop ? -1 : 1);
256  Int_t end = (fFromTop ? 0 : nX+1);
257  Int_t start = (fFromTop ? nX : 1);
258  if (fVerbose)
259  Info("Teminate", "Integrating %s from bin %d to 1",
260  h->GetName(), nX);
261  for (Int_t i = start; i != end; i += dBin) {
262  Double_t curInt = (fFromTop ?
263  h->Integral(i, start) :
264  h->Integral(start,i));
265  if (curInt < 0) continue;
266  Double_t curCent = curInt / total * 100;
267  cent->SetBinContent(i, curCent);
268  if (fVerbose)
269  Info("Terminate", "Bin %3d -> %9f/%9f -> %5.1f%%",
270  i, curInt, total, curCent);
271  }
272  }
273  virtual void Print(Option_t* option="nh") const
274  {
275  TString opt(option); opt.ToLower();
276  if (opt.Contains("n"))
277  Printf("1D Estimator: %s (%screasing)",
278  GetName(),(fFromTop ? "de" : "in"));
279  if (opt.Contains("h") && fHistogram) {
280  Int_t nBin = fHistogram->GetNbinsX();
281  Printf(" %d bins between %f and %f",nBin,
282  fHistogram->GetXaxis()->GetXmin(),
283  fHistogram->GetXaxis()->GetXmax());
284  }
285  }
287 };
288 
289 //____________________________________________________________________
291 {
292  // TH2* fBvsC;
295  : Fast1DCentEstimator("B"),
296  // fBvsC(0),
297  fkFactor(1000)
298  {}
308  virtual void Setup(TCollection* out, TTree* tree, UShort_t sNN,
309  Bool_t tgtA, Bool_t projA)
310  {
311  fHistogram = MakeHistogram(sNN, tgtA, projA);
312 #if 0
313  if (fHistogram &&
314  fHistogram->GetXaxis()) {
315  if (fHistogram->GetXaxis()->GetXbins() &&
316  fHistogram->GetXaxis()->GetXbins()->GetArray()) {
317  fBvsC = new TH2D("bVsC", "Impact parameter vs Centrality",
318  fHistogram->GetXaxis()->GetNbins(),
319  fHistogram->GetXaxis()->GetXbins()->GetArray(),
320  20, 0, 100);
321  }
322  else {
323  fBvsC = new TH2D("bVsC", "Impact parameter vs Centrality",
324  fHistogram->GetXaxis()->GetNbins(),
325  fHistogram->GetXaxis()->GetXmin(),
326  fHistogram->GetXaxis()->GetXmax(),
327  20, 0, 100);
328  }
329  }
330  if (fBvsC) {
331  fBvsC->SetDirectory(0);
332  fBvsC->SetXTitle("b [fm]");
333  fBvsC->SetYTitle("Centrality [%]");
334  out->Add(fBvsC);
335  }
336  else
337  Warning("BCentEstimator", "Couldn't make bVsC histogram");
338 #endif
339  Fast1DCentEstimator::Setup(out, tree, sNN, tgtA, projA);
340  // if (tree) tree->Branch(GetName(), &fB, "value/D");
341  }
348  {
349  return static_cast<TH1*>(l->FindObject(Form("raw%s",GetName())));
350  }
352  {
353  TArrayD cents;
354  TArrayD bins; // In B
355  if (tgtA && projA) { // Pb-Pb
356  if (sNN == 2760) {
357  // PbPb @ 2.76TeV only
358  // Updated 4th of November 2014 from
359  // cern.ch/twiki/bin/view/ALICE/CentStudies
360  // #Tables_with_centrality_bins_AN1
361  Double_t bs[] = { 0, 1.57, 2.22, 2.71, 3.13,
362  3.50, 4.94, 6.05, 6.98, 7.81,
363  8.55, 9.23, 9.88, 10.47, 11.04,
364  11.58, 12.09, 12.58, 13.05, 13.52,
365  13.97, 14.43, 14.96, 15.67, 20.00 };
366  Double_t cs[] = { 0.5, 1.5, 2.5, 3.5, 4.5,
367  7.5, 12.5, 17.5, 22.5, 27.5,
368  32.5, 37.5, 42.5, 47.5, 52.5,
369  57.5, 62.5, 67.5, 72.5, 77.5,
370  82.5, 87.5, 92.5, 97.5 };
371  cents.Set(24,cs);
372  bins.Set(25,bs);
373  }
374  else if (sNN == 5023) {
375  // PbPb @ 5.02TeV only
376  // https://twiki.cern.ch/twiki/bin/viewauth/ALICE/CentralityCodeSnippets
377  Double_t bs[] = { 0.00, 1.56, 2.22, 2.71, 3.13,
378  3.51, 3.84, 4.15, 4.43, 4.71,
379  4.96, 6.08, 7.01, 7.84, 8.59,
380  9.27, 9.92, 10.5, 11.1, 11.6,
381  12.1, 12.6, 13.1, 13.6, 14.0,
382  14.5, 15.0, 15.7, 19.6 }; // 29
383  Double_t cs[] = { 0.5, 1.5, 2.5, 3.5, 4.5,
384  5.5, 6.5, 7.5, 8.5, 9.5,
385  12.5, 17.5, 22.5, 27.5, 32.5,
386  37.5, 42.5, 47.5, 52.5, 57.5,
387  62.5, 67.5, 72.5, 77.5, 82.5,
388  87.5, 92.5, 97.5 };
389  cents.Set(28,cs);
390  bins.Set(29,bs);
391  }
392  }
393  else if (tgtA || projA) { // p-Pb or Pb-p
394  if (sNN == 5023) {
395  Double_t cs[] = { 2.5, 7.5, 15., 30.,
396  50., 70., 90. };
397  Double_t bs[] = { 0, 1.83675, 2.59375, 3.66875,
398  5.18625, 6.35475, 7.40225, 13.8577 };
399  cents.Set(7, cs);
400  bins.Set(8, bs);
401  }
402  }
403  if (bins.GetSize() <= 0 || cents.GetSize() <= 0 ) // Nothing defined
404  return 0;
405  printf("b bins: ");
406  for (Int_t i = 0; i < bins.GetSize(); i++) {
407  bins[i] *= fkFactor; // Scale to 1/1000 fm
408  printf("%s%7.1f", i!=0 ? "-" : "", bins[i]);
409  }
410  Printf("");
411  TH1* h = new TH1D(Form("raw%s",GetName()), "B to Centrality",
412  bins.GetSize()-1, bins.GetArray());
413  h->SetDirectory(0);
414  h->SetXTitle("b\\hbox{ [10^{3}fm]}");
415  h->SetYTitle("c\\hbox{ [\\%]}");
416  h->SetBinContent(0,1);
417  for (Int_t i = 1; i <= cents.GetSize(); i++) {
418  h->SetBinContent(i, cents[i-1]);
419  }
420  return h;
421  }
426  void PreEvent() { fCache = 100*fkFactor; };
433  {
434  // In Ap (EPOS) we have many spec in the target, meaning they will
435  // be detected on the C side (p is the projectile, A is the target)
436  // if (!fSpectators) return;
437  fCache = h.fB * fkFactor;
438  Info("", "Cache=%ld", fCache);
439  // if (fBvsC) fBvsC->Fill(fCache, h.fC);
440  }
444  virtual void Process(const TParticle*) {};
449  void PostEvent() {};
455  virtual void Terminate(TCollection* out)
456  {
457  TH1* h = GetHistogram(out);
458  TH1* cent = static_cast<TH1*>(h->Clone(GetName()));
459  cent->SetDirectory(0);
460  cent->SetYTitle("Centrality [%]");
461  cent->SetTitle(Form("%s mapping", GetName()));
462  out->Add(cent);
463  // Scale to number of workers
464  Double_t scale = cent->GetBinContent(0);
465  cent->Scale(1./scale);
466  }
474  {
475  Double_t ret = (fHistogram ?
476  fHistogram->GetBinContent(fHistogram->FindBin(b*fkFactor)) :
477  200);
478  // Info("", "Look-up of %f (%f) -> %5.1f%%", b*fkFactor, b, ret);
479  return ret;
480  }
482 };
483 
484 //____________________________________________________________________
489 {
495  FastNchCentEstimator(const char* name="")
496  : Fast1DCentEstimator(name)
497  {}
507  virtual void Process(const TParticle* p)
508  {
509  if (!IsCharged(p)) return;
510  if (!Accept(p)) return;
511  fCache++;
512  }
520  virtual Bool_t Accept(const TParticle* p) = 0;
521  virtual void Print(Option_t* option="nh") const
522  {
523  TString opt(option); opt.ToLower();
524  if (opt.Contains("n"))
525  Printf("1D Nch Estimator: %s (%screasing)",
526  GetName(),(fFromTop ? "de" : "in"));
527  opt.ReplaceAll("n", "");
529  }
531 };
532 
533 //____________________________________________________________________
538 {
552  V0CentEstimator(Short_t mode=0, Bool_t onlyPrimary=false)
553  : FastNchCentEstimator(Form("%s%s",
554  (mode < 0 ? "V0C" :
555  mode > 0 ? "V0A" : "V0M"),
556  (onlyPrimary ? "P" : ""))),
557  fMode(mode),
558  fOnlyPrimary(onlyPrimary)
559  {
560  fAEtaMin = +2.8;
561  fAEtaMax = +5.1;
562  fCEtaMin = -3.7;
563  fCEtaMax = -1.7;
564  }
565  void Flip(Bool_t onlySign=true)
566  {
567  if (!onlySign) {
568  std::swap(fAEtaMin, fCEtaMax); // a=[-1.7,*] c=[*,+2.8]
569  std::swap(fAEtaMax, fCEtaMin); // a=[*,-3.7] c=[+5.1,*]
570  }
571  else {
572  std::swap(fAEtaMin, fAEtaMax); // a=[+5.1,+2.8]
573  std::swap(fCEtaMin, fCEtaMax); // c=[-1.7,-3.7]
574  }
575  // onlySign !onlySign
576  fAEtaMin *= -1; // a=[-5.1,*] a=[+1.7,*]
577  fAEtaMax *= -1; // a=[*,-2.8] a=[*,+3.7]
578  fCEtaMin *= -1; // c=[+1.7,*] c=[-5.1,*]
579  fCEtaMax *= -1; // c=[*,+3.7] c=[*,-2.8]
580  Printf("Flipped %s", (onlySign ? "sign" : "acceptance"));
581  }
592  void Setup(TCollection* l, TTree* tree, UShort_t sNN,
593  Bool_t tgtA, Bool_t projA)
594  {
595  Bool_t isAA = (tgtA && projA);
596  Bool_t isPA = (tgtA ^ projA); // XOR
597  UInt_t max = (isAA ? 13000 : isPA ? 800 : 300);
598  UInt_t dBin = (isAA ? 10 : isPA ? 1 : 1);
599  Color_t color = (fMode < 0 ? kRed : fMode > 0 ? kBlue : kGreen)+2;
600  fHistogram = new TH1D(Form("raw%s",GetName()),
601  Form("%s #it{N}_{ch} distribution", GetName()),
602  max/dBin, 0, (fMode == 0 ? 2 : 1)*max);
603  fHistogram->SetXTitle("#it{N}_{ch}");
604  fHistogram->SetYTitle("Raw #it{P}(#it{N}_{ch})");
605  fHistogram->SetDirectory(0);
606  fHistogram->SetLineColor(color);
607  fHistogram->SetFillColor(color);
608  fHistogram->SetMarkerColor(color);
609  fHistogram->SetMarkerStyle(20);
610  fHistogram->SetFillStyle(3002);
611 
612  Fast1DCentEstimator::Setup(l, tree, sNN, tgtA, projA);
613  }
622  Bool_t Accept(const TParticle* p)
623  {
624  if (fOnlyPrimary && !IsPrimary(p)) return false;
625  Double_t eta = Eta(p);
626  Bool_t v0A = ((eta >= fAEtaMin) && (eta <= fAEtaMax));
627  Bool_t v0C = ((eta >= fCEtaMin) && (eta <= fCEtaMax));
628  if (fMode < 0) return v0C;
629  if (fMode > 0) return v0A;
630  return v0A || v0C;
631  }
638  {
639  return static_cast<TH1*>(l->FindObject(Form("raw%s",GetName())));
640  }
641  virtual void Print(Option_t* option="nah") const
642  {
643  TString opt(option); opt.ToLower();
644  if (opt.Contains("n"))
645  Printf("V0 Estimator: %s (%screasing) mode=%d",
646  GetName(),(fFromTop ? "de" : "in"), fMode);
647  if (opt.Contains("a"))
648  Printf(" A: eta=[%f,%f], C: eta=[%f,%f]",
650  opt.ReplaceAll("n", "");
652  }
654 };
655 //____________________________________________________________________
660 {
669  : FastNchCentEstimator(Form("RefMult%02dd%02d",
670  Int_t(etaCut), Int_t(100*etaCut)%100)),
671  fEtaCut(etaCut)
672  {}
683  void Setup(TCollection* l, TTree* tree, UShort_t sNN,
684  Bool_t tgtA, Bool_t projA)
685  {
686  Bool_t isAA = (tgtA && projA);
687  Bool_t isPA = (tgtA ^ projA); // XOR
688  UInt_t max = (isAA ? 15000 : isPA ? 900 : 200);
689  UInt_t dBin = (isAA ? 10 : isPA ? 1 : 1);
690  Color_t color = kMagenta;
691  fHistogram = new TH1D(Form("raw%s",GetName()),
692  Form("#it{N}_{ch} |#it{#eta}|<%5.2f distribution",
693  fEtaCut), max/dBin, 0, max);
694  fHistogram->SetXTitle("#it{N}_{ch}");
695  fHistogram->SetYTitle("Raw #it{P}(#it{N}_{ch})");
696  fHistogram->SetDirectory(0);
697  fHistogram->SetLineColor(color);
698  fHistogram->SetFillColor(color);
699  fHistogram->SetMarkerColor(color);
700  fHistogram->SetMarkerStyle(20);
701  fHistogram->SetFillStyle(3002);
702 
703  Fast1DCentEstimator::Setup(l, tree, sNN, tgtA, projA);
704  }
713  Bool_t Accept(const TParticle* p)
714  {
715  Double_t eta = Eta(p);
716  if (!IsPrimary(p)) return false;
717  if (TMath::Abs(eta) > fEtaCut) return false;
718  return true;
719  }
726  {
727  return static_cast<TH1*>(l->FindObject(Form("raw%s",GetName())));
728  }
730 };
731 //____________________________________________________________________
736 {
737  // What to put in
746  ULong64_t fNSpec;
756  Bool_t neutrons=true,
757  Bool_t spectators=false,
758  Bool_t primary=false)
759  : Fast1DCentEstimator(Form("Z%c%c%c%c",
760  (neutrons ? 'N' : 'P'),
761  (mode < 0 ? 'C' : (mode > 0 ? 'A' : 'M')),
762  (spectators ? 'S' : 'E'),
763  (primary ? 'P' : 'A'))),
764  fMode(mode),
765  fNeutrons(neutrons),
766  fSpectators(spectators),
767  fPrimary(primary),
768  fMinEta(0),
769  fMaxEta(0),
770  fMaxPhi(-1)
771  {
772  fFromTop = (fSpectators ? false : true);
773  if (fNeutrons) {
774  const Double_t zN = 1161.3; // distance to ZN cm
775  const Double_t dN = 7; // Area of ZN
776  const Double_t rN = TMath::Sqrt(2*dN*dN); // Radius of ZN
777  const Double_t tN = TMath::ATan2(rN,zN); // Largest angle ZN
778  fMinEta = -TMath::Log(TMath::Tan(tN/2));
779  fMaxEta = TMath::Infinity();
780  fMaxPhi = -1;
781  }
782  else {
783  const Double_t zP = 1156.3; // distance to ZC cm
784  const Double_t wP = 20.8; // Width of cal
785  const Double_t hP = 12; // Width of cal
786  const Double_t oP = 19;
787  const Double_t rO = TMath::Sqrt(TMath::Power(oP+wP/2,2)+hP*hP/4);
788  const Double_t rI = TMath::Sqrt(TMath::Power(oP-wP/2,2)+hP*hP/4);
789  const Double_t tO = TMath::ATan2(rO, zP);
790  const Double_t tI = TMath::ATan2(rI, zP);
791  fMinEta = -TMath::Log(TMath::Tan(tO/2));
792  fMaxEta = -TMath::Log(TMath::Tan(tI/2));
793  fMaxPhi = TMath::ATan2(hP/2,oP);
794  }
795  }
806  void Setup(TCollection* l, TTree* tree, UShort_t sNN,
807  Bool_t tgtA, Bool_t projA)
808  {
809  Bool_t isAA = (tgtA && projA);
810  Bool_t isPA = (tgtA ^ projA); // XOR
811  UInt_t max = (isAA ? 2000 : isPA ? 300 : 30);
812  UInt_t dBin = (isAA ? 10 : isPA ? 1 : 1);
813  Color_t color = (fMode < 0 ? kRed : fMode > 0 ? kBlue : kGreen)+2;
814 
815  TString nTxt;
816  nTxt.Form("#it{N}_{%s%c} (%s)",
817  fSpectators ? "spec," : "", fNeutrons ? 'n' : 'p',
818  fPrimary ? "primary" : "all");
819  fHistogram = new TH1D(Form("raw%s",GetName()),
820  Form("%s %s distribution", GetName(), nTxt.Data()),
821  max/dBin, 0, (fMode == 0 ? 2 : 1)*max);
822  fHistogram->SetXTitle(nTxt);
823  fHistogram->SetYTitle(Form("Raw #it{P}(%s)", nTxt.Data()));
824  fHistogram->SetDirectory(0);
825  fHistogram->SetLineColor(color);
826  fHistogram->SetFillColor(color);
827  fHistogram->SetMarkerColor(color);
828  fHistogram->SetMarkerStyle(20);
829  fHistogram->SetFillStyle(3002);
830 
831  Fast1DCentEstimator::Setup(l, tree, sNN, tgtA, projA);
832  }
838  virtual void Process(const TParticle* p)
839  {
840  // Check if we should count emitted
841  if (fSpectators) return;
842 
843  // Check if we should count only primaries
844  if (fPrimary && !IsPrimary(p)) return;
845  else if (!fPrimary && p->GetStatusCode()==4) return;
846  // if ((fFlags & kPrimary) != 0 && !IsPrimary(p)) return;
847 
848  // Check if this is a spectator
849  // Int_t status = p->GetStatusCode();
850  // if ((fFlags & kSpectators) == 0 && (status==13 || status==14)) return;
851 
852  // Check particle type
853  Int_t aPdg = TMath::Abs(p->GetPdgCode());
854  if (fNeutrons && aPdg != kNeutron)
855  // Looking for neutrons
856  return;
857  if (!fNeutrons && aPdg != kProton)
858  // Looking for protons
859  return;
860 
861  // Check side
862  Double_t eta = Eta(p);
863  if (fMode < 0 && eta > 0) return; // Wrong side
864  if (fMode > 0 && eta < 0) return; // Wrong side
865 
866  // Check acceptance
867  Double_t aeta = TMath::Abs(eta);
868  if (aeta > fMaxEta || aeta < fMinEta) return; // Not acceptance
869 
870  if (fMaxPhi > 0) {
871  Double_t phi = TMath::Abs(Phi(p) - (eta < 0 ? 0 : TMath::Pi()));
872  if (phi > fMaxPhi) return; // Not acceptance
873  }
874 
875  // Decrement
876  fCache++;
877  // Printf("%s: increment from %d (%f) -> %lld",
878  // GetName(), aPdg, eta, fCache);
879  }
881  {
882  // In Ap (EPOS) we have many spec in the target, meaning they will
883  // be detected on the C side (p is the projectile, A is the target)
884  // if (!fSpectators) return;
885  fCache = 0;
886  fNSpec = 0;
887  Int_t nC = 0;
888  Int_t nA = 0;
889  if (fNeutrons) {
890  nC = h.fNSpecNtgt;
891  nA = h.fNSpecNproj;
892  }
893  else {
894  nC = h.fNSpecPtgt;
895  nA = h.fNSpecPproj;
896  }
897  if (fMode < 0) fNSpec = nC;
898  else if (fMode > 0) fNSpec = nA;
899  else fNSpec = nC + nA;
900  // Printf("%s: Initial cache value: %d", GetName(), fNSpec);
901  }
902  virtual void PostEvent()
903  {
904  // Info(GetName(), "Nspec=%lld Nneu=%lld", fNSpec, fCache);
905  if (fSpectators) fCache = fNSpec;
906  else {
907  if (fNSpec > fCache) {
908  // Warning("PostEvent", "Nspec (%d) > Nneu (%d)", fNSpec, fCache);
909  fNSpec = fCache;
910  }
911  if (!fPrimary) fCache -= fNSpec;
912  }
914  }
921  {
922  return static_cast<TH1*>(l->FindObject(Form("raw%s",GetName())));
923  }
924  virtual void Print(Option_t* option="nah") const
925  {
926  TString opt(option); opt.ToLower();
927  if (opt.Contains("n"))
928  Printf("ZN Estimator: %s (%screasing) mode=%d - %s %s",
929  GetName(),(fFromTop ? "de" : "in"), fMode,
930  (fSpectators ? "spectator" : "emitted"),
931  (fPrimary ? "primary" : "all"),
932  (fNeutrons ? "neutrons" : "protons"));
933  if (opt.Contains("a"))
934  Printf(" eta=[%f,%f], phi=%f", fMinEta, fMaxEta, fMaxPhi);
935  opt.ReplaceAll("n", "");
937  }
939 };
940 
941 //
942 // EOF
943 //
Int_t color[]
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 TH1 * GetHistogram(TCollection *l)=0
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
ClassDef(V0CentEstimator, 1)
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:228
Definition: External.C:212
virtual void Print(Option_t *option="nah") const
Bool_t Accept(const TParticle *p)
virtual TH1 * GetHistogram(TCollection *l)
void Setup(TCollection *l, TTree *tree, UShort_t, Bool_t, Bool_t)
V0CentEstimator(Short_t mode=0, Bool_t onlyPrimary=false)
Int_t mode
Definition: anaM.C:40
virtual Bool_t Accept(const TParticle *p)=0
virtual TH1 * GetHistogram(TCollection *l)
ClassDef(FastNchCentEstimator, 1)
short Short_t
Definition: External.C:23
ClassDef(RefMultEstimator, 1)
virtual void PostEvent()
void Setup(TCollection *l, TTree *tree, UShort_t sNN, Bool_t tgtA, Bool_t projA)
void SetVerbose(Bool_t verb)
ClassDef(FastCentEstimator, 1)
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)
ClassDef(Fast1DCentEstimator, 1)
ClassDef(BCentEstimator, 2)
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)
ClassDef(ZNCentEstimator, 1)
bool Bool_t
Definition: External.C:53
FastNchCentEstimator(const char *name="")
virtual void PreEvent()
Bool_t Accept(const TParticle *p)
virtual void Process(const TParticle *p)
virtual void Print(Option_t *option="nh") const
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="")