]> git.uio.no Git - u/mrichter/AliRoot.git/blob - AliGeant4/AliTrackingAction.cxx
bab74ea1f3055387b6ab23b11a038817730370d5
[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 "AliSensitiveDetector.h"
9 #include "AliRun.h"
10 #include "AliGlobals.h"  
11 #include "TG4StepManager.h"
12
13 #include <G4TrackingManager.hh>
14 #include <G4Track.hh>
15 #include <G4Event.hh>
16 #include <G4SDManager.hh>
17 #include <G4VSensitiveDetector.hh>
18 #include <G4VHitsCollection.hh>
19
20 #include <TTree.h>
21 #include <TParticle.h>
22
23 // static data members
24 AliTrackingAction* AliTrackingAction::fgInstance = 0;
25
26 // static methods
27 AliTrackingAction* AliTrackingAction::Instance() {
28 // 
29   return fgInstance; 
30 }
31
32 AliTrackingAction::AliTrackingAction()
33   : fParticles(0),
34     fPrimaryTrackID(0),
35     fVerboseLevel(2),
36     fSavePrimaries(true),
37     fPrimariesCounter(0),
38     fParticlesCounter(0),
39     fLastParticleIndex(-1)
40 {
41 //
42   if (fgInstance) { 
43     AliGlobals::Exception("AliTrackingAction constructed twice."); 
44   }
45
46   fMessenger = new AliTrackingActionMessenger(this);
47   fgInstance = this;
48 }
49
50 AliTrackingAction::AliTrackingAction(const AliTrackingAction& right) {
51 //
52   AliGlobals::Exception("AliTrackingAction is protected from copying.");
53 }
54
55 AliTrackingAction::~AliTrackingAction() {
56 //
57   delete fMessenger;
58 }
59
60 // operators
61
62 AliTrackingAction& 
63 AliTrackingAction::operator=(const AliTrackingAction &right)
64 {
65   // check assignement to self
66   if (this == &right) return *this;
67   
68   AliGlobals::Exception("AliTrackingAction is protected from assigning.");
69
70   return *this;
71 }
72
73 // public methods
74
75 void AliTrackingAction::PrepareNewEvent()
76 {
77 // Called by G4 kernel at the beginning of event.
78 // ---
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     // improve this later with retrieving the generation process
127     // (primary particles are stored 
128     //  by AlStackingAction in ClassifyNewTrack() method)
129     G4String origin = "secondary"; 
130     SaveParticle(aTrack, origin);
131   }
132 }
133
134 void AliTrackingAction::PostTrackingAction(const G4Track* aTrack)
135 {
136 // Called by G4 kernel after finishing tracking.
137 // ---
138
139   G4String particleName 
140     = aTrack->GetDynamicParticle()->GetDefinition()->GetParticleName();
141   if (particleName == "opticalphoton") {
142     G4cout << "$$$ Track " <<  aTrack->GetTrackID()
143            << " is optical photon." << G4endl;
144   }                     
145 }
146
147 void AliTrackingAction::SaveAndDestroyTrack()
148 {
149 // Calls AliRun::PurifyKine and fills trees of hits
150 // after finishing tracking of each primary track.
151 // !! This method has to be also called from AlEventAction::EndOfEventAction() 
152 // for storing the last primary track of the current event.
153 // --- 
154
155   if (fPrimaryTrackID>0)
156   {
157      if (fVerboseLevel == 3) { 
158        G4cout << "$$$ Primary track " << fPrimaryTrackID << G4endl;
159      } 
160      else if ( fVerboseLevel == 2 &&  fPrimaryTrackID % 10 == 0 ) {
161          G4cout << "$$$ Primary track " << fPrimaryTrackID  << G4endl;
162      } 
163      else if ( fVerboseLevel == 1 &&  fPrimaryTrackID % 100 == 0 ) {
164          G4cout << "$$$ Primary track " << fPrimaryTrackID  << G4endl;
165      } 
166      
167      // aliroot
168      G4int lastSavedTrack 
169        = gAlice->PurifyKine(fLastParticleIndex, fParticlesCounter);
170      G4int nofPuredSavedTracks 
171        = gAlice->GetNtrack();
172      fLastParticleIndex = lastSavedTrack;
173      fParticlesCounter = nofPuredSavedTracks;
174
175      if(gAlice->TreeH()) gAlice->TreeH()->Fill();
176      gAlice->ResetHits();
177    }
178    else { 
179      fLastParticleIndex = fPrimariesCounter-1; 
180    }
181    fPrimaryTrackID = 0;
182 }  
183
184 void AliTrackingAction::SaveParticle(const G4Track* track, 
185                                      G4String processName)
186 {
187 // Converts G4track to TParticle and saves it in AliRun::fParticles
188 // array.
189 // ----
190
191   fParticlesCounter++;
192
193   // track history
194   G4int firstDaughter = -1; 
195   G4int lastDaughter = -1;      
196   G4int parentID = track->GetParentID();
197   G4int motherIndex1;
198   G4int motherIndex2 = -1;
199   if (parentID == 0) { 
200     motherIndex1 = -1; 
201     fPrimariesCounter++;    
202   }
203   else {
204     // set first/last child for already defined particles
205     motherIndex1 = GetParticleIndex(parentID);
206     // aliroot
207     TParticle* parent 
208       = dynamic_cast<TParticle*>(fParticles->UncheckedAt(motherIndex1));
209     if (parent) {  
210       if (parent->GetFirstDaughter()<0) 
211         parent->SetFirstDaughter(fParticlesCounter);
212       parent->SetLastDaughter(fParticlesCounter);
213     }
214     else {
215       AliGlobals::Exception(
216         "AliTrackingAction::SaveParticle: Unknown particle type");
217     }           
218   };
219      
220   // ?? status of particle 
221   // temporarily used for storing trackID from G4
222   G4int ks = track->GetTrackID();
223    
224   // particle type
225   // is this G3 ekvivalent ?? - check
226   G4int pdg = track->GetDefinition()->GetPDGEncoding();
227
228   // track kinematics
229   G4ThreeVector momentum = track->GetMomentum(); 
230   G4double px = momentum.x()/GeV;
231   G4double py = momentum.y()/GeV;
232   G4double pz = momentum.z()/GeV;
233   G4double e = track->GetTotalEnergy()/GeV;
234
235   G4ThreeVector position = track->GetPosition(); 
236   G4double vx = position.x()/cm;
237   G4double vy = position.y()/cm;
238   G4double vz = position.z()/cm;
239   // time of production - check if ekvivalent with G4
240   G4double t = track->GetGlobalTime();
241
242   G4ThreeVector polarization = track->GetPolarization(); 
243   G4double polX = polarization.x();
244   G4double polY = polarization.y();
245   G4double polZ = polarization.z();
246
247   // aliroot
248   TClonesArray& theCollectionRef = *fParticles;
249   G4int nofParticles = theCollectionRef.GetEntriesFast();
250   TParticle* particle 
251    = new(theCollectionRef[nofParticles]) 
252        TParticle(pdg, ks, motherIndex1, motherIndex1, 
253          firstDaughter, lastDaughter, px, py, pz, e, vx, vy, vz, t);
254   particle->SetPolarisation(polX, polY, polZ);
255   particle->SetBit(kKeepBit, false); 
256 }
257
258 G4int AliTrackingAction::GetParticleIndex(G4int trackID)
259 {
260 // Converts the trackID into the index of the particle
261 // in AliRun::fParticles array.
262 // ---
263
264   if (trackID <= fPrimariesCounter) { 
265     return trackID-1; 
266   }
267   else
268     for (G4int i=fLastParticleIndex+1; i<fParticlesCounter; i++) {
269       // aliroot
270       TParticle* particle
271         = (TParticle*)fParticles->UncheckedAt(i);
272       if (trackID == particle->GetStatusCode()) return i;
273     }
274
275   AliGlobals::Exception("AliTrackingAction::GetParticleIndex() has failed.");
276   return 0;   
277 }