]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EVE/alice-macros/trd_clusters.C
Adding more information to the debug output
[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
27   AliEveEventManager::AssertGeometry();
28
29   AliRunLoader *rl = AliEveEventManager::AssertRunLoader();
30   rl->LoadRecPoints("TRD");
31
32   TTree *recPoints = rl->GetTreeR("TRD", kFALSE);
33   if (recPoints == 0)
34     return 0;
35
36   TObjArray *TRDcluster = 0x0;
37   recPoints->SetBranchAddress("TRDcluster", &TRDcluster);
38
39   TEvePointSet *clusters = new TEvePointSet(kMaxClusters);
40   clusters->SetOwnIds(kTRUE);
41
42
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);
51       Float_t g[3]; //global coordinates
52       c->GetGlobalXYZ(g);
53       clusters->SetNextPoint(g[0], g[1], g[2]);
54       AliCluster *atp = new AliCluster(*c);
55       clusters->SetPointId(atp);
56     }
57     TRDcluster->Clear();
58   }
59
60   if(clusters->Size() == 0 && gEve->GetKeepEmptyCont() == kFALSE) {
61     Warning("trd_clusters.C", "No TRD clusters");
62     delete clusters;
63     return 0;
64   }
65
66   char form[1000];
67   sprintf(form,"TRD Clusters");
68   clusters->SetName(form);
69
70   char tip[1000];
71   sprintf(tip,"N=%d", clusters->Size());
72   clusters->SetTitle(tip);
73
74   const TString viz_tag("TRD Clusters");
75   if (gEve->FindVizDBEntry(viz_tag) == 0)
76   {
77     TEvePointSet* m = new TEvePointSet();
78     m->SetMarkerColor(4);
79     m->SetMarkerSize(0.2);
80     m->SetMarkerStyle(2);
81     gEve->InsertVizDBEntry(viz_tag, m);
82   }
83   // The above can be removed when going to new root - then call:
84   // clusters->ApplyVizTag(viz_tag, "Clusters");
85   clusters->ApplyVizTag(viz_tag);
86
87   gEve->AddElement(clusters, cont);
88
89   gEve->Redraw3D();
90
91   return clusters;
92 }