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