]> git.uio.no Git - u/mrichter/AliRoot.git/blob - AliGeant4/AliTrackingAction.cxx
updated commands description
[u/mrichter/AliRoot.git] / AliGeant4 / AliTrackingAction.cxx
1 // $Id$
2 // Category: event
3 //
4 // Author: I. Hrivnacova
5 //
6 // Class AliTrackingAction
7 // -----------------------
8 // See the class description in the header file.
9
10 #include "AliTrackingAction.h"
11 #include "AliTrackInformation.h"
12 #include "AliRun.h"
13 #include "AliGlobals.h"  
14 #include "TG4StepManager.h"
15 #include "TG4PhysicsManager.h"
16 #include "TG4ParticlesManager.h"
17
18 #include <G4Track.hh>
19 #include <G4TrackVector.hh>
20 #include <G4VUserTrackInformation.hh>
21 #include <G4TrackingManager.hh>
22 #include <G4SteppingManager.hh>
23 #include <G4UImanager.hh>
24
25 // static data members
26 AliTrackingAction* AliTrackingAction::fgInstance = 0;
27
28 //_____________________________________________________________________________
29 AliTrackingAction::AliTrackingAction()
30   : fPrimaryTrackID(0),
31     fVerboseLevel(2),
32     fNewVerboseLevel(0),
33     fNewVerboseTrackID(-1),
34     fSavePrimaries(true),
35     fTrackCounter(0),
36     fMessenger(this)
37 {
38 //
39   if (fgInstance) { 
40     AliGlobals::Exception("AliTrackingAction constructed twice."); 
41   }
42
43   fgInstance = this;
44 }
45
46 //_____________________________________________________________________________
47 AliTrackingAction::AliTrackingAction(const AliTrackingAction& right) 
48   : fMessenger(this) {
49 //
50   AliGlobals::Exception("AliTrackingAction is protected from copying.");
51 }
52
53 //_____________________________________________________________________________
54 AliTrackingAction::~AliTrackingAction() {
55 //
56 }
57
58 // operators
59
60 //_____________________________________________________________________________
61 AliTrackingAction& 
62 AliTrackingAction::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 // private methods
73
74 //_____________________________________________________________________________
75 AliTrackInformation* AliTrackingAction::GetTrackInformation(
76                                            const G4Track* track,
77                                            const G4String& method) const
78 {
79 // Returns user track information.
80 // ---
81  
82 #ifdef TGEANT4_DEBUG
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;
95 #else  
96   return (AliTrackInformation*)track->GetUserInformation();
97 #endif  
98 }    
99   
100 // public methods
101
102 //_____________________________________________________________________________
103 void AliTrackingAction::PrepareNewEvent()
104 {
105 // Called by G4 kernel at the beginning of event.
106 // ---
107
108   fTrackCounter = 0;
109
110   // set g4 stepping manager pointer
111   TG4StepManager* stepManager = TG4StepManager::Instance();
112   stepManager->SetSteppingManager(fpTrackingManager->GetSteppingManager());
113 }
114
115 //_____________________________________________________________________________
116 void AliTrackingAction::PreTrackingAction(const G4Track* aTrack)
117 {
118 // Called by G4 kernel before starting tracking.
119 // ---
120
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) { 
126     // in AliRoot (from V3.0) track numbering starts from 0
127     trackIndex = trackID-1; 
128   } 
129   else { 
130     trackIndex = gAlice->GetNtrack();
131   }
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);  
146
147   // set current track number
148   gAlice->SetCurrentTrack(trackIndex);
149
150   if (parentID == 0) {  
151     // finish previous primary track
152     FinishPrimaryTrack();
153     fPrimaryTrackID = aTrack->GetTrackID();
154     
155     // begin this primary track
156     gAlice->BeginPrimary();
157   }
158   else { 
159     // save secondary particles info 
160     SaveTrack(aTrack);
161   }
162   
163   // verbose
164   if (trackID == fNewVerboseTrackID) {
165       G4String command = "/tracking/verbose ";
166       AliGlobals::AppendNumberToString(command, fNewVerboseLevel);
167       G4UImanager::GetUIpointer()->ApplyCommand(command);
168   }    
169
170   // aliroot pre track actions
171   gAlice->PreTrack();
172 }
173
174 //_____________________________________________________________________________
175 void AliTrackingAction::PostTrackingAction(const G4Track* aTrack)
176 {
177 // Called by G4 kernel after finishing tracking.
178 // ---
179
180   fTrackCounter++;
181   
182   // set parent track particle index to all secondary tracks 
183   G4TrackVector* secondaryTracks 
184     = fpTrackingManager->GetSteppingManager()->GetSecondary();
185   if (secondaryTracks){
186     G4int i;
187     for (i=0; i<secondaryTracks->size(); i++) {
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
213   gAlice->PostTrack();
214 }
215
216 //_____________________________________________________________________________
217 void AliTrackingAction::FinishPrimaryTrack()
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
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;
242 }  
243
244 //_____________________________________________________________________________
245 void AliTrackingAction::SaveTrack(const G4Track* track)
246 {
247 // Get all needed parameters from G4track and pass them
248 // to AliRun::SetTrack() that creates corresponding TParticle
249 // in the AliRun::fParticles array.
250 // ----
251
252   // parent particle index 
253   G4int parentID = track->GetParentID();
254   G4int motherIndex;
255   if (parentID == 0) { 
256     motherIndex = -1; 
257   }
258   else {
259     motherIndex 
260       = GetTrackInformation(track,"SaveTrack")->GetParentParticleID();
261   }
262      
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
270   G4int pdg 
271     = TG4ParticlesManager::Instance()
272       ->GetPDGEncodingFast(track->GetDefinition());
273
274   // track kinematics  
275   G4ThreeVector momentum = track->GetMomentum(); 
276   
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
294   // production process
295   AliMCProcess mcProcess;  
296   const G4VProcess* kpProcess = track->GetCreatorProcess();
297   if (!kpProcess) {
298     mcProcess = kPPrimary;
299   }
300   else {  
301     mcProcess = TG4PhysicsManager::Instance()->GetMCProcess(kpProcess);  
302     // distinguish kPDeltaRay from kPEnergyLoss  
303     if (mcProcess == kPEnergyLoss) mcProcess = kPDeltaRay;
304   }  
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                    
311 }
312
313 //_____________________________________________________________________________
314 void 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
323 //_____________________________________________________________________________
324 void 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 }