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