]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EVE/alice-macros/trd_clusters.C
cosmetics
[u/mrichter/AliRoot.git] / EVE / alice-macros / trd_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 <TObjArray.h>
12 #include <TTree.h>
13 #include <TEveManager.h>
14 #include <TEveElement.h>
15 #include <TEvePointSet.h>
16 #include <TGeoMatrix.h>
17
18 #include <AliCluster.h>
19 #include <AliGeomManager.h>
20 #include <AliRunLoader.h>
21 #include <AliTRDcluster.h>
22 #include <AliEveEventManager.h>
23 #else
24 class TEvePointSet;
25 class TEveElement;
26 #endif
27
28 TEvePointSet* trd_clusters(TEveElement *cont = 0)
29 {
30   const Int_t kMaxClusters = 18 * 6 * 24 *10;
31
32   AliEveEventManager::AssertGeometry();
33
34   AliRunLoader *rl = AliEveEventManager::AssertRunLoader();
35   rl->LoadRecPoints("TRD");
36
37   TTree *recPoints = rl->GetTreeR("TRD", kFALSE);
38   if (recPoints == 0)
39     return 0;
40
41   TObjArray *TRDcluster = 0x0;
42   recPoints->SetBranchAddress("TRDcluster", &TRDcluster);
43
44   TEvePointSet *clusters = new TEvePointSet(kMaxClusters);
45   clusters->SetOwnIds(kTRUE);
46
47   Int_t nentr=(Int_t)recPoints->GetEntries();
48   for (Int_t i=0; i<nentr; i++) {
49     if (!recPoints->GetEvent(i)) continue;
50
51     Int_t ncl=TRDcluster->GetEntriesFast();
52
53     while (ncl--) {
54       AliTRDcluster *c = (AliTRDcluster*)TRDcluster->UncheckedAt(ncl);
55       //Float_t g[3]; //global coordinates
56       //c->GetGlobalXYZ(g);
57       
58       Int_t fVolumeId = c->GetVolumeId();
59       const TGeoHMatrix *mt =AliGeomManager::GetTracking2LocalMatrix(fVolumeId);;
60       Double_t txyz[3] = {c->GetX(), c->GetY(), c->GetZ()};
61       Double_t lxyz[3] = {0, 0, 0};
62       mt->LocalToMaster(txyz,lxyz);
63    
64       TGeoHMatrix *mlIdeal = AliGeomManager::GetOrigGlobalMatrix(fVolumeId);
65       Double_t gxyzIdeal[3] = {0, 0, 0};
66       mlIdeal->LocalToMaster(lxyz,gxyzIdeal);
67       
68       clusters->SetNextPoint(gxyzIdeal[0], gxyzIdeal[1], gxyzIdeal[2]);
69       
70       AliCluster *atp = new AliCluster(*c);
71       clusters->SetPointId(atp);
72     }
73     TRDcluster->Clear();
74   }
75
76   rl->UnloadRecPoints("TRD");
77
78   if(clusters->Size() == 0 && gEve->GetKeepEmptyCont() == kFALSE) {
79     Warning("trd_clusters.C", "No TRD clusters");
80     delete clusters;
81     return 0;
82   }
83
84   clusters->SetName("TRD Clusters");
85
86   clusters->SetTitle(Form("N=%d", clusters->Size()));
87
88   const TString viz_tag("REC Clusters TRD");
89
90   clusters->ApplyVizTag(viz_tag, "Clusters");
91
92   gEve->AddElement(clusters, cont);
93
94   gEve->Redraw3D();
95
96   return clusters;
97 }