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