AliRoot Core  edcc906 (edcc906)
PoissonHit.C
Go to the documentation of this file.
1 //____________________________________________________________________
2 //
3 // $Id$
4 //
5 // Script that contains a class to draw hits, using the
6 // AliFMDInputHits class in the util library.
7 //
8 // It draws the energy loss versus the p/(mq^2). It can be overlayed
9 // with the Bethe-Bloc curve to show how the simulation behaves
10 // relative to the expected.
11 //
12 // Use the script `Compile.C' to compile this class using ACLic.
13 //
14 #include "Poisson.C"
15 #include <TMath.h>
16 #include <TCanvas.h>
17 #include <AliFMDHit.h>
18 #include <AliFMDGeometry.h>
19 
31 class PoissonHit : public Poisson
32 {
33 protected:
34  TH2D* fHits; // Histogram
35  TH2D* fDiff; // Histogram
36 public:
45  PoissonHit(Double_t threshold=.3,
46  Int_t nEta=120, Float_t minEta=-6, Float_t maxEta=6,
47  Int_t nPhi=4, Float_t minPhi=0, Float_t maxPhi=2*TMath::Pi())
48  : Poisson(threshold, nEta, minEta, maxEta, nPhi, minPhi, maxPhi)
49  {
50  AddLoad(kHits);
52  fHits = new TH2D(*fEmpty);
53  fHits->SetName("hits");
54  fHits->SetTitle("# of hits");
55  fDiff = new TH2D(*fEmpty);
56  fDiff->SetName("diff");
57  fDiff->SetTitle("Difference between poisson and hits");
58  fHits->SetXTitle("#eta");
59  fHits->SetYTitle("#phi");
60  fHits->SetZTitle("N");
61  fDiff->SetXTitle("#eta");
62  fDiff->SetYTitle("#phi");
63  fDiff->SetZTitle("#frac{N_{hit}-N_{poisson}}{N_{hit}}");
64  }
67  virtual Bool_t Init()
68  {
69  if (!Poisson::Init()) return kFALSE;
72  return kTRUE;
73  }
82  void PhysicalCoordinates(Double_t x, Double_t y, Double_t z,
83  Double_t& eta, Double_t& phi)
84  {
85  Double_t r, theta;
86  phi = TMath::ATan2(y, x);
87  r = TMath::Sqrt(y * y + x * x);
88  theta = TMath::ATan2(r, z);
89  eta = -TMath::Log(TMath::Tan(theta / 2));
90  if (phi < 0) phi += 2 * TMath::Pi();
91  }
95  virtual Bool_t ProcessHit(AliFMDHit* hit, TParticle*)
96  {
97  Double_t x, y, z;
98 #if 0
100  geom->Detector2XYZ(hit->Detector(),hit->Ring(),hit->Sector(),
101  hit->Strip(),x,y,z);
102 #else
103  x = hit->X();
104  y = hit->Y();
105  z = hit->Z();
106 #endif
107  Double_t eta, phi;
108  PhysicalCoordinates(x, y, z, eta, phi);
109  fHits->Fill(eta, phi);
110  return kTRUE;
111  }
115  virtual Bool_t Begin(Int_t event)
116  {
117  if (!Poisson::Begin(event)) return kFALSE;
118  fHits->Clear();
119  fDiff->Clear();
120  return kTRUE;
121  }
125  virtual Bool_t End()
126  {
127  if (!Poisson::End()) return kFALSE;
128  fDiff->Add(fMult,fHits,-1.,1.);
129  fDiff->Divide(fHits);
130  if (!gROOT->IsBatch()) {
131  gStyle->SetPalette(1);
132  TCanvas* c1 = new TCanvas("hits", "Hit multiplicity");
133  c1->SetFillColor(0);
134  fHits->Draw("colz");
135  TCanvas* c2 = new TCanvas("diff", "Difference between Hit and poisson");
136  c2->SetFillColor(0);
137  fDiff->Draw("colz");
138  TCanvas* c3 = new TCanvas("empty", "# of Empty strips");
139  c3->SetFillColor(0);
140  fEmpty->Draw("colz");
141  TCanvas* c4 = new TCanvas("total", "Total # of strips");
142  c4->SetFillColor(0);
143  fTotal->Draw("colz");
144  }
145  return kTRUE;
146  }
147 
148  ClassDef(PoissonHit,0);
149 };
150 
151 //____________________________________________________________________
152 //
153 // EOF
154 //
virtual Bool_t End()
Definition: PoissonHit.C:125
Geometry mananger for the FMD.
void PhysicalCoordinates(Double_t x, Double_t y, Double_t z, Double_t &eta, Double_t &phi)
Definition: PoissonHit.C:82
virtual Bool_t Begin(Int_t event)
Definition: Poisson.C:89
TStyle * gStyle
Hit in the FMD.
TH2D * fHits
Definition: PoissonHit.C:34
virtual void InitTransformations(Bool_t force=kFALSE)
virtual Bool_t ProcessHit(AliFMDHit *hit, TParticle *)
Definition: PoissonHit.C:95
TH2D * fTotal
Definition: Poisson.C:39
TROOT * gROOT
Float_t Z() const
Definition: AliHit.h:23
Make a poisson reconstruction.
Definition: Poisson.C:35
Float_t X() const
Definition: AliHit.h:21
TH2D * fDiff
Definition: PoissonHit.C:35
Float_t Y() const
Definition: AliHit.h:22
virtual Bool_t Init()
Definition: PoissonHit.C:67
Char_t Ring() const
Definition: AliFMDHit.h:76
TH2D * fEmpty
Definition: Poisson.C:38
Singleton object of FMD geometry descriptions and parameters. This class is a singleton that handles ...
UShort_t Detector() const
Definition: AliFMDHit.h:74
virtual void AddLoad(ETrees tree)
Definition: AliFMDInput.h:134
void Detector2XYZ(UShort_t detector, Char_t ring, UShort_t sector, UShort_t strip, Double_t &x, Double_t &y, Double_t &z) const
UShort_t Sector() const
Definition: AliFMDHit.h:78
virtual void Init()
UShort_t Strip() const
Definition: AliFMDHit.h:80
TH2D * fMult
Definition: Poisson.C:40
AliFMDhit is the hit class for the FMD. Hits are the information that comes from a Monte Carlo at eac...
Definition: AliFMDHit.h:30
virtual Bool_t End()
Definition: Poisson.C:137
PoissonHit(Double_t threshold=.3, Int_t nEta=120, Float_t minEta=-6, Float_t maxEta=6, Int_t nPhi=4, Float_t minPhi=0, Float_t maxPhi=2 *TMath::Pi())
Definition: PoissonHit.C:45
TEveGeoShape * geom
Definition: tpc_tracks.C:10
static AliFMDGeometry * Instance()
Make a poisson reconstruction and compare to simulated hits.
Definition: PoissonHit.C:31
virtual Bool_t Init()
Definition: Poisson.C:79
virtual Bool_t Begin(Int_t event)
Definition: PoissonHit.C:115