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