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