]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliStack.cxx
Fixes for bug #52499: Field polarities inconsistiency
[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 = GetParticleMapEntry(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=GetParticleMapEntry(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 = GetParticleMapEntry(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 = GetParticleMapEntry(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 = GetParticleMapEntry(i);
348       parent = part->GetFirstMother();
349       if(parent>=0) {
350           father = GetParticleMapEntry(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 = GetParticleMapEntry(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] = GetParticleMapEntry(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 = GetParticleMapEntry(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 = GetParticleMapEntry(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   TParticle *part;
537   for(Int_t i=0; i<fHgwmk+1; ++i) {
538     if((part=GetParticleMapEntry(i))) {
539       fParticleBuffer = 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
562 void AliStack::FlagTrack(Int_t track)
563 {
564   //
565   // Flags a track and all its family tree to be kept
566   //
567   
568   TParticle *particle;
569
570   Int_t curr=track;
571   while(1) {
572     particle = GetParticleMapEntry(curr);
573     
574     // If the particle is flagged the three from here upward is saved already
575     if(particle->TestBit(kKeepBit)) return;
576     
577     // Save this particle
578     particle->SetBit(kKeepBit);
579     
580     // Move to father if any
581     if((curr=particle->GetFirstMother())==-1) return;
582   }
583 }
584  
585 //_____________________________________________________________________________
586 void AliStack::KeepTrack(Int_t track)
587
588   //
589   // Flags a track to be kept
590   //
591   
592   fParticleMap.At(track)->SetBit(kKeepBit);
593 }
594
595 //_____________________________________________________________________________
596 void  AliStack::Clean(Int_t size) 
597 {
598   //
599   // Reset stack data except for fTreeK
600   //
601   
602   fNtrack=0;
603   fNprimary=0;
604   fHgwmk=0;
605   fLoadPoint=0;
606   fCurrent = -1;
607   ResetArrays(size);
608 }
609
610 //_____________________________________________________________________________
611 void  AliStack::Reset(Int_t size) 
612 {
613   //
614   // Reset stack data including fTreeK
615   //
616
617   Clean(size);
618   delete fParticleBuffer; fParticleBuffer = 0;
619   fTreeK = 0x0;
620 }
621
622 //_____________________________________________________________________________
623 void  AliStack::ResetArrays(Int_t size) 
624 {
625   //
626   // Resets stack arrays
627   //
628   fParticles.Clear();
629   fParticleMap.Clear();
630   if (size>0) fParticleMap.Expand(size);
631 }
632
633 //_____________________________________________________________________________
634 void AliStack::SetHighWaterMark(Int_t)
635 {
636   //
637   // Set high water mark for last track in event
638   //
639     
640     fHgwmk = fNtrack-1;
641     fCurrentPrimary=fHgwmk;
642     // Set also number of primary tracks
643     fNprimary = fHgwmk+1;
644 }
645
646 //_____________________________________________________________________________
647 TParticle* AliStack::Particle(Int_t i)
648 {
649   //
650   // Return particle with specified ID
651
652   if(!fParticleMap.At(i)) {
653     Int_t nentries = fParticles.GetEntriesFast();
654     // algorithmic way of getting entry index
655     // (primary particles are filled after secondaries)
656     Int_t entry = TreeKEntry(i);
657     // check whether algorithmic way and 
658     // and the fParticleFileMap[i] give the same;
659     // give the fatal error if not
660     if (entry != fParticleFileMap[i]) {
661       AliFatal(Form(
662         "!! The algorithmic way and map are different: !!\n entry: %d map: %d",
663         entry, fParticleFileMap[i])); 
664     } 
665     // Load particle at entry into fParticleBuffer
666     TreeK()->GetEntry(entry);
667     // Add to the TClonesarray
668     new (fParticles[nentries]) TParticle(*fParticleBuffer);
669     // Store a pointer in the TObjArray
670     fParticleMap.AddAt(fParticles[nentries],i);
671   }
672   return GetParticleMapEntry(i);
673 }
674
675 //_____________________________________________________________________________
676 TParticle* AliStack::ParticleFromTreeK(Int_t id) const
677 {
678 // 
679 // return pointer to TParticle with label id
680 //
681   Int_t entry;
682   if ((entry = TreeKEntry(id)) < 0) return 0;
683   if (fTreeK->GetEntry(entry)<=0) return 0;
684   return fParticleBuffer;
685 }
686
687 //_____________________________________________________________________________
688 Int_t AliStack::TreeKEntry(Int_t id) const 
689 {
690 //
691 // Return entry number in the TreeK for particle with label id
692 // Return negative number if label>fNtrack
693 //
694 // The order of particles in TreeK reflects the order of the transport of primaries and production of secondaries:
695 //
696 // Before transport there are fNprimary particles on the stack.
697 // They are transported one by one and secondaries (fNtrack - fNprimary) are produced. 
698 // After the transport of each particles secondaries are written to the TreeK
699 // They occupy the entries 0 ... fNtrack - fNprimary - 1
700 // The primaries are written after they have been transported and occupy 
701 // fNtrack - fNprimary .. fNtrack - 1
702
703   Int_t entry;
704   if (id<fNprimary)
705     entry = id+fNtrack-fNprimary;
706   else 
707     entry = id-fNprimary;
708   return entry;
709 }
710
711 //_____________________________________________________________________________
712 Int_t AliStack::GetCurrentParentTrackNumber() const
713 {
714   //
715   // Return number of the parent of the current track
716   //
717   
718   TParticle* current = GetParticleMapEntry(fCurrent);
719
720   if (current) 
721     return current->GetFirstMother();
722   else {
723     AliWarning("Current track not found in the stack");
724     return -1;
725   }  
726 }
727  
728 //_____________________________________________________________________________
729 Int_t AliStack::GetPrimary(Int_t id)
730 {
731   //
732   // Return number of primary that has generated track
733   //
734   
735   int current, parent;
736   //
737   parent=id;
738   while (1) {
739     current=parent;
740     parent=Particle(current)->GetFirstMother();
741     if(parent<0) return current;
742   }
743 }
744  
745 //_____________________________________________________________________________
746 void AliStack::DumpPart (Int_t i) const
747 {
748   //
749   // Dumps particle i in the stack
750   //
751   GetParticleMapEntry(i)->Print();
752 }
753
754 //_____________________________________________________________________________
755 void AliStack::DumpPStack ()
756 {
757   //
758   // Dumps the particle stack
759   //
760
761   Int_t i;
762
763   printf("\n\n=======================================================================\n");
764   for (i=0;i<fNtrack;i++) 
765     {
766       TParticle* particle = Particle(i);
767       if (particle) {
768         printf("-> %d ",i); particle->Print();
769         printf("--------------------------------------------------------------\n");
770       }
771       else 
772         Warning("DumpPStack", "No particle with id %d.", i); 
773     }    
774
775   printf("\n=======================================================================\n\n");
776   
777   // print  particle file map
778   // printf("\nParticle file map: \n");
779   // for (i=0; i<fNtrack; i++) 
780   //     printf("   %d th entry: %d \n",i,fParticleFileMap[i]);
781 }
782
783
784 //_____________________________________________________________________________
785 void AliStack::DumpLoadedStack() const
786 {
787   //
788   // Dumps the particle in the stack
789   // that are loaded in memory.
790   //
791
792   printf(
793          "\n\n=======================================================================\n");
794   for (Int_t i=0;i<fNtrack;i++) 
795     {
796       TParticle* particle = GetParticleMapEntry(i);
797       if (particle) {
798         printf("-> %d ",i); particle->Print();
799         printf("--------------------------------------------------------------\n");
800       }
801       else {    
802         printf("-> %d  Particle not loaded.\n",i);
803         printf("--------------------------------------------------------------\n");
804       } 
805     }
806   printf(
807          "\n=======================================================================\n\n");
808 }
809
810 //
811 // protected methods
812 //
813
814 //_____________________________________________________________________________
815 void AliStack::CleanParents()
816 {
817   //
818   // Clean particles stack
819   // Set parent/daughter relations
820   //
821   
822   TParticle *part;
823   int i;
824   for(i=0; i<fHgwmk+1; i++) {
825     part = GetParticleMapEntry(i);
826     if(part) if(!part->TestBit(kDaughtersBit)) {
827       part->SetFirstDaughter(-1);
828       part->SetLastDaughter(-1);
829     }
830   }
831 }
832
833 //_____________________________________________________________________________
834 TParticle* AliStack::GetNextParticle()
835 {
836   //
837   // Return next particle from stack of particles
838   //
839   
840   TParticle* particle = 0;
841   
842   // search secondaries
843   //for(Int_t i=fNtrack-1; i>=0; i--) {
844   for(Int_t i=fNtrack-1; i>fHgwmk; i--) {
845       particle = GetParticleMapEntry(i);
846       if ((particle) && (!particle->TestBit(kDoneBit))) {
847           fCurrent=i;    
848           return particle;
849       }   
850   }    
851
852   // take next primary if all secondaries were done
853   while (fCurrentPrimary>=0) {
854       fCurrent = fCurrentPrimary;    
855       particle = GetParticleMapEntry(fCurrentPrimary--);
856       if ((particle) && (!particle->TestBit(kDoneBit))) {
857           return particle;
858       } 
859   }
860   
861   // nothing to be tracked
862   fCurrent = -1;
863  
864   
865   return particle;  
866 }
867 //__________________________________________________________________________________________
868
869 void AliStack::ConnectTree(TTree* tree)
870 {
871 //
872 //  Creates branch for writing particles
873 //
874
875   fTreeK = tree;
876     
877   AliDebug(1, "Connecting TreeK");
878   if (fTreeK == 0x0)
879    {
880     if (TreeK() == 0x0)
881      {
882       AliFatal("Parameter is NULL");//we don't like such a jokes
883       return;
884      }
885     return;//in this case TreeK() calls back this method (ConnectTree) 
886            //tree after setting fTreeK, the rest was already executed
887            //it is safe to return now
888    }
889
890  //  Create a branch for particles   
891   
892   AliDebug(2, Form("Tree name is %s",fTreeK->GetName()));
893    
894   if (fTreeK->GetDirectory())
895    {
896      AliDebug(2, Form("and dir is %s",fTreeK->GetDirectory()->GetName()));
897    }    
898   else
899     AliWarning("DIR IS NOT SET !!!");
900   
901   TBranch *branch=fTreeK->GetBranch("Particles");
902   if(branch == 0x0)
903    {
904     branch = fTreeK->Branch("Particles", "TParticle", &fParticleBuffer, 4000);
905     AliDebug(2, "Creating Branch in Tree");
906    }  
907   else
908    {
909     AliDebug(2, "Branch Found in Tree");
910     branch->SetAddress(&fParticleBuffer);
911    }
912   if (branch->GetDirectory())
913    {
914     AliDebug(1, Form("Branch Dir Name is %s",branch->GetDirectory()->GetName()));
915    } 
916   else
917     AliWarning("Branch Dir is NOT SET");
918 }
919
920 //_____________________________________________________________________________
921
922 Bool_t AliStack::GetEvent()
923 {
924 //
925 // Get new event from TreeK
926
927     // Reset/Create the particle stack
928     Int_t size = (Int_t)TreeK()->GetEntries();
929     ResetArrays(size);
930     return kTRUE;
931 }
932 //_____________________________________________________________________________
933
934 Bool_t AliStack::IsStable(Int_t pdg) const
935 {
936   //
937   // Decide whether particle (pdg) is stable
938   //
939   
940   
941   // All ions/nucleons are considered as stable
942   // Nuclear code is 10LZZZAAAI
943   if(pdg>1000000000)return kTRUE;
944
945   const Int_t kNstable = 15;
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     kK0Short,           // K0s
955     kK0Long,            // K0l
956     kProton,            // Proton 
957     kNeutron,           // Neutron
958     kLambda0,           // Lambda_0
959     kSigmaMinus,        // Sigma Minus
960     kSigmaPlus,         // Sigma Plus
961     3312,               // Xsi Minus 
962     3322,               // Xsi 
963     3334,               // Omega
964   };
965     
966   Bool_t isStable = kFALSE;
967   for (i = 0; i < kNstable; i++) {
968     if (pdg == TMath::Abs(pdgStable[i])) {
969       isStable = kTRUE;
970       break;
971     }
972   }
973   
974   return isStable;
975 }
976
977 Bool_t AliStack::IsPhysicalPrimary(Int_t index)
978 {
979     //
980     // Test if a particle is a physical primary according to the following definition:
981     // Particles produced in the collision including products of strong and
982     // electromagnetic decay and excluding feed-down from weak decays of strange
983     // particles.
984     //
985     TParticle* p = Particle(index);
986     Int_t ist = p->GetStatusCode();
987     
988     //
989     // Initial state particle
990     if (ist > 1) return kFALSE;
991     
992     Int_t pdg = TMath::Abs(p->GetPdgCode());
993     
994     if (!IsStable(pdg)) return kFALSE;
995     
996     if (index < GetNprimary()) {
997 //
998 // Particle produced by generator
999         return kTRUE;
1000     } else {
1001 //
1002 // Particle produced during transport
1003 //
1004
1005         Int_t imo =  p->GetFirstMother();
1006         TParticle* pm  = Particle(imo);
1007         Int_t mpdg = TMath::Abs(pm->GetPdgCode());
1008 // Check if it comes from a pi0 decay
1009 //
1010 // What about the pi0 Dalitz ??
1011 //      if ((mpdg == kPi0) && (imo < GetNprimary())) return kTRUE; 
1012
1013 // Check if this is a heavy flavor decay product
1014         Int_t mfl  = Int_t (mpdg / TMath::Power(10, Int_t(TMath::Log10(mpdg))));
1015         //
1016         // Light hadron
1017         if (mfl < 4) return kFALSE;
1018         
1019         //
1020         // Heavy flavor hadron produced by generator
1021         if (imo <  GetNprimary()) {
1022             return kTRUE;
1023         }
1024         
1025         // To be sure that heavy flavor has not been produced in a secondary interaction
1026         // Loop back to the generated mother
1027         while (imo >=  GetNprimary()) {
1028             imo = pm->GetFirstMother();
1029             pm  =  Particle(imo);
1030         }
1031         mpdg = TMath::Abs(pm->GetPdgCode());
1032         mfl  = Int_t (mpdg / TMath::Power(10, Int_t(TMath::Log10(mpdg))));
1033
1034         if (mfl < 4) {
1035             return kFALSE;
1036         } else {
1037             return kTRUE;
1038         } 
1039     } // produced by generator ?
1040