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