]> git.uio.no Git - u/mrichter/AliRoot.git/blob - FMD/scripts/DrawHits.C
Added new library libFMDutil. This library contains utility classes that
[u/mrichter/AliRoot.git] / FMD / scripts / DrawHits.C
1 //
2 // $Id$
3 //
4 // Script that contains a class to draw hits, using the
5 // AliFMDInputHits class in the util library. 
6 //
7 // It draws the energy loss versus the p/(mq^2).  It can be overlayed
8 // with the Bethe-Bloc curve to show how the simulation behaves
9 // relative to the expected. 
10 //
11 // Use the script `Compile.C' to compile this class using ACLic. 
12 //
13 #include <TH2D.h>
14 #include <AliFMDHit.h>
15 #include <AliFMDInput.h>
16 #include <iostream>
17 #include <TStyle.h>
18
19 class DrawHits : public AliFMDInputHits
20 {
21 private:
22   TH2D* fElossVsPMQ; // Histogram 
23 public:
24   DrawHits() 
25   { 
26     fElossVsPMQ = new TH2D("bad", "#Delta E vs. p/(mq^{2})>1GeV", 
27                      1000, 1, 100, 50, 0.00001, 10);
28     fElossVsPMQ->SetXTitle("p/(mq^{2}) [GeV/GeV]");
29     fElossVsPMQ->SetYTitle("#Delta E [MeV]");
30   }
31   Bool_t ProcessHit(AliFMDHit* hit, TParticle*) 
32   {
33     if (!hit) {
34       std::cout << "No hit" << std::endl;
35       return kFALSE;
36     }
37     Float_t pmq = 0;
38     if (hit->M() != 0 && hit->Q() != 0) 
39       pmq = hit->P() / hit->M() / TMath::Power(hit->Q()/3, 2);
40     fElossVsPMQ->Fill(pmq, hit->Edep());
41     return kTRUE;
42   }
43   Bool_t Finish()
44   {
45     gStyle->SetPalette(1);
46     fElossVsPMQ->SetStats(kFALSE);
47     fElossVsPMQ->Draw("COLZ");
48     return kTRUE;
49   }
50 };
51
52 //____________________________________________________________________
53 //
54 // EOF
55 //