]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliStack.cxx
Classes for CTP time parameters added (Plamen)
[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   fParticleMap.AddAtAndExpand(particle, fNtrack);//CHECK!!
230
231   if(parent>=0) {
232       particle = GetParticleMapEntry(parent);
233       if (particle) {
234           particle->SetLastDaughter(fNtrack);
235           if(particle->GetFirstDaughter()<0) particle->SetFirstDaughter(fNtrack);
236       }
237       else {
238           AliError(Form("Parent %d does not exist",parent));
239       }
240   } else { 
241       //
242       // This is a primary track. Set high water mark for this event
243       fHgwmk = fNtrack;
244       //
245       // Set also number if primary tracks
246       fNprimary = fHgwmk+1;
247       fCurrentPrimary++;
248   }
249   ntr = fNtrack++;
250 }
251
252 //_____________________________________________________________________________
253 TParticle*  AliStack::PopNextTrack(Int_t& itrack)
254 {
255   //
256   // Returns next track from stack of particles
257   //
258   
259
260   TParticle* track = GetNextParticle();
261
262   if (track) {
263     itrack = fCurrent;
264     track->SetBit(kDoneBit);
265   }
266   else
267     itrack = -1;
268   
269   fCurrentTrack = track;
270   return track;
271 }
272
273 //_____________________________________________________________________________
274 TParticle*  AliStack::PopPrimaryForTracking(Int_t i)
275 {
276   //
277   // Returns i-th primary particle if it is flagged to be tracked,
278   // 0 otherwise
279   //
280   
281   TParticle* particle = Particle(i);
282   
283   if (!particle->TestBit(kDoneBit)) {
284     fCurrentTrack = particle;
285     return particle;
286   }
287   else
288     return 0;
289 }      
290
291 //_____________________________________________________________________________
292 Bool_t AliStack::PurifyKine()
293 {
294   //
295   // Compress kinematic tree keeping only flagged particles
296   // and renaming the particle id's in all the hits
297   //
298
299   int nkeep = fHgwmk + 1, parent, i;
300   TParticle *part, *father;
301   fTrackLabelMap.Set(fParticleMap.GetLast()+1);
302
303   // Save in Header total number of tracks before compression
304   // If no tracks generated return now
305   if(fHgwmk+1 == fNtrack) return kFALSE;
306
307   // First pass, invalid Daughter information
308   for(i=0; i<fNtrack; i++) {
309       // Preset map, to be removed later
310       if(i<=fHgwmk) fTrackLabelMap[i]=i ; 
311       else {
312           fTrackLabelMap[i] = -99;
313           if((part=GetParticleMapEntry(i))) {
314 //
315 //        Check of this track should be kept for physics reasons 
316               if (KeepPhysics(part)) KeepTrack(i);
317 //
318               part->ResetBit(kDaughtersBit);
319               part->SetFirstDaughter(-1);
320               part->SetLastDaughter(-1);
321           }
322       }
323   }
324   // Invalid daughter information for the parent of the first particle
325   // generated. This may or may not be the current primary according to
326   // whether decays have been recorded among the primaries
327   part = GetParticleMapEntry(fHgwmk+1);
328   fParticleMap.At(part->GetFirstMother())->ResetBit(kDaughtersBit);
329   // Second pass, build map between old and new numbering
330   for(i=fHgwmk+1; i<fNtrack; i++) {
331       if(fParticleMap.At(i)->TestBit(kKeepBit)) {
332           // This particle has to be kept
333           fTrackLabelMap[i]=nkeep;
334           // If old and new are different, have to move the pointer
335           if(i!=nkeep) fParticleMap[nkeep]=fParticleMap.At(i);
336           part = GetParticleMapEntry(nkeep);
337           // as the parent is always *before*, it must be already
338           // in place. This is what we are checking anyway!
339           if((parent=part->GetFirstMother())>fHgwmk) {
340               if(fTrackLabelMap[parent]==-99) Fatal("PurifyKine","fTrackLabelMap[%d] = -99!\n",parent);
341               else part->SetFirstMother(fTrackLabelMap[parent]);}
342           nkeep++;
343       }
344   }
345   
346   // Fix daughters information
347   for (i=fHgwmk+1; i<nkeep; i++) {
348       part = GetParticleMapEntry(i);
349       parent = part->GetFirstMother();
350       if(parent>=0) {
351           father = GetParticleMapEntry(parent);
352           if(father->TestBit(kDaughtersBit)) {
353               
354               if(i<father->GetFirstDaughter()) father->SetFirstDaughter(i);
355               if(i>father->GetLastDaughter())  father->SetLastDaughter(i);
356           } else {
357               // Initialise daughters info for first pass
358               father->SetFirstDaughter(i);
359               father->SetLastDaughter(i);
360               father->SetBit(kDaughtersBit);
361           }
362       }
363   }
364   //
365   // Now the output bit, from fHgwmk to nkeep we write everything and we erase
366   if(nkeep > fParticleFileMap.GetSize()) fParticleFileMap.Set(Int_t (nkeep*1.5));
367   for (i=fHgwmk+1; i<nkeep; ++i) {
368       fParticleBuffer = GetParticleMapEntry(i);
369       fParticleFileMap[i]=static_cast<Int_t>(TreeK()->GetEntries());
370       TreeK()->Fill();
371       fParticleMap[i]=fParticleBuffer=0;
372   }
373   
374   for (i = nkeep; i < fNtrack; ++i) fParticleMap[i]=0;
375   
376   Int_t toshrink = fNtrack-fHgwmk-1;
377   fLoadPoint-=toshrink;
378   
379   for(i=fLoadPoint; i<fLoadPoint+toshrink; ++i) fParticles.RemoveAt(i);
380   fNtrack=nkeep;
381   fHgwmk=nkeep-1;
382   return kTRUE;
383 }
384
385
386 Bool_t AliStack::ReorderKine()
387 {
388 //
389 // In some transport code children might not come in a continuous sequence.
390 // In this case the stack  has  to  be reordered in order to establish the 
391 // mother daughter relation using index ranges.
392 //    
393   if(fHgwmk+1 == fNtrack) return kFALSE;
394
395   //
396   // Howmany secondaries have been produced ?
397   Int_t nNew = fNtrack - fHgwmk - 1;
398     
399   if (nNew > 0) {
400       Int_t i, j;
401       TArrayI map1(nNew);
402       //
403       // Copy pointers to temporary array
404       TParticle** tmp = new TParticle*[nNew];
405       
406       for (i = 0; i < nNew; i++) {
407           if (fParticleMap.At(fHgwmk + 1 + i)) {
408               tmp[i] = GetParticleMapEntry(fHgwmk + 1 + i);
409           } else {
410               tmp[i] = 0x0;
411           }
412           map1[i] = -99;
413       }
414   
415       
416       //
417       // Reset  LoadPoint 
418       // 
419       Int_t loadPoint = fHgwmk + 1;
420       //
421       // Re-Push particles into stack 
422       // The outer loop is over parents, the inner over children.
423       // -1 refers to the primary particle
424       //
425       for (i = -1; i < nNew-1; i++) {
426           Int_t ipa;
427           TParticle* parP;
428           if (i == -1) {
429               ipa  = tmp[0]->GetFirstMother();
430               parP = GetParticleMapEntry(ipa);
431           } else {
432               ipa = (fHgwmk + 1 + i);
433               // Skip deleted particles
434               if (!tmp[i])                          continue;
435               // Skip particles without children
436               if (tmp[i]->GetFirstDaughter() == -1) continue;
437               parP = tmp[i];
438           }
439           // Reset daughter information
440
441           Int_t idaumin = parP->GetFirstDaughter() - fHgwmk - 1;
442           Int_t idaumax = parP->GetLastDaughter()  - fHgwmk - 1;
443           parP->SetFirstDaughter(-1);
444           parP->SetLastDaughter(-1);
445           for (j = idaumin; j <= idaumax; j++) {
446               // Skip deleted particles
447               if (!tmp[j])        continue;
448               // Skip particles already handled
449               if (map1[j] != -99) continue;
450               Int_t jpa = tmp[j]->GetFirstMother();
451               // Check if daughter of current parent
452               if (jpa == ipa) {
453                   fParticleMap[loadPoint] = tmp[j];
454                   // Re-establish daughter information
455                   parP->SetLastDaughter(loadPoint);
456                   if (parP->GetFirstDaughter() == -1) parP->SetFirstDaughter(loadPoint);
457                   // Set Mother information
458                   if (i != -1) {
459                       tmp[j]->SetFirstMother(map1[i]);
460                   } 
461                   // Build the map
462                   map1[j] = loadPoint;
463                   // Increase load point
464                   loadPoint++;
465               }
466           } // children
467       } // parents
468
469       delete[] tmp;
470
471       //
472       // Build map for remapping of hits
473       // 
474       fTrackLabelMap.Set(fNtrack);
475       for (i = 0; i < fNtrack; i ++) {
476           if (i <= fHgwmk) {
477               fTrackLabelMap[i] = i;
478           } else{
479               fTrackLabelMap[i] = map1[i - fHgwmk -1];
480           }
481       }
482   } // new particles poduced
483   
484   return kTRUE;
485 }
486
487 Bool_t AliStack::KeepPhysics(const TParticle* part)
488 {
489     //
490     // Some particles have to kept on the stack for reasons motivated
491     // by physics analysis. Decision is put here.
492     //
493     Bool_t keep = kFALSE;
494     //
495     // Keep first-generation daughter from primaries with heavy flavor 
496     //
497     Int_t parent = part->GetFirstMother();
498     if (parent >= 0 && parent <= fHgwmk) {
499         TParticle* father = GetParticleMapEntry(parent);
500         Int_t kf = father->GetPdgCode();
501         kf = TMath::Abs(kf);
502         Int_t kfl = kf;
503         // meson ?
504         if  (kfl > 10) kfl/=100;
505         // baryon
506         if (kfl > 10) kfl/=10;
507         if (kfl > 10) kfl/=10;
508         if (kfl >= 4) {
509             keep = kTRUE;
510         }
511     }
512     return keep;
513 }
514
515 //_____________________________________________________________________________
516 void AliStack::FinishEvent()
517 {
518 //
519 // Write out the kinematics that was not yet filled
520 //
521   
522 // Update event header
523
524   if (!TreeK()) {
525 //    Fatal("FinishEvent", "No kinematics tree is defined.");
526 //    Don't panic this is a probably a lego run
527       return;
528   }  
529   
530   CleanParents();
531    if(TreeK()->GetEntries() ==0) {
532     // set the fParticleFileMap size for the first time
533     fParticleFileMap.Set(fHgwmk+1);
534   }
535
536   Bool_t allFilled = kFALSE;
537   TParticle *part;
538   for(Int_t i=0; i<fHgwmk+1; ++i) {
539     if((part=GetParticleMapEntry(i))) {
540       fParticleBuffer = part;
541       fParticleFileMap[i]= static_cast<Int_t>(TreeK()->GetEntries());
542       TreeK()->Fill();
543       fParticleBuffer=0;      
544       fParticleMap.AddAt(0,i);      
545       
546       // When all primaries were filled no particle!=0
547       // should be left => to be removed later.
548       if (allFilled) AliWarning(Form("Why != 0 part # %d?\n",i));
549     }
550     else 
551     {
552       // // printf("Why = 0 part # %d?\n",i); => We know.
553       // break;
554       // we don't break now in order to be sure there is no
555       // particle !=0 left.
556       // To be removed later and replaced with break.
557        if(!allFilled) allFilled = kTRUE;
558     }
559   }
560
561 //_____________________________________________________________________________
562
563 void AliStack::FlagTrack(Int_t track)
564 {
565   //
566   // Flags a track and all its family tree to be kept
567   //
568   
569   TParticle *particle;
570
571   Int_t curr=track;
572   while(1) {
573     particle = GetParticleMapEntry(curr);
574     
575     // If the particle is flagged the three from here upward is saved already
576     if(particle->TestBit(kKeepBit)) return;
577     
578     // Save this particle
579     particle->SetBit(kKeepBit);
580     
581     // Move to father if any
582     if((curr=particle->GetFirstMother())==-1) return;
583   }
584 }
585  
586 //_____________________________________________________________________________
587 void AliStack::KeepTrack(Int_t track)
588
589   //
590   // Flags a track to be kept
591   //
592   
593   fParticleMap.At(track)->SetBit(kKeepBit);
594 }
595
596 //_____________________________________________________________________________
597 void  AliStack::Clean(Int_t size) 
598 {
599   //
600   // Reset stack data except for fTreeK
601   //
602   
603   fNtrack=0;
604   fNprimary=0;
605   fHgwmk=0;
606   fLoadPoint=0;
607   fCurrent = -1;
608   ResetArrays(size);
609 }
610
611 //_____________________________________________________________________________
612 void  AliStack::Reset(Int_t size) 
613 {
614   //
615   // Reset stack data including fTreeK
616   //
617
618   Clean(size);
619   delete fParticleBuffer; fParticleBuffer = 0;
620   fTreeK = 0x0;
621 }
622
623 //_____________________________________________________________________________
624 void  AliStack::ResetArrays(Int_t size) 
625 {
626   //
627   // Resets stack arrays
628   //
629   fParticles.Clear();
630   fParticleMap.Clear();
631   if (size>0) fParticleMap.Expand(size);
632 }
633
634 //_____________________________________________________________________________
635 void AliStack::SetHighWaterMark(Int_t)
636 {
637   //
638   // Set high water mark for last track in event
639   //
640     
641     fHgwmk = fNtrack-1;
642     fCurrentPrimary=fHgwmk;
643     // Set also number of primary tracks
644     fNprimary = fHgwmk+1;
645 }
646
647 //_____________________________________________________________________________
648 TParticle* AliStack::Particle(Int_t i)
649 {
650   //
651   // Return particle with specified ID
652
653   if(!fParticleMap.At(i)) {
654     Int_t nentries = fParticles.GetEntriesFast();
655     // algorithmic way of getting entry index
656     // (primary particles are filled after secondaries)
657     Int_t entry = TreeKEntry(i);
658     // check whether algorithmic way and 
659     // and the fParticleFileMap[i] give the same;
660     // give the fatal error if not
661     if (entry != fParticleFileMap[i]) {
662       AliFatal(Form(
663         "!! The algorithmic way and map are different: !!\n entry: %d map: %d",
664         entry, fParticleFileMap[i])); 
665     } 
666     // Load particle at entry into fParticleBuffer
667     TreeK()->GetEntry(entry);
668     // Add to the TClonesarray
669     new (fParticles[nentries]) TParticle(*fParticleBuffer);
670     // Store a pointer in the TObjArray
671     fParticleMap.AddAt(fParticles[nentries],i);
672   }
673   return GetParticleMapEntry(i);
674 }
675
676 //_____________________________________________________________________________
677 TParticle* AliStack::ParticleFromTreeK(Int_t id) const
678 {
679 // 
680 // return pointer to TParticle with label id
681 //
682   Int_t entry;
683   if ((entry = TreeKEntry(id)) < 0) return 0;
684   if (fTreeK->GetEntry(entry)<=0) return 0;
685   return fParticleBuffer;
686 }
687
688 //_____________________________________________________________________________
689 Int_t AliStack::TreeKEntry(Int_t id) const 
690 {
691 //
692 // Return entry number in the TreeK for particle with label id
693 // Return negative number if label>fNtrack
694 //
695 // The order of particles in TreeK reflects the order of the transport of primaries and production of secondaries:
696 //
697 // Before transport there are fNprimary particles on the stack.
698 // They are transported one by one and secondaries (fNtrack - fNprimary) are produced. 
699 // After the transport of each particles secondaries are written to the TreeK
700 // They occupy the entries 0 ... fNtrack - fNprimary - 1
701 // The primaries are written after they have been transported and occupy 
702 // fNtrack - fNprimary .. fNtrack - 1
703
704   Int_t entry;
705   if (id<fNprimary)
706     entry = id+fNtrack-fNprimary;
707   else 
708     entry = id-fNprimary;
709   return entry;
710 }
711
712 //_____________________________________________________________________________
713 Int_t AliStack::GetCurrentParentTrackNumber() const
714 {
715   //
716   // Return number of the parent of the current track
717   //
718   
719   TParticle* current = GetParticleMapEntry(fCurrent);
720
721   if (current) 
722     return current->GetFirstMother();
723   else {
724     AliWarning("Current track not found in the stack");
725     return -1;
726   }  
727 }
728  
729 //_____________________________________________________________________________
730 Int_t AliStack::GetPrimary(Int_t id)
731 {
732   //
733   // Return number of primary that has generated track
734   //
735   
736   int current, parent;
737   //
738   parent=id;
739   while (1) {
740     current=parent;
741     parent=Particle(current)->GetFirstMother();
742     if(parent<0) return current;
743   }
744 }
745  
746 //_____________________________________________________________________________
747 void AliStack::DumpPart (Int_t i) const
748 {
749   //
750   // Dumps particle i in the stack
751   //
752   GetParticleMapEntry(i)->Print();
753 }
754
755 //_____________________________________________________________________________
756 void AliStack::DumpPStack ()
757 {
758   //
759   // Dumps the particle stack
760   //
761
762   Int_t i;
763
764   printf("\n\n=======================================================================\n");
765   for (i=0;i<fNtrack;i++) 
766     {
767       TParticle* particle = Particle(i);
768       if (particle) {
769         printf("-> %d ",i); particle->Print();
770         printf("--------------------------------------------------------------\n");
771       }
772       else 
773         Warning("DumpPStack", "No particle with id %d.", i); 
774     }    
775
776   printf("\n=======================================================================\n\n");
777   
778   // print  particle file map
779   // printf("\nParticle file map: \n");
780   // for (i=0; i<fNtrack; i++) 
781   //     printf("   %d th entry: %d \n",i,fParticleFileMap[i]);
782 }
783
784
785 //_____________________________________________________________________________
786 void AliStack::DumpLoadedStack() const
787 {
788   //
789   // Dumps the particle in the stack
790   // that are loaded in memory.
791   //
792
793   printf(
794          "\n\n=======================================================================\n");
795   for (Int_t i=0;i<fNtrack;i++) 
796     {
797       TParticle* particle = GetParticleMapEntry(i);
798       if (particle) {
799         printf("-> %d ",i); particle->Print();
800         printf("--------------------------------------------------------------\n");
801       }
802       else {    
803         printf("-> %d  Particle not loaded.\n",i);
804         printf("--------------------------------------------------------------\n");
805       } 
806     }
807   printf(
808          "\n=======================================================================\n\n");
809 }
810
811 //_____________________________________________________________________________
812 void  AliStack::SetCurrentTrack(Int_t track)
813
814   fCurrent = track; 
815   if (fCurrent < fNprimary) fCurrentTrack = Particle(track);
816 }
817
818
819 //_____________________________________________________________________________
820 //
821 // protected methods
822 //
823
824 //_____________________________________________________________________________
825 void AliStack::CleanParents()
826 {
827   //
828   // Clean particles stack
829   // Set parent/daughter relations
830   //
831   
832   TParticle *part;
833   int i;
834   for(i=0; i<fHgwmk+1; i++) {
835     part = GetParticleMapEntry(i);
836     if(part) if(!part->TestBit(kDaughtersBit)) {
837       part->SetFirstDaughter(-1);
838       part->SetLastDaughter(-1);
839     }
840   }
841 }
842
843 //_____________________________________________________________________________
844 TParticle* AliStack::GetNextParticle()
845 {
846   //
847   // Return next particle from stack of particles
848   //
849   
850   TParticle* particle = 0;
851   
852   // search secondaries
853   //for(Int_t i=fNtrack-1; i>=0; i--) {
854   for(Int_t i=fNtrack-1; i>fHgwmk; i--) {
855       particle = GetParticleMapEntry(i);
856       if ((particle) && (!particle->TestBit(kDoneBit))) {
857           fCurrent=i;    
858           return particle;
859       }   
860   }    
861
862   // take next primary if all secondaries were done
863   while (fCurrentPrimary>=0) {
864       fCurrent = fCurrentPrimary;    
865       particle = GetParticleMapEntry(fCurrentPrimary--);
866       if ((particle) && (!particle->TestBit(kDoneBit))) {
867           return particle;
868       } 
869   }
870   
871   // nothing to be tracked
872   fCurrent = -1;
873  
874   
875   return particle;  
876 }
877 //__________________________________________________________________________________________
878
879 void AliStack::ConnectTree(TTree* tree)
880 {
881 //
882 //  Creates branch for writing particles
883 //
884
885   fTreeK = tree;
886     
887   AliDebug(1, "Connecting TreeK");
888   if (fTreeK == 0x0)
889    {
890     if (TreeK() == 0x0)
891      {
892       AliFatal("Parameter is NULL");//we don't like such a jokes
893       return;
894      }
895     return;//in this case TreeK() calls back this method (ConnectTree) 
896            //tree after setting fTreeK, the rest was already executed
897            //it is safe to return now
898    }
899
900  //  Create a branch for particles   
901   
902   AliDebug(2, Form("Tree name is %s",fTreeK->GetName()));
903    
904   if (fTreeK->GetDirectory())
905    {
906      AliDebug(2, Form("and dir is %s",fTreeK->GetDirectory()->GetName()));
907    }    
908   else
909     AliWarning("DIR IS NOT SET !!!");
910   
911   TBranch *branch=fTreeK->GetBranch("Particles");
912   if(branch == 0x0)
913    {
914     branch = fTreeK->Branch("Particles", "TParticle", &fParticleBuffer, 4000);
915     AliDebug(2, "Creating Branch in Tree");
916    }  
917   else
918    {
919     AliDebug(2, "Branch Found in Tree");
920     branch->SetAddress(&fParticleBuffer);
921    }
922   if (branch->GetDirectory())
923    {
924     AliDebug(1, Form("Branch Dir Name is %s",branch->GetDirectory()->GetName()));
925    } 
926   else
927     AliWarning("Branch Dir is NOT SET");
928 }
929
930 //_____________________________________________________________________________
931
932 Bool_t AliStack::GetEvent()
933 {
934 //
935 // Get new event from TreeK
936
937     // Reset/Create the particle stack
938     Int_t size = (Int_t)TreeK()->GetEntries();
939     ResetArrays(size);
940     return kTRUE;
941 }
942 //_____________________________________________________________________________
943
944 Bool_t AliStack::IsStable(Int_t pdg) const
945 {
946   //
947   // Decide whether particle (pdg) is stable
948   //
949   
950   
951   // All ions/nucleons are considered as stable
952   // Nuclear code is 10LZZZAAAI
953   if(pdg>1000000000)return kTRUE;
954
955   const Int_t kNstable = 15;
956   Int_t i;
957   
958   Int_t pdgStable[kNstable] = {
959     kGamma,             // Photon
960     kElectron,          // Electron
961     kMuonPlus,          // Muon 
962     kPiPlus,            // Pion
963     kKPlus,             // Kaon
964     kK0Short,           // K0s
965     kK0Long,            // K0l
966     kProton,            // Proton 
967     kNeutron,           // Neutron
968     kLambda0,           // Lambda_0
969     kSigmaMinus,        // Sigma Minus
970     kSigmaPlus,         // Sigma Plus
971     3312,               // Xsi Minus 
972     3322,               // Xsi 
973     3334,               // Omega
974   };
975     
976   Bool_t isStable = kFALSE;
977   for (i = 0; i < kNstable; i++) {
978     if (pdg == TMath::Abs(pdgStable[i])) {
979       isStable = kTRUE;
980       break;
981     }
982   }
983   
984   return isStable;
985 }
986
987 //_____________________________________________________________________________
988 Bool_t AliStack::IsPhysicalPrimary(Int_t index)
989 {
990     //
991     // Test if a particle is a physical primary according to the following definition:
992     // Particles produced in the collision including products of strong and
993     // electromagnetic decay and excluding feed-down from weak decays of strange
994     // particles.
995     //
996     TParticle* p = Particle(index);
997     Int_t ist = p->GetStatusCode();
998     
999     //
1000     // Initial state particle
1001     if (ist > 1) return kFALSE;
1002     
1003     Int_t pdg = TMath::Abs(p->GetPdgCode());
1004     
1005     if (!IsStable(pdg)) return kFALSE;
1006     
1007     if (index < GetNprimary()) {
1008 //
1009 // Particle produced by generator
1010         return kTRUE;
1011     } else {
1012 //
1013 // Particle produced during transport
1014 //
1015
1016         Int_t imo =  p->GetFirstMother();
1017         TParticle* pm  = Particle(imo);
1018         Int_t mpdg = TMath::Abs(pm->GetPdgCode());
1019 // Check if it comes from a pi0 decay
1020 //
1021 // What about the pi0 Dalitz ??
1022 //      if ((mpdg == kPi0) && (imo < GetNprimary())) return kTRUE; 
1023
1024 // Check if this is a heavy flavor decay product
1025         Int_t mfl  = Int_t (mpdg / TMath::Power(10, Int_t(TMath::Log10(mpdg))));
1026         //
1027         // Light hadron
1028         if (mfl < 4) return kFALSE;
1029         
1030         //
1031         // Heavy flavor hadron produced by generator
1032         if (imo <  GetNprimary()) {
1033             return kTRUE;
1034         }
1035         
1036         // To be sure that heavy flavor has not been produced in a secondary interaction
1037         // Loop back to the generated mother
1038         while (imo >=  GetNprimary()) {
1039             imo = pm->GetFirstMother();
1040             pm  =  Particle(imo);
1041         }
1042         mpdg = TMath::Abs(pm->GetPdgCode());
1043         mfl  = Int_t (mpdg / TMath::Power(10, Int_t(TMath::Log10(mpdg))));
1044
1045         if (mfl < 4) {
1046             return kFALSE;
1047         } else {
1048             return kTRUE;
1049         } 
1050     } // produced by generator ?
1051