]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EVE/alice-macros/itsU_digits.C
Overlaps corrected, new shape of sectors
[u/mrichter/AliRoot.git] / EVE / alice-macros / itsU_digits.C
1 // $Id: its_digits.C 30728 2009-01-22 18:14:34Z 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 #if !defined(__CINT__) || defined(__MAKECINT__)
11 #include <vector> // bypass a cint problem in root-5.20
12 #include <TSystem.h>
13 #include <TEveManager.h>
14 #include <TGeoManager.h>
15 #include <TEveElement.h>
16 #include <TEvePointSet.h>
17 #include <TObjArray.h>
18 #include <TClonesArray.h>
19 #include <TTree.h>
20
21 #include "AliEveEventManager.h"
22 #include "AliGeomManager.h"
23 #include "../ITS/UPGRADE/AliITSUGeomTGeo.h"
24 #include "../ITS/UPGRADE/AliITSUSegmentationPix.h"
25 #include "../ITS/UPGRADE/AliITSUDigitPix.h"
26 #include "AliRunLoader.h"
27 #include "AliEveITSUModule.h"
28
29 #else
30
31 class TEveElement;
32 class TEvePointSet;
33
34 #endif
35
36
37 void itsU_digits() 
38   // Int_t mode            = 63,
39   // Bool_t check_empty    = kTRUE,
40   // Bool_t scaled_modules = kFALSE)
41 {
42  
43   gSystem->Load("libITSUpgradeBase");
44   gSystem->Load("libITSUpgradeSim");
45
46   gGeoManager = gEve->GetGeometry("geometry.root");
47   //  if (AliGeomManager::GetGeometry() == 0) AliGeomManager::LoadGeometry("geometry.root");
48
49   AliITSUGeomTGeo* gm = new AliITSUGeomTGeo(kTRUE);
50   TObjArray segmArr;
51   AliITSUSegmentationPix::LoadSegmentations(&segmArr, AliITSUGeomTGeo::GetITSsegmentationFileName());
52
53   //
54   //  Int_t nLayers = gm->GetNLayers();
55   Int_t nModules = gm->GetNModules();
56
57   AliEveEventManager::AssertGeometry();
58   AliRunLoader* rl =  AliEveEventManager::AssertRunLoader();
59   rl->LoadDigits("ITS");
60   TTree* digTree = rl->GetTreeD("ITS", false);
61
62   //  printf(" -------------------- event %d - nModules %d \n",rl->GetEventNumber(), nModules);
63
64   //DIGITS INIT
65
66   TClonesArray *digArr= new TClonesArray("AliITSUDigitPix");
67   digTree->SetBranchAddress("ITSDigitsPix",&digArr);
68
69   TEveElementList* evITSdig = new TEveElementList("ITSU Digits - as modules");
70   TEveElementList* layITSdig = 0;
71
72   TEveElementList* evITSdigPoints = new TEveElementList("ITSU Digits - as points");
73   TEvePointSet* layITSdigPoints = 0;
74   
75   Int_t layOld =-1;
76
77   for (int imod=0;imod<nModules;imod++) {
78
79     digTree->GetEntry(imod);     
80     // clone array, because digArr is reused for the next module
81     // and the digit pointers for prev. modules are lost
82     TClonesArray *digArrClone =  (TClonesArray*)digArr->Clone("digits clone");
83  
84     int ndig  = digArrClone->GetEntries();
85     if (ndig<1) continue; 
86
87     int detType = gm->GetModuleDetTypeID(imod);
88     AliITSUSegmentationPix* segm = (AliITSUSegmentationPix*)segmArr.At(detType);
89     int lay,lad,det;
90     gm->GetModuleId(imod, lay,lad,det);
91     //   printf("\nModule %3d: (det %2d in ladder %2d of Layer %d) | NDigits: %4d\n",imod,det,lad,lay,ndig);
92     //
93     
94     if (lay!=layOld) { // assumes order in the digits !
95       layITSdigPoints = new TEvePointSet(Form("ITSU Digits (points) - Layer %d",lay),10000);
96       //  layITSdigPoints->ApplyVizTag(viz_tag, "Clusters");
97       layITSdigPoints->SetMarkerColor(kBlack);
98       layITSdigPoints->SetMarkerStyle(kFullDotMedium);
99       layITSdigPoints->SetRnrSelf(kTRUE);
100       layITSdigPoints->SetOwnIds(kTRUE);
101
102       evITSdigPoints->AddElement(layITSdigPoints);   
103
104     }
105
106     AliEveITSUModule *evMod  = new AliEveITSUModule(gm,imod,lay,lad,det);
107  
108     for (int idig=0;idig<ndig;idig++) {
109    
110       AliITSUDigitPix *pDig = (AliITSUDigitPix*)digArrClone->At(idig);
111       evMod->SetDigitInQuad(pDig);
112   
113       /* printf("#%3d digit (0x%lx), col:%4d/row:%3d signal: %5d e-,  generated by tracks ",
114              idig,(ULong_t)pDig,pDig->GetCoord1(),pDig->GetCoord2(),pDig->GetSignalPix()); 
115       for (int itr=0;itr<pDig->GetNTracks();itr++) if (pDig->GetTrack(itr)>=0) printf(" %5d",pDig->GetTrack(itr)); printf("\n");
116       */
117    
118       // extract global coordinates
119       Float_t x,z;
120       segm->DetToLocal(pDig->GetCoord2(),pDig->GetCoord1(),x,z);
121
122       Double_t loc[3]={x,0,z}; 
123       Double_t glob[3]; 
124       gm->LocalToGlobal(imod,loc,glob);
125    
126       layITSdigPoints->SetNextPoint(glob[0], glob[1], glob[2]);
127       layITSdigPoints->SetPointId(pDig);
128    
129       //   printf("  loc(%.3lf,%.3lf,%.3lf)->glob(%.3lf,%.3lf,%.3lf) \n",loc[0],loc[1],loc[2],glob[0],glob[1],glob[2]);
130     }
131
132  
133     if (lay!=layOld) { // assumes order in the digits !
134       if (layOld>=0) {
135         evITSdig->AddElement(layITSdig);   
136       }
137       layITSdig = new TEveElementList(Form("ITSU Digits - Layer %d",lay));  
138       layOld=lay;
139
140     }
141
142     layITSdig->AddElement(evMod);
143   
144   }
145   evITSdig->AddElement(layITSdig);     // add list of digits of last layer
146
147   //  gEve->AddElement(evITSdigPoints);
148   gEve->AddElement(evITSdig);
149
150   gEve->Redraw3D();
151 }