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