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