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