]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EVE/alice-macros/tpc_clusters.C
Coverity
[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
10 #ifdef __CINT__
11
12 class TEveElement;
13 class TEvePointSet;
14
15 #else
16
17 #include <TEveManager.h>
18 #include <TColor.h>
19 #include <TEvePointSet.h>
20 #include <EveBase/AliEveEventManager.h>
21
22 #include <AliRunLoader.h>
23 #include <AliCluster.h>
24 #include <TPC/AliTPCClustersRow.h>
25 #include <TPC/AliTPCclusterMI.h>
26
27 #endif
28
29 TEvePointSet* tpc_clusters(TEveElement* cont=0, Float_t maxR=270)
30 {
31
32   const Int_t kMaxCl=100*160;
33
34   Int_t fNColorBins = 5;
35
36   AliEveEventManager::AssertGeometry();
37
38   AliRunLoader* rl = AliEveEventManager::AssertRunLoader();
39   rl->LoadRecPoints("TPC");
40
41   TTree *cTree = rl->GetTreeR("TPC", false);
42   if (cTree == 0)
43     return 0;
44
45   AliTPCClustersRow *clrow = new AliTPCClustersRow();
46   clrow->SetClass("AliTPCclusterMI");
47   clrow->SetArray(kMaxCl);
48   cTree->SetBranchAddress("Segment", &clrow);
49
50   TEvePointSet* clusters = new TEvePointSet(kMaxCl);
51   clusters->SetOwnIds(kTRUE);
52
53   TEvePointSetArray * cc = new TEvePointSetArray("TPC Clusters Colorized");
54   cc->SetMainColor(kRed);
55   cc->SetMarkerStyle(4);
56   cc->SetMarkerSize(0.4);
57   cc->InitBins("Cluster Charge", fNColorBins, 0., fNColorBins*60.);
58   
59   cc->GetBin(0)->SetMainColor(kGray);
60   cc->GetBin(0)->SetMarkerSize(0.4);
61   cc->GetBin(1)->SetMainColor(kBlue);
62   cc->GetBin(1)->SetMarkerSize(0.42);
63   cc->GetBin(2)->SetMainColor(kCyan);
64   cc->GetBin(2)->SetMarkerSize(0.44);
65   cc->GetBin(3)->SetMainColor(kGreen);
66   cc->GetBin(3)->SetMarkerSize(0.46);
67   cc->GetBin(4)->SetMainColor(kYellow);
68   cc->GetBin(4)->SetMarkerSize(0.48);
69   cc->GetBin(5)->SetMainColor(kRed);
70   cc->GetBin(5)->SetMarkerSize(0.50);
71   cc->GetBin(6)->SetMainColor(kMagenta);
72   cc->GetBin(6)->SetMarkerSize(0.52);
73
74   Float_t maxRsqr = maxR*maxR;
75   Int_t nentr=(Int_t)cTree->GetEntries();
76   for (Int_t i=0; i<nentr; i++)
77   {
78     if (!cTree->GetEvent(i)) continue;
79
80     TClonesArray *cl = clrow->GetArray();
81     Int_t ncl = cl->GetEntriesFast();
82
83     while (ncl--)
84     {
85
86       AliTPCclusterMI* clusterMi = (AliTPCclusterMI*) cl->At(ncl);
87
88       AliCluster *c = (AliCluster*) cl->UncheckedAt(ncl);
89       Float_t g[3]; //global coordinates
90       c->GetGlobalXYZ(g);
91       if (g[0]*g[0]+g[1]*g[1] < maxRsqr)
92       {
93         cc->Fill(g[0], g[1], g[2], clusterMi->GetQ());
94         clusters->SetNextPoint(g[0], g[1], g[2]);
95         AliCluster *atp = new AliCluster(*c);
96         clusters->SetPointId(atp);
97       }
98     }
99     cl->Clear();
100   }
101
102   delete clrow;
103
104   rl->UnloadRecPoints("TPC");
105
106   if (clusters->Size() == 0 && gEve->GetKeepEmptyCont() == kFALSE)
107   {
108     Warning("tpc_clusters.C", "No TPC clusters");
109     delete clusters;
110     return 0;
111   }
112
113   clusters->SetName("TPC Clusters");
114
115   clusters->SetTitle(Form("N=%d", clusters->Size()));
116
117   const TString viz_tag("REC Clusters TPC");
118
119   clusters->ApplyVizTag(viz_tag, "Clusters");
120
121 //  clusters->SetRnrSelf(kFALSE);
122 //  clusters->SetRnrChildren(kFALSE);    
123
124 //  gEve->AddElement(clusters, cont);
125
126   cc->SetRnrSelf(kTRUE);
127
128   gEve->AddElement(cc);
129
130   gEve->Redraw3D();
131
132   return clusters;
133 }