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