]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOSvImpacts.cxx
Timing added to reconstruction
[u/mrichter/AliRoot.git] / PHOS / AliPHOSvImpacts.cxx
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 vImpacts of PHOS Manager class.
20 // This class inherits from v1 and adds impacts storing.
21 // Impacts stands for exact values of track coming to the detectors
22 // EMC, CPV or PPSD.
23 // Impacts are written to the same tree as hits are
24 // but in separate branches.
25 //---
26 //*-- Author: Yuri Kharlov (IHEP, Protvino/SUBATECH, Nantes)
27
28
29 // --- ROOT system ---
30
31 #include "TTree.h"
32
33 // --- Standard library ---
34
35 // --- AliRoot header files ---
36
37 #include "AliRun.h"
38 #include "AliMC.h"
39 #include "AliPHOSvImpacts.h"
40 #include "AliPHOSGeometry.h"
41 #include "AliPHOSImpact.h"
42
43 ClassImp(AliPHOSvImpacts)
44
45 //____________________________________________________________________________
46 AliPHOSvImpacts::AliPHOSvImpacts():AliPHOSv1()
47 {
48   // ctor
49 }
50
51 //____________________________________________________________________________
52 AliPHOSvImpacts::AliPHOSvImpacts(const char *name, const char *title):
53 AliPHOSv1(name,title) 
54 {
55   // ctor : title is used to identify the layout
56   //
57   // We store hits :
58   //   - fHits (the "normal" one), which retains the hits associated with
59   //     the current primary particle being tracked
60   //     (this array is reset after each primary has been tracked).
61   //     This part inherits from AliPHOSv1
62   //
63   // We store impacts :
64   //  - fEMCImpacts, fCPVImpacts which are
65   //    TList of EMC and CPV modules respectively, each
66   //    modules contains TClonesArray of AliPHOSImpacts
67   
68   fEMCImpacts  = new TList();
69   fCPVImpacts  = new TList();
70
71   Int_t nPHOSModules = GetGeometry()->GetNModules();
72   Int_t nCPVModules  = GetGeometry()->GetNModules();
73
74   Int_t iPHOSModule;
75   TClonesArray * impacts;
76   for (iPHOSModule=0; iPHOSModule<nPHOSModules; iPHOSModule++) {
77     fEMCImpacts->Add(new TClonesArray("AliPHOSImpact",200)) ;
78     fNEMCImpacts[iPHOSModule] = 0;
79     impacts = (TClonesArray *)fEMCImpacts->At(iPHOSModule);
80   }
81   for (iPHOSModule=0; iPHOSModule<nCPVModules; iPHOSModule++) {
82     fCPVImpacts->Add(new TClonesArray("AliPHOSImpact",200)) ;
83     fNCPVImpacts[iPHOSModule] = 0;
84     impacts = (TClonesArray *)fCPVImpacts->At(iPHOSModule);
85   }
86
87 }
88
89 //____________________________________________________________________________
90 AliPHOSvImpacts::~AliPHOSvImpacts()
91 {
92   // dtor
93
94   // Delete hits
95   if ( fHits ) {
96     fHits->Delete() ; 
97     delete fHits ;
98     fHits = 0 ; 
99   }
100
101   // Delete impacts in EMC, CPV
102   if ( fEMCImpacts ) {
103     fEMCImpacts->Delete() ; 
104     delete fEMCImpacts ;
105     fEMCImpacts = 0 ; 
106   }
107   if ( fCPVImpacts ) {
108     fCPVImpacts->Delete() ; 
109     delete fCPVImpacts ;
110     fCPVImpacts = 0 ; 
111   }
112 }
113
114 //____________________________________________________________________________
115 void AliPHOSvImpacts::AddImpact( char* det, Int_t shunt, Int_t primary, Int_t track, Int_t module,
116                            Int_t pid, TLorentzVector p, Float_t *xyz)
117 {
118   // Add an impact to the impact list.
119
120   TClonesArray * impacts = 0;
121   Int_t         nImpacts = 0;
122
123   if (strcmp(det,"EMC ")==0) {
124     impacts = (TClonesArray *)fEMCImpacts->At(module);
125     nImpacts= fNEMCImpacts[module];
126     fNEMCImpacts[module]++ ;
127   }
128   else if (strcmp(det,"CPV ")==0) {
129     impacts = (TClonesArray *)fCPVImpacts->At(module);
130     nImpacts= fNCPVImpacts[module];
131     fNCPVImpacts[module]++ ;
132   }
133
134   new((*impacts)[nImpacts]) AliPHOSImpact(shunt,primary,track,pid,p,xyz) ;
135
136   if (fDebug==1) {
137     printf("Module %d %s: ",module,det);
138     ((AliPHOSImpact*)(impacts->At(nImpacts)))->Print();
139   }
140 }
141
142 //____________________________________________________________________________
143 void AliPHOSvImpacts::MakeBranch(Option_t *opt, const char *file)
144 {  
145   // Create new branch in the current Hits Root Tree containing
146   // a list of PHOS impacts (exact values of track coming to detector)
147
148   AliDetector::MakeBranch(opt,file);
149   
150   Int_t bufferSize = 32000 ;
151   Int_t splitlevel = 0 ;
152   gAlice->TreeH()->Branch("PHOSEmcImpacts" , "TList", &fEMCImpacts , bufferSize, splitlevel);
153   gAlice->TreeH()->Branch("PHOSCpvImpacts" , "TList", &fCPVImpacts , bufferSize, splitlevel);
154   
155 }
156
157 //____________________________________________________________________________
158 void AliPHOSvImpacts::ResetHits() 
159 {              
160   // Reset impact branches for EMC, CPV and PPSD
161
162   AliDetector::ResetHits();
163
164   Int_t i;
165   for (i=0; i<GetGeometry()->GetNModules(); i++) {
166     ((TClonesArray*)fEMCImpacts->At(i)) -> Clear();
167     fNEMCImpacts[i] = 0 ;
168   }
169
170   for (i=0; i<GetGeometry()->GetNModules(); i++) {
171     ((TClonesArray*)fCPVImpacts->At(i)) -> Clear();
172     fNCPVImpacts[i] = 0 ;
173   }
174   
175 }
176
177 //_____________________________________________________________________________
178 void AliPHOSvImpacts::StepManager(void)
179 {
180   // Find impacts (tracks which enter the EMC, CPV)
181   // and add them to to the impact lists
182
183   AliPHOSv1::StepManager();
184
185   Float_t xyzm[3], xyzd[3], pm[3], pd[3];
186   TLorentzVector pmom     ;           // Lorentz momentum of the particle initiated hit
187   TLorentzVector pos      ;           // Lorentz vector of the track current position
188   Int_t          copy     ;
189
190   Int_t tracknumber =  gAlice->CurrentTrack() ; 
191   Int_t primary     =  gAlice->GetPrimary( gAlice->CurrentTrack() ); 
192   TString name      =  GetGeometry()->GetName() ; 
193
194   // Add impact to EMC
195
196   if( gMC->CurrentVolID(copy) == gMC->VolId("PXTL") &&
197       gMC->IsTrackEntering() ) {
198     gMC->TrackMomentum(pmom);
199     gMC->TrackPosition(pos) ;
200
201     Int_t i;
202     for (i=0; i<3; i++) xyzm[i] = pos[i];
203
204     for (i=0; i<3; i++) {
205       xyzm[i] = pos[i] ;
206       pm[i]   = pmom[i];
207     }
208     gMC -> Gmtod (xyzm, xyzd, 1);    // transform coordinate from master to daughter system
209     gMC -> Gmtod (pm,   pd,   2);    // transform 3-momentum from master to daughter system
210
211     // Select tracks coming to the crystal from up or down sides
212     if (pd[1]<0 && xyzd[1] >  GetGeometry()->GetCrystalSize(1)/2-0.001 ||
213         pd[1]>0 && xyzd[1] < -GetGeometry()->GetCrystalSize(1)/2+0.001) {
214       Int_t pid = gMC->TrackPid();
215       Int_t module;
216       gMC->CurrentVolOffID(10,module);
217       module--;
218       AddImpact("EMC ",fIshunt, primary,tracknumber, module, pid, pmom, xyzm);
219     }
220   }
221
222   // Add impact to CPV
223
224   if( gMC->CurrentVolID(copy) == gMC->VolId("PCPQ") &&
225       gMC->IsTrackEntering() ) {
226     gMC->TrackMomentum(pmom);
227     gMC->TrackPosition(pos) ;
228
229     Int_t i;
230     for (i=0; i<3; i++) xyzm[i] = pos[i];
231
232     for (i=0; i<3; i++) {
233       xyzm[i] = pos[i] ;
234       pm[i]   = pmom[i];
235     }
236     Int_t pid = gMC->TrackPid();
237     Int_t module;
238     gMC->CurrentVolOffID(3,module);
239     module--;
240     AddImpact("CPV ",fIshunt, primary,tracknumber, module, pid, pmom, xyzm);
241   }
242   
243 }