]> git.uio.no Git - u/mrichter/AliRoot.git/blame - AliGeant4/AliTrackingAction.cxx
added inheritance from AliVerbose;
[u/mrichter/AliRoot.git] / AliGeant4 / AliTrackingAction.cxx
CommitLineData
676fb573 1// $Id$
2// Category: event
3//
7005154f 4// Author: I. Hrivnacova
5//
6// Class AliTrackingAction
7// -----------------------
676fb573 8// See the class description in the header file.
9
10#include "AliTrackingAction.h"
4240c1e8 11#include "AliTrackInformation.h"
676fb573 12#include "AliRun.h"
13#include "AliGlobals.h"
14#include "TG4StepManager.h"
c8bc1274 15#include "TG4PhysicsManager.h"
e7057a42 16#include "TG4ParticlesManager.h"
676fb573 17
676fb573 18#include <G4Track.hh>
4240c1e8 19#include <G4TrackVector.hh>
20#include <G4VUserTrackInformation.hh>
21#include <G4TrackingManager.hh>
22#include <G4SteppingManager.hh>
78ca1e9c 23#include <G4UImanager.hh>
676fb573 24
676fb573 25// static data members
26AliTrackingAction* AliTrackingAction::fgInstance = 0;
27
78ca1e9c 28//_____________________________________________________________________________
676fb573 29AliTrackingAction::AliTrackingAction()
4240c1e8 30 : fPrimaryTrackID(0),
676fb573 31 fVerboseLevel(2),
aafc96be 32 fNewVerboseLevel(0),
33 fNewVerboseTrackID(-1),
676fb573 34 fSavePrimaries(true),
7005154f 35 fTrackCounter(0),
36 fMessenger(this)
676fb573 37{
38//
39 if (fgInstance) {
40 AliGlobals::Exception("AliTrackingAction constructed twice.");
41 }
42
676fb573 43 fgInstance = this;
44}
45
78ca1e9c 46//_____________________________________________________________________________
7005154f 47AliTrackingAction::AliTrackingAction(const AliTrackingAction& right)
48 : fMessenger(this) {
676fb573 49//
50 AliGlobals::Exception("AliTrackingAction is protected from copying.");
51}
52
78ca1e9c 53//_____________________________________________________________________________
676fb573 54AliTrackingAction::~AliTrackingAction() {
55//
676fb573 56}
57
58// operators
59
78ca1e9c 60//_____________________________________________________________________________
676fb573 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
4240c1e8 72// private methods
73
78ca1e9c 74//_____________________________________________________________________________
4240c1e8 75AliTrackInformation* AliTrackingAction::GetTrackInformation(
76 const G4Track* track,
77 const G4String& method) const
78{
79// Returns user track information.
80// ---
81
7005154f 82#ifdef TGEANT4_DEBUG
4240c1e8 83 G4VUserTrackInformation* trackInfo = track->GetUserInformation();
84 if (!trackInfo) return 0;
85
86 AliTrackInformation* aliTrackInfo
87 = dynamic_cast<AliTrackInformation*>(trackInfo);
88 if (!aliTrackInfo) {
89 G4String text = "AliTrackingAction::" + method + ":\n";
90 text = text + " Unknown track information type";
91 AliGlobals::Exception(text);
92 }
93
94 return aliTrackInfo;
7005154f 95#else
96 return (AliTrackInformation*)track->GetUserInformation();
97#endif
4240c1e8 98}
99
676fb573 100// public methods
101
78ca1e9c 102//_____________________________________________________________________________
676fb573 103void AliTrackingAction::PrepareNewEvent()
104{
105// Called by G4 kernel at the beginning of event.
106// ---
107
3bfabcfc 108 fTrackCounter = 0;
109
676fb573 110 // set g4 stepping manager pointer
4240c1e8 111 TG4StepManager* stepManager = TG4StepManager::Instance();
112 stepManager->SetSteppingManager(fpTrackingManager->GetSteppingManager());
676fb573 113}
114
78ca1e9c 115//_____________________________________________________________________________
9bcb6317 116void AliTrackingAction::PreTrackingAction(const G4Track* aTrack)
676fb573 117{
118// Called by G4 kernel before starting tracking.
119// ---
120
676fb573 121 // track index in the particles array
122 G4int trackID = aTrack->GetTrackID();
123 G4int parentID = aTrack->GetParentID();
124 Int_t trackIndex;
125 if (parentID==0) {
4240c1e8 126 // in AliRoot (from V3.0) track numbering starts from 0
127 trackIndex = trackID-1;
676fb573 128 }
129 else {
4240c1e8 130 trackIndex = gAlice->GetNtrack();
676fb573 131 }
4240c1e8 132
133 // set track index to track information
134 AliTrackInformation* trackInfo
135 = GetTrackInformation(aTrack, "PreTrackingAction");
136 if (!trackInfo) {
137 // create track information and set it to G4Track
138 // if it does not yet exist
139 trackInfo = new AliTrackInformation(trackIndex);
140 fpTrackingManager->SetUserTrackInformation(trackInfo);
141 // the track information is deleted together with its
142 // G4Track object
143 }
144 else
145 trackInfo->SetTrackParticleID(trackIndex);
676fb573 146
4240c1e8 147 // set current track number
148 gAlice->SetCurrentTrack(trackIndex);
676fb573 149
150 if (parentID == 0) {
4240c1e8 151 // finish previous primary track
152 FinishPrimaryTrack();
676fb573 153 fPrimaryTrackID = aTrack->GetTrackID();
887f9ee9 154
155 // begin this primary track
156 gAlice->BeginPrimary();
676fb573 157 }
158 else {
159 // save secondary particles info
4240c1e8 160 SaveTrack(aTrack);
9bcb6317 161 }
b17c5755 162
aafc96be 163 // verbose
164 if (trackID == fNewVerboseTrackID) {
165 G4String command = "/tracking/verbose ";
166 AliGlobals::AppendNumberToString(command, fNewVerboseLevel);
167 G4UImanager::GetUIpointer()->ApplyCommand(command);
168 }
169
4240c1e8 170 // aliroot pre track actions
b17c5755 171 gAlice->PreTrack();
676fb573 172}
173
78ca1e9c 174//_____________________________________________________________________________
9bcb6317 175void AliTrackingAction::PostTrackingAction(const G4Track* aTrack)
676fb573 176{
177// Called by G4 kernel after finishing tracking.
178// ---
179
3bfabcfc 180 fTrackCounter++;
4240c1e8 181
182 // set parent track particle index to all secondary tracks
183 G4TrackVector* secondaryTracks
184 = fpTrackingManager->GetSteppingManager()->GetSecondary();
185 if (secondaryTracks){
186 G4int i;
ff849423 187 for (i=0; i<secondaryTracks->size(); i++) {
4240c1e8 188 G4Track* track = (*secondaryTracks)[i];
189
190 if (track->GetUserInformation()) {
191 // this should never happen
192 G4String text = "AliTrackingAction::PostTrackingAction:\n";
193 text = text + " Inconsistent track information.";
194 AliGlobals::Exception(text);
195 }
196
197 // get parent track index
198 AliTrackInformation* aliParentInfo
199 = GetTrackInformation(aTrack, "PostTrackingAction");
200 G4int parentParticleID
201 = aliParentInfo->GetTrackParticleID();
202
203 // create track information and set it to the G4Track
204 AliTrackInformation* trackInfo
205 = new AliTrackInformation(-1, parentParticleID);
206 track->SetUserInformation(trackInfo);
207 // the track information is deleted together with its
208 // G4Track object
209 }
210 }
211
212 // aliroot post track actions
b17c5755 213 gAlice->PostTrack();
676fb573 214}
215
78ca1e9c 216//_____________________________________________________________________________
4240c1e8 217void AliTrackingAction::FinishPrimaryTrack()
676fb573 218{
219// Calls AliRun::PurifyKine and fills trees of hits
220// after finishing tracking of each primary track.
221// !! This method has to be also called from AlEventAction::EndOfEventAction()
222// for storing the last primary track of the current event.
223// ---
224
4240c1e8 225 if (fPrimaryTrackID>0) {
226
227 // verbose
228 if (fVerboseLevel == 3) {
229 G4cout << "$$$ Primary track " << fPrimaryTrackID << G4endl;
230 }
231 else if ( fVerboseLevel == 2 && fPrimaryTrackID % 10 == 0 ) {
232 G4cout << "$$$ Primary track " << fPrimaryTrackID << G4endl;
233 }
234 else if ( fVerboseLevel == 1 && fPrimaryTrackID % 100 == 0 ) {
235 G4cout << "$$$ Primary track " << fPrimaryTrackID << G4endl;
236 }
237
238 // aliroot finish primary track
239 gAlice->FinishPrimary();
240 }
241 fPrimaryTrackID = 0;
676fb573 242}
243
78ca1e9c 244//_____________________________________________________________________________
4240c1e8 245void AliTrackingAction::SaveTrack(const G4Track* track)
676fb573 246{
4240c1e8 247// Get all needed parameters from G4track and pass them
248// to AliRun::SetTrack() that creates corresponding TParticle
249// in the AliRun::fParticles array.
676fb573 250// ----
251
4240c1e8 252 // parent particle index
676fb573 253 G4int parentID = track->GetParentID();
4240c1e8 254 G4int motherIndex;
676fb573 255 if (parentID == 0) {
4240c1e8 256 motherIndex = -1;
676fb573 257 }
258 else {
4240c1e8 259 motherIndex
260 = GetTrackInformation(track,"SaveTrack")->GetParentParticleID();
261 }
676fb573 262
4240c1e8 263 //G4cout << "SaveTrack: TrackID = " << track->GetTrackID()
264 // << " Parent ID = " << track->GetParentID()
265 // << " Index = " << gAlice->CurrentTrack()
266 // << " Parent Index = " << motherIndex
267 // << G4endl;
268
269 // PDG code
e7057a42 270 G4int pdg
271 = TG4ParticlesManager::Instance()
272 ->GetPDGEncodingFast(track->GetDefinition());
676fb573 273
4240c1e8 274 // track kinematics
676fb573 275 G4ThreeVector momentum = track->GetMomentum();
4240c1e8 276
676fb573 277 G4double px = momentum.x()/GeV;
278 G4double py = momentum.y()/GeV;
279 G4double pz = momentum.z()/GeV;
280 G4double e = track->GetTotalEnergy()/GeV;
281
282 G4ThreeVector position = track->GetPosition();
283 G4double vx = position.x()/cm;
284 G4double vy = position.y()/cm;
285 G4double vz = position.z()/cm;
286 // time of production - check if ekvivalent with G4
287 G4double t = track->GetGlobalTime();
288
289 G4ThreeVector polarization = track->GetPolarization();
290 G4double polX = polarization.x();
291 G4double polY = polarization.y();
292 G4double polZ = polarization.z();
293
4240c1e8 294 // production process
c8bc1274 295 AliMCProcess mcProcess;
296 const G4VProcess* kpProcess = track->GetCreatorProcess();
297 if (!kpProcess) {
298 mcProcess = kPPrimary;
299 }
300 else {
e7057a42 301 mcProcess = TG4PhysicsManager::Instance()->GetMCProcess(kpProcess);
c8bc1274 302 // distinguish kPDeltaRay from kPEnergyLoss
303 if (mcProcess == kPEnergyLoss) mcProcess = kPDeltaRay;
304 }
4240c1e8 305
306 G4int ntr;
307 // create particle
308 gAlice->SetTrack(1, motherIndex, pdg, px, py, pz, e, vx, vy, vz, t,
309 polX, polY, polZ, mcProcess, ntr);
310
676fb573 311}
312
78ca1e9c 313//_____________________________________________________________________________
aafc96be 314void AliTrackingAction::SetNewVerboseLevel(G4int level)
315{
316// Set the new verbose level that will be set when the track with
317// specified track ID (fNewVerboseTrackID) starts tracking.
318// ---
319
320 fNewVerboseLevel = level;
321}
322
78ca1e9c 323//_____________________________________________________________________________
aafc96be 324void AliTrackingAction::SetNewVerboseTrackID(G4int trackID)
325{
326// Set the trackID for which the new verbose level (fNewVerboseLevel)
327// will be applied.
328// ---
329
330 fNewVerboseTrackID = trackID;
331}