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"
18 #include <G4TrackVector.hh>
19 #include <G4VUserTrackInformation.hh>
20 #include <G4TrackingManager.hh>
21 #include <G4SteppingManager.hh>
22 #include <G4UImanager.hh>
24 // static data members
25 AliTrackingAction* AliTrackingAction::fgInstance = 0;
27 //_____________________________________________________________________________
28 AliTrackingAction::AliTrackingAction()
32 fNewVerboseTrackID(-1),
39 AliGlobals::Exception("AliTrackingAction constructed twice.");
45 //_____________________________________________________________________________
46 AliTrackingAction::AliTrackingAction(const AliTrackingAction& right)
49 AliGlobals::Exception("AliTrackingAction is protected from copying.");
52 //_____________________________________________________________________________
53 AliTrackingAction::~AliTrackingAction() {
59 //_____________________________________________________________________________
61 AliTrackingAction::operator=(const AliTrackingAction &right)
63 // check assignement to self
64 if (this == &right) return *this;
66 AliGlobals::Exception("AliTrackingAction is protected from assigning.");
73 //_____________________________________________________________________________
74 AliTrackInformation* AliTrackingAction::GetTrackInformation(
76 const G4String& method) const
78 // Returns user track information.
82 G4VUserTrackInformation* trackInfo = track->GetUserInformation();
83 if (!trackInfo) return 0;
85 AliTrackInformation* aliTrackInfo
86 = dynamic_cast<AliTrackInformation*>(trackInfo);
88 G4String text = "AliTrackingAction::" + method + ":\n";
89 text = text + " Unknown track information type";
90 AliGlobals::Exception(text);
95 return (AliTrackInformation*)track->GetUserInformation();
101 //_____________________________________________________________________________
102 void AliTrackingAction::PrepareNewEvent()
104 // Called by G4 kernel at the beginning of event.
109 // set g4 stepping manager pointer
110 TG4StepManager* stepManager = TG4StepManager::Instance();
111 stepManager->SetSteppingManager(fpTrackingManager->GetSteppingManager());
114 //_____________________________________________________________________________
115 void AliTrackingAction::PreTrackingAction(const G4Track* aTrack)
117 // Called by G4 kernel before starting tracking.
120 // track index in the particles array
121 G4int trackID = aTrack->GetTrackID();
122 G4int parentID = aTrack->GetParentID();
125 // in AliRoot (from V3.0) track numbering starts from 0
126 trackIndex = trackID-1;
129 trackIndex = gAlice->GetNtrack();
132 // set track index to track information
133 AliTrackInformation* trackInfo
134 = GetTrackInformation(aTrack, "PreTrackingAction");
136 // create track information and set it to G4Track
137 // if it does not yet exist
138 trackInfo = new AliTrackInformation(trackIndex);
139 fpTrackingManager->SetUserTrackInformation(trackInfo);
140 // the track information is deleted together with its
144 trackInfo->SetTrackParticleID(trackIndex);
146 // set current track number
147 gAlice->SetCurrentTrack(trackIndex);
150 // finish previous primary track
151 FinishPrimaryTrack();
152 fPrimaryTrackID = aTrack->GetTrackID();
154 // begin this primary track
155 gAlice->BeginPrimary();
158 // save secondary particles info
163 if (trackID == fNewVerboseTrackID) {
164 G4String command = "/tracking/verbose ";
165 AliGlobals::AppendNumberToString(command, fNewVerboseLevel);
166 G4UImanager::GetUIpointer()->ApplyCommand(command);
169 // aliroot pre track actions
173 //_____________________________________________________________________________
174 void AliTrackingAction::PostTrackingAction(const G4Track* aTrack)
176 // Called by G4 kernel after finishing tracking.
181 // set parent track particle index to all secondary tracks
182 G4TrackVector* secondaryTracks
183 = fpTrackingManager->GetSteppingManager()->GetSecondary();
184 if (secondaryTracks){
186 for (i=0; i<secondaryTracks->size(); i++) {
187 G4Track* track = (*secondaryTracks)[i];
189 if (track->GetUserInformation()) {
190 // this should never happen
191 G4String text = "AliTrackingAction::PostTrackingAction:\n";
192 text = text + " Inconsistent track information.";
193 AliGlobals::Exception(text);
196 // get parent track index
197 AliTrackInformation* aliParentInfo
198 = GetTrackInformation(aTrack, "PostTrackingAction");
199 G4int parentParticleID
200 = aliParentInfo->GetTrackParticleID();
202 // create track information and set it to the G4Track
203 AliTrackInformation* trackInfo
204 = new AliTrackInformation(-1, parentParticleID);
205 track->SetUserInformation(trackInfo);
206 // the track information is deleted together with its
211 // aliroot post track actions
215 //_____________________________________________________________________________
216 void AliTrackingAction::FinishPrimaryTrack()
218 // Calls AliRun::PurifyKine and fills trees of hits
219 // after finishing tracking of each primary track.
220 // !! This method has to be also called from AlEventAction::EndOfEventAction()
221 // for storing the last primary track of the current event.
224 if (fPrimaryTrackID>0) {
227 if (fVerboseLevel == 3) {
228 G4cout << "$$$ Primary track " << fPrimaryTrackID << G4endl;
230 else if ( fVerboseLevel == 2 && fPrimaryTrackID % 10 == 0 ) {
231 G4cout << "$$$ Primary track " << fPrimaryTrackID << G4endl;
233 else if ( fVerboseLevel == 1 && fPrimaryTrackID % 100 == 0 ) {
234 G4cout << "$$$ Primary track " << fPrimaryTrackID << G4endl;
237 // aliroot finish primary track
238 gAlice->FinishPrimary();
243 //_____________________________________________________________________________
244 void AliTrackingAction::SaveTrack(const G4Track* track)
246 // Get all needed parameters from G4track and pass them
247 // to AliRun::SetTrack() that creates corresponding TParticle
248 // in the AliRun::fParticles array.
251 // parent particle index
252 G4int parentID = track->GetParentID();
259 = GetTrackInformation(track,"SaveTrack")->GetParentParticleID();
262 //G4cout << "SaveTrack: TrackID = " << track->GetTrackID()
263 // << " Parent ID = " << track->GetParentID()
264 // << " Index = " << gAlice->CurrentTrack()
265 // << " Parent Index = " << motherIndex
269 G4int pdg = track->GetDefinition()->GetPDGEncoding();
272 G4ThreeVector momentum = track->GetMomentum();
274 G4double px = momentum.x()/GeV;
275 G4double py = momentum.y()/GeV;
276 G4double pz = momentum.z()/GeV;
277 G4double e = track->GetTotalEnergy()/GeV;
279 G4ThreeVector position = track->GetPosition();
280 G4double vx = position.x()/cm;
281 G4double vy = position.y()/cm;
282 G4double vz = position.z()/cm;
283 // time of production - check if ekvivalent with G4
284 G4double t = track->GetGlobalTime();
286 G4ThreeVector polarization = track->GetPolarization();
287 G4double polX = polarization.x();
288 G4double polY = polarization.y();
289 G4double polZ = polarization.z();
291 // production process
292 AliMCProcess mcProcess;
293 const G4VProcess* kpProcess = track->GetCreatorProcess();
295 mcProcess = kPPrimary;
298 TG4PhysicsManager* pPhysicsManager = TG4PhysicsManager::Instance();
299 mcProcess = pPhysicsManager->GetMCProcess(kpProcess);
300 // distinguish kPDeltaRay from kPEnergyLoss
301 if (mcProcess == kPEnergyLoss) mcProcess = kPDeltaRay;
306 gAlice->SetTrack(1, motherIndex, pdg, px, py, pz, e, vx, vy, vz, t,
307 polX, polY, polZ, mcProcess, ntr);
311 //_____________________________________________________________________________
312 void AliTrackingAction::SetNewVerboseLevel(G4int level)
314 // Set the new verbose level that will be set when the track with
315 // specified track ID (fNewVerboseTrackID) starts tracking.
318 fNewVerboseLevel = level;
321 //_____________________________________________________________________________
322 void AliTrackingAction::SetNewVerboseTrackID(G4int trackID)
324 // Set the trackID for which the new verbose level (fNewVerboseLevel)
328 fNewVerboseTrackID = trackID;