]>
Commit | Line | Data |
---|---|---|
9e1a0ddb | 1 | /************************************************************************** |
2 | * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. * | |
3 | * * | |
4 | * Author: The ALICE Off-line Project. * | |
5 | * Contributors are mentioned in the code where appropriate. * | |
6 | * * | |
7 | * Permission to use, copy, modify and distribute this software and its * | |
8 | * documentation strictly for non-commercial purposes is hereby granted * | |
9 | * without fee, provided that the above copyright notice appears in all * | |
10 | * copies and that both the copyright notice and this permission notice * | |
11 | * appear in the supporting documentation. The authors make no claims * | |
12 | * about the suitability of this software for any purpose. It is * | |
13 | * provided "as is" without express or implied warranty. * | |
14 | **************************************************************************/ | |
15 | ||
acd84897 | 16 | /* $Id$ */ |
9e1a0ddb | 17 | |
18 | /////////////////////////////////////////////////////////////////////////////// | |
19 | // // | |
af7ba10c | 20 | // Particles stack class // |
21 | // Implements the TMCVirtualStack of the Virtual Monte Carlo // | |
22 | // Holds the particles transported during simulation // | |
23 | // Is used to compare results of reconstruction with simulation // | |
24 | // Author A.Morsch // | |
9e1a0ddb | 25 | // // |
26 | /////////////////////////////////////////////////////////////////////////////// | |
27 | ||
9e1a0ddb | 28 | |
7ca4655f | 29 | #include <TClonesArray.h> |
9e1a0ddb | 30 | #include <TObjArray.h> |
7ca4655f | 31 | #include <TPDGCode.h> |
9e1a0ddb | 32 | #include <TParticle.h> |
fe046ade | 33 | #include <TParticlePDG.h> |
9e1a0ddb | 34 | #include <TTree.h> |
942a9039 | 35 | #include <TDirectory.h> |
9e1a0ddb | 36 | |
21bf7095 | 37 | #include "AliLog.h" |
88cb7938 | 38 | #include "AliStack.h" |
9e1a0ddb | 39 | |
40 | ClassImp(AliStack) | |
41 | ||
e2afb3b6 | 42 | //_______________________________________________________________________ |
43 | AliStack::AliStack(): | |
304b9ce3 | 44 | fParticles("TParticle", 1000), |
d94af0c1 | 45 | fParticleMap(), |
e2afb3b6 | 46 | fParticleFileMap(0), |
47 | fParticleBuffer(0), | |
90e48c0c | 48 | fCurrentTrack(0), |
e2afb3b6 | 49 | fTreeK(0), |
50 | fNtrack(0), | |
51 | fNprimary(0), | |
52 | fCurrent(-1), | |
53 | fCurrentPrimary(-1), | |
54 | fHgwmk(0), | |
88cb7938 | 55 | fLoadPoint(0), |
942a9039 | 56 | fTrackLabelMap(0) |
9e1a0ddb | 57 | { |
58 | // | |
e2afb3b6 | 59 | // Default constructor |
9e1a0ddb | 60 | // |
9e1a0ddb | 61 | } |
62 | ||
e2afb3b6 | 63 | //_______________________________________________________________________ |
942a9039 | 64 | AliStack::AliStack(Int_t size, const char* /*evfoldname*/): |
d94af0c1 | 65 | fParticles("TParticle",1000), |
66 | fParticleMap(size), | |
e2afb3b6 | 67 | fParticleFileMap(0), |
68 | fParticleBuffer(0), | |
90e48c0c | 69 | fCurrentTrack(0), |
e2afb3b6 | 70 | fTreeK(0), |
71 | fNtrack(0), | |
72 | fNprimary(0), | |
73 | fCurrent(-1), | |
74 | fCurrentPrimary(-1), | |
75 | fHgwmk(0), | |
88cb7938 | 76 | fLoadPoint(0), |
942a9039 | 77 | fTrackLabelMap(0) |
e2afb3b6 | 78 | { |
79 | // | |
80 | // Constructor | |
81 | // | |
82 | } | |
9e1a0ddb | 83 | |
e2afb3b6 | 84 | //_______________________________________________________________________ |
85 | AliStack::AliStack(const AliStack& st): | |
9938d01c | 86 | TVirtualMCStack(st), |
d94af0c1 | 87 | fParticles("TParticle",1000), |
88 | fParticleMap(*(st.Particles())), | |
9938d01c | 89 | fParticleFileMap(st.fParticleFileMap), |
90 | fParticleBuffer(0), | |
91 | fCurrentTrack(0), | |
92 | fTreeK((TTree*)(st.fTreeK->Clone())), | |
93 | fNtrack(st.GetNtrack()), | |
94 | fNprimary(st.GetNprimary()), | |
95 | fCurrent(-1), | |
96 | fCurrentPrimary(-1), | |
97 | fHgwmk(0), | |
98 | fLoadPoint(0), | |
942a9039 | 99 | fTrackLabelMap(0) |
9e1a0ddb | 100 | { |
9938d01c | 101 | // Copy constructor |
9e1a0ddb | 102 | } |
103 | ||
9938d01c | 104 | |
e2afb3b6 | 105 | //_______________________________________________________________________ |
6c4904c2 | 106 | void AliStack::Copy(TObject&) const |
e2afb3b6 | 107 | { |
21bf7095 | 108 | AliFatal("Not implemented!"); |
e2afb3b6 | 109 | } |
9e1a0ddb | 110 | |
e2afb3b6 | 111 | //_______________________________________________________________________ |
9e1a0ddb | 112 | AliStack::~AliStack() |
113 | { | |
114 | // | |
115 | // Destructor | |
116 | // | |
117 | ||
d94af0c1 | 118 | fParticles.Clear(); |
9e1a0ddb | 119 | } |
120 | ||
121 | // | |
122 | // public methods | |
123 | // | |
124 | ||
125 | //_____________________________________________________________________________ | |
642f15cf | 126 | void AliStack::PushTrack(Int_t done, Int_t parent, Int_t pdg, Float_t *pmom, |
88cb7938 | 127 | Float_t *vpos, Float_t *polar, Float_t tof, |
128 | TMCProcess mech, Int_t &ntr, Float_t weight, Int_t is) | |
9e1a0ddb | 129 | { |
130 | // | |
131 | // Load a track on the stack | |
132 | // | |
b4bb4536 | 133 | // done 1 if the track has to be transported |
134 | // 0 if not | |
9e1a0ddb | 135 | // parent identifier of the parent track. -1 for a primary |
136 | // pdg particle code | |
137 | // pmom momentum GeV/c | |
138 | // vpos position | |
139 | // polar polarisation | |
140 | // tof time of flight in seconds | |
141 | // mecha production mechanism | |
142 | // ntr on output the number of the track stored | |
143 | // | |
144 | ||
9e1a0ddb | 145 | // const Float_t tlife=0; |
146 | ||
147 | // | |
148 | // Here we get the static mass | |
149 | // For MC is ok, but a more sophisticated method could be necessary | |
150 | // if the calculated mass is required | |
151 | // also, this method is potentially dangerous if the mass | |
152 | // used in the MC is not the same of the PDG database | |
153 | // | |
fe046ade | 154 | TParticlePDG* pmc = TDatabasePDG::Instance()->GetParticle(pdg); |
155 | if (pmc) { | |
156 | Float_t mass = TDatabasePDG::Instance()->GetParticle(pdg)->Mass(); | |
157 | Float_t e=TMath::Sqrt(mass*mass+pmom[0]*pmom[0]+ | |
158 | pmom[1]*pmom[1]+pmom[2]*pmom[2]); | |
159 | ||
ccf7a81f | 160 | // printf("Loading mass %f ene %f No %d ip %d parent %d done %d pos %f %f %f mom %f %f %f kS %d m \n", |
161 | // mass,e,fNtrack,pdg,parent,done,vpos[0],vpos[1],vpos[2],pmom[0],pmom[1],pmom[2],kS); | |
162 | ||
9e1a0ddb | 163 | |
642f15cf | 164 | PushTrack(done, parent, pdg, pmom[0], pmom[1], pmom[2], e, |
fe046ade | 165 | vpos[0], vpos[1], vpos[2], tof, polar[0], polar[1], polar[2], |
166 | mech, ntr, weight, is); | |
167 | } else { | |
21bf7095 | 168 | AliWarning(Form("Particle type %d not defined in PDG Database !", pdg)); |
169 | AliWarning("Particle skipped !"); | |
fe046ade | 170 | } |
9e1a0ddb | 171 | } |
172 | ||
173 | //_____________________________________________________________________________ | |
642f15cf | 174 | void AliStack::PushTrack(Int_t done, Int_t parent, Int_t pdg, |
9e1a0ddb | 175 | Double_t px, Double_t py, Double_t pz, Double_t e, |
176 | Double_t vx, Double_t vy, Double_t vz, Double_t tof, | |
177 | Double_t polx, Double_t poly, Double_t polz, | |
98490ea9 | 178 | TMCProcess mech, Int_t &ntr, Double_t weight, Int_t is) |
9e1a0ddb | 179 | { |
180 | // | |
181 | // Load a track on the stack | |
182 | // | |
b4bb4536 | 183 | // done 1 if the track has to be transported |
184 | // 0 if not | |
9e1a0ddb | 185 | // parent identifier of the parent track. -1 for a primary |
186 | // pdg particle code | |
187 | // kS generation status code | |
188 | // px, py, pz momentum GeV/c | |
189 | // vx, vy, vz position | |
190 | // polar polarisation | |
191 | // tof time of flight in seconds | |
192 | // mech production mechanism | |
193 | // ntr on output the number of the track stored | |
194 | // | |
195 | // New method interface: | |
196 | // arguments were changed to be in correspondence with TParticle | |
197 | // constructor. | |
198 | // Note: the energy is not calculated from the static mass but | |
199 | // it is passed by argument e. | |
200 | ||
9e1a0ddb | 201 | const Int_t kFirstDaughter=-1; |
202 | const Int_t kLastDaughter=-1; | |
6308b375 | 203 | |
204 | ||
9e1a0ddb | 205 | TParticle* particle |
d94af0c1 | 206 | = new(fParticles[fLoadPoint++]) |
47c8bcbe | 207 | TParticle(pdg, is, parent, -1, kFirstDaughter, kLastDaughter, |
9e1a0ddb | 208 | px, py, pz, e, vx, vy, vz, tof); |
209 | ||
210 | particle->SetPolarisation(polx, poly, polz); | |
211 | particle->SetWeight(weight); | |
212 | particle->SetUniqueID(mech); | |
213 | ||
f4414561 | 214 | |
215 | ||
216 | if(!done) { | |
217 | particle->SetBit(kDoneBit); | |
218 | } else { | |
219 | particle->SetBit(kTransportBit); | |
220 | } | |
221 | ||
222 | ||
9e1a0ddb | 223 | |
224 | // Declare that the daughter information is valid | |
225 | particle->SetBit(kDaughtersBit); | |
226 | // Add the particle to the stack | |
6308b375 | 227 | |
228 | ||
d94af0c1 | 229 | fParticleMap.AddAtAndExpand(particle, fNtrack);//CHECK!! |
9e1a0ddb | 230 | |
231 | if(parent>=0) { | |
d94af0c1 | 232 | particle = dynamic_cast<TParticle*>(fParticleMap.At(parent)); |
942a9039 | 233 | if (particle) { |
234 | particle->SetLastDaughter(fNtrack); | |
235 | if(particle->GetFirstDaughter()<0) particle->SetFirstDaughter(fNtrack); | |
236 | } | |
237 | else { | |
238 | AliError(Form("Parent %d does not exist",parent)); | |
239 | } | |
240 | } else { | |
241 | // | |
242 | // This is a primary track. Set high water mark for this event | |
243 | fHgwmk = fNtrack; | |
244 | // | |
245 | // Set also number if primary tracks | |
246 | fNprimary = fHgwmk+1; | |
247 | fCurrentPrimary++; | |
9e1a0ddb | 248 | } |
249 | ntr = fNtrack++; | |
250 | } | |
251 | ||
252 | //_____________________________________________________________________________ | |
642f15cf | 253 | TParticle* AliStack::PopNextTrack(Int_t& itrack) |
b9d0a01d | 254 | { |
255 | // | |
256 | // Returns next track from stack of particles | |
257 | // | |
258 | ||
259 | ||
260 | TParticle* track = GetNextParticle(); | |
261 | ||
262 | if (track) { | |
263 | itrack = fCurrent; | |
264 | track->SetBit(kDoneBit); | |
265 | } | |
fe046ade | 266 | else |
b9d0a01d | 267 | itrack = -1; |
fe046ade | 268 | |
269 | fCurrentTrack = track; | |
b9d0a01d | 270 | return track; |
271 | } | |
272 | ||
b9d0a01d | 273 | //_____________________________________________________________________________ |
642f15cf | 274 | TParticle* AliStack::PopPrimaryForTracking(Int_t i) |
b9d0a01d | 275 | { |
276 | // | |
277 | // Returns i-th primary particle if it is flagged to be tracked, | |
278 | // 0 otherwise | |
279 | // | |
280 | ||
281 | TParticle* particle = Particle(i); | |
282 | ||
283 | if (!particle->TestBit(kDoneBit)) | |
284 | return particle; | |
285 | else | |
286 | return 0; | |
287 | } | |
288 | ||
9e1a0ddb | 289 | //_____________________________________________________________________________ |
77db8cbf | 290 | Bool_t AliStack::PurifyKine() |
9e1a0ddb | 291 | { |
292 | // | |
293 | // Compress kinematic tree keeping only flagged particles | |
294 | // and renaming the particle id's in all the hits | |
295 | // | |
296 | ||
77db8cbf | 297 | int nkeep = fHgwmk + 1, parent, i; |
9e1a0ddb | 298 | TParticle *part, *father; |
d94af0c1 | 299 | fTrackLabelMap.Set(fParticleMap.GetLast()+1); |
9e1a0ddb | 300 | |
301 | // Save in Header total number of tracks before compression | |
9e1a0ddb | 302 | // If no tracks generated return now |
77db8cbf | 303 | if(fHgwmk+1 == fNtrack) return kFALSE; |
9e1a0ddb | 304 | |
9e1a0ddb | 305 | // First pass, invalid Daughter information |
942a9039 | 306 | for(i=0; i<fNtrack; i++) { |
307 | // Preset map, to be removed later | |
308 | if(i<=fHgwmk) fTrackLabelMap[i]=i ; | |
309 | else { | |
310 | fTrackLabelMap[i] = -99; | |
d94af0c1 | 311 | if((part=dynamic_cast<TParticle*>(fParticleMap.At(i)))) { |
e44c5340 | 312 | // |
313 | // Check of this track should be kept for physics reasons | |
942a9039 | 314 | if (KeepPhysics(part)) KeepTrack(i); |
e44c5340 | 315 | // |
942a9039 | 316 | part->ResetBit(kDaughtersBit); |
317 | part->SetFirstDaughter(-1); | |
318 | part->SetLastDaughter(-1); | |
319 | } | |
76969085 | 320 | } |
9e1a0ddb | 321 | } |
322 | // Invalid daughter information for the parent of the first particle | |
323 | // generated. This may or may not be the current primary according to | |
324 | // whether decays have been recorded among the primaries | |
d94af0c1 | 325 | part = dynamic_cast<TParticle*>(fParticleMap.At(fHgwmk+1)); |
326 | fParticleMap.At(part->GetFirstMother())->ResetBit(kDaughtersBit); | |
9e1a0ddb | 327 | // Second pass, build map between old and new numbering |
328 | for(i=fHgwmk+1; i<fNtrack; i++) { | |
d94af0c1 | 329 | if(fParticleMap.At(i)->TestBit(kKeepBit)) { |
942a9039 | 330 | // This particle has to be kept |
331 | fTrackLabelMap[i]=nkeep; | |
332 | // If old and new are different, have to move the pointer | |
d94af0c1 | 333 | if(i!=nkeep) fParticleMap[nkeep]=fParticleMap.At(i); |
334 | part = dynamic_cast<TParticle*>(fParticleMap.At(nkeep)); | |
942a9039 | 335 | // as the parent is always *before*, it must be already |
336 | // in place. This is what we are checking anyway! | |
337 | if((parent=part->GetFirstMother())>fHgwmk) | |
338 | if(fTrackLabelMap[parent]==-99) Fatal("PurifyKine","fTrackLabelMap[%d] = -99!\n",parent); | |
339 | else part->SetFirstMother(fTrackLabelMap[parent]); | |
942a9039 | 340 | nkeep++; |
341 | } | |
9e1a0ddb | 342 | } |
343 | ||
344 | // Fix daughters information | |
345 | for (i=fHgwmk+1; i<nkeep; i++) { | |
d94af0c1 | 346 | part = dynamic_cast<TParticle*>(fParticleMap.At(i)); |
942a9039 | 347 | parent = part->GetFirstMother(); |
348 | if(parent>=0) { | |
d94af0c1 | 349 | father = dynamic_cast<TParticle*>(fParticleMap.At(parent)); |
942a9039 | 350 | if(father->TestBit(kDaughtersBit)) { |
351 | ||
352 | if(i<father->GetFirstDaughter()) father->SetFirstDaughter(i); | |
353 | if(i>father->GetLastDaughter()) father->SetLastDaughter(i); | |
354 | } else { | |
355 | // Initialise daughters info for first pass | |
356 | father->SetFirstDaughter(i); | |
357 | father->SetLastDaughter(i); | |
358 | father->SetBit(kDaughtersBit); | |
359 | } | |
9e1a0ddb | 360 | } |
9e1a0ddb | 361 | } |
942a9039 | 362 | // |
363 | // Now the output bit, from fHgwmk to nkeep we write everything and we erase | |
364 | if(nkeep > fParticleFileMap.GetSize()) fParticleFileMap.Set(Int_t (nkeep*1.5)); | |
365 | for (i=fHgwmk+1; i<nkeep; ++i) { | |
d94af0c1 | 366 | fParticleBuffer = dynamic_cast<TParticle*>(fParticleMap.At(i)); |
942a9039 | 367 | fParticleFileMap[i]=static_cast<Int_t>(TreeK()->GetEntries()); |
368 | TreeK()->Fill(); | |
d94af0c1 | 369 | fParticleMap[i]=fParticleBuffer=0; |
942a9039 | 370 | } |
9e1a0ddb | 371 | |
d94af0c1 | 372 | for (i = nkeep; i < fNtrack; ++i) fParticleMap[i]=0; |
6308b375 | 373 | |
942a9039 | 374 | Int_t toshrink = fNtrack-fHgwmk-1; |
375 | fLoadPoint-=toshrink; | |
2b22f272 | 376 | |
d94af0c1 | 377 | for(i=fLoadPoint; i<fLoadPoint+toshrink; ++i) fParticles.RemoveAt(i); |
942a9039 | 378 | fNtrack=nkeep; |
379 | fHgwmk=nkeep-1; | |
77db8cbf | 380 | return kTRUE; |
6308b375 | 381 | } |
382 | ||
1166db55 | 383 | |
77db8cbf | 384 | Bool_t AliStack::ReorderKine() |
6308b375 | 385 | { |
386 | // | |
387 | // In some transport code children might not come in a continuous sequence. | |
388 | // In this case the stack has to be reordered in order to establish the | |
389 | // mother daughter relation using index ranges. | |
390 | // | |
77db8cbf | 391 | if(fHgwmk+1 == fNtrack) return kFALSE; |
6308b375 | 392 | |
393 | // | |
394 | // Howmany secondaries have been produced ? | |
395 | Int_t nNew = fNtrack - fHgwmk - 1; | |
a1d92a3f | 396 | |
6308b375 | 397 | if (nNew > 0) { |
398 | Int_t i, j; | |
6308b375 | 399 | TArrayI map1(nNew); |
400 | // | |
401 | // Copy pointers to temporary array | |
402 | TParticle** tmp = new TParticle*[nNew]; | |
403 | ||
404 | for (i = 0; i < nNew; i++) { | |
d94af0c1 | 405 | if (fParticleMap.At(fHgwmk + 1 + i)) { |
406 | tmp[i] = (TParticle*) (fParticleMap.At(fHgwmk + 1 + i)); | |
6308b375 | 407 | } else { |
408 | tmp[i] = 0x0; | |
409 | } | |
410 | map1[i] = -99; | |
411 | } | |
412 | ||
413 | ||
414 | // | |
415 | // Reset LoadPoint | |
416 | // | |
1166db55 | 417 | Int_t loadPoint = fHgwmk + 1; |
6308b375 | 418 | // |
419 | // Re-Push particles into stack | |
420 | // The outer loop is over parents, the inner over children. | |
421 | // -1 refers to the primary particle | |
422 | // | |
a1d92a3f | 423 | for (i = -1; i < nNew-1; i++) { |
6308b375 | 424 | Int_t ipa; |
425 | TParticle* parP; | |
426 | if (i == -1) { | |
427 | ipa = tmp[0]->GetFirstMother(); | |
d94af0c1 | 428 | parP =dynamic_cast<TParticle*>(fParticleMap.At(ipa)); |
6308b375 | 429 | } else { |
430 | ipa = (fHgwmk + 1 + i); | |
431 | // Skip deleted particles | |
432 | if (!tmp[i]) continue; | |
433 | // Skip particles without children | |
434 | if (tmp[i]->GetFirstDaughter() == -1) continue; | |
435 | parP = tmp[i]; | |
436 | } | |
437 | // Reset daughter information | |
a1d92a3f | 438 | |
439 | Int_t idaumin = parP->GetFirstDaughter() - fHgwmk - 1; | |
440 | Int_t idaumax = parP->GetLastDaughter() - fHgwmk - 1; | |
6308b375 | 441 | parP->SetFirstDaughter(-1); |
442 | parP->SetLastDaughter(-1); | |
a1d92a3f | 443 | for (j = idaumin; j <= idaumax; j++) { |
6308b375 | 444 | // Skip deleted particles |
445 | if (!tmp[j]) continue; | |
446 | // Skip particles already handled | |
447 | if (map1[j] != -99) continue; | |
448 | Int_t jpa = tmp[j]->GetFirstMother(); | |
449 | // Check if daughter of current parent | |
450 | if (jpa == ipa) { | |
d94af0c1 | 451 | fParticleMap[loadPoint] = tmp[j]; |
6308b375 | 452 | // Re-establish daughter information |
1166db55 | 453 | parP->SetLastDaughter(loadPoint); |
454 | if (parP->GetFirstDaughter() == -1) parP->SetFirstDaughter(loadPoint); | |
6308b375 | 455 | // Set Mother information |
456 | if (i != -1) { | |
457 | tmp[j]->SetFirstMother(map1[i]); | |
458 | } | |
459 | // Build the map | |
1166db55 | 460 | map1[j] = loadPoint; |
6308b375 | 461 | // Increase load point |
1166db55 | 462 | loadPoint++; |
6308b375 | 463 | } |
464 | } // children | |
465 | } // parents | |
466 | ||
467 | delete[] tmp; | |
468 | ||
469 | // | |
470 | // Build map for remapping of hits | |
471 | // | |
942a9039 | 472 | fTrackLabelMap.Set(fNtrack); |
6308b375 | 473 | for (i = 0; i < fNtrack; i ++) { |
474 | if (i <= fHgwmk) { | |
942a9039 | 475 | fTrackLabelMap[i] = i; |
6308b375 | 476 | } else{ |
942a9039 | 477 | fTrackLabelMap[i] = map1[i - fHgwmk -1]; |
6308b375 | 478 | } |
479 | } | |
6308b375 | 480 | } // new particles poduced |
77db8cbf | 481 | |
482 | return kTRUE; | |
9e1a0ddb | 483 | } |
484 | ||
e44c5340 | 485 | Bool_t AliStack::KeepPhysics(TParticle* part) |
486 | { | |
487 | // | |
488 | // Some particles have to kept on the stack for reasons motivated | |
489 | // by physics analysis. Decision is put here. | |
490 | // | |
491 | Bool_t keep = kFALSE; | |
492 | // | |
493 | // Keep first-generation daughter from primaries with heavy flavor | |
494 | // | |
495 | Int_t parent = part->GetFirstMother(); | |
496 | if (parent >= 0 && parent <= fHgwmk) { | |
e2afb3b6 | 497 | TParticle* father = dynamic_cast<TParticle*>(Particles()->At(parent)); |
e44c5340 | 498 | Int_t kf = father->GetPdgCode(); |
499 | kf = TMath::Abs(kf); | |
500 | Int_t kfl = kf; | |
501 | // meson ? | |
502 | if (kfl > 10) kfl/=100; | |
503 | // baryon | |
504 | if (kfl > 10) kfl/=10; | |
505 | if (kfl > 10) kfl/=10; | |
506 | if (kfl >= 4) { | |
507 | keep = kTRUE; | |
508 | } | |
509 | } | |
510 | return keep; | |
511 | } | |
512 | ||
9e1a0ddb | 513 | //_____________________________________________________________________________ |
514 | void AliStack::FinishEvent() | |
515 | { | |
88cb7938 | 516 | // |
517 | // Write out the kinematics that was not yet filled | |
518 | // | |
9e1a0ddb | 519 | |
88cb7938 | 520 | // Update event header |
9e1a0ddb | 521 | |
88cb7938 | 522 | if (!TreeK()) { |
9e1a0ddb | 523 | // Fatal("FinishEvent", "No kinematics tree is defined."); |
524 | // Don't panic this is a probably a lego run | |
525 | return; | |
9e1a0ddb | 526 | } |
527 | ||
528 | CleanParents(); | |
942a9039 | 529 | if(TreeK()->GetEntries() ==0) { |
9e1a0ddb | 530 | // set the fParticleFileMap size for the first time |
531 | fParticleFileMap.Set(fHgwmk+1); | |
532 | } | |
533 | ||
534 | Bool_t allFilled = kFALSE; | |
535 | TObject *part; | |
536 | for(Int_t i=0; i<fHgwmk+1; ++i) | |
d94af0c1 | 537 | if((part=fParticleMap.At(i))) { |
e2afb3b6 | 538 | fParticleBuffer = dynamic_cast<TParticle*>(part); |
88cb7938 | 539 | fParticleFileMap[i]= static_cast<Int_t>(TreeK()->GetEntries()); |
540 | TreeK()->Fill(); | |
da2af24c | 541 | fParticleBuffer=0; |
d94af0c1 | 542 | fParticleMap.AddAt(0,i); |
9e1a0ddb | 543 | |
544 | // When all primaries were filled no particle!=0 | |
545 | // should be left => to be removed later. | |
21bf7095 | 546 | if (allFilled) AliWarning(Form("Why != 0 part # %d?\n",i)); |
9e1a0ddb | 547 | } |
88cb7938 | 548 | else |
549 | { | |
9e1a0ddb | 550 | // // printf("Why = 0 part # %d?\n",i); => We know. |
551 | // break; | |
88cb7938 | 552 | // we don't break now in order to be sure there is no |
553 | // particle !=0 left. | |
554 | // To be removed later and replaced with break. | |
555 | if(!allFilled) allFilled = kTRUE; | |
556 | } | |
9e1a0ddb | 557 | } |
9e1a0ddb | 558 | //_____________________________________________________________________________ |
88cb7938 | 559 | |
9e1a0ddb | 560 | void AliStack::FlagTrack(Int_t track) |
561 | { | |
562 | // | |
563 | // Flags a track and all its family tree to be kept | |
564 | // | |
565 | ||
566 | TParticle *particle; | |
567 | ||
568 | Int_t curr=track; | |
569 | while(1) { | |
d94af0c1 | 570 | particle=dynamic_cast<TParticle*>(fParticleMap.At(curr)); |
9e1a0ddb | 571 | |
572 | // If the particle is flagged the three from here upward is saved already | |
573 | if(particle->TestBit(kKeepBit)) return; | |
574 | ||
575 | // Save this particle | |
576 | particle->SetBit(kKeepBit); | |
577 | ||
578 | // Move to father if any | |
579 | if((curr=particle->GetFirstMother())==-1) return; | |
580 | } | |
581 | } | |
582 | ||
583 | //_____________________________________________________________________________ | |
ed54bf31 | 584 | void AliStack::KeepTrack(Int_t track) |
9e1a0ddb | 585 | { |
586 | // | |
587 | // Flags a track to be kept | |
588 | // | |
589 | ||
d94af0c1 | 590 | fParticleMap.At(track)->SetBit(kKeepBit); |
9e1a0ddb | 591 | } |
592 | ||
593 | //_____________________________________________________________________________ | |
ae96e41f | 594 | void AliStack::Clean(Int_t size) |
9e1a0ddb | 595 | { |
596 | // | |
ae96e41f | 597 | // Reset stack data except for fTreeK |
9e1a0ddb | 598 | // |
88cb7938 | 599 | |
9e1a0ddb | 600 | fNtrack=0; |
601 | fNprimary=0; | |
602 | fHgwmk=0; | |
603 | fLoadPoint=0; | |
604 | fCurrent = -1; | |
605 | ResetArrays(size); | |
606 | } | |
607 | ||
ae96e41f | 608 | //_____________________________________________________________________________ |
609 | void AliStack::Reset(Int_t size) | |
610 | { | |
611 | // | |
612 | // Reset stack data including fTreeK | |
613 | // | |
614 | ||
615 | Clean(size); | |
c38db470 | 616 | delete fParticleBuffer; fParticleBuffer = 0; |
ae96e41f | 617 | fTreeK = 0x0; |
618 | } | |
619 | ||
9e1a0ddb | 620 | //_____________________________________________________________________________ |
621 | void AliStack::ResetArrays(Int_t size) | |
622 | { | |
623 | // | |
624 | // Resets stack arrays | |
625 | // | |
d94af0c1 | 626 | fParticles.Clear(); |
627 | fParticleMap.Clear(); | |
628 | if (size>0) fParticleMap.Expand(size); | |
9e1a0ddb | 629 | } |
630 | ||
631 | //_____________________________________________________________________________ | |
e2afb3b6 | 632 | void AliStack::SetHighWaterMark(Int_t) |
9e1a0ddb | 633 | { |
634 | // | |
635 | // Set high water mark for last track in event | |
636 | // | |
6308b375 | 637 | |
942a9039 | 638 | fHgwmk = fNtrack-1; |
639 | fCurrentPrimary=fHgwmk; | |
640 | // Set also number of primary tracks | |
641 | fNprimary = fHgwmk+1; | |
9e1a0ddb | 642 | } |
643 | ||
644 | //_____________________________________________________________________________ | |
645 | TParticle* AliStack::Particle(Int_t i) | |
646 | { | |
647 | // | |
648 | // Return particle with specified ID | |
942a9039 | 649 | |
d94af0c1 | 650 | if(!fParticleMap.At(i)) { |
651 | Int_t nentries = fParticles.GetEntriesFast(); | |
9e1a0ddb | 652 | // algorithmic way of getting entry index |
653 | // (primary particles are filled after secondaries) | |
e94530da | 654 | Int_t entry = TreeKEntry(i); |
9e1a0ddb | 655 | // check whether algorithmic way and |
656 | // and the fParticleFileMap[i] give the same; | |
657 | // give the fatal error if not | |
658 | if (entry != fParticleFileMap[i]) { | |
21bf7095 | 659 | AliFatal(Form( |
9e1a0ddb | 660 | "!! The algorithmic way and map are different: !!\n entry: %d map: %d", |
21bf7095 | 661 | entry, fParticleFileMap[i])); |
9e1a0ddb | 662 | } |
942a9039 | 663 | // Load particle at entry into fParticleBuffer |
88cb7938 | 664 | TreeK()->GetEntry(entry); |
942a9039 | 665 | // Add to the TClonesarray |
d94af0c1 | 666 | new (fParticles[nentries]) TParticle(*fParticleBuffer); |
942a9039 | 667 | // Store a pointer in the TObjArray |
d94af0c1 | 668 | fParticleMap.AddAt(fParticles[nentries],i); |
9e1a0ddb | 669 | } |
d94af0c1 | 670 | return dynamic_cast<TParticle*>(fParticleMap.At(i)); |
9e1a0ddb | 671 | } |
672 | ||
e94530da | 673 | //_____________________________________________________________________________ |
674 | TParticle* AliStack::ParticleFromTreeK(Int_t id) const | |
675 | { | |
676 | // | |
677 | // return pointer to TParticle with label id | |
678 | // | |
679 | Int_t entry; | |
680 | if ((entry = TreeKEntry(id)) < 0) return 0; | |
681 | if (fTreeK->GetEntry(entry)<=0) return 0; | |
682 | return fParticleBuffer; | |
683 | } | |
684 | ||
685 | //_____________________________________________________________________________ | |
686 | Int_t AliStack::TreeKEntry(Int_t id) const | |
687 | { | |
688 | // | |
942a9039 | 689 | // Return entry number in the TreeK for particle with label id |
690 | // Return negative number if label>fNtrack | |
691 | // | |
692 | // The order of particles in TreeK reflects the order of the transport of primaries and production of secondaries: | |
e94530da | 693 | // |
942a9039 | 694 | // Before transport there are fNprimary particles on the stack. |
695 | // They are transported one by one and secondaries (fNtrack - fNprimary) are produced. | |
696 | // After the transport of each particles secondaries are written to the TreeK | |
697 | // They occupy the entries 0 ... fNtrack - fNprimary - 1 | |
698 | // The primaries are written after they have been transported and occupy | |
699 | // fNtrack - fNprimary .. fNtrack - 1 | |
700 | ||
e94530da | 701 | Int_t entry; |
702 | if (id<fNprimary) | |
703 | entry = id+fNtrack-fNprimary; | |
704 | else | |
705 | entry = id-fNprimary; | |
706 | return entry; | |
707 | } | |
708 | ||
5d8718b8 | 709 | //_____________________________________________________________________________ |
642f15cf | 710 | Int_t AliStack::GetCurrentParentTrackNumber() const |
5d8718b8 | 711 | { |
712 | // | |
713 | // Return number of the parent of the current track | |
714 | // | |
715 | ||
d94af0c1 | 716 | TParticle* current = (TParticle*)fParticleMap.At(fCurrent); |
5d8718b8 | 717 | |
718 | if (current) | |
719 | return current->GetFirstMother(); | |
720 | else { | |
21bf7095 | 721 | AliWarning("Current track not found in the stack"); |
5d8718b8 | 722 | return -1; |
723 | } | |
724 | } | |
725 | ||
9e1a0ddb | 726 | //_____________________________________________________________________________ |
f5f55563 | 727 | Int_t AliStack::GetPrimary(Int_t id) |
9e1a0ddb | 728 | { |
729 | // | |
730 | // Return number of primary that has generated track | |
731 | // | |
732 | ||
733 | int current, parent; | |
9e1a0ddb | 734 | // |
735 | parent=id; | |
736 | while (1) { | |
737 | current=parent; | |
f5f55563 | 738 | parent=Particle(current)->GetFirstMother(); |
9e1a0ddb | 739 | if(parent<0) return current; |
740 | } | |
741 | } | |
742 | ||
743 | //_____________________________________________________________________________ | |
744 | void AliStack::DumpPart (Int_t i) const | |
745 | { | |
746 | // | |
747 | // Dumps particle i in the stack | |
748 | // | |
d94af0c1 | 749 | dynamic_cast<TParticle*>(fParticleMap.At(i))->Print(); |
9e1a0ddb | 750 | } |
751 | ||
752 | //_____________________________________________________________________________ | |
753 | void AliStack::DumpPStack () | |
754 | { | |
755 | // | |
756 | // Dumps the particle stack | |
757 | // | |
758 | ||
3ab6f951 | 759 | Int_t i; |
760 | ||
88cb7938 | 761 | printf("\n\n=======================================================================\n"); |
3ab6f951 | 762 | for (i=0;i<fNtrack;i++) |
9e1a0ddb | 763 | { |
764 | TParticle* particle = Particle(i); | |
765 | if (particle) { | |
766 | printf("-> %d ",i); particle->Print(); | |
767 | printf("--------------------------------------------------------------\n"); | |
768 | } | |
769 | else | |
770 | Warning("DumpPStack", "No particle with id %d.", i); | |
771 | } | |
772 | ||
88cb7938 | 773 | printf("\n=======================================================================\n\n"); |
9e1a0ddb | 774 | |
775 | // print particle file map | |
776 | printf("\nParticle file map: \n"); | |
3ab6f951 | 777 | for (i=0; i<fNtrack; i++) |
9e1a0ddb | 778 | printf(" %d th entry: %d \n",i,fParticleFileMap[i]); |
779 | } | |
780 | ||
781 | ||
782 | //_____________________________________________________________________________ | |
783 | void AliStack::DumpLoadedStack() const | |
784 | { | |
785 | // | |
786 | // Dumps the particle in the stack | |
787 | // that are loaded in memory. | |
788 | // | |
789 | ||
9e1a0ddb | 790 | printf( |
791 | "\n\n=======================================================================\n"); | |
792 | for (Int_t i=0;i<fNtrack;i++) | |
793 | { | |
d94af0c1 | 794 | TParticle* particle = dynamic_cast<TParticle*>(fParticleMap[i]); |
9e1a0ddb | 795 | if (particle) { |
796 | printf("-> %d ",i); particle->Print(); | |
797 | printf("--------------------------------------------------------------\n"); | |
798 | } | |
799 | else { | |
800 | printf("-> %d Particle not loaded.\n",i); | |
801 | printf("--------------------------------------------------------------\n"); | |
802 | } | |
803 | } | |
804 | printf( | |
805 | "\n=======================================================================\n\n"); | |
806 | } | |
807 | ||
808 | // | |
809 | // protected methods | |
810 | // | |
811 | ||
812 | //_____________________________________________________________________________ | |
813 | void AliStack::CleanParents() | |
814 | { | |
815 | // | |
816 | // Clean particles stack | |
817 | // Set parent/daughter relations | |
818 | // | |
819 | ||
9e1a0ddb | 820 | TParticle *part; |
821 | int i; | |
822 | for(i=0; i<fHgwmk+1; i++) { | |
d94af0c1 | 823 | part = dynamic_cast<TParticle*>(fParticleMap.At(i)); |
9e1a0ddb | 824 | if(part) if(!part->TestBit(kDaughtersBit)) { |
825 | part->SetFirstDaughter(-1); | |
826 | part->SetLastDaughter(-1); | |
827 | } | |
828 | } | |
829 | } | |
830 | ||
831 | //_____________________________________________________________________________ | |
832 | TParticle* AliStack::GetNextParticle() | |
833 | { | |
834 | // | |
835 | // Return next particle from stack of particles | |
836 | // | |
837 | ||
838 | TParticle* particle = 0; | |
839 | ||
840 | // search secondaries | |
841 | //for(Int_t i=fNtrack-1; i>=0; i--) { | |
842 | for(Int_t i=fNtrack-1; i>fHgwmk; i--) { | |
d94af0c1 | 843 | particle = dynamic_cast<TParticle*>(fParticleMap.At(i)); |
9e1a0ddb | 844 | if ((particle) && (!particle->TestBit(kDoneBit))) { |
845 | fCurrent=i; | |
9e1a0ddb | 846 | return particle; |
847 | } | |
848 | } | |
849 | ||
850 | // take next primary if all secondaries were done | |
851 | while (fCurrentPrimary>=0) { | |
852 | fCurrent = fCurrentPrimary; | |
d94af0c1 | 853 | particle = dynamic_cast<TParticle*>(fParticleMap.At(fCurrentPrimary--)); |
9e1a0ddb | 854 | if ((particle) && (!particle->TestBit(kDoneBit))) { |
9e1a0ddb | 855 | return particle; |
856 | } | |
857 | } | |
858 | ||
859 | // nothing to be tracked | |
860 | fCurrent = -1; | |
6308b375 | 861 | |
862 | ||
9e1a0ddb | 863 | return particle; |
864 | } | |
88cb7938 | 865 | //__________________________________________________________________________________________ |
9e1a0ddb | 866 | |
88cb7938 | 867 | TTree* AliStack::TreeK() |
868 | { | |
869 | //returns TreeK | |
942a9039 | 870 | return fTreeK; |
88cb7938 | 871 | } |
9e1a0ddb | 872 | //__________________________________________________________________________________________ |
88cb7938 | 873 | |
942a9039 | 874 | void AliStack::ConnectTree(TTree* tree) |
9e1a0ddb | 875 | { |
876 | // | |
88cb7938 | 877 | // Creates branch for writing particles |
878 | // | |
ae96e41f | 879 | |
942a9039 | 880 | fTreeK = tree; |
881 | ||
21bf7095 | 882 | AliDebug(1, "Connecting TreeK"); |
88cb7938 | 883 | if (fTreeK == 0x0) |
884 | { | |
885 | if (TreeK() == 0x0) | |
886 | { | |
21bf7095 | 887 | AliFatal("Parameter is NULL");//we don't like such a jokes |
88cb7938 | 888 | return; |
889 | } | |
890 | return;//in this case TreeK() calls back this method (ConnectTree) | |
891 | //tree after setting fTreeK, the rest was already executed | |
892 | //it is safe to return now | |
893 | } | |
894 | ||
895 | // Create a branch for particles | |
896 | ||
21bf7095 | 897 | AliDebug(2, Form("Tree name is %s",fTreeK->GetName())); |
88cb7938 | 898 | |
899 | if (fTreeK->GetDirectory()) | |
900 | { | |
21bf7095 | 901 | AliDebug(2, Form("and dir is %s",fTreeK->GetDirectory()->GetName())); |
88cb7938 | 902 | } |
903 | else | |
21bf7095 | 904 | AliWarning("DIR IS NOT SET !!!"); |
88cb7938 | 905 | |
942a9039 | 906 | TBranch *branch=fTreeK->GetBranch("Particles"); |
88cb7938 | 907 | if(branch == 0x0) |
908 | { | |
942a9039 | 909 | branch = fTreeK->Branch("Particles", "TParticle", &fParticleBuffer, 4000); |
21bf7095 | 910 | AliDebug(2, "Creating Branch in Tree"); |
88cb7938 | 911 | } |
912 | else | |
913 | { | |
21bf7095 | 914 | AliDebug(2, "Branch Found in Tree"); |
88cb7938 | 915 | branch->SetAddress(&fParticleBuffer); |
916 | } | |
917 | if (branch->GetDirectory()) | |
918 | { | |
21bf7095 | 919 | AliDebug(1, Form("Branch Dir Name is %s",branch->GetDirectory()->GetName())); |
88cb7938 | 920 | } |
921 | else | |
21bf7095 | 922 | AliWarning("Branch Dir is NOT SET"); |
9e1a0ddb | 923 | } |
88cb7938 | 924 | |
5d8718b8 | 925 | //_____________________________________________________________________________ |
88cb7938 | 926 | |
927 | Bool_t AliStack::GetEvent() | |
9e1a0ddb | 928 | { |
929 | // | |
930 | // Get new event from TreeK | |
ccf7a81f | 931 | |
9e1a0ddb | 932 | // Reset/Create the particle stack |
88cb7938 | 933 | Int_t size = (Int_t)TreeK()->GetEntries(); |
9e1a0ddb | 934 | ResetArrays(size); |
ccf7a81f | 935 | return kTRUE; |
9e1a0ddb | 936 | } |
88cb7938 | 937 | //_____________________________________________________________________________ |
938 | ||
1fed0e8b | 939 | Bool_t AliStack::IsStable(Int_t pdg) const |
940 | { | |
941 | // | |
942 | // Decide whether particle (pdg) is stable | |
943 | // | |
944 | ||
945 | const Int_t kNstable = 14; | |
946 | Int_t i; | |
947 | ||
948 | Int_t pdgStable[kNstable] = { | |
949 | kGamma, // Photon | |
950 | kElectron, // Electron | |
951 | kMuonPlus, // Muon | |
952 | kPiPlus, // Pion | |
953 | kKPlus, // Kaon | |
954 | kProton, // Proton | |
955 | kNeutron, // Neutron | |
956 | kLambda0, // Lambda_0 | |
957 | kSigmaMinus, // Sigma Minus | |
958 | kSigma0, // Sigma_0 | |
959 | kSigmaPlus, // Sigma Plus | |
960 | 3312, // Xsi Minus | |
961 | 3322, // Xsi | |
962 | 3334, // Omega | |
963 | }; | |
964 | ||
965 | Bool_t isStable = kFALSE; | |
966 | for (i = 0; i < kNstable; i++) { | |
967 | if (pdg == TMath::Abs(pdgStable[i])) { | |
968 | isStable = kTRUE; | |
969 | break; | |
970 | } | |
971 | } | |
972 | ||
973 | return isStable; | |
974 | } | |
975 | ||
976 | Bool_t AliStack::IsPhysicalPrimary(Int_t index) | |
977 | { | |
978 | // | |
979 | // Test if a particle is a physical primary according to the following definition: | |
980 | // Particles produced in the collision including products of strong and | |
981 | // electromagnetic decay and excluding feed-down from weak decays of strange | |
982 | // particles. | |
983 | // | |
984 | TParticle* p = Particle(index); | |
985 | Int_t ist = p->GetStatusCode(); | |
986 | ||
987 | // | |
988 | // Initial state particle | |
989 | if (ist > 20) return kFALSE; | |
990 | ||
991 | Int_t pdg = TMath::Abs(p->GetPdgCode()); | |
992 | ||
993 | if (!IsStable(pdg)) return kFALSE; | |
994 | ||
995 | if (index < GetNprimary()) { | |
996 | // | |
997 | // Particle produced by generator | |
998 | return kTRUE; | |
999 | } else { | |
1000 | // | |
1001 | // Particle produced during transport | |
1002 | // | |
1b13a70f | 1003 | |
1fed0e8b | 1004 | Int_t imo = p->GetFirstMother(); |
1005 | TParticle* pm = Particle(imo); | |
1006 | Int_t mpdg = TMath::Abs(pm->GetPdgCode()); | |
1b13a70f | 1007 | // Check if it comes from a pi0 decay |
1008 | // | |
1009 | // What about the pi0 Dalitz ?? | |
1010 | // if ((mpdg == kPi0) && (imo < GetNprimary())) return kTRUE; | |
1011 | ||
1012 | // Check if this is a heavy flavor decay product | |
1fed0e8b | 1013 | Int_t mfl = Int_t (mpdg / TMath::Power(10, Int_t(TMath::Log10(mpdg)))); |
1014 | // | |
1015 | // Light hadron | |
1016 | if (mfl < 4) return kFALSE; | |
1017 | ||
1018 | // | |
1019 | // Heavy flavor hadron produced by generator | |
1020 | if (imo < GetNprimary()) { | |
1021 | return kTRUE; | |
1022 | } | |
1023 | ||
1024 | // To be sure that heavy flavor has not been produced in a secondary interaction | |
1025 | // Loop back to the generated mother | |
1026 | while (imo >= GetNprimary()) { | |
1b13a70f | 1027 | imo = pm->GetFirstMother(); |
1fed0e8b | 1028 | pm = Particle(imo); |
1029 | } | |
1030 | mpdg = TMath::Abs(pm->GetPdgCode()); | |
1031 | mfl = Int_t (mpdg / TMath::Power(10, Int_t(TMath::Log10(mpdg)))); | |
1032 | ||
1033 | if (mfl < 4) { | |
1034 | return kFALSE; | |
1035 | } else { | |
1036 | return kTRUE; | |
1037 | } | |
1038 | } // produced by generator ? | |
1039 | } |