]> git.uio.no Git - u/mrichter/AliRoot.git/blob - FMD/scripts/DrawDigitsRecs.C
Updated reconstruction to (optionally) make angle correction.
[u/mrichter/AliRoot.git] / FMD / scripts / DrawDigitsRecs.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 <AliFMDUShortMap.h>
19 #include <AliFMDFloatMap.h>
20 #include <AliFMDRecPoint.h>
21 #include <AliESDFMD.h>
22 #include <AliLog.h>
23 #include <iostream>
24 #include <TStyle.h>
25 #include <TArrayF.h>
26 #include <TCanvas.h>
27
28 /** @class DrawDigitsRecs
29     @brief Draw digit ADC versus Rec point mult
30     @code 
31     Root> .L Compile.C
32     Root> Compile("DrawDigitsRecs.C")
33     Root> DrawDigitsRecs c
34     Root> c.Run();
35     @endcode
36     @ingroup FMD_script
37  */
38 class DrawDigitsRecs : public AliFMDInput
39 {
40 private:
41   TH2D* fAdcVsSingle; // Histogram 
42   AliFMDUShortMap fMap;
43 public:
44   //__________________________________________________________________
45   TArrayF MakeLogScale(Int_t n, Double_t min, Double_t max) 
46   {
47     TArrayF bins(n+1);
48     Float_t dp   = n / TMath::Log10(max / min);
49     Float_t pmin = TMath::Log10(min);
50     bins[0]      = min;
51     for (Int_t i = 1; i < n+1; i++) {
52       Float_t p = pmin + i / dp;
53       bins[i]   = TMath::Power(10, p);
54     }
55     return bins;
56   }
57   //__________________________________________________________________
58   DrawDigitsRecs(Int_t m=1100, Double_t amin=-0.5, Double_t amax=1099.5,
59                  Int_t n=105, Double_t mmin=-0.5, Double_t mmax=20.5) 
60   { 
61     AddLoad(kDigits);
62     AddLoad(kRecPoints);
63     fAdcVsSingle = new TH2D("adcVsSingle", "ADC vs. Multiplicity (strip)", 
64                             m, amin, amax, n, mmin, mmax);
65     fAdcVsSingle->SetXTitle("ADC value");
66     fAdcVsSingle->SetYTitle("Strip Multiplicity");
67   }
68   //__________________________________________________________________
69   /** Begining of event
70       @param ev Event number
71       @return @c false on error */
72   Bool_t Begin(Int_t ev) 
73   {
74     fMap.Reset();
75     return AliFMDInput::Begin(ev);
76   }
77   //__________________________________________________________________
78   Bool_t ProcessDigit(AliFMDDigit* digit) 
79   {
80     // Cache the energy loss 
81     if (!digit) return kFALSE;
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     if (str > 511) {
87       AliWarning(Form("Bad strip number %d in digit", str));
88       return kTRUE;
89     }
90     fMap(det, rng, sec, str) = digit->Counts();
91     return kTRUE;
92   }
93   //__________________________________________________________________
94   Bool_t ProcessRecPoint(AliFMDRecPoint* single)
95   {
96     if (!single) return kFALSE;
97     UShort_t det = single->Detector();
98     Char_t   rng = single->Ring();
99     UShort_t sec = single->Sector();
100     UShort_t str = single->Strip();
101     if (str > 511) {
102       AliWarning(Form("Bad strip number %d in single", str));
103       return kFALSE;
104     }
105     fAdcVsSingle->Fill(fMap(det, rng, sec, str), single->Particles());
106     return kTRUE;
107   }
108   //__________________________________________________________________
109   Bool_t Finish()
110   {
111     gStyle->SetPalette(1);
112     gStyle->SetOptTitle(0);
113     gStyle->SetCanvasColor(0);
114     gStyle->SetCanvasBorderSize(0);
115     gStyle->SetPadColor(0);
116     gStyle->SetPadBorderSize(0);
117
118     fAdcVsSingle->SetStats(kFALSE);
119     fAdcVsSingle->Draw("COLZ");
120     return kTRUE;
121   }
122
123   ClassDef(DrawDigitsRecs,0);
124   
125 };
126
127 //____________________________________________________________________
128 //
129 // EOF
130 //