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