]> git.uio.no Git - u/mrichter/AliRoot.git/blob - FMD/scripts/DrawDigitsRecs.C
0. General code clean-up, including messages, and the like.
[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 : public AliFMDInputDigits
28 {
29 private:
30   TH2D* fAdcVsSingle; // Histogram 
31   AliFMDUShortMap fMap;
32 public:
33   //__________________________________________________________________
34   TArrayF MakeLogScale(Int_t n, Double_t min, Double_t max) 
35   {
36     TArrayF bins(n+1);
37     Float_t dp   = n / TMath::Log10(max / min);
38     Float_t pmin = TMath::Log10(min);
39     bins[0]      = min;
40     for (Int_t i = 1; i < n+1; i++) {
41       Float_t p = pmin + i / dp;
42       bins[i]   = TMath::Power(10, p);
43     }
44     return bins;
45   }
46   //__________________________________________________________________
47   DrawDigitsRecs(Int_t m=1100, Double_t amin=-0.5, Double_t amax=1099.5,
48                  Int_t n=105, Double_t mmin=-0.5, Double_t mmax=20.5) 
49   { 
50     AddLoad(kRecPoints);
51     fAdcVsSingle = new TH2D("adcVsSingle", "ADC vs. Multiplicity (strip)", 
52                             m, amin, amax, n, mmin, mmax);
53     fAdcVsSingle->SetXTitle("ADC value");
54     fAdcVsSingle->SetYTitle("Strip Multiplicity");
55   }
56   //__________________________________________________________________
57   Bool_t Begin(Int_t ev) 
58   {
59     fMap.Reset();
60     return AliFMDInputDigits::Begin(ev);
61   }
62   //__________________________________________________________________
63   Bool_t ProcessDigit(AliFMDDigit* digit) 
64   {
65     // Cache the energy loss 
66     if (!digit) return kFALSE;
67     UShort_t det = digit->Detector();
68     Char_t   rng = digit->Ring();
69     UShort_t sec = digit->Sector();
70     UShort_t str = digit->Strip();
71     if (str > 511) {
72       AliWarning(Form("Bad strip number %d in digit", str));
73       return kTRUE;
74     }
75     fMap(det, rng, sec, str) = digit->Counts();
76     return kTRUE;
77   }
78   //__________________________________________________________________
79   Bool_t ProcessRecPoint(AliFMDRecPoint* single)
80   {
81     if (!single) continue;
82     UShort_t det = single->Detector();
83     Char_t   rng = single->Ring();
84     UShort_t sec = single->Sector();
85     UShort_t str = single->Strip();
86     if (str > 511) {
87       AliWarning(Form("Bad strip number %d in single", str));
88       continue;
89     }
90     fAdcVsSingle->Fill(fMap(det, rng, sec, str), single->Particles());
91     return kTRUE;
92   }
93   //__________________________________________________________________
94   Bool_t Finish()
95   {
96     gStyle->SetPalette(1);
97     gStyle->SetOptTitle(0);
98     gStyle->SetCanvasColor(0);
99     gStyle->SetCanvasBorderSize(0);
100     gStyle->SetPadColor(0);
101     gStyle->SetPadBorderSize(0);
102
103     fAdcVsSingle->SetStats(kFALSE);
104     fAdcVsSingle->Draw("COLZ");
105     return kTRUE;
106   }
107
108   ClassDef(DrawDigitsRecs,0);
109   
110 };
111
112 //____________________________________________________________________
113 //
114 // EOF
115 //