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