]> git.uio.no Git - u/mrichter/AliRoot.git/blame - AliGeant4/AliTrackingAction.cxx
added calls to AliRun::Pre/PostTrack
[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"
8#include "AliRun.h"
9#include "AliGlobals.h"
10#include "TG4StepManager.h"
11
12#include <G4TrackingManager.hh>
13#include <G4Track.hh>
676fb573 14
71c03be6 15#include <TTree.h>
676fb573 16#include <TParticle.h>
c97337f9 17#include <TClonesArray.h>
676fb573 18
19// static data members
20AliTrackingAction* AliTrackingAction::fgInstance = 0;
21
22// static methods
23AliTrackingAction* AliTrackingAction::Instance() {
24//
25 return fgInstance;
26}
27
28AliTrackingAction::AliTrackingAction()
29 : fParticles(0),
30 fPrimaryTrackID(0),
31 fVerboseLevel(2),
32 fSavePrimaries(true),
33 fPrimariesCounter(0),
34 fParticlesCounter(0),
3bfabcfc 35 fTrackCounter(0),
676fb573 36 fLastParticleIndex(-1)
37{
38//
39 if (fgInstance) {
40 AliGlobals::Exception("AliTrackingAction constructed twice.");
41 }
42
43 fMessenger = new AliTrackingActionMessenger(this);
44 fgInstance = this;
45}
46
47AliTrackingAction::AliTrackingAction(const AliTrackingAction& right) {
48//
49 AliGlobals::Exception("AliTrackingAction is protected from copying.");
50}
51
52AliTrackingAction::~AliTrackingAction() {
53//
54 delete fMessenger;
55}
56
57// operators
58
59AliTrackingAction&
60AliTrackingAction::operator=(const AliTrackingAction &right)
61{
62 // check assignement to self
63 if (this == &right) return *this;
64
65 AliGlobals::Exception("AliTrackingAction is protected from assigning.");
66
67 return *this;
68}
69
70// public methods
71
72void AliTrackingAction::PrepareNewEvent()
73{
74// Called by G4 kernel at the beginning of event.
75// ---
76
3bfabcfc 77 fTrackCounter = 0;
78
676fb573 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 }
b17c5755 131
132 gAlice->PreTrack();
676fb573 133}
134
9bcb6317 135void AliTrackingAction::PostTrackingAction(const G4Track* aTrack)
676fb573 136{
137// Called by G4 kernel after finishing tracking.
138// ---
139
3bfabcfc 140 fTrackCounter++;
b17c5755 141
142 gAlice->PostTrack();
676fb573 143}
144
145void AliTrackingAction::SaveAndDestroyTrack()
146{
147// Calls AliRun::PurifyKine and fills trees of hits
148// after finishing tracking of each primary track.
149// !! This method has to be also called from AlEventAction::EndOfEventAction()
150// for storing the last primary track of the current event.
151// ---
152
153 if (fPrimaryTrackID>0)
154 {
155 if (fVerboseLevel == 3) {
5f1d09c5 156 G4cout << "$$$ Primary track " << fPrimaryTrackID << G4endl;
676fb573 157 }
158 else if ( fVerboseLevel == 2 && fPrimaryTrackID % 10 == 0 ) {
5f1d09c5 159 G4cout << "$$$ Primary track " << fPrimaryTrackID << G4endl;
676fb573 160 }
161 else if ( fVerboseLevel == 1 && fPrimaryTrackID % 100 == 0 ) {
5f1d09c5 162 G4cout << "$$$ Primary track " << fPrimaryTrackID << G4endl;
676fb573 163 }
164
165 // aliroot
166 G4int lastSavedTrack
167 = gAlice->PurifyKine(fLastParticleIndex, fParticlesCounter);
168 G4int nofPuredSavedTracks
169 = gAlice->GetNtrack();
170 fLastParticleIndex = lastSavedTrack;
171 fParticlesCounter = nofPuredSavedTracks;
172
173 if(gAlice->TreeH()) gAlice->TreeH()->Fill();
174 gAlice->ResetHits();
175 }
176 else {
177 fLastParticleIndex = fPrimariesCounter-1;
178 }
179 fPrimaryTrackID = 0;
180}
181
182void AliTrackingAction::SaveParticle(const G4Track* track,
183 G4String processName)
184{
185// Converts G4track to TParticle and saves it in AliRun::fParticles
186// array.
187// ----
188
189 fParticlesCounter++;
190
191 // track history
192 G4int firstDaughter = -1;
193 G4int lastDaughter = -1;
194 G4int parentID = track->GetParentID();
195 G4int motherIndex1;
196 G4int motherIndex2 = -1;
197 if (parentID == 0) {
198 motherIndex1 = -1;
199 fPrimariesCounter++;
200 }
201 else {
202 // set first/last child for already defined particles
203 motherIndex1 = GetParticleIndex(parentID);
204 // aliroot
205 TParticle* parent
206 = dynamic_cast<TParticle*>(fParticles->UncheckedAt(motherIndex1));
207 if (parent) {
208 if (parent->GetFirstDaughter()<0)
209 parent->SetFirstDaughter(fParticlesCounter);
210 parent->SetLastDaughter(fParticlesCounter);
211 }
212 else {
213 AliGlobals::Exception(
214 "AliTrackingAction::SaveParticle: Unknown particle type");
215 }
216 };
217
218 // ?? status of particle
219 // temporarily used for storing trackID from G4
220 G4int ks = track->GetTrackID();
221
222 // particle type
223 // is this G3 ekvivalent ?? - check
224 G4int pdg = track->GetDefinition()->GetPDGEncoding();
225
226 // track kinematics
227 G4ThreeVector momentum = track->GetMomentum();
228 G4double px = momentum.x()/GeV;
229 G4double py = momentum.y()/GeV;
230 G4double pz = momentum.z()/GeV;
231 G4double e = track->GetTotalEnergy()/GeV;
232
233 G4ThreeVector position = track->GetPosition();
234 G4double vx = position.x()/cm;
235 G4double vy = position.y()/cm;
236 G4double vz = position.z()/cm;
237 // time of production - check if ekvivalent with G4
238 G4double t = track->GetGlobalTime();
239
240 G4ThreeVector polarization = track->GetPolarization();
241 G4double polX = polarization.x();
242 G4double polY = polarization.y();
243 G4double polZ = polarization.z();
244
245 // aliroot
246 TClonesArray& theCollectionRef = *fParticles;
247 G4int nofParticles = theCollectionRef.GetEntriesFast();
248 TParticle* particle
249 = new(theCollectionRef[nofParticles])
250 TParticle(pdg, ks, motherIndex1, motherIndex1,
251 firstDaughter, lastDaughter, px, py, pz, e, vx, vy, vz, t);
252 particle->SetPolarisation(polX, polY, polZ);
b65f8990 253 particle->SetBit(kKeepBit, false);
676fb573 254}
255
256G4int AliTrackingAction::GetParticleIndex(G4int trackID)
257{
258// Converts the trackID into the index of the particle
259// in AliRun::fParticles array.
260// ---
261
262 if (trackID <= fPrimariesCounter) {
263 return trackID-1;
264 }
265 else
266 for (G4int i=fLastParticleIndex+1; i<fParticlesCounter; i++) {
267 // aliroot
268 TParticle* particle
269 = (TParticle*)fParticles->UncheckedAt(i);
270 if (trackID == particle->GetStatusCode()) return i;
271 }
272
273 AliGlobals::Exception("AliTrackingAction::GetParticleIndex() has failed.");
274 return 0;
275}