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