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
22 ///////////////////////////////////////////////////////////////////////////////
27 #include <TObjArray.h>
28 #include <TParticle.h>
29 #include <TParticlePDG.h>
33 #include "AliModule.h"
36 #include "AliRunLoader.h"
41 //_______________________________________________________________________
54 fEventFolderName(AliConfig::fgkDefaultEventFolderName)
57 // Default constructor
61 //_______________________________________________________________________
62 AliStack::AliStack(Int_t size, const char* evfoldname):
63 fParticles(new TClonesArray("TParticle",1000)),
64 fParticleMap(new TObjArray(size)),
74 fEventFolderName(evfoldname)
81 //_______________________________________________________________________
82 AliStack::AliStack(const AliStack& st):
102 //_______________________________________________________________________
103 void AliStack::Copy(TObject&) const
105 Fatal("Copy","Not implemented!\n");
108 //_______________________________________________________________________
109 AliStack::~AliStack()
116 fParticles->Delete();
120 //PH??? delete fTreeK;
127 //_____________________________________________________________________________
128 void AliStack::PushTrack(Int_t done, Int_t parent, Int_t pdg, Float_t *pmom,
129 Float_t *vpos, Float_t *polar, Float_t tof,
130 TMCProcess mech, Int_t &ntr, Float_t weight, Int_t is)
133 // Load a track on the stack
135 // done 0 if the track has to be transported
137 // parent identifier of the parent track. -1 for a primary
139 // pmom momentum GeV/c
141 // polar polarisation
142 // tof time of flight in seconds
143 // mecha production mechanism
144 // ntr on output the number of the track stored
147 // const Float_t tlife=0;
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
156 TParticlePDG* pmc = TDatabasePDG::Instance()->GetParticle(pdg);
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]);
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);
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);
170 Warning("PushTrack", "Particle type %d not defined in PDG Database !\n", pdg);
171 Warning("PushTrack", "Particle skipped !\n");
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)
183 // Load a track on the stack
185 // done 0 if the track has to be transported
187 // parent identifier of the parent track. -1 for a primary
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
197 // New method interface:
198 // arguments were changed to be in correspondence with TParticle
200 // Note: the energy is not calculated from the static mass but
201 // it is passed by argument e.
204 const Int_t kFirstDaughter=-1;
205 const Int_t kLastDaughter=-1;
207 TClonesArray &particles = *fParticles;
209 = new(particles[fLoadPoint++])
210 TParticle(pdg, is, parent, -1, kFirstDaughter, kLastDaughter,
211 px, py, pz, e, vx, vy, vz, tof);
213 particle->SetPolarisation(polx, poly, polz);
214 particle->SetWeight(weight);
215 particle->SetUniqueID(mech);
217 if(!done) particle->SetBit(kDoneBit);
219 // Declare that the daughter information is valid
220 particle->SetBit(kDaughtersBit);
221 // Add the particle to the stack
222 fParticleMap->AddAtAndExpand(particle, fNtrack);//CHECK!!
225 particle = dynamic_cast<TParticle*>(fParticleMap->At(parent));
227 particle->SetLastDaughter(fNtrack);
228 if(particle->GetFirstDaughter()<0) particle->SetFirstDaughter(fNtrack);
231 printf("Error in AliStack::PushTrack: Parent %d does not exist\n",parent);
236 // This is a primary track. Set high water mark for this event
239 // Set also number if primary tracks
240 fNprimary = fHgwmk+1;
246 //_____________________________________________________________________________
247 TParticle* AliStack::PopNextTrack(Int_t& itrack)
250 // Returns next track from stack of particles
254 TParticle* track = GetNextParticle();
258 track->SetBit(kDoneBit);
263 fCurrentTrack = track;
267 //_____________________________________________________________________________
268 TParticle* AliStack::PopPrimaryForTracking(Int_t i)
271 // Returns i-th primary particle if it is flagged to be tracked,
275 TParticle* particle = Particle(i);
277 if (!particle->TestBit(kDoneBit))
283 //_____________________________________________________________________________
284 void AliStack::PurifyKine()
287 // Compress kinematic tree keeping only flagged particles
288 // and renaming the particle id's in all the hits
291 TObjArray &particles = *fParticleMap;
292 int nkeep=fHgwmk+1, parent, i;
293 TParticle *part, *father;
294 TArrayI map(particles.GetLast()+1);
296 // Save in Header total number of tracks before compression
298 // If no tracks generated return now
299 if(fHgwmk+1 == fNtrack) return;
301 // First pass, invalid Daughter information
302 for(i=0; i<fNtrack; i++) {
303 // Preset map, to be removed later
304 if(i<=fHgwmk) map[i]=i ;
307 if((part=dynamic_cast<TParticle*>(particles.At(i)))) {
309 // Check of this track should be kept for physics reasons
310 if (KeepPhysics(part)) KeepTrack(i);
312 part->ResetBit(kDaughtersBit);
313 part->SetFirstDaughter(-1);
314 part->SetLastDaughter(-1);
318 // Invalid daughter information for the parent of the first particle
319 // generated. This may or may not be the current primary according to
320 // whether decays have been recorded among the primaries
321 part = dynamic_cast<TParticle*>(particles.At(fHgwmk+1));
322 particles.At(part->GetFirstMother())->ResetBit(kDaughtersBit);
323 // Second pass, build map between old and new numbering
324 for(i=fHgwmk+1; i<fNtrack; i++) {
325 if(particles.At(i)->TestBit(kKeepBit)) {
327 // This particle has to be kept
329 // If old and new are different, have to move the pointer
330 if(i!=nkeep) particles[nkeep]=particles.At(i);
331 part = dynamic_cast<TParticle*>(particles.At(nkeep));
333 // as the parent is always *before*, it must be already
334 // in place. This is what we are checking anyway!
335 if((parent=part->GetFirstMother())>fHgwmk)
336 if(map[parent]==-99) Fatal("PurifyKine","map[%d] = -99!\n",parent);
337 else part->SetFirstMother(map[parent]);
343 // Fix daughters information
344 for (i=fHgwmk+1; i<nkeep; i++) {
345 part = dynamic_cast<TParticle*>(particles.At(i));
346 parent = part->GetFirstMother();
348 father = dynamic_cast<TParticle*>(particles.At(parent));
349 if(father->TestBit(kDaughtersBit)) {
351 if(i<father->GetFirstDaughter()) father->SetFirstDaughter(i);
352 if(i>father->GetLastDaughter()) father->SetLastDaughter(i);
354 // Initialise daughters info for first pass
355 father->SetFirstDaughter(i);
356 father->SetLastDaughter(i);
357 father->SetBit(kDaughtersBit);
362 // Now loop on all registered hit lists
363 TList* hitLists = gAlice->GetMCApp()->GetHitLists();
364 TIter next(hitLists);
365 TCollection *hitList;
366 while((hitList = dynamic_cast<TCollection*>(next()))) {
367 TIter nexthit(hitList);
369 while((hit = dynamic_cast<AliHit*>(nexthit()))) {
370 hit->SetTrack(map[hit->GetTrack()]);
375 // This for detectors which have a special mapping mechanism
376 // for hits, such as TPC and TRD
379 TObjArray* modules = gAlice->Modules();
380 TIter nextmod(modules);
382 while((detector = dynamic_cast<AliModule*>(nextmod()))) {
383 detector->RemapTrackHitIDs(map.GetArray());
384 detector->RemapTrackReferencesIDs(map.GetArray());
387 gAlice->GetMCApp()->RemapTrackReferencesIDs(map.GetArray());
389 // Now the output bit, from fHgwmk to nkeep we write everything and we erase
390 if(nkeep>fParticleFileMap.GetSize()) fParticleFileMap.Set(Int_t (nkeep*1.5));
392 for (i=fHgwmk+1; i<nkeep; ++i) {
393 fParticleBuffer = dynamic_cast<TParticle*>(particles.At(i));
394 fParticleFileMap[i]=static_cast<Int_t>(TreeK()->GetEntries());
396 particles[i]=fParticleBuffer=0;
399 for (i=nkeep; i<fNtrack; ++i) particles[i]=0;
401 Int_t toshrink = fNtrack-fHgwmk-1;
402 fLoadPoint-=toshrink;
403 for(i=fLoadPoint; i<fLoadPoint+toshrink; ++i) fParticles->RemoveAt(i);
410 Bool_t AliStack::KeepPhysics(TParticle* part)
413 // Some particles have to kept on the stack for reasons motivated
414 // by physics analysis. Decision is put here.
416 Bool_t keep = kFALSE;
418 // Keep first-generation daughter from primaries with heavy flavor
420 Int_t parent = part->GetFirstMother();
421 if (parent >= 0 && parent <= fHgwmk) {
422 TParticle* father = dynamic_cast<TParticle*>(Particles()->At(parent));
423 Int_t kf = father->GetPdgCode();
427 if (kfl > 10) kfl/=100;
429 if (kfl > 10) kfl/=10;
430 if (kfl > 10) kfl/=10;
438 //_____________________________________________________________________________
439 void AliStack::FinishEvent()
442 // Write out the kinematics that was not yet filled
445 // Update event header
449 // Fatal("FinishEvent", "No kinematics tree is defined.");
450 // Don't panic this is a probably a lego run
455 if(TreeK()->GetEntries() ==0) {
456 // set the fParticleFileMap size for the first time
457 fParticleFileMap.Set(fHgwmk+1);
460 Bool_t allFilled = kFALSE;
462 for(Int_t i=0; i<fHgwmk+1; ++i)
463 if((part=fParticleMap->At(i))) {
464 fParticleBuffer = dynamic_cast<TParticle*>(part);
465 fParticleFileMap[i]= static_cast<Int_t>(TreeK()->GetEntries());
467 //PH (*fParticleMap)[i]=fParticleBuffer=0;
469 fParticleMap->AddAt(0,i);
471 // When all primaries were filled no particle!=0
472 // should be left => to be removed later.
473 if (allFilled) printf("Why != 0 part # %d?\n",i);
477 // // printf("Why = 0 part # %d?\n",i); => We know.
479 // we don't break now in order to be sure there is no
480 // particle !=0 left.
481 // To be removed later and replaced with break.
482 if(!allFilled) allFilled = kTRUE;
485 //_____________________________________________________________________________
487 void AliStack::FlagTrack(Int_t track)
490 // Flags a track and all its family tree to be kept
497 particle=dynamic_cast<TParticle*>(fParticleMap->At(curr));
499 // If the particle is flagged the three from here upward is saved already
500 if(particle->TestBit(kKeepBit)) return;
502 // Save this particle
503 particle->SetBit(kKeepBit);
505 // Move to father if any
506 if((curr=particle->GetFirstMother())==-1) return;
510 //_____________________________________________________________________________
511 void AliStack::KeepTrack(Int_t track)
514 // Flags a track to be kept
517 fParticleMap->At(track)->SetBit(kKeepBit);
520 //_____________________________________________________________________________
521 void AliStack::Reset(Int_t size)
536 //_____________________________________________________________________________
537 void AliStack::ResetArrays(Int_t size)
540 // Resets stack arrays
546 fParticles = new TClonesArray("TParticle",1000);
548 fParticleMap->Clear();
549 if (size>0) fParticleMap->Expand(size);}
551 fParticleMap = new TObjArray(size);
554 //_____________________________________________________________________________
555 void AliStack::SetHighWaterMark(Int_t)
558 // Set high water mark for last track in event
562 fCurrentPrimary=fHgwmk;
564 // Set also number of primary tracks
565 fNprimary = fHgwmk+1;
569 //_____________________________________________________________________________
570 TParticle* AliStack::Particle(Int_t i)
573 // Return particle with specified ID
575 //PH if(!(*fParticleMap)[i]) {
576 if(!fParticleMap->At(i)) {
577 Int_t nentries = fParticles->GetEntriesFast();
578 // algorithmic way of getting entry index
579 // (primary particles are filled after secondaries)
580 Int_t entry = TreeKEntry(i);
581 // check whether algorithmic way and
582 // and the fParticleFileMap[i] give the same;
583 // give the fatal error if not
584 if (entry != fParticleFileMap[i]) {
586 "!! The algorithmic way and map are different: !!\n entry: %d map: %d",
587 entry, fParticleFileMap[i]);
590 TreeK()->GetEntry(entry);
591 new ((*fParticles)[nentries]) TParticle(*fParticleBuffer);
592 fParticleMap->AddAt((*fParticles)[nentries],i);
594 //PH return dynamic_cast<TParticle *>((*fParticleMap)[i]);
595 return dynamic_cast<TParticle*>(fParticleMap->At(i));
598 //_____________________________________________________________________________
599 TParticle* AliStack::ParticleFromTreeK(Int_t id) const
602 // return pointer to TParticle with label id
605 if ((entry = TreeKEntry(id)) < 0) return 0;
606 if (fTreeK->GetEntry(entry)<=0) return 0;
607 return fParticleBuffer;
610 //_____________________________________________________________________________
611 Int_t AliStack::TreeKEntry(Int_t id) const
614 // return entry number in the TreeK for particle with label id
615 // return negative number if label>fNtrack
619 entry = id+fNtrack-fNprimary;
621 entry = id-fNprimary;
625 //_____________________________________________________________________________
626 Int_t AliStack::GetCurrentParentTrackNumber() const
629 // Return number of the parent of the current track
632 TParticle* current = (TParticle*)fParticleMap->At(fCurrent);
635 return current->GetFirstMother();
637 Warning("GetCurrentParentTrackNumber", "Current track not found in the stack");
642 //_____________________________________________________________________________
643 Int_t AliStack::GetPrimary(Int_t id)
646 // Return number of primary that has generated track
654 parent=Particle(current)->GetFirstMother();
655 if(parent<0) return current;
659 //_____________________________________________________________________________
660 void AliStack::DumpPart (Int_t i) const
663 // Dumps particle i in the stack
666 //PH dynamic_cast<TParticle*>((*fParticleMap)[i])->Print();
667 dynamic_cast<TParticle*>(fParticleMap->At(i))->Print();
670 //_____________________________________________________________________________
671 void AliStack::DumpPStack ()
674 // Dumps the particle stack
679 printf("\n\n=======================================================================\n");
680 for (i=0;i<fNtrack;i++)
682 TParticle* particle = Particle(i);
684 printf("-> %d ",i); particle->Print();
685 printf("--------------------------------------------------------------\n");
688 Warning("DumpPStack", "No particle with id %d.", i);
691 printf("\n=======================================================================\n\n");
693 // print particle file map
694 printf("\nParticle file map: \n");
695 for (i=0; i<fNtrack; i++)
696 printf(" %d th entry: %d \n",i,fParticleFileMap[i]);
700 //_____________________________________________________________________________
701 void AliStack::DumpLoadedStack() const
704 // Dumps the particle in the stack
705 // that are loaded in memory.
708 TObjArray &particles = *fParticleMap;
710 "\n\n=======================================================================\n");
711 for (Int_t i=0;i<fNtrack;i++)
713 TParticle* particle = dynamic_cast<TParticle*>(particles[i]);
715 printf("-> %d ",i); particle->Print();
716 printf("--------------------------------------------------------------\n");
719 printf("-> %d Particle not loaded.\n",i);
720 printf("--------------------------------------------------------------\n");
724 "\n=======================================================================\n\n");
731 //_____________________________________________________________________________
732 void AliStack::CleanParents()
735 // Clean particles stack
736 // Set parent/daughter relations
739 TObjArray &particles = *fParticleMap;
742 for(i=0; i<fHgwmk+1; i++) {
743 part = dynamic_cast<TParticle*>(particles.At(i));
744 if(part) if(!part->TestBit(kDaughtersBit)) {
745 part->SetFirstDaughter(-1);
746 part->SetLastDaughter(-1);
751 //_____________________________________________________________________________
752 TParticle* AliStack::GetNextParticle()
755 // Return next particle from stack of particles
758 TParticle* particle = 0;
760 // search secondaries
761 //for(Int_t i=fNtrack-1; i>=0; i--) {
762 for(Int_t i=fNtrack-1; i>fHgwmk; i--) {
763 particle = dynamic_cast<TParticle*>(fParticleMap->At(i));
764 if ((particle) && (!particle->TestBit(kDoneBit))) {
770 // take next primary if all secondaries were done
771 while (fCurrentPrimary>=0) {
772 fCurrent = fCurrentPrimary;
773 particle = dynamic_cast<TParticle*>(fParticleMap->At(fCurrentPrimary--));
774 if ((particle) && (!particle->TestBit(kDoneBit))) {
779 // nothing to be tracked
783 //__________________________________________________________________________________________
785 TTree* AliStack::TreeK()
794 AliRunLoader *rl = AliRunLoader::GetRunLoader(fEventFolderName);
797 Fatal("TreeK","Can not get RunLoader from event folder named %s",fEventFolderName.Data());
798 return 0x0;//pro forma
800 fTreeK = rl->TreeK();
807 //don't panic - could be Lego
808 if (AliLoader::GetDebug())
810 Warning("TreeK","Can not get TreeK from RL. Ev. Folder is %s",fEventFolderName.Data());
814 return fTreeK;//never reached
816 //__________________________________________________________________________________________
818 void AliStack::ConnectTree()
821 // Creates branch for writing particles
823 if (AliLoader::GetDebug()) Info("ConnectTree","Connecting TreeK");
828 Fatal("ConnectTree","Parameter is NULL");//we don't like such a jokes
831 return;//in this case TreeK() calls back this method (ConnectTree)
832 //tree after setting fTreeK, the rest was already executed
833 //it is safe to return now
836 // Create a branch for particles
838 if (AliLoader::GetDebug())
839 Info("ConnectTree","Tree name is %s",fTreeK->GetName());
841 if (fTreeK->GetDirectory())
843 if (AliLoader::GetDebug())
844 Info("ConnectTree","and dir is %s",fTreeK->GetDirectory()->GetName());
847 Warning("ConnectTree","DIR IS NOT SET !!!");
849 TBranch *branch=fTreeK->GetBranch(AliRunLoader::fgkKineBranchName);
852 branch = fTreeK->Branch(AliRunLoader::fgkKineBranchName, "TParticle", &fParticleBuffer, 4000);
853 if (AliLoader::GetDebug()) Info("ConnectTree","Creating Branch in Tree");
857 if (AliLoader::GetDebug()) Info("ConnectTree","Branch Found in Tree");
858 branch->SetAddress(&fParticleBuffer);
860 if (branch->GetDirectory())
862 if (AliLoader::GetDebug())
863 Info("ConnectTree","Branch Dir Name is %s",branch->GetDirectory()->GetName());
866 Warning("ConnectTree","Branch Dir is NOT SET");
868 //__________________________________________________________________________________________
871 void AliStack::BeginEvent()
877 //_____________________________________________________________________________
878 void AliStack::FinishRun()
880 // Clean TreeK information
882 //_____________________________________________________________________________
884 Bool_t AliStack::GetEvent()
887 // Get new event from TreeK
889 // Reset/Create the particle stack
892 if (TreeK() == 0x0) //forces connecting
894 Error("GetEvent","cannot find Kine Tree for current event\n");
898 Int_t size = (Int_t)TreeK()->GetEntries();
902 //_____________________________________________________________________________
904 void AliStack::SetEventFolderName(const char* foldname)
906 //Sets event folder name
907 fEventFolderName = foldname;