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