]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EVE/alice-macros/tof_clusters.C
Change to support compilation when ALICE_INSTALL is not same as ALICE_ROOT
[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
9 #if !defined(__CINT__) || defined(__MAKECINT__)
10 #include <TObjArray.h>
11 #include <TMath.h>
12 #include <TTree.h>
13 #include <TString.h>
14 #include <TEveManager.h>
15 #include <TEvePointSet.h>
16
17 #include <AliCluster.h>
18 #include <AliRunLoader.h>
19 #include <AliEveEventManager.h>
20 #else
21 class TEveElement;
22 class TEvePointSet;
23 #endif
24
25 TEvePointSet* tof_clusters(TEveElement* cont=0, Float_t maxR=390)
26 {
27   AliEveEventManager::AssertGeometry();
28
29   AliRunLoader* rl = AliEveEventManager::AssertRunLoader();
30   rl->LoadRecPoints("TOF");
31
32   TTree *cTree = rl->GetTreeR("TOF", false);
33   if (cTree == 0)
34     return 0;
35
36   TObjArray *cl = 0x0;
37   cTree->SetBranchAddress("TOF", &cl);
38
39   TEvePointSet* clusters = new TEvePointSet(10000);
40   clusters->SetOwnIds(kTRUE);
41
42   Float_t maxRsqr = maxR*maxR;
43
44   Int_t nentr=(Int_t)cTree->GetEntries();
45   for (Int_t i=0; i<nentr; i++)
46   {
47     if (!cTree->GetEvent(i)) continue;
48
49     Int_t ncl=cl->GetEntriesFast();
50     while (ncl--)
51     {
52       AliCluster *c=(AliCluster*)cl->UncheckedAt(ncl);
53       Float_t g[3]; //global coordinates
54       c->GetGlobalXYZ(g);
55
56       if (g[0]*g[0]+g[1]*g[1] < maxRsqr) {
57         clusters->SetNextPoint(g[0], g[1], g[2]);
58         AliCluster *atp = new AliCluster(*c);
59         clusters->SetPointId(atp);
60       }
61     }
62   }
63
64   rl->UnloadRecPoints("TOF");
65
66   if (clusters->Size() == 0 && gEve->GetKeepEmptyCont() == kFALSE)
67   {
68     Warning("tof_clusters.C", "No TOF clusters");
69     delete clusters;
70     return 0;
71   }
72
73   clusters->SetName("TOF Clusters");
74   clusters->SetTitle(Form("N=%d", clusters->Size()));
75
76   const TString viz_tag("REC Clusters TOF");
77   clusters->ApplyVizTag(viz_tag, "Clusters");
78
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   rl->UnloadRecPoints("TOF");
135
136   if (clusters->Size() == 0 && gEve->GetKeepEmptyCont() == kFALSE) {
137     Warning("tof_clusters.C", "No TOF clusters");
138     delete clusters;
139     return 0;
140   }
141
142   clusters->SetName(Form("REC Clusters TOF"));
143   clusters->SetTitle(Form("N=%d", clusters->Size()));
144
145   const TString viz_tag("REC Clusters TOF");
146   clusters->ApplyVizTag(viz_tag, "Clusters");
147   
148   gEve->AddElement(clusters, cont);
149   gEve->Redraw3D();
150
151   return clusters;
152 }