1 //____________________________________________________________________
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.
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.
12 // Use the script `Compile.C' to compile this class using ACLic.
15 #include <AliFMDHit.h>
16 #include <AliFMDDigit.h>
17 #include <AliFMDInput.h>
18 #include <AliFMDEdepMap.h>
24 /** @class DrawHitsDigits
25 @brief Draw hit energy loss versus digit ADC
28 Root> Compile("DrawHitsDigits.C")
29 Root> DrawHitsDigits c
34 class DrawHitsDigits : public AliFMDInput
37 TH2D* fElossVsAdc; // Histogram
40 //__________________________________________________________________
41 TArrayF MakeLogScale(Int_t n, Double_t min, Double_t max)
44 Float_t dp = n / TMath::Log10(max / min);
45 Float_t pmin = TMath::Log10(min);
47 for (Int_t i = 1; i < n+1; i++) {
48 Float_t p = pmin + i / dp;
49 bins[i] = TMath::Power(10, p);
53 //__________________________________________________________________
54 DrawHitsDigits(Int_t n=900, Double_t emin=1e-3, Double_t emax=10,
55 Int_t m=1100, Double_t amin=-0.5, Double_t amax=1099.5)
59 TArrayF eloss(MakeLogScale(n, emin, emax));
62 for (Int_t i = 1; i < m+1; i++) adcs[i] = adcs[i-1] + (amax-amin)/m;
63 fElossVsAdc = new TH2D("bad", "#Delta E vs. ADC",
64 eloss.fN-1, eloss.fArray, adcs.fN-1, adcs.fArray);
65 fElossVsAdc->SetXTitle("#Delta E/#Delta x [MeV/cm]");
66 fElossVsAdc->SetYTitle("ADC value");
68 //__________________________________________________________________
70 @param ev Event number
71 @return @c false on error */
72 Bool_t Begin(Int_t ev)
75 return AliFMDInput::Begin(ev);
77 //__________________________________________________________________
78 Bool_t ProcessHit(AliFMDHit* hit, TParticle*)
80 // Cache the energy loss
81 if (!hit) return kFALSE;
82 UShort_t det = hit->Detector();
83 Char_t rng = hit->Ring();
84 UShort_t sec = hit->Sector();
85 UShort_t str = hit->Strip();
87 AliWarning(Form("Bad strip number %d in hit", str));
90 fMap(det, rng, sec, str).fEdep += hit->Edep();
91 fMap(det, rng, sec, str).fN++;
94 //__________________________________________________________________
95 Bool_t ProcessDigit(AliFMDDigit* digit)
97 if (!digit) return kFALSE;
98 UShort_t det = digit->Detector();
99 Char_t rng = digit->Ring();
100 UShort_t sec = digit->Sector();
101 UShort_t str = digit->Strip();
103 AliWarning(Form("Bad strip number %d in digit", str));
106 fElossVsAdc->Fill(fMap(det, rng, sec, str).fEdep, digit->Counts());
109 //__________________________________________________________________
112 gStyle->SetPalette(1);
113 gStyle->SetOptTitle(0);
114 gStyle->SetCanvasColor(0);
115 gStyle->SetCanvasBorderSize(0);
116 gStyle->SetPadColor(0);
117 gStyle->SetPadBorderSize(0);
118 fElossVsAdc->SetStats(kFALSE);
119 fElossVsAdc->Draw("COLZ");
123 ClassDef(DrawHitsDigits,0);
126 //____________________________________________________________________