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