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