]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOSvImpacts.cxx
8d888eda690638ac94e741c7936796a458b83344
[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()
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   //        GPS2 = 5 modules (EMC + PPSD)
57   //        IHEP = 5 modules (EMC + CPV )
58   //        MIXT = 4 modules (EMC + CPV ) and 1 module (EMC + PPSD)
59   //
60   // We store hits :
61   //   - fHits (the "normal" one), which retains the hits associated with
62   //     the current primary particle being tracked
63   //     (this array is reset after each primary has been tracked).
64   //     This part inherits from AliPHOSv1
65   //
66   // We store impacts :
67   //  - fEMCImpacts, fCPVImpacts, fPPSDImpacts which are
68   //    TList of EMC, CPV and PPSD modules respectively, each
69   //    modules contains TClonesArray of AliPHOSImpacts
70   
71   fEMCImpacts  = new TList();
72   fCPVImpacts  = new TList();
73   fPPSDImpacts = new TList();
74
75   Int_t nPHOSModules = GetGeometry()->GetNModules();
76   Int_t nCPVModules  = GetGeometry()->GetNCPVModules();
77   Int_t nPPSDModules = GetGeometry()->GetNPPSDModules();
78
79   Int_t iPHOSModule;
80   TClonesArray * impacts;
81   for (iPHOSModule=0; iPHOSModule<nPHOSModules; iPHOSModule++) {
82     fEMCImpacts->Add(new TClonesArray("AliPHOSImpact",200)) ;
83     fNEMCImpacts[iPHOSModule] = 0;
84     impacts = (TClonesArray *)fEMCImpacts->At(iPHOSModule);
85   }
86   for (iPHOSModule=0; iPHOSModule<nCPVModules; iPHOSModule++) {
87     fCPVImpacts->Add(new TClonesArray("AliPHOSImpact",200)) ;
88     fNCPVImpacts[iPHOSModule] = 0;
89     impacts = (TClonesArray *)fCPVImpacts->At(iPHOSModule);
90   }
91   for (iPHOSModule=0; iPHOSModule<nPPSDModules; iPHOSModule++) {
92     fPPSDImpacts->Add(new TClonesArray("AliPHOSImpact",200)) ;
93     fNPPSDImpacts[iPHOSModule] = 0;
94     impacts = (TClonesArray *)fPPSDImpacts->At(iPHOSModule);
95   }
96
97 }
98
99 //____________________________________________________________________________
100 AliPHOSvImpacts::~AliPHOSvImpacts()
101 {
102   // dtor
103
104   // Delete hits
105   if ( fHits ) {
106     fHits->Delete() ; 
107     delete fHits ;
108     fHits = 0 ; 
109   }
110
111   // Delete impacts in EMC, CPV and PPSD
112   if ( fEMCImpacts ) {
113     fEMCImpacts->Delete() ; 
114     delete fEMCImpacts ;
115     fEMCImpacts = 0 ; 
116   }
117   if ( fCPVImpacts ) {
118     fCPVImpacts->Delete() ; 
119     delete fCPVImpacts ;
120     fCPVImpacts = 0 ; 
121   }
122   if ( fPPSDImpacts ) {
123     fPPSDImpacts->Delete() ; 
124     delete fPPSDImpacts ;
125     fPPSDImpacts = 0 ; 
126   }
127 }
128
129 //____________________________________________________________________________
130 void AliPHOSvImpacts::AddImpact( char* det, Int_t shunt, Int_t primary, Int_t track, Int_t module,
131                            Int_t pid, TLorentzVector p, Float_t *xyz)
132 {
133   // Add an impact to the impact list.
134
135   TClonesArray * impacts = 0;
136   Int_t         nImpacts = 0;
137
138   if (strcmp(det,"EMC ")==0) {
139     impacts = (TClonesArray *)fEMCImpacts->At(module);
140     nImpacts= fNEMCImpacts[module];
141     fNEMCImpacts[module]++ ;
142   }
143   else if (strcmp(det,"CPV ")==0) {
144     impacts = (TClonesArray *)fCPVImpacts->At(module);
145     nImpacts= fNCPVImpacts[module];
146     fNCPVImpacts[module]++ ;
147   }
148   else if (strcmp(det,"PPSD")==0) {
149     impacts = (TClonesArray *)fPPSDImpacts->At(module);
150     nImpacts= fNPPSDImpacts[module];
151     fNPPSDImpacts[module]++ ;
152   }
153
154   new((*impacts)[nImpacts]) AliPHOSImpact(shunt,primary,track,pid,p,xyz) ;
155
156   if (fDebug==1) {
157     printf("Module %d %s: ",module,det);
158     ((AliPHOSImpact*)(impacts->At(nImpacts)))->Print();
159   }
160 }
161
162 //____________________________________________________________________________
163 void AliPHOSvImpacts::MakeBranch(Option_t *opt, const char *file)
164 {  
165   // Create new branch in the current Hits Root Tree containing
166   // a list of PHOS impacts (exact values of track coming to detector)
167
168   AliDetector::MakeBranch(opt,file);
169   
170   Int_t bufferSize = 32000 ;
171   Int_t splitlevel = 0 ;
172   gAlice->TreeH()->Branch("PHOSEmcImpacts" , "TList", &fEMCImpacts , bufferSize, splitlevel);
173   gAlice->TreeH()->Branch("PHOSCpvImpacts" , "TList", &fCPVImpacts , bufferSize, splitlevel);
174   gAlice->TreeH()->Branch("PHOSPpsdImpacts", "TList", &fPPSDImpacts, bufferSize, splitlevel);
175   
176 }
177
178 //____________________________________________________________________________
179 void AliPHOSvImpacts::ResetHits() 
180 {              
181   // Reset impact branches for EMC, CPV and PPSD
182
183   AliDetector::ResetHits();
184
185   Int_t i;
186   for (i=0; i<GetGeometry()->GetNModules(); i++) {
187     ((TClonesArray*)fEMCImpacts->At(i)) -> Clear();
188     fNEMCImpacts[i] = 0 ;
189   }
190
191   if ( strcmp(GetGeometry()->GetName(),"IHEP") == 0 || strcmp(GetGeometry()->GetName(),"MIXT") == 0 ) {
192     for (i=0; i<GetGeometry()->GetNCPVModules(); i++) {
193       ((TClonesArray*)fCPVImpacts->At(i)) -> Clear();
194       fNCPVImpacts[i] = 0 ;
195     }
196   }
197
198   if ( strcmp(GetGeometry()->GetName(),"GPS2") == 0 || strcmp(GetGeometry()->GetName(),"MIXT") == 0 ) {
199     for (i=0; i<GetGeometry()->GetNPPSDModules(); i++) {
200       ((TClonesArray*)fPPSDImpacts->At(i)) -> Clear();
201       fNPPSDImpacts[i] = 0 ;
202     }
203   }
204   
205 }
206
207 //_____________________________________________________________________________
208 void AliPHOSvImpacts::StepManager(void)
209 {
210   // Find impacts (tracks which enter the EMC, CPV or PPSD)
211   // and add them to to the impact lists
212
213   AliPHOSv1::StepManager();
214
215   Float_t xyzm[3], xyzd[3], pm[3], pd[3];
216   TLorentzVector pmom     ;           // Lorentz momentum of the particle initiated hit
217   TLorentzVector pos      ;           // Lorentz vector of the track current position
218   Int_t          copy     ;
219
220   Int_t tracknumber =  gAlice->CurrentTrack() ; 
221   Int_t primary     =  gAlice->GetPrimary( gAlice->CurrentTrack() ); 
222   TString name      =  GetGeometry()->GetName() ; 
223
224   // Add impact to EMC
225
226   if( gMC->CurrentVolID(copy) == gMC->VolId("PXTL") &&
227       gMC->IsTrackEntering() ) {
228     gMC->TrackMomentum(pmom);
229     gMC->TrackPosition(pos) ;
230
231     Int_t i;
232     for (i=0; i<3; i++) xyzm[i] = pos[i];
233
234     for (i=0; i<3; i++) {
235       xyzm[i] = pos[i] ;
236       pm[i]   = pmom[i];
237     }
238     gMC -> Gmtod (xyzm, xyzd, 1);    // transform coordinate from master to daughter system
239     gMC -> Gmtod (pm,   pd,   2);    // transform 3-momentum from master to daughter system
240
241     // Select tracks coming to the crystal from up or down sides
242     if (pd[1]<0 && xyzd[1] >  GetGeometry()->GetCrystalSize(1)/2-0.001 ||
243         pd[1]>0 && xyzd[1] < -GetGeometry()->GetCrystalSize(1)/2+0.001) {
244       Int_t pid = gMC->TrackPid();
245       Int_t module;
246       gMC->CurrentVolOffID(10,module);
247       if ( name == "MIXT" && strcmp(gMC->CurrentVolOffName(10),"PHO1") == 0 )
248         module += GetGeometry()->GetNModules() - GetGeometry()->GetNPPSDModules();
249       module--;
250       AddImpact("EMC ",fIshunt, primary,tracknumber, module, pid, pmom, xyzm);
251     }
252   }
253
254   // Add impact to CPV
255
256   if( (name == "IHEP" || name == "MIXT") &&
257       gMC->CurrentVolID(copy) == gMC->VolId("PCPQ") &&
258       gMC->IsTrackEntering() ) {
259     gMC->TrackMomentum(pmom);
260     gMC->TrackPosition(pos) ;
261
262     Int_t i;
263     for (i=0; i<3; i++) xyzm[i] = pos[i];
264
265     for (i=0; i<3; i++) {
266       xyzm[i] = pos[i] ;
267       pm[i]   = pmom[i];
268     }
269     Int_t pid = gMC->TrackPid();
270     Int_t module;
271     gMC->CurrentVolOffID(3,module);
272     module--;
273     AddImpact("CPV ",fIshunt, primary,tracknumber, module, pid, pmom, xyzm);
274   }
275
276   // Add impact to PPSD
277
278   if( (name == "GPS2" || name == "MIXT") &&
279       gMC->CurrentVolID(copy) == gMC->VolId("PPCE") &&
280       gMC->IsTrackEntering() ) {
281     gMC->TrackMomentum(pmom);
282     gMC->TrackPosition(pos) ;
283
284     Int_t i;
285     for (i=0; i<3; i++) xyzm[i] = pos[i];
286
287     for (i=0; i<3; i++) {
288       xyzm[i] = pos[i] ;
289       pm[i]   = pmom[i];
290     }
291     gMC -> Gmtod (xyzm, xyzd, 1);    // transform coordinate from master to daughter system
292     gMC -> Gmtod (pm,   pd,   2);    // transform 3-momentum from master to daughter system
293
294     // Select tracks coming to the crystal from up or down sides
295     if (pd[1]<0 && xyzd[1] >  (GetGeometry()->GetConversionGap() +  GetGeometry()->GetAvalancheGap())/2-0.001 ||
296         pd[1]>0 && xyzd[1] < -(GetGeometry()->GetConversionGap() +  GetGeometry()->GetAvalancheGap())/2+0.001) {
297       Int_t pid = gMC->TrackPid();
298       Int_t module;
299       gMC->CurrentVolOffID(5,module);
300       module--;
301       AddImpact("PPSD",fIshunt, primary,tracknumber, module, pid, pmom, xyzm);
302     }
303   }
304 }