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