]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/AliEMCALv2.cxx
AddAlignableVolumes() method added for EMCAL
[u/mrichter/AliRoot.git] / EMCAL / AliEMCALv2.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-2004, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16 /* $Id$ */
17
18 //_________________________________________________________________________
19 //*-- Implementation version v2 of EMCAL Manager class; SHASHLYK version
20 //*-- An object of this class does not produce digits
21 //*-- It is the one to use if you do want to produce outputs in TREEH 
22 //*--                  
23 //*-- Author : Aleksei Pavlinov (WSU)
24
25 // This Class not stores information on all particles prior to EMCAL entry - in order to facilitate analysis.
26 // This is done by setting fIShunt =2, and flagging all parents of particles entering the EMCAL.
27
28 // --- ROOT system ---
29 #include <TParticle.h>
30 #include <TVirtualMC.h>
31 #include <TBrowser.h>
32 #include <TH1.h>
33 #include <TH2.h>
34
35 // --- Standard library ---
36
37 // --- AliRoot header files ---
38 #include "AliEMCALv2.h"
39 #include "AliEMCALHit.h"
40 #include "AliEMCALGeometry.h"
41 #include "AliRun.h"
42 #include "AliHeader.h"
43 #include "AliMC.h"
44 #include "AliPoints.h"
45 // for TRD1 case only; May 31,2006
46
47 ClassImp(AliEMCALv2)
48
49 //______________________________________________________________________
50 AliEMCALv2::AliEMCALv2():AliEMCALv1(), fGeometry(0){
51   // ctor
52 }
53
54 //______________________________________________________________________
55 AliEMCALv2::AliEMCALv2(const char *name, const char *title): AliEMCALv1(name,title) {
56     // Standard Creator.
57
58     fHits= new TClonesArray("AliEMCALHit",1000);
59     gAlice->GetMCApp()->AddHitList(fHits);
60
61     fNhits    = 0;
62     fIshunt   = 2; // All hits are associated with particles entering the calorimeter
63     fTimeCut  = 30e-09;
64
65     fGeometry = GetGeometry(); 
66     fHDe = fHNhits = 0;
67     //    if (gDebug>0){
68     if (1){
69       TH1::AddDirectory(0);
70       fHDe = new TH1F("fHDe","De in EMCAL", 1000, 0., 1.);
71       fHNhits = new TH1F("fHNhits","#hits in EMCAL", 2001, -0.5, 2000.5);
72       fHistograms->Add(fHDe);
73       fHistograms->Add(fHNhits);
74       TH1::AddDirectory(1);
75     }
76 }
77
78 //______________________________________________________________________
79 AliEMCALv2::AliEMCALv2(const AliEMCALv2 & emcal):AliEMCALv1(emcal)
80 {
81   //copy ctor
82   fGeometry = emcal.fGeometry;
83   fHDe = emcal.fHDe;
84   fHNhits = emcal.fHNhits;
85
86 }
87
88 //______________________________________________________________________
89 AliEMCALv2::~AliEMCALv2(){
90     // dtor
91
92     if ( fHits) {
93         fHits->Delete();
94         delete fHits;
95         fHits = 0;
96     }
97 }
98
99 //______________________________________________________________________
100 void AliEMCALv2::AddHit(Int_t shunt, Int_t primary, Int_t tracknumber, Int_t iparent, Float_t ienergy, 
101                         Int_t id, Float_t * hits,Float_t * p){
102     // Add a hit to the hit list.
103     // An EMCAL hit is the sum of all hits in a tower section
104     //   originating from the same entering particle 
105     static Int_t hitCounter;
106     static AliEMCALHit *newHit, *curHit;
107     static Bool_t deja;
108
109     deja = kFALSE;
110
111     newHit = new AliEMCALHit(shunt, primary, tracknumber, iparent, ienergy, id, hits, p);
112     for ( hitCounter = fNhits-1; hitCounter >= 0 && !deja; hitCounter-- ) {
113         curHit = (AliEMCALHit*) (*fHits)[hitCounter];
114         // We add hits with the same tracknumber, while GEANT treats
115         // primaries succesively
116         if(curHit->GetPrimary() != primary) 
117           break;
118         if( *curHit == *newHit ) {
119             *curHit = *curHit + *newHit;
120             deja = kTRUE;
121             //            break; // 30-aug-04 by PAI 
122         } // end if
123     } // end for hitCounter
124     
125     if ( !deja ) {
126         new((*fHits)[fNhits]) AliEMCALHit(*newHit);
127         fNhits++;
128     }
129     //    printf(" fNhits %i \n", fNhits); 
130     delete newHit;
131 }
132 //______________________________________________________________________
133 void AliEMCALv2::StepManager(void){
134   // Accumulates hits as long as the track stays in a tower
135
136   // position wrt MRS and energy deposited
137   static Float_t        xyzte[5]={0.,0.,0.,0.,0.};// position wrt MRS, time and energy deposited
138   static Float_t        pmom[4]={0.,0.,0.,0.};
139   static TLorentzVector pos; // Lorentz vector of the track current position.
140   static TLorentzVector mom; // Lorentz vector of the track current momentum.
141   static Float_t ienergy = 0;
142   static TString curVolName;
143   static int supModuleNumber, moduleNumber, yNumber, xNumber, absid;
144   static int keyGeom=1;
145   static char *vn = "SCMX"; // Apr 13, 2006 - only TRD1 case now
146   //  static char *vn = "SX"; // 15-mar-05
147   static int nSMOP[7]={1,3,5,7,9,11}; // 30-mar-05
148   static int nSMON[7]={2,4,6,8,10,12};
149   static Float_t depositedEnergy=0.0; 
150
151   if(keyGeom == 0) {
152     keyGeom = 2;
153     if(gMC->VolId("PBMO")==0 || gMC->VolId("WSUC")==1) {
154       vn      = "SCMX";   // old TRD2(TRD1) or WSUC
155       keyGeom = 1;
156     }    
157     printf("AliEMCALv2::StepManager():  keyGeom %i : Sensetive volume %s \n", 
158     keyGeom, vn); 
159     if(gMC->VolId("WSUC")==1) printf(" WSUC - cosmic ray stand geometry \n");
160   }
161   Int_t tracknumber =  gAlice->GetMCApp()->GetCurrentTrackNumber();
162
163   curVolName = gMC->CurrentVolName();
164   if(curVolName.Contains(vn)) { // We are in a scintillator layer
165     //    printf(" keyGeom %i : Sensetive volume %s (%s) \n", keyGeom, curVolName.Data(), vn); 
166     
167     if( ((depositedEnergy = gMC->Edep()) > 0.)  && (gMC->TrackTime() < fTimeCut)){// Track is inside a scintillator and deposits some energy
168       //       Info("StepManager "," entry %i DE %f",++ientry, depositedEnergy); // for testing
169        if (fCurPrimary==-1) 
170         fCurPrimary=gAlice->GetMCApp()->GetPrimary(tracknumber);
171
172       if (fCurParent==-1 || tracknumber != fCurTrack) {
173         // Check parentage
174         Int_t parent=tracknumber;
175         if (fCurParent != -1) {
176           while (parent != fCurParent && parent != -1) {
177             TParticle *part=gAlice->GetMCApp()->Particle(parent);
178             parent=part->GetFirstMother();
179           }
180         }
181         if (fCurParent==-1 || parent==-1) {
182           Int_t parent=tracknumber;
183           TParticle *part=gAlice->GetMCApp()->Particle(parent);
184           while (parent != -1 && fGeometry->IsInEMCAL(part->Vx(),part->Vy(),part->Vz())) {
185             parent=part->GetFirstMother();
186             if (parent!=-1) 
187               part=gAlice->GetMCApp()->Particle(parent);
188           } 
189           fCurParent=parent;
190           if (fCurParent==-1)
191             Error("StepManager","Cannot find parent");
192           else {
193             TParticle *part=gAlice->GetMCApp()->Particle(fCurParent);
194             ienergy = part->Energy(); 
195           }
196           while (parent != -1) {
197             part=gAlice->GetMCApp()->Particle(parent);
198             part->SetBit(kKeepBit);
199             parent=part->GetFirstMother();
200           }
201         }
202         fCurTrack=tracknumber;
203       }    
204       gMC->TrackPosition(pos);
205       xyzte[0] = pos[0];
206       xyzte[1] = pos[1];
207       xyzte[2] = pos[2];
208       xyzte[3] = gMC->TrackTime() ;       
209       
210       gMC->TrackMomentum(mom);
211       pmom[0] = mom[0];
212       pmom[1] = mom[1];
213       pmom[2] = mom[2];
214       pmom[3] = mom[3];
215       
216       //      if(ientry%200 > 0) return; // testing
217       supModuleNumber = moduleNumber = yNumber = xNumber = absid = 0;
218       if(keyGeom >= 1) { // old style
219         gMC->CurrentVolOffID(4, supModuleNumber);
220         gMC->CurrentVolOffID(3, moduleNumber);
221         gMC->CurrentVolOffID(1, yNumber);
222         gMC->CurrentVolOffID(0, xNumber); // really x number now
223         if(strcmp(gMC->CurrentVolOffName(4),"SM10")==0) supModuleNumber += 10; // 13-oct-05
224       } else {
225         gMC->CurrentVolOffID(5, supModuleNumber);
226         gMC->CurrentVolOffID(4, moduleNumber);
227         gMC->CurrentVolOffID(1, yNumber);
228         gMC->CurrentVolOffID(0, xNumber);
229         if     (strcmp(gMC->CurrentVolOffName(5),"SMOP")==0) supModuleNumber = nSMOP[supModuleNumber-1];
230         else if(strcmp(gMC->CurrentVolOffName(5),"SMON")==0) supModuleNumber = nSMON[supModuleNumber-1];
231         else   assert(0); // something wrong
232       }
233       absid = fGeometry->GetAbsCellId(supModuleNumber-1, moduleNumber-1, yNumber-1, xNumber-1);
234     
235       if (absid < 0) Fatal("StepManager()", "Wrong id ") ;
236
237       Float_t lightYield =  depositedEnergy ;
238       // Apply Birk's law (copied from G3BIRK)
239
240       if (gMC->TrackCharge()!=0) { // Check
241           Float_t birkC1Mod = 0;
242         if (fBirkC0==1){ // Apply correction for higher charge states
243           if (TMath::Abs(gMC->TrackCharge())>=2) birkC1Mod = fBirkC1*7.2/12.6;
244           else                                   birkC1Mod = fBirkC1;
245         }
246
247         Float_t dedxcm;
248         if (gMC->TrackStep()>0)  dedxcm=1000.*gMC->Edep()/gMC->TrackStep();
249         else                     dedxcm=0;
250         lightYield=lightYield/(1.+birkC1Mod*dedxcm+fBirkC2*dedxcm*dedxcm);
251       } 
252
253       // use sampling fraction to get original energy --HG
254       //      xyzte[4] = lightYield * fGeometry->GetSampling();
255       xyzte[4] = lightYield; // 15-dec-04
256         
257       if (gDebug == -2) 
258       printf("#sm %2i #m %3i #x %1i #z %1i -> absid %i : xyzte[4] = %f\n",
259       supModuleNumber,moduleNumber,yNumber,xNumber,absid, xyzte[4]);
260
261       AddHit(fIshunt, fCurPrimary,tracknumber, fCurParent, ienergy, absid,  xyzte, pmom);
262     } // there is deposited energy
263   }
264 }
265
266 void AliEMCALv2::FinishEvent()
267
268   // Calculate deposit energy and fill control histogram; 26-may-05
269   static double de=0.;
270   fHNhits->Fill(double(fHits->GetEntries()));
271   de = GetDepositEnergy(0);
272   if(fHDe) fHDe->Fill(de);
273 }
274
275 Double_t AliEMCALv2::GetDepositEnergy(int print)
276
277   // 23-mar-05 - for testing
278   if(fHits == 0) return 0.;
279   AliEMCALHit  *hit=0;
280   Double_t de=0.;
281   for(int ih=0; ih<fHits->GetEntries(); ih++) {
282     hit = (AliEMCALHit*)fHits->UncheckedAt(ih);
283     de += hit->GetEnergy();
284   }
285   if(print>0) {
286     cout<<"AliEMCALv2::GetDepositEnergy() : fHits "<<fHits<<endl; 
287     printf(" #hits %i de %f \n", fHits->GetEntries(), de);
288     if(print>1) {
289       printf(" #primary particles %i\n", gAlice->GetHeader()->GetNprimary()); 
290     }
291   }
292   return de;
293 }
294
295 void AliEMCALv2::Browse(TBrowser* b)
296 {
297   TObject::Browse(b);
298 }
299
300 void AliEMCALv2::DrawCalorimeterCut(const char *name, int axis, double dcut)
301
302   // Size of tower is 5.6x5.6x24.8 (25.0); cut on Z axiz
303   TString g(fGeometry->GetName());
304   g.ToUpper();
305   gMC->Gsatt("*", "seen", 0);
306
307   int fill = 1;
308
309   if(axis!=1) SetVolumeAttributes("STPL", 1, 1, fill);
310
311   int colo=4;
312   TString sn(name);
313   if(sn.Contains("SCM")) colo=5;
314   SetVolumeAttributes(name, 1, colo, fill);
315   if(g.Contains("110DEG") && sn=="SMOD") SetVolumeAttributes("SM10", 1, colo, fill);
316
317   TString st(GetTitle());
318   st += ", zcut, ";
319   st += name;
320
321   char *optShad = "on", *optHide="on";
322   double cxy=0.02;
323   if     (axis==1) {
324     dcut = 0.;
325     //    optHide = optShad = "off";
326   } else if(axis==2){
327     if(dcut==0.) SetVolumeAttributes("SCM0", 1, 5, fill); // yellow    
328   }
329
330   gMC->Gdopt("hide", optHide);
331   gMC->Gdopt("shad", optShad);
332
333   gROOT->ProcessLine("TGeant3 *g3 = (TGeant3*)gMC");
334   char cmd[200];
335   sprintf(cmd,"g3->Gdrawc(\"XEN1\", %i, %5.1f, 12., 8., %3.2f, %3.2f)\n", axis, dcut, cxy,cxy);
336   printf("\n%s\n",cmd); gROOT->ProcessLine(cmd);
337
338   sprintf(cmd,"gMC->Gdhead(1111, \"%s\")\n", st.Data());
339   printf("%s\n",cmd); gROOT->ProcessLine(cmd);
340 }
341
342 void AliEMCALv2::DrawSuperModuleCut(const char *name, int axis, double dcut, int fill)
343
344  // Size of tower is 5.6x5.6x24.8 (25.0); cut on Z axiz
345   TString sn(GetGeometry()->GetName());
346   sn.ToUpper();
347   char *tit[3]={"xcut", "ycut", "zcut"};
348   if(axis<1) axis=1; if(axis>3) axis=3;
349
350   gMC->Gsatt("*", "seen", 0);
351   //  int fill = 6;
352
353   // SetVolumeAttributes("SMOD", 1, 1, fill);  // black
354   SetVolumeAttributes(name, 1, 5, fill);    // yellow 
355
356   double cxy=0.055, x0=10., y0=10.;
357   char *optShad = "on", *optHide="on";
358   if(sn.Contains("TRD1")) {
359     SetVolumeAttributes("STPL", 1, 3, fill);  // green 
360     if     (axis==1) {
361       gMC->Gsatt("STPL", "seen", 0);
362       dcut = 0.;
363       optHide = "off";
364       optShad = "off";
365     } else if(axis==3) cxy = 0.1;
366   } else if(sn.Contains("TRD2")) {
367     y0 = -10.;
368     if (axis==2) cxy=0.06;
369   }
370   gMC->Gdopt("hide", optHide);
371   gMC->Gdopt("shad", optShad);
372
373   TString st("Shish-Kebab, Compact, SMOD, ");
374   if(sn.Contains("TWIST")) st = "Shish-Kebab, Twist, SMOD, ";
375   st += tit[axis-1];;
376
377   gROOT->ProcessLine("TGeant3 *g3 = (TGeant3*)gMC");
378   char cmd[200];
379   sprintf(cmd,"g3->Gdrawc(\"SMOD\", %i, %6.3f, %5.1f, %5.1f, %5.3f, %5.3f)\n",axis,dcut,x0,y0,cxy,cxy);
380   printf("\n%s\n",cmd); gROOT->ProcessLine(cmd);
381
382   sprintf(cmd,"gMC->Gdhead(1111, \"%s\")\n", st.Data());
383   printf("%s\n",cmd); gROOT->ProcessLine(cmd);
384   // hint for testing
385   if(sn.Contains("TRD1")){
386     printf("Begin of super module\n");
387     printf("g3->Gdrawc(\"SMOD\", 2,  0.300, 89., 10., 0.5, 0.5)\n");
388     printf("Center of super modules\n");
389     printf("g3->Gdrawc(\"SMOD\", 2,  0.300, 0., 10., 0.5, 0.5)\n");
390     printf("end of super modules\n");
391     printf("g3->Gdrawc(\"SMOD\", 2,  0.300, -70., 10., 0.5, 0.5)\n");
392   } else if(sn.Contains("TRD2")){
393     printf("Begin of super module  ** TRD2 ** \n");
394     printf("g3->Gdrawc(\"SMOD\", 2,  0.00,  40., -80, 0.2, 0.2)\n");
395     printf("end of super modules\n");
396     printf("g3->Gdrawc(\"SMOD\", 2,  0.00, -20., -80, 0.2, 0.2)\n");
397
398     printf(" ***  Z cut (Y|X plane)\n scale 0.4\n");
399     printf("g3->Gdrawc(\"SMOD\", 3,  -170.,-165.,10.,0.4,0.4)\n");
400     printf(" scale 0.2\n");
401     printf("g3->Gdrawc(\"SMOD\", 3,  -170.,-80.,10.,0.2,0.2)\n");
402     printf(" scale 0.12\n");
403     printf("g3->Gdrawc(\"SMOD\", 3,   -170.,-45.,10.,0.12,0.12)\n");
404   }
405 }
406
407 void AliEMCALv2::DrawTowerCut(const char *name, int axis, double dcut, int fill, char *optShad)
408
409   // Size of tower is 5.6x5.6x24.8 (25.0); cut on Z axiz
410   if(axis<1) axis=1; if(axis>3) axis=3;
411   TString mn(name); mn.ToUpper();
412   TString sn(GetGeometry()->GetName());
413
414   gMC->Gsatt("*", "seen", 0);
415   gMC->Gsatt("*", "fill", fill);
416
417   // int mainColo=5; // yellow
418   if(mn == "EMOD") {
419     SetVolumeAttributes(mn.Data(), 1, 1, fill);
420     SetVolumeAttributes("SCM0", 1, 5, fill); // yellow
421     SetVolumeAttributes("SCPA", 1, 3, fill); // green - 13-may-05
422   } else if(mn == "SCM0") { // first division 
423     SetVolumeAttributes(mn.Data(), 1, 1, fill);
424     SetVolumeAttributes("SCMY", 1, 5, fill); // yellow
425   } else if(mn == "SCMY") { // first division 
426     SetVolumeAttributes(mn.Data(), 1, 1, fill); 
427     if(sn.Contains("TEST") && sn.Contains("3X3")) {
428       SetVolumeAttributes("SCX1", 1, 5, fill); // yellow
429       SetVolumeAttributes("SCX2", 1, 2, fill); // red
430       SetVolumeAttributes("SCX3", 1, 3, fill); // green
431     } else {
432       SetVolumeAttributes("SCMX", 1, 5, fill); // yellow
433     }
434   } else if(mn == "SCMX" || mn.Contains("SCX")) {
435     SetVolumeAttributes(mn.Data(), 1, 5, fill);// yellow
436     SetVolumeAttributes("PBTI", 1, 4, fill);
437   } else {
438     printf("<W> for volume |%s| did not defined volume attributes\n", mn.Data());
439   }
440
441   //  TString st("Shish-Kebab, 2x2 mm sampling, 62 layers, ");
442   TString st("Shashlyk, 2x2 mm sampling, 62 layers, ");
443   if    (sn.Contains("25"))   st = "Shish-Kebab, 5x5 mm sampling, 25 layers, ";
444   else if(sn.Contains("MAY05"))   st = "Shish-Kebab, 5x5 mm sampling, 77 layers, ";
445   if(sn.Contains("TRD1")) st += " TRD1, ";
446   if(sn.Contains("3X3"))  st += " 3x3, ";
447   st += name;
448
449   gROOT->ProcessLine("TGeant3 *g3 = (TGeant3*)gMC");
450   double cx=0.78, cy=2.;
451   if     (axis==1 && sn.Contains("TRD")==0) {
452    cx = cy = 1.;
453    gMC->Gsatt("PBTI", "seen", 0);
454   } else if (sn.Contains("TEST") && sn.Contains("3X3")) {
455     cx = cy = 0.7;
456     if (axis==3) cx = cy = 1.;
457   } else if (sn.Contains("TRD2")) {
458     if (axis==3) cx = cy = 2.;
459   } else if (mn.Contains("EMOD")) {
460     cx = cy = 0.5;
461   }
462   char cmd[200];
463
464   gMC->Gdopt("hide", optShad);
465   gMC->Gdopt("shad", optShad);
466
467   sprintf(cmd,"g3->Gdrawc(\"%s\", %i, %6.2f, 10., 10., %5.3f, %5.3f)\n",name, axis,dcut,cx,cy);
468   printf("\n%s\n",cmd); gROOT->ProcessLine(cmd);
469
470   sprintf(cmd,"gMC->Gdhead(1111, \"%s\")\n", st.Data());
471   printf("%s\n",cmd); gROOT->ProcessLine(cmd);
472 }
473   
474 void AliEMCALv2::DrawAlicWithHits(int mode)
475
476  // 20-sep-04; does not work now
477   static TH2F *h2;
478   if(h2==0) h2 = new TH2F("h2","test fo hits", 60,0.5,60.5, 28,0.5,28.5);
479   else      h2->Reset();
480   if(mode==0) {
481     gROOT->ProcessLine("TGeant3 *g3 = (TGeant3*)gMC");
482     gMC->Gsatt("*","seen",0);
483     gMC->Gsatt("scm0","seen",5);
484
485     gROOT->ProcessLine("g3->Gdrawc(\"ALIC\", 1,   2.0, 12., -130, 0.3, 0.3)");
486     // g3->Gdrawc("ALIC", 1,   2.0, 10., -14, 0.05, 0.05)
487   }
488
489   TClonesArray *hits = Hits();
490   Int_t nhits = hits->GetEntries(), absId, nSupMod, nTower, nIphi, nIeta, iphi, ieta;
491   AliEMCALHit *hit = 0;
492   Double_t de, des=0.;
493   for(Int_t i=0; i<nhits; i++) {
494     hit   = (AliEMCALHit*)hits->UncheckedAt(i);
495     absId = hit->GetId();
496     de    = hit->GetEnergy();
497     des += de;
498     if(fGeometry->GetCellIndex(absId, nSupMod, nTower, nIphi, nIeta)){
499       //      printf(" de %f abs id %i smod %i tower %i | cell iphi %i : ieta %i\n",
500       // de, absId, nSupMod, nTower, nIphi, nIeta); 
501       if(nSupMod==3) {
502         fGeometry->GetCellPhiEtaIndexInSModule(nSupMod,nTower,nIphi,nIeta, iphi,ieta);
503         // printf(" iphi %i : ieta %i\n", iphi,ieta);
504         h2->Fill(double(ieta),double(iphi), de);
505       }
506     }
507     //    hit->Dump();
508   }
509   printf(" #hits %i : sum de %f -> %f \n", nhits, des, h2->Integral());
510   if(mode>0 && h2->Integral()>0.) h2->Draw();
511 }
512
513 void AliEMCALv2::SetVolumeAttributes(const char *name, int seen, int color, int fill)
514 {
515  /* seen=-2:volume is visible but none of its descendants;
516          -1:volume is not visible together with all its descendants;
517           0:volume is not visible;
518           1:volume is     visible. */
519   gMC->Gsatt(name, "seen", seen);
520   gMC->Gsatt(name, "colo", color);
521   gMC->Gsatt(name, "fill", fill); 
522   printf(" %s : seen %i color %i fill %i \n", name, seen, color, fill);
523
524
525 void AliEMCALv2::TestIndexTransition(int pri, int idmax)
526
527  // Test for EMCAL_SHISH geometry
528   TString sn(fGeometry->GetName());
529   if(!sn.Contains("SHISH")) {
530     printf("Wrong geometry |%s| ! Bye \n", sn.Data());
531     return; 
532   }
533
534   Int_t nSupMod, nTower, nIphi, nIeta, idNew, nGood=0;
535   if(idmax==0) idmax = fGeometry->GetNCells();
536   for(Int_t id=1; id<=idmax; id++) {
537     if(!fGeometry->GetCellIndex(id, nSupMod, nTower, nIphi, nIeta)){
538       printf(" Wrong abs ID %i : #cells %i\n", id, fGeometry->GetNCells());
539       break;  
540     }
541     idNew = fGeometry->GetAbsCellId(nSupMod, nTower, nIphi, nIeta);
542     if(id != idNew || pri>0) {
543       printf(" ID %i : %i <- new id\n", id, idNew);
544       printf(" nSupMod   %i ",  nSupMod);
545       printf(" nTower    %i ", nTower);
546       printf(" nIphi     %i ", nIphi);
547       printf(" nIeta     %i \n", nIeta);
548      
549     } else nGood++;
550   }
551   printf(" Good decoding %i : %i <- #cells \n", nGood, fGeometry->GetNCells());
552 }