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