]> git.uio.no Git - u/mrichter/AliRoot.git/blob - FMD/scripts/CheckAlign.C
add calculation and histograms for MC cross section
[u/mrichter/AliRoot.git] / FMD / scripts / CheckAlign.C
1 //____________________________________________________________________
2 //
3 // $Id$
4 //
5 // Script that contains a class to compare the raw data written to the
6 // digits it's created from.
7 //
8 // Use the script `Compile.C' to compile this class using ACLic. 
9 //
10 #include <AliLog.h>
11 #include <AliFMDHit.h>
12 #include <AliFMDDigit.h>
13 #include <AliFMDInput.h>
14 #include <AliFMDUShortMap.h>
15 #include <AliFMDParameters.h>
16 #include <AliFMDGeometry.h>
17 #include <AliFMDRing.h>
18 #include <AliFMDDetector.h>
19 #include <iostream>
20 #include <TH3D.h>
21 #include <TStyle.h>
22 #include <TArrayF.h>
23 #include <TParticle.h>
24 #include <TCanvas.h>
25
26 /** @class CheckAlign
27     @brief Check alignment 
28     @code 
29     Root> .L Compile.C
30     Root> Compile("CheckAlign.C")
31     Root> CheckAlign c
32     Root> c.Run();
33     @endcode
34     @ingroup FMD_script
35  */
36 class CheckAlign : public AliFMDInput
37 {
38 public:
39   CheckAlign()
40   {
41     AddLoad(kHits);
42     AddLoad(kDigits);
43     AddLoad(kGeometry);
44     AliFMDGeometry* geom = AliFMDGeometry::Instance();
45     geom->Init();
46     // geom->InitTransformations();
47     Double_t iR  = geom->GetRing('I')->GetHighR()+5;
48     Double_t oR  = geom->GetRing('O')->GetHighR()+5;
49     Double_t z1l = geom->GetDetector(1)->GetInnerZ()-5;
50     Double_t z1h = z1l+10;
51     Double_t z2l = geom->GetDetector(2)->GetOuterZ()-5;
52     Double_t z2h = geom->GetDetector(2)->GetInnerZ()+5;
53     Double_t z3l = geom->GetDetector(3)->GetOuterZ()-5;
54     Double_t z3h = geom->GetDetector(3)->GetInnerZ()+5;
55     
56     f1Hits   = new TH3D("hits1", "FMD1 hits", 
57                         300,-iR,iR,300,-iR,iR,100,z1l,z1h);
58     f1Hits->SetMarkerColor(2);
59     f1Hits->SetMarkerStyle(3);
60       
61     f2Hits   = new TH3D("hits2", "FMD2 hits", 
62                         300,-oR,oR,300,-oR,oR,100,z2l,z2h);
63     f2Hits->SetMarkerColor(2);
64     f2Hits->SetMarkerStyle(3);
65
66     f3Hits   = new TH3D("hits3", "FMD3 hits", 
67                         300,-oR,oR,300,-oR,oR,100,z3l,z3h);
68     f3Hits->SetMarkerColor(2);
69     f3Hits->SetMarkerStyle(3);
70
71     f1Digits = new TH3D("digits1", "FMD1 digits", 
72                         300,-iR,iR,300,-iR,iR,100,z1l,z1h); 
73     f1Digits->SetMarkerColor(4);
74     f1Digits->SetMarkerStyle(2);
75
76     f2Digits = new TH3D("digits2", "FMD2 digits", 
77                         300,-oR,oR,300,-oR,oR,100,z2l,z2h);
78     f2Digits->SetMarkerColor(4);
79     f2Digits->SetMarkerStyle(2);
80
81     f3Digits = new TH3D("digits3", "FMD3 hits", 
82                         300,-oR,oR,300,-oR,oR,100,z3l,z3h);
83     f3Digits->SetMarkerColor(4);
84     f3Digits->SetMarkerStyle(2);
85   }
86   Bool_t Init() 
87   {
88     Bool_t ret = AliFMDInput::Init();
89     AliFMDGeometry* geom = AliFMDGeometry::Instance();
90     geom->Init();
91     geom->InitTransformations();
92     AliFMDParameters* param = AliFMDParameters::Instance();
93     param->Init();
94     return ret;
95   }
96   
97   Bool_t ProcessHit(AliFMDHit* hit, TParticle*)
98   {
99     // Cache the energy loss 
100     if (!hit) return kFALSE;
101     UShort_t det = hit->Detector();
102     UShort_t str = hit->Strip();
103     if (str > 511) {
104       AliWarning(Form("Bad strip number %d in hit", str));
105       return kTRUE;
106     }
107     switch (det) {
108     case 1: f1Hits->Fill(hit->X(), hit->Y(), hit->Z()); break;
109     case 2: f2Hits->Fill(hit->X(), hit->Y(), hit->Z()); break;
110     case 3: f3Hits->Fill(hit->X(), hit->Y(), hit->Z()); break;
111     }
112     return kTRUE;
113   }
114   Bool_t ProcessDigit(AliFMDDigit* digit)
115   {
116     // Cache the energy loss 
117     if (!digit) return kFALSE;
118     UShort_t det = digit->Detector();
119     Char_t   rng = digit->Ring();
120     UShort_t sec = digit->Sector();
121     UShort_t str = digit->Strip();
122     if (str > 511) {
123       AliWarning(Form("Bad strip number %d in digit", str));
124       return kTRUE;
125     }
126     AliFMDParameters* param = AliFMDParameters::Instance();
127     if (digit->Counts() < (param->GetPedestal(det, rng, sec, str) + 4 *
128                            param->GetPedestalWidth(det, rng, sec, str)))
129       return kTRUE;
130     AliFMDGeometry* geom = AliFMDGeometry::Instance();
131     Double_t x, y, z;
132     geom->Detector2XYZ(det, rng, sec, str, x, y, z);
133     switch (det) {
134     case 1: f1Digits->Fill(x, y , z); break;
135     case 2: f2Digits->Fill(x, y , z); break;
136     case 3: f3Digits->Fill(x, y , z); break;
137     }
138     return kTRUE;
139   }
140   //__________________________________________________________________
141   Bool_t Finish()
142   {
143     gStyle->SetPalette(1);
144     gStyle->SetOptTitle(0);
145     gStyle->SetCanvasColor(0);
146     gStyle->SetCanvasBorderSize(0);
147     gStyle->SetPadColor(0);
148     gStyle->SetPadBorderSize(0);
149
150     TCanvas* c1 = new TCanvas("FMD1","FMD1");
151     c1->cd();
152     f1Hits->Draw();
153     f1Digits->Draw("same");
154
155     TCanvas* c2 = new TCanvas("FMD2","FMD2");
156     c2->cd();
157     f2Hits->Draw();
158     f2Digits->Draw("same");
159
160     TCanvas* c3 = new TCanvas("FMD3","FMD3");
161     c3->cd();
162     f3Hits->Draw();
163     f3Digits->Draw("same");
164     return kTRUE;
165   }
166 protected:
167   TH3D* f1Hits;
168   TH3D* f2Hits;
169   TH3D* f3Hits;
170   TH3D* f1Digits;
171   TH3D* f2Digits;
172   TH3D* f3Digits;
173 };
174
175
176 //____________________________________________________________________
177 //
178 // EOF
179 //
180
181
182   
183