4 // Author: I. Hrivnacova
6 // Class AliTrackingAction
7 // -----------------------
8 // See the class description in the header file.
10 #include "AliTrackingAction.h"
11 #include "AliTrackInformation.h"
13 #include "AliGlobals.h"
14 #include "TG4StepManager.h"
15 #include "TG4PhysicsManager.h"
16 #include "TG4ParticlesManager.h"
19 #include <G4TrackVector.hh>
20 #include <G4VUserTrackInformation.hh>
21 #include <G4TrackingManager.hh>
22 #include <G4SteppingManager.hh>
23 #include <G4UImanager.hh>
25 // static data members
26 AliTrackingAction* AliTrackingAction::fgInstance = 0;
28 //_____________________________________________________________________________
29 AliTrackingAction::AliTrackingAction()
33 fNewVerboseTrackID(-1),
40 AliGlobals::Exception("AliTrackingAction constructed twice.");
46 //_____________________________________________________________________________
47 AliTrackingAction::AliTrackingAction(const AliTrackingAction& right)
50 AliGlobals::Exception("AliTrackingAction is protected from copying.");
53 //_____________________________________________________________________________
54 AliTrackingAction::~AliTrackingAction() {
60 //_____________________________________________________________________________
62 AliTrackingAction::operator=(const AliTrackingAction &right)
64 // check assignement to self
65 if (this == &right) return *this;
67 AliGlobals::Exception("AliTrackingAction is protected from assigning.");
74 //_____________________________________________________________________________
75 AliTrackInformation* AliTrackingAction::GetTrackInformation(
77 const G4String& method) const
79 // Returns user track information.
83 G4VUserTrackInformation* trackInfo = track->GetUserInformation();
84 if (!trackInfo) return 0;
86 AliTrackInformation* aliTrackInfo
87 = dynamic_cast<AliTrackInformation*>(trackInfo);
89 G4String text = "AliTrackingAction::" + method + ":\n";
90 text = text + " Unknown track information type";
91 AliGlobals::Exception(text);
96 return (AliTrackInformation*)track->GetUserInformation();
102 //_____________________________________________________________________________
103 void AliTrackingAction::PrepareNewEvent()
105 // Called by G4 kernel at the beginning of event.
110 // set g4 stepping manager pointer
111 TG4StepManager* stepManager = TG4StepManager::Instance();
112 stepManager->SetSteppingManager(fpTrackingManager->GetSteppingManager());
115 //_____________________________________________________________________________
116 void AliTrackingAction::PreTrackingAction(const G4Track* aTrack)
118 // Called by G4 kernel before starting tracking.
121 // track index in the particles array
122 G4int trackID = aTrack->GetTrackID();
123 G4int parentID = aTrack->GetParentID();
126 // in AliRoot (from V3.0) track numbering starts from 0
127 trackIndex = trackID-1;
130 trackIndex = gAlice->GetNtrack();
133 // set track index to track information
134 AliTrackInformation* trackInfo
135 = GetTrackInformation(aTrack, "PreTrackingAction");
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
145 trackInfo->SetTrackParticleID(trackIndex);
147 // set current track number
148 gAlice->SetCurrentTrack(trackIndex);
151 // finish previous primary track
152 FinishPrimaryTrack();
153 fPrimaryTrackID = aTrack->GetTrackID();
155 // begin this primary track
156 gAlice->BeginPrimary();
159 // save secondary particles info
164 if (trackID == fNewVerboseTrackID) {
165 G4String command = "/tracking/verbose ";
166 AliGlobals::AppendNumberToString(command, fNewVerboseLevel);
167 G4UImanager::GetUIpointer()->ApplyCommand(command);
170 // aliroot pre track actions
174 //_____________________________________________________________________________
175 void AliTrackingAction::PostTrackingAction(const G4Track* aTrack)
177 // Called by G4 kernel after finishing tracking.
182 // set parent track particle index to all secondary tracks
183 G4TrackVector* secondaryTracks
184 = fpTrackingManager->GetSteppingManager()->GetSecondary();
185 if (secondaryTracks){
187 for (i=0; i<secondaryTracks->size(); i++) {
188 G4Track* track = (*secondaryTracks)[i];
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);
197 // get parent track index
198 AliTrackInformation* aliParentInfo
199 = GetTrackInformation(aTrack, "PostTrackingAction");
200 G4int parentParticleID
201 = aliParentInfo->GetTrackParticleID();
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
212 // aliroot post track actions
216 //_____________________________________________________________________________
217 void AliTrackingAction::FinishPrimaryTrack()
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.
225 if (fPrimaryTrackID>0) {
228 if (fVerboseLevel == 3) {
229 G4cout << "$$$ Primary track " << fPrimaryTrackID << G4endl;
231 else if ( fVerboseLevel == 2 && fPrimaryTrackID % 10 == 0 ) {
232 G4cout << "$$$ Primary track " << fPrimaryTrackID << G4endl;
234 else if ( fVerboseLevel == 1 && fPrimaryTrackID % 100 == 0 ) {
235 G4cout << "$$$ Primary track " << fPrimaryTrackID << G4endl;
238 // aliroot finish primary track
239 gAlice->FinishPrimary();
244 //_____________________________________________________________________________
245 void AliTrackingAction::SaveTrack(const G4Track* track)
247 // Get all needed parameters from G4track and pass them
248 // to AliRun::SetTrack() that creates corresponding TParticle
249 // in the AliRun::fParticles array.
252 // parent particle index
253 G4int parentID = track->GetParentID();
260 = GetTrackInformation(track,"SaveTrack")->GetParentParticleID();
263 //G4cout << "SaveTrack: TrackID = " << track->GetTrackID()
264 // << " Parent ID = " << track->GetParentID()
265 // << " Index = " << gAlice->CurrentTrack()
266 // << " Parent Index = " << motherIndex
271 = TG4ParticlesManager::Instance()
272 ->GetPDGEncodingFast(track->GetDefinition());
275 G4ThreeVector momentum = track->GetMomentum();
277 G4double px = momentum.x()/GeV;
278 G4double py = momentum.y()/GeV;
279 G4double pz = momentum.z()/GeV;
280 G4double e = track->GetTotalEnergy()/GeV;
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();
289 G4ThreeVector polarization = track->GetPolarization();
290 G4double polX = polarization.x();
291 G4double polY = polarization.y();
292 G4double polZ = polarization.z();
294 // production process
295 AliMCProcess mcProcess;
296 const G4VProcess* kpProcess = track->GetCreatorProcess();
298 mcProcess = kPPrimary;
301 mcProcess = TG4PhysicsManager::Instance()->GetMCProcess(kpProcess);
302 // distinguish kPDeltaRay from kPEnergyLoss
303 if (mcProcess == kPEnergyLoss) mcProcess = kPDeltaRay;
308 gAlice->SetTrack(1, motherIndex, pdg, px, py, pz, e, vx, vy, vz, t,
309 polX, polY, polZ, mcProcess, ntr);
313 //_____________________________________________________________________________
314 void AliTrackingAction::SetNewVerboseLevel(G4int level)
316 // Set the new verbose level that will be set when the track with
317 // specified track ID (fNewVerboseTrackID) starts tracking.
320 fNewVerboseLevel = level;
323 //_____________________________________________________________________________
324 void AliTrackingAction::SetNewVerboseTrackID(G4int trackID)
326 // Set the trackID for which the new verbose level (fNewVerboseLevel)
330 fNewVerboseTrackID = trackID;