]>
Commit | Line | Data |
---|---|---|
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__ | |
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" | |
503bfbc8 | 20 | #include "TRD/AliTRDcluster.h" |
af13c843 | 21 | #endif |
22 | ||
503bfbc8 | 23 | TEvePointSet* trd_clusters(TEveElement *cont = 0) |
24 | { | |
25 | const Int_t kMaxClusters = 18 * 6 * 24 *10; | |
af13c843 | 26 | |
8ebd7df4 | 27 | AliEveEventManager::AssertGeometry(); |
af13c843 | 28 | |
503bfbc8 | 29 | AliRunLoader *rl = AliEveEventManager::AssertRunLoader(); |
30 | rl->LoadRecPoints("TRD"); | |
af13c843 | 31 | |
503bfbc8 | 32 | TTree *recPoints = rl->GetTreeR("TRD", kFALSE); |
8ebd7df4 | 33 | if (recPoints == 0) |
34 | return 0; | |
35 | ||
36 | TObjArray *TRDcluster = 0x0; | |
503bfbc8 | 37 | recPoints->SetBranchAddress("TRDcluster", &TRDcluster); |
af13c843 | 38 | |
8ebd7df4 | 39 | TEvePointSet *clusters = new TEvePointSet(kMaxClusters); |
40 | clusters->SetOwnIds(kTRUE); | |
41 | ||
42 | ||
af13c843 | 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); | |
503bfbc8 | 53 | clusters->SetNextPoint(g[0], g[1], g[2]); |
54 | AliCluster *atp = new AliCluster(*c); | |
55 | clusters->SetPointId(atp); | |
af13c843 | 56 | } |
57 | TRDcluster->Clear(); | |
58 | } | |
59 | ||
9fb2a3ed | 60 | rl->UnloadRecPoints("TRD"); |
61 | ||
af13c843 | 62 | if(clusters->Size() == 0 && gEve->GetKeepEmptyCont() == kFALSE) { |
8ebd7df4 | 63 | Warning("trd_clusters.C", "No TRD clusters"); |
af13c843 | 64 | delete clusters; |
65 | return 0; | |
66 | } | |
67 | ||
30650838 | 68 | clusters->SetName("TRD Clusters"); |
af13c843 | 69 | |
9dcd42ea | 70 | clusters->SetTitle(Form("N=%d", clusters->Size())); |
f6afd0e1 | 71 | |
9dcd42ea | 72 | const TString viz_tag("REC Clusters TRD"); |
30650838 | 73 | |
68ca2fe7 | 74 | clusters->ApplyVizTag(viz_tag, "Clusters"); |
f6afd0e1 | 75 | |
af13c843 | 76 | gEve->AddElement(clusters, cont); |
f6afd0e1 | 77 | |
af13c843 | 78 | gEve->Redraw3D(); |
79 | ||
80 | return clusters; | |
81 | } |