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