]> git.uio.no Git - u/mrichter/AliRoot.git/blob - FMD/scripts/DrawDigitsRecs.C
Some changes to AliFMDInput and scripts
[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   DrawDigitsRecs(Int_t m=1100, Double_t amin=-0.5, Double_t amax=1099.5,
46                  Int_t n=105, Double_t mmin=-0.5, Double_t mmax=20.5) 
47   { 
48     AddLoad(kDigits);
49     AddLoad(kRecPoints);
50     fAdcVsSingle = new TH2D("adcVsSingle", "ADC vs. Multiplicity (strip)", 
51                             m, amin, amax, n, mmin, mmax);
52     fAdcVsSingle->Sumw2();
53     fAdcVsSingle->SetXTitle("ADC value");
54     fAdcVsSingle->SetYTitle("Strip Multiplicity");
55   }
56   //__________________________________________________________________
57   /** Begining of event
58       @param ev Event number
59       @return @c false on error */
60   Bool_t Begin(Int_t ev) 
61   {
62     fMap.Reset();
63     return AliFMDInput::Begin(ev);
64   }
65   //__________________________________________________________________
66   Bool_t ProcessDigit(AliFMDDigit* digit) 
67   {
68     // Cache the energy loss 
69     if (!digit) return kFALSE;
70     UShort_t det = digit->Detector();
71     Char_t   rng = digit->Ring();
72     UShort_t sec = digit->Sector();
73     UShort_t str = digit->Strip();
74     if (str > 511) {
75       AliWarning(Form("Bad strip number %d in digit", str));
76       return kTRUE;
77     }
78     fMap(det, rng, sec, str) = digit->Counts();
79     return kTRUE;
80   }
81   //__________________________________________________________________
82   Bool_t ProcessRecPoint(AliFMDRecPoint* single)
83   {
84     if (!single) return kFALSE;
85     UShort_t det = single->Detector();
86     Char_t   rng = single->Ring();
87     UShort_t sec = single->Sector();
88     UShort_t str = single->Strip();
89     if (str > 511) {
90       AliWarning(Form("Bad strip number %d in single", str));
91       return kFALSE;
92     }
93     fAdcVsSingle->Fill(fMap(det, rng, sec, str), single->Particles());
94     return kTRUE;
95   }
96   //__________________________________________________________________
97   Bool_t Finish()
98   {
99     gStyle->SetPalette(1);
100     gStyle->SetOptTitle(0);
101     gStyle->SetCanvasColor(0);
102     gStyle->SetCanvasBorderSize(0);
103     gStyle->SetPadColor(0);
104     gStyle->SetPadBorderSize(0);
105     gStyle->SetOptLogz(kTRUE);
106     fAdcVsSingle->SetStats(0);
107     fAdcVsSingle->Draw("COLZ");
108     return kTRUE;
109   }
110
111   ClassDef(DrawDigitsRecs,0);
112   
113 };
114
115 //____________________________________________________________________
116 //
117 // EOF
118 //