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