]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EVE/alice-macros/tof_clusters.C
a1cee3214d982c52d4679bd00e333194ee5cce55
[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 <TMath.h>
16
17 #include <TEveManager.h>
18 #include <TEvePointSet.h>
19 #include <EveBase/AliEveEventManager.h>
20
21 #include <AliRunLoader.h>
22 #include <AliCluster.h>
23
24 #endif
25
26 TEvePointSet* tof_clusters(TEveElement* cont=0, Float_t maxR=390)
27 {
28   AliEveEventManager::AssertGeometry();
29
30   AliRunLoader* rl = AliEveEventManager::AssertRunLoader();
31   rl->LoadRecPoints("TOF");
32
33   TTree *cTree = rl->GetTreeR("TOF", false);
34   if (cTree == 0)
35     return 0;
36
37   TObjArray *cl = 0x0;
38   cTree->SetBranchAddress("TOF", &cl);
39
40   TEvePointSet* clusters = new TEvePointSet(10000);
41   clusters->SetOwnIds(kTRUE);
42
43   Float_t maxRsqr = maxR*maxR;
44
45   Int_t nentr=(Int_t)cTree->GetEntries();
46   for (Int_t i=0; i<nentr; i++)
47   {
48     if (!cTree->GetEvent(i)) continue;
49
50     Int_t ncl=cl->GetEntriesFast();
51     while (ncl--)
52     {
53       AliCluster *c=(AliCluster*)cl->UncheckedAt(ncl);
54       Float_t g[3]; //global coordinates
55       c->GetGlobalXYZ(g);
56
57       if (g[0]*g[0]+g[1]*g[1] < maxRsqr) {
58         clusters->SetNextPoint(g[0], g[1], g[2]);
59         AliCluster *atp = new AliCluster(*c);
60         clusters->SetPointId(atp);
61       }
62     }
63   }
64
65   rl->UnloadRecPoints("TOF");
66
67   if (clusters->Size() == 0 && gEve->GetKeepEmptyCont() == kFALSE)
68   {
69     Warning("tof_clusters.C", "No TOF clusters");
70     delete clusters;
71     return 0;
72   }
73
74   clusters->SetName("TOF Clusters");
75   clusters->SetTitle(Form("N=%d", clusters->Size()));
76
77   const TString viz_tag("REC Clusters TOF");
78   clusters->ApplyVizTag(viz_tag, "Clusters");
79
80   gEve->AddElement(clusters, cont);
81   gEve->Redraw3D();
82
83   return clusters;
84 }
85
86 TEvePointSet* tof_clusters_sec(Int_t selectedSector,
87                                TEveElement* cont=0, Float_t maxR=390)
88 {
89   AliEveEventManager::AssertGeometry();
90
91   AliRunLoader* rl = AliEveEventManager::AssertRunLoader();
92   rl->LoadRecPoints("TOF");
93
94   TTree *cTree = rl->GetTreeR("TOF", false);
95   if (cTree == 0)
96     return 0;
97
98   TObjArray *cl = 0x0;
99   cTree->SetBranchAddress("TOF", &cl);
100
101   TEvePointSet* clusters = new TEvePointSet(10000);
102   clusters->SetOwnIds(kTRUE);
103
104   Float_t maxRsqr = maxR*maxR;
105   Double_t phiAngle = maxR*maxR;
106
107   Int_t nentr=(Int_t)cTree->GetEntries();
108   for (Int_t i=0; i<nentr; i++) {
109     if (!cTree->GetEvent(i)) continue;
110
111     Int_t ncl=cl->GetEntriesFast();
112     while (ncl--) {
113       AliCluster *c=(AliCluster*)cl->UncheckedAt(ncl);
114       Float_t g[3]; //global coordinates
115       c->GetGlobalXYZ(g);
116
117       phiAngle = 180./TMath::Pi()*(TMath::ATan2(-g[1],-g[0])+TMath::Pi());
118
119       if (g[0]*g[0]+g[1]*g[1] < maxRsqr) {
120         if (
121             (selectedSector!=-1 &&
122              phiAngle>=selectedSector*20. && phiAngle<(selectedSector+1)*20.)
123             ||
124             selectedSector==-1)
125           {
126             clusters->SetNextPoint(g[0], g[1], g[2]);
127             AliCluster *atp = new AliCluster(*c);
128             clusters->SetPointId(atp);
129           }
130       }
131
132     }
133   }
134
135   rl->UnloadRecPoints("TOF");
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   clusters->SetName(Form("REC Clusters TOF"));
144   clusters->SetTitle(Form("N=%d", clusters->Size()));
145
146   const TString viz_tag("REC Clusters TOF");
147   clusters->ApplyVizTag(viz_tag, "Clusters");
148   
149   gEve->AddElement(clusters, cont);
150   gEve->Redraw3D();
151
152   return clusters;
153 }