AliRoot Core  edcc906 (edcc906)
CheckAlign.C
Go to the documentation of this file.
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 
36 class CheckAlign : public AliFMDInput
37 {
38 public:
40  {
41  AddLoad(kHits);
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();
90  geom->Init();
91  geom->InitTransformations();
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  }
127  if (digit->Counts() < (param->GetPedestal(det, rng, sec, str) + 4 *
128  param->GetPedestalWidth(det, rng, sec, str)))
129  return kTRUE;
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 
class for digits
Definition: AliFMDDigit.h:28
Double_t GetOuterZ() const
Geometry mananger for the FMD.
Bool_t ProcessHit(AliFMDHit *hit, TParticle *)
Definition: CheckAlign.C:97
TStyle * gStyle
AliFMDDetector * GetDetector(Int_t i) const
Hit in the FMD.
Double_t GetHighR() const
Definition: AliFMDRing.h:194
Double_t GetInnerZ() const
virtual void InitTransformations(Bool_t force=kFALSE)
Digits for the FMD.
FMD utility classes for reading FMD data.
UShort_t Sector() const
Char_t Ring() const
This class is a singleton that handles various parameters of the FMD detectors. This class reads from...
Manager of FMD parameters.
Bool_t ProcessDigit(AliFMDDigit *digit)
Definition: CheckAlign.C:114
FMD ring geometry parameters.
Check alignment.
Definition: CheckAlign.C:36
Float_t Z() const
Definition: AliHit.h:23
Per strip of unisgned shorts (16 bit) data.
#define AliWarning(message)
Definition: AliLog.h:541
Float_t X() const
Definition: AliHit.h:21
Float_t Y() const
Definition: AliHit.h:22
virtual Bool_t Init()
Bool_t Init()
Definition: CheckAlign.C:86
TH3D * f3Hits
Definition: CheckAlign.C:169
TH3D * f2Hits
Definition: CheckAlign.C:168
Singleton object of FMD geometry descriptions and parameters. This class is a singleton that handles ...
UShort_t Detector() const
Definition: AliFMDHit.h:74
UShort_t Strip() const
TH3D * f3Digits
Definition: CheckAlign.C:172
virtual void AddLoad(ETrees tree)
Definition: AliFMDInput.h:134
TH3D * f1Digits
Definition: CheckAlign.C:170
static AliFMDParameters * Instance()
void Detector2XYZ(UShort_t detector, Char_t ring, UShort_t sector, UShort_t strip, Double_t &x, Double_t &y, Double_t &z) const
virtual void Init()
Base class for reading in various FMD data. The class loops over all found events. For each event the specified data is read in. The class then loops over all elements of the read data, and process these with user defined code.
Definition: AliFMDInput.h:106
UShort_t Strip() const
Definition: AliFMDHit.h:80
AliFMDhit is the hit class for the FMD. Hits are the information that comes from a Monte Carlo at eac...
Definition: AliFMDHit.h:30
TH3D * f2Digits
Definition: CheckAlign.C:171
UShort_t Detector() const
AliFMDRing * GetRing(Char_t i) const
Bool_t Finish()
Definition: CheckAlign.C:141
UShort_t Counts() const
Definition: AliFMDDigit.h:147
TH3D * f1Hits
Definition: CheckAlign.C:167
TEveGeoShape * geom
Definition: tpc_tracks.C:10
static AliFMDGeometry * Instance()
Float_t GetPedestal(UShort_t detector, Char_t ring, UShort_t sector, UShort_t strip) const
Float_t GetPedestalWidth(UShort_t detector, Char_t ring, UShort_t sector, UShort_t strip) const
Sub-detector base class declaration.
UShort_t Init(Bool_t forceReInit=kFALSE, UInt_t what=kAll)