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