]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EVE/Alieve/AliEvePMDModule.cxx
Put black-listed classes out of Alieve namespace.
[u/mrichter/AliRoot.git] / EVE / Alieve / AliEvePMDModule.cxx
1 // $Id$
2 // Main authors: Matevz Tadel & Alja Mrak-Tadel: 2006, 2007
3
4 /**************************************************************************
5  * Copyright(c) 1998-2008, ALICE Experiment at CERN, all rights reserved. *
6  * See http://aliceinfo.cern.ch/Offline/AliRoot/License.html for          *
7  * full copyright notice.                                                 * 
8  **************************************************************************/
9
10 #include "AliEvePMDModule.h"
11
12 #include "AliPMDdigit.h"
13 #include "AliPMDddldata.h"
14
15 #include <TClonesArray.h>
16
17
18 const Float_t AliEvePMDModule::fgkRad        = 0.25;
19 const Float_t AliEvePMDModule::fgkSqRoot3    = 1.732050808;
20 const Float_t AliEvePMDModule::fgkZpos       = 0.;
21 Int_t         AliEvePMDModule::fPreTotPads   = 0;
22 Int_t         AliEvePMDModule::fCpvTotPads   = 0;
23 Int_t         AliEvePMDModule::fPreTotAdc    = 0;
24 Int_t         AliEvePMDModule::fCpvTotAdc    = 0;
25
26
27 //______________________________________________________________________
28 // AliEvePMDModule
29 //
30
31
32 ClassImp(AliEvePMDModule)
33
34 AliEvePMDModule::AliEvePMDModule():
35   fH1(0),
36   fX(0.),
37   fY(0.),
38   fZ(0.),
39   fNPads(0),
40   fAdc(0)
41 {
42
43 }
44 // -------------------------------------------------------------------- //
45 void AliEvePMDModule::DisplayInit(Int_t ism)
46 {
47   TString smodule = "Module";
48   smodule+= ism;
49   
50   Float_t xism =0, yism = 0;
51   Float_t dxism =0, dyism = 0;
52
53   GenerateBox(ism,xism,yism,dxism,dyism);
54
55   TEveFrameBox *pmdModBox = new TEveFrameBox();
56   pmdModBox->SetAAQuadXY(xism, yism, 0, dxism, dyism);
57   pmdModBox->SetFrameColor((Color_t) 31);
58   pmdModBox->SetFrameFill(kTRUE);
59   SetFrame(pmdModBox);
60
61   SetName(smodule.Data());
62   SetOwnIds(kTRUE);
63   Reset(TEveQuadSet::kQT_HexagonXY, kFALSE, 32);
64
65   fH1 = new TH1F("fH1", smodule.Data(), 100, 0., 1000.);
66   fH1->SetDirectory(0);
67   fH1->GetXaxis()->SetTitle("Single Cell Edep (adc)");
68 }
69
70 // -------------------------------------------------------------------- //
71
72 void AliEvePMDModule::DisplayDigitsData(Int_t ism, TTree *pmdt)
73 {
74   DisplayInit(ism);
75
76   Int_t det, smn, irow, icol, adc;
77   Int_t xpad = 0, ypad = 0;
78   Float_t xpos, ypos;
79
80   TClonesArray *digits = new TClonesArray("AliPMDdigit", 0); 
81
82   TBranch *branch = pmdt->GetBranch("PMDDigit");
83   branch->SetAddress(&digits);
84
85   AliPMDdigit  *pmddigit;
86
87   branch->GetEntry(ism); 
88   Int_t nentries = digits->GetLast();
89   //printf("%d\n", nentries);
90
91   for (Int_t ient = 0; ient < nentries+1; ient++)
92     {
93       pmddigit = (AliPMDdigit*)digits->UncheckedAt(ient);
94       
95       det    = pmddigit->GetDetector();
96       smn    = pmddigit->GetSMNumber();
97       irow   = pmddigit->GetRow();
98       icol   = pmddigit->GetColumn();
99       adc    = (Int_t) pmddigit->GetADC();
100
101       if(smn <12)
102         {
103           xpad = icol;
104           ypad = irow;
105         }
106       else if(smn >=12 && smn < 24)
107         {
108           xpad = irow;
109           ypad = icol;
110         }
111               
112       RectGeomCellPos(smn, xpad, ypad, xpos, ypos);
113
114       AddHexagon(xpos, ypos, fgkZpos, fgkRad);
115       
116       QuadValue(adc);
117           
118       QuadId(new AliPMDdigit(*pmddigit));
119       // new TNamed(Form("Quad with idx=%d", ient),
120       //  "This title is not confusing."));
121
122       ++fNPads;
123       fAdc += adc;
124
125       if (det == 0)
126         {
127           fPreTotAdc += (Int_t) adc;
128           ++fPreTotPads;
129         }
130       if (det == 1)
131         {
132           fCpvTotAdc += (Int_t) adc;
133           ++fCpvTotPads;
134         }
135       fH1->Fill((Float_t)adc);
136
137     }
138
139   RefitPlex();
140
141   fHMTrans.SetPos(fX, fY, fZ);
142
143   delete digits;
144 }
145
146 // -------------------------------------------------------------------- //
147
148 void AliEvePMDModule::DisplayRawData(Int_t ism, TObjArray *ddlcont)
149 {
150   DisplayInit(ism);
151
152   if (ism > 23) ism -= 24;
153
154   Int_t det, smn, irow, icol, adc;
155   Int_t xpad = 0, ypad = 0;
156   Float_t xpos, ypos;
157
158   Int_t nentries = ddlcont->GetEntries();
159   //printf("%d\n", nentries);
160
161   for (Int_t ient = 0; ient < nentries; ient++)
162     {
163       AliPMDddldata *pmdddl = (AliPMDddldata*)ddlcont->UncheckedAt(ient);
164       
165       det    = pmdddl->GetDetector();
166       smn    = pmdddl->GetSMN();
167       if (smn != ism) continue;
168       irow   = pmdddl->GetRow();
169       icol   = pmdddl->GetColumn();
170       adc    = pmdddl->GetSignal();
171
172       if(smn <12)
173         {
174           xpad = icol;
175           ypad = irow;
176         }
177       else if(smn >=12 && smn < 24)
178         {
179           xpad = irow;
180           ypad = icol;
181         }
182               
183       RectGeomCellPos(smn, xpad, ypad, xpos, ypos);
184
185       AddHexagon(xpos, ypos, fgkZpos, fgkRad);
186       
187       QuadValue(adc);
188           
189       QuadId(new AliPMDddldata(*pmdddl));
190       //new TNamed(Form("Quad with idx=%d", ient),
191       //  "This title is not confusing."));
192
193       ++fNPads;
194       fAdc += adc;
195
196       if (det == 0)
197         {
198           fPreTotAdc += (Int_t) adc;
199           ++fPreTotPads;
200         }
201       if (det == 1)
202         {
203           fCpvTotAdc += (Int_t) adc;
204           ++fCpvTotPads;
205         }
206
207       fH1->Fill((Float_t) adc);
208     }
209
210   RefitPlex();
211
212   fHMTrans.SetPos(fX, fY, fZ);
213   
214 }
215
216 // -------------------------------------------------------------------- //
217
218 void AliEvePMDModule::RectGeomCellPos(Int_t ism, Int_t xpad, Int_t ypad,
219                                 Float_t &xpos, Float_t &ypos)
220 {
221   // This routine finds the cell eta,phi for the new PMD rectangular 
222   // geometry in ALICE
223   // Authors : Bedanga Mohanty and Dipak Mishra - 29.4.2003
224   // modified by B. K. Nandi for change of coordinate sys
225   //
226   // SMA  ---> Supermodule Type A           ( SM - 0)
227   // SMAR ---> Supermodule Type A ROTATED   ( SM - 1)
228   // SMB  ---> Supermodule Type B           ( SM - 2)
229   // SMBR ---> Supermodule Type B ROTATED   ( SM - 3)
230   //
231   // ism   : Serial module number from 0 to 23 for each plane
232
233  
234   // Corner positions (x,y) of the 24 unit moudles in ALICE PMD
235
236   const Double_t kXcorner[24] =
237     {
238       74.8833,  53.0045, 31.1255,    //Type-A
239       74.8833,  53.0045, 31.1255,    //Type-A
240       -74.8833, -53.0044, -31.1255,  //Type-AR
241       -74.8833, -53.0044, -31.1255,  //Type-AR
242       8.9165, -33.7471,            //Type-B
243       8.9165, -33.7471,            //Type-B
244       8.9165, -33.7471,            //Type-B
245       -8.9165, 33.7471,            //Type-BR
246       -8.9165, 33.7471,            //Type-BR
247       -8.9165, 33.7471,            //Type-BR
248     };
249
250   
251   const Double_t kYcorner[24] =
252     {
253       86.225,  86.225,  86.225,      //Type-A
254       37.075,  37.075,  37.075,      //Type-A
255       -86.225, -86.225, -86.225,     //Type-AR
256       -37.075, -37.075, -37.075,     //Type-AR
257       86.225,  86.225,               //Type-B
258       61.075,  61.075,               //Type-B
259       35.925,  35.925,               //Type-B
260       -86.225, -86.225,              //Type-BR
261       -61.075, -61.075,              //Type-BR
262       -35.925, -35.925               //Type-BR
263     };
264
265   
266   //  const Float_t kSqroot3      = 1.732050808;  // sqrt(3.);
267   //  const Float_t kCellRadius   = 0.25;
268   
269   //
270   //Every even row of cells is shifted and placed
271   //in geant so this condition
272   //
273   Float_t shift = 0.0;
274   if(ypad%2 == 0)
275     {
276       shift = 0.0;
277     }
278   else
279     {
280       shift = 0.25;
281     }
282
283
284   if(ism < 6)
285     {
286       ypos = kYcorner[ism] - (Float_t) xpad*fgkRad*2.0 - shift;
287       xpos = kXcorner[ism] - (Float_t) ypad*fgkSqRoot3*fgkRad;
288     }
289   else if(ism >=6 && ism < 12)
290     {
291       ypos = kYcorner[ism] + (Float_t) xpad*fgkRad*2.0 + shift;
292       xpos = kXcorner[ism] + (Float_t) ypad*fgkSqRoot3*fgkRad;
293     }
294   else if(ism >= 12 && ism < 18)
295     {
296       ypos = kYcorner[ism] - (Float_t) xpad*fgkRad*2.0 - shift;
297       xpos = kXcorner[ism] - (Float_t) ypad*fgkSqRoot3*fgkRad;
298     }
299   else if(ism >= 18 && ism < 24)
300     {
301       ypos = kYcorner[ism] + (Float_t) xpad*fgkRad*2.0 + shift;
302       xpos = kXcorner[ism] + (Float_t) ypad*fgkSqRoot3*fgkRad;
303     }
304
305 }
306
307 // -------------------------------------------------------------------- //
308
309 void AliEvePMDModule::GenerateBox(Int_t ism, Float_t &xism, Float_t &yism,
310                             Float_t &dxism, Float_t &dyism)
311 {
312   const Float_t kDia     = 0.50;
313
314   const Double_t kXcorner[24] =
315     {
316       74.8833,  53.0045, 31.1255,    //Type-A
317       74.8833,  53.0045, 31.1255,    //Type-A
318       -74.8833, -53.0044, -31.1255,  //Type-AR
319       -74.8833, -53.0044, -31.1255,  //Type-AR
320       8.9165, -33.7471,            //Type-B
321       8.9165, -33.7471,            //Type-B
322       8.9165, -33.7471,            //Type-B
323       -8.9165, 33.7471,            //Type-BR
324       -8.9165, 33.7471,            //Type-BR
325       -8.9165, 33.7471,            //Type-BR
326     };
327
328   
329   const Double_t kYcorner[24] =
330     {
331       86.225,  86.225,  86.225,      //Type-A
332       37.075,  37.075,  37.075,      //Type-A
333       -86.225, -86.225, -86.225,     //Type-AR
334       -37.075, -37.075, -37.075,     //Type-AR
335       86.225,  86.225,               //Type-B
336       61.075,  61.075,               //Type-B
337       35.925,  35.925,               //Type-B
338       -86.225, -86.225,              //Type-BR
339       -61.075, -61.075,              //Type-BR
340       -35.925, -35.925               //Type-BR
341     };
342
343
344   if (ism > 23) ism -= 24;
345
346
347   if (ism < 6)
348     {
349       xism  = kXcorner[ism] + fgkRad;
350       yism  = kYcorner[ism] + fgkRad;
351       dxism = -fgkRad*fgkSqRoot3*48.;
352       dyism = -kDia*96. - fgkRad;
353   }
354   if (ism >= 6 && ism < 12)
355     {
356       xism  = kXcorner[ism] - fgkRad;
357       yism  = kYcorner[ism] - fgkRad;
358       dxism = fgkRad*fgkSqRoot3*48.;
359       dyism = kDia*96. + fgkRad;
360   }
361   if (ism >= 12 && ism < 18)
362     {
363       xism  = kXcorner[ism] + fgkRad;
364       yism  = kYcorner[ism] + fgkRad;
365       dxism = -fgkRad*fgkSqRoot3*96.;
366       dyism = -kDia*48. - fgkRad;
367   }
368   if (ism >= 18 && ism < 24)
369     {
370       xism  = kXcorner[ism] - fgkRad;
371       yism  = kYcorner[ism] - fgkRad;
372       dxism = fgkRad*fgkSqRoot3*96.;
373       dyism = kDia*48. + fgkRad;
374   }
375   
376 }
377
378 // -------------------------------------------------------------------- //
379
380 void AliEvePMDModule::SetPosition(Float_t x, Float_t y, Float_t z)
381 {
382   fX = x;
383   fY = y;
384   fZ = z;
385 }
386
387 // -------------------------------------------------------------------- //