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