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 cout << "########## AliMC::AddParticles" << endl;
234 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);
236 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);
239 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);
240 //Anti-Hyper hydrogen 4
241 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);
244 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);
245 //Anti-Hyper helium 4
246 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);
249 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);
251 //Anti-Lambda-Neutron
252 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);
255 TVirtualMC::GetMC()->DefineParticle(1020000020, "Hdibaryon", kPTNeutron, 2.23, 0.0, 2.632e-10,"Hadron", 0.0, 0, 1, 0, 0, 0, 0, 0, 2, kFALSE);
257 TVirtualMC::GetMC()->DefineParticle(-1020000020, "AntiHdibaryon", kPTNeutron, 2.23, 0.0, 2.632e-10,"Hadron", 0.0, 0, 1, 0, 0, 0, 0, 0, 2, kFALSE);
260 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);
263 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);
265 //Lambda-Neutron-Neutron
266 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);
268 //Anti-Lambda-Neutron-Neutron
269 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);
271 //Resonances not in Generators
272 // f0(980) assume 70 MeV as width (PDG: 40 to 100 MeV)
273 TVirtualMC::GetMC()->DefineParticle(9010221, "f0_980", kPTNeutron, 0.98 , 0.0, 9.403e-24,"Hadron", 7e-2, 0, 0, 0, 0, 0, 0, 0, 0, kTRUE);
275 // f2(1270) (PDG: width = 185 MeV)
276 TVirtualMC::GetMC()->DefineParticle(225, "f2_1270", kPTNeutron, 1.275 , 0.0, 3.558e-24,"Hadron", 0.185, 0, 0, 0, 0, 0, 0, 0, 0, kTRUE);
278 // Define the 2- and 3-body phase space decay for the Hyper-Triton
282 for (Int_t kz = 0; kz < 6; kz++) {
289 mode[0][0] = 1000020030; // Helium3
290 mode[0][1] = -211; // negative pion
293 mode[1][0] = 1000010020; // deuteron
294 mode[1][1] = 2212; // proton
295 mode[1][2] = -211; // negative pion
297 TVirtualMC::GetMC()->SetDecayMode(1010010030,bratio,mode);
301 // Define the 2- and 3-body phase space decay for the Anti-Hyper-Triton
305 for (Int_t kz = 0; kz < 6; kz++) {
312 amode[0][0] = -1000020030; // anti- Helium3
313 amode[0][1] = 211; // positive pion
315 amode[1][0] = -1000010020; // anti-deuteron
316 amode[1][1] = -2212; // anti-proton
317 amode[1][2] = 211; // positive pion
319 TVirtualMC::GetMC()->SetDecayMode(-1010010030,abratio,amode);
321 ////// ----------Hypernuclei with Mass=4 ----------- //////////
323 // Define the 2- and 3-body phase space decay for the Hyper Hydrogen 4
328 for (Int_t kz = 0; kz < 6; kz++) {
335 mode3[0][0] = 1000020040; // Helium4
336 mode3[0][1] = -211; // negative pion
339 mode3[1][0] = 1000010030; // tritium
340 mode3[1][1] = 2212; // proton
341 mode3[1][2] = -211; // negative pion
343 TVirtualMC::GetMC()->SetDecayMode(1010010040,bratio3,mode3);
346 // Define the 2- and 3-body phase space decay for the Hyper Hydrogen 4
350 for (Int_t kz = 0; kz < 6; kz++) {
357 amode3[0][0] = -1000020040; // anti- Helium4
358 amode3[0][1] = 211; // positive pion
360 amode3[1][0] = -1000010030; // anti-tritium
361 amode3[1][1] = -2212; // anti-proton
362 amode3[1][2] = 211; // positive pion
364 TVirtualMC::GetMC()->SetDecayMode(-1010010040,abratio3,amode3);
367 // Define the 3-body phase space decay for the Hyper Helium 4
371 for (Int_t kz = 0; kz < 6; kz++) {
378 mode4[0][0] = 1000020030; // Helium3
379 mode4[0][1] = -211; // negative pion
380 mode4[0][2] = 2212; // proton
382 TVirtualMC::GetMC()->SetDecayMode(1010020040,bratio4,mode4);
385 // Define the 2-body phase space decay for the Anti-Hyper Helium 4
389 for (Int_t kz = 0; kz < 6; kz++) {
396 amode4[0][0] = -1000020030; // anti-Helium 3
397 amode4[0][1] = 211; // positive pion
398 amode4[0][2] = -2212; // anti proton
400 TVirtualMC::GetMC()->SetDecayMode(-1010020040,abratio4,amode4);
403 // Define the 2-body phase space decay for the Lambda-neutron boundstate
407 for (Int_t kz = 0; kz < 6; kz++) {
414 mode1[0][0] = 1000010020; // deuteron
415 mode1[0][1] = -211; // negative pion
417 TVirtualMC::GetMC()->SetDecayMode(1010000020,bratio1,mode1);
420 // Define the 2-body phase space decay for the Anti-Lambda-neutron boundstate
424 for (Int_t kz = 0; kz < 6; kz++) {
431 amode1[0][0] = -1000010020; // anti-deuteron
432 amode1[0][1] = 211; // positive pion
434 TVirtualMC::GetMC()->SetDecayMode(-1010000020,abratio1,amode1);
436 // Define the 2-body phase space decay for the H-Dibaryon
440 for (Int_t kz = 0; kz < 6; kz++) {
447 mode2[0][0] = 3122; // Lambda
448 mode2[0][1] = 2212; // proton
449 mode2[0][2] = -211; // negative pion
451 TVirtualMC::GetMC()->SetDecayMode(1020000020,bratio2,mode2);
453 // Define the 2-body phase space decay for the Anti-H-Dibaryon
457 for (Int_t kz = 0; kz < 6; kz++) {
464 amode2[0][0] = -3122; // anti-deuteron
465 amode2[0][1] = -2212; // anti-proton
466 amode2[0][2] = 211; // positive pion
468 TVirtualMC::GetMC()->SetDecayMode(-1020000020,abratio2,amode2);
470 // Define the 2-body phase space decay for the Xi0P
474 for (Int_t kz = 0; kz < 6; kz++) {
481 mode5[0][0] = 3122; // Lambda
482 mode5[0][1] = 2212; // proton
484 TVirtualMC::GetMC()->SetDecayMode(1030000020,bratio5,mode5);
486 // Define the 2-body phase space decay for the Anti-Xi0P
490 for (Int_t kz = 0; kz < 6; kz++) {
497 amode5[0][0] = -3122; // anti-Lambda
498 amode5[0][1] = -2212; // anti-proton
500 TVirtualMC::GetMC()->SetDecayMode(-1030000020,abratio5,amode5);
502 // Define the 2-body phase space decay for the Lambda-Neutron-Neutron
506 for (Int_t kz = 0; kz < 6; kz++) {
513 mode6[0][0] = 1000010030; // triton
514 mode6[0][1] = -211; // pion
516 TVirtualMC::GetMC()->SetDecayMode(1010000030,bratio6,mode6);
518 // Define the 2-body phase space decay for the Anti-Lambda-Neutron-Neutron
522 for (Int_t kz = 0; kz < 6; kz++) {
529 amode6[0][0] = -1000010030; // anti-triton
530 amode6[0][1] = 211; // pion
532 TVirtualMC::GetMC()->SetDecayMode(-1010000030,abratio6,amode6);
534 ///////////////////////////////////////////////////////////////////
536 // Define the 2-body phase space decay for the f0(980)
538 // Float_t bratio[6];
540 for (Int_t kz = 0; kz < 6; kz++) {
547 mode[0][0] = 211; // pion
548 mode[0][1] = -211; // pion
550 TVirtualMC::GetMC()->SetDecayMode(9010221,bratio,mode);
552 // Define the 2-body phase space decay for the f2(1270)
554 // Float_t bratio[6];
556 for (Int_t kz = 0; kz < 6; kz++) {
563 mode[0][0] = 211; // pion
564 mode[0][1] = -211; // pion
566 TVirtualMC::GetMC()->SetDecayMode(225,bratio,mode);
568 // Lambda1520/Lambda1520bar
570 TVirtualMC::GetMC()->DefineParticle(3124, "Lambda1520", kPTNeutron, 1.5195 , 0.0, 4.22e-23,"Hadron", 0.0156, 3, -1, 0, 0, 0, 0, 3, 0, kTRUE);
571 TVirtualMC::GetMC()->DefineParticle(-3124, "Lambda1520bar", kPTNeutron, 1.5195 , 0.0, 4.22e-23,"Hadron", 0.0156, 3, -1, 0, 0, 0, 0, 3, 0, kTRUE);
573 // Lambda1520 decay modes
576 bratio[0] = 0.223547;
581 bratio[1] = 0.223547;
585 // L(1520) -> Sigma+ pi-
586 bratio[2] = 0.139096;
590 // L(1520) -> Sigma0 pi0
591 bratio[3] = 0.139096;
595 // L(1520) -> Sigma- pi+
596 bratio[4] = 0.139096;
600 // The other decay modes are neglected
605 TVirtualMC::GetMC()->SetDecayMode(3124,bratio,mode);
607 // Lambda1520bar decay modes
609 // L(1520)bar -> p- K+
610 bratio[0] = 0.223547;
614 // L(1520)bar -> nbar K0bar
615 bratio[1] = 0.223547;
619 // L(1520)bar -> Sigmabar- pi+
620 bratio[2] = 0.139096;
624 // L(1520)bar -> Sigma0bar pi0
625 bratio[3] = 0.139096;
629 // L(1520)bar -> Sigmabar+ pi-
630 bratio[4] = 0.139096;
634 // The other decay modes are neglected
639 TVirtualMC::GetMC()->SetDecayMode(-3124,bratio,mode);
641 // --------------------------------------------------------------------
644 //_______________________________________________________________________
645 void AliMC::InitGeometry()
648 // Initialize detectors
651 AliInfo("Initialisation:");
653 TIter next(gAlice->Modules());
655 while((detector = dynamic_cast<AliModule*>(next()))) {
658 AliInfo(Form("%10s R:%.2fs C:%.2fs",
659 detector->GetName(),stw.RealTime(),stw.CpuTime()));
663 //_______________________________________________________________________
664 void AliMC::SetGeometryFromCDB()
666 // Set the loading of geometry from cdb instead of creating it
667 // A default CDB storage needs to be set before this method is called
668 if(AliCDBManager::Instance()->IsDefaultStorageSet() &&
669 AliCDBManager::Instance()->GetRun() >= 0)
670 AliSimulation::Instance()->SetGeometryFile("*OCDB*");
672 AliError("Loading of geometry from CDB ignored. First set a default CDB storage!");
675 //_______________________________________________________________________
676 Bool_t AliMC::IsGeometryFromCDB() const
678 return (strcmp(AliSimulation::Instance()->GetGeometryFile(),"*OCDB*")==0);
681 //_______________________________________________________________________
682 void AliMC::SetAllAlignableVolumes()
685 // Add alignable volumes (TGeoPNEntries) looping on all
689 AliInfo(Form("Setting entries for all alignable volumes of active detectors"));
691 TIter next(gAlice->Modules());
692 while((detector = dynamic_cast<AliModule*>(next()))) {
693 detector->AddAlignableVolumes();
697 //_______________________________________________________________________
698 void AliMC::GeneratePrimaries()
701 // Generate primary particles and fill them in the stack.
704 Generator()->Generate();
707 //_______________________________________________________________________
708 void AliMC::SetGenerator(AliGenerator *generator)
711 // Load the event generator
713 if(!fGenerator) fGenerator = generator;
716 //_______________________________________________________________________
717 void AliMC::ResetGenerator(AliGenerator *generator)
720 // Load the event generator
724 AliWarning(Form("Replacing generator %s with %s",
725 fGenerator->GetName(),generator->GetName()));
728 AliWarning(Form("Replacing generator %s with NULL",
729 fGenerator->GetName()));
732 fGenerator = generator;
735 //_______________________________________________________________________
736 void AliMC::FinishRun()
738 // Clean generator information
739 AliDebug(1, "fGenerator->FinishRun()");
740 fGenerator->FinishRun();
742 // Monitoring information
745 fMonitor->Export("timing.root");
748 //Output energy summary tables
749 AliDebug(1, "EnergySummary()");
750 ToAliDebug(1, EnergySummary());
753 //_______________________________________________________________________
754 void AliMC::BeginPrimary()
757 // Called at the beginning of each primary track
762 ResetTrackReferences();
765 //_______________________________________________________________________
766 void AliMC::PreTrack()
768 // Actions before the track's transport
770 //verbose.PreTrack();
772 TObjArray &dets = *gAlice->Modules();
775 for(Int_t i=0; i<=gAlice->GetNdets(); i++)
776 if((module = dynamic_cast<AliModule*>(dets[i])))
780 //_______________________________________________________________________
781 void AliMC::Stepping()
784 // Called at every step during transport
786 //verbose.Stepping();
788 Int_t id = DetFromMate(TVirtualMC::GetMC()->CurrentMedium());
792 if ( TVirtualMC::GetMC()->IsNewTrack() &&
793 TVirtualMC::GetMC()->TrackTime() == 0. &&
795 fRDecayMax > fRDecayMin &&
796 TVirtualMC::GetMC()->TrackPid() == fDecayPdg )
798 FixParticleDecaytime();
801 // --- If monitoring timing was requested, monitor the step
802 if (fUseMonitoring) {
804 fMonitor = new AliTransportMonitor(TVirtualMC::GetMC()->NofVolumes()+1);
807 if (TVirtualMC::GetMC()->IsNewTrack() || TVirtualMC::GetMC()->TrackTime() == 0. || TVirtualMC::GetMC()->TrackStep()<1.1E-10) {
808 fMonitor->DummyStep();
812 Int_t volId = TVirtualMC::GetMC()->CurrentVolID(copy);
813 Int_t pdg = TVirtualMC::GetMC()->TrackPid();
814 TLorentzVector xyz, pxpypz;
815 TVirtualMC::GetMC()->TrackPosition(xyz);
816 TVirtualMC::GetMC()->TrackMomentum(pxpypz);
817 fMonitor->StepInfo(volId, pdg, pxpypz.E(), xyz.X(), xyz.Y(), xyz.Z());
821 // --- If lego option, do it and leave
822 if (AliSimulation::Instance()->Lego())
823 AliSimulation::Instance()->Lego()->StepManager();
826 //Update energy deposition tables
827 AddEnergyDeposit(TVirtualMC::GetMC()->CurrentVolID(copy),TVirtualMC::GetMC()->Edep());
829 // write tracke reference for track which is dissapearing - MI
831 if (TVirtualMC::GetMC()->IsTrackDisappeared() && !(TVirtualMC::GetMC()->IsTrackAlive())) {
832 if (TVirtualMC::GetMC()->Etot() > 0.05) AddTrackReference(GetCurrentTrackNumber(),
833 AliTrackReference::kDisappeared);
838 //Call the appropriate stepping routine;
839 AliModule *det = dynamic_cast<AliModule*>(gAlice->Modules()->At(id));
840 if(det && det->StepManagerIsEnabled()) {
846 //_______________________________________________________________________
847 void AliMC::EnergySummary()
850 // Print summary of deposited energy
856 Int_t kn, i, left, j, id;
857 const Float_t kzero=0;
858 Int_t ievent=AliRunLoader::Instance()->GetHeader()->GetEvent()+1;
860 // Energy loss information
862 printf("***************** Energy Loss Information per event (GEV) *****************\n");
863 for(kn=1;kn<fEventEnergy.GetSize();kn++) {
866 fEventEnergy[ndep]=kn;
871 ed2=100*TMath::Sqrt(TMath::Max(ed2-ed*ed,kzero))/ed;
874 fSummEnergy[ndep]=ed;
875 fSum2Energy[ndep]=TMath::Min(static_cast<Float_t>(99.),TMath::Max(ed2,kzero));
880 for(kn=0;kn<(ndep-1)/3+1;kn++) {
882 for(i=0;i<(3<left?3:left);i++) {
884 id=Int_t (fEventEnergy[j]+0.1);
885 printf(" %s %10.3f +- %10.3f%%;",TVirtualMC::GetMC()->VolName(id),fSummEnergy[j],fSum2Energy[j]);
890 // Relative energy loss in different detectors
891 printf("******************** Relative Energy Loss per event ********************\n");
892 printf("Total energy loss per event %10.3f GeV\n",edtot);
893 for(kn=0;kn<(ndep-1)/5+1;kn++) {
895 for(i=0;i<(5<left?5:left);i++) {
897 id=Int_t (fEventEnergy[j]+0.1);
898 printf(" %s %10.3f%%;",TVirtualMC::GetMC()->VolName(id),100*fSummEnergy[j]/edtot);
902 for(kn=0;kn<75;kn++) printf("*");
906 // Reset the TArray's
907 // fEventEnergy.Set(0);
908 // fSummEnergy.Set(0);
909 // fSum2Energy.Set(0);
912 //_____________________________________________________________________________
913 void AliMC::BeginEvent()
916 // Clean-up previous event
918 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
919 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
920 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
921 AliDebug(1, " BEGINNING EVENT ");
922 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
923 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
924 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
926 AliRunLoader *runloader=AliRunLoader::Instance();
928 /*******************************/
929 /* Clean after eventual */
931 /*******************************/
934 //Set the next event in Run Loader -> Cleans trees (TreeK and all trees in detectors),
935 gAlice->SetEventNrInRun(gAlice->GetEventNrInRun()+1);
936 runloader->SetEventNumber(gAlice->GetEventNrInRun());// sets new files, cleans the previous event stuff, if necessary, etc.,
937 AliDebug(1, Form("EventNr is %d",gAlice->GetEventNrInRun()));
939 fEventEnergy.Reset();
940 // Clean detector information
942 if (runloader->Stack())
943 runloader->Stack()->Reset();//clean stack -> tree is unloaded
945 runloader->MakeStack();//or make a new one
947 // Random engine status
950 if ( fSaveRndmStatus || fSaveRndmEventStatus) {
951 TString fileName="random";
952 if ( fSaveRndmEventStatus ) {
954 fileName += gAlice->GetEventNrInRun();
958 // write ROOT random engine status
959 cout << "Saving random engine status in " << fileName.Data() << endl;
960 TFile f(fileName.Data(),"RECREATE");
961 gRandom->Write(fileName.Data());
964 if ( fReadRndmStatus ) {
965 //read ROOT random engine status
966 cout << "Reading random engine status from " << fRndmFileName.Data() << endl;
967 TFile f(fRndmFileName.Data());
968 gRandom->Read(fRndmFileName.Data());
971 if(AliSimulation::Instance()->Lego() == 0x0)
973 AliDebug(1, "fRunLoader->MakeTree(K)");
974 runloader->MakeTree("K");
977 AliDebug(1, "TVirtualMC::GetMC()->SetStack(fRunLoader->Stack())");
978 TVirtualMC::GetMC()->SetStack(runloader->Stack());//Was in InitMC - but was moved here
979 //because we don't have guarantee that
980 //stack pointer is not going to change from event to event
981 //since it bellobgs to header and is obtained via RunLoader
983 // Reset all Detectors & kinematics & make/reset trees
986 runloader->GetHeader()->Reset(AliCDBManager::Instance()->GetRun(),gAlice->GetEvNumber(),
987 gAlice->GetEventNrInRun());
988 // fRunLoader->WriteKinematics("OVERWRITE"); is there any reason to rewrite here since MakeTree does so
990 if(AliSimulation::Instance()->Lego())
992 AliSimulation::Instance()->Lego()->BeginEvent();
997 AliDebug(1, "ResetHits()");
1000 AliDebug(1, "fRunLoader->MakeTree(H)");
1001 runloader->MakeTree("H");
1005 MakeTmpTrackRefsTree();
1006 //create new branches and SetAdresses
1007 TIter next(gAlice->Modules());
1008 AliModule *detector;
1009 while((detector = (AliModule*)next()))
1011 AliDebug(2, Form("%s->MakeBranch(H)",detector->GetName()));
1012 detector->MakeBranch("H");
1016 //_______________________________________________________________________
1017 void AliMC::ResetHits()
1020 // Reset all Detectors hits
1022 TIter next(gAlice->Modules());
1023 AliModule *detector;
1024 while((detector = dynamic_cast<AliModule*>(next()))) {
1025 detector->ResetHits();
1029 //_______________________________________________________________________
1030 void AliMC::ResetDigits()
1033 // Reset all Detectors digits
1035 TIter next(gAlice->Modules());
1036 AliModule *detector;
1037 while((detector = dynamic_cast<AliModule*>(next()))) {
1038 detector->ResetDigits();
1042 //_______________________________________________________________________
1043 void AliMC::ResetSDigits()
1046 // Reset all Detectors digits
1048 TIter next(gAlice->Modules());
1049 AliModule *detector;
1050 while((detector = dynamic_cast<AliModule*>(next()))) {
1051 detector->ResetSDigits();
1055 //_______________________________________________________________________
1056 void AliMC::PostTrack()
1058 // Posts tracks for each module
1060 TObjArray &dets = *gAlice->Modules();
1063 for(Int_t i=0; i<=gAlice->GetNdets(); i++)
1064 if((module = dynamic_cast<AliModule*>(dets[i])))
1065 module->PostTrack();
1068 //_______________________________________________________________________
1069 void AliMC::FinishPrimary()
1072 // Called at the end of each primary track
1075 AliRunLoader *runloader=AliRunLoader::Instance();
1076 // static Int_t count=0;
1077 // const Int_t times=10;
1078 // This primary is finished, purify stack
1079 #if ROOT_VERSION_CODE > 262152
1080 if (!(TVirtualMC::GetMC()->SecondariesAreOrdered())) {
1081 if (runloader->Stack()->ReorderKine()) RemapHits();
1084 if (runloader->Stack()->PurifyKine()) RemapHits();
1086 TIter next(gAlice->Modules());
1087 AliModule *detector;
1088 while((detector = dynamic_cast<AliModule*>(next()))) {
1089 detector->FinishPrimary();
1090 AliLoader* loader = detector->GetLoader();
1093 TTree* treeH = loader->TreeH();
1094 if (treeH) treeH->Fill(); //can be Lego run and treeH can not exist
1098 // Write out track references if any
1099 if (fTmpTreeTR) fTmpTreeTR->Fill();
1102 void AliMC::RemapHits()
1105 // Remaps the track labels of the hits
1106 AliRunLoader *runloader=AliRunLoader::Instance();
1107 AliStack* stack = runloader->Stack();
1108 TList* hitLists = GetHitLists();
1109 TIter next(hitLists);
1110 TCollection *hitList;
1112 while((hitList = dynamic_cast<TCollection*>(next()))) {
1113 TIter nexthit(hitList);
1115 while((hit = dynamic_cast<AliHit*>(nexthit()))) {
1116 hit->SetTrack(stack->TrackLabel(hit->GetTrack()));
1121 // This for detectors which have a special mapping mechanism
1122 // for hits, such as TPC and TRD
1126 TObjArray* modules = gAlice->Modules();
1127 TIter nextmod(modules);
1129 while((module = (AliModule*) nextmod())) {
1130 AliDetector* det = dynamic_cast<AliDetector*> (module);
1131 if (det) det->RemapTrackHitIDs(stack->TrackLabelMap());
1134 RemapTrackReferencesIDs(stack->TrackLabelMap());
1137 //_______________________________________________________________________
1138 void AliMC::FinishEvent()
1141 // Called at the end of the event.
1144 if(AliSimulation::Instance()->Lego()) AliSimulation::Instance()->Lego()->FinishEvent();
1146 TIter next(gAlice->Modules());
1147 AliModule *detector;
1148 while((detector = dynamic_cast<AliModule*>(next()))) {
1149 detector->FinishEvent();
1152 //Update the energy deposit tables
1154 for(i=0;i<fEventEnergy.GetSize();i++)
1156 fSummEnergy[i]+=fEventEnergy[i];
1157 fSum2Energy[i]+=fEventEnergy[i]*fEventEnergy[i];
1160 AliRunLoader *runloader=AliRunLoader::Instance();
1162 AliHeader* header = runloader->GetHeader();
1163 AliStack* stack = runloader->Stack();
1164 if ( (header == 0x0) || (stack == 0x0) )
1165 {//check if we got header and stack. If not cry and exit aliroot
1166 AliFatal("Can not get the stack or header from LOADER");
1167 return;//never reached
1169 // Update Header information
1170 header->SetNprimary(stack->GetNprimary());
1171 header->SetNtrack(stack->GetNtrack());
1172 header->SetTimeStamp(AliSimulation::Instance()->GenerateTimeStamp());
1174 // Write out the kinematics
1175 if (!AliSimulation::Instance()->Lego()) stack->FinishEvent();
1177 // Synchronize the TreeTR with TreeK
1178 if (fTmpTreeTR) ReorderAndExpandTreeTR();
1180 // Write out the event Header information
1181 TTree* treeE = runloader->TreeE();
1184 header->SetStack(stack);
1189 AliError("Can not get TreeE from RL");
1192 if(AliSimulation::Instance()->Lego() == 0x0)
1194 runloader->WriteKinematics("OVERWRITE");
1195 runloader->WriteTrackRefs("OVERWRITE");
1196 runloader->WriteHits("OVERWRITE");
1199 AliDebug(1, "<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<");
1200 AliDebug(1, "<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<");
1201 AliDebug(1, "<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<");
1202 AliDebug(1, " FINISHING EVENT ");
1203 AliDebug(1, "<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<");
1204 AliDebug(1, "<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<");
1205 AliDebug(1, "<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<");
1208 //_______________________________________________________________________
1211 // MC initialization
1214 //=================Create Materials and geometry
1215 TVirtualMC::GetMC()->Init();
1216 // Set alignable volumes for the whole geometry (with old root)
1217 #if ROOT_VERSION_CODE < 331527
1218 SetAllAlignableVolumes();
1220 //Read the cuts for all materials
1222 //Build the special IMEDIA table
1225 //Compute cross-sections
1226 TVirtualMC::GetMC()->BuildPhysics();
1228 //Initialise geometry deposition table
1229 fEventEnergy.Set(TVirtualMC::GetMC()->NofVolumes()+1);
1230 fSummEnergy.Set(TVirtualMC::GetMC()->NofVolumes()+1);
1231 fSum2Energy.Set(TVirtualMC::GetMC()->NofVolumes()+1);
1233 // Register MC in configuration
1234 AliConfig::Instance()->Add(TVirtualMC::GetMC());
1237 //_______________________________________________________________________
1238 void AliMC::MediaTable()
1241 // Built media table to get from the media number to
1245 Int_t kz, nz, idt, lz, i, k, ind;
1247 TObjArray &dets = *gAlice->Detectors();
1249 Int_t ndets=gAlice->GetNdets();
1251 // For all detectors
1252 for (kz=0;kz<ndets;kz++) {
1253 // If detector is defined
1254 if((det=dynamic_cast<AliModule*>(dets[kz]))) {
1255 TArrayI &idtmed = *(det->GetIdtmed());
1256 for(nz=0;nz<100;nz++) {
1258 // Find max and min material number
1259 if((idt=idtmed[nz])) {
1260 det->LoMedium() = det->LoMedium() < idt ? det->LoMedium() : idt;
1261 det->HiMedium() = det->HiMedium() > idt ? det->HiMedium() : idt;
1264 if(det->LoMedium() > det->HiMedium()) {
1265 det->LoMedium() = 0;
1266 det->HiMedium() = 0;
1268 if(det->HiMedium() > fImedia->GetSize()) {
1269 AliError(Form("Increase fImedia from %d to %d",
1270 fImedia->GetSize(),det->HiMedium()));
1273 // Tag all materials in rage as belonging to detector kz
1274 for(lz=det->LoMedium(); lz<= det->HiMedium(); lz++) {
1281 // Print summary table
1282 AliInfo("Tracking media ranges:");
1284 for(i=0;i<(ndets-1)/6+1;i++) {
1285 for(k=0;k< (6<ndets-i*6?6:ndets-i*6);k++) {
1287 det=dynamic_cast<AliModule*>(dets[ind]);
1289 printf(" %6s: %3d -> %3d;",det->GetName(),det->LoMedium(),
1292 printf(" %6s: %3d -> %3d;","NULL",0,0);
1299 //_______________________________________________________________________
1300 void AliMC::ReadTransPar()
1303 // Read filename to set the transport parameters
1307 const Int_t kncuts=10;
1308 const Int_t knflags=12;
1309 const Int_t knpars=kncuts+knflags;
1310 const char kpars[knpars][7] = {"CUTGAM" ,"CUTELE","CUTNEU","CUTHAD","CUTMUO",
1311 "BCUTE","BCUTM","DCUTE","DCUTM","PPCUTM","ANNI",
1312 "BREM","COMP","DCAY","DRAY","HADR","LOSS",
1313 "MULS","PAIR","PHOT","RAYL","STRA"};
1317 Float_t cut[kncuts];
1318 Int_t flag[knflags];
1319 Int_t i, itmed, iret, jret, ktmed, kz;
1322 // See whether the file is there
1323 filtmp=gSystem->ExpandPathName(fTransParName.Data());
1324 lun=fopen(filtmp,"r");
1327 AliWarning(Form("File %s does not exist!",fTransParName.Data()));
1332 // Initialise cuts and flags
1333 for(i=0;i<kncuts;i++) cut[i]=-99;
1334 for(i=0;i<knflags;i++) flag[i]=-99;
1337 // Read up to the end of line excluded
1338 iret=fscanf(lun,"%255[^\n]",line);
1344 // Read the end of line
1345 jret = fscanf(lun,"%*c");
1347 if(line[0]=='*') continue;
1349 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",
1350 detName,&itmed,&cut[0],&cut[1],&cut[2],&cut[3],&cut[4],&cut[5],&cut[6],&cut[7],&cut[8],
1351 &cut[9],&flag[0],&flag[1],&flag[2],&flag[3],&flag[4],&flag[5],&flag[6],&flag[7],
1352 &flag[8],&flag[9],&flag[10],&flag[11]);
1356 AliWarning(Form("Error reading file %s",fTransParName.Data()));
1359 // Check that the module exist
1360 AliModule *mod = gAlice->GetModule(detName);
1362 // Get the array of media numbers
1363 TArrayI &idtmed = *mod->GetIdtmed();
1364 // Check that the tracking medium code is valid
1365 if(0<=itmed && itmed < 100) {
1366 ktmed=idtmed[itmed];
1368 AliWarning(Form("Invalid tracking medium code %d for %s",itmed,mod->GetName()));
1371 // Set energy thresholds
1372 for(kz=0;kz<kncuts;kz++) {
1374 AliDebug(2, Form("%-6s set to %10.3E for tracking medium code %4d for %s",
1375 kpars[kz],cut[kz],itmed,mod->GetName()));
1376 TVirtualMC::GetMC()->Gstpar(ktmed,kpars[kz],cut[kz]);
1379 // Set transport mechanisms
1380 for(kz=0;kz<knflags;kz++) {
1382 AliDebug(2, Form("%-6s set to %10d for tracking medium code %4d for %s",
1383 kpars[kncuts+kz],flag[kz],itmed,mod->GetName()));
1384 TVirtualMC::GetMC()->Gstpar(ktmed,kpars[kncuts+kz],Float_t(flag[kz]));
1388 AliWarning(Form("Invalid medium code %d",itmed));
1392 AliDebug(1, Form("%s not present",detName));
1398 //_______________________________________________________________________
1399 void AliMC::SetTransPar(const char *filename)
1402 // Sets the file name for transport parameters
1404 fTransParName = filename;
1407 //_______________________________________________________________________
1408 void AliMC::AddHit(Int_t id, Int_t track, Int_t *vol, Float_t *hits) const
1411 // Add a hit to detector id
1413 TObjArray &dets = *gAlice->Modules();
1414 if(dets[id]) static_cast<AliModule*>(dets[id])->AddHit(track,vol,hits);
1417 //_______________________________________________________________________
1418 void AliMC::AddDigit(Int_t id, Int_t *tracks, Int_t *digits) const
1421 // Add digit to detector id
1423 TObjArray &dets = *gAlice->Modules();
1424 if(dets[id]) static_cast<AliModule*>(dets[id])->AddDigit(tracks,digits);
1427 //_______________________________________________________________________
1428 Int_t AliMC::GetCurrentTrackNumber() const {
1430 // Returns current track
1432 return AliRunLoader::Instance()->Stack()->GetCurrentTrackNumber();
1435 //_______________________________________________________________________
1436 void AliMC::DumpPart (Int_t i) const
1439 // Dumps particle i in the stack
1441 AliRunLoader * runloader = AliRunLoader::Instance();
1442 if (runloader->Stack())
1443 runloader->Stack()->DumpPart(i);
1446 //_______________________________________________________________________
1447 void AliMC::DumpPStack () const
1450 // Dumps the particle stack
1452 AliRunLoader * runloader = AliRunLoader::Instance();
1453 if (runloader->Stack())
1454 runloader->Stack()->DumpPStack();
1457 //_______________________________________________________________________
1458 Int_t AliMC::GetNtrack() const {
1460 // Returns number of tracks in stack
1463 AliRunLoader * runloader = AliRunLoader::Instance();
1464 if (runloader->Stack())
1465 ntracks = runloader->Stack()->GetNtrack();
1469 //_______________________________________________________________________
1470 Int_t AliMC::GetPrimary(Int_t track) const
1473 // return number of primary that has generated track
1475 Int_t nprimary = -999;
1476 AliRunLoader * runloader = AliRunLoader::Instance();
1477 if (runloader->Stack())
1478 nprimary = runloader->Stack()->GetPrimary(track);
1482 //_______________________________________________________________________
1483 TParticle* AliMC::Particle(Int_t i) const
1485 // Returns the i-th particle from the stack taking into account
1486 // the remaping done by PurifyKine
1487 AliRunLoader * runloader = AliRunLoader::Instance();
1489 if (runloader->Stack())
1490 return runloader->Stack()->Particle(i);
1494 //_______________________________________________________________________
1495 const TObjArray* AliMC::Particles() const {
1497 // Returns pointer to Particles array
1499 AliRunLoader * runloader = AliRunLoader::Instance();
1501 if (runloader->Stack())
1502 return runloader->Stack()->Particles();
1506 //_______________________________________________________________________
1507 void AliMC::PushTrack(Int_t done, Int_t parent, Int_t pdg, const Float_t *pmom,
1508 const Float_t *vpos, const Float_t *polar, Float_t tof,
1509 TMCProcess mech, Int_t &ntr, Float_t weight, Int_t is) const
1511 // Delegate to stack
1513 AliRunLoader * runloader = AliRunLoader::Instance();
1515 if (runloader->Stack())
1516 runloader->Stack()->PushTrack(done, parent, pdg, pmom, vpos, polar, tof,
1517 mech, ntr, weight, is);
1520 //_______________________________________________________________________
1521 void AliMC::PushTrack(Int_t done, Int_t parent, Int_t pdg,
1522 Double_t px, Double_t py, Double_t pz, Double_t e,
1523 Double_t vx, Double_t vy, Double_t vz, Double_t tof,
1524 Double_t polx, Double_t poly, Double_t polz,
1525 TMCProcess mech, Int_t &ntr, Float_t weight, Int_t is) const
1527 // Delegate to stack
1529 AliRunLoader * runloader = AliRunLoader::Instance();
1531 if (runloader->Stack())
1532 runloader->Stack()->PushTrack(done, parent, pdg, px, py, pz, e, vx, vy, vz, tof,
1533 polx, poly, polz, mech, ntr, weight, is);
1536 //_______________________________________________________________________
1537 void AliMC::SetHighWaterMark(Int_t nt) const
1540 // Set high water mark for last track in event
1541 AliRunLoader * runloader = AliRunLoader::Instance();
1543 if (runloader->Stack())
1544 runloader->Stack()->SetHighWaterMark(nt);
1547 //_______________________________________________________________________
1548 void AliMC::KeepTrack(Int_t track) const
1551 // Delegate to stack
1553 AliRunLoader * runloader = AliRunLoader::Instance();
1555 if (runloader->Stack())
1556 runloader->Stack()->KeepTrack(track);
1559 //_______________________________________________________________________
1560 void AliMC::FlagTrack(Int_t track) const
1562 // Delegate to stack
1564 AliRunLoader * runloader = AliRunLoader::Instance();
1566 if (runloader->Stack())
1567 runloader->Stack()->FlagTrack(track);
1570 //_______________________________________________________________________
1571 void AliMC::SetCurrentTrack(Int_t track) const
1574 // Set current track number
1576 AliRunLoader * runloader = AliRunLoader::Instance();
1578 if (runloader->Stack())
1579 runloader->Stack()->SetCurrentTrack(track);
1582 //_______________________________________________________________________
1583 AliTrackReference* AliMC::AddTrackReference(Int_t label, Int_t id)
1586 // add a trackrefernce to the list
1587 Int_t primary = GetPrimary(label);
1588 Particle(primary)->SetBit(kKeepBit);
1590 Int_t nref = fTmpTrackReferences.GetEntriesFast();
1591 return new(fTmpTrackReferences[nref]) AliTrackReference(label, id);
1596 //_______________________________________________________________________
1597 void AliMC::ResetTrackReferences()
1600 // Reset all references
1602 fTmpTrackReferences.Clear();
1605 //_______________________________________________________________________
1606 void AliMC::RemapTrackReferencesIDs(const Int_t *map)
1609 // Remapping track reference
1610 // Called at finish primary
1613 Int_t nEntries = fTmpTrackReferences.GetEntries();
1614 for (Int_t i=0; i < nEntries; i++){
1615 AliTrackReference * ref = dynamic_cast<AliTrackReference*>(fTmpTrackReferences.UncheckedAt(i));
1617 Int_t newID = map[ref->GetTrack()];
1618 if (newID>=0) ref->SetTrack(newID);
1620 ref->SetBit(kNotDeleted,kFALSE);
1621 fTmpTrackReferences.RemoveAt(i);
1625 fTmpTrackReferences.Compress();
1628 //_______________________________________________________________________
1629 void AliMC::FixParticleDecaytime()
1632 // Fix the particle decay time according to rmin and rmax for decays
1636 TVirtualMC::GetMC()->TrackMomentum(p);
1637 Double_t tmin, tmax;
1640 // Transverse velocity
1641 Double_t vt = p.Pt() / p.E();
1643 if ((b = ((AliMagF*)TGeoGlobalMagField::Instance()->GetField())->SolenoidField()) > 0.) { // [kG]
1647 Double_t rho = p.Pt() / 0.0003 / b; // [cm]
1649 // Revolution frequency
1651 Double_t omega = vt / rho;
1653 // Maximum and minimum decay time
1655 // Check for curlers first
1656 const Double_t kOvRhoSqr2 = 1./(rho*TMath::Sqrt(2.));
1657 if (fRDecayMax * kOvRhoSqr2 > 1.) return;
1661 tmax = TMath::ACos((1.-fRDecayMax*kOvRhoSqr2)*(1.+fRDecayMax*kOvRhoSqr2)) / omega; // [ct]
1662 tmin = TMath::ACos((1.-fRDecayMin*kOvRhoSqr2)*(1.+fRDecayMin*kOvRhoSqr2)) / omega; // [ct]
1664 tmax = fRDecayMax / vt; // [ct]
1665 tmin = fRDecayMin / vt; // [ct]
1669 // Dial t using the two limits
1670 Double_t t = tmin + (tmax - tmin) * gRandom->Rndm(); // [ct]
1673 // Force decay time in transport code
1675 TVirtualMC::GetMC()->ForceDecayTime(t / 2.99792458e10);
1678 void AliMC::MakeTmpTrackRefsTree()
1680 // Make the temporary track reference tree
1681 fTmpFileTR = new TFile("TrackRefsTmp.root", "recreate");
1682 fTmpTreeTR = new TTree("TreeTR", "Track References");
1683 TClonesArray* pRef = &fTmpTrackReferences;
1684 fTmpTreeTR->Branch("TrackReferences", &pRef, 4000);
1687 //_______________________________________________________________________
1688 void AliMC::ReorderAndExpandTreeTR()
1691 // Reorder and expand the temporary track reference tree in order to match the kinematics tree
1694 AliRunLoader *rl = AliRunLoader::Instance();
1697 AliDebug(1, "fRunLoader->MakeTrackRefsContainer()");
1698 rl->MakeTrackRefsContainer();
1699 TTree * treeTR = rl->TreeTR();
1700 // make branch for central track references
1701 TClonesArray* pRef = &fTrackReferences;
1702 treeTR->Branch("TrackReferences", &pRef);
1704 AliStack* stack = rl->Stack();
1705 Int_t np = stack->GetNprimary();
1706 Int_t nt = fTmpTreeTR->GetEntries();
1708 // Loop over tracks and find the secondaries with the help of the kine tree
1711 for (Int_t ip = np - 1; ip > -1; ip--) {
1712 TParticle *part = stack->Particle(ip);
1713 //printf("Particle %5d %5d %5d %5d %5d \n", ip, part->GetPdgCode(), part->GetFirstMother(), part->GetFirstDaughter(), part->GetLastDaughter());
1715 // Skip primaries that have not been transported
1716 Int_t dau1 = part->GetFirstDaughter();
1718 if (!part->TestBit(kTransportBit)) continue;
1720 fTmpTreeTR->GetEntry(it++);
1721 Int_t nh = fTmpTrackReferences.GetEntries();
1722 // Determine range of secondaries produced by this primary
1724 Int_t inext = ip - 1;
1727 part = stack->Particle(inext);
1728 dau2 = part->GetFirstDaughter();
1729 if (!(part->TestBit(kTransportBit)) || dau2 == -1 || dau2 < np) {
1730 // if (dau2 == -1 || dau2 < np) {
1736 dau2 = stack->GetNtrack() - 1;
1739 } // find upper bound
1741 // printf("Check (1) %5d %5d %5d %5d %5d \n", ip, np, it, dau1, dau2);
1743 // Loop over reference hits and find secondary label
1744 for (Int_t id = dau1; (id <= dau2) && (dau1 > -1); id++) {
1745 for (Int_t ih = 0; ih < nh; ih++) {
1746 AliTrackReference* tr = (AliTrackReference*) fTmpTrackReferences.At(ih);
1747 Int_t label = tr->Label();
1749 if (label == ip) continue;
1750 if (label > dau2 || label < dau1)
1751 AliWarning(Form("Track Reference Label out of range !: %5d %5d %5d \n", label, dau1, dau2));
1754 Int_t nref = fTrackReferences.GetEntriesFast();
1755 new(fTrackReferences[nref]) AliTrackReference(*tr);
1759 fTrackReferences.Clear();
1764 // Now loop again and write the primaries
1766 for (Int_t ip = 0; ip < np; ip++) {
1767 TParticle* part = stack->Particle(ip);
1768 // if ((part->GetFirstDaughter() == -1 && part->GetStatusCode() <= 1) || part->GetFirstDaughter() >= np)
1769 if (part->TestBit(kTransportBit))
1771 // Skip particles that have not been transported
1772 fTmpTreeTR->GetEntry(it--);
1773 Int_t nh = fTmpTrackReferences.GetEntries();
1775 // Loop over reference hits and find primary labels
1776 for (Int_t ih = 0; ih < nh; ih++) {
1777 AliTrackReference* tr = (AliTrackReference*) fTmpTrackReferences.At(ih);
1778 Int_t label = tr->Label();
1780 Int_t nref = fTrackReferences.GetEntriesFast();
1781 new(fTrackReferences[nref]) AliTrackReference(*tr);
1786 fTrackReferences.Clear();
1790 if (ifills != stack->GetNtrack())
1791 AliWarning(Form("Number of entries in TreeTR (%5d) unequal to TreeK (%5d) \n", ifills, stack->GetNtrack()));
1795 fTmpFileTR->Close();
1797 fTmpTrackReferences.Clear();
1798 gSystem->Exec("rm -rf TrackRefsTmp.root");