]> git.uio.no Git - u/mrichter/AliRoot.git/blame - EVE/alice-macros/tpc_clusters.C
Adding more information to the debug output
[u/mrichter/AliRoot.git] / EVE / alice-macros / tpc_clusters.C
CommitLineData
d810d0de 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 *
51346b82 7 * full copyright notice. *
d810d0de 8 **************************************************************************/
cb4245bb 9
32e219c2 10#ifdef __CINT__
ae1afd97 11
84aff7a4 12class TEveElement;
13class TEvePointSet;
32e219c2 14
15#else
16
84aff7a4 17#include <TEveManager.h>
18#include <TEvePointSet.h>
ea5d9491 19#include <EveBase/AliEveEventManager.h>
32e219c2 20
21#include <AliRunLoader.h>
22#include <AliCluster.h>
23#include <TPC/AliTPCClustersRow.h>
24
25#endif
26
84aff7a4 27TEvePointSet* tpc_clusters(TEveElement* cont=0, Float_t maxR=270)
ae1afd97 28{
29 const Int_t kMaxCl=100*160;
30
d810d0de 31 AliEveEventManager::AssertGeometry();
ae1afd97 32
d810d0de 33 AliRunLoader* rl = AliEveEventManager::AssertRunLoader();
ae1afd97 34 rl->LoadRecPoints("TPC");
35
787c89c5 36 TTree *cTree = rl->GetTreeR("TPC", false);
37 if (cTree == 0)
38 return 0;
39
40 AliTPCClustersRow *clrow = new AliTPCClustersRow();
ae1afd97 41 clrow->SetClass("AliTPCclusterMI");
42 clrow->SetArray(kMaxCl);
32e219c2 43 cTree->SetBranchAddress("Segment", &clrow);
ae1afd97 44
787c89c5 45 TEvePointSet* clusters = new TEvePointSet(kMaxCl);
46 clusters->SetOwnIds(kTRUE);
47
48
ae1afd97 49 Float_t maxRsqr = maxR*maxR;
ae1afd97 50 Int_t nentr=(Int_t)cTree->GetEntries();
0abe4b3a 51 for (Int_t i=0; i<nentr; i++)
52 {
ae1afd97 53 if (!cTree->GetEvent(i)) continue;
54
0abe4b3a 55 TClonesArray *cl = clrow->GetArray();
56 Int_t ncl = cl->GetEntriesFast();
ae1afd97 57
0abe4b3a 58 while (ncl--)
59 {
60 AliCluster *c = (AliCluster*) cl->UncheckedAt(ncl);
ae1afd97 61 Float_t g[3]; //global coordinates
62 c->GetGlobalXYZ(g);
63 if (g[0]*g[0]+g[1]*g[1] < maxRsqr)
64 {
65 clusters->SetNextPoint(g[0], g[1], g[2]);
66 AliCluster *atp = new AliCluster(*c);
67 clusters->SetPointId(atp);
68 }
69 }
70 cl->Clear();
71 }
72
73 delete clrow;
74
442b2a2f 75 rl->UnloadRecPoints("TPC");
6b6f47ac 76
0abe4b3a 77 if (clusters->Size() == 0 && gEve->GetKeepEmptyCont() == kFALSE)
78 {
8ebd7df4 79 Warning("tpc_clusters.C", "No TPC clusters");
ae1afd97 80 delete clusters;
81 return 0;
82 }
83
ae1afd97 84 char form[1000];
85 sprintf(form,"TPC Clusters");
86 clusters->SetName(form);
87
88 char tip[1000];
89 sprintf(tip,"N=%d", clusters->Size());
90 clusters->SetTitle(tip);
f6afd0e1 91
92 const TString viz_tag("TPC Clusters");
93 if (gEve->FindVizDBEntry(viz_tag) == 0)
94 {
95 TEvePointSet* m = new TEvePointSet();
96 m->SetMarkerColor(4);
97 m->SetMarkerSize(0.2);
98 m->SetMarkerStyle(2);
99 gEve->InsertVizDBEntry(viz_tag, m);
100 }
101 // The above can be removed when going to new root - then call:
102 // clusters->ApplyVizTag(viz_tag, "Clusters");
103 clusters->ApplyVizTag(viz_tag);
104
84aff7a4 105 gEve->AddElement(clusters, cont);
f6afd0e1 106
84aff7a4 107 gEve->Redraw3D();
ae1afd97 108
109 return clusters;
110}