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 // This class is extracted from the AliRun class
19 // and contains all the MC-related functionality
20 // The number of dependencies has to be reduced...
21 // Author: F.Carminati
22 // Federico.Carminati@cern.ch
28 #include <TClonesArray.h>
30 #include <TGeoGlobalMagField.h>
31 #include <TGeoManager.h>
32 #include <TParticle.h>
34 #include <TStopwatch.h>
36 #include <TVirtualMC.h>
37 #include <TMCVerbose.h>
40 #include "AliCDBEntry.h"
41 #include "AliCDBManager.h"
42 #include "AliCDBStorage.h"
43 #include "AliDetector.h"
44 #include "AliGenerator.h"
45 #include "AliGeomManager.h"
46 #include "AliHeader.h"
53 #include "AliSimulation.h"
55 #include "AliTrackReference.h"
56 #include "AliTransportMonitor.h"
62 //_______________________________________________________________________
65 fSaveRndmStatus(kFALSE),
66 fSaveRndmEventStatus(kFALSE),
67 fReadRndmStatus(kFALSE),
68 fUseMonitoring(kFALSE),
69 fRndmFileName("random.root"),
92 //_______________________________________________________________________
93 AliMC::AliMC(const char *name, const char *title) :
94 TVirtualMCApplication(name, title),
96 fSaveRndmStatus(kFALSE),
97 fSaveRndmEventStatus(kFALSE),
98 fReadRndmStatus(kFALSE),
99 fUseMonitoring(kFALSE),
100 fRndmFileName("random.root"),
109 fImedia(new TArrayI(1000)),
112 fHitLists(new TList()),
115 fTrackReferences("AliTrackReference", 100),
116 fTmpTrackReferences("AliTrackReference", 100)
119 // Set transport parameters
122 // Prepare the tracking medium lists
123 for(Int_t i=0;i<1000;i++) (*fImedia)[i]=-99;
126 //_______________________________________________________________________
134 // Delete track references
137 //_______________________________________________________________________
138 void AliMC::ConstructGeometry()
141 // Either load geometry from file or create it through usual
142 // loop on detectors. In the first case the method
143 // AliModule::CreateMaterials() only builds fIdtmed and is postponed
144 // at InitGeometry().
147 if(AliSimulation::Instance()->IsGeometryFromFile()){ //load geometry either from CDB or from file
148 if(IsGeometryFromCDB()){
149 AliInfo("Loading geometry from CDB default storage");
150 AliCDBPath path("GRP","Geometry","Data");
151 AliCDBEntry *entry=AliCDBManager::Instance()->Get(path.GetPath());
152 if(!entry) AliFatal("Unable to load geometry from CDB!");
154 gGeoManager = (TGeoManager*) entry->GetObject();
155 if (!gGeoManager) AliFatal("TGeoManager object not found in the specified CDB entry!");
158 const char *geomfilename = AliSimulation::Instance()->GetGeometryFile();
159 if(gSystem->ExpandPathName(geomfilename)){
160 AliInfo(Form("Loading geometry from file:\n %40s",geomfilename));
161 TGeoManager::Import(geomfilename);
163 AliInfo(Form("Geometry file %40s not found!\n",geomfilename));
167 TVirtualMC::GetMC()->SetRootGeometry();
169 // Create modules, materials, geometry
170 if (!gGeoManager) new TGeoManager("ALICE", "ALICE geometry");
172 TIter next(gAlice->Modules());
174 AliDebug(1, "Geometry creation:");
175 while((detector = dynamic_cast<AliModule*>(next()))) {
177 // Initialise detector materials and geometry
178 detector->CreateMaterials();
179 detector->CreateGeometry();
180 AliInfo(Form("%10s R:%.2fs C:%.2fs",
181 detector->GetName(),stw.RealTime(),stw.CpuTime()));
187 //_______________________________________________________________________
188 Bool_t AliMC::MisalignGeometry()
190 // Call misalignment code if AliSimulation object was defined.
192 if(!AliSimulation::Instance()->IsGeometryFromFile()){
193 //Set alignable volumes for the whole geometry
194 SetAllAlignableVolumes();
196 // Misalign geometry via AliSimulation instance
197 if (!AliSimulation::Instance()) return kFALSE;
198 AliGeomManager::SetGeometry(gGeoManager);
199 if(!AliGeomManager::CheckSymNamesLUT("ALL"))
200 AliFatal("Current loaded geometry differs in the definition of symbolic names!");
202 return AliSimulation::Instance()->MisalignGeometry(AliRunLoader::Instance());
205 //_______________________________________________________________________
206 void AliMC::ConstructOpGeometry()
209 // Loop all detector modules and call DefineOpticalProperties() method
212 TIter next(gAlice->Modules());
214 AliInfo("Optical properties definition");
215 while((detector = dynamic_cast<AliModule*>(next()))) {
216 // Initialise detector geometry
217 if(AliSimulation::Instance()->IsGeometryFromFile()) detector->CreateMaterials();
218 // Initialise detector optical properties
219 detector->DefineOpticalProperties();
223 #include <TPDGCode.h>
224 //_______________________________________________________________________
225 void AliMC::AddParticles()
228 // Add particles (not present in Geant3 or Geant4)
231 // --------------------------------------------------------------------
232 // An example of adding a particle He5 with defined decay mode
233 // (TO BE REMOVED after a useful code is added)
235 cout << "########## AliMC::AddParticles" << endl;
238 TVirtualMC::GetMC()->DefineParticle(1010010030, "HyperTriton", kPTHadron, 2.99131 , 1.0, 2.632e-10,"Ion", 0.0, 0, 1, 0, 0, 0, 0, 0, 3, kFALSE);
240 TVirtualMC::GetMC()->DefineParticle(-1010010030, "AntiHyperTriton", kPTHadron, 2.99131 , 1.0, 2.632e-10,"Ion", 0.0, 0, 1, 0, 0, 0, 0, 0, 3, kFALSE);
243 TVirtualMC::GetMC()->DefineParticle(1010010040, "Hyperhydrog4", kPTHadron, 3.931 , 1.0, 2.632e-10,"Ion", 0.0, 0, 1, 0, 0, 0, 0, 0, 4, kFALSE);
244 //Anti-Hyper hydrogen 4
245 TVirtualMC::GetMC()->DefineParticle(-1010010040, "AntiHyperhydrog4", kPTHadron, 3.931 , 1.0, 2.632e-10,"Ion", 0.0, 0, 1, 0, 0, 0, 0, 0, 4, kFALSE);
248 TVirtualMC::GetMC()->DefineParticle(1010020040, "Hyperhelium4", kPTHadron, 3.929 , 2.0, 2.632e-10,"Ion", 0.0, 0, 1, 0, 0, 0, 0, 0, 4, kFALSE);
249 //Anti-Hyper helium 4
250 TVirtualMC::GetMC()->DefineParticle(-1010020040, "AntiHyperhelium4", kPTHadron, 3.929 , 2.0, 2.632e-10,"Ion", 0.0, 0, 1, 0, 0, 0, 0, 0, 4, kFALSE);
254 TVirtualMC::GetMC()->DefineParticle(1010000020, "LambdaNeutron", kPTNeutron, 2.054 , 0.0, 2.632e-10,"Hadron", 0.0, 0, 1, 0, 0, 0, 0, 0, 2, kFALSE);
256 //Anti-Lambda-Neutron
257 TVirtualMC::GetMC()->DefineParticle(-1010000020, "AntiLambdaNeutron", kPTNeutron, 2.054 , 0.0, 2.632e-10,"Hadron", 0.0, 0, 1, 0, 0, 0, 0, 0, 2, kFALSE);
260 TVirtualMC::GetMC()->DefineParticle(1020000020, "Hdibaryon", kPTNeutron, 2.21, 0.0, 2.632e-10,"Hadron", 0.0, 0, 1, 0, 0, 0, 0, 0, 2, kFALSE);
262 TVirtualMC::GetMC()->DefineParticle(-1020000020, "AntiHdibaryon", kPTNeutron, 2.21 , 0.0, 2.632e-10,"Hadron", 0.0, 0, 1, 0, 0, 0, 0, 0, 2, kFALSE);
265 TVirtualMC::GetMC()->DefineParticle(1030000020, "Xi0Proton", kPTHadron, 2.248 , 1.0, 1.333e-10,"Ion", 0.0, 0, 1, 0, 0, 0, 0, 0, 2, kFALSE);
268 TVirtualMC::GetMC()->DefineParticle(-1030000020, "AntiXi0Proton", kPTHadron, 2.248 , 1.0, 1.333e-10,"Ion", 0.0, 0, 1, 0, 0, 0, 0, 0, 2, kFALSE);
270 //Lambda-Neutron-Neutron
271 TVirtualMC::GetMC()->DefineParticle(1010000030, "LambdaNeutronNeutron", kPTNeutron, 2.982 , 0.0, 2.632e-10,"Hadron", 0.0, 0, 1, 0, 0, 0, 0, 0, 2, kFALSE);
273 //Anti-Lambda-Neutron-Neutron
274 TVirtualMC::GetMC()->DefineParticle(-1010000030, "AntiLambdaNeutronNeutron", kPTNeutron, 2.982 , 0.0, 2.632e-10,"Hadron", 0.0, 0, 1, 0, 0, 0, 0, 0, 2, kFALSE);
277 // Define the 2- and 3-body phase space decay for the Hyper-Triton
281 for (Int_t kz = 0; kz < 6; kz++) {
288 mode[0][0] = 1000020030; // Helium3
289 mode[0][1] = -211; // negative pion
292 mode[1][0] = 1000010020; // deuteron
293 mode[1][1] = 2212; // proton
294 mode[1][2] = -211; // negative pion
296 TVirtualMC::GetMC()->SetDecayMode(1010010030,bratio,mode);
300 // Define the 2- and 3-body phase space decay for the Anti-Hyper-Triton
304 for (Int_t kz = 0; kz < 6; kz++) {
311 amode[0][0] = -1000020030; // anti- Helium3
312 amode[0][1] = 211; // positive pion
314 amode[1][0] = -1000010020; // anti-deuteron
315 amode[1][1] = -2212; // anti-proton
316 amode[1][2] = 211; // positive pion
318 TVirtualMC::GetMC()->SetDecayMode(-1010010030,abratio,amode);
320 ////// ----------Hypernuclei with Mass=4 ----------- //////////
322 // Define the 2- and 3-body phase space decay for the Hyper Hydrogen 4
327 for (Int_t kz = 0; kz < 6; kz++) {
334 mode3[0][0] = 1000020040; // Helium4
335 mode3[0][1] = -211; // negative pion
338 mode3[1][0] = 1000010030; // tritium
339 mode3[1][1] = 2212; // proton
340 mode3[1][2] = -211; // negative pion
342 TVirtualMC::GetMC()->SetDecayMode(1010010040,bratio3,mode3);
345 // Define the 2- and 3-body phase space decay for the Hyper Hydrogen 4
349 for (Int_t kz = 0; kz < 6; kz++) {
356 amode3[0][0] = -1000020040; // anti- Helium4
357 amode3[0][1] = 211; // positive pion
359 amode3[1][0] = -1000010030; // anti-tritium
360 amode3[1][1] = -2212; // anti-proton
361 amode3[1][2] = 211; // positive pion
363 TVirtualMC::GetMC()->SetDecayMode(-1010010040,abratio3,amode3);
366 // Define the 3-body phase space decay for the Hyper Helium 4
370 for (Int_t kz = 0; kz < 6; kz++) {
377 mode4[0][0] = 1000020030; // Helium3
378 mode4[0][1] = -211; // negative pion
379 mode4[0][2] = 2212; // proton
381 TVirtualMC::GetMC()->SetDecayMode(1010020040,bratio4,mode4);
384 // Define the 2-body phase space decay for the Anti-Hyper Helium 4
388 for (Int_t kz = 0; kz < 6; kz++) {
395 amode4[0][0] = -1000020030; // anti-Helium 3
396 amode4[0][1] = 211; // positive pion
397 amode4[0][2] = -2212; // anti proton
399 TVirtualMC::GetMC()->SetDecayMode(-1010020040,abratio4,amode4);
402 // Define the 2-body phase space decay for the Lambda-neutron boundstate
406 for (Int_t kz = 0; kz < 6; kz++) {
413 mode1[0][0] = 1000010020; // deuteron
414 mode1[0][1] = -211; // negative pion
416 TVirtualMC::GetMC()->SetDecayMode(1010000020,bratio1,mode1);
419 // Define the 2-body phase space decay for the Anti-Lambda-neutron boundstate
423 for (Int_t kz = 0; kz < 6; kz++) {
430 amode1[0][0] = -1000010020; // anti-deuteron
431 amode1[0][1] = 211; // positive pion
433 TVirtualMC::GetMC()->SetDecayMode(-1010000020,abratio1,amode1);
435 // Define the 2-body phase space decay for the H-Dibaryon
439 for (Int_t kz = 0; kz < 6; kz++) {
446 mode2[0][0] = 3122; // Lambda
447 mode2[0][1] = 2212; // proton
448 mode2[0][2] = -211; // negative pion
450 TVirtualMC::GetMC()->SetDecayMode(1020000020,bratio2,mode2);
452 // Define the 2-body phase space decay for the Anti-H-Dibaryon
456 for (Int_t kz = 0; kz < 6; kz++) {
463 amode2[0][0] = -3122; // anti-deuteron
464 amode2[0][1] = -2212; // anti-proton
465 amode2[0][2] = 211; // positive pion
467 TVirtualMC::GetMC()->SetDecayMode(-1020000020,abratio2,amode2);
469 // Define the 2-body phase space decay for the Xi0P
473 for (Int_t kz = 0; kz < 6; kz++) {
480 mode5[0][0] = 3122; // Lambda
481 mode5[0][1] = 2212; // proton
483 TVirtualMC::GetMC()->SetDecayMode(1030000020,bratio5,mode5);
485 // Define the 2-body phase space decay for the Anti-Xi0P
489 for (Int_t kz = 0; kz < 6; kz++) {
496 amode5[0][0] = -3122; // anti-Lambda
497 amode5[0][1] = -2212; // anti-proton
499 TVirtualMC::GetMC()->SetDecayMode(-1030000020,abratio5,amode5);
501 // Define the 2-body phase space decay for the Lambda-Neutron-Neutron
505 for (Int_t kz = 0; kz < 6; kz++) {
512 mode6[0][0] = 1000010030; // triton
513 mode6[0][1] = -211; // pion
515 TVirtualMC::GetMC()->SetDecayMode(1010000030,bratio6,mode6);
517 // Define the 2-body phase space decay for the Anti-Lambda-Neutron-Neutron
521 for (Int_t kz = 0; kz < 6; kz++) {
528 amode6[0][0] = -1000010030; // anti-triton
529 amode6[0][1] = 211; // pion
531 TVirtualMC::GetMC()->SetDecayMode(-1010000030,abratio6,amode6);
533 // end of the example
534 // --------------------------------------------------------------------
537 //_______________________________________________________________________
538 void AliMC::InitGeometry()
541 // Initialize detectors
544 AliInfo("Initialisation:");
546 TIter next(gAlice->Modules());
548 while((detector = dynamic_cast<AliModule*>(next()))) {
551 AliInfo(Form("%10s R:%.2fs C:%.2fs",
552 detector->GetName(),stw.RealTime(),stw.CpuTime()));
556 //_______________________________________________________________________
557 void AliMC::SetGeometryFromCDB()
559 // Set the loading of geometry from cdb instead of creating it
560 // A default CDB storage needs to be set before this method is called
561 if(AliCDBManager::Instance()->IsDefaultStorageSet() &&
562 AliCDBManager::Instance()->GetRun() >= 0)
563 AliSimulation::Instance()->SetGeometryFile("*OCDB*");
565 AliError("Loading of geometry from CDB ignored. First set a default CDB storage!");
568 //_______________________________________________________________________
569 Bool_t AliMC::IsGeometryFromCDB() const
571 return (strcmp(AliSimulation::Instance()->GetGeometryFile(),"*OCDB*")==0);
574 //_______________________________________________________________________
575 void AliMC::SetAllAlignableVolumes()
578 // Add alignable volumes (TGeoPNEntries) looping on all
582 AliInfo(Form("Setting entries for all alignable volumes of active detectors"));
584 TIter next(gAlice->Modules());
585 while((detector = dynamic_cast<AliModule*>(next()))) {
586 detector->AddAlignableVolumes();
590 //_______________________________________________________________________
591 void AliMC::GeneratePrimaries()
594 // Generate primary particles and fill them in the stack.
597 Generator()->Generate();
600 //_______________________________________________________________________
601 void AliMC::SetGenerator(AliGenerator *generator)
604 // Load the event generator
606 if(!fGenerator) fGenerator = generator;
609 //_______________________________________________________________________
610 void AliMC::ResetGenerator(AliGenerator *generator)
613 // Load the event generator
617 AliWarning(Form("Replacing generator %s with %s",
618 fGenerator->GetName(),generator->GetName()));
621 AliWarning(Form("Replacing generator %s with NULL",
622 fGenerator->GetName()));
625 fGenerator = generator;
628 //_______________________________________________________________________
629 void AliMC::FinishRun()
631 // Clean generator information
632 AliDebug(1, "fGenerator->FinishRun()");
633 fGenerator->FinishRun();
635 // Monitoring information
638 fMonitor->Export("timing.root");
641 //Output energy summary tables
642 AliDebug(1, "EnergySummary()");
643 ToAliDebug(1, EnergySummary());
646 //_______________________________________________________________________
647 void AliMC::BeginPrimary()
650 // Called at the beginning of each primary track
655 ResetTrackReferences();
658 //_______________________________________________________________________
659 void AliMC::PreTrack()
661 // Actions before the track's transport
663 //verbose.PreTrack();
665 TObjArray &dets = *gAlice->Modules();
668 for(Int_t i=0; i<=gAlice->GetNdets(); i++)
669 if((module = dynamic_cast<AliModule*>(dets[i])))
673 //_______________________________________________________________________
674 void AliMC::Stepping()
677 // Called at every step during transport
679 //verbose.Stepping();
681 Int_t id = DetFromMate(TVirtualMC::GetMC()->CurrentMedium());
685 if ( TVirtualMC::GetMC()->IsNewTrack() &&
686 TVirtualMC::GetMC()->TrackTime() == 0. &&
688 fRDecayMax > fRDecayMin &&
689 TVirtualMC::GetMC()->TrackPid() == fDecayPdg )
691 FixParticleDecaytime();
694 // --- If monitoring timing was requested, monitor the step
695 if (fUseMonitoring) {
697 fMonitor = new AliTransportMonitor(TVirtualMC::GetMC()->NofVolumes()+1);
700 if (TVirtualMC::GetMC()->IsNewTrack() || TVirtualMC::GetMC()->TrackTime() == 0. || TVirtualMC::GetMC()->TrackStep()<1.1E-10) {
701 fMonitor->DummyStep();
705 Int_t volId = TVirtualMC::GetMC()->CurrentVolID(copy);
706 Int_t pdg = TVirtualMC::GetMC()->TrackPid();
707 TLorentzVector xyz, pxpypz;
708 TVirtualMC::GetMC()->TrackPosition(xyz);
709 TVirtualMC::GetMC()->TrackMomentum(pxpypz);
710 fMonitor->StepInfo(volId, pdg, pxpypz.E(), xyz.X(), xyz.Y(), xyz.Z());
714 // --- If lego option, do it and leave
715 if (AliSimulation::Instance()->Lego())
716 AliSimulation::Instance()->Lego()->StepManager();
719 //Update energy deposition tables
720 AddEnergyDeposit(TVirtualMC::GetMC()->CurrentVolID(copy),TVirtualMC::GetMC()->Edep());
722 // write tracke reference for track which is dissapearing - MI
724 if (TVirtualMC::GetMC()->IsTrackDisappeared() && !(TVirtualMC::GetMC()->IsTrackAlive())) {
725 if (TVirtualMC::GetMC()->Etot() > 0.05) AddTrackReference(GetCurrentTrackNumber(),
726 AliTrackReference::kDisappeared);
731 //Call the appropriate stepping routine;
732 AliModule *det = dynamic_cast<AliModule*>(gAlice->Modules()->At(id));
733 if(det && det->StepManagerIsEnabled()) {
739 //_______________________________________________________________________
740 void AliMC::EnergySummary()
743 // Print summary of deposited energy
749 Int_t kn, i, left, j, id;
750 const Float_t kzero=0;
751 Int_t ievent=AliRunLoader::Instance()->GetHeader()->GetEvent()+1;
753 // Energy loss information
755 printf("***************** Energy Loss Information per event (GEV) *****************\n");
756 for(kn=1;kn<fEventEnergy.GetSize();kn++) {
759 fEventEnergy[ndep]=kn;
764 ed2=100*TMath::Sqrt(TMath::Max(ed2-ed*ed,kzero))/ed;
767 fSummEnergy[ndep]=ed;
768 fSum2Energy[ndep]=TMath::Min(static_cast<Float_t>(99.),TMath::Max(ed2,kzero));
773 for(kn=0;kn<(ndep-1)/3+1;kn++) {
775 for(i=0;i<(3<left?3:left);i++) {
777 id=Int_t (fEventEnergy[j]+0.1);
778 printf(" %s %10.3f +- %10.3f%%;",TVirtualMC::GetMC()->VolName(id),fSummEnergy[j],fSum2Energy[j]);
783 // Relative energy loss in different detectors
784 printf("******************** Relative Energy Loss per event ********************\n");
785 printf("Total energy loss per event %10.3f GeV\n",edtot);
786 for(kn=0;kn<(ndep-1)/5+1;kn++) {
788 for(i=0;i<(5<left?5:left);i++) {
790 id=Int_t (fEventEnergy[j]+0.1);
791 printf(" %s %10.3f%%;",TVirtualMC::GetMC()->VolName(id),100*fSummEnergy[j]/edtot);
795 for(kn=0;kn<75;kn++) printf("*");
799 // Reset the TArray's
800 // fEventEnergy.Set(0);
801 // fSummEnergy.Set(0);
802 // fSum2Energy.Set(0);
805 //_____________________________________________________________________________
806 void AliMC::BeginEvent()
809 // Clean-up previous event
811 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
812 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
813 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
814 AliDebug(1, " BEGINNING EVENT ");
815 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
816 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
817 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
819 AliRunLoader *runloader=AliRunLoader::Instance();
821 /*******************************/
822 /* Clean after eventual */
824 /*******************************/
827 //Set the next event in Run Loader -> Cleans trees (TreeK and all trees in detectors),
828 gAlice->SetEventNrInRun(gAlice->GetEventNrInRun()+1);
829 runloader->SetEventNumber(gAlice->GetEventNrInRun());// sets new files, cleans the previous event stuff, if necessary, etc.,
830 AliDebug(1, Form("EventNr is %d",gAlice->GetEventNrInRun()));
832 fEventEnergy.Reset();
833 // Clean detector information
835 if (runloader->Stack())
836 runloader->Stack()->Reset();//clean stack -> tree is unloaded
838 runloader->MakeStack();//or make a new one
840 // Random engine status
843 if ( fSaveRndmStatus || fSaveRndmEventStatus) {
844 TString fileName="random";
845 if ( fSaveRndmEventStatus ) {
847 fileName += gAlice->GetEventNrInRun();
851 // write ROOT random engine status
852 cout << "Saving random engine status in " << fileName.Data() << endl;
853 TFile f(fileName.Data(),"RECREATE");
854 gRandom->Write(fileName.Data());
857 if ( fReadRndmStatus ) {
858 //read ROOT random engine status
859 cout << "Reading random engine status from " << fRndmFileName.Data() << endl;
860 TFile f(fRndmFileName.Data());
861 gRandom->Read(fRndmFileName.Data());
864 if(AliSimulation::Instance()->Lego() == 0x0)
866 AliDebug(1, "fRunLoader->MakeTree(K)");
867 runloader->MakeTree("K");
870 AliDebug(1, "TVirtualMC::GetMC()->SetStack(fRunLoader->Stack())");
871 TVirtualMC::GetMC()->SetStack(runloader->Stack());//Was in InitMC - but was moved here
872 //because we don't have guarantee that
873 //stack pointer is not going to change from event to event
874 //since it bellobgs to header and is obtained via RunLoader
876 // Reset all Detectors & kinematics & make/reset trees
879 runloader->GetHeader()->Reset(AliCDBManager::Instance()->GetRun(),gAlice->GetEvNumber(),
880 gAlice->GetEventNrInRun());
881 // fRunLoader->WriteKinematics("OVERWRITE"); is there any reason to rewrite here since MakeTree does so
883 if(AliSimulation::Instance()->Lego())
885 AliSimulation::Instance()->Lego()->BeginEvent();
890 AliDebug(1, "ResetHits()");
893 AliDebug(1, "fRunLoader->MakeTree(H)");
894 runloader->MakeTree("H");
898 MakeTmpTrackRefsTree();
899 //create new branches and SetAdresses
900 TIter next(gAlice->Modules());
902 while((detector = (AliModule*)next()))
904 AliDebug(2, Form("%s->MakeBranch(H)",detector->GetName()));
905 detector->MakeBranch("H");
909 //_______________________________________________________________________
910 void AliMC::ResetHits()
913 // Reset all Detectors hits
915 TIter next(gAlice->Modules());
917 while((detector = dynamic_cast<AliModule*>(next()))) {
918 detector->ResetHits();
922 //_______________________________________________________________________
923 void AliMC::ResetDigits()
926 // Reset all Detectors digits
928 TIter next(gAlice->Modules());
930 while((detector = dynamic_cast<AliModule*>(next()))) {
931 detector->ResetDigits();
935 //_______________________________________________________________________
936 void AliMC::ResetSDigits()
939 // Reset all Detectors digits
941 TIter next(gAlice->Modules());
943 while((detector = dynamic_cast<AliModule*>(next()))) {
944 detector->ResetSDigits();
948 //_______________________________________________________________________
949 void AliMC::PostTrack()
951 // Posts tracks for each module
953 TObjArray &dets = *gAlice->Modules();
956 for(Int_t i=0; i<=gAlice->GetNdets(); i++)
957 if((module = dynamic_cast<AliModule*>(dets[i])))
961 //_______________________________________________________________________
962 void AliMC::FinishPrimary()
965 // Called at the end of each primary track
968 AliRunLoader *runloader=AliRunLoader::Instance();
969 // static Int_t count=0;
970 // const Int_t times=10;
971 // This primary is finished, purify stack
972 #if ROOT_VERSION_CODE > 262152
973 if (!(TVirtualMC::GetMC()->SecondariesAreOrdered())) {
974 if (runloader->Stack()->ReorderKine()) RemapHits();
977 if (runloader->Stack()->PurifyKine()) RemapHits();
979 TIter next(gAlice->Modules());
981 while((detector = dynamic_cast<AliModule*>(next()))) {
982 detector->FinishPrimary();
983 AliLoader* loader = detector->GetLoader();
986 TTree* treeH = loader->TreeH();
987 if (treeH) treeH->Fill(); //can be Lego run and treeH can not exist
991 // Write out track references if any
992 if (fTmpTreeTR) fTmpTreeTR->Fill();
995 void AliMC::RemapHits()
998 // Remaps the track labels of the hits
999 AliRunLoader *runloader=AliRunLoader::Instance();
1000 AliStack* stack = runloader->Stack();
1001 TList* hitLists = GetHitLists();
1002 TIter next(hitLists);
1003 TCollection *hitList;
1005 while((hitList = dynamic_cast<TCollection*>(next()))) {
1006 TIter nexthit(hitList);
1008 while((hit = dynamic_cast<AliHit*>(nexthit()))) {
1009 hit->SetTrack(stack->TrackLabel(hit->GetTrack()));
1014 // This for detectors which have a special mapping mechanism
1015 // for hits, such as TPC and TRD
1019 TObjArray* modules = gAlice->Modules();
1020 TIter nextmod(modules);
1022 while((module = (AliModule*) nextmod())) {
1023 AliDetector* det = dynamic_cast<AliDetector*> (module);
1024 if (det) det->RemapTrackHitIDs(stack->TrackLabelMap());
1027 RemapTrackReferencesIDs(stack->TrackLabelMap());
1030 //_______________________________________________________________________
1031 void AliMC::FinishEvent()
1034 // Called at the end of the event.
1037 if(AliSimulation::Instance()->Lego()) AliSimulation::Instance()->Lego()->FinishEvent();
1039 TIter next(gAlice->Modules());
1040 AliModule *detector;
1041 while((detector = dynamic_cast<AliModule*>(next()))) {
1042 detector->FinishEvent();
1045 //Update the energy deposit tables
1047 for(i=0;i<fEventEnergy.GetSize();i++)
1049 fSummEnergy[i]+=fEventEnergy[i];
1050 fSum2Energy[i]+=fEventEnergy[i]*fEventEnergy[i];
1053 AliRunLoader *runloader=AliRunLoader::Instance();
1055 AliHeader* header = runloader->GetHeader();
1056 AliStack* stack = runloader->Stack();
1057 if ( (header == 0x0) || (stack == 0x0) )
1058 {//check if we got header and stack. If not cry and exit aliroot
1059 AliFatal("Can not get the stack or header from LOADER");
1060 return;//never reached
1062 // Update Header information
1063 header->SetNprimary(stack->GetNprimary());
1064 header->SetNtrack(stack->GetNtrack());
1065 header->SetTimeStamp(AliSimulation::Instance()->GenerateTimeStamp());
1067 // Write out the kinematics
1068 if (!AliSimulation::Instance()->Lego()) stack->FinishEvent();
1070 // Synchronize the TreeTR with TreeK
1071 if (fTmpTreeTR) ReorderAndExpandTreeTR();
1073 // Write out the event Header information
1074 TTree* treeE = runloader->TreeE();
1077 header->SetStack(stack);
1082 AliError("Can not get TreeE from RL");
1085 if(AliSimulation::Instance()->Lego() == 0x0)
1087 runloader->WriteKinematics("OVERWRITE");
1088 runloader->WriteTrackRefs("OVERWRITE");
1089 runloader->WriteHits("OVERWRITE");
1092 AliDebug(1, "<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<");
1093 AliDebug(1, "<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<");
1094 AliDebug(1, "<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<");
1095 AliDebug(1, " FINISHING EVENT ");
1096 AliDebug(1, "<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<");
1097 AliDebug(1, "<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<");
1098 AliDebug(1, "<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<");
1101 //_______________________________________________________________________
1104 // MC initialization
1107 //=================Create Materials and geometry
1108 TVirtualMC::GetMC()->Init();
1109 // Set alignable volumes for the whole geometry (with old root)
1110 #if ROOT_VERSION_CODE < 331527
1111 SetAllAlignableVolumes();
1113 //Read the cuts for all materials
1115 //Build the special IMEDIA table
1118 //Compute cross-sections
1119 TVirtualMC::GetMC()->BuildPhysics();
1121 //Initialise geometry deposition table
1122 fEventEnergy.Set(TVirtualMC::GetMC()->NofVolumes()+1);
1123 fSummEnergy.Set(TVirtualMC::GetMC()->NofVolumes()+1);
1124 fSum2Energy.Set(TVirtualMC::GetMC()->NofVolumes()+1);
1126 // Register MC in configuration
1127 AliConfig::Instance()->Add(TVirtualMC::GetMC());
1130 //_______________________________________________________________________
1131 void AliMC::MediaTable()
1134 // Built media table to get from the media number to
1138 Int_t kz, nz, idt, lz, i, k, ind;
1140 TObjArray &dets = *gAlice->Detectors();
1142 Int_t ndets=gAlice->GetNdets();
1144 // For all detectors
1145 for (kz=0;kz<ndets;kz++) {
1146 // If detector is defined
1147 if((det=dynamic_cast<AliModule*>(dets[kz]))) {
1148 TArrayI &idtmed = *(det->GetIdtmed());
1149 for(nz=0;nz<100;nz++) {
1151 // Find max and min material number
1152 if((idt=idtmed[nz])) {
1153 det->LoMedium() = det->LoMedium() < idt ? det->LoMedium() : idt;
1154 det->HiMedium() = det->HiMedium() > idt ? det->HiMedium() : idt;
1157 if(det->LoMedium() > det->HiMedium()) {
1158 det->LoMedium() = 0;
1159 det->HiMedium() = 0;
1161 if(det->HiMedium() > fImedia->GetSize()) {
1162 AliError(Form("Increase fImedia from %d to %d",
1163 fImedia->GetSize(),det->HiMedium()));
1166 // Tag all materials in rage as belonging to detector kz
1167 for(lz=det->LoMedium(); lz<= det->HiMedium(); lz++) {
1174 // Print summary table
1175 AliInfo("Tracking media ranges:");
1177 for(i=0;i<(ndets-1)/6+1;i++) {
1178 for(k=0;k< (6<ndets-i*6?6:ndets-i*6);k++) {
1180 det=dynamic_cast<AliModule*>(dets[ind]);
1182 printf(" %6s: %3d -> %3d;",det->GetName(),det->LoMedium(),
1185 printf(" %6s: %3d -> %3d;","NULL",0,0);
1192 //_______________________________________________________________________
1193 void AliMC::ReadTransPar()
1196 // Read filename to set the transport parameters
1200 const Int_t kncuts=10;
1201 const Int_t knflags=12;
1202 const Int_t knpars=kncuts+knflags;
1203 const char kpars[knpars][7] = {"CUTGAM" ,"CUTELE","CUTNEU","CUTHAD","CUTMUO",
1204 "BCUTE","BCUTM","DCUTE","DCUTM","PPCUTM","ANNI",
1205 "BREM","COMP","DCAY","DRAY","HADR","LOSS",
1206 "MULS","PAIR","PHOT","RAYL","STRA"};
1210 Float_t cut[kncuts];
1211 Int_t flag[knflags];
1212 Int_t i, itmed, iret, jret, ktmed, kz;
1215 // See whether the file is there
1216 filtmp=gSystem->ExpandPathName(fTransParName.Data());
1217 lun=fopen(filtmp,"r");
1220 AliWarning(Form("File %s does not exist!",fTransParName.Data()));
1225 // Initialise cuts and flags
1226 for(i=0;i<kncuts;i++) cut[i]=-99;
1227 for(i=0;i<knflags;i++) flag[i]=-99;
1230 // Read up to the end of line excluded
1231 iret=fscanf(lun,"%255[^\n]",line);
1237 // Read the end of line
1238 jret = fscanf(lun,"%*c");
1240 if(line[0]=='*') continue;
1242 iret=sscanf(line,"%6s %d %f %f %f %f %f %f %f %f %f %f %d %d %d %d %d %d %d %d %d %d %d %d",
1243 detName,&itmed,&cut[0],&cut[1],&cut[2],&cut[3],&cut[4],&cut[5],&cut[6],&cut[7],&cut[8],
1244 &cut[9],&flag[0],&flag[1],&flag[2],&flag[3],&flag[4],&flag[5],&flag[6],&flag[7],
1245 &flag[8],&flag[9],&flag[10],&flag[11]);
1249 AliWarning(Form("Error reading file %s",fTransParName.Data()));
1252 // Check that the module exist
1253 AliModule *mod = gAlice->GetModule(detName);
1255 // Get the array of media numbers
1256 TArrayI &idtmed = *mod->GetIdtmed();
1257 // Check that the tracking medium code is valid
1258 if(0<=itmed && itmed < 100) {
1259 ktmed=idtmed[itmed];
1261 AliWarning(Form("Invalid tracking medium code %d for %s",itmed,mod->GetName()));
1264 // Set energy thresholds
1265 for(kz=0;kz<kncuts;kz++) {
1267 AliDebug(2, Form("%-6s set to %10.3E for tracking medium code %4d for %s",
1268 kpars[kz],cut[kz],itmed,mod->GetName()));
1269 TVirtualMC::GetMC()->Gstpar(ktmed,kpars[kz],cut[kz]);
1272 // Set transport mechanisms
1273 for(kz=0;kz<knflags;kz++) {
1275 AliDebug(2, Form("%-6s set to %10d for tracking medium code %4d for %s",
1276 kpars[kncuts+kz],flag[kz],itmed,mod->GetName()));
1277 TVirtualMC::GetMC()->Gstpar(ktmed,kpars[kncuts+kz],Float_t(flag[kz]));
1281 AliWarning(Form("Invalid medium code %d",itmed));
1285 AliDebug(1, Form("%s not present",detName));
1291 //_______________________________________________________________________
1292 void AliMC::SetTransPar(const char *filename)
1295 // Sets the file name for transport parameters
1297 fTransParName = filename;
1300 //_______________________________________________________________________
1301 void AliMC::AddHit(Int_t id, Int_t track, Int_t *vol, Float_t *hits) const
1304 // Add a hit to detector id
1306 TObjArray &dets = *gAlice->Modules();
1307 if(dets[id]) static_cast<AliModule*>(dets[id])->AddHit(track,vol,hits);
1310 //_______________________________________________________________________
1311 void AliMC::AddDigit(Int_t id, Int_t *tracks, Int_t *digits) const
1314 // Add digit to detector id
1316 TObjArray &dets = *gAlice->Modules();
1317 if(dets[id]) static_cast<AliModule*>(dets[id])->AddDigit(tracks,digits);
1320 //_______________________________________________________________________
1321 Int_t AliMC::GetCurrentTrackNumber() const {
1323 // Returns current track
1325 return AliRunLoader::Instance()->Stack()->GetCurrentTrackNumber();
1328 //_______________________________________________________________________
1329 void AliMC::DumpPart (Int_t i) const
1332 // Dumps particle i in the stack
1334 AliRunLoader * runloader = AliRunLoader::Instance();
1335 if (runloader->Stack())
1336 runloader->Stack()->DumpPart(i);
1339 //_______________________________________________________________________
1340 void AliMC::DumpPStack () const
1343 // Dumps the particle stack
1345 AliRunLoader * runloader = AliRunLoader::Instance();
1346 if (runloader->Stack())
1347 runloader->Stack()->DumpPStack();
1350 //_______________________________________________________________________
1351 Int_t AliMC::GetNtrack() const {
1353 // Returns number of tracks in stack
1356 AliRunLoader * runloader = AliRunLoader::Instance();
1357 if (runloader->Stack())
1358 ntracks = runloader->Stack()->GetNtrack();
1362 //_______________________________________________________________________
1363 Int_t AliMC::GetPrimary(Int_t track) const
1366 // return number of primary that has generated track
1368 Int_t nprimary = -999;
1369 AliRunLoader * runloader = AliRunLoader::Instance();
1370 if (runloader->Stack())
1371 nprimary = runloader->Stack()->GetPrimary(track);
1375 //_______________________________________________________________________
1376 TParticle* AliMC::Particle(Int_t i) const
1378 // Returns the i-th particle from the stack taking into account
1379 // the remaping done by PurifyKine
1380 AliRunLoader * runloader = AliRunLoader::Instance();
1382 if (runloader->Stack())
1383 return runloader->Stack()->Particle(i);
1387 //_______________________________________________________________________
1388 const TObjArray* AliMC::Particles() const {
1390 // Returns pointer to Particles array
1392 AliRunLoader * runloader = AliRunLoader::Instance();
1394 if (runloader->Stack())
1395 return runloader->Stack()->Particles();
1399 //_______________________________________________________________________
1400 void AliMC::PushTrack(Int_t done, Int_t parent, Int_t pdg, const Float_t *pmom,
1401 const Float_t *vpos, const Float_t *polar, Float_t tof,
1402 TMCProcess mech, Int_t &ntr, Float_t weight, Int_t is) const
1404 // Delegate to stack
1406 AliRunLoader * runloader = AliRunLoader::Instance();
1408 if (runloader->Stack())
1409 runloader->Stack()->PushTrack(done, parent, pdg, pmom, vpos, polar, tof,
1410 mech, ntr, weight, is);
1413 //_______________________________________________________________________
1414 void AliMC::PushTrack(Int_t done, Int_t parent, Int_t pdg,
1415 Double_t px, Double_t py, Double_t pz, Double_t e,
1416 Double_t vx, Double_t vy, Double_t vz, Double_t tof,
1417 Double_t polx, Double_t poly, Double_t polz,
1418 TMCProcess mech, Int_t &ntr, Float_t weight, Int_t is) const
1420 // Delegate to stack
1422 AliRunLoader * runloader = AliRunLoader::Instance();
1424 if (runloader->Stack())
1425 runloader->Stack()->PushTrack(done, parent, pdg, px, py, pz, e, vx, vy, vz, tof,
1426 polx, poly, polz, mech, ntr, weight, is);
1429 //_______________________________________________________________________
1430 void AliMC::SetHighWaterMark(Int_t nt) const
1433 // Set high water mark for last track in event
1434 AliRunLoader * runloader = AliRunLoader::Instance();
1436 if (runloader->Stack())
1437 runloader->Stack()->SetHighWaterMark(nt);
1440 //_______________________________________________________________________
1441 void AliMC::KeepTrack(Int_t track) const
1444 // Delegate to stack
1446 AliRunLoader * runloader = AliRunLoader::Instance();
1448 if (runloader->Stack())
1449 runloader->Stack()->KeepTrack(track);
1452 //_______________________________________________________________________
1453 void AliMC::FlagTrack(Int_t track) const
1455 // Delegate to stack
1457 AliRunLoader * runloader = AliRunLoader::Instance();
1459 if (runloader->Stack())
1460 runloader->Stack()->FlagTrack(track);
1463 //_______________________________________________________________________
1464 void AliMC::SetCurrentTrack(Int_t track) const
1467 // Set current track number
1469 AliRunLoader * runloader = AliRunLoader::Instance();
1471 if (runloader->Stack())
1472 runloader->Stack()->SetCurrentTrack(track);
1475 //_______________________________________________________________________
1476 AliTrackReference* AliMC::AddTrackReference(Int_t label, Int_t id)
1479 // add a trackrefernce to the list
1480 Int_t primary = GetPrimary(label);
1481 Particle(primary)->SetBit(kKeepBit);
1483 Int_t nref = fTmpTrackReferences.GetEntriesFast();
1484 return new(fTmpTrackReferences[nref]) AliTrackReference(label, id);
1489 //_______________________________________________________________________
1490 void AliMC::ResetTrackReferences()
1493 // Reset all references
1495 fTmpTrackReferences.Clear();
1498 //_______________________________________________________________________
1499 void AliMC::RemapTrackReferencesIDs(const Int_t *map)
1502 // Remapping track reference
1503 // Called at finish primary
1506 Int_t nEntries = fTmpTrackReferences.GetEntries();
1507 for (Int_t i=0; i < nEntries; i++){
1508 AliTrackReference * ref = dynamic_cast<AliTrackReference*>(fTmpTrackReferences.UncheckedAt(i));
1510 Int_t newID = map[ref->GetTrack()];
1511 if (newID>=0) ref->SetTrack(newID);
1513 ref->SetBit(kNotDeleted,kFALSE);
1514 fTmpTrackReferences.RemoveAt(i);
1518 fTmpTrackReferences.Compress();
1521 //_______________________________________________________________________
1522 void AliMC::FixParticleDecaytime()
1525 // Fix the particle decay time according to rmin and rmax for decays
1529 TVirtualMC::GetMC()->TrackMomentum(p);
1530 Double_t tmin, tmax;
1533 // Transverse velocity
1534 Double_t vt = p.Pt() / p.E();
1536 if ((b = ((AliMagF*)TGeoGlobalMagField::Instance()->GetField())->SolenoidField()) > 0.) { // [kG]
1540 Double_t rho = p.Pt() / 0.0003 / b; // [cm]
1542 // Revolution frequency
1544 Double_t omega = vt / rho;
1546 // Maximum and minimum decay time
1548 // Check for curlers first
1549 const Double_t kOvRhoSqr2 = 1./(rho*TMath::Sqrt(2.));
1550 if (fRDecayMax * kOvRhoSqr2 > 1.) return;
1554 tmax = TMath::ACos((1.-fRDecayMax*kOvRhoSqr2)*(1.+fRDecayMax*kOvRhoSqr2)) / omega; // [ct]
1555 tmin = TMath::ACos((1.-fRDecayMin*kOvRhoSqr2)*(1.+fRDecayMin*kOvRhoSqr2)) / omega; // [ct]
1557 tmax = fRDecayMax / vt; // [ct]
1558 tmin = fRDecayMin / vt; // [ct]
1562 // Dial t using the two limits
1563 Double_t t = tmin + (tmax - tmin) * gRandom->Rndm(); // [ct]
1566 // Force decay time in transport code
1568 TVirtualMC::GetMC()->ForceDecayTime(t / 2.99792458e10);
1571 void AliMC::MakeTmpTrackRefsTree()
1573 // Make the temporary track reference tree
1574 fTmpFileTR = new TFile("TrackRefsTmp.root", "recreate");
1575 fTmpTreeTR = new TTree("TreeTR", "Track References");
1576 TClonesArray* pRef = &fTmpTrackReferences;
1577 fTmpTreeTR->Branch("TrackReferences", &pRef, 4000);
1580 //_______________________________________________________________________
1581 void AliMC::ReorderAndExpandTreeTR()
1584 // Reorder and expand the temporary track reference tree in order to match the kinematics tree
1587 AliRunLoader *rl = AliRunLoader::Instance();
1590 AliDebug(1, "fRunLoader->MakeTrackRefsContainer()");
1591 rl->MakeTrackRefsContainer();
1592 TTree * treeTR = rl->TreeTR();
1593 // make branch for central track references
1594 TClonesArray* pRef = &fTrackReferences;
1595 treeTR->Branch("TrackReferences", &pRef);
1597 AliStack* stack = rl->Stack();
1598 Int_t np = stack->GetNprimary();
1599 Int_t nt = fTmpTreeTR->GetEntries();
1601 // Loop over tracks and find the secondaries with the help of the kine tree
1604 for (Int_t ip = np - 1; ip > -1; ip--) {
1605 TParticle *part = stack->Particle(ip);
1606 //printf("Particle %5d %5d %5d %5d %5d \n", ip, part->GetPdgCode(), part->GetFirstMother(), part->GetFirstDaughter(), part->GetLastDaughter());
1608 // Skip primaries that have not been transported
1609 Int_t dau1 = part->GetFirstDaughter();
1611 if (!part->TestBit(kTransportBit)) continue;
1613 fTmpTreeTR->GetEntry(it++);
1614 Int_t nh = fTmpTrackReferences.GetEntries();
1615 // Determine range of secondaries produced by this primary
1617 Int_t inext = ip - 1;
1620 part = stack->Particle(inext);
1621 dau2 = part->GetFirstDaughter();
1622 if (!(part->TestBit(kTransportBit)) || dau2 == -1 || dau2 < np) {
1623 // if (dau2 == -1 || dau2 < np) {
1629 dau2 = stack->GetNtrack() - 1;
1632 } // find upper bound
1634 // printf("Check (1) %5d %5d %5d %5d %5d \n", ip, np, it, dau1, dau2);
1636 // Loop over reference hits and find secondary label
1637 for (Int_t id = dau1; (id <= dau2) && (dau1 > -1); id++) {
1638 for (Int_t ih = 0; ih < nh; ih++) {
1639 AliTrackReference* tr = (AliTrackReference*) fTmpTrackReferences.At(ih);
1640 Int_t label = tr->Label();
1642 if (label == ip) continue;
1643 if (label > dau2 || label < dau1)
1644 AliWarning(Form("Track Reference Label out of range !: %5d %5d %5d \n", label, dau1, dau2));
1647 Int_t nref = fTrackReferences.GetEntriesFast();
1648 new(fTrackReferences[nref]) AliTrackReference(*tr);
1652 fTrackReferences.Clear();
1657 // Now loop again and write the primaries
1659 for (Int_t ip = 0; ip < np; ip++) {
1660 TParticle* part = stack->Particle(ip);
1661 // if ((part->GetFirstDaughter() == -1 && part->GetStatusCode() <= 1) || part->GetFirstDaughter() >= np)
1662 if (part->TestBit(kTransportBit))
1664 // Skip particles that have not been transported
1665 fTmpTreeTR->GetEntry(it--);
1666 Int_t nh = fTmpTrackReferences.GetEntries();
1668 // Loop over reference hits and find primary labels
1669 for (Int_t ih = 0; ih < nh; ih++) {
1670 AliTrackReference* tr = (AliTrackReference*) fTmpTrackReferences.At(ih);
1671 Int_t label = tr->Label();
1673 Int_t nref = fTrackReferences.GetEntriesFast();
1674 new(fTrackReferences[nref]) AliTrackReference(*tr);
1679 fTrackReferences.Clear();
1683 if (ifills != stack->GetNtrack())
1684 AliWarning(Form("Number of entries in TreeTR (%5d) unequal to TreeK (%5d) \n", ifills, stack->GetNtrack()));
1688 fTmpFileTR->Close();
1690 fTmpTrackReferences.Clear();
1691 gSystem->Exec("rm -rf TrackRefsTmp.root");