]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EVE/alice-macros/tpc_clusters.C
Remove trailing whitespace.
[u/mrichter/AliRoot.git] / EVE / alice-macros / tpc_clusters.C
1 // $Id$
2 // Main authors: Matevz Tadel & Alja Mrak-Tadel: 2006, 2007
3
4 /**************************************************************************
5  * Copyright(c) 1998-2008, ALICE Experiment at CERN, all rights reserved. *
6  * See http://aliceinfo.cern.ch/Offline/AliRoot/License.html for          *
7  * full copyright notice.                                                 *
8  **************************************************************************/
9 #ifdef __CINT__
10
11 namespace TEveUtil
12 {
13 class TEveElement;
14 class TEvePointSet;
15 }
16
17 #else
18
19 #include <TEveManager.h>
20 #include <TEvePointSet.h>
21 #include <Alieve/AliEveEventManager.h>
22
23 #include <AliRunLoader.h>
24 #include <AliCluster.h>
25 #include <TPC/AliTPCClustersRow.h>
26
27 #endif
28
29 TEvePointSet* tpc_clusters(TEveElement* cont=0, Float_t maxR=270)
30 {
31   const Int_t kMaxCl=100*160;
32
33   AliEveEventManager::AssertGeometry();
34
35   TEvePointSet* clusters = new TEvePointSet(kMaxCl);
36   clusters->SetOwnIds(kTRUE);
37
38   AliRunLoader* rl = AliEveEventManager::AssertRunLoader();
39   rl->LoadRecPoints("TPC");
40
41   AliTPCClustersRow *clrow=new AliTPCClustersRow();
42   clrow->SetClass("AliTPCclusterMI");
43   clrow->SetArray(kMaxCl);
44
45   TTree *cTree = rl->GetTreeR("TPC", false);
46   cTree->SetBranchAddress("Segment", &clrow);
47
48   Float_t maxRsqr = maxR*maxR;
49   TClonesArray *cl=clrow->GetArray();
50   Int_t nentr=(Int_t)cTree->GetEntries();
51   for (Int_t i=0; i<nentr; i++) {
52     if (!cTree->GetEvent(i)) continue;
53
54     Int_t ncl=cl->GetEntriesFast();
55
56     while (ncl--) {
57       AliCluster *c=(AliCluster*)cl->UncheckedAt(ncl);
58       Float_t g[3]; //global coordinates
59       c->GetGlobalXYZ(g);
60       if (g[0]*g[0]+g[1]*g[1] < maxRsqr)
61       {
62         clusters->SetNextPoint(g[0], g[1], g[2]);
63         AliCluster *atp = new AliCluster(*c);
64         clusters->SetPointId(atp);
65       }
66     }
67     cl->Clear();
68   }
69
70   delete clrow;
71
72   if(clusters->Size() == 0 && gEve->GetKeepEmptyCont() == kFALSE) {
73     Warning("tpc_clusters", "No TPC clusters");
74     delete clusters;
75     return 0;
76   }
77
78   clusters->SetMarkerStyle(2);
79   clusters->SetMarkerSize(0.2);
80   clusters->SetMarkerColor(4);
81
82   char form[1000];
83   sprintf(form,"TPC Clusters");
84   clusters->SetName(form);
85
86   char tip[1000];
87   sprintf(tip,"N=%d", clusters->Size());
88   clusters->SetTitle(tip);
89   gEve->AddElement(clusters, cont);
90   gEve->Redraw3D();
91
92   return clusters;
93 }