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