]> git.uio.no Git - u/mrichter/AliRoot.git/blob - AliGeant4/AliTrackingAction.cxx
fixed date (in header) replaced with java script
[u/mrichter/AliRoot.git] / AliGeant4 / AliTrackingAction.cxx
1 // $Id$
2 // Category: event
3 //
4 // See the class description in the header file.
5
6 #include "AliTrackingAction.h"
7 #include "AliTrackingActionMessenger.h"
8 #include "AliRun.h"
9 #include "AliGlobals.h"  
10 #include "TG4StepManager.h"
11 #include "TG4PhysicsManager.h"
12
13 #include <G4TrackingManager.hh>
14 #include <G4Track.hh>
15
16 #include <TTree.h>
17 #include <TParticle.h>
18 #include <TClonesArray.h>
19
20 // static data members
21 AliTrackingAction* AliTrackingAction::fgInstance = 0;
22
23 // static methods
24 AliTrackingAction* AliTrackingAction::Instance() {
25 // 
26   return fgInstance; 
27 }
28
29 AliTrackingAction::AliTrackingAction()
30   : fParticles(0),
31     fPrimaryTrackID(0),
32     fVerboseLevel(2),
33     fSavePrimaries(true),
34     fPrimariesCounter(0),
35     fParticlesCounter(0),
36     fTrackCounter(0),
37     fLastParticleIndex(-1)
38 {
39 //
40   if (fgInstance) { 
41     AliGlobals::Exception("AliTrackingAction constructed twice."); 
42   }
43
44   fMessenger = new AliTrackingActionMessenger(this);
45   fgInstance = this;
46 }
47
48 AliTrackingAction::AliTrackingAction(const AliTrackingAction& right) {
49 //
50   AliGlobals::Exception("AliTrackingAction is protected from copying.");
51 }
52
53 AliTrackingAction::~AliTrackingAction() {
54 //
55   delete fMessenger;
56 }
57
58 // operators
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 // public methods
72
73 void AliTrackingAction::PrepareNewEvent()
74 {
75 // Called by G4 kernel at the beginning of event.
76 // ---
77
78   fTrackCounter = 0;
79
80   // aliroot
81   if (!fParticles) fParticles = gAlice->Particles();
82
83   // set the particles and primaries counter to the
84   // number of tracks already stored in AliRun::fParticles
85   // (fParticles can already contain primary particles
86   //  saved by AliGenerator)
87   G4int nofTracks = gAlice->GetNtrack();
88   fPrimariesCounter = nofTracks;
89   fParticlesCounter = nofTracks;
90   fLastParticleIndex = nofTracks-1;
91
92   // set g4 stepping manager pointer
93   G4SteppingManager* pG4StepManager 
94     = fpTrackingManager->GetSteppingManager();
95   TG4StepManager* pStepManager = TG4StepManager::Instance();
96   pStepManager->SetSteppingManager(pG4StepManager);
97 }
98
99 void AliTrackingAction::PreTrackingAction(const G4Track* aTrack)
100 {
101 // Called by G4 kernel before starting tracking.
102 // ---
103
104   // aliroot
105   // track index in the particles array
106   G4int trackID = aTrack->GetTrackID();
107   G4int parentID = aTrack->GetParentID();
108   Int_t trackIndex;
109   if (parentID==0) { 
110     trackIndex = trackID; 
111   } 
112   else { 
113     trackIndex = GetNofSavedTracks(); 
114   }
115
116   // in AliRoot (from V3.0) track numbering starts from 0
117   gAlice->SetCurrentTrack(--trackIndex);
118
119   if (parentID == 0) {  
120     // save primary track
121     SaveAndDestroyTrack();
122     fPrimaryTrackID = aTrack->GetTrackID();
123   }
124   else { 
125     // save secondary particles info 
126     SaveParticle(aTrack);
127   }
128   
129   gAlice->PreTrack();
130 }
131
132 void AliTrackingAction::PostTrackingAction(const G4Track* aTrack)
133 {
134 // Called by G4 kernel after finishing tracking.
135 // ---
136
137   fTrackCounter++;
138
139   gAlice->PostTrack();
140 }
141
142 void AliTrackingAction::SaveAndDestroyTrack()
143 {
144 // Calls AliRun::PurifyKine and fills trees of hits
145 // after finishing tracking of each primary track.
146 // !! This method has to be also called from AlEventAction::EndOfEventAction() 
147 // for storing the last primary track of the current event.
148 // --- 
149
150   if (fPrimaryTrackID>0)
151   {
152      if (fVerboseLevel == 3) { 
153        G4cout << "$$$ Primary track " << fPrimaryTrackID << G4endl;
154      } 
155      else if ( fVerboseLevel == 2 &&  fPrimaryTrackID % 10 == 0 ) {
156          G4cout << "$$$ Primary track " << fPrimaryTrackID  << G4endl;
157      } 
158      else if ( fVerboseLevel == 1 &&  fPrimaryTrackID % 100 == 0 ) {
159          G4cout << "$$$ Primary track " << fPrimaryTrackID  << G4endl;
160      } 
161      
162      // aliroot
163      G4int lastSavedTrack 
164        = gAlice->PurifyKine(fLastParticleIndex, fParticlesCounter);
165      G4int nofPuredSavedTracks 
166        = gAlice->GetNtrack();
167      fLastParticleIndex = lastSavedTrack;
168      fParticlesCounter = nofPuredSavedTracks;
169
170      if(gAlice->TreeH()) gAlice->TreeH()->Fill();
171      gAlice->ResetHits();
172    }
173    else { 
174      fLastParticleIndex = fPrimariesCounter-1; 
175    }
176    fPrimaryTrackID = 0;
177 }  
178
179 void AliTrackingAction::SaveParticle(const G4Track* track)
180 {
181 // Converts G4track to TParticle and saves it in AliRun::fParticles
182 // array.
183 // ----
184
185   fParticlesCounter++;
186
187   // track history
188   G4int firstDaughter = -1; 
189   G4int lastDaughter = -1;      
190   G4int parentID = track->GetParentID();
191   G4int motherIndex1;
192   G4int motherIndex2 = -1;
193   if (parentID == 0) { 
194     motherIndex1 = -1; 
195     fPrimariesCounter++;    
196   }
197   else {
198     // set first/last child for already defined particles
199     motherIndex1 = GetParticleIndex(parentID);
200     // aliroot
201     TParticle* parent 
202       = dynamic_cast<TParticle*>(fParticles->UncheckedAt(motherIndex1));
203     if (parent) {  
204       if (parent->GetFirstDaughter()<0) 
205         parent->SetFirstDaughter(fParticlesCounter);
206       parent->SetLastDaughter(fParticlesCounter);
207     }
208     else {
209       AliGlobals::Exception(
210         "AliTrackingAction::SaveParticle: Unknown particle type");
211     }           
212   };
213      
214   // ?? status of particle 
215   // temporarily used for storing trackID from G4
216   G4int ks = track->GetTrackID();
217    
218   // particle type
219   // is this G3 ekvivalent ?? - check
220   G4int pdg = track->GetDefinition()->GetPDGEncoding();
221
222   // track kinematics
223   G4ThreeVector momentum = track->GetMomentum(); 
224   G4double px = momentum.x()/GeV;
225   G4double py = momentum.y()/GeV;
226   G4double pz = momentum.z()/GeV;
227   G4double e = track->GetTotalEnergy()/GeV;
228
229   G4ThreeVector position = track->GetPosition(); 
230   G4double vx = position.x()/cm;
231   G4double vy = position.y()/cm;
232   G4double vz = position.z()/cm;
233   // time of production - check if ekvivalent with G4
234   G4double t = track->GetGlobalTime();
235
236   G4ThreeVector polarization = track->GetPolarization(); 
237   G4double polX = polarization.x();
238   G4double polY = polarization.y();
239   G4double polZ = polarization.z();
240
241   // create TParticle
242   TClonesArray& theCollectionRef = *fParticles;
243   G4int nofParticles = theCollectionRef.GetEntriesFast();
244   TParticle* particle 
245    = new(theCollectionRef[nofParticles]) 
246        TParticle(pdg, ks, motherIndex1, motherIndex1, 
247          firstDaughter, lastDaughter, px, py, pz, e, vx, vy, vz, t);
248   particle->SetPolarisation(polX, polY, polZ);
249   particle->SetBit(kKeepBit, false); 
250   
251   // set production process
252   AliMCProcess mcProcess;  
253   const G4VProcess* kpProcess = track->GetCreatorProcess();
254   if (!kpProcess) {
255     mcProcess = kPPrimary;
256   }
257   else {  
258     TG4PhysicsManager* pPhysicsManager = TG4PhysicsManager::Instance();
259     mcProcess = pPhysicsManager->GetMCProcess(kpProcess);  
260     // distinguish kPDeltaRay from kPEnergyLoss  
261     if (mcProcess == kPEnergyLoss) mcProcess = kPDeltaRay;
262   }  
263   particle->SetUniqueID(mcProcess);  
264 }
265
266 G4int AliTrackingAction::GetParticleIndex(G4int trackID)
267 {
268 // Converts the trackID into the index of the particle
269 // in AliRun::fParticles array.
270 // ---
271
272   if (trackID <= fPrimariesCounter) { 
273     return trackID-1; 
274   }
275   else
276     for (G4int i=fLastParticleIndex+1; i<fParticlesCounter; i++) {
277       // aliroot
278       TParticle* particle
279         = (TParticle*)fParticles->UncheckedAt(i);
280       if (trackID == particle->GetStatusCode()) return i;
281     }
282
283   AliGlobals::Exception("AliTrackingAction::GetParticleIndex() has failed.");
284   return 0;   
285 }