]> git.uio.no Git - u/mrichter/AliRoot.git/blob - FMD/scripts/DrawHitsDigits.C
Process codes corrected/completed.
[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   DrawHitsDigits(Int_t n=900, Double_t emin=1e-3, Double_t emax=10, 
43                  Int_t m=1100, Double_t amin=-0.5, Double_t amax=1099.5) 
44   { 
45     AddLoad(kDigits);
46     AddLoad(kHits);
47     TArrayF eloss(MakeLogScale(n, emin, emax));
48     TArrayF adcs(m+1);
49     adcs[0] = amin;
50     for (Int_t i = 1; i < m+1; i++) adcs[i] = adcs[i-1] + (amax-amin)/m;
51     fElossVsAdc = new TH2D("bad", "#Delta E vs. ADC", 
52                            eloss.fN-1, eloss.fArray, adcs.fN-1, adcs.fArray);
53     fElossVsAdc->SetXTitle("#Delta E/#Delta x [MeV/cm]");
54     fElossVsAdc->SetYTitle("ADC value");
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 ProcessHit(AliFMDHit* hit, TParticle*) 
67   {
68     // Cache the energy loss 
69     if (!hit) return kFALSE;
70     UShort_t det = hit->Detector();
71     Char_t   rng = hit->Ring();
72     UShort_t sec = hit->Sector();
73     UShort_t str = hit->Strip();
74     if (str > 511) {
75       AliWarning(Form("Bad strip number %d in hit", str));
76       return kTRUE;
77     }
78     fMap(det, rng, sec, str).fEdep += hit->Edep();
79     fMap(det, rng, sec, str).fN++;
80     return kTRUE;
81   }
82   //__________________________________________________________________
83   Bool_t ProcessDigit(AliFMDDigit* digit)
84   {
85     if (!digit) return kFALSE;
86     UShort_t det = digit->Detector();
87     Char_t   rng = digit->Ring();
88     UShort_t sec = digit->Sector();
89     UShort_t str = digit->Strip();
90     if (str > 511) {
91       AliWarning(Form("Bad strip number %d in digit", str));
92       return kFALSE;
93     }
94     fElossVsAdc->Fill(fMap(det, rng, sec, str).fEdep, digit->Counts());
95     return kTRUE;
96   }
97   //__________________________________________________________________
98   Bool_t Finish()
99   {
100     gStyle->SetPalette(1);
101     gStyle->SetOptTitle(0);
102     gStyle->SetCanvasColor(0);
103     gStyle->SetCanvasBorderSize(0);
104     gStyle->SetPadColor(0);
105     gStyle->SetPadBorderSize(0);
106     fElossVsAdc->SetStats(kFALSE);
107     fElossVsAdc->Draw("COLZ");
108     return kTRUE;
109   }
110
111   ClassDef(DrawHitsDigits,0);
112 };
113
114 //____________________________________________________________________
115 //
116 // EOF
117 //