AliPhysics  ad6828d (ad6828d)
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
dNdetaAnalysis.C
Go to the documentation of this file.
1 #ifndef __CINT__
2 # include <TH1.h>
3 # include <TH2.h>
4 # include <TMath.h>
5 # include <TParticle.h>
6 # include <TObjArray.h>
7 # include <TString.h>
8 #else
9 class TH1;
10 class TParticle;
11 #endif
12 #include "FastAnalysis.C"
13 #include "FastCentHelper.C"
14 namespace dNdeta {
15  //====================================================================
19  struct Base : public FastAnalysis
20  {
22 
23  Base(Bool_t verbose=false, Int_t monitor=0)
24  : FastAnalysis(verbose,monitor), fdNdeta(0)
25  {}
31  static TH1D* CreatedNdeta()
32  {
33  Double_t maxEta = 9; // 10;
34  Int_t nEta = Int_t(200 * maxEta/5+.5);
35  TH1D* eta = new TH1D("dNdeta",
36  "Charged particle pseudo-rapidity density",
37  nEta, -maxEta, +maxEta);
38  eta->Sumw2();
39  eta->SetXTitle("#it{#eta}");
40  eta->SetYTitle("1/N d#it{N}_{ch}/d#it{#eta}");
41  eta->SetMarkerColor(kRed+2);
42  eta->SetMarkerStyle(20);
43  eta->SetDirectory(0);
44 
45  return eta;
46  }
50  virtual void SlaveBegin(TTree*)
51  {
52  Info("SlaveBegin", "Making dN/deta histogram");
54  fdNdeta->SetMarkerColor(kBlack);
55  fdNdeta->SetMarkerStyle(21);
56  fdNdeta->SetTitle(GetName());
57  fOutput->Add(fdNdeta);
58  }
64  virtual Bool_t ProcessParticle(const TParticle* p)
65  {
66  Double_t pT = p->Pt();
67  Double_t pZ = p->Pz();
68  Double_t theta = TMath::ATan2(pT, pZ);
69  Double_t eta = (pT < 1e-10 ? 1024 : -TMath::Log(TMath::Tan(theta/2)));
70  if (TMath::Abs(eta) > 1000) return false;
71 
72  Fill(eta);
73  return true;
74  }
75  virtual void Fill(Double_t eta) { fdNdeta->Fill(eta); }
80  virtual void Terminate()
81  {
82  fOK = GetEventCount();
83  if (fOK <= 0) {
84  SetStatus(-1);
85  Warning("Terminate", "No events selected");
86  return;
87  }
88  fdNdeta = static_cast<TH1*>(GetOutputObject("dNdeta", TH1::Class()));
89  if (!fdNdeta) {
90  SetStatus(-1);
91  Warning("Terminate", "No dN/deta histogram found");
92  return;
93  }
94  Printf("A total of %ld events", fOK);
95  TH1* copy = static_cast<TH1*>(fdNdeta->Clone("before"));
96  fOutput->Add(copy);
97 
98  fdNdeta->Scale(1. / fOK, "width");
99  fdNdeta->Draw();
100  SetStatus(0);
101  }
108  {
109  TObject* m1 = new TNamed("dNdeta", "");
110  m1->SetUniqueID(0x8); // Scale
111 
112  TList* ret = new TList;
113  ret->Add(m1);
114 
115  return ret;
116  }
117  ClassDef(Base,1);
118  };
119 
120  //====================================================================
125  struct NSD : public Base
126  {
133  NSD(Bool_t verbose=false, Int_t monitor=0)
134  : Base(verbose, monitor)
135  {}
143  {
144  if (fHeader->fType & 0x1) return false;
145  return true;
146  }
147  ClassDef(NSD,1);
148  };
149 
150  //====================================================================
155  struct INEL : public Base
156  {
163  INEL(Bool_t verbose=false, Int_t monitor=0)
164  : Base(verbose,monitor)
165  {}
171  virtual Bool_t ProcessHeader() { return true; }
172  ClassDef(INEL,1);
173  };
174 
175  //====================================================================
180  struct INELGt0 : public Base
181  {
188  INELGt0(Bool_t verbose=false, Int_t monitor=0)
189  : Base(verbose,monitor)
190  {}
197  {
198  if (fHeader->fType & 0x4) return true;
199  return false;
200  }
201  ClassDef(INELGt0,1);
202  };
203 
204  //====================================================================
205  struct Cent : public Base
206  {
208  const Long_t fMinEvents;
217  Cent(const char* method="V0M", Bool_t verbose=true, Int_t monitor=0)
218  : Base(verbose,monitor),
219  fHelper(method),
220  fMinEvents(100),
221  fCentBin(0)
222  {
224  }
230  virtual TH1* CreateSingle() { return Base::CreatedNdeta(); }
236  virtual void SlaveBegin(TTree* t)
237  {
238  Base::SlaveBegin(t);
239  // We use the parent fdNdeta histogram as a cache
240  fdNdeta->SetName("cache");
241  fdNdeta->SetTitle("Cache");
242 
244  fOutput->ls();
245  }
247  {
248  if (!Base::SetupEstimator()) return false;
249 
251  return true;
252  }
253 
257  virtual void Clear(Option_t* option="")
258  {
259  Base::Clear(option);
260  fdNdeta->Reset(); // Clear the cache
261  }
268  {
269  Double_t cent = GetCentrality();
271  fEventMult,
272  fHeader->fB,
274  fHeader->fNbin);
275  return fCentBin >= 0;
276  }
285  virtual void ProcessParticles()
286  {
287  // Check we got a bin
288  if (fCentBin < 0) return;
289 
290  // Find the histogram to update
291  TH1* out = static_cast<TH1*>(fHelper.fCentList->At(fCentBin-1));
292  // If we still have no histogram, return immediately
293  if (!out) return;
294 
295  // Use parent function to fill cache
297  if (fVerbose) {
298  Double_t n0 = fdNdeta->GetBinContent(fdNdeta->GetNbinsX()/2);
299  Double_t d0 = fdNdeta->GetXaxis()->GetBinWidth(fdNdeta->GetNbinsX()
300  /2);
301  Double_t eta0 = (n0 / d0);
302  Printf("Centrality %6.2f-%6.2f (bin %4d) "
303  "Nch=%8.1f/%4.2f=%8.1f out=%s (%d)",
304  fHelper.fCentAxis->GetBinLowEdge(fCentBin),
305  fHelper.fCentAxis->GetBinUpEdge(fCentBin),
306  fCentBin,
307  n0, d0, eta0,
308  out->GetName(),
309  Int_t(out->GetBinContent(0)));
310  }
311 
312  // Make sure we have nothing in the underflow bin.
313  fdNdeta->SetBinContent(0,0);
314 
315  // Add our cache to the appropriate bin
316  out->Add(fdNdeta);
317  // Increment underflow bin for the event count
318  out->SetBinContent(0, out->GetBinContent(0)+1);
319  }
326  virtual void Terminate()
327  {
328  fOK = GetEventCount();
329  if (fOK <= 0) {
330  SetStatus(-1);
331  Warning("Terminate", "No events selected");
332  return;
333  }
334 
335  if (!fHelper.Finalize(fOutput, fMinEvents))
336  SetStatus(-1);
337  }
344  {
345  TObject* m1 = new TNamed("cent", "hist text30");
346  TObject* m2 = new TNamed("centAcc", "hist text30");
347  TObject* m3 = new TNamed("byCent", "e");
348 
349  m3->SetUniqueID(0x8); // Scale
350  TList* ret = new TList;
351  ret->Add(m1);
352  ret->Add(m2);
353  ret->Add(m3);
354 
355  return ret;
356  }
357  ClassDef(Cent,1);
358  };
359 
360  //====================================================================
366  struct Mult : public Cent
367  {
375  Mult(const char* method="RefMult00d80",
376  Bool_t verbose=false,
377  Int_t monitor=0)
378  : Cent(method, verbose,monitor)
379  {
380  // +1 +2 +3 +3 +5, +5, +5, +5,+10,+10,+10,+10,+10,+10,+10
381  Double_t bins[]={ 0, 3, 6, 9, 14, 19, 24, 29, 39, 49, 59, 69, 79, 89, 99};
382  fHelper.SetCentAxis(14, bins);
383  }
384  ClassDef(Mult,1);
385  };
386 } // End of namespace
387 //====================================================================
388 /*
389  * The function to make our analyser
390  */
392 {
396  dNdetaMaker() : FastAnalysis::Maker("dNdeta") {}
397 
408  FastAnalysis* Make(const TString& subtype,
409  Int_t monitor,
410  Bool_t verbose,
411  TString& uopt)
412  {
413  TString t(subtype);
414  if (t.EqualTo("INEL")) return new dNdeta::INEL(verbose,monitor);
415  else if (t.EqualTo("NSD")) return new dNdeta::NSD(verbose,monitor);
416  else if (t.EqualTo("INELGt0")) return new dNdeta::INELGt0(verbose,monitor);
417  else if (t.BeginsWith("MULT") || t.BeginsWith("CENT")) {
418  TString w(t(4, t.Length()-4));
419  if (!(w.BeginsWith("RefMult") ||
420  w.BeginsWith("ZNA") ||
421  w.BeginsWith("ZNC") ||
422  w.BeginsWith("ZPA") ||
423  w.BeginsWith("ZPC") ||
424  w.BeginsWith("V0M") ||
425  w.BeginsWith("V0A") ||
426  w.BeginsWith("V0C") ||
427  w.BeginsWith("B") ||
428  w.IsNull())) {
429  Printf("Warning: dNdetaMaker::Make: Unknown estimator: %s",
430  w.Data());
431  return 0;
432  }
433  if (t.BeginsWith("MULT"))
434  return new dNdeta::Mult(w, verbose, monitor);
435  else
436  return new dNdeta::Cent(w, verbose, monitor);
437  }
438  Printf("Error: dNdetaMaker::Run: Invalid spec: %s", t.Data());
439  return 0;
440  }
445  void List() const
446  {
447  Printf(" INEL - inelastic");
448  Printf(" INELGt0 - inelastic with at least 1 particle in |eta|<1");
449  Printf(" NSD - Non-single diffractive");
450  Printf(" CENT<est> - Centrality classes. <est> is one of ");
451  Printf(" ZNA - ZNA signal");
452  Printf(" ZNC - ZNC signal");
453  Printf(" ZPA - ZPA signal");
454  Printf(" ZPC - ZPC signal");
455  Printf(" V0M - V0-A + -C");
456  Printf(" V0A - V0-A");
457  Printf(" V0C - V0-C");
458  }
462  const char* Script() const { return __FILE__; }
463 };
464 
465 // ------------------------------------------------------------------
466 // Create instance of maker
468 
469 #ifdef __MAKECINT__
470 #pragma link C++ nestedclasses;
471 #pragma link C++ namespace dNdeta;
472 #endif
473 
474 //
475 // EOF
476 //
TObject * GetOutputObject(const char *name, TClass *cls)
Definition: FastAnalysis.C:404
virtual Bool_t ProcessHeader()
virtual Bool_t ProcessParticle(const TParticle *p)
Cent(const char *method="V0M", Bool_t verbose=true, Int_t monitor=0)
double Double_t
Definition: External.C:58
FastShortHeader * fHeader
Definition: FastAnalysis.C:67
virtual Bool_t SetupEstimator()
Definition: FastAnalysis.C:181
TH1 * fCentHist
Definition: FastAnalysis.C:81
This script defines classes for looping over the data produced by FastSim.C.
NSD(Bool_t verbose=false, Int_t monitor=0)
virtual Bool_t ProcessHeader()
virtual void Fill(Double_t eta)
void CreateHistos(TCollection *output, TH1D *(*callback)())
FastCentHelper fHelper
ClassDef(INELGt0, 1)
dNdetaMaker * _dNdetaMaker
const char * Script() const
virtual void Clear(Option_t *option="")
virtual TList * GetMonitorObjects()
void SetCentAxis(Int_t n, Double_t low, Double_t high)
ClassDef(Cent, 1)
int Int_t
Definition: External.C:63
virtual Bool_t SetupEstimator()
virtual void SlaveBegin(TTree *t)
Int_t method
static TH1D * CreatedNdeta()
virtual Bool_t ProcessHeader()
Definition: External.C:212
virtual void Terminate()
Base(Bool_t verbose=false, Int_t monitor=0)
Long_t GetEventCount()
Definition: FastAnalysis.C:386
Double_t GetCentrality() const
Definition: FastAnalysis.C:356
virtual void ProcessParticles()
Definition: FastAnalysis.C:475
Bool_t fVerbose
Definition: FastAnalysis.C:71
Mult(const char *method="RefMult00d80", Bool_t verbose=false, Int_t monitor=0)
INEL(Bool_t verbose=false, Int_t monitor=0)
ClassDef(Base, 1)
virtual Bool_t ProcessHeader()
ULong64_t fEventMult
Definition: FastAnalysis.C:77
ULong_t fOK
Definition: FastAnalysis.C:83
void CreateDiagnostics(TCollection *output, TH1 *centHist)
Bool_t Finalize(TCollection *output, Long_t minEvents)
ClassDef(NSD, 1)
Int_t CheckCentrality(Double_t cent, Double_t mult, Double_t b, Double_t nPart, Double_t nBin)
virtual void ProcessParticles()
INELGt0(Bool_t verbose=false, Int_t monitor=0)
FastAnalysis * Make(const TString &subtype, Int_t monitor, Bool_t verbose, TString &uopt)
void List() const
virtual void Terminate()
const char Option_t
Definition: External.C:48
bool Bool_t
Definition: External.C:53
virtual TList * GetMonitorObjects()
Maker(const char *type="")
Definition: FastAnalysis.C:944
const Long_t fMinEvents
TString fCentMethod
Definition: FastAnalysis.C:79
virtual TH1 * CreateSingle()
virtual void SlaveBegin(TTree *)
Definition: External.C:196
ClassDef(Mult, 1)
virtual void Clear(Option_t *option="")
Definition: FastAnalysis.C:460
ClassDef(INEL, 1)