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