1 #ifndef FASTCENTESTIMATOR_C 2 #define FASTCENTESTIMATOR_C 8 # include <TCollection.h> 10 # include <TParticle.h> 42 :
TObject(), fName(name), fVerbose(false)
53 const char*
GetName()
const {
return fName.Data(); }
80 virtual void Process(
const TParticle* p) = 0;
93 Printf(
"%s: %s", ClassName(),
GetName());
106 Double_t theta = TMath::ATan2(pT, pZ);
119 Double_t tanth = TMath::Tan(theta/2);
120 if (tanth < 1e-100)
return +1e10;
121 if (tanth > 1e100)
return -1e10;
136 Double_t phi = TMath::ATan2(px, py);
148 return p->TestBit(BIT(14));
159 return p->TestBit(BIT(15));
170 return p->TestBit(BIT(16));
204 virtual const char*
MultSpec()
const {
return "l"; }
215 if (fHistogram && l) {
216 Info(
"Setup",
"Adding histogram %s to output",
217 fHistogram->GetName());
220 TString leaves; leaves.Form(
"value/%s", MultSpec());
221 if (tree) tree->Branch(
GetName(), &fCache, leaves.Data());
235 if (
fVerbose) Info(
"PostEvent",
" Got %lld %s particles",
237 fHistogram->Fill(fCache);
251 TH1* h = GetHistogram(out);
253 cent->SetDirectory(0);
254 cent->SetYTitle(
"Centrality [%]");
255 cent->SetTitle(Form(
"%s mapping",
GetName()));
259 Int_t nX = h->GetNbinsX();
261 Int_t dBin = (fFromTop ? -1 : 1);
262 Int_t end = (fFromTop ? 0 : nX+1);
263 Int_t start = (fFromTop ? nX : 1);
265 Info(
"Teminate",
"Integrating %s from bin %d to 1",
267 for (
Int_t i = start; i != end; i += dBin) {
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);
275 Info(
"Terminate",
"Bin %3d -> %9f/%9f -> %5.1f%%",
276 i, curInt, total, curCent);
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());
317 fHistogram = MakeHistogram(sNN, tgtA, projA);
328 return static_cast<TH1*
>(l->FindObject(Form(
"raw%s",
GetName())));
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 };
353 else if (sNN == 5023) {
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 };
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,
371 else if (sNN == 5440) {
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
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,
387 else if (tgtA || projA) {
389 Double_t cs[] = { 2.5, 7.5, 15., 30.,
391 Double_t bs[] = { 0, 1.83675, 2.59375, 3.66875,
392 5.18625, 6.35475, 7.40225, 13.8577 };
397 if (bins.GetSize() <= 0 || cents.GetSize() <= 0 ) {
399 Warning(
"MakeHistogram",
"No bins defined for sNN=%d (%c-%c)",
400 sNN, tgtA ?
'A' :
'p', projA ?
'A' :
'p');
404 for (
Int_t i = 0; i < bins.GetSize(); i++) {
406 printf(
"%s%7.1f", i!=0 ?
"-" :
"", bins[i]);
410 bins.GetSize()-1, bins.GetArray());
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]);
435 fCache = h.
fB * fkFactor;
455 TH1* h = GetHistogram(out);
457 Warning(
"Terminate",
"No histogram on input");
462 cent->SetDirectory(0);
463 cent->SetYTitle(
"Centrality [%]");
464 cent->SetTitle(Form(
"%s mapping",
GetName()));
467 Double_t scale = cent->GetBinContent(0);
468 cent->Scale(1./scale);
479 fHistogram->GetBinContent(fHistogram->FindBin(b*fkFactor)) :
513 if (!Accept(p))
return;
523 virtual Bool_t Accept(
const TParticle* p) = 0;
527 if (opt.Contains(
"n"))
528 Printf(
"1D Nch Estimator: %s (%screasing)",
529 GetName(),(fFromTop ?
"de" :
"in"));
530 opt.ReplaceAll(
"n",
"");
558 mode > 0 ?
"V0A" :
"V0M"),
559 (onlyPrimary ?
"P" :
""))),
561 fOnlyPrimary(onlyPrimary)
583 Printf(
"Flipped %s", (onlySign ?
"sign" :
"acceptance"));
598 Bool_t isAA = (tgtA && projA);
599 Bool_t isPA = (tgtA ^ projA);
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;
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);
627 if (fOnlyPrimary && !
IsPrimary(p))
return false;
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;
642 return static_cast<TH1*
>(l->FindObject(Form(
"raw%s",
GetName())));
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",
"");
689 Bool_t isAA = (tgtA && projA);
690 Bool_t isPA = (tgtA ^ projA);
691 UInt_t max = (isAA ? 15000 : isPA ? 900 : 200);
692 UInt_t dBin = (isAA ? 10 : isPA ? 1 : 1);
693 Color_t
color = kMagenta;
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);
720 if (TMath::Abs(eta) > fEtaCut)
return false;
730 return static_cast<TH1*
>(l->FindObject(Form(
"raw%s",
GetName())));
763 (neutrons ?
'N' :
'P'),
764 (
mode < 0 ?
'C' : (
mode > 0 ?
'A' :
'M')),
765 (spectators ?
'S' :
'E'),
766 (primary ?
'P' :
'A'))),
769 fSpectators(spectators),
775 fFromTop = (fSpectators ?
false :
true);
779 const Double_t rN = TMath::Sqrt(2*dN*dN);
780 const Double_t tN = TMath::ATan2(rN,zN);
781 fMinEta = -TMath::Log(TMath::Tan(tN/2));
782 fMaxEta = TMath::Infinity();
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);
812 Bool_t isAA = (tgtA && projA);
813 Bool_t isPA = (tgtA ^ projA);
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;
819 nTxt.Form(
"#it{N}_{%s%c} (%s)",
820 fSpectators ?
"spec," :
"", fNeutrons ?
'n' :
'p',
821 fPrimary ?
"primary" :
"all");
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);
844 if (fSpectators)
return;
848 else if (!fPrimary && p->GetStatusCode()==4)
return;
856 Int_t aPdg = TMath::Abs(p->GetPdgCode());
860 if (!fNeutrons && aPdg !=
kProton)
866 if (fMode < 0 && eta > 0)
return;
867 if (fMode > 0 && eta < 0)
return;
871 if (aeta > fMaxEta || aeta < fMinEta)
return;
874 Double_t phi = TMath::Abs(
Phi(p) - (eta < 0 ? 0 : TMath::Pi()));
875 if (phi > fMaxPhi)
return;
900 if (fMode < 0) fNSpec = nC;
901 else if (fMode > 0) fNSpec = nA;
902 else fNSpec = nC + nA;
908 if (fSpectators) fCache = fNSpec;
910 if (fNSpec > fCache) {
914 if (!fPrimary) fCache -= fNSpec;
925 return static_cast<TH1*
>(l->FindObject(Form(
"raw%s",
GetName())));
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",
"");
Int_t color[]
print message on plot with ok/not ok
virtual void Terminate(TCollection *out)
virtual void Setup(TCollection *out, TTree *tree, UShort_t sNN, Bool_t tgtA, Bool_t projA)
static Bool_t IsPrimary(const TParticle *p)
virtual ~Fast1DCentEstimator()
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
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)
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
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)
virtual TH1 * GetHistogram(TCollection *l)
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 ~FastNchCentEstimator()
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()
void Setup(TCollection *l, TTree *tree, UShort_t sNN, Bool_t tgtA, Bool_t projA)
void ProcessHeader(FastShortHeader &h)
FastNchCentEstimator(const char *name="")
virtual void Setup(TCollection *l, TTree *tree, UShort_t, Bool_t, Bool_t)
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 *)
virtual TH1 * GetHistogram(TCollection *l)
Fast1DCentEstimator(const char *name="")