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