]> git.uio.no Git - u/mrichter/AliRoot.git/blob - FMD/scripts/DrawHits.C
Fixes and removals
[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 #include <TArrayF.h>
19
20 class DrawHits : public AliFMDInputHits
21 {
22 private:
23   TH2D* fElossVsPMQ; // Histogram 
24 public:
25   DrawHits() 
26   { 
27     Double_t emax = 1e4;
28     Double_t emin = 1e-5;
29     Int_t    n    = 901;
30     TArrayF tkine(n);
31     Float_t dp   = 1/TMath::Log10(emax/emin)/10;
32     Float_t pmin = TMath::Log10(emin);
33     tkine[0]     = emin;
34     for (Int_t i=1; i < tkine.fN; i++) {
35       Float_t el = pmin + i * dp;
36       tkine[i]   = TMath::Power(10, el);
37     }
38     Double_t dmin = .1;
39     Double_t dmax = 10;
40     Int_t    nd = 1000;
41     TArrayF eloss(nd+1);
42     eloss[0] = dmin;
43     for (Int_t i = 1; i < eloss.fN; i++){
44       eloss[i] = dmin + i * (dmax-dmin)/(eloss.fN-1);
45     }
46     fElossVsPMQ = new TH2D("bad", "#Delta E vs. p/(mq^{2})", 
47                            tkine.fN-1, tkine.fArray, eloss.fN-1, eloss.fArray);
48     fElossVsPMQ->SetXTitle("p/(mq^{2})");
49     fElossVsPMQ->SetYTitle("#Delta E [MeV]");
50   }
51   Bool_t ProcessHit(AliFMDHit* hit, TParticle*) 
52   {
53     if (!hit) {
54       std::cout << "No hit" << std::endl;
55       return kFALSE;
56     }
57     if (TMath::Abs(hit->Pdg()) != 211) return kTRUE;
58     // Float_t pmq = 0;
59     // Float_t q = hit->Q() / 3.;
60     // if (hit->M() != 0 && hit->Q() != 0) 
61     //  pmq = hit->P() / (hit->M()*q*q);
62     Float_t pmq = hit->P();
63     fElossVsPMQ->Fill(pmq, hit->Edep());
64     return kTRUE;
65   }
66   Bool_t Finish()
67   {
68     gStyle->SetPalette(1);
69     gStyle->SetOptTitle(0);
70     fElossVsPMQ->SetStats(kFALSE);
71     fElossVsPMQ->Draw("COLZ");
72     return kTRUE;
73   }
74 };
75
76 //____________________________________________________________________
77 //
78 // EOF
79 //