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