]> git.uio.no Git - u/mrichter/AliRoot.git/blame - EVE/alice-macros/trd_clusters.C
Restored compilation and functionallity (Ruben)
[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
10#ifdef __CINT__
11class TEvePointSet;
12class TEveElement;
13#else
14#include <TEveManager.h>
15#include <TEvePointSet.h>
176fc00c 16#include <TTree.h>
af13c843 17#include <EveBase/AliEveEventManager.h>
18
19#include "AliRunLoader.h"
20#include "AliCluster.h"
503bfbc8 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
43
af13c843 44 Int_t nentr=(Int_t)recPoints->GetEntries();
45 for (Int_t i=0; i<nentr; i++) {
46 if (!recPoints->GetEvent(i)) continue;
47
48 Int_t ncl=TRDcluster->GetEntriesFast();
49
50 while (ncl--) {
51 AliTRDcluster *c = (AliTRDcluster*)TRDcluster->UncheckedAt(ncl);
52 Float_t g[3]; //global coordinates
53 c->GetGlobalXYZ(g);
503bfbc8 54 clusters->SetNextPoint(g[0], g[1], g[2]);
55 AliCluster *atp = new AliCluster(*c);
56 clusters->SetPointId(atp);
af13c843 57 }
58 TRDcluster->Clear();
59 }
60
9fb2a3ed 61 rl->UnloadRecPoints("TRD");
62
af13c843 63 if(clusters->Size() == 0 && gEve->GetKeepEmptyCont() == kFALSE) {
8ebd7df4 64 Warning("trd_clusters.C", "No TRD clusters");
af13c843 65 delete clusters;
66 return 0;
67 }
68
30650838 69 clusters->SetName("TRD Clusters");
af13c843 70
9dcd42ea 71 clusters->SetTitle(Form("N=%d", clusters->Size()));
f6afd0e1 72
9dcd42ea 73 const TString viz_tag("REC Clusters TRD");
30650838 74
68ca2fe7 75 clusters->ApplyVizTag(viz_tag, "Clusters");
f6afd0e1 76
af13c843 77 gEve->AddElement(clusters, cont);
f6afd0e1 78
af13c843 79 gEve->Redraw3D();
80
81 return clusters;
82}