1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
18 ///////////////////////////////////////////////////////////////////////////////
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 //
26 ///////////////////////////////////////////////////////////////////////////////
29 #include <TObjArray.h>
30 #include <TParticle.h>
31 #include <TParticlePDG.h>
35 #include "AliModule.h"
38 #include "AliRunLoader.h"
43 //_______________________________________________________________________
56 fEventFolderName(AliConfig::GetDefaultEventFolderName())
59 // Default constructor
63 //_______________________________________________________________________
64 AliStack::AliStack(Int_t size, const char* evfoldname):
65 fParticles(new TClonesArray("TParticle",1000)),
66 fParticleMap(new TObjArray(size)),
76 fEventFolderName(evfoldname)
83 //_______________________________________________________________________
84 AliStack::AliStack(const AliStack& st):
104 //_______________________________________________________________________
105 void AliStack::Copy(TObject&) const
107 Fatal("Copy","Not implemented!\n");
110 //_______________________________________________________________________
111 AliStack::~AliStack()
118 fParticles->Delete();
122 //PH??? delete fTreeK;
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)
135 // Load a track on the stack
137 // done 0 if the track has to be transported
139 // parent identifier of the parent track. -1 for a primary
141 // pmom momentum GeV/c
143 // polar polarisation
144 // tof time of flight in seconds
145 // mecha production mechanism
146 // ntr on output the number of the track stored
149 // const Float_t tlife=0;
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
158 TParticlePDG* pmc = TDatabasePDG::Instance()->GetParticle(pdg);
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]);
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);
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);
172 Warning("PushTrack", "Particle type %d not defined in PDG Database !\n", pdg);
173 Warning("PushTrack", "Particle skipped !\n");
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)
185 // Load a track on the stack
187 // done 0 if the track has to be transported
189 // parent identifier of the parent track. -1 for a primary
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
199 // New method interface:
200 // arguments were changed to be in correspondence with TParticle
202 // Note: the energy is not calculated from the static mass but
203 // it is passed by argument e.
206 const Int_t kFirstDaughter=-1;
207 const Int_t kLastDaughter=-1;
210 TClonesArray &particles = *fParticles;
212 = new(particles[fLoadPoint++])
213 TParticle(pdg, is, parent, -1, kFirstDaughter, kLastDaughter,
214 px, py, pz, e, vx, vy, vz, tof);
216 particle->SetPolarisation(polx, poly, polz);
217 particle->SetWeight(weight);
218 particle->SetUniqueID(mech);
220 if(!done) particle->SetBit(kDoneBit);
222 // Declare that the daughter information is valid
223 particle->SetBit(kDaughtersBit);
224 // Add the particle to the stack
227 fParticleMap->AddAtAndExpand(particle, fNtrack);//CHECK!!
230 particle = dynamic_cast<TParticle*>(fParticleMap->At(parent));
232 particle->SetLastDaughter(fNtrack);
233 if(particle->GetFirstDaughter()<0) particle->SetFirstDaughter(fNtrack);
236 printf("Error in AliStack::PushTrack: Parent %d does not exist\n",parent);
241 // This is a primary track. Set high water mark for this event
244 // Set also number if primary tracks
245 fNprimary = fHgwmk+1;
251 //_____________________________________________________________________________
252 TParticle* AliStack::PopNextTrack(Int_t& itrack)
255 // Returns next track from stack of particles
259 TParticle* track = GetNextParticle();
263 track->SetBit(kDoneBit);
268 fCurrentTrack = track;
272 //_____________________________________________________________________________
273 TParticle* AliStack::PopPrimaryForTracking(Int_t i)
276 // Returns i-th primary particle if it is flagged to be tracked,
280 TParticle* particle = Particle(i);
282 if (!particle->TestBit(kDoneBit))
288 //_____________________________________________________________________________
289 void AliStack::PurifyKine()
292 // Compress kinematic tree keeping only flagged particles
293 // and renaming the particle id's in all the hits
296 TObjArray &particles = *fParticleMap;
297 int nkeep=fHgwmk+1, parent, i;
298 TParticle *part, *father;
299 TArrayI map(particles.GetLast()+1);
301 // Save in Header total number of tracks before compression
302 // If no tracks generated return now
303 if(fHgwmk+1 == fNtrack) return;
305 // First pass, invalid Daughter information
307 for(i=0; i<fNtrack; i++) {
308 // Preset map, to be removed later
309 if(i<=fHgwmk) map[i]=i ;
312 if((part=dynamic_cast<TParticle*>(particles.At(i)))) {
314 // Check of this track should be kept for physics reasons
315 if (KeepPhysics(part)) KeepTrack(i);
317 part->ResetBit(kDaughtersBit);
318 part->SetFirstDaughter(-1);
319 part->SetLastDaughter(-1);
323 // Invalid daughter information for the parent of the first particle
324 // generated. This may or may not be the current primary according to
325 // whether decays have been recorded among the primaries
326 part = dynamic_cast<TParticle*>(particles.At(fHgwmk+1));
327 particles.At(part->GetFirstMother())->ResetBit(kDaughtersBit);
328 // Second pass, build map between old and new numbering
329 for(i=fHgwmk+1; i<fNtrack; i++) {
330 if(particles.At(i)->TestBit(kKeepBit)) {
332 // This particle has to be kept
334 // If old and new are different, have to move the pointer
335 if(i!=nkeep) particles[nkeep]=particles.At(i);
336 part = dynamic_cast<TParticle*>(particles.At(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(map[parent]==-99) Fatal("PurifyKine","map[%d] = -99!\n",parent);
342 else part->SetFirstMother(map[parent]);
348 // Fix daughters information
349 for (i=fHgwmk+1; i<nkeep; i++) {
350 part = dynamic_cast<TParticle*>(particles.At(i));
351 parent = part->GetFirstMother();
353 father = dynamic_cast<TParticle*>(particles.At(parent));
354 if(father->TestBit(kDaughtersBit)) {
356 if(i<father->GetFirstDaughter()) father->SetFirstDaughter(i);
357 if(i>father->GetLastDaughter()) father->SetLastDaughter(i);
359 // Initialise daughters info for first pass
360 father->SetFirstDaughter(i);
361 father->SetLastDaughter(i);
362 father->SetBit(kDaughtersBit);
367 // Now loop on all registered hit lists
368 TList* hitLists = gAlice->GetMCApp()->GetHitLists();
369 TIter next(hitLists);
370 TCollection *hitList;
372 while((hitList = dynamic_cast<TCollection*>(next()))) {
373 TIter nexthit(hitList);
375 while((hit = dynamic_cast<AliHit*>(nexthit()))) {
376 hit->SetTrack(map[hit->GetTrack()]);
381 // This for detectors which have a special mapping mechanism
382 // for hits, such as TPC and TRD
385 TObjArray* modules = gAlice->Modules();
386 TIter nextmod(modules);
388 while((detector = dynamic_cast<AliModule*>(nextmod()))) {
389 detector->RemapTrackHitIDs(map.GetArray());
390 detector->RemapTrackReferencesIDs(map.GetArray());
393 gAlice->GetMCApp()->RemapTrackReferencesIDs(map.GetArray());
395 // Now the output bit, from fHgwmk to nkeep we write everything and we erase
396 if(nkeep>fParticleFileMap.GetSize()) fParticleFileMap.Set(Int_t (nkeep*1.5));
398 for (i=fHgwmk+1; i<nkeep; ++i) {
399 fParticleBuffer = dynamic_cast<TParticle*>(particles.At(i));
400 fParticleFileMap[i]=static_cast<Int_t>(TreeK()->GetEntries());
402 particles[i]=fParticleBuffer=0;
405 for (i=nkeep; i<fNtrack; ++i) particles[i]=0;
407 Int_t toshrink = fNtrack-fHgwmk-1;
408 fLoadPoint-=toshrink;
411 for(i=fLoadPoint; i<fLoadPoint+toshrink; ++i) fParticles->RemoveAt(i);
417 void AliStack::ReorderKine()
420 // In some transport code children might not come in a continuous sequence.
421 // In this case the stack has to be reordered in order to establish the
422 // mother daughter relation using index ranges.
424 if(fHgwmk+1 == fNtrack) return;
427 // Howmany secondaries have been produced ?
428 Int_t nNew = fNtrack - fHgwmk - 1;
432 TObjArray &particles = *fParticleMap;
435 // Copy pointers to temporary array
436 TParticle** tmp = new TParticle*[nNew];
438 for (i = 0; i < nNew; i++) {
439 if (particles.At(fHgwmk + 1 + i)) {
440 // tmp[i] = (TParticle*) (particles.At(fHgwmk + 1 + i))->Clone();
441 tmp[i] = (TParticle*) (particles.At(fHgwmk + 1 + i));
443 // if (((TParticle*) (particles.At(fHgwmk + 1 + i)))->TestBit(kKeepBit))
444 // tmp[i]->SetBit(kKeepBit);
445 // if (((TParticle*) (particles.At(fHgwmk + 1 + i)))->TestBit(kDoneBit))
446 // tmp[i]->SetBit(kDoneBit);
457 fLoadPoint = fHgwmk + 1;
459 // Re-Push particles into stack
460 // The outer loop is over parents, the inner over children.
461 // -1 refers to the primary particle
463 for (i = -1; i < nNew - 1; i++) {
467 ipa = tmp[0]->GetFirstMother();
468 parP =dynamic_cast<TParticle*>(particles.At(ipa));
470 ipa = (fHgwmk + 1 + i);
471 // Skip deleted particles
472 if (!tmp[i]) continue;
473 // Skip particles without children
474 if (tmp[i]->GetFirstDaughter() == -1) continue;
477 // Reset daughter information
478 parP->SetFirstDaughter(-1);
479 parP->SetLastDaughter(-1);
480 for (j = i + 1; j < nNew; j++) {
481 // Skip deleted particles
482 if (!tmp[j]) continue;
483 // Skip particles already handled
484 if (map1[j] != -99) continue;
485 Int_t jpa = tmp[j]->GetFirstMother();
486 // Check if daughter of current parent
488 particles[fLoadPoint] = tmp[j];
489 // Re-establish daughter information
490 parP->SetLastDaughter(fLoadPoint);
491 if (parP->GetFirstDaughter() == -1) parP->SetFirstDaughter(fLoadPoint);
492 // Set Mother information
494 tmp[j]->SetFirstMother(map1[i]);
497 map1[j] = fLoadPoint;
498 // Increase load point
507 // Build map for remapping of hits
509 TArrayI map(fNtrack);
510 for (i = 0; i < fNtrack; i ++) {
514 map[i] = map1[i - fHgwmk -1];
518 // Now loop on all registered hit lists
520 TList* hitLists = gAlice->GetMCApp()->GetHitLists();
521 TIter next(hitLists);
522 TCollection *hitList;
524 while((hitList = dynamic_cast<TCollection*>(next()))) {
525 TIter nexthit(hitList);
527 while((hit = dynamic_cast<AliHit*>(nexthit()))) {
528 hit->SetTrack(map[hit->GetTrack()]);
533 // This for detectors which have a special mapping mechanism
534 // for hits, such as TPC and TRD
537 TObjArray* modules = gAlice->Modules();
538 TIter nextmod(modules);
540 while((detector = dynamic_cast<AliModule*>(nextmod()))) {
541 detector->RemapTrackHitIDs(map.GetArray());
542 detector->RemapTrackReferencesIDs(map.GetArray());
544 gAlice->GetMCApp()->RemapTrackReferencesIDs(map.GetArray());
545 } // new particles poduced
548 Bool_t AliStack::KeepPhysics(TParticle* part)
551 // Some particles have to kept on the stack for reasons motivated
552 // by physics analysis. Decision is put here.
554 Bool_t keep = kFALSE;
556 // Keep first-generation daughter from primaries with heavy flavor
558 Int_t parent = part->GetFirstMother();
559 if (parent >= 0 && parent <= fHgwmk) {
560 TParticle* father = dynamic_cast<TParticle*>(Particles()->At(parent));
561 Int_t kf = father->GetPdgCode();
565 if (kfl > 10) kfl/=100;
567 if (kfl > 10) kfl/=10;
568 if (kfl > 10) kfl/=10;
576 //_____________________________________________________________________________
577 void AliStack::FinishEvent()
580 // Write out the kinematics that was not yet filled
583 // Update event header
587 // Fatal("FinishEvent", "No kinematics tree is defined.");
588 // Don't panic this is a probably a lego run
593 if(TreeK()->GetEntries() ==0) {
594 // set the fParticleFileMap size for the first time
595 fParticleFileMap.Set(fHgwmk+1);
598 Bool_t allFilled = kFALSE;
600 for(Int_t i=0; i<fHgwmk+1; ++i)
601 if((part=fParticleMap->At(i))) {
602 fParticleBuffer = dynamic_cast<TParticle*>(part);
603 fParticleFileMap[i]= static_cast<Int_t>(TreeK()->GetEntries());
605 //PH (*fParticleMap)[i]=fParticleBuffer=0;
607 fParticleMap->AddAt(0,i);
609 // When all primaries were filled no particle!=0
610 // should be left => to be removed later.
611 if (allFilled) printf("Why != 0 part # %d?\n",i);
615 // // printf("Why = 0 part # %d?\n",i); => We know.
617 // we don't break now in order to be sure there is no
618 // particle !=0 left.
619 // To be removed later and replaced with break.
620 if(!allFilled) allFilled = kTRUE;
623 //_____________________________________________________________________________
625 void AliStack::FlagTrack(Int_t track)
628 // Flags a track and all its family tree to be kept
635 particle=dynamic_cast<TParticle*>(fParticleMap->At(curr));
637 // If the particle is flagged the three from here upward is saved already
638 if(particle->TestBit(kKeepBit)) return;
640 // Save this particle
641 particle->SetBit(kKeepBit);
643 // Move to father if any
644 if((curr=particle->GetFirstMother())==-1) return;
648 //_____________________________________________________________________________
649 void AliStack::KeepTrack(Int_t track)
652 // Flags a track to be kept
655 fParticleMap->At(track)->SetBit(kKeepBit);
658 //_____________________________________________________________________________
659 void AliStack::Reset(Int_t size)
674 //_____________________________________________________________________________
675 void AliStack::ResetArrays(Int_t size)
678 // Resets stack arrays
684 fParticles = new TClonesArray("TParticle",1000);
686 fParticleMap->Clear();
687 if (size>0) fParticleMap->Expand(size);}
689 fParticleMap = new TObjArray(size);
692 //_____________________________________________________________________________
693 void AliStack::SetHighWaterMark(Int_t)
696 // Set high water mark for last track in event
701 fCurrentPrimary=fHgwmk;
702 // Set also number of primary tracks
703 fNprimary = fHgwmk+1;
707 //_____________________________________________________________________________
708 TParticle* AliStack::Particle(Int_t i)
711 // Return particle with specified ID
713 //PH if(!(*fParticleMap)[i]) {
714 if(!fParticleMap->At(i)) {
715 Int_t nentries = fParticles->GetEntriesFast();
716 // algorithmic way of getting entry index
717 // (primary particles are filled after secondaries)
718 Int_t entry = TreeKEntry(i);
719 // check whether algorithmic way and
720 // and the fParticleFileMap[i] give the same;
721 // give the fatal error if not
722 if (entry != fParticleFileMap[i]) {
724 "!! The algorithmic way and map are different: !!\n entry: %d map: %d",
725 entry, fParticleFileMap[i]);
728 TreeK()->GetEntry(entry);
729 new ((*fParticles)[nentries]) TParticle(*fParticleBuffer);
730 fParticleMap->AddAt((*fParticles)[nentries],i);
732 //PH return dynamic_cast<TParticle *>((*fParticleMap)[i]);
733 return dynamic_cast<TParticle*>(fParticleMap->At(i));
736 //_____________________________________________________________________________
737 TParticle* AliStack::ParticleFromTreeK(Int_t id) const
740 // return pointer to TParticle with label id
743 if ((entry = TreeKEntry(id)) < 0) return 0;
744 if (fTreeK->GetEntry(entry)<=0) return 0;
745 return fParticleBuffer;
748 //_____________________________________________________________________________
749 Int_t AliStack::TreeKEntry(Int_t id) const
752 // return entry number in the TreeK for particle with label id
753 // return negative number if label>fNtrack
757 entry = id+fNtrack-fNprimary;
759 entry = id-fNprimary;
763 //_____________________________________________________________________________
764 Int_t AliStack::GetCurrentParentTrackNumber() const
767 // Return number of the parent of the current track
770 TParticle* current = (TParticle*)fParticleMap->At(fCurrent);
773 return current->GetFirstMother();
775 Warning("GetCurrentParentTrackNumber", "Current track not found in the stack");
780 //_____________________________________________________________________________
781 Int_t AliStack::GetPrimary(Int_t id)
784 // Return number of primary that has generated track
792 parent=Particle(current)->GetFirstMother();
793 if(parent<0) return current;
797 //_____________________________________________________________________________
798 void AliStack::DumpPart (Int_t i) const
801 // Dumps particle i in the stack
804 //PH dynamic_cast<TParticle*>((*fParticleMap)[i])->Print();
805 dynamic_cast<TParticle*>(fParticleMap->At(i))->Print();
808 //_____________________________________________________________________________
809 void AliStack::DumpPStack ()
812 // Dumps the particle stack
817 printf("\n\n=======================================================================\n");
818 for (i=0;i<fNtrack;i++)
820 TParticle* particle = Particle(i);
822 printf("-> %d ",i); particle->Print();
823 printf("--------------------------------------------------------------\n");
826 Warning("DumpPStack", "No particle with id %d.", i);
829 printf("\n=======================================================================\n\n");
831 // print particle file map
832 printf("\nParticle file map: \n");
833 for (i=0; i<fNtrack; i++)
834 printf(" %d th entry: %d \n",i,fParticleFileMap[i]);
838 //_____________________________________________________________________________
839 void AliStack::DumpLoadedStack() const
842 // Dumps the particle in the stack
843 // that are loaded in memory.
846 TObjArray &particles = *fParticleMap;
848 "\n\n=======================================================================\n");
849 for (Int_t i=0;i<fNtrack;i++)
851 TParticle* particle = dynamic_cast<TParticle*>(particles[i]);
853 printf("-> %d ",i); particle->Print();
854 printf("--------------------------------------------------------------\n");
857 printf("-> %d Particle not loaded.\n",i);
858 printf("--------------------------------------------------------------\n");
862 "\n=======================================================================\n\n");
869 //_____________________________________________________________________________
870 void AliStack::CleanParents()
873 // Clean particles stack
874 // Set parent/daughter relations
877 TObjArray &particles = *fParticleMap;
880 for(i=0; i<fHgwmk+1; i++) {
881 part = dynamic_cast<TParticle*>(particles.At(i));
882 if(part) if(!part->TestBit(kDaughtersBit)) {
883 part->SetFirstDaughter(-1);
884 part->SetLastDaughter(-1);
889 //_____________________________________________________________________________
890 TParticle* AliStack::GetNextParticle()
893 // Return next particle from stack of particles
896 TParticle* particle = 0;
898 // search secondaries
899 //for(Int_t i=fNtrack-1; i>=0; i--) {
900 for(Int_t i=fNtrack-1; i>fHgwmk; i--) {
901 particle = dynamic_cast<TParticle*>(fParticleMap->At(i));
902 if ((particle) && (!particle->TestBit(kDoneBit))) {
908 // take next primary if all secondaries were done
909 while (fCurrentPrimary>=0) {
910 fCurrent = fCurrentPrimary;
911 particle = dynamic_cast<TParticle*>(fParticleMap->At(fCurrentPrimary--));
912 if ((particle) && (!particle->TestBit(kDoneBit))) {
917 // nothing to be tracked
923 //__________________________________________________________________________________________
925 TTree* AliStack::TreeK()
934 AliRunLoader *rl = AliRunLoader::GetRunLoader(fEventFolderName);
937 Fatal("TreeK","Can not get RunLoader from event folder named %s",fEventFolderName.Data());
938 return 0x0;//pro forma
940 fTreeK = rl->TreeK();
947 //don't panic - could be Lego
948 if (AliLoader::GetDebug())
950 Warning("TreeK","Can not get TreeK from RL. Ev. Folder is %s",fEventFolderName.Data());
954 return fTreeK;//never reached
956 //__________________________________________________________________________________________
958 void AliStack::ConnectTree()
961 // Creates branch for writing particles
963 if (AliLoader::GetDebug()) Info("ConnectTree","Connecting TreeK");
968 Fatal("ConnectTree","Parameter is NULL");//we don't like such a jokes
971 return;//in this case TreeK() calls back this method (ConnectTree)
972 //tree after setting fTreeK, the rest was already executed
973 //it is safe to return now
976 // Create a branch for particles
978 if (AliLoader::GetDebug())
979 Info("ConnectTree","Tree name is %s",fTreeK->GetName());
981 if (fTreeK->GetDirectory())
983 if (AliLoader::GetDebug())
984 Info("ConnectTree","and dir is %s",fTreeK->GetDirectory()->GetName());
987 Warning("ConnectTree","DIR IS NOT SET !!!");
989 TBranch *branch=fTreeK->GetBranch(AliRunLoader::GetKineBranchName());
992 branch = fTreeK->Branch(AliRunLoader::GetKineBranchName(), "TParticle", &fParticleBuffer, 4000);
993 if (AliLoader::GetDebug()) Info("ConnectTree","Creating Branch in Tree");
997 if (AliLoader::GetDebug()) Info("ConnectTree","Branch Found in Tree");
998 branch->SetAddress(&fParticleBuffer);
1000 if (branch->GetDirectory())
1002 if (AliLoader::GetDebug())
1003 Info("ConnectTree","Branch Dir Name is %s",branch->GetDirectory()->GetName());
1006 Warning("ConnectTree","Branch Dir is NOT SET");
1008 //__________________________________________________________________________________________
1011 void AliStack::BeginEvent()
1013 // start a new event
1017 //_____________________________________________________________________________
1018 void AliStack::FinishRun()
1020 // Clean TreeK information
1022 //_____________________________________________________________________________
1024 Bool_t AliStack::GetEvent()
1027 // Get new event from TreeK
1029 // Reset/Create the particle stack
1032 if (TreeK() == 0x0) //forces connecting
1034 Error("GetEvent","cannot find Kine Tree for current event\n");
1038 Int_t size = (Int_t)TreeK()->GetEntries();
1042 //_____________________________________________________________________________
1044 void AliStack::SetEventFolderName(const char* foldname)
1046 //Sets event folder name
1047 fEventFolderName = foldname;