]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EVE/alice-macros/trd_clusters.C
Update to new method names in AliESDv0, load them in visscan_init.C.
[u/mrichter/AliRoot.git] / EVE / alice-macros / trd_clusters.C
1 // $Id: tpc_clusters.C 23497 2008-01-23 20:43:14Z mtadel $
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 <EveBase/AliEveEventManager.h>
17
18 #include "AliRunLoader.h"
19 #include "AliCluster.h"
20 #include "TRD/AliTRDcluster.h"
21 #endif
22
23 TEvePointSet* trd_clusters(TEveElement *cont = 0)
24 {
25   const Int_t kMaxClusters = 18 * 6 * 24 *10;
26   AliEveEventManager::AssertGeometry();
27
28   TEvePointSet *clusters = new TEvePointSet(kMaxClusters);
29   clusters->SetOwnIds(kTRUE);
30
31   AliRunLoader *rl = AliEveEventManager::AssertRunLoader();
32   rl->LoadRecPoints("TRD");
33
34   TObjArray *TRDcluster = 0x0;
35   TTree *recPoints = rl->GetTreeR("TRD", kFALSE);
36   recPoints->SetBranchAddress("TRDcluster", &TRDcluster);
37
38   Int_t nentr=(Int_t)recPoints->GetEntries();
39   for (Int_t i=0; i<nentr; i++) {
40     if (!recPoints->GetEvent(i)) continue;
41
42     Int_t ncl=TRDcluster->GetEntriesFast();
43
44     while (ncl--) {
45       AliTRDcluster *c = (AliTRDcluster*)TRDcluster->UncheckedAt(ncl);
46       Float_t g[3]; //global coordinates
47       c->GetGlobalXYZ(g);
48       clusters->SetNextPoint(g[0], g[1], g[2]);
49       AliCluster *atp = new AliCluster(*c);
50       clusters->SetPointId(atp);
51     }
52     TRDcluster->Clear();
53   }
54
55   if(clusters->Size() == 0 && gEve->GetKeepEmptyCont() == kFALSE) {
56     Warning("trd_clusters", "No TRD clusters");
57     delete clusters;
58     return 0;
59   }
60
61   clusters->SetMarkerStyle(2);
62   clusters->SetMarkerSize(0.2);
63   clusters->SetMarkerColor(4);
64
65   char form[1000];
66   sprintf(form,"TRD Clusters");
67   clusters->SetName(form);
68
69   char tip[1000];
70   sprintf(tip,"N=%d", clusters->Size());
71   clusters->SetTitle(tip);
72   gEve->AddElement(clusters, cont);
73   gEve->Redraw3D();
74
75   return clusters;
76 }