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