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