]> git.uio.no Git - u/mrichter/AliRoot.git/blame - AliGeant4/AliTrackingAction.cxx
creating mag. field moved from constructor to SetMagField()
[u/mrichter/AliRoot.git] / AliGeant4 / AliTrackingAction.cxx
CommitLineData
676fb573 1// $Id$
2// Category: event
3//
4// See the class description in the header file.
5
6#include "AliTrackingAction.h"
7#include "AliTrackingActionMessenger.h"
9bcb6317 8#include "AliSensitiveDetector.h"
676fb573 9#include "AliRun.h"
10#include "AliGlobals.h"
11#include "TG4StepManager.h"
12
13#include <G4TrackingManager.hh>
14#include <G4Track.hh>
15#include <G4Event.hh>
16#include <G4SDManager.hh>
9bcb6317 17#include <G4VSensitiveDetector.hh>
676fb573 18#include <G4VHitsCollection.hh>
19
20#include <TParticle.h>
21
22// static data members
23AliTrackingAction* AliTrackingAction::fgInstance = 0;
24
25// static methods
26AliTrackingAction* AliTrackingAction::Instance() {
27//
28 return fgInstance;
29}
30
31AliTrackingAction::AliTrackingAction()
32 : fParticles(0),
33 fPrimaryTrackID(0),
34 fVerboseLevel(2),
35 fSavePrimaries(true),
36 fPrimariesCounter(0),
37 fParticlesCounter(0),
38 fLastParticleIndex(-1)
39{
40//
41 if (fgInstance) {
42 AliGlobals::Exception("AliTrackingAction constructed twice.");
43 }
44
45 fMessenger = new AliTrackingActionMessenger(this);
46 fgInstance = this;
47}
48
49AliTrackingAction::AliTrackingAction(const AliTrackingAction& right) {
50//
51 AliGlobals::Exception("AliTrackingAction is protected from copying.");
52}
53
54AliTrackingAction::~AliTrackingAction() {
55//
56 delete fMessenger;
57}
58
59// operators
60
61AliTrackingAction&
62AliTrackingAction::operator=(const AliTrackingAction &right)
63{
64 // check assignement to self
65 if (this == &right) return *this;
66
67 AliGlobals::Exception("AliTrackingAction is protected from assigning.");
68
69 return *this;
70}
71
72// public methods
73
74void AliTrackingAction::PrepareNewEvent()
75{
76// Called by G4 kernel at the beginning of event.
77// ---
78
79 // aliroot
80 if (!fParticles) fParticles = gAlice->Particles();
81
82 // set the particles and primaries counter to the
83 // number of tracks already stored in AliRun::fParticles
84 // (fParticles can already contain primary particles
85 // saved by AliGenerator)
86 G4int nofTracks = gAlice->GetNtrack();
87 fPrimariesCounter = nofTracks;
88 fParticlesCounter = nofTracks;
89 fLastParticleIndex = nofTracks-1;
90
91 // set g4 stepping manager pointer
92 G4SteppingManager* pG4StepManager
93 = fpTrackingManager->GetSteppingManager();
94 TG4StepManager* pStepManager = TG4StepManager::Instance();
95 pStepManager->SetSteppingManager(pG4StepManager);
96}
97
9bcb6317 98void AliTrackingAction::PreTrackingAction(const G4Track* aTrack)
676fb573 99{
100// Called by G4 kernel before starting tracking.
101// ---
102
103 // aliroot
104 // track index in the particles array
105 G4int trackID = aTrack->GetTrackID();
106 G4int parentID = aTrack->GetParentID();
107 Int_t trackIndex;
108 if (parentID==0) {
109 trackIndex = trackID;
110 }
111 else {
112 trackIndex = GetNofSavedTracks();
113 }
114
115 // in AliRoot (from V3.0) track numbering starts from 0
116 gAlice->SetCurrentTrack(--trackIndex);
117
118 if (parentID == 0) {
119 // save primary track
120 SaveAndDestroyTrack();
121 fPrimaryTrackID = aTrack->GetTrackID();
122 }
123 else {
124 // save secondary particles info
125 // improve this later with retrieving the generation process
126 // (primary particles are stored
127 // by AlStackingAction in ClassifyNewTrack() method)
128 G4String origin = "secondary";
129 SaveParticle(aTrack, origin);
9bcb6317 130 }
676fb573 131}
132
9bcb6317 133void AliTrackingAction::PostTrackingAction(const G4Track* aTrack)
676fb573 134{
135// Called by G4 kernel after finishing tracking.
136// ---
137
138 G4String particleName
139 = aTrack->GetDynamicParticle()->GetDefinition()->GetParticleName();
140 if (particleName == "opticalphoton") {
141 G4cout << "$$$ Track " << aTrack->GetTrackID()
5f1d09c5 142 << " is optical photon." << G4endl;
676fb573 143 }
144}
145
146void AliTrackingAction::SaveAndDestroyTrack()
147{
148// Calls AliRun::PurifyKine and fills trees of hits
149// after finishing tracking of each primary track.
150// !! This method has to be also called from AlEventAction::EndOfEventAction()
151// for storing the last primary track of the current event.
152// ---
153
154 if (fPrimaryTrackID>0)
155 {
156 if (fVerboseLevel == 3) {
5f1d09c5 157 G4cout << "$$$ Primary track " << fPrimaryTrackID << G4endl;
676fb573 158 }
159 else if ( fVerboseLevel == 2 && fPrimaryTrackID % 10 == 0 ) {
5f1d09c5 160 G4cout << "$$$ Primary track " << fPrimaryTrackID << G4endl;
676fb573 161 }
162 else if ( fVerboseLevel == 1 && fPrimaryTrackID % 100 == 0 ) {
5f1d09c5 163 G4cout << "$$$ Primary track " << fPrimaryTrackID << G4endl;
676fb573 164 }
165
166 // aliroot
167 G4int lastSavedTrack
168 = gAlice->PurifyKine(fLastParticleIndex, fParticlesCounter);
169 G4int nofPuredSavedTracks
170 = gAlice->GetNtrack();
171 fLastParticleIndex = lastSavedTrack;
172 fParticlesCounter = nofPuredSavedTracks;
173
174 if(gAlice->TreeH()) gAlice->TreeH()->Fill();
175 gAlice->ResetHits();
176 }
177 else {
178 fLastParticleIndex = fPrimariesCounter-1;
179 }
180 fPrimaryTrackID = 0;
181}
182
183void AliTrackingAction::SaveParticle(const G4Track* track,
184 G4String processName)
185{
186// Converts G4track to TParticle and saves it in AliRun::fParticles
187// array.
188// ----
189
190 fParticlesCounter++;
191
192 // track history
193 G4int firstDaughter = -1;
194 G4int lastDaughter = -1;
195 G4int parentID = track->GetParentID();
196 G4int motherIndex1;
197 G4int motherIndex2 = -1;
198 if (parentID == 0) {
199 motherIndex1 = -1;
200 fPrimariesCounter++;
201 }
202 else {
203 // set first/last child for already defined particles
204 motherIndex1 = GetParticleIndex(parentID);
205 // aliroot
206 TParticle* parent
207 = dynamic_cast<TParticle*>(fParticles->UncheckedAt(motherIndex1));
208 if (parent) {
209 if (parent->GetFirstDaughter()<0)
210 parent->SetFirstDaughter(fParticlesCounter);
211 parent->SetLastDaughter(fParticlesCounter);
212 }
213 else {
214 AliGlobals::Exception(
215 "AliTrackingAction::SaveParticle: Unknown particle type");
216 }
217 };
218
219 // ?? status of particle
220 // temporarily used for storing trackID from G4
221 G4int ks = track->GetTrackID();
222
223 // particle type
224 // is this G3 ekvivalent ?? - check
225 G4int pdg = track->GetDefinition()->GetPDGEncoding();
226
227 // track kinematics
228 G4ThreeVector momentum = track->GetMomentum();
229 G4double px = momentum.x()/GeV;
230 G4double py = momentum.y()/GeV;
231 G4double pz = momentum.z()/GeV;
232 G4double e = track->GetTotalEnergy()/GeV;
233
234 G4ThreeVector position = track->GetPosition();
235 G4double vx = position.x()/cm;
236 G4double vy = position.y()/cm;
237 G4double vz = position.z()/cm;
238 // time of production - check if ekvivalent with G4
239 G4double t = track->GetGlobalTime();
240
241 G4ThreeVector polarization = track->GetPolarization();
242 G4double polX = polarization.x();
243 G4double polY = polarization.y();
244 G4double polZ = polarization.z();
245
246 // aliroot
247 TClonesArray& theCollectionRef = *fParticles;
248 G4int nofParticles = theCollectionRef.GetEntriesFast();
249 TParticle* particle
250 = new(theCollectionRef[nofParticles])
251 TParticle(pdg, ks, motherIndex1, motherIndex1,
252 firstDaughter, lastDaughter, px, py, pz, e, vx, vy, vz, t);
253 particle->SetPolarisation(polX, polY, polZ);
b65f8990 254 particle->SetBit(kKeepBit, false);
676fb573 255}
256
257G4int AliTrackingAction::GetParticleIndex(G4int trackID)
258{
259// Converts the trackID into the index of the particle
260// in AliRun::fParticles array.
261// ---
262
263 if (trackID <= fPrimariesCounter) {
264 return trackID-1;
265 }
266 else
267 for (G4int i=fLastParticleIndex+1; i<fParticlesCounter; i++) {
268 // aliroot
269 TParticle* particle
270 = (TParticle*)fParticles->UncheckedAt(i);
271 if (trackID == particle->GetStatusCode()) return i;
272 }
273
274 AliGlobals::Exception("AliTrackingAction::GetParticleIndex() has failed.");
275 return 0;
276}