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