]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EVE/alice-macros/itsU_clusters.C
SetSeed dummy implementations
[u/mrichter/AliRoot.git] / EVE / alice-macros / itsU_clusters.C
1 // $Id: its_clusters.C 55060 2012-03-09 18:13:17Z quark $
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 <TClonesArray.h>
12 #include <TBranch.h>
13 #include <TTree.h>
14 #include <TEveManager.h>
15 #include <TEveElement.h>
16 #include <TEvePointSet.h>
17
18 #include "AliRunLoader.h"
19 #include "AliCluster.h"
20 #include "AliEveEventManager.h"
21 #include "AliGeomManager.h"
22 #include "../ITS/UPGRADE/AliITSUGeomTGeo.h"
23 #include "../ITS/UPGRADE/AliITSUSegmentationPix.h"
24 #include "../ITS/UPGRADE/AliITSUClusterPix.h"
25 #include "TGeoManager.h"
26
27 #else
28 class TEveElement;
29 class TEvePointSet;
30 class TTree;
31 class TBranch;
32 #endif
33
34 void itsU_clusters(TEveElement* cont=0, Float_t maxR=50)
35 {
36   AliEveEventManager::AssertGeometry();
37
38   AliRunLoader* rl = AliEveEventManager::AssertRunLoader();
39   rl->LoadRecPoints("ITS");
40
41   gGeoManager = gEve->GetGeometry("geometry.root");
42   AliITSUGeomTGeo* gm = new AliITSUGeomTGeo(kTRUE,kTRUE);
43   AliITSUClusterPix::SetGeom(gm);
44   TTree *cTree = rl->GetTreeR("ITS", false);
45   if (cTree == 0)
46     return ;
47
48   TObjArray layerClus;
49   int nlr=0;
50   while(1) {
51     TBranch* br = cTree->GetBranch(Form("ITSRecPoints%d",nlr));
52     if (!br) break;
53     TClonesArray* clr = 0;
54     br->SetAddress(&clr);
55     layerClus.AddLast(clr);
56     nlr++;
57   }
58
59   TEveElementList* evClusters = new TEveElementList("ITS clusters");
60   TEvePointSet* layClusters = 0;
61   //  layClusters->SetOwnIds(kTRUE);
62
63   Int_t layOld =-1;
64   cTree->GetEntry(0);
65   for (int ilr=0;ilr<nlr;ilr++) {
66     TClonesArray* clr = (TClonesArray*)layerClus.At(ilr);
67     int ncl = clr->GetEntries();
68     //      AliITSUSegmentationPix* segm = (AliITSUSegmentationPix*)fGM->GetSegmentation(ilr);
69     Float_t maxRsqr = maxR*maxR;
70     for (Int_t icl = 0; icl < ncl; ++icl) {
71       AliITSUClusterPix *c = (AliITSUClusterPix*) clr->UncheckedAt(icl);
72       Int_t id = c->GetVolumeId();
73       int lay,sta,ssta,mod,chip;
74       gm->GetChipId(id, lay,sta,ssta,mod,chip);
75       
76       if (lay!=layOld) { // assumes order in the digits !
77         layOld=lay;
78         layClusters = new TEvePointSet(Form("ITS clusters - Layer %d",lay),10000);
79         //  layClusters->ApplyVizTag(viz_tag, "Clusters");
80         layClusters->SetMarkerColor(kBlue);
81         layClusters->SetMarkerStyle(21);//kFullDotMedium);
82         layClusters->SetRnrSelf(kTRUE);
83         layClusters->SetOwnIds(kTRUE);
84         
85         evClusters->AddElement(layClusters);   
86       }
87       
88       float g[3]; 
89       c->GetGlobalXYZ(g);
90       //      gm->LocalToGlobal(mod,l,g);
91       if (g[0]*g[0] + g[1]*g[1] < maxRsqr) {
92         layClusters->SetNextPoint(g[0], g[1], g[2]);
93         AliITSUClusterPix *atp = new AliITSUClusterPix(*c);
94         layClusters->SetPointId(atp);
95 //      printf("%d: mod %d: loc(%.4lf,%.4lf,%.4lf); glob(%.4lf,%.4lf,%.4lf); \n",
96 //                     icl,atp->GetVolumeId(), l[0],l[1],l[2], g[0],g[1],g[2]);
97       }
98     }
99   }
100   
101   rl->UnloadRecPoints("ITS");
102   
103   //  clusters->SetTitle(Form("N=%d", clusters->Size()));
104   
105   const TString viz_tag("REC Clusters ITS");
106   
107   evClusters->ApplyVizTag(viz_tag, "Clusters");
108   
109   gEve->AddElement(evClusters, cont);
110   
111   gEve->Redraw3D();
112
113 }