]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EVE/alice-macros/emcal_digits.C
Remove trailing whitespace.
[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 void emcal_digits()
10 {
11   AliRunLoader* rl =  AliEveEventManager::AssertRunLoader();
12
13   rl->LoadgAlice();
14   AliEMCAL         * emcal = (AliEMCAL*) rl->GetAliRun()->GetDetector("EMCAL");
15   AliEMCALGeometry * geom  = emcal->GetGeometry();
16
17   rl->LoadDigits("EMCAL");
18   TTree* dt = rl->GetTreeD("EMCAL", kFALSE);
19
20   gGeoManager = gEve->GetGeometry("$REVESYS/alice-data/alice_fullgeo.root");
21   TGeoNode* node = gGeoManager->GetTopVolume()->FindNode("XEN1_1");
22
23   TGeoBBox* bbbox = (TGeoBBox*) node->GetDaughter(0) ->GetVolume()->GetShape();
24   bbbox->Dump();
25   TGeoBBox* sbbox = (TGeoBBox*) node->GetDaughter(10)->GetVolume()->GetShape();
26   sbbox->Dump();
27
28   TEveElementList* l = new TEveElementList("EMCAL");
29   l->SetTitle("Tooltip");
30   gEve->AddElement(l);
31
32   TEveFrameBox* frame_big = new TEveFrameBox();
33   frame_big->SetAABoxCenterHalfSize(0, 0, 0, bbbox->GetDX(), bbbox->GetDY(), bbbox->GetDZ());
34
35   TEveFrameBox* frame_sml = new TEveFrameBox();
36   frame_sml->SetAABoxCenterHalfSize(0, 0, 0, sbbox->GetDX(), sbbox->GetDY(), sbbox->GetDZ());
37
38   gStyle->SetPalette(1, 0);
39   TEveRGBAPalette* pal = new TEveRGBAPalette(0, 512);
40   pal->SetLimits(0, 1024);
41
42   TEveQuadSet* smodules[12];
43
44   for (Int_t sm=0; sm<12; ++sm)
45   {
46     TEveQuadSet* q = new TEveQuadSet(Form("SM %d", sm+1));
47     q->SetOwnIds(kTRUE);
48     q->Reset(TEveQuadSet::kQT_RectangleYZFixedDimX, kFALSE, 32);
49     q->SetDefWidth (geom->GetPhiTileSize());
50     q->SetDefHeight(geom->GetEtaTileSize());
51
52     q->RefHMTrans().SetFrom(*node->GetDaughter(sm)->GetMatrix());
53
54     q->SetFrame(sm < 10 ? frame_big : frame_sml);
55     q->SetPalette(pal);
56
57     gEve->AddElement(q, l);
58     smodules[sm] = q;
59   }
60
61   TClonesArray *digits = 0;
62   dt->SetBranchAddress("EMCAL", &digits);
63   dt->GetEntry(0);
64   Int_t nEnt = digits->GetEntriesFast();
65   AliEMCALDigit * dig;
66
67   Int_t iEvent  = -1 ;
68   Float_t amp   = -1 ;
69   Float_t time  = -1 ;
70   Int_t id      = -1 ;
71   Int_t iSupMod =  0 ;
72   Int_t iTower  =  0 ;
73   Int_t iIphi   =  0 ;
74   Int_t iIeta   =  0 ;
75   Int_t iphi    =  0 ;
76   Int_t ieta    =  0 ;
77   Double_t x, y, z;
78
79   for (Int_t idig = 0; idig < nEnt; ++idig)
80   {
81     dig = static_cast<AliEMCALDigit *>(digits->At(idig));
82
83     if(dig != 0) {
84       id   = dig->GetId() ; //cell (digit) label
85       amp  = dig->GetAmp(); //amplitude in cell (digit)
86       time = dig->GetTime();//time of creation of digit after collision
87
88       cout<<"Cell ID "<<id<<" Amp "<<amp<<endl;//" time "<<time<<endl;
89
90       //Geometry methods
91       geom->GetCellIndex(id,iSupMod,iTower,iIphi,iIeta);
92       //Gives SuperModule and Tower numbers
93       geom->GetCellPhiEtaIndexInSModule(iSupMod,iTower,
94                                         iIphi, iIeta,iphi,ieta);
95       //Gives label of cell in eta-phi position per each supermodule
96
97       cout<< "SModule "<<iSupMod<<"; Tower "<<iTower
98           <<"; Eta "<<iIeta<<"; Phi "<<iIphi
99           <<"; Cell Eta "<<ieta<<"; Cell Phi "<<iphi<<endl;
100
101       geom->RelPosCellInSModule(id, x, y, z);
102       cout << x <<" "<< y <<" "<< z <<endl;
103
104       TEveQuadSet* q = smodules[iSupMod];
105       q->AddQuad(y, z);
106       q->QuadValue(amp);
107       q->QuadId(dig);
108     } else {
109       cout<<"Digit pointer 0x0"<<endl;
110     }
111   }
112
113   for (Int_t sm = 0; sm < 12; ++sm)
114   {
115     smodules[iSupMod]->RefitPlex();
116   }
117
118   gEve->Redraw3D();
119 }