]> git.uio.no Git - u/mrichter/AliRoot.git/blame_incremental - EMCAL/AliEMCALv1.cxx
DeltaPhi function corrected; Additional checks in GetKStar
[u/mrichter/AliRoot.git] / EMCAL / AliEMCALv1.cxx
... / ...
CommitLineData
1/**************************************************************************
2 * Copyright(c) 1998-1999, 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 v1 of EMCAL Manager class
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: Sahal Yacoob (LBL /UCT)
24//*-- : Jennifer Klay (LBL)
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// 15/02/2002 .... Yves Schutz
29// 1. fSamplingFraction and fLayerToPreshowerRatio have been removed
30// 2. Timing signal is collected and added to hit
31
32// --- ROOT system ---
33#include "TParticle.h"
34#include "TVirtualMC.h"
35
36// --- Standard library ---
37
38// --- AliRoot header files ---
39#include "AliEMCALv1.h"
40#include "AliEMCALHit.h"
41#include "AliEMCALGeometry.h"
42#include "AliRun.h"
43#include "AliMC.h"
44
45ClassImp(AliEMCALv1)
46
47
48//______________________________________________________________________
49AliEMCALv1::AliEMCALv1():AliEMCALv0(){
50 // ctor
51
52}
53
54//______________________________________________________________________
55AliEMCALv1::AliEMCALv1(const char *name, const char *title):
56 AliEMCALv0(name,title){
57 // Standard Creator.
58
59 fHits= new TClonesArray("AliEMCALHit",1000);
60 gAlice->GetMCApp()->AddHitList(fHits);
61
62 fNhits = 0;
63 fIshunt = 2; // All hits are associated with particles entering the calorimeter
64}
65
66//______________________________________________________________________
67AliEMCALv1::~AliEMCALv1(){
68 // dtor
69
70 if ( fHits) {
71 fHits->Delete();
72 delete fHits;
73 fHits = 0;
74 }
75}
76//______________________________________________________________________
77void AliEMCALv1::AddHit(Int_t shunt, Int_t primary, Int_t tracknumber, Int_t iparent, Float_t ienergy,
78 Int_t id, Float_t * hits,Float_t * p){
79 // Add a hit to the hit list.
80 // An EMCAL hit is the sum of all hits in a tower section
81 // originating from the same entering particle
82 Int_t hitCounter;
83
84 AliEMCALHit *newHit;
85 AliEMCALHit *curHit;
86 Bool_t deja = kFALSE;
87
88 newHit = new AliEMCALHit(shunt, primary, tracknumber, iparent, ienergy, id, hits, p);
89 for ( hitCounter = fNhits-1; hitCounter >= 0 && !deja; hitCounter-- ) {
90 curHit = (AliEMCALHit*) (*fHits)[hitCounter];
91 // We add hits with the same tracknumber, while GEANT treats
92 // primaries succesively
93 if(curHit->GetPrimary() != primary)
94 break;
95 if( *curHit == *newHit ) {
96 *curHit = *curHit + *newHit;
97 deja = kTRUE;
98 } // end if
99 } // end for hitCounter
100
101 if ( !deja ) {
102 new((*fHits)[fNhits]) AliEMCALHit(*newHit);
103 fNhits++;
104 } // end if
105
106 delete newHit;
107}
108//______________________________________________________________________
109void AliEMCALv1::StepManager(void){
110 // Accumulates hits as long as the track stays in a tower
111
112 Int_t id[2]; // (phi, Eta) indices
113 // position wrt MRS and energy deposited
114 Float_t xyzte[5]={0.,0.,0.,0.,0.};// position wrt MRS, time and energy deposited
115 Float_t pmom[4]={0.,0.,0.,0.};
116 TLorentzVector pos; // Lorentz vector of the track current position.
117 TLorentzVector mom; // Lorentz vector of the track current momentum.
118 Int_t tracknumber = gAlice->GetMCApp()->GetCurrentTrackNumber();
119 Int_t primary = 0;
120 static Int_t iparent = 0;
121 static Float_t ienergy = 0;
122 Int_t copy = 0;
123
124 AliEMCALGeometry * geom = GetGeometry() ;
125
126 if(gMC->IsTrackEntering() && (strcmp(gMC->CurrentVolName(),"XALU") == 0)){ // This Particle in enterring the Calorimeter
127 gMC->TrackPosition(pos) ;
128 xyzte[0] = pos[0] ;
129 xyzte[1] = pos[1] ;
130 xyzte[2] = pos[2] ;
131 if ( (xyzte[0]*xyzte[0] + xyzte[1]*xyzte[1])
132 < (geom->GetEnvelop(0)+geom->GetGap2Active()+1.5 )*(geom->GetEnvelop(0)+geom->GetGap2Active()+1.5 ) ) {
133 iparent = tracknumber;
134 gMC->TrackMomentum(mom);
135 ienergy = mom[3];
136 TParticle * part = 0 ;
137 Int_t parent = iparent ;
138 while ( parent != -1 ) { // <------------- flags this particle to be kept and
139 //all the ancestors of this particle
140 part = gAlice->GetMCApp()->Particle(parent) ;
141 part->SetBit(kKeepBit);
142 parent = part->GetFirstMother() ;
143 }
144 }
145 }
146 if(gMC->CurrentVolID(copy) == gMC->VolId("XPHI") ) { // We are in a Scintillator Layer
147
148 Float_t depositedEnergy ;
149
150 if( (depositedEnergy = gMC->Edep()) > 0.){// Track is inside a scintillator and deposits some energy
151
152 // use sampling fraction to get original energy --HG
153 depositedEnergy = depositedEnergy * geom->GetSampling();
154
155 gMC->TrackPosition(pos);
156 xyzte[0] = pos[0];
157 xyzte[1] = pos[1];
158 xyzte[2] = pos[2];
159 xyzte[3] = gMC->TrackTime() ;
160
161 gMC->TrackMomentum(mom);
162 pmom[0] = mom[0];
163 pmom[1] = mom[1];
164 pmom[2] = mom[2];
165 pmom[3] = mom[3];
166
167 gMC->CurrentVolOffID(1, id[0]); // get the POLY copy number;
168 gMC->CurrentVolID(id[1]); // get the phi number inside the layer
169
170 Int_t tower = (id[0]-1) % geom->GetNZ() + 1 + (id[1] - 1) * geom->GetNZ() ;
171 Int_t layer = static_cast<Int_t>((id[0]-1)/(geom->GetNZ())) + 1 ;
172 Int_t absid = tower;
173 Int_t nlayers = geom->GetNECLayers();
174 if ((layer > nlayers)||(layer<1))
175 Fatal("StepManager", "Wrong calculation of layer number: layer = %d > %d\n", layer, nlayers) ;
176
177
178 Float_t lightYield = depositedEnergy ;
179 xyzte[4] = lightYield ;
180
181 primary = gAlice->GetMCApp()->GetPrimary(tracknumber);
182
183 if (gDebug == 2)
184 printf("StepManager: id0 = %d, id1 = %d, absid = %d tower = %d layer = %d energy = %f\n", id[0], id[1], absid, tower, layer, xyzte[4]) ;
185
186 AddHit(fIshunt, primary,tracknumber, iparent, ienergy, absid, xyzte, pmom);
187 } // there is deposited energy
188 }
189}