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