]> git.uio.no Git - u/mrichter/AliRoot.git/blob - FMD/scripts/DrawHitsDigits.C
Check-in to fix a few problems:
[u/mrichter/AliRoot.git] / FMD / scripts / DrawHitsDigits.C
1 //____________________________________________________________________
2 //
3 // $Id$
4 //
5 // Script that contains a class to draw eloss from hits, versus ADC
6 // counts from digits, using the 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 <TH2D.h>
15 #include <AliFMDHit.h>
16 #include <AliFMDDigit.h>
17 #include <AliFMDInput.h>
18 #include <AliFMDEdepMap.h>
19 #include <iostream>
20 #include <TStyle.h>
21 #include <TArrayF.h>
22 #include <AliLog.h>
23 #include <TMath.h>
24
25 /** @class DrawHitsDigits
26     @brief Draw hit energy loss versus digit ADC
27     @code 
28     Root> .L Compile.C
29     Root> Compile("DrawHitsDigits.C")
30     Root> DrawHitsDigits c
31     Root> c.Run();
32     @endcode
33     @ingroup FMD_script
34  */
35 class DrawHitsDigits : public AliFMDInput
36 {
37 private:
38   TH2D* fElossVsAdc; // Histogram 
39   AliFMDEdepMap fMap;
40 public:
41   //__________________________________________________________________
42   TArrayF MakeLogScale(Int_t n, Double_t min, Double_t max) 
43   {
44     TArrayF bins(n+1);
45     Float_t dp   = n / TMath::Log10(max / min);
46     Float_t pmin = TMath::Log10(min);
47     bins[0]      = min;
48     for (Int_t i = 1; i < n+1; i++) {
49       Float_t p = pmin + i / dp;
50       bins[i]   = TMath::Power(10, p);
51     }
52     return bins;
53   }
54   //__________________________________________________________________
55   DrawHitsDigits(Int_t n=900, Double_t emin=1e-3, Double_t emax=10, 
56                  Int_t m=1100, Double_t amin=-0.5, Double_t amax=1099.5) 
57   { 
58     AddLoad(kDigits);
59     AddLoad(kHits);
60     TArrayF eloss(MakeLogScale(n, emin, emax));
61     TArrayF adcs(m+1);
62     adcs[0] = amin;
63     for (Int_t i = 1; i < m+1; i++) adcs[i] = adcs[i-1] + (amax-amin)/m;
64     fElossVsAdc = new TH2D("bad", "#Delta E vs. ADC", 
65                            eloss.fN-1, eloss.fArray, adcs.fN-1, adcs.fArray);
66     fElossVsAdc->SetXTitle("#Delta E/#Delta x [MeV/cm]");
67     fElossVsAdc->SetYTitle("ADC value");
68   }
69   //__________________________________________________________________
70   /** Begining of event
71       @param ev Event number
72       @return @c false on error */
73   Bool_t Begin(Int_t ev) 
74   {
75     fMap.Reset();
76     return AliFMDInput::Begin(ev);
77   }
78   //__________________________________________________________________
79   Bool_t ProcessHit(AliFMDHit* hit, TParticle*) 
80   {
81     // Cache the energy loss 
82     if (!hit) return kFALSE;
83     UShort_t det = hit->Detector();
84     Char_t   rng = hit->Ring();
85     UShort_t sec = hit->Sector();
86     UShort_t str = hit->Strip();
87     if (str > 511) {
88       AliWarning(Form("Bad strip number %d in hit", str));
89       return kTRUE;
90     }
91     fMap(det, rng, sec, str).fEdep += hit->Edep();
92     fMap(det, rng, sec, str).fN++;
93     return kTRUE;
94   }
95   //__________________________________________________________________
96   Bool_t ProcessDigit(AliFMDDigit* digit)
97   {
98     if (!digit) return kFALSE;
99     UShort_t det = digit->Detector();
100     Char_t   rng = digit->Ring();
101     UShort_t sec = digit->Sector();
102     UShort_t str = digit->Strip();
103     if (str > 511) {
104       AliWarning(Form("Bad strip number %d in digit", str));
105       return kFALSE;
106     }
107     fElossVsAdc->Fill(fMap(det, rng, sec, str).fEdep, digit->Counts());
108     return kTRUE;
109   }
110   //__________________________________________________________________
111   Bool_t Finish()
112   {
113     gStyle->SetPalette(1);
114     gStyle->SetOptTitle(0);
115     gStyle->SetCanvasColor(0);
116     gStyle->SetCanvasBorderSize(0);
117     gStyle->SetPadColor(0);
118     gStyle->SetPadBorderSize(0);
119     fElossVsAdc->SetStats(kFALSE);
120     fElossVsAdc->Draw("COLZ");
121     return kTRUE;
122   }
123
124   ClassDef(DrawHitsDigits,0);
125 };
126
127 //____________________________________________________________________
128 //
129 // EOF
130 //