]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EVE/alice-macros/tpc_clusters.C
Prototype for visualization-macro manager and gui.
[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 <TEvePointSet.h>
19 #include <EveBase/AliEveEventManager.h>
20
21 #include <AliRunLoader.h>
22 #include <AliCluster.h>
23 #include <TPC/AliTPCClustersRow.h>
24
25 #endif
26
27 TEvePointSet* tpc_clusters(TEveElement* cont=0, Float_t maxR=270)
28 {
29   const Int_t kMaxCl=100*160;
30
31   AliEveEventManager::AssertGeometry();
32
33   AliRunLoader* rl = AliEveEventManager::AssertRunLoader();
34   rl->LoadRecPoints("TPC");
35
36   TTree *cTree = rl->GetTreeR("TPC", false);
37   if (cTree == 0)
38     return 0;
39
40   AliTPCClustersRow *clrow = new AliTPCClustersRow();
41   clrow->SetClass("AliTPCclusterMI");
42   clrow->SetArray(kMaxCl);
43   cTree->SetBranchAddress("Segment", &clrow);
44
45   TEvePointSet* clusters = new TEvePointSet(kMaxCl);
46   clusters->SetOwnIds(kTRUE);
47
48
49   Float_t maxRsqr = maxR*maxR;
50   Int_t nentr=(Int_t)cTree->GetEntries();
51   for (Int_t i=0; i<nentr; i++)
52   {
53     if (!cTree->GetEvent(i)) continue;
54
55     TClonesArray *cl = clrow->GetArray();
56     Int_t ncl = cl->GetEntriesFast();
57
58     while (ncl--)
59     {
60       AliCluster *c = (AliCluster*) cl->UncheckedAt(ncl);
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
75   rl->UnloadRecPoints("TPC");
76
77   if (clusters->Size() == 0 && gEve->GetKeepEmptyCont() == kFALSE)
78   {
79     Warning("tpc_clusters.C", "No TPC clusters");
80     delete clusters;
81     return 0;
82   }
83
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);
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
105   gEve->AddElement(clusters, cont);
106
107   gEve->Redraw3D();
108
109   return clusters;
110 }