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