SDigits has become a copy of hits and threshold to associate primary particles has...
[u/mrichter/AliRoot.git] / EMCAL / AliEMCALv1.cxx
CommitLineData
b13bbe81 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
fceb09ed 26// This Class not stores information on all particles prior to EMCAL entry - in order to facilitate analysis.
27// This is done by setting fIShunt =2, and flagging all parents of particles entering the EMCAL.
28
1f4d29d2 29// 15/02/2002 .... Yves Schutz
30// 1. fSamplingFraction and fLayerToPreshowerRatio have been removed
31// 2. Timing signal is collected and added to hit
b13bbe81 32
33// --- ROOT system ---
34#include "TPGON.h"
35#include "TTUBS.h"
36#include "TNode.h"
37#include "TRandom.h"
38#include "TTree.h"
39#include "TGeometry.h"
fceb09ed 40#include "TParticle.h"
b13bbe81 41
42// --- Standard library ---
43
44#include <stdio.h>
45#include <string.h>
46#include <stdlib.h>
47#include <strstream.h>
48#include <iostream.h>
49#include <math.h>
50// --- AliRoot header files ---
51
52#include "AliEMCALv1.h"
53#include "AliEMCALHit.h"
54#include "AliEMCALGeometry.h"
55#include "AliConst.h"
56#include "AliRun.h"
57#include "AliMC.h"
58
59ClassImp(AliEMCALv1)
60
61
62//______________________________________________________________________
63AliEMCALv1::AliEMCALv1():AliEMCALv0(){
64 // ctor
65}
1f4d29d2 66
b13bbe81 67//______________________________________________________________________
68AliEMCALv1::AliEMCALv1(const char *name, const char *title):
69 AliEMCALv0(name,title){
70 // Standard Creator.
71
72 fHits= new TClonesArray("AliEMCALHit",1000);
73 gAlice->AddHitList(fHits);
74
75 fNhits = 0;
fceb09ed 76 fIshunt = 2; // All hits are associated with particles entering the calorimeter
b13bbe81 77}
1f4d29d2 78
b13bbe81 79//______________________________________________________________________
80AliEMCALv1::~AliEMCALv1(){
81 // dtor
82
83 if ( fHits) {
84 fHits->Delete();
85 delete fHits;
86 fHits = 0;
87 }
88}
89//______________________________________________________________________
61e0abb5 90void AliEMCALv1::AddHit(Int_t shunt, Int_t primary, Int_t tracknumber, Int_t iparent, Float_t ienergy,
91 Int_t id, Float_t * hits,Float_t * p){
b13bbe81 92 // Add a hit to the hit list.
ffa6d63b 93 // An EMCAL hit is the sum of all hits in a single segment
94 // originating from the same enterring particle
b13bbe81 95 Int_t hitCounter;
96
97 AliEMCALHit *newHit;
98 AliEMCALHit *curHit;
99 Bool_t deja = kFALSE;
100
61e0abb5 101 newHit = new AliEMCALHit(shunt, primary, tracknumber, iparent, ienergy, id, hits, p);
b13bbe81 102 for ( hitCounter = fNhits-1; hitCounter >= 0 && !deja; hitCounter-- ) {
103 curHit = (AliEMCALHit*) (*fHits)[hitCounter];
61e0abb5 104 // We add hits with the same tracknumber, while GEANT treats
b13bbe81 105 // primaries succesively
106 if(curHit->GetPrimary() != primary) break;
107 if( *curHit == *newHit ) {
108 *curHit = *curHit + *newHit;
109 deja = kTRUE;
110 } // end if
111 } // end for hitCounter
112
113 if ( !deja ) {
114 new((*fHits)[fNhits]) AliEMCALHit(*newHit);
115 fNhits++;
116 } // end if
117
118 delete newHit;
119}
120//______________________________________________________________________
121void AliEMCALv1::StepManager(void){
1f4d29d2 122 // Accumulates hits as long as the track stays in a single
123 // crystal or PPSD gas Cell
124
125 Int_t id[2]; // (layer, phi, Eta) indices
126 Int_t absid;
127 // position wrt MRS and energy deposited
128 Float_t xyzte[5]={0.,0.,0.,0.,0.};// position wrt MRS, time and energy deposited
129 Float_t pmom[4]={0.,0.,0.,0.};
130 TLorentzVector pos; // Lorentz vector of the track current position.
131 TLorentzVector mom; // Lorentz vector of the track current momentum.
132 Int_t tracknumber = gAlice->CurrentTrack();
133 Int_t primary = 0;
134 static Int_t iparent = 0;
135 static Float_t ienergy = 0;
136 Int_t copy = 0;
137
138
139 if(gMC->IsTrackEntering() && (strcmp(gMC->CurrentVolName(),"XALU") == 0)){ // This Particle in enterring the Calorimeter
140 gMC->TrackPosition(pos) ;
141 xyzte[0] = pos[0] ;
142 xyzte[1] = pos[1] ;
143 xyzte[2] = pos[2] ;
144 if ( (xyzte[0]*xyzte[0] + xyzte[1]*xyzte[1])
145 < (fGeom->GetEnvelop(0)+fGeom->GetGap2Active()+1.5 )*(fGeom->GetEnvelop(0)+fGeom->GetGap2Active()+1.5 ) ) {
146 iparent = tracknumber;
147 gMC->TrackMomentum(mom);
148 ienergy = mom[3];
149 TParticle * part = 0 ;
150 Int_t parent = iparent ;
151 while ( parent != -1 ) { // <------------- flags this particle to be kept and
152 //all the ancestors of this particle
153 part = gAlice->Particle(parent) ;
154 part->SetBit(kKeepBit);
155 parent = part->GetFirstMother() ;
fceb09ed 156 }
1f4d29d2 157 }
158 }
159 if(gMC->CurrentVolID(copy) == gMC->VolId("XPHI") ) { // We are in a Scintillator Layer
160
b13bbe81 161 gMC->CurrentVolOffID(1, id[0]); // get the POLY copy number;
162 gMC->CurrentVolID(id[1]); // get the phi number inside the layer
163 primary = gAlice->GetPrimary(tracknumber);
b13bbe81 164 gMC->TrackPosition(pos);
165 gMC->TrackMomentum(mom);
1f4d29d2 166 xyzte[0] = pos[0];
167 xyzte[1] = pos[1];
168 xyzte[2] = pos[2];
169 xyzte[3] = gMC->TrackTime() ;
170 xyzte[4] = gMC->Edep();
61e0abb5 171 pmom[0] = mom[0];
172 pmom[1] = mom[1];
173 pmom[2] = mom[2];
174 pmom[3] = mom[3];
175
1f4d29d2 176 if(xyzte[4] > 0.){// Track is inside a scintillator and deposits some energy
177 absid = (id[0]-1)*(fGeom->GetNPhi()) + id[1];
178 AddHit(fIshunt, primary,tracknumber, iparent, ienergy, absid, xyzte, pmom);
b13bbe81 179 } // there is deposited energy
1f4d29d2 180 }
b13bbe81 181}