]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EVE/alice-macros/emcal_digits.C
Increase transparency
[u/mrichter/AliRoot.git] / EVE / alice-macros / emcal_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 #ifndef __CINT__
10
11 #include <TEveManager.h>
12 #include <TEveQuadSet.h>
13 #include <TGeoNode.h>
14 #include <TGeoBBox.h>
15 #include <TGeoManager.h>
16 #include <TStyle.h>
17 #include <TEveTrans.h>
18 #include <TClonesArray.h>
19
20 #include <EveBase/AliEveEventManager.h>
21
22 #include <AliRunLoader.h>
23 #include <AliCluster.h>
24 #include <EMCAL/AliEMCALGeometry.h>
25 #include <EMCAL/AliEMCALDigit.h>
26
27 #include <Riostream.h>
28 #endif
29
30 void emcal_digits()
31 {
32   AliEveEventManager::AssertGeometry();
33
34   TGeoNode* node = gGeoManager->GetTopVolume()->FindNode("XEN1_1");
35   if (!node) return;
36
37   Int_t nModules = node->GetNdaughters();
38
39   TEveElementList* l = new TEveElementList("EMCAL");
40   l->SetTitle("Tooltip");
41   gEve->AddElement(l);
42
43   TGeoBBox* bbbox = (TGeoBBox*) node->GetDaughter(0) ->GetVolume()->GetShape();
44   TEveFrameBox* frame_big = new TEveFrameBox();
45   frame_big->SetFrameColorRGBA(200,200,0,50);
46   frame_big->SetAABoxCenterHalfSize(0, 0, 0, bbbox->GetDX(), bbbox->GetDY(), bbbox->GetDZ());
47
48   TEveFrameBox* frame_sml = 0x0;
49
50   if (nModules==12) {
51     TGeoBBox* sbbox = (TGeoBBox*) node->GetDaughter(10)->GetVolume()->GetShape();
52     frame_sml = new TEveFrameBox();
53     frame_sml->SetFrameColorRGBA(200,200,0,50);
54     frame_sml->SetAABoxCenterHalfSize(0, 0, 0, sbbox->GetDX(), sbbox->GetDY(), sbbox->GetDZ());
55   }
56
57   gStyle->SetPalette(1, 0);
58   TEveRGBAPalette* pal = new TEveRGBAPalette(0, 512);
59   pal->SetLimits(0, 1024);
60
61   TEveQuadSet* smodules[12];
62   memset(smodules,0,12*sizeof(TEveQuadSet*));
63
64
65   AliEMCALGeometry * geom  = AliEMCALGeometry::GetInstance();  
66   if (!geom) geom = AliEMCALGeometry::GetInstance("","");
67
68   for (Int_t sm=0; sm<nModules; ++sm)
69   {
70     TEveQuadSet* q = new TEveQuadSet(Form("SM %d", sm+1));
71     q->SetOwnIds(kTRUE);
72     q->Reset(TEveQuadSet::kQT_RectangleYZFixedDimX, kFALSE, 32);
73     q->SetDefWidth (geom->GetPhiTileSize());
74     q->SetDefHeight(geom->GetEtaTileSize());
75
76     q->RefMainTrans().SetFrom(*node->GetDaughter(sm)->GetMatrix());
77
78     q->SetFrame(sm < 10 ? frame_big : frame_sml);
79     q->SetPalette(pal);
80
81     gEve->AddElement(q, l);
82     smodules[sm] = q;
83   }
84
85   AliRunLoader* rl =  AliEveEventManager::AssertRunLoader();
86
87   rl->LoadDigits("EMCAL");
88   TTree* dt = rl->GetTreeD("EMCAL", kFALSE);
89   if (!dt) return;
90
91   TClonesArray *digits = 0;
92   dt->SetBranchAddress("EMCAL", &digits);
93   dt->GetEntry(0);
94   Int_t nEnt = digits->GetEntriesFast();
95   AliEMCALDigit * dig;
96
97   Float_t amp   = -1 ;
98   Float_t time  = -1 ;
99   Int_t id      = -1 ;
100   Int_t iSupMod =  0 ;
101   Int_t iTower  =  0 ;
102   Int_t iIphi   =  0 ;
103   Int_t iIeta   =  0 ;
104   Int_t iphi    =  0 ;
105   Int_t ieta    =  0 ;
106   Double_t x, y, z;
107
108   for (Int_t idig = 0; idig < nEnt; ++idig)
109   {
110     dig = static_cast<AliEMCALDigit *>(digits->At(idig));
111
112     if(dig != 0) {
113       id   = dig->GetId() ; //cell (digit) label
114       amp  = dig->GetAmp(); //amplitude in cell (digit)
115       time = dig->GetTime();//time of creation of digit after collision
116
117       cout<<"Cell ID "<<id<<" Amp "<<amp<<endl;//" time "<<time<<endl;
118
119       //Geometry methods
120       geom->GetCellIndex(id,iSupMod,iTower,iIphi,iIeta);
121       //Gives SuperModule and Tower numbers
122       geom->GetCellPhiEtaIndexInSModule(iSupMod,iTower,
123                                         iIphi, iIeta,iphi,ieta);
124       //Gives label of cell in eta-phi position per each supermodule
125
126       cout<< "SModule "<<iSupMod<<"; Tower "<<iTower
127           <<"; Eta "<<iIeta<<"; Phi "<<iIphi
128           <<"; Cell Eta "<<ieta<<"; Cell Phi "<<iphi<<endl;
129
130       geom->RelPosCellInSModule(id, x, y, z);
131       cout << x <<" "<< y <<" "<< z <<endl;
132
133       TEveQuadSet* q = smodules[iSupMod];
134       if (q) {
135         q->AddQuad(y, z);
136         q->QuadValue(TMath::Nint(amp));
137         q->QuadId(new AliEMCALDigit(*dig));
138       }
139     } else {
140       cout<<"Digit pointer 0x0"<<endl;
141     }
142   }
143
144   rl->UnloadDigits("EMCAL");
145
146
147   for (Int_t sm = 0; sm < nModules; ++sm)
148   {
149     smodules[iSupMod]->RefitPlex();
150   }
151
152   gEve->Redraw3D();
153 }