]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/AliEMCALv3.cxx
Include cassert (Christian)
[u/mrichter/AliRoot.git] / EMCAL / AliEMCALv3.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 "AliEMCALv3.h"
30 #include "AliEMCALHitv1.h"
31
32 #include "TParticle.h"
33 #include "TVirtualMC.h"
34 #include "TBrowser.h"
35 #include "TH1.h"
36 #include "TH2.h"
37 #include <cassert>
38
39 // --- Standard library ---
40
41 // --- AliRoot header files ---
42 #include "AliEMCALGeometry.h"
43 #include "AliRun.h"
44 #include "AliHeader.h"
45 #include "AliMC.h"
46 #include "AliPoints.h"
47 // for TRD1,2 case
48
49 ClassImp(AliEMCALv3)
50
51   //extern "C" void gdaxis_(float &x0, float &y0, float &z0, float &axsiz);
52
53
54 //______________________________________________________________________
55 AliEMCALv3::AliEMCALv3():AliEMCALv1(){
56   // ctor
57
58 }
59
60 //______________________________________________________________________
61 AliEMCALv3::AliEMCALv3(const char *name, const char *title): AliEMCALv1(name,title) {
62     // Standard Creator.
63
64     fHits= new TClonesArray("AliEMCALHitv1",2000);
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     fHDe = fHNhits = 0;
73     //    if (gDebug>0){
74     if (1){
75       TH1::AddDirectory(0);
76       fHDe = new TH1F("fHDe","De in EMCAL", 1000, 0., 1.);
77       fHNhits = new TH1F("fHNhits","#hits in EMCAL", 1001, -0.5, 1000.5);
78       int nbin = 77*2;
79       fHDeDz = new TH1F("fHDeDz","Longitudinal shower profile", 6*nbin, 0.0, 0.16*nbin);
80       
81       fHistograms->Add(fHDe);
82       fHistograms->Add(fHNhits);
83       fHistograms->Add(fHDeDz);
84       TH1::AddDirectory(1);
85     }
86 }
87
88 //______________________________________________________________________
89 AliEMCALv3::AliEMCALv3(const AliEMCALv3 & emcal):AliEMCALv1(emcal)
90 {
91   fGeometry = emcal.fGeometry;
92   fHDe = emcal.fHDe;
93   fHNhits = emcal.fHNhits;
94   fHDeDz = emcal.fHDeDz;
95
96 }
97
98 //______________________________________________________________________
99 AliEMCALv3::~AliEMCALv3(){
100     // dtor
101
102     if ( fHits) {
103         fHits->Delete();
104         delete fHits;
105         fHits = 0;
106     }
107 }
108
109 //______________________________________________________________________
110 void AliEMCALv3::AddHit(Int_t shunt, Int_t primary, Int_t tracknumber, Int_t iparent, Float_t ienergy, 
111                         Int_t id, Float_t * hits,Float_t * p){
112     // Add a hit to the hit list.
113     // An EMCAL hit is the sum of all hits in a tower section
114     //   originating from the same entering particle 
115     static Int_t hitCounter;
116     static AliEMCALHitv1 *newHit, *curHit;
117     static Bool_t deja;
118
119     deja = kFALSE;
120
121     newHit = new AliEMCALHitv1(shunt, primary, tracknumber, iparent, ienergy, id, hits, p);
122     for ( hitCounter = fNhits-1; hitCounter >= 0 && !deja; hitCounter-- ) {
123         curHit = (AliEMCALHitv1*) (*fHits)[hitCounter];
124         // We add hits with the same tracknumber, while GEANT treats
125         // primaries succesively
126         if(curHit->GetPrimary() != primary) 
127           break;
128         if( *curHit == *newHit ) {
129             *curHit = *curHit + *newHit;
130             deja = kTRUE;
131             //            break; // 30-aug-04 by PAI 
132         } // end if
133     } // end for hitCounter
134     
135     if ( !deja ) {
136         new((*fHits)[fNhits]) AliEMCALHitv1(*newHit);
137         fNhits++;
138     }
139     //    printf(" fNhits %i \n", fNhits); 
140     delete newHit;
141 }
142 //______________________________________________________________________
143 void AliEMCALv3::StepManager(void){
144   // Accumulates hits as long as the track stays in a tower
145
146   // position wrt MRS and energy deposited
147   static Float_t        xyzte[5]={0.,0.,0.,0.,0.};// position wrt MRS, time and energy deposited
148   static Float_t        pmom[4]={0.,0.,0.,0.};
149   static TLorentzVector pos; // Lorentz vector of the track current position.
150   static TLorentzVector mom; // Lorentz vector of the track current momentum.
151   static Float_t ienergy = 0;
152   static TString curVolName;
153   static int supModuleNumber, moduleNumber, yNumber, xNumber, absid;
154   static int keyGeom=0;
155   static char *vn = "SX"; // 15-mar-05
156   static int nSMOP[7]={1,3,5,7,9,11}; // 30-mar-05
157   static int nSMON[7]={2,4,6,8,10,12};
158   static Float_t depositedEnergy=0.0; 
159
160   if(keyGeom == 0) {
161     keyGeom = 2;
162     if(gMC->VolId("PBMO")==0 || gMC->VolId("WSUC")==1) {
163       vn      = "SCMX";   // old TRD2(TRD1) or WSUC
164       keyGeom = 1;
165     }    
166     printf("AliEMCALv3::StepManager():  keyGeom %i : Sensetive volume %s \n", 
167     keyGeom, vn); 
168     if(gMC->VolId("WSUC")==1) printf(" WSUC - cosmic ray stand geometry \n");
169   }
170   Int_t tracknumber =  gAlice->GetMCApp()->GetCurrentTrackNumber();
171
172   curVolName = gMC->CurrentVolName();
173   if(curVolName.Contains(vn)) { // We are in a scintillator layer
174     //    printf(" keyGeom %i : Sensetive volume %s (%s) \n", keyGeom, curVolName.Data(), vn); 
175     
176     if( ((depositedEnergy = gMC->Edep()) > 0.)  && (gMC->TrackTime() < fTimeCut)){// Track is inside a scintillator and deposits some energy
177       //       Info("StepManager "," entry %i DE %f",++ientry, depositedEnergy); // for testing
178        if (fCurPrimary==-1) 
179         fCurPrimary=gAlice->GetMCApp()->GetPrimary(tracknumber);
180
181       if (fCurParent==-1 || tracknumber != fCurTrack) {
182         // Check parentage
183         Int_t parent=tracknumber;
184         if (fCurParent != -1) {
185           while (parent != fCurParent && parent != -1) {
186             TParticle *part=gAlice->GetMCApp()->Particle(parent);
187             parent=part->GetFirstMother();
188           }
189         }
190         if (fCurParent==-1 || parent==-1) {
191           Int_t parent=tracknumber;
192           TParticle *part=gAlice->GetMCApp()->Particle(parent);
193           while (parent != -1 && fGeometry->IsInEMCAL(part->Vx(),part->Vy(),part->Vz())) {
194             parent=part->GetFirstMother();
195             if (parent!=-1) 
196               part=gAlice->GetMCApp()->Particle(parent);
197           } 
198           fCurParent=parent;
199           if (fCurParent==-1)
200             Error("StepManager","Cannot find parent");
201           else {
202             TParticle *part=gAlice->GetMCApp()->Particle(fCurParent);
203             ienergy = part->Energy(); 
204           }
205           while (parent != -1) {
206             part=gAlice->GetMCApp()->Particle(parent);
207             part->SetBit(kKeepBit);
208             parent=part->GetFirstMother();
209           }
210         }
211         fCurTrack=tracknumber;
212       }    
213       gMC->TrackPosition(pos);
214       xyzte[0] = pos[0];
215       xyzte[1] = pos[1];
216       xyzte[2] = pos[2];
217       xyzte[3] = gMC->TrackTime() ;       
218       
219       gMC->TrackMomentum(mom);
220       pmom[0] = mom[0];
221       pmom[1] = mom[1];
222       pmom[2] = mom[2];
223       pmom[3] = mom[3];
224       
225       //      if(ientry%200 > 0) return; // testing
226       supModuleNumber = moduleNumber = yNumber = xNumber = absid = 0;
227       if(keyGeom >= 1) { // old style
228         gMC->CurrentVolOffID(4, supModuleNumber);
229         gMC->CurrentVolOffID(3, moduleNumber);
230         gMC->CurrentVolOffID(1, yNumber);
231         gMC->CurrentVolOffID(0, xNumber); // really x number now
232         if(strcmp(gMC->CurrentVolOffName(4),"SM10")==0) supModuleNumber += 10; // 13-oct-05
233       } else {
234         gMC->CurrentVolOffID(5, supModuleNumber);
235         gMC->CurrentVolOffID(4, moduleNumber);
236         gMC->CurrentVolOffID(1, yNumber);
237         gMC->CurrentVolOffID(0, xNumber);
238         if     (strcmp(gMC->CurrentVolOffName(5),"SMOP")==0) supModuleNumber = nSMOP[supModuleNumber-1];
239         else if(strcmp(gMC->CurrentVolOffName(5),"SMON")==0) supModuleNumber = nSMON[supModuleNumber-1];
240         else   assert(0); // something wrong
241       }
242       absid = fGeometry->GetAbsCellId(supModuleNumber, moduleNumber, yNumber, xNumber);
243     
244       if (absid == -1) Fatal("StepManager()", "Wrong id ") ;
245
246       Float_t lightYield =  depositedEnergy ;
247       // Apply Birk's law (copied from G3BIRK)
248
249       if (gMC->TrackCharge()!=0) { // Check
250           Float_t BirkC1_mod = 0;
251         if (fBirkC0==1){ // Apply correction for higher charge states
252           if (TMath::Abs(gMC->TrackCharge())>=2) BirkC1_mod=fBirkC1*7.2/12.6;
253           else                                    BirkC1_mod=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.+BirkC1_mod*dedxcm+fBirkC2*dedxcm*dedxcm);
260       } 
261
262       xyzte[4] = lightYield; // 15-dec-04
263       Float_t xd[5]; // see Gmtod
264       gMC->Gmtod(xyzte, xd, 1);
265       xd[3] = xyzte[3];
266       xd[4] = xyzte[4];
267       //      printf(" x %7.2f y %7.2f z %7.2f de %f\n", xd[0], xd[1], xd[2], xd[4]);
268         
269       if (gDebug == -2) 
270       printf("#sm %2i #m %3i #x %1i #z %1i -> absid %i : xyzte[4] = %f\n",
271       supModuleNumber,moduleNumber,yNumber,xNumber,absid, xd[4]);
272
273       AddHit(fIshunt, fCurPrimary,tracknumber, fCurParent, ienergy, absid,  xd, pmom);
274     } // there is deposited energy
275   }
276 }
277
278 void AliEMCALv3::FinishEvent()
279 { // 26-may-05
280   static double de=0.;
281   fHNhits->Fill(double(fHits->GetEntries()));
282   de = GetDepositEnergy(0);
283   if(fHDe) fHDe->Fill(de);
284 }
285
286 Double_t AliEMCALv3::GetDepositEnergy(int print)
287 { // 23-mar-05 - for testing
288   //  cout<<"AliEMCALv3::GetDepositEnergy() : fHits "<<fHits<<endl; 
289   if(fHits == 0) return 0.;
290   static AliEMCALHitv1  *hit=0;
291   static Double_t de=0., deHit=0., zhl=0.0, zShift = 0.16*77;
292   for(int ih=0; ih<fHits->GetEntries(); ih++) {
293     hit = (AliEMCALHitv1*)fHits->UncheckedAt(ih);
294     deHit  = (double)hit->GetEnergy();
295     zhl    = (double)hit->Z() + zShift;
296     if(zhl>2.*zShift) zhl = 2*zShift - 0.002;
297     de    += deHit;
298     if(fHDeDz) fHDeDz->Fill(zhl, deHit);
299   }
300   if(print>0) {
301     printf(" #hits %i de %f \n", fHits->GetEntries(), de);
302     if(print>1) {
303       printf(" #primary particles %i\n", gAlice->GetHeader()->GetNprimary()); 
304     }
305   }
306   return de;
307 }
308
309 void AliEMCALv3::Browse(TBrowser* b)
310 {
311   TObject::Browse(b);
312 }