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