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