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