]> git.uio.no Git - u/mrichter/AliRoot.git/blob - FMD/scripts/DrawHitsDigits.C
Fixed some coding style violations.
[u/mrichter/AliRoot.git] / FMD / scripts / DrawHitsDigits.C
1 //
2 // $Id$
3 //
4 // Script that contains a class to draw eloss from hits, versus ADC
5 // counts from digits, using the 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 <AliFMDDigit.h>
16 #include <AliFMDInput.h>
17 #include <AliFMDEdepMap.h>
18 #include <iostream>
19 #include <TStyle.h>
20 #include <TArrayF.h>
21
22 class DrawHitsDigits : public AliFMDInputHits
23 {
24 private:
25   TH2D* fElossVsAdc; // Histogram 
26   AliFMDEdepMap fMap;
27 public:
28   TArrayF MakeLogScale(Int_t n, Double_t min, Double_t max) 
29   {
30     TArrayF bins(n+1);
31     Float_t dp   = n / TMath::Log10(max / min);
32     Float_t pmin = TMath::Log10(min);
33     bins[0]      = min;
34     for (Int_t i = 1; i < n+1; i++) {
35       Float_t p = pmin + i / dp;
36       bins[i]   = TMath::Power(10, p);
37     }
38     return bins;
39   }
40   DrawHitsDigits(Int_t n=900, Double_t emin=1e-3, Double_t emax=10, 
41                  Int_t m=1100, Double_t amin=-0.5, Double_t amax=1099.5) 
42   { 
43     AddLoad(kDigits);
44     TArrayF eloss(MakeLogScale(n, emin, emax));
45     TArrayF adcs(m+1);
46     adcs[0] = amin;
47     for (Int_t i = 1; i < m+1; i++) adcs[i] = adcs[i-1] + (amax-amin)/m;
48     fElossVsAdc = new TH2D("bad", "#Delta E vs. ADC", 
49                            eloss.fN-1, eloss.fArray, adcs.fN-1, adcs.fArray);
50     fElossVsAdc->SetXTitle("#Delta E/#Delta x [MeV/cm]");
51     fElossVsAdc->SetYTitle("ADC value");
52   }
53   Bool_t Begin(Int_t ev) 
54   {
55     fMap.Reset();
56     return AliFMDInputHits::Begin(ev);
57   }
58   Bool_t ProcessHit(AliFMDHit* hit, TParticle*) 
59   {
60     // Cache the energy loss 
61     if (!hit) return kFALSE;
62     UShort_t det = hit->Detector();
63     Char_t   rng = hit->Ring();
64     UShort_t sec = hit->Sector();
65     UShort_t str = hit->Strip();
66     fMap(det, rng, sec, str).fEdep += hit->Edep();
67     fMap(det, rng, sec, str).fN++;
68     return kTRUE;
69   }
70   Bool_t Event() 
71   {
72     if (!AliFMDInputHits::Event()) return kFALSE;
73     Int_t nEv = fTreeD->GetEntries();
74     for (Int_t i = 0; i < nEv; i++) {
75       Int_t digitRead  = fTreeD->GetEntry(i);
76       if (digitRead <= 0) continue;
77       Int_t nDigit = fArrayD->GetEntries();
78       if (nDigit <= 0) continue;
79       for (Int_t j = 0; j < nDigit; j++) {
80         AliFMDDigit* digit = static_cast<AliFMDDigit*>(fArrayD->At(j));
81         if (!digit) continue;
82         UShort_t det = digit->Detector();
83         Char_t   rng = digit->Ring();
84         UShort_t sec = digit->Sector();
85         UShort_t str = digit->Strip();
86         fElossVsAdc->Fill(fMap(det, rng, sec, str).fEdep, digit->Counts());
87       }    
88     }
89     return kTRUE;
90   }
91   Bool_t Finish()
92   {
93     gStyle->SetPalette(1);
94     gStyle->SetOptTitle(0);
95     gStyle->SetCanvasColor(0);
96     gStyle->SetCanvasBorderSize(0);
97     gStyle->SetPadColor(0);
98     gStyle->SetPadBorderSize(0);
99     fElossVsAdc->SetStats(kFALSE);
100     fElossVsAdc->Draw("COLZ");
101     return kTRUE;
102   }
103 };
104
105 //____________________________________________________________________
106 //
107 // EOF
108 //