]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EVE/alice-macros/tof_clusters.C
fcb9e26c8eb159a8568ab3cff00a8365c8d1fd0c
[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   char form[1000];
69   sprintf(form,"TOF Clusters");
70   clusters->SetName(form);
71
72   char tip[1000];
73   sprintf(tip,"N=%d", clusters->Size());
74   clusters->SetTitle(tip);
75
76   const TString viz_tag("TOF Clusters");
77   // when going to new root call:
78   // clusters->ApplyVizTag(viz_tag, "Clusters");
79   clusters->ApplyVizTag(viz_tag);
80
81   gEve->AddElement(clusters, cont);
82
83   gEve->Redraw3D();
84
85   return clusters;
86 }
87
88 TEvePointSet* tof_clusters_sec(Int_t selectedSector,
89                                TEveElement* cont=0, Float_t maxR=390)
90 {
91   AliEveEventManager::AssertGeometry();
92
93   AliRunLoader* rl = AliEveEventManager::AssertRunLoader();
94   rl->LoadRecPoints("TOF");
95
96   TTree *cTree = rl->GetTreeR("TOF", false);
97   if (cTree == 0)
98     return 0;
99
100   TObjArray *cl = 0x0;
101   cTree->SetBranchAddress("TOF", &cl);
102
103   TEvePointSet* clusters = new TEvePointSet(10000);
104   clusters->SetOwnIds(kTRUE);
105
106   Float_t maxRsqr = maxR*maxR;
107   Double_t phiAngle = maxR*maxR;
108
109   Int_t nentr=(Int_t)cTree->GetEntries();
110   for (Int_t i=0; i<nentr; i++) {
111     if (!cTree->GetEvent(i)) continue;
112
113     Int_t ncl=cl->GetEntriesFast();
114     while (ncl--) {
115       AliCluster *c=(AliCluster*)cl->UncheckedAt(ncl);
116       Float_t g[3]; //global coordinates
117       c->GetGlobalXYZ(g);
118
119       phiAngle = 180./TMath::Pi()*(TMath::ATan2(-g[1],-g[0])+TMath::Pi());
120
121       if (g[0]*g[0]+g[1]*g[1] < maxRsqr) {
122         if (
123             (selectedSector!=-1 &&
124              phiAngle>=selectedSector*20. && phiAngle<(selectedSector+1)*20.)
125             ||
126             selectedSector==-1)
127           {
128             clusters->SetNextPoint(g[0], g[1], g[2]);
129             AliCluster *atp = new AliCluster(*c);
130             clusters->SetPointId(atp);
131           }
132       }
133
134     }
135   }
136
137   if (clusters->Size() == 0 && gEve->GetKeepEmptyCont() == kFALSE) {
138     Warning("tof_clusters.C", "No TOF clusters");
139     delete clusters;
140     return 0;
141   }
142
143   char form[1000];
144   sprintf(form,"TOF Clusters");
145   clusters->SetName(form);
146
147   char tip[1000];
148   sprintf(tip,"N=%d", clusters->Size());
149   clusters->SetTitle(tip);
150
151   const TString viz_tag("TOF Clusters");
152   // when going to new root call:
153   // clusters->ApplyVizTag(viz_tag, "Clusters");
154   clusters->ApplyVizTag(viz_tag);
155
156   gEve->AddElement(clusters, cont);
157
158   gEve->Redraw3D();
159
160   return clusters;
161 }