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