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