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 **************************************************************************/
8 #include <AliEveITSUModule.h>
9 #include <TGeoMatrix.h>
10 #include <TEveTrans.h>
11 #include <TClonesArray.h>
15 #include <TEveManager.h>
16 #include <TGeoManager.h>
17 #include <AliEveEventManager.h>
18 #include <AliGeomManager.h>
20 //#include "UPGRADE/AliITSUDigitPix.h"
22 #include <AliITSUGeomTGeo.h>
23 #include <AliITSsegmentation.h>
24 #include <AliITSUDigitPix.h>
25 //______________________________________________________________________________
27 // Visualization of an ITS Upgrade module.
29 ClassImp(AliEveITSUModule)
31 Bool_t AliEveITSUModule::fgStaticInitDone = 0;
33 TEveFrameBox* AliEveITSUModule::fgITSUFrameBox = 0;
34 TEveFrameBox* AliEveITSUModule::fgITSUFrameBoxDead = 0;
35 TEveRGBAPalette* AliEveITSUModule::fgITSUPalette = 0;
37 AliITSUGeomTGeo* fGM = 0;
38 const AliITSsegmentation* fSegm = 0;
40 /******************************************************************************/
42 AliEveITSUModule::AliEveITSUModule(const Text_t* n, const Text_t* t) :
49 fAtLeastOneDigit(kFALSE)
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)),
62 fAtLeastOneDigit(kFALSE)
68 fGM = gm; // ITSU Geometry Manager
69 fgStaticInitDone = kFALSE;
70 fSegm = fGM->GetSegmentation(layer);
71 fDpx = fSegm->Dpx(0); // pixel pitch in x
72 fDpz = fSegm->Dpz(0); // pixel pitch in z
77 AliEveITSUModule::~AliEveITSUModule()
84 /******************************************************************************/
86 void AliEveITSUModule::InitStatics()
88 // Initialize static variables.
90 // Warning all sensor sizes are cm
91 // In Eve half-lengths/widths are used, hence 1/2.
93 if (fgStaticInitDone) return;
94 fgStaticInitDone = kTRUE;
96 Float_t dx = fSegm->Dx(); // dimension in x in cm
97 Float_t dz = fSegm->Dz(); // dimension in y in cm
98 Float_t dy = 0;// ? eventuelly a few 100 micron, right?
101 fgITSUFrameBox = new TEveFrameBox();
102 fgITSUFrameBox->SetAAQuadXZ(-dx/2, dy, -dz/2, dx, dz);
103 fgITSUFrameBox->SetFrameColor(kBlue-4);
104 fgITSUFrameBox->SetFrameFill(kTRUE);
105 fgITSUFrameBox->IncRefCount();
107 fgITSUPalette = new TEveRGBAPalette(0,1);
108 fgITSUPalette->IncRefCount();
110 fgITSUFrameBoxDead = new TEveFrameBox();
111 fgITSUFrameBoxDead->SetAAQuadXZ(-dx/2, dy, -dz/2, dx, dz);
112 fgITSUFrameBoxDead->SetFrameColor(kRed);
113 fgITSUFrameBoxDead->SetFrameFill(kTRUE);
114 fgITSUFrameBoxDead->IncRefCount();
121 /******************************************************************************/
123 void AliEveITSUModule::SetID(Int_t gid, Bool_t trans)
127 static const TEveException kEH("AliEveITSUModule::SetID ");
130 if (!fgStaticInitDone)
135 SetFrame(fgITSUFrameBox);
136 SetPalette(fgITSUPalette);
146 /******************************************************************************/
148 void AliEveITSUModule::SetDigitInQuad(AliITSUDigitPix *pDig)
150 // Sets a digit from source in a visualization structure - called quads.
153 if (!fAtLeastOneDigit) {
154 Reset(kQT_RectangleXZFixedY, kFALSE, 32);
155 fAtLeastOneDigit = kTRUE;
159 fSegm->DetToLocal(pDig->GetCoord2(),pDig->GetCoord1(),x,z);
161 AddQuad(x-fDpx/2, z-fDpz/2, fDpx, fDpz);
164 Int_t intSignal = pDig->GetSignalPix();
165 QuadValue(intSignal);
166 if (fgITSUPalette->GetMaxVal()<intSignal) {
167 fgITSUPalette->SetMax(intSignal);
168 fgITSUPalette->MinMaxValChanged();
174 /******************************************************************************/
176 void AliEveITSUModule::SetTrans()
178 // Set transformation matrix
180 const TGeoHMatrix *mat = fGM->GetMatrixSens(fID);
181 fMainTrans->SetFrom(*mat);
186 /******************************************************************************/
188 void AliEveITSUModule::Print(Option_t* ) const
190 // Print object summary information.
192 printf("AliEveITSUModule: ModuleId: %d, layer %d, ladder %d, detector %d\n",
193 fID, fkLayer, fkLadder, fkDetector);
197 /******************************************************************************/
199 void AliEveITSUModule::DigitSelected(Int_t idx)
201 // Override secondary select (alt-click) from TEveQuadSet.
203 // for (Int_t i=0;i<7;i++) {
205 DigitBase_t *qb = GetDigit(idx);
206 TObject *obj = GetId(idx);
207 AliITSUDigitPix *pDig = dynamic_cast<AliITSUDigitPix*>(obj);
208 printf("AliEveITSUModule::QuadSelected ");
209 printf(" idx=%d, value=%d, obj=0x%lx, digit=0x%lx\n",
210 idx, qb->fValue, (ULong_t)obj, (ULong_t)pDig);
213 fSegm->DetToLocal(pDig->GetCoord2(),pDig->GetCoord1(),x,z);
214 printf(" Digit info: mod|lay/lad/det=%d|%d/%d/%d; row/col=%3d/%4d; \n",
215 fID,fkLayer,fkLadder,fkDetector,
216 pDig->GetCoord2(),pDig->GetCoord1());
217 printf(" local (x,z)=(%.4lf,%.4lf)cm; signal:%5d e-; generated by tracks ",
218 x,z,pDig->GetSignalPix());
219 for (int itr=0;itr<pDig->GetNTracks();itr++)
220 if (pDig->GetTrack(itr)>=0) printf(" %5d",pDig->GetTrack(itr)); printf("\n");