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