]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EVE/alice-macros/vzero_digits.C
Put all classes in Alieve/ sub-module out of the Alieve namespace and
[u/mrichter/AliRoot.git] / EVE / alice-macros / vzero_digits.C
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 //   fV0CHeight1         =    2.5; // height of cell 1
10 //   fV0CHeight2         =    4.4; // height of cell 2
11 //   fV0CHeight3         =    7.4; // height of cell 3
12 //   fV0CHeight4         =   12.5; // height of cell 4
13 //   fV0CRMin            =    4.6; // inner radius of box
14
15 //   fV0AR0     =  4.2;  // Radius of hole
16 //   fV0AR1     =  7.6;  // Maximun radius of 1st cell
17 //   fV0AR2     = 13.8; // Maximun radius of 2nd cell
18 //   fV0AR3     = 22.7; // Maximun radius of 3rd cell
19 //   fV0AR4     = 41.3; // Maximun radius of 4th cell
20
21 void vzero_digits()
22 {
23   static const Float_t RadC[] = { 4.6, 7.1, 11.5, 18.9, 31.4 };
24   static const Float_t RadA[] = { 4.2, 7.6, 13.8, 22.7, 41.4 };
25   static const Float_t RadEps = 0.4;
26   static const Float_t PhiEps = 0.025;
27   static const Float_t PhiStp = TMath::TwoPi()/8.0;
28
29   gStyle->SetPalette(1, 0);
30
31   AliRunLoader* rl =  AliEveEventManager::AssertRunLoader();
32   rl->LoadDigits("VZERO");
33
34   TTree* dt = rl->GetTreeD("VZERO", false);
35   TClonesArray* dca = 0;
36   dt->SetBranchAddress("VZERODigit", &dca);
37   dt->GetEntry(0);
38
39   Float_t v[12];
40
41   TEveQuadSet* modR = new TEveQuadSet("V0R");
42   modR->Reset(TEveQuadSet::kQT_FreeQuad, kFALSE, 32);
43
44   TEveQuadSet* modL = new TEveQuadSet("V0L");
45   modL->Reset(TEveQuadSet::kQT_FreeQuad, kFALSE, 32);
46  
47   Int_t numEntr = dca->GetEntriesFast();
48   for (Int_t entr=0; entr<numEntr; ++entr)
49   {
50     AliVZEROdigit* d = (AliVZEROdigit*) dca->UncheckedAt(entr);
51     Int_t i = d->PMNumber();
52
53     if (i < 32)   // AliEveV0 Right
54     {
55       TEveQuadSet* module = modR;
56       Int_t ri = i / 8;
57       Int_t pi = i % 8;
58       Float_t minR = RadC[ri]  + RadEps, maxR = RadC[ri+1] - RadEps;
59       Float_t minP = pi*PhiStp + PhiEps, maxP = (pi+1)*PhiStp - PhiEps;
60
61       v[ 0] = minR*TMath::Cos(minP); v[ 1] = minR*TMath::Sin(minP); v[ 2] = 0;
62       v[ 3] = maxR*TMath::Cos(minP); v[ 4] = maxR*TMath::Sin(minP); v[ 5] = 0;
63       v[ 6] = maxR*TMath::Cos(maxP); v[ 7] = maxR*TMath::Sin(maxP); v[ 8] = 0;
64       v[ 9] = minR*TMath::Cos(maxP); v[10] = minR*TMath::Sin(maxP); v[11] = 0;
65
66       module->AddQuad(v);
67       module->QuadValue(d->ADC());
68       module->QuadId(d);
69     }
70     else          // AliEveV0 Left
71     {
72       TEveQuadSet* module = modL;
73       Int_t ri = (i-32) / 8;
74       Int_t pi = i % 8;
75       Float_t minR = RadA[ri]  + RadEps, maxR = RadA[ri+1] - RadEps;
76       Float_t minP = pi*PhiStp + PhiEps, maxP = (pi+1)*PhiStp - PhiEps;
77
78       v[ 0] = minR*TMath::Cos(minP); v[ 1] = minR*TMath::Sin(minP); v[ 2] = 0;
79       v[ 3] = maxR*TMath::Cos(minP); v[ 4] = maxR*TMath::Sin(minP); v[ 5] = 0;
80       v[ 6] = maxR*TMath::Cos(maxP); v[ 7] = maxR*TMath::Sin(maxP); v[ 8] = 0;
81       v[ 9] = minR*TMath::Cos(maxP); v[10] = minR*TMath::Sin(maxP); v[11] = 0;
82
83       module->AddQuad(v);
84       module->QuadValue(d->ADC());
85       module->QuadId(d);
86     }
87   }
88
89   modL->RefHMTrans().SetPos(0, 0, 324);
90   modR->RefHMTrans().SetPos(0, 0, -84);
91
92   gEve->AddElement(modL);
93   gEve->AddElement(modR);
94
95   gEve->Redraw3D();
96 }