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