]>
Commit | Line | Data |
---|---|---|
676fb573 | 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> | |
676fb573 | 14 | |
71c03be6 | 15 | #include <TTree.h> |
676fb573 | 16 | #include <TParticle.h> |
c97337f9 | 17 | #include <TClonesArray.h> |
676fb573 | 18 | |
19 | // static data members | |
20 | AliTrackingAction* AliTrackingAction::fgInstance = 0; | |
21 | ||
22 | // static methods | |
23 | AliTrackingAction* AliTrackingAction::Instance() { | |
24 | // | |
25 | return fgInstance; | |
26 | } | |
27 | ||
28 | AliTrackingAction::AliTrackingAction() | |
29 | : fParticles(0), | |
30 | fPrimaryTrackID(0), | |
31 | fVerboseLevel(2), | |
32 | fSavePrimaries(true), | |
33 | fPrimariesCounter(0), | |
34 | fParticlesCounter(0), | |
3bfabcfc | 35 | fTrackCounter(0), |
676fb573 | 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 | ||
3bfabcfc | 77 | fTrackCounter = 0; |
78 | ||
676fb573 | 79 | // aliroot |
80 | if (!fParticles) fParticles = gAlice->Particles(); | |
81 | ||
82 | // set the particles and primaries counter to the | |
83 | // number of tracks already stored in AliRun::fParticles | |
84 | // (fParticles can already contain primary particles | |
85 | // saved by AliGenerator) | |
86 | G4int nofTracks = gAlice->GetNtrack(); | |
87 | fPrimariesCounter = nofTracks; | |
88 | fParticlesCounter = nofTracks; | |
89 | fLastParticleIndex = nofTracks-1; | |
90 | ||
91 | // set g4 stepping manager pointer | |
92 | G4SteppingManager* pG4StepManager | |
93 | = fpTrackingManager->GetSteppingManager(); | |
94 | TG4StepManager* pStepManager = TG4StepManager::Instance(); | |
95 | pStepManager->SetSteppingManager(pG4StepManager); | |
96 | } | |
97 | ||
9bcb6317 | 98 | void AliTrackingAction::PreTrackingAction(const G4Track* aTrack) |
676fb573 | 99 | { |
100 | // Called by G4 kernel before starting tracking. | |
101 | // --- | |
102 | ||
103 | // aliroot | |
104 | // track index in the particles array | |
105 | G4int trackID = aTrack->GetTrackID(); | |
106 | G4int parentID = aTrack->GetParentID(); | |
107 | Int_t trackIndex; | |
108 | if (parentID==0) { | |
109 | trackIndex = trackID; | |
110 | } | |
111 | else { | |
112 | trackIndex = GetNofSavedTracks(); | |
113 | } | |
114 | ||
115 | // in AliRoot (from V3.0) track numbering starts from 0 | |
116 | gAlice->SetCurrentTrack(--trackIndex); | |
117 | ||
118 | if (parentID == 0) { | |
119 | // save primary track | |
120 | SaveAndDestroyTrack(); | |
121 | fPrimaryTrackID = aTrack->GetTrackID(); | |
122 | } | |
123 | else { | |
124 | // save secondary particles info | |
125 | // improve this later with retrieving the generation process | |
126 | // (primary particles are stored | |
127 | // by AlStackingAction in ClassifyNewTrack() method) | |
128 | G4String origin = "secondary"; | |
129 | SaveParticle(aTrack, origin); | |
9bcb6317 | 130 | } |
676fb573 | 131 | } |
132 | ||
9bcb6317 | 133 | void AliTrackingAction::PostTrackingAction(const G4Track* aTrack) |
676fb573 | 134 | { |
135 | // Called by G4 kernel after finishing tracking. | |
136 | // --- | |
137 | ||
3bfabcfc | 138 | fTrackCounter++; |
676fb573 | 139 | } |
140 | ||
141 | void AliTrackingAction::SaveAndDestroyTrack() | |
142 | { | |
143 | // Calls AliRun::PurifyKine and fills trees of hits | |
144 | // after finishing tracking of each primary track. | |
145 | // !! This method has to be also called from AlEventAction::EndOfEventAction() | |
146 | // for storing the last primary track of the current event. | |
147 | // --- | |
148 | ||
149 | if (fPrimaryTrackID>0) | |
150 | { | |
151 | if (fVerboseLevel == 3) { | |
5f1d09c5 | 152 | G4cout << "$$$ Primary track " << fPrimaryTrackID << G4endl; |
676fb573 | 153 | } |
154 | else if ( fVerboseLevel == 2 && fPrimaryTrackID % 10 == 0 ) { | |
5f1d09c5 | 155 | G4cout << "$$$ Primary track " << fPrimaryTrackID << G4endl; |
676fb573 | 156 | } |
157 | else if ( fVerboseLevel == 1 && fPrimaryTrackID % 100 == 0 ) { | |
5f1d09c5 | 158 | G4cout << "$$$ Primary track " << fPrimaryTrackID << G4endl; |
676fb573 | 159 | } |
160 | ||
161 | // aliroot | |
162 | G4int lastSavedTrack | |
163 | = gAlice->PurifyKine(fLastParticleIndex, fParticlesCounter); | |
164 | G4int nofPuredSavedTracks | |
165 | = gAlice->GetNtrack(); | |
166 | fLastParticleIndex = lastSavedTrack; | |
167 | fParticlesCounter = nofPuredSavedTracks; | |
168 | ||
169 | if(gAlice->TreeH()) gAlice->TreeH()->Fill(); | |
170 | gAlice->ResetHits(); | |
171 | } | |
172 | else { | |
173 | fLastParticleIndex = fPrimariesCounter-1; | |
174 | } | |
175 | fPrimaryTrackID = 0; | |
176 | } | |
177 | ||
178 | void AliTrackingAction::SaveParticle(const G4Track* track, | |
179 | G4String processName) | |
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 | // aliroot | |
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); | |
b65f8990 | 249 | particle->SetBit(kKeepBit, false); |
676fb573 | 250 | } |
251 | ||
252 | G4int AliTrackingAction::GetParticleIndex(G4int trackID) | |
253 | { | |
254 | // Converts the trackID into the index of the particle | |
255 | // in AliRun::fParticles array. | |
256 | // --- | |
257 | ||
258 | if (trackID <= fPrimariesCounter) { | |
259 | return trackID-1; | |
260 | } | |
261 | else | |
262 | for (G4int i=fLastParticleIndex+1; i<fParticlesCounter; i++) { | |
263 | // aliroot | |
264 | TParticle* particle | |
265 | = (TParticle*)fParticles->UncheckedAt(i); | |
266 | if (trackID == particle->GetStatusCode()) return i; | |
267 | } | |
268 | ||
269 | AliGlobals::Exception("AliTrackingAction::GetParticleIndex() has failed."); | |
270 | return 0; | |
271 | } |