]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EVE/alice-macros/tof_clusters.C
AliEveEventManager
[u/mrichter/AliRoot.git] / EVE / alice-macros / tof_clusters.C
1 // Main authors: Matevz Tadel & Alja Mrak-Tadel: 2006, 2007
2
3 /**************************************************************************
4  * Copyright(c) 1998-2008, ALICE Experiment at CERN, all rights reserved. *
5  * See http://aliceinfo.cern.ch/Offline/AliRoot/License.html for          *
6  * full copyright notice.                                                 *
7  **************************************************************************/
8 #ifdef __CINT__
9
10 class TEveElement;
11 class TEvePointSet;
12
13 #else
14
15 #include <TEveManager.h>
16 #include <TEvePointSet.h>
17 #include <EveBase/AliEveEventManager.h>
18
19 #include <AliRunLoader.h>
20 #include <AliCluster.h>
21
22 #endif
23
24 TEvePointSet* tof_clusters(TEveElement* cont=0, Float_t maxR=390)
25 {
26   AliEveEventManager::AssertGeometry();
27
28   AliRunLoader* rl = AliEveEventManager::AssertRunLoader();
29   rl->LoadRecPoints("TOF");
30
31   TTree *cTree = rl->GetTreeR("TOF", false);
32   if (cTree == 0)
33     return 0;
34
35   TObjArray *cl = 0x0;
36   cTree->SetBranchAddress("TOF", &cl);
37
38   TEvePointSet* clusters = new TEvePointSet(10000);
39   clusters->SetOwnIds(kTRUE);
40
41   Float_t maxRsqr = maxR*maxR;
42
43   Int_t nentr=(Int_t)cTree->GetEntries();
44   for (Int_t i=0; i<nentr; i++) {
45     if (!cTree->GetEvent(i)) continue;
46
47     Int_t ncl=cl->GetEntriesFast();
48     while (ncl--) {
49       AliCluster *c=(AliCluster*)cl->UncheckedAt(ncl);
50       Float_t g[3]; //global coordinates
51       c->GetGlobalXYZ(g);
52
53       if (g[0]*g[0]+g[1]*g[1] < maxRsqr) {
54         clusters->SetNextPoint(g[0], g[1], g[2]);
55         AliCluster *atp = new AliCluster(*c);
56         clusters->SetPointId(atp);
57       }
58
59     }
60   }
61
62   if (clusters->Size() == 0 && gEve->GetKeepEmptyCont() == kFALSE) {
63     Warning("tof_clusters.C", "No TOF clusters");
64     delete clusters;
65     return 0;
66   }
67
68   clusters->SetMarkerStyle(2);
69   clusters->SetMarkerSize(0.2);
70   clusters->SetMarkerColor(4);
71
72   char form[1000];
73   sprintf(form,"TOF Clusters");
74   clusters->SetName(form);
75
76   char tip[1000];
77   sprintf(tip,"N=%d", clusters->Size());
78   clusters->SetTitle(tip);
79   gEve->AddElement(clusters, cont);
80   gEve->Redraw3D();
81
82   return clusters;
83 }
84
85 TEvePointSet* tof_clusters_sec(Int_t selectedSector,
86                                TEveElement* cont=0, Float_t maxR=390)
87 {
88   AliEveEventManager::AssertGeometry();
89
90   AliRunLoader* rl = AliEveEventManager::AssertRunLoader();
91   rl->LoadRecPoints("TOF");
92
93   TTree *cTree = rl->GetTreeR("TOF", false);
94   if (cTree == 0)
95     return 0;
96
97   TObjArray *cl = 0x0;
98   cTree->SetBranchAddress("TOF", &cl);
99
100   TEvePointSet* clusters = new TEvePointSet(10000);
101   clusters->SetOwnIds(kTRUE);
102
103   Float_t maxRsqr = maxR*maxR;
104   Double_t phiAngle = maxR*maxR;
105
106   Int_t nentr=(Int_t)cTree->GetEntries();
107   for (Int_t i=0; i<nentr; i++) {
108     if (!cTree->GetEvent(i)) continue;
109
110     Int_t ncl=cl->GetEntriesFast();
111     while (ncl--) {
112       AliCluster *c=(AliCluster*)cl->UncheckedAt(ncl);
113       Float_t g[3]; //global coordinates
114       c->GetGlobalXYZ(g);
115
116       phiAngle = 180./TMath::Pi()*(TMath::ATan2(-g[1],-g[0])+TMath::Pi());
117
118       if (g[0]*g[0]+g[1]*g[1] < maxRsqr) {
119         if (
120             (selectedSector!=-1 &&
121              phiAngle>=selectedSector*20. && phiAngle<(selectedSector+1)*20.)
122             ||
123             selectedSector==-1)
124           {
125             clusters->SetNextPoint(g[0], g[1], g[2]);
126             AliCluster *atp = new AliCluster(*c);
127             clusters->SetPointId(atp);
128           }
129       }
130
131     }
132   }
133
134   if (clusters->Size() == 0 && gEve->GetKeepEmptyCont() == kFALSE) {
135     Warning("tof_clusters.C", "No TOF clusters");
136     delete clusters;
137     return 0;
138   }
139
140   clusters->SetMarkerStyle(2);
141   clusters->SetMarkerSize(0.2);
142   clusters->SetMarkerColor(4);
143
144   char form[1000];
145   sprintf(form,"TOF Clusters");
146   clusters->SetName(form);
147
148   char tip[1000];
149   sprintf(tip,"N=%d", clusters->Size());
150   clusters->SetTitle(tip);
151   gEve->AddElement(clusters, cont);
152   gEve->Redraw3D();
153
154   return clusters;
155 }