]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EVE/EveDet/AliEveITSUModule.cxx
Added ITSUModule (Stefan R.)
[u/mrichter/AliRoot.git] / EVE / EveDet / AliEveITSUModule.cxx
1
2 /**************************************************************************
3  * Copyright(c) 1998-2008, ALICE Experiment at CERN, all rights reserved. *
4  * See http://aliceinfo.cern.ch/Offline/AliRoot/License.html for          *
5  * full copyright notice.                                                 *
6  **************************************************************************/
7
8 #include <AliEveITSUModule.h>
9 #include <TGeoMatrix.h>
10 #include <TEveTrans.h>
11 #include <TClonesArray.h>
12 #include <TStyle.h>
13 #include "AliLog.h"
14
15 #include <TEveManager.h>
16 #include <TGeoManager.h>
17 #include <AliEveEventManager.h>
18 #include <AliGeomManager.h>
19
20 //#include "UPGRADE/AliITSUDigitPix.h"
21
22 #include <AliITSUGeomTGeo.h>
23 #include <AliITSUSegmentationPix.h>
24 #include <AliITSUDigitPix.h>
25 //______________________________________________________________________________
26 //
27 // Visualization of an ITS Upgrade module.
28
29 ClassImp(AliEveITSUModule)
30
31 Bool_t AliEveITSUModule::fgStaticInitDone = 0;
32
33 TEveFrameBox*    AliEveITSUModule::fgITSUFrameBox     = 0;
34 TEveFrameBox*    AliEveITSUModule::fgITSUFrameBoxDead = 0;
35 TEveRGBAPalette* AliEveITSUModule::fgITSUPalette  = 0;
36
37 AliITSUGeomTGeo* fGM                 = 0;
38 AliITSUSegmentationPix * fSegm       = 0;
39
40 /******************************************************************************/
41
42 AliEveITSUModule::AliEveITSUModule(const Text_t* n, const Text_t* t) :
43   TEveQuadSet(n, t),
44   fID(0),
45   fkLayer(0),
46   fkLadder(0),
47   fkDetector(0),
48   fDpx(0), fDpz(0),
49   fAtLeastOneDigit(kFALSE)
50 {
51   // Constructor.
52   
53 }
54
55 AliEveITSUModule::AliEveITSUModule(AliITSUGeomTGeo *gm, Int_t id, Int_t layer, Int_t ladder, Int_t detector) :
56   TEveQuadSet(Form("ITSU module %d; (lay,lad,det)=(%d,%d,%d)", id,layer,ladder,detector),Form("%d",id)),
57   fID(id),
58   fkLayer(layer),
59   fkLadder(ladder),
60   fkDetector(detector),
61   fDpx(0), fDpz(0),
62   fAtLeastOneDigit(kFALSE)
63 {
64
65   // 
66   // constructor
67   //
68
69   fGM = gm; // ITSU Geometry Manager
70
71   fgStaticInitDone = kFALSE; 
72
73   TObjArray segmArr;
74   AliITSUSegmentationPix::LoadSegmentations(&segmArr, AliITSUGeomTGeo::GetITSsegmentationFileName());
75  
76   int detType = fGM->GetModuleDetTypeID(id);
77   fSegm = (AliITSUSegmentationPix*)segmArr.At(detType);
78   
79   fDpx = fSegm->Dpx();  // pixel pitch in x
80   fDpz = fSegm->Dpz(0); // pixel pitch in z
81
82   SetID(id);
83   
84
85 }
86
87 AliEveITSUModule::~AliEveITSUModule()
88 {
89   // Destructor.
90
91  
92 }
93
94 /******************************************************************************/
95
96 void AliEveITSUModule::InitStatics()
97 {
98   // Initialize static variables.
99   //
100   // Warning all sensor sizes are cm
101   // In Eve half-lengths/widths are used, hence 1/2.
102
103   if (fgStaticInitDone) return;
104   fgStaticInitDone = kTRUE;
105
106   Float_t dx =  fSegm->Dpx() * fSegm->GetNRow(); // dimension in x in cm
107   Float_t dz =  fSegm->Dpz(0)* fSegm->GetNCol(); // dimension in y in cm
108   Float_t dy =  0;// ? eventuelly a few 100 micron, right?
109
110   {
111     fgITSUFrameBox = new TEveFrameBox();
112     fgITSUFrameBox->SetAAQuadXZ(-dx/2, dy, -dz/2, dx, dz);
113     fgITSUFrameBox->SetFrameColor(kBlue-4);
114     fgITSUFrameBox->SetFrameFill(kTRUE);
115     fgITSUFrameBox->IncRefCount();
116
117     fgITSUPalette  = new TEveRGBAPalette(0,1);
118     fgITSUPalette->IncRefCount();
119
120     fgITSUFrameBoxDead = new TEveFrameBox();
121     fgITSUFrameBoxDead->SetAAQuadXZ(-dx/2, dy, -dz/2, dx, dz);
122     fgITSUFrameBoxDead->SetFrameColor(kRed);
123     fgITSUFrameBoxDead->SetFrameFill(kTRUE);
124     fgITSUFrameBoxDead->IncRefCount();
125   }
126
127
128 }
129
130
131 /******************************************************************************/
132
133 void AliEveITSUModule::SetID(Int_t gid, Bool_t trans)
134 {
135   // Set detector id.
136
137   static const TEveException kEH("AliEveITSUModule::SetID ");
138
139   fID = gid;
140   if (!fgStaticInitDone)
141   {
142     InitStatics();
143   }
144
145   SetFrame(fgITSUFrameBox);
146   SetPalette(fgITSUPalette);
147
148   RefitPlex();
149   ComputeBBox();
150   InitMainTrans();
151   if (trans)
152     SetTrans();
153
154 }
155
156 /******************************************************************************/
157
158 void AliEveITSUModule::SetDigitInQuad(AliITSUDigitPix *pDig)
159 {
160   // Sets a digit from source in a visualization structure - called quads.
161
162
163   if (!fAtLeastOneDigit) {
164     Reset(kQT_RectangleXZFixedY, kFALSE, 32);
165     fAtLeastOneDigit = kTRUE;
166   }
167
168   Float_t x,z;
169   fSegm->DetToLocal(pDig->GetCoord2(),pDig->GetCoord1(),x,z);
170
171   AddQuad(x-fDpx/2, z-fDpz/2, fDpx, fDpz);
172   QuadId(pDig);
173
174   Int_t intSignal = pDig->GetSignalPix();
175   QuadValue(intSignal);
176   if (fgITSUPalette->GetMaxVal()<intSignal) {
177     fgITSUPalette->SetMax(intSignal);
178     fgITSUPalette->MinMaxValChanged();
179   }
180
181
182 }
183
184 /******************************************************************************/
185
186 void AliEveITSUModule::SetTrans()
187 {
188   // Set transformation matrix 
189    
190   const TGeoHMatrix *mat = fGM->GetMatrixSens(fID);
191   fMainTrans->SetFrom(*mat);
192
193 }
194
195
196 /******************************************************************************/
197
198 void AliEveITSUModule::Print(Option_t* ) const
199 {
200   // Print object summary information.
201
202   printf("AliEveITSUModule: ModuleId: %d, layer %d, ladder %d, detector %d\n",
203          fID, fkLayer, fkLadder, fkDetector);
204
205 }
206
207 /******************************************************************************/
208
209 void AliEveITSUModule::DigitSelected(Int_t idx)
210 {
211   // Override secondary select (alt-click) from TEveQuadSet.
212  
213   //  for (Int_t i=0;i<7;i++) {
214   //    idx=i;
215   DigitBase_t *qb  = GetDigit(idx);
216   TObject     *obj = GetId(idx);
217   AliITSUDigitPix *pDig  = dynamic_cast<AliITSUDigitPix*>(obj);
218   printf("AliEveITSUModule::QuadSelected "); 
219   printf("  idx=%d, value=%d, obj=0x%lx, digit=0x%lx\n",
220          idx, qb->fValue, (ULong_t)obj, (ULong_t)pDig);
221   if (pDig) { 
222     Float_t x,z;
223     fSegm->DetToLocal(pDig->GetCoord2(),pDig->GetCoord1(),x,z);
224     printf(" Digit info: mod|lay/lad/det=%d|%d/%d/%d; row/col=%3d/%4d; \n",
225            fID,fkLayer,fkLadder,fkDetector,
226            pDig->GetCoord2(),pDig->GetCoord1());
227     printf(" local (x,z)=(%.4lf,%.4lf)cm; signal:%5d e-;  generated by tracks ",
228            x,z,pDig->GetSignalPix()); 
229     for (int itr=0;itr<pDig->GetNTracks();itr++)  
230       if (pDig->GetTrack(itr)>=0) printf(" %5d",pDig->GetTrack(itr)); printf("\n");
231   }
232   
233