]> git.uio.no Git - u/mrichter/AliRoot.git/blame - FMD/scripts/DrawHitsSDigits.C
bug fix
[u/mrichter/AliRoot.git] / FMD / scripts / DrawHitsSDigits.C
CommitLineData
faf80567 1//____________________________________________________________________
2//
3// $Id: DrawHitsSDigits.C 20907 2007-09-25 08:44:03Z cholm $
4//
5// Script that contains a class to draw eloss from hits, versus ADC
6// counts from sdigits, 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 <AliFMDSDigit.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 DrawHitsSDigits
26 @brief Draw hit energy loss versus sdigit ADC
27 @code
28 Root> .L Compile.C
29 Root> Compile("DrawHitsSDigits.C")
30 Root> DrawHitsSDigits c
31 Root> c.Run();
32 @endcode
33 @ingroup FMD_script
34 */
35class DrawHitsSDigits : public AliFMDInput
36{
37private:
38 TH2D* fElossVsAdc; // Histogram
39 AliFMDEdepMap fMap;
40public:
41 //__________________________________________________________________
42 DrawHitsSDigits(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(kSDigits);
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 Info("ProcessHit", "Processing hit %p", hit);
69 // Cache the energy loss
70 if (!hit) return kFALSE;
71 UShort_t det = hit->Detector();
72 Char_t rng = hit->Ring();
73 UShort_t sec = hit->Sector();
74 UShort_t str = hit->Strip();
75 if (str > 511) {
76 AliWarning(Form("Bad strip number %d in hit", str));
77 return kTRUE;
78 }
79 fMap(det, rng, sec, str).fEdep += hit->Edep();
80 fMap(det, rng, sec, str).fN++;
81 hit->Print();
82 return kTRUE;
83 }
84 //__________________________________________________________________
85 Bool_t ProcessSDigit(AliFMDSDigit* sdigit)
86 {
87 Info("ProcessSDigit", "Processing sdigit %p", sdigit);
88 if (!sdigit) return kFALSE;
89 UShort_t det = sdigit->Detector();
90 Char_t rng = sdigit->Ring();
91 UShort_t sec = sdigit->Sector();
92 UShort_t str = sdigit->Strip();
93 if (str > 511) {
94 AliWarning(Form("Bad strip number %d in sdigit", str));
95 return kFALSE;
96 }
97 fElossVsAdc->Fill(fMap(det, rng, sec, str).fEdep, sdigit->Counts());
98 sdigit->Print();
99 return kTRUE;
100 }
101 //__________________________________________________________________
102 Bool_t Finish()
103 {
104 gStyle->SetPalette(1);
105 gStyle->SetOptTitle(0);
106 gStyle->SetCanvasColor(0);
107 gStyle->SetCanvasBorderSize(0);
108 gStyle->SetPadColor(0);
109 gStyle->SetPadBorderSize(0);
110 fElossVsAdc->SetStats(kFALSE);
111 fElossVsAdc->Draw("COLZ");
112 return kTRUE;
113 }
114
115 ClassDef(DrawHitsSDigits,0);
116};
117
118//____________________________________________________________________
119//
120// EOF
121//