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