]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/AliEMCALv2.cxx
Removing meaningless const declarations (linuxicc)
[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 #include "AliEMCALGetter.h"
46 // for TRD1,2 case
47
48 ClassImp(AliEMCALv2)
49
50   //extern "C" void gdaxis_(float &x0, float &y0, float &z0, float &axsiz);
51
52
53 //______________________________________________________________________
54 AliEMCALv2::AliEMCALv2():AliEMCALv1(), fGeometry(0){
55   // ctor
56
57 }
58
59 //______________________________________________________________________
60 AliEMCALv2::AliEMCALv2(const char *name, const char *title): AliEMCALv1(name,title) {
61     // Standard Creator.
62
63   //    fGeant3 = (TGeant3*)gMC;
64     fHits= new TClonesArray("AliEMCALHit",1000);
65     gAlice->GetMCApp()->AddHitList(fHits);
66
67     fNhits    = 0;
68     fIshunt   = 2; // All hits are associated with particles entering the calorimeter
69     fTimeCut  = 30e-09;
70
71     fGeometry = GetGeometry(); 
72     fHDe = fHNhits = 0;
73     //    if (gDebug>0){
74     if (1){
75       TH1::AddDirectory(0);
76       fHDe = new TH1F("fHDe","De in EMCAL", 1000, 0., 1.);
77       fHNhits = new TH1F("fHNhits","#hits in EMCAL", 1001, -0.5, 1000.5);
78       fHistograms->Add(fHDe);
79       fHistograms->Add(fHNhits);
80       TH1::AddDirectory(1);
81     }
82 }
83
84 //______________________________________________________________________
85 AliEMCALv2::~AliEMCALv2(){
86     // dtor
87
88     if ( fHits) {
89         fHits->Delete();
90         delete fHits;
91         fHits = 0;
92     }
93 }
94
95 //______________________________________________________________________
96 void AliEMCALv2::AddHit(Int_t shunt, Int_t primary, Int_t tracknumber, Int_t iparent, Float_t ienergy, 
97                         Int_t id, Float_t * hits,Float_t * p){
98     // Add a hit to the hit list.
99     // An EMCAL hit is the sum of all hits in a tower section
100     //   originating from the same entering particle 
101     static Int_t hitCounter;
102     static AliEMCALHit *newHit, *curHit;
103     static Bool_t deja;
104
105     deja = kFALSE;
106
107     newHit = new AliEMCALHit(shunt, primary, tracknumber, iparent, ienergy, id, hits, p);
108     for ( hitCounter = fNhits-1; hitCounter >= 0 && !deja; hitCounter-- ) {
109         curHit = (AliEMCALHit*) (*fHits)[hitCounter];
110         // We add hits with the same tracknumber, while GEANT treats
111         // primaries succesively
112         if(curHit->GetPrimary() != primary) 
113           break;
114         if( *curHit == *newHit ) {
115             *curHit = *curHit + *newHit;
116             deja = kTRUE;
117             //            break; // 30-aug-04 by PAI 
118         } // end if
119     } // end for hitCounter
120     
121     if ( !deja ) {
122         new((*fHits)[fNhits]) AliEMCALHit(*newHit);
123         fNhits++;
124     }
125     //    printf(" fNhits %i \n", fNhits); 
126     delete newHit;
127 }
128 //______________________________________________________________________
129 void AliEMCALv2::StepManager(void){
130   // Accumulates hits as long as the track stays in a tower
131
132   // position wrt MRS and energy deposited
133   static Float_t        xyzte[5]={0.,0.,0.,0.,0.};// position wrt MRS, time and energy deposited
134   static Float_t        pmom[4]={0.,0.,0.,0.};
135   static TLorentzVector pos; // Lorentz vector of the track current position.
136   static TLorentzVector mom; // Lorentz vector of the track current momentum.
137   static Float_t ienergy = 0;
138   static TString curVolName;
139   static int supModuleNumber, moduleNumber, yNumber, xNumber, absid;
140   static int keyGeom=0;
141   static char *vn = "SX"; // 15-mar-05
142   static int nSMOP[7]={1,3,5,7,9,11}; // 30-mar-05
143   static int nSMON[7]={2,4,6,8,10,12};
144   static Float_t depositedEnergy=0.0; 
145
146   if(keyGeom == 0) {
147     keyGeom = 2;
148     if(gMC->VolId("PBMO")==0 || gMC->VolId("WSUC")==1) {
149       vn      = "SCMX";   // old TRD2(TRD1) or WSUC
150       keyGeom = 1;
151     }    
152     printf("AliEMCALv2::StepManager():  keyGeom %i : Sensetive volume %s \n", 
153     keyGeom, vn); 
154     if(gMC->VolId("WSUC")==1) printf(" WSUC - cosmic ray stand geometry \n");
155   }
156   Int_t tracknumber =  gAlice->GetMCApp()->GetCurrentTrackNumber();
157
158   curVolName = gMC->CurrentVolName();
159   if(curVolName.Contains(vn)) { // We are in a scintillator layer
160     //    printf(" keyGeom %i : Sensetive volume %s (%s) \n", keyGeom, curVolName.Data(), vn); 
161     
162     if( ((depositedEnergy = gMC->Edep()) > 0.)  && (gMC->TrackTime() < fTimeCut)){// Track is inside a scintillator and deposits some energy
163       //       Info("StepManager "," entry %i DE %f",++ientry, depositedEnergy); // for testing
164        if (fCurPrimary==-1) 
165         fCurPrimary=gAlice->GetMCApp()->GetPrimary(tracknumber);
166
167       if (fCurParent==-1 || tracknumber != fCurTrack) {
168         // Check parentage
169         Int_t parent=tracknumber;
170         if (fCurParent != -1) {
171           while (parent != fCurParent && parent != -1) {
172             TParticle *part=gAlice->GetMCApp()->Particle(parent);
173             parent=part->GetFirstMother();
174           }
175         }
176         if (fCurParent==-1 || parent==-1) {
177           Int_t parent=tracknumber;
178           TParticle *part=gAlice->GetMCApp()->Particle(parent);
179           while (parent != -1 && fGeometry->IsInEMCAL(part->Vx(),part->Vy(),part->Vz())) {
180             parent=part->GetFirstMother();
181             if (parent!=-1) 
182               part=gAlice->GetMCApp()->Particle(parent);
183           } 
184           fCurParent=parent;
185           if (fCurParent==-1)
186             Error("StepManager","Cannot find parent");
187           else {
188             TParticle *part=gAlice->GetMCApp()->Particle(fCurParent);
189             ienergy = part->Energy(); 
190           }
191           while (parent != -1) {
192             part=gAlice->GetMCApp()->Particle(parent);
193             part->SetBit(kKeepBit);
194             parent=part->GetFirstMother();
195           }
196         }
197         fCurTrack=tracknumber;
198       }    
199       gMC->TrackPosition(pos);
200       xyzte[0] = pos[0];
201       xyzte[1] = pos[1];
202       xyzte[2] = pos[2];
203       xyzte[3] = gMC->TrackTime() ;       
204       
205       gMC->TrackMomentum(mom);
206       pmom[0] = mom[0];
207       pmom[1] = mom[1];
208       pmom[2] = mom[2];
209       pmom[3] = mom[3];
210       
211       //      if(ientry%200 > 0) return; // testing
212       supModuleNumber = moduleNumber = yNumber = xNumber = absid = 0;
213       if(keyGeom >= 1) { // old style
214         gMC->CurrentVolOffID(4, supModuleNumber);
215         gMC->CurrentVolOffID(3, moduleNumber);
216         gMC->CurrentVolOffID(1, yNumber);
217         gMC->CurrentVolOffID(0, xNumber); // really x number now
218         if(strcmp(gMC->CurrentVolOffName(4),"SM10")==0) supModuleNumber += 10; // 13-oct-05
219       } else {
220         gMC->CurrentVolOffID(5, supModuleNumber);
221         gMC->CurrentVolOffID(4, moduleNumber);
222         gMC->CurrentVolOffID(1, yNumber);
223         gMC->CurrentVolOffID(0, xNumber);
224         if     (strcmp(gMC->CurrentVolOffName(5),"SMOP")==0) supModuleNumber = nSMOP[supModuleNumber-1];
225         else if(strcmp(gMC->CurrentVolOffName(5),"SMON")==0) supModuleNumber = nSMON[supModuleNumber-1];
226         else   assert(0); // something wrong
227       }
228       absid = fGeometry->GetAbsCellId(supModuleNumber, moduleNumber, yNumber, xNumber);
229     
230       if (absid == -1) Fatal("StepManager()", "Wrong id ") ;
231
232       Float_t lightYield =  depositedEnergy ;
233       // Apply Birk's law (copied from G3BIRK)
234
235       if (gMC->TrackCharge()!=0) { // Check
236           Float_t BirkC1_mod = 0;
237         if (fBirkC0==1){ // Apply correction for higher charge states
238           if (TMath::Abs(gMC->TrackCharge())>=2) BirkC1_mod=fBirkC1*7.2/12.6;
239           else                                    BirkC1_mod=fBirkC1;
240         }
241
242         Float_t dedxcm;
243         if (gMC->TrackStep()>0)  dedxcm=1000.*gMC->Edep()/gMC->TrackStep();
244         else                     dedxcm=0;
245         lightYield=lightYield/(1.+BirkC1_mod*dedxcm+fBirkC2*dedxcm*dedxcm);
246       } 
247
248       // use sampling fraction to get original energy --HG
249       //      xyzte[4] = lightYield * fGeometry->GetSampling();
250       xyzte[4] = lightYield; // 15-dec-04
251         
252       if (gDebug == -2) 
253       printf("#sm %2i #m %3i #x %1i #z %1i -> absid %i : xyzte[4] = %f\n",
254       supModuleNumber,moduleNumber,yNumber,xNumber,absid, xyzte[4]);
255
256       AddHit(fIshunt, fCurPrimary,tracknumber, fCurParent, ienergy, absid,  xyzte, pmom);
257     } // there is deposited energy
258   }
259 }
260
261 void AliEMCALv2::FinishEvent()
262 { // 26-may-05
263   static double de=0.;
264   fHNhits->Fill(double(fHits->GetEntries()));
265   de = GetDepositEnergy(0);
266   if(fHDe) fHDe->Fill(de);
267 }
268
269 Double_t AliEMCALv2::GetDepositEnergy(int print)
270 { // 23-mar-05 - for testing
271   if(fHits == 0) return 0.;
272   AliEMCALHit  *hit=0;
273   Double_t de=0.;
274   for(int ih=0; ih<fHits->GetEntries(); ih++) {
275     hit = (AliEMCALHit*)fHits->UncheckedAt(ih);
276     de += hit->GetEnergy();
277   }
278   if(print>0) {
279     printf(" #hits %i de %f \n", fHits->GetEntries(), de);
280     if(print>1) {
281       printf(" #primary particles %i\n", gAlice->GetHeader()->GetNprimary()); 
282     }
283   }
284   return de;
285 }
286
287 void AliEMCALv2::Browse(TBrowser* b)
288 {
289   TObject::Browse(b);
290 }
291
292 void AliEMCALv2::DrawCalorimeterCut(const char *name, int axis, double dcut)
293 { // Size of tower is 5.6x5.6x24.8 (25.0); cut on Z axiz
294   TString g(fGeometry->GetName());
295   g.ToUpper();
296   gMC->Gsatt("*", "seen", 0);
297
298   int fill = 1;
299
300   if(axis!=1) SetVolumeAttributes("STPL", 1, 1, fill);
301
302   int colo=4;
303   TString sn(name);
304   if(sn.Contains("SCM")) colo=5;
305   SetVolumeAttributes(name, 1, colo, fill);
306   if(g.Contains("110DEG") && sn=="SMOD") SetVolumeAttributes("SM10", 1, colo, fill);
307
308   TString st(GetTitle());
309   st += ", zcut, ";
310   st += name;
311
312   char *optShad = "on", *optHide="on";
313   double cxy=0.02;
314   if     (axis==1) {
315     dcut = 0.;
316     //    optHide = optShad = "off";
317   } else if(axis==2){
318     if(dcut==0.) SetVolumeAttributes("SCM0", 1, 5, fill); // yellow    
319   }
320
321   gMC->Gdopt("hide", optHide);
322   gMC->Gdopt("shad", optShad);
323
324   gROOT->ProcessLine("TGeant3 *g3 = (TGeant3*)gMC");
325   char cmd[200];
326   sprintf(cmd,"g3->Gdrawc(\"XEN1\", %i, %5.1f, 12., 8., %3.2f, %3.2f)\n", axis, dcut, cxy,cxy);
327   printf("\n%s\n",cmd); gROOT->ProcessLine(cmd);
328
329   sprintf(cmd,"gMC->Gdhead(1111, \"%s\")\n", st.Data());
330   printf("%s\n",cmd); gROOT->ProcessLine(cmd);
331 }
332
333 void AliEMCALv2::DrawSuperModuleCut(const char *name, int axis, double dcut, int fill)
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 { // Size of tower is 5.6x5.6x24.8 (25.0); cut on Z axiz
399   if(axis<1) axis=1; if(axis>3) axis=3;
400   TString mn(name); mn.ToUpper();
401   TString sn(GetGeometry()->GetName());
402
403   gMC->Gsatt("*", "seen", 0);
404   gMC->Gsatt("*", "fill", fill);
405
406   // int mainColo=5; // yellow
407   if(mn == "EMOD") {
408     SetVolumeAttributes(mn.Data(), 1, 1, fill);
409     SetVolumeAttributes("SCM0", 1, 5, fill); // yellow
410     SetVolumeAttributes("SCPA", 1, 3, fill); // green - 13-may-05
411   } else if(mn == "SCM0") { // first division 
412     SetVolumeAttributes(mn.Data(), 1, 1, fill);
413     SetVolumeAttributes("SCMY", 1, 5, fill); // yellow
414   } else if(mn == "SCMY") { // first division 
415     SetVolumeAttributes(mn.Data(), 1, 1, fill); 
416     if(sn.Contains("TEST") && sn.Contains("3X3")) {
417       SetVolumeAttributes("SCX1", 1, 5, fill); // yellow
418       SetVolumeAttributes("SCX2", 1, 2, fill); // red
419       SetVolumeAttributes("SCX3", 1, 3, fill); // green
420     } else {
421       SetVolumeAttributes("SCMX", 1, 5, fill); // yellow
422     }
423   } else if(mn == "SCMX" || mn.Contains("SCX")) {
424     SetVolumeAttributes(mn.Data(), 1, 5, fill);// yellow
425     SetVolumeAttributes("PBTI", 1, 4, fill);
426   } else {
427     printf("<W> for volume |%s| did not defined volume attributes\n", mn.Data());
428   }
429
430   //  TString st("Shish-Kebab, 2x2 mm sampling, 62 layers, ");
431   TString st("Shashlyk, 2x2 mm sampling, 62 layers, ");
432   if    (sn.Contains("25"))   st = "Shish-Kebab, 5x5 mm sampling, 25 layers, ";
433   else if(sn.Contains("MAY05"))   st = "Shish-Kebab, 5x5 mm sampling, 77 layers, ";
434   if(sn.Contains("TRD1")) st += " TRD1, ";
435   if(sn.Contains("3X3"))  st += " 3x3, ";
436   st += name;
437
438   gROOT->ProcessLine("TGeant3 *g3 = (TGeant3*)gMC");
439   double cx=0.78, cy=2.;
440   if     (axis==1 && sn.Contains("TRD")==0) {
441    cx = cy = 1.;
442    gMC->Gsatt("PBTI", "seen", 0);
443   } else if (sn.Contains("TEST") && sn.Contains("3X3")) {
444     cx = cy = 0.7;
445     if (axis==3) cx = cy = 1.;
446   } else if (sn.Contains("TRD2")) {
447     if (axis==3) cx = cy = 2.;
448   } else if (mn.Contains("EMOD")) {
449     cx = cy = 0.5;
450   }
451   char cmd[200];
452
453   gMC->Gdopt("hide", optShad);
454   gMC->Gdopt("shad", optShad);
455
456   sprintf(cmd,"g3->Gdrawc(\"%s\", %i, %6.2f, 10., 10., %5.3f, %5.3f)\n",name, axis,dcut,cx,cy);
457   printf("\n%s\n",cmd); gROOT->ProcessLine(cmd);
458
459   sprintf(cmd,"gMC->Gdhead(1111, \"%s\")\n", st.Data());
460   printf("%s\n",cmd); gROOT->ProcessLine(cmd);
461 }
462   
463 void AliEMCALv2::DrawAlicWithHits(int mode)
464 { // 20-sep-04; does not work now
465   static TH2F *h2;
466   if(h2==0) h2 = new TH2F("h2","test fo hits", 60,0.5,60.5, 28,0.5,28.5);
467   else      h2->Reset();
468   if(mode==0) {
469     gROOT->ProcessLine("TGeant3 *g3 = (TGeant3*)gMC");
470     gMC->Gsatt("*","seen",0);
471     gMC->Gsatt("scm0","seen",5);
472
473     gROOT->ProcessLine("g3->Gdrawc(\"ALIC\", 1,   2.0, 12., -130, 0.3, 0.3)");
474     // g3->Gdrawc("ALIC", 1,   2.0, 10., -14, 0.05, 0.05)
475   }
476
477   TClonesArray *hits = Hits();
478   Int_t nhits = hits->GetEntries(), absId, nSupMod, nTower, nIphi, nIeta, iphi, ieta;
479   AliEMCALHit *hit = 0;
480   Double_t de, des=0.;
481   for(Int_t i=0; i<nhits; i++) {
482     hit   = (AliEMCALHit*)hits->UncheckedAt(i);
483     absId = hit->GetId();
484     de    = hit->GetEnergy();
485     des += de;
486     if(fGeometry->GetCellIndex(absId, nSupMod, nTower, nIphi, nIeta)){
487       //      printf(" de %f abs id %i smod %i tower %i | cell iphi %i : ieta %i\n",
488       // de, absId, nSupMod, nTower, nIphi, nIeta); 
489       if(nSupMod==3) {
490         fGeometry->GetCellPhiEtaIndexInSModule(nSupMod,nTower,nIphi,nIeta, iphi,ieta);
491         // printf(" iphi %i : ieta %i\n", iphi,ieta);
492         h2->Fill(double(ieta),double(iphi), de);
493       }
494     }
495     //    hit->Dump();
496   }
497   printf(" #hits %i : sum de %f -> %f \n", nhits, des, h2->Integral());
498   if(mode>0 && h2->Integral()>0.) h2->Draw();
499 }
500
501 void AliEMCALv2::SetVolumeAttributes(const char *name, int seen, int color, int fill)
502 {
503  /* seen=-2:volume is visible but none of its descendants;
504          -1:volume is not visible together with all its descendants;
505           0:volume is not visible;
506           1:volume is     visible. */
507   gMC->Gsatt(name, "seen", seen);
508   gMC->Gsatt(name, "colo", color);
509   gMC->Gsatt(name, "fill", fill); 
510   printf(" %s : seen %i color %i fill %i \n", name, seen, color, fill);
511
512
513 void AliEMCALv2::TestIndexTransition(int pri, int idmax)
514 { // for EMCAL_SHISH geometry
515   TString sn(fGeometry->GetName());
516   if(!sn.Contains("SHISH")) {
517     printf("Wrong geometry |%s| ! Bye \n", sn.Data());
518     return; 
519   }
520
521   Int_t nSupMod, nTower, nIphi, nIeta, idNew, nGood=0;
522   if(idmax==0) idmax = fGeometry->GetNCells();
523   for(Int_t id=1; id<=idmax; id++) {
524     if(!fGeometry->GetCellIndex(id, nSupMod, nTower, nIphi, nIeta)){
525       printf(" Wrong abs ID %i : #cells %i\n", id, fGeometry->GetNCells());
526       break;  
527     }
528     idNew = fGeometry->GetAbsCellId(nSupMod, nTower, nIphi, nIeta);
529     if(id != idNew || pri>0) {
530       printf(" ID %i : %i <- new id\n", id, idNew);
531       printf(" nSupMod   %i ",  nSupMod);
532       printf(" nTower    %i ", nTower);
533       printf(" nIphi     %i ", nIphi);
534       printf(" nIeta     %i \n", nIeta);
535      
536     } else nGood++;
537   }
538   printf(" Good decoding %i : %i <- #cells \n", nGood, fGeometry->GetNCells());
539 }
540
541 /* try to draw hits
542  gMC->Gsatt("*","seen",0);
543  gMC->Gsatt("scm0","seen",5);
544  g3->Gdrawc("alic", 1,   2.0, 12., -125, 0.3, 0.3);
545 */
546 //void AliEMCALv2::Gdaxis(float x0, float y0, float z0, float <axsiz) {gdaxis_(x0,y0,z0,axsiz);}
547
548 /*
549
550 void AliEMCALv2::RemapTrackHitIDs(Int_t *map) {
551   // remap track index numbers for primary and parent indices
552   // (Called by AliStack::PurifyKine)
553   if (Hits()==0)
554     return;
555   TIter hit_it(Hits());
556   Int_t i_hit=0;
557   while (AliEMCALHit *hit=dynamic_cast<AliEMCALHit*>(hit_it()) ) {
558     if (map[hit->GetIparent()]==-99)
559       cout << "Remapping, found -99 for parent id " << hit->GetIparent() << ", " << map[hit->GetIparent()] << ", i_hit " << i_hit << endl;
560     hit->SetIparent(map[hit->GetIparent()]);
561     if (map[hit->GetPrimary()]==-99)
562       cout << "Remapping, found -99 for primary id " << hit->GetPrimary() << ", " << map[hit->GetPrimary()] << ", i_hit " << i_hit << endl;
563     hit->SetPrimary(map[hit->GetPrimary()]);
564     i_hit++;
565   }
566 }
567
568 void AliEMCALv2::FinishPrimary() {
569   fCurPrimary=-1;
570   fCurParent=-1;
571   fCurTrack=-1;
572 }
573 */