]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PHOS/AliPHOSv0hits.cxx
Put back the two bellows in front of the absorber.
[u/mrichter/AliRoot.git] / PHOS / AliPHOSv0hits.cxx
CommitLineData
994235dd 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// Version of AliPHOSv0 which allows for keeping all hits in TreeH
19//
20//*-- Author: Gines MARTINEZ (SUBATECH)
21
22
23// --- ROOT system ---
24
25#include "TBRIK.h"
26#include "TNode.h"
27#include "TRandom.h"
28
29// --- Standard library ---
30
31#include <stdio.h>
32#include <string.h>
33#include <stdlib.h>
34#include <strstream.h>
35
36// --- AliRoot header files ---
37
84630675 38#include "AliPHOSv0hits.h"
994235dd 39#include "AliPHOSHit.h"
40#include "AliPHOSDigit.h"
41#include "AliPHOSReconstructioner.h"
42#include "AliRun.h"
43#include "AliConst.h"
44
45ClassImp(AliPHOSv0hits)
46
47//____________________________________________________________________________
48AliPHOSv0hits::AliPHOSv0hits()
49{
88714635 50 // default ctor
994235dd 51 fNTmpHits = 0 ;
52 fTmpHits = 0 ;
53}
54
55//____________________________________________________________________________
56AliPHOSv0hits::AliPHOSv0hits(const char *name, const char *title):
57 AliPHOSv0(name,title)
58{
88714635 59 // ctor
994235dd 60 fHits= new TClonesArray("AliPHOSHit",1000) ;
61}
62
63//____________________________________________________________________________
64AliPHOSv0hits::~AliPHOSv0hits()
65{
66 // dtor
67 fTmpHits->Delete() ;
68 delete fTmpHits ;
69 fTmpHits = 0 ;
70
ba81c761 71 fEmcRecPoints->Delete() ;
72 delete fEmcRecPoints ;
73 fEmcRecPoints = 0 ;
994235dd 74
ba81c761 75 fPpsdRecPoints->Delete() ;
76 delete fPpsdRecPoints ;
77 fPpsdRecPoints = 0 ;
994235dd 78
79 fTrackSegments->Delete() ;
80 delete fTrackSegments ;
81 fTrackSegments = 0 ;
ba81c761 82
994235dd 83}
84
85//____________________________________________________________________________
66380251 86void AliPHOSv0hits::AddHit(Int_t shunt, Int_t primary, Int_t tracknumber, Int_t Id, Float_t * hits)
994235dd 87{
88 // Add a hit to the hit list.
89 // In this version of AliPHOSv0, a PHOS hit is real geant
90 // hits in a single crystal or in a single PPSD gas cell
91
c6553021 92 // cout << "Primary is " << primary << endl;
93 //cout << "Tracknumber is " << tracknumber << endl;
94 //cout << "Vol Id is " << Id << endl;
95 //cout << "hits is " << hits[0] << " " << hits[1] << " " << hits[2] << " " << hits[3] <<endl;
994235dd 96
55ad9645 97 // cout << " Adding a hit number " << fNhits << endl ;
994235dd 98
99 TClonesArray &ltmphits = *fHits ;
100 AliPHOSHit *newHit ;
101
102 // fHits->Print("");
103
66380251 104 newHit = new AliPHOSHit(shunt, primary, tracknumber, Id, hits) ;
994235dd 105
106 // We DO want to save in TreeH the raw hits
107 // TClonesArray &lhits = *fHits;
55ad9645 108 // cout << " Adding a hit to fHits TCloneArray number " << fNhits << endl ;
994235dd 109 new(ltmphits[fNhits]) AliPHOSHit(*newHit) ;
110 fNhits++ ;
111
55ad9645 112 // cout << " Added a hit to fHits TCloneArray number " << fNhits << endl ;
994235dd 113 // We do not want to save in TreeH the raw hits
114 // new(lhits[fNhits]) AliPHOSHit(*newHit) ;
115 // fNhits++ ;
116
117 // Please note that the fTmpHits array must survive up to the
118 // end of the events, so it does not appear e.g. in ResetHits() (
119 // which is called at the end of each primary).
120
121 delete newHit;
122
123}
124
125
126
127
128//___________________________________________________________________________
129void AliPHOSv0hits::FinishEvent()
130{
131 // Makes the digits from the sum of summed hit in a single crystal or PPSD gas cell
132 // Adds to the energy the electronic noise
133 // Keeps digits with energy above fDigitThreshold
134
135 // Save the cumulated hits instead of raw hits (need to create the branch myself)
136 // It is put in the Digit Tree because the TreeH is filled after each primary
137 // and the TreeD at the end of the event.
138
139 Int_t i ;
140 Int_t relid[4];
141 Int_t j ;
142 TClonesArray &lDigits = *fDigits ;
143 AliPHOSHit * hit ;
144 AliPHOSDigit * newdigit ;
145 AliPHOSDigit * curdigit ;
146 Bool_t deja = kFALSE ;
147
148 for ( i = 0 ; i < fNhits ; i++ ) {
149 hit = (AliPHOSHit*)fHits->At(i) ;
150 newdigit = new AliPHOSDigit( hit->GetPrimary(), hit->GetId(), Digitize( hit->GetEnergy() ) ) ;
151 deja =kFALSE ;
152 for ( j = 0 ; j < fNdigits ; j++) {
153 curdigit = (AliPHOSDigit*) lDigits[j] ;
154 if ( *curdigit == *newdigit) {
155 *curdigit = *curdigit + *newdigit ;
156 deja = kTRUE ;
157 }
158 }
159 if ( !deja ) {
160 new(lDigits[fNdigits]) AliPHOSDigit(* newdigit) ;
161 fNdigits++ ;
162 }
163
164 delete newdigit ;
165 }
166
167 // Noise induced by the PIN diode of the PbWO crystals
168
169 Float_t energyandnoise ;
170 for ( i = 0 ; i < fNdigits ; i++ ) {
171 newdigit = (AliPHOSDigit * ) fDigits->At(i) ;
172 fGeom->AbsToRelNumbering(newdigit->GetId(), relid) ;
173
174 if (relid[1]==0){ // Digits belong to EMC (PbW0_4 crystals)
175 energyandnoise = newdigit->GetAmp() + Digitize(gRandom->Gaus(0., fPinElectronicNoise)) ;
176
177 if (energyandnoise < 0 )
178 energyandnoise = 0 ;
179
180 if ( newdigit->GetAmp() < fDigitThreshold ) // if threshold not surpassed, remove digit from list
181 fDigits->RemoveAt(i) ;
182 }
183 }
184
185 fDigits->Compress() ;
186
187 fNdigits = fDigits->GetEntries() ;
188 for (i = 0 ; i < fNdigits ; i++) {
189 newdigit = (AliPHOSDigit *) fDigits->At(i) ;
190 newdigit->SetIndexInList(i) ;
191 }
192
193}
194
677e29a3 195void AliPHOSv0hits::StepManager(void)
196{
197 // Accumulates hits as long as the track stays in a single crystal or PPSD gas Cell
198
199 Int_t relid[4] ; // (box, layer, row, column) indices
200 Float_t xyze[4] ; // position wrt MRS and energy deposited
201 TLorentzVector pos ;
202 Int_t copy ;
203
204 Int_t tracknumber = gAlice->CurrentTrack() ;
20427a27 205 Int_t primary = gAlice->GetPrimary( gAlice->CurrentTrack() );
677e29a3 206
207 TString name = fGeom->GetName() ;
208 if ( name == "GPS2" ) { // the CPV is a PPSD
209 if( gMC->CurrentVolID(copy) == gMC->VolId("GCEL") ) // We are inside a gas cell
210 {
211 gMC->TrackPosition(pos) ;
212 xyze[0] = pos[0] ;
213 xyze[1] = pos[1] ;
214 xyze[2] = pos[2] ;
215 xyze[3] = gMC->Edep() ;
216
217 if ( xyze[3] != 0 ) { // there is deposited energy
218 gMC->CurrentVolOffID(5, relid[0]) ; // get the PHOS Module number
219 gMC->CurrentVolOffID(3, relid[1]) ; // get the Micromegas Module number
220 // 1-> Geom->GetNumberOfModulesPhi() * fGeom->GetNumberOfModulesZ() upper
221 // > fGeom->GetNumberOfModulesPhi() * fGeom->GetNumberOfModulesZ() lower
222 gMC->CurrentVolOffID(1, relid[2]) ; // get the row number of the cell
223 gMC->CurrentVolID(relid[3]) ; // get the column number
224
225 // get the absolute Id number
226
227 Int_t absid ;
228 fGeom->RelToAbsNumbering(relid, absid) ;
229
230 // add current hit to the hit list
66380251 231 AddHit(fIshunt, primary, tracknumber, absid, xyze);
677e29a3 232
233 } // there is deposited energy
234 } // We are inside the gas of the CPV
235 } // GPS2 configuration
236
237 if(gMC->CurrentVolID(copy) == gMC->VolId("PXTL") ) // We are inside a PBWO crystal
238 {
239 gMC->TrackPosition(pos) ;
240 xyze[0] = pos[0] ;
241 xyze[1] = pos[1] ;
242 xyze[2] = pos[2] ;
243 xyze[3] = gMC->Edep() ;
244
245 if ( xyze[3] != 0 ) {
246 gMC->CurrentVolOffID(10, relid[0]) ; // get the PHOS module number ;
247 relid[1] = 0 ; // means PBW04
248 gMC->CurrentVolOffID(4, relid[2]) ; // get the row number inside the module
249 gMC->CurrentVolOffID(3, relid[3]) ; // get the cell number inside the module
250
251 // get the absolute Id number
252
253 Int_t absid ;
254 fGeom->RelToAbsNumbering(relid, absid) ;
255
256 // add current hit to the hit list
257
66380251 258 AddHit(fIshunt, primary,tracknumber, absid, xyze);
677e29a3 259
260 } // there is deposited energy
261 } // we are inside a PHOS Xtal
262}
263