]> git.uio.no Git - u/mrichter/AliRoot.git/blob - FMD/scripts/DrawDigitsRecs.C
Various small fixes. Make sure Emacs knows it's C++ mode, and the like.
[u/mrichter/AliRoot.git] / FMD / scripts / DrawDigitsRecs.C
1 //
2 // $Id$
3 //
4 // Script that contains a class to draw eloss from hits, versus ADC
5 // counts from digits, using the AliFMDInputHits class in the util library. 
6 //
7 // It draws the energy loss versus the p/(mq^2).  It can be overlayed
8 // with the Bethe-Bloc curve to show how the simulation behaves
9 // relative to the expected. 
10 //
11 // Use the script `Compile.C' to compile this class using ACLic. 
12 //
13 #include <TH2D.h>
14 #include <AliFMDHit.h>
15 #include <AliFMDDigit.h>
16 #include <AliFMDInput.h>
17 #include <AliFMDUShortMap.h>
18 #include <AliFMDFloatMap.h>
19 #include <AliFMDMultStrip.h>
20 #include <AliFMDMultRegion.h>
21 #include <iostream>
22 #include <TStyle.h>
23 #include <TArrayF.h>
24 #include <TCanvas.h>
25
26 class DrawDigitsRecs : public AliFMDInputDigits
27 {
28 private:
29   TH2D* fAdcVsSingle; // Histogram 
30   TH2D* fAdcVsRegion; // Histogram 
31   TH2D* fSingleVsRegion; // Histogram 
32   AliFMDUShortMap fMap;
33   AliFMDFloatMap fEta;
34   AliFMDFloatMap fPhi;
35   AliFMDFloatMap fMult;
36 public:
37   TArrayF MakeLogScale(Int_t n, Double_t min, Double_t max) 
38   {
39     TArrayF bins(n+1);
40     Float_t dp   = n / TMath::Log10(max / min);
41     Float_t pmin = TMath::Log10(min);
42     bins[0]      = min;
43     for (Int_t i = 1; i < n+1; i++) {
44       Float_t p = pmin + i / dp;
45       bins[i]   = TMath::Power(10, p);
46     }
47     return bins;
48   }
49   DrawDigitsRecs(Int_t m=1100, Double_t amin=-0.5, Double_t amax=1099.5,
50                  Int_t n=105, Double_t mmin=-0.5, Double_t mmax=20.5) 
51   { 
52     AddLoad(kRecPoints);
53     fAdcVsSingle = new TH2D("adcVsSingle", "ADC vs. Multiplicity (strip)", 
54                             m, amin, amax, n, mmin, mmax);
55     fAdcVsSingle->SetXTitle("ADC value");
56     fAdcVsSingle->SetYTitle("Strip Multiplicity");
57
58     mmin *= 20;
59     mmax *= 20;
60     amin *= 50;
61     amax *= 50;
62     fAdcVsRegion = new TH2D("adcVsRegion", "ADC vs Muliplcity (region)", 
63                             m, amin, amax, n, mmin, mmax);
64     fAdcVsRegion->SetXTitle("ADC value");
65     fAdcVsRegion->SetYTitle("Region Multiplicity");
66
67     fSingleVsRegion = new TH2D("singleVsRegion", "Single vs Region", 
68                                n, mmin, mmax, n, mmin, mmax);
69     fSingleVsRegion->SetXTitle("Strip Multiplicity");
70     fSingleVsRegion->SetYTitle("Region Multiplicity");
71   }
72   Bool_t Begin(Int_t ev) 
73   {
74     fMap.Reset();
75     fEta.Reset();
76     fPhi.Reset();
77     return AliFMDInputDigits::Begin(ev);
78   }
79   Bool_t ProcessDigit(AliFMDDigit* digit) 
80   {
81     // Cache the energy loss 
82     if (!digit) return kFALSE;
83     UShort_t det = digit->Detector();
84     Char_t   rng = digit->Ring();
85     UShort_t sec = digit->Sector();
86     UShort_t str = digit->Strip();
87     if (str > 511) {
88       AliWarning(Form("Bad strip number %d in digit", str));
89       return kTRUE;
90     }
91     fMap(det, rng, sec, str) = digit->Counts();
92     return kTRUE;
93   }
94   Bool_t Event() 
95   {
96     if (!AliFMDInputDigits::Event()) return kFALSE;
97     Int_t nEv = fTreeR->GetEntries();
98     for (Int_t i = 0; i < nEv; i++) {
99       Int_t recoRead  = fTreeR->GetEntry(i);
100       if (recoRead <= 0) continue;
101       Int_t nSingle = fArrayN->GetEntries();
102       if (nSingle <= 0) continue;
103       for (Int_t j = 0; j < nSingle; j++) {
104         AliFMDMultStrip* single = 
105           static_cast<AliFMDMultStrip*>(fArrayN->At(j));
106         if (!single) continue;
107         UShort_t det = single->Detector();
108         Char_t   rng = single->Ring();
109         UShort_t sec = single->Sector();
110         UShort_t str = single->Strip();
111         if (str > 511) {
112           AliWarning(Form("Bad strip number %d in single", str));
113           continue;
114         }
115         fAdcVsSingle->Fill(fMap(det, rng, sec, str), single->Particles());
116         fEta(det, rng, sec, str)  = single->Eta();
117         fPhi(det, rng, sec, str)  = single->Phi() * 180 / TMath::Pi();
118         fMult(det, rng, sec, str) = single->Particles();
119       }    
120
121       Int_t nRegion = fArrayP->GetEntries();
122       if (nRegion <= 0) continue;
123       for (Int_t j = 0; j < nRegion; j++) {
124         AliFMDMultRegion* region = 
125           static_cast<AliFMDMultRegion*>(fArrayP->At(j));
126         if (!region) continue;
127         UShort_t det  = region->Detector();
128         Char_t   rng  = region->Ring();
129         UInt_t   suma = 0;
130         Float_t  summ = 0;
131         
132         for (size_t sec = 0; sec < AliFMDMap::kMaxSectors; sec++) {
133           Float_t phi = fPhi(det, rng, sec, 0);     
134           Bool_t ok   = (phi <= region->MaxPhi() && phi >= region->MinPhi());
135           if (!ok) continue;
136           for (size_t str = 0; str < AliFMDMap::kMaxStrips; str++) {
137             Float_t eta  = fEta(det, rng, sec, str);
138             Float_t sign = eta < 0 ? -1. : 1.;
139             ok           = (sign * eta <= sign * region->MaxEta() 
140                             && sign * eta >= sign * region->MinEta());
141             if (!ok) continue;
142             suma += fMap(det, rng, sec, str);
143             summ += fMult(det, rng, sec, str);
144           }
145         }
146         fAdcVsRegion->Fill(suma, region->Particles()); 
147         fSingleVsRegion->Fill(summ, region->Particles());
148      }    
149     }
150     return kTRUE;
151   }
152   Bool_t Finish()
153   {
154     gStyle->SetPalette(1);
155     gStyle->SetOptTitle(0);
156     gStyle->SetCanvasColor(0);
157     gStyle->SetCanvasBorderSize(0);
158     gStyle->SetPadColor(0);
159     gStyle->SetPadBorderSize(0);
160
161     new TCanvas("c1", "ADC vs. Strip Multiplicity");
162     fAdcVsSingle->SetStats(kFALSE);
163     fAdcVsSingle->Draw("COLZ");
164
165     new TCanvas("c2", "ADC vs. Region Multiplicity");
166     fAdcVsRegion->SetStats(kFALSE);
167     fAdcVsRegion->Draw("COLZ");
168
169     new TCanvas("c3", "Single vs. Region Multiplicity");
170     fSingleVsRegion->SetStats(kFALSE);
171     fSingleVsRegion->Draw("COLZ");
172
173     return kTRUE;
174   }
175
176   ClassDef(DrawDigitsRecs,0);
177   
178 };
179
180 //____________________________________________________________________
181 //
182 // EOF
183 //