]> git.uio.no Git - u/mrichter/AliRoot.git/blame - EVE/alice-macros/muon_hits.C
Moved point-count variable declaration within the track loop.
[u/mrichter/AliRoot.git] / EVE / alice-macros / muon_hits.C
CommitLineData
34ea2f0a 1// $Id$
2
3//void
4//muon_hits(const char *varexp = "fX:fY:fZ",
5// const char *selection = "",
6// Option_t *option = "goff")
7void muon_hits(Int_t iShowCha = 1)
8{
9 char name[128];
10 char title[128];
11
12 AliRunLoader* rl = Alieve::Event::AssertRunLoader();
13 rl->LoadHits("MUON");
14
15 TTree* ht = rl->GetTreeH("MUON", false);
16 //ht->Draw(varexp, selection, option);
17 //ht->Print();
18
19 for (Int_t iSta = 1; iSta <= 7; iSta++) {
20
21 for (Int_t iCha = 1; iCha <= 2; iCha++) {
22
23 Int_t iChamber = (iSta-1) * 2 + iCha;
24
25 if (iChamber != iShowCha) continue;
26
27 sprintf(name,"M-ST%1dCH%1d/Hits",iSta,iCha);
28 Reve::RenderElementList* l = new Reve::RenderElementList(name);
29 if (iSta <= 5) {
30 sprintf(title,"Station %1d chamber %1d (tracking)",iSta,iCha);
31 } else {
32 sprintf(title,"Station %1d chamber %1d (trigger)",iSta,iCha);
33 }
34
35 l->SetTitle(title);
36 l->SetMainColor((Color_t)4);
5b96ea20 37 gReve->AddRenderElement(l);
34ea2f0a 38
39 AliMpDEIterator ite;
40 for ( ite.First(iChamber-1); ! ite.IsDone(); ite.Next() ) {
41
42 Int_t detElemId = ite.CurrentDE();
43
44 Int_t nTracks = ht->GetEntries();
45 printf("Found %d tracks. \n",nTracks);
46 for (Int_t it = 0; it < nTracks; it++) {
47
48 TClonesArray *hits = 0;
49 ht->SetBranchAddress("MUONHits",&hits);
50
51 ht->GetEvent(it);
52
53 Int_t nHits = hits->GetEntriesFast();
54 //printf("Found %d hits in track %d. \n",nHits,it);
55
56 Int_t nstep = 10;
57 Reve::PointSet* points = new Reve::PointSet(nHits*nstep);
58 sprintf(name,"Track %d: hits",it);
59 points->SetName(name);
60 sprintf(title,"TreeH entry %d",it);
61 points->SetTitle(title);
62 Int_t np = 0;
63
64 AliMUONHit *hit;
65 Float_t x, y, z;
66 for (Int_t ih = 0; ih < nHits; ih++) {
67
68 hit = (AliMUONHit*)hits->UncheckedAt(ih);
69 if (detElemId != hit->DetElemId()) continue;
70
71 x = hit->X();
72 y = hit->Y();
73 z = hit->Z();
74
75 Float_t dstep = 2.0*TMath::Pi() / (Float_t)nstep;
76 Float_t d;
77 Float_t r = 1.0;
78 for (Int_t istep = 1; istep < (nstep+1); istep++) {
79
80 d = istep * dstep;
81 Float_t xc = x + r * TMath::Cos(d);
82 Float_t yc = y + r * TMath::Sin(d);
83
84 points->SetPoint(np,z,yc,xc);
85 np++;
86
87 }
88
89 }
90
91 if (np > 0) {
5b96ea20 92 gReve->AddRenderElement(l, points);
93 } else {
94 delete points;
34ea2f0a 95 }
34ea2f0a 96
97 }
98
34ea2f0a 99 } // end detElemId loop
100
101 } // end chamber loop
102
103 } // end station loop
104
5b96ea20 105 gReve->Redraw3D();
106
34ea2f0a 107 return;
108}