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