]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EVE/alice-macros/trd_clusters.C
1b39c3d36904cf46a712733c6465acae9d3aa9ed
[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 #ifdef __CINT__
11 class TEvePointSet;
12 class TEveElement;
13 #else
14 #include <TEveManager.h>
15 #include <TEvePointSet.h>
16 #include <TTree.h>
17 #include <EveBase/AliEveEventManager.h>
18
19 #include "AliRunLoader.h"
20 #include "AliCluster.h"
21 #include "TRD/AliTRDcluster.h"
22 #endif
23
24 TEvePointSet* trd_clusters(TEveElement *cont = 0)
25 {
26   const Int_t kMaxClusters = 18 * 6 * 24 *10;
27
28   AliEveEventManager::AssertGeometry();
29
30   AliRunLoader *rl = AliEveEventManager::AssertRunLoader();
31   rl->LoadRecPoints("TRD");
32
33   TTree *recPoints = rl->GetTreeR("TRD", kFALSE);
34   if (recPoints == 0)
35     return 0;
36
37   TObjArray *TRDcluster = 0x0;
38   recPoints->SetBranchAddress("TRDcluster", &TRDcluster);
39
40   TEvePointSet *clusters = new TEvePointSet(kMaxClusters);
41   clusters->SetOwnIds(kTRUE);
42
43
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);
54       clusters->SetNextPoint(g[0], g[1], g[2]);
55       AliCluster *atp = new AliCluster(*c);
56       clusters->SetPointId(atp);
57     }
58     TRDcluster->Clear();
59   }
60
61   rl->UnloadRecPoints("TRD");
62
63   if(clusters->Size() == 0 && gEve->GetKeepEmptyCont() == kFALSE) {
64     Warning("trd_clusters.C", "No TRD clusters");
65     delete clusters;
66     return 0;
67   }
68
69   clusters->SetName("TRD Clusters");
70
71   clusters->SetTitle(Form("N=%d", clusters->Size()));
72
73   const TString viz_tag("REC Clusters TRD");
74
75   clusters->ApplyVizTag(viz_tag, "Clusters");
76
77   gEve->AddElement(clusters, cont);
78
79   gEve->Redraw3D();
80
81   return clusters;
82 }