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