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