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