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