]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOSvImpacts.cxx
Obsolete - removed
[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::AddImpact( char* det, Int_t shunt, Int_t primary, Int_t track, Int_t module,
121                            Int_t pid, TLorentzVector p, Float_t *xyz)
122 {
123   // Add an impact to the impact list.
124
125   TClonesArray * impacts = 0;
126   Int_t         nImpacts = 0;
127
128   if (strcmp(det,"EMC ")==0) {
129     impacts = dynamic_cast<TClonesArray *>(fEMCImpacts->At(module));
130     nImpacts= fNEMCImpacts[module];
131     fNEMCImpacts[module]++ ;
132   }
133   else if (strcmp(det,"CPV ")==0) {
134     impacts = dynamic_cast<TClonesArray *>(fCPVImpacts->At(module));
135     nImpacts= fNCPVImpacts[module];
136     fNCPVImpacts[module]++ ;
137   }
138
139   new((*impacts)[nImpacts]) AliPHOSImpact(shunt,primary,track,pid,p,xyz) ;
140
141   if (fDebug==1) {
142     printf("Module %d %s: ",module,det);
143     (dynamic_cast<AliPHOSImpact*>((impacts->At(nImpacts))))->Print();
144   }
145 }
146
147 //____________________________________________________________________________
148 void AliPHOSvImpacts::MakeBranch(Option_t *opt)
149 {  
150   // Create new branch in the current Hits Root Tree containing
151   // a list of PHOS impacts (exact values of track coming to detector)
152
153   AliDetector::MakeBranch(opt);
154   
155   Int_t bufferSize = 32000 ;
156   Int_t splitlevel = 0 ;
157   TreeH()->Branch("PHOSEmcImpacts" , "TList", &fEMCImpacts , bufferSize, splitlevel);
158   TreeH()->Branch("PHOSCpvImpacts" , "TList", &fCPVImpacts , bufferSize, splitlevel);
159   
160 }
161
162 //____________________________________________________________________________
163 void AliPHOSvImpacts::ResetHits() 
164 {              
165   // Reset impact branches for EMC, CPV and PPSD
166
167   AliDetector::ResetHits();
168
169   Int_t i;
170   for (i=0; i<GetGeometry()->GetNModules(); i++) {
171     (dynamic_cast<TClonesArray*>(fEMCImpacts->At(i))) -> Clear();
172     fNEMCImpacts[i] = 0 ;
173   }
174
175   for (i=0; i<GetGeometry()->GetNModules(); i++) {
176     (dynamic_cast<TClonesArray*>(fCPVImpacts->At(i))) -> Clear();
177     fNCPVImpacts[i] = 0 ;
178   }
179   
180 }
181
182 //_____________________________________________________________________________
183 void AliPHOSvImpacts::StepManager(void)
184 {
185   // Find impacts (tracks which enter the EMC, CPV)
186   // and add them to to the impact lists
187
188   AliPHOSv1::StepManager();
189
190   Float_t xyzm[3], xyzd[3], pm[3], pd[3];
191   TLorentzVector pmom     ;           // Lorentz momentum of the particle initiated hit
192   TLorentzVector pos      ;           // Lorentz vector of the track current position
193   Int_t          copy     ;
194
195   Int_t tracknumber =  gAlice->GetMCApp()->GetCurrentTrackNumber() ; 
196   Int_t primary     =  gAlice->GetMCApp()->GetPrimary( gAlice->GetMCApp()->GetCurrentTrackNumber() ); 
197   TString name      =  GetGeometry()->GetName() ; 
198
199   // Add impact to EMC
200
201   if( gMC->CurrentVolID(copy) == gMC->VolId("PXTL") &&
202       gMC->IsTrackEntering() ) {
203     gMC->TrackMomentum(pmom);
204     gMC->TrackPosition(pos) ;
205
206     Int_t i;
207     for (i=0; i<3; i++) xyzm[i] = pos[i];
208
209     for (i=0; i<3; i++) {
210       xyzm[i] = pos[i] ;
211       pm[i]   = pmom[i];
212     }
213     gMC -> Gmtod (xyzm, xyzd, 1);    // transform coordinate from master to daughter system
214     gMC -> Gmtod (pm,   pd,   2);    // transform 3-momentum from master to daughter system
215
216     // Select tracks coming to the crystal from up or down sides
217     if (pd[1]<0 && xyzd[1] >  GetGeometry()->GetCrystalSize(1)/2-0.1 ||
218         pd[1]>0 && xyzd[1] < -GetGeometry()->GetCrystalSize(1)/2+0.1) {
219     // Select tracks coming to the crystal from up or down sides
220       Int_t pid = gMC->TrackPid();
221       Int_t module;
222       gMC->CurrentVolOffID(10,module);
223       module--;
224       AddImpact("EMC ",fIshunt, primary,tracknumber, module, pid, pmom, xyzm);
225     }
226   }
227
228   // Add impact to CPV
229
230   if( gMC->CurrentVolID(copy) == gMC->VolId("PCPQ") &&
231       gMC->IsTrackEntering() ) {
232     gMC->TrackMomentum(pmom);
233     gMC->TrackPosition(pos) ;
234
235     Int_t i;
236     for (i=0; i<3; i++) xyzm[i] = pos[i];
237
238     for (i=0; i<3; i++) {
239       xyzm[i] = pos[i] ;
240       pm[i]   = pmom[i];
241     }
242     Int_t pid = gMC->TrackPid();
243     Int_t module;
244     gMC->CurrentVolOffID(3,module);
245     module--;
246     AddImpact("CPV ",fIshunt, primary,tracknumber, module, pid, pmom, xyzm);
247   }
248   
249 }