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
25 #include <TStopwatch.h>
27 #include <TVirtualMC.h>
29 #include "AliDetector.h"
30 #include "AliGenerator.h"
31 #include "AliHeader.h"
37 #include "AliTrackReference.h"
42 //_______________________________________________________________________
61 //_______________________________________________________________________
62 AliMC::AliMC(const char *name, const char *title) :
63 TVirtualMCApplication(name, title),
71 fImedia(new TArrayI(1000)),
74 fHitLists(new TList()),
75 fTrackReferences(new TClonesArray("AliTrackReference", 100))
78 // Set transport parameters
81 // Prepare the tracking medium lists
82 for(Int_t i=0;i<1000;i++) (*fImedia)[i]=-99;
86 //_______________________________________________________________________
87 AliMC::AliMC(const AliMC &mc) :
88 TVirtualMCApplication(mc),
103 // Copy constructor for AliRun
108 //_______________________________________________________________________
116 // Delete track references
117 if (fTrackReferences) {
118 fTrackReferences->Delete();
119 delete fTrackReferences;
120 fTrackReferences = 0;
125 //_______________________________________________________________________
126 void AliMC::Copy(TObject &) const
128 //dummy Copy function
129 Fatal("Copy","Not implemented!\n");
132 //_______________________________________________________________________
133 void AliMC::ConstructGeometry()
136 // Create modules, materials, geometry
140 TIter next(gAlice->Modules());
142 if (GetDebug()) Info("ConstructGeometry","Geometry creation:");
143 while((detector = dynamic_cast<AliModule*>(next()))) {
145 // Initialise detector materials and geometry
146 detector->CreateMaterials();
147 detector->CreateGeometry();
148 printf("%10s R:%.2fs C:%.2fs\n",
149 detector->GetName(),stw.RealTime(),stw.CpuTime());
153 //_______________________________________________________________________
154 void AliMC::InitGeometry()
157 // Initialize detectors and display geometry
160 printf("Initialisation:\n");
162 TIter next(gAlice->Modules());
164 while((detector = dynamic_cast<AliModule*>(next()))) {
166 // Initialise detector and display geometry
168 detector->BuildGeometry();
169 printf("%10s R:%.2fs C:%.2fs\n",
170 detector->GetName(),stw.RealTime(),stw.CpuTime());
175 //_______________________________________________________________________
176 void AliMC::GeneratePrimaries()
179 // Generate primary particles and fill them in the stack.
182 Generator()->Generate();
185 //_______________________________________________________________________
186 void AliMC::SetGenerator(AliGenerator *generator)
189 // Load the event generator
191 if(!fGenerator) fGenerator = generator;
194 //_______________________________________________________________________
195 void AliMC::ResetGenerator(AliGenerator *generator)
198 // Load the event generator
202 Warning("ResetGenerator","Replacing generator %s with %s\n",
203 fGenerator->GetName(),generator->GetName());
205 Warning("ResetGenerator","Replacing generator %s with NULL\n",
206 fGenerator->GetName());
207 fGenerator = generator;
210 //_______________________________________________________________________
211 void AliMC::FinishRun()
213 // Clean generator information
214 if (GetDebug()) Info("FinishRun"," fGenerator->FinishRun()");
215 fGenerator->FinishRun();
217 //Output energy summary tables
219 Info("FinishRun"," EnergySummary()");
225 //_______________________________________________________________________
226 void AliMC::BeginPrimary()
229 // Called at the beginning of each primary track
234 ResetTrackReferences();
238 //_______________________________________________________________________
239 void AliMC::PreTrack()
241 // Actions before the track's transport
242 TObjArray &dets = *gAlice->Modules();
245 for(Int_t i=0; i<=gAlice->GetNdets(); i++)
246 if((module = dynamic_cast<AliModule*>(dets[i])))
252 //_______________________________________________________________________
253 void AliMC::Stepping()
256 // Called at every step during transport
259 Int_t id = DetFromMate(gMC->GetMedium());
263 // --- If lego option, do it and leave
265 gAlice->Lego()->StepManager();
268 //Update energy deposition tables
269 AddEnergyDeposit(gMC->CurrentVolID(copy),gMC->Edep());
271 // write tracke reference for track which is dissapearing - MI
272 if (gMC->IsTrackDisappeared()) {
273 if (gMC->Etot()>0.05) AddTrackReference(GetCurrentTrackNumber());
276 //Call the appropriate stepping routine;
277 AliModule *det = dynamic_cast<AliModule*>(gAlice->Modules()->At(id));
278 if(det && det->StepManagerIsEnabled()) {
279 fMCQA->StepManager(id);
285 //_______________________________________________________________________
286 void AliMC::EnergySummary()
289 // Print summary of deposited energy
295 Int_t kn, i, left, j, id;
296 const Float_t kzero=0;
297 Int_t ievent=gAlice->GetRunLoader()->GetHeader()->GetEvent()+1;
299 // Energy loss information
301 printf("***************** Energy Loss Information per event (GEV) *****************\n");
302 for(kn=1;kn<fEventEnergy.GetSize();kn++) {
305 fEventEnergy[ndep]=kn;
310 ed2=100*TMath::Sqrt(TMath::Max(ed2-ed*ed,kzero))/ed;
313 fSummEnergy[ndep]=ed;
314 fSum2Energy[ndep]=TMath::Min(static_cast<Float_t>(99.),TMath::Max(ed2,kzero));
319 for(kn=0;kn<(ndep-1)/3+1;kn++) {
321 for(i=0;i<(3<left?3:left);i++) {
323 id=Int_t (fEventEnergy[j]+0.1);
324 printf(" %s %10.3f +- %10.3f%%;",gMC->VolName(id),fSummEnergy[j],fSum2Energy[j]);
329 // Relative energy loss in different detectors
330 printf("******************** Relative Energy Loss per event ********************\n");
331 printf("Total energy loss per event %10.3f GeV\n",edtot);
332 for(kn=0;kn<(ndep-1)/5+1;kn++) {
334 for(i=0;i<(5<left?5:left);i++) {
336 id=Int_t (fEventEnergy[j]+0.1);
337 printf(" %s %10.3f%%;",gMC->VolName(id),100*fSummEnergy[j]/edtot);
341 for(kn=0;kn<75;kn++) printf("*");
345 // Reset the TArray's
346 // fEventEnergy.Set(0);
347 // fSummEnergy.Set(0);
348 // fSum2Energy.Set(0);
351 //_____________________________________________________________________________
352 void AliMC::BeginEvent()
355 // Clean-up previous event
359 Info("BeginEvent",">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
360 Info("BeginEvent",">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
361 Info("BeginEvent",">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
362 Info("BeginEvent"," BEGINNING EVENT ");
363 Info("BeginEvent",">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
364 Info("BeginEvent",">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
365 Info("BeginEvent",">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
368 AliRunLoader *runloader=gAlice->GetRunLoader();
370 /*******************************/
371 /* Clean after eventual */
373 /*******************************/
376 //Set the next event in Run Loader -> Cleans trees (TreeK and all trees in detectors),
377 gAlice->SetEventNrInRun(gAlice->GetEventNrInRun()+1);
378 runloader->SetEventNumber(gAlice->GetEventNrInRun());// sets new files, cleans the previous event stuff, if necessary, etc.,
379 if (GetDebug()) Info("BeginEvent","EventNr is %d",gAlice->GetEventNrInRun());
381 fEventEnergy.Reset();
382 // Clean detector information
384 if (runloader->Stack())
385 runloader->Stack()->Reset();//clean stack -> tree is unloaded
387 runloader->MakeStack();//or make a new one
389 if(gAlice->Lego() == 0x0)
391 if (GetDebug()) Info("BeginEvent"," fRunLoader->MakeTree(K)");
392 runloader->MakeTree("K");
395 if (GetDebug()) Info("BeginEvent"," gMC->SetStack(fRunLoader->Stack())");
396 gMC->SetStack(gAlice->GetRunLoader()->Stack());//Was in InitMC - but was moved here
397 //because we don't have guarantee that
398 //stack pointer is not going to change from event to event
399 //since it bellobgs to header and is obtained via RunLoader
401 // Reset all Detectors & kinematics & make/reset trees
404 runloader->GetHeader()->Reset(gAlice->GetRunNumber(),gAlice->GetEvNumber(),
405 gAlice->GetEventNrInRun());
406 // fRunLoader->WriteKinematics("OVERWRITE"); is there any reason to rewrite here since MakeTree does so
410 gAlice->Lego()->BeginEvent();
415 if (GetDebug()) Info("BeginEvent"," ResetHits()");
418 if (GetDebug()) Info("BeginEvent"," fRunLoader->MakeTree(H)");
419 runloader->MakeTree("H");
421 if (GetDebug()) Info("BeginEvent"," fRunLoader->MakeTrackRefsContainer()");
422 runloader->MakeTrackRefsContainer();//for insurance
425 //create new branches and SetAdresses
426 TIter next(gAlice->Modules());
428 while((detector = (AliModule*)next()))
430 if (GetDebug()) Info("BeginEvent"," %s->MakeBranch(H)",detector->GetName());
431 detector->MakeBranch("H");
432 if (GetDebug()) Info("BeginEvent"," %s->MakeBranchTR()",detector->GetName());
433 detector->MakeBranchTR();
434 if (GetDebug()) Info("BeginEvent"," %s->SetTreeAddress()",detector->GetName());
435 detector->SetTreeAddress();
437 // make branch for AliRun track References
438 TTree * treeTR = runloader->TreeTR();
440 // make branch for central track references
441 if (!fTrackReferences) fTrackReferences = new TClonesArray("AliTrackReference",0);
443 branch = treeTR->Branch("AliRun",&fTrackReferences);
444 branch->SetAddress(&fTrackReferences);
449 //_______________________________________________________________________
450 void AliMC::ResetHits()
453 // Reset all Detectors hits
455 TIter next(gAlice->Modules());
457 while((detector = dynamic_cast<AliModule*>(next()))) {
458 detector->ResetHits();
462 //_______________________________________________________________________
463 void AliMC::PostTrack()
465 // Posts tracks for each module
466 TObjArray &dets = *gAlice->Modules();
469 for(Int_t i=0; i<=gAlice->GetNdets(); i++)
470 if((module = dynamic_cast<AliModule*>(dets[i])))
474 //_______________________________________________________________________
475 void AliMC::FinishPrimary()
478 // Called at the end of each primary track
480 AliRunLoader *runloader=gAlice->GetRunLoader();
481 // static Int_t count=0;
482 // const Int_t times=10;
483 // This primary is finished, purify stack
484 runloader->Stack()->PurifyKine();
486 TIter next(gAlice->Modules());
488 while((detector = dynamic_cast<AliModule*>(next()))) {
489 detector->FinishPrimary();
490 AliLoader* loader = detector->GetLoader();
493 TTree* treeH = loader->TreeH();
494 if (treeH) treeH->Fill(); //can be Lego run and treeH can not exist
498 // Write out track references if any
499 if (runloader->TreeTR()) runloader->TreeTR()->Fill();
502 //_______________________________________________________________________
503 void AliMC::FinishEvent()
506 // Called at the end of the event.
510 if(gAlice->Lego()) gAlice->Lego()->FinishEvent();
512 TIter next(gAlice->Modules());
514 while((detector = dynamic_cast<AliModule*>(next()))) {
515 detector->FinishEvent();
518 //Update the energy deposit tables
520 for(i=0;i<fEventEnergy.GetSize();i++)
522 fSummEnergy[i]+=fEventEnergy[i];
523 fSum2Energy[i]+=fEventEnergy[i]*fEventEnergy[i];
526 AliRunLoader *runloader=gAlice->GetRunLoader();
528 AliHeader* header = runloader->GetHeader();
529 AliStack* stack = runloader->Stack();
530 if ( (header == 0x0) || (stack == 0x0) )
531 {//check if we got header and stack. If not cry and exit aliroot
532 Fatal("AliRun","Can not get the stack or header from LOADER");
533 return;//never reached
535 // Update Header information
536 header->SetNprimary(stack->GetNprimary());
537 header->SetNtrack(stack->GetNtrack());
540 // Write out the kinematics
541 stack->FinishEvent();
543 // Write out the event Header information
544 TTree* treeE = runloader->TreeE();
547 header->SetStack(stack);
552 Error("FinishEvent","Can not get TreeE from RL");
555 if(gAlice->Lego() == 0x0)
557 runloader->WriteKinematics("OVERWRITE");
558 runloader->WriteTrackRefs("OVERWRITE");
559 runloader->WriteHits("OVERWRITE");
564 Info("FinishEvent","<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<");
565 Info("FinishEvent","<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<");
566 Info("FinishEvent","<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<");
567 Info("FinishEvent"," FINISHING EVENT ");
568 Info("FinishEvent","<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<");
569 Info("FinishEvent","<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<");
570 Info("FinishEvent","<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<");
574 //_______________________________________________________________________
575 void AliMC::Field(const Double_t* x, Double_t* b) const
577 // Calculates field "b" at point "x"
581 //_______________________________________________________________________
586 //=================Create Materials and geometry
588 //Read the cuts for all materials
590 //Build the special IMEDIA table
593 //Compute cross-sections
596 //Write Geometry object to current file.
597 gAlice->GetRunLoader()->WriteGeometry();
599 //Initialise geometry deposition table
600 fEventEnergy.Set(gMC->NofVolumes()+1);
601 fSummEnergy.Set(gMC->NofVolumes()+1);
602 fSum2Energy.Set(gMC->NofVolumes()+1);
605 fMCQA = new AliMCQA(gAlice->GetNdets());
607 // Register MC in configuration
608 AliConfig::Instance()->Add(gMC);
612 //_______________________________________________________________________
613 void AliMC::MediaTable()
616 // Built media table to get from the media number to
620 Int_t kz, nz, idt, lz, i, k, ind;
622 TObjArray &dets = *gAlice->Detectors();
624 Int_t ndets=gAlice->GetNdets();
627 for (kz=0;kz<ndets;kz++) {
628 // If detector is defined
629 if((det=dynamic_cast<AliModule*>(dets[kz]))) {
630 TArrayI &idtmed = *(det->GetIdtmed());
631 for(nz=0;nz<100;nz++) {
632 // Find max and min material number
633 if((idt=idtmed[nz])) {
634 det->LoMedium() = det->LoMedium() < idt ? det->LoMedium() : idt;
635 det->HiMedium() = det->HiMedium() > idt ? det->HiMedium() : idt;
638 if(det->LoMedium() > det->HiMedium()) {
642 if(det->HiMedium() > fImedia->GetSize()) {
643 Error("MediaTable","Increase fImedia from %d to %d",
644 fImedia->GetSize(),det->HiMedium());
647 // Tag all materials in rage as belonging to detector kz
648 for(lz=det->LoMedium(); lz<= det->HiMedium(); lz++) {
655 // Print summary table
656 printf(" Traking media ranges:\n");
657 for(i=0;i<(ndets-1)/6+1;i++) {
658 for(k=0;k< (6<ndets-i*6?6:ndets-i*6);k++) {
660 det=dynamic_cast<AliModule*>(dets[ind]);
662 printf(" %6s: %3d -> %3d;",det->GetName(),det->LoMedium(),
665 printf(" %6s: %3d -> %3d;","NULL",0,0);
671 //_______________________________________________________________________
672 void AliMC::ReadTransPar()
675 // Read filename to set the transport parameters
679 const Int_t kncuts=10;
680 const Int_t knflags=11;
681 const Int_t knpars=kncuts+knflags;
682 const char kpars[knpars][7] = {"CUTGAM" ,"CUTELE","CUTNEU","CUTHAD","CUTMUO",
683 "BCUTE","BCUTM","DCUTE","DCUTM","PPCUTM","ANNI",
684 "BREM","COMP","DCAY","DRAY","HADR","LOSS",
685 "MULS","PAIR","PHOT","RAYL"};
691 Int_t i, itmed, iret, ktmed, kz;
694 // See whether the file is there
695 filtmp=gSystem->ExpandPathName(fTransParName.Data());
696 lun=fopen(filtmp,"r");
699 Warning("ReadTransPar","File %s does not exist!\n",fTransParName.Data());
704 printf(" "); for(i=0;i<60;i++) printf("*"); printf("\n");
705 printf(" *%59s\n","*");
706 printf(" * Please check carefully what you are doing!%10s\n","*");
707 printf(" *%59s\n","*");
711 // Initialise cuts and flags
712 for(i=0;i<kncuts;i++) cut[i]=-99;
713 for(i=0;i<knflags;i++) flag[i]=-99;
715 for(i=0;i<256;i++) line[i]='\0';
716 // Read up to the end of line excluded
717 iret=fscanf(lun,"%[^\n]",line);
722 printf(" *%59s\n","*");
723 printf(" "); for(i=0;i<60;i++) printf("*"); printf("\n");
727 // Read the end of line
730 if(line[0]=='*') continue;
732 iret=sscanf(line,"%s %d %f %f %f %f %f %f %f %f %f %f %d %d %d %d %d %d %d %d %d %d %d",
733 detName,&itmed,&cut[0],&cut[1],&cut[2],&cut[3],&cut[4],&cut[5],&cut[6],&cut[7],&cut[8],
734 &cut[9],&flag[0],&flag[1],&flag[2],&flag[3],&flag[4],&flag[5],&flag[6],&flag[7],
735 &flag[8],&flag[9],&flag[10]);
739 Warning("ReadTransPar","Error reading file %s\n",fTransParName.Data());
742 // Check that the module exist
743 AliModule *mod = gAlice->GetModule(detName);
745 // Get the array of media numbers
746 TArrayI &idtmed = *mod->GetIdtmed();
747 // Check that the tracking medium code is valid
748 if(0<=itmed && itmed < 100) {
751 Warning("ReadTransPar","Invalid tracking medium code %d for %s\n",itmed,mod->GetName());
754 // Set energy thresholds
755 for(kz=0;kz<kncuts;kz++) {
757 if(fDebug) printf(" * %-6s set to %10.3E for tracking medium code %4d for %s\n",
758 kpars[kz],cut[kz],itmed,mod->GetName());
759 gMC->Gstpar(ktmed,kpars[kz],cut[kz]);
762 // Set transport mechanisms
763 for(kz=0;kz<knflags;kz++) {
765 if(fDebug) printf(" * %-6s set to %10d for tracking medium code %4d for %s\n",
766 kpars[kncuts+kz],flag[kz],itmed,mod->GetName());
767 gMC->Gstpar(ktmed,kpars[kncuts+kz],Float_t(flag[kz]));
771 Warning("ReadTransPar","Invalid medium code %d *\n",itmed);
775 if(fDebug) printf("%s::ReadTransParModule: %s not present\n",ClassName(),detName);
781 //_______________________________________________________________________
782 void AliMC::SetTransPar(const char *filename)
785 // Sets the file name for transport parameters
787 fTransParName = filename;
790 //_______________________________________________________________________
791 void AliMC::Browse(TBrowser *b) const
794 // Called when the item "Run" is clicked on the left pane
795 // of the Root browser.
796 // It displays the Root Trees and all detectors.
798 //detectors are in folders anyway
799 b->Add(fMCQA,"AliMCQA");
803 //_______________________________________________________________________
804 void AliMC::AddHit(Int_t id, Int_t track, Int_t *vol, Float_t *hits) const
807 // Add a hit to detector id
809 TObjArray &dets = *gAlice->Modules();
810 if(dets[id]) dynamic_cast<AliModule*>(dets[id])->AddHit(track,vol,hits);
813 //_______________________________________________________________________
814 void AliMC::AddDigit(Int_t id, Int_t *tracks, Int_t *digits) const
817 // Add digit to detector id
819 TObjArray &dets = *gAlice->Modules();
820 if(dets[id]) dynamic_cast<AliModule*>(dets[id])->AddDigit(tracks,digits);
823 //_______________________________________________________________________
824 Int_t AliMC::GetCurrentTrackNumber() const {
826 // Returns current track
828 return gAlice->GetRunLoader()->Stack()->GetCurrentTrackNumber();
831 //_______________________________________________________________________
832 void AliMC::DumpPart (Int_t i) const
835 // Dumps particle i in the stack
837 AliRunLoader * runloader = gAlice->GetRunLoader();
838 if (runloader->Stack())
839 runloader->Stack()->DumpPart(i);
842 //_______________________________________________________________________
843 void AliMC::DumpPStack () const
846 // Dumps the particle stack
848 AliRunLoader * runloader = gAlice->GetRunLoader();
849 if (runloader->Stack())
850 runloader->Stack()->DumpPStack();
853 //_______________________________________________________________________
854 Int_t AliMC::GetNtrack() const {
856 // Returns number of tracks in stack
859 AliRunLoader * runloader = gAlice->GetRunLoader();
860 if (runloader->Stack())
861 ntracks = runloader->Stack()->GetNtrack();
865 //_______________________________________________________________________
866 Int_t AliMC::GetPrimary(Int_t track) const
869 // return number of primary that has generated track
871 Int_t nprimary = -999;
872 AliRunLoader * runloader = gAlice->GetRunLoader();
873 if (runloader->Stack())
874 nprimary = runloader->Stack()->GetPrimary(track);
878 //_______________________________________________________________________
879 TParticle* AliMC::Particle(Int_t i) const
881 // Returns the i-th particle from the stack taking into account
882 // the remaping done by PurifyKine
883 AliRunLoader * runloader = gAlice->GetRunLoader();
885 if (runloader->Stack())
886 return runloader->Stack()->Particle(i);
890 //_______________________________________________________________________
891 TObjArray* AliMC::Particles() const {
893 // Returns pointer to Particles array
895 AliRunLoader * runloader = gAlice->GetRunLoader();
897 if (runloader->Stack())
898 return runloader->Stack()->Particles();
902 //_______________________________________________________________________
903 void AliMC::PushTrack(Int_t done, Int_t parent, Int_t pdg, Float_t *pmom,
904 Float_t *vpos, Float_t *polar, Float_t tof,
905 TMCProcess mech, Int_t &ntr, Float_t weight, Int_t is) const
909 AliRunLoader * runloader = gAlice->GetRunLoader();
911 if (runloader->Stack())
912 runloader->Stack()->PushTrack(done, parent, pdg, pmom, vpos, polar, tof,
913 mech, ntr, weight, is);
916 //_______________________________________________________________________
917 void AliMC::PushTrack(Int_t done, Int_t parent, Int_t pdg,
918 Double_t px, Double_t py, Double_t pz, Double_t e,
919 Double_t vx, Double_t vy, Double_t vz, Double_t tof,
920 Double_t polx, Double_t poly, Double_t polz,
921 TMCProcess mech, Int_t &ntr, Float_t weight, Int_t is) const
925 AliRunLoader * runloader = gAlice->GetRunLoader();
927 if (runloader->Stack())
928 runloader->Stack()->PushTrack(done, parent, pdg, px, py, pz, e, vx, vy, vz, tof,
929 polx, poly, polz, mech, ntr, weight, is);
932 //_______________________________________________________________________
933 void AliMC::SetHighWaterMark(Int_t nt) const
936 // Set high water mark for last track in event
937 AliRunLoader * runloader = gAlice->GetRunLoader();
939 if (runloader->Stack())
940 runloader->Stack()->SetHighWaterMark(nt);
943 //_______________________________________________________________________
944 void AliMC::KeepTrack(Int_t track) const
949 AliRunLoader * runloader = gAlice->GetRunLoader();
951 if (runloader->Stack())
952 runloader->Stack()->KeepTrack(track);
955 //_______________________________________________________________________
956 void AliMC::FlagTrack(Int_t track) const
960 AliRunLoader * runloader = gAlice->GetRunLoader();
962 if (runloader->Stack())
963 runloader->Stack()->FlagTrack(track);
966 //_______________________________________________________________________
967 void AliMC::SetCurrentTrack(Int_t track) const
970 // Set current track number
972 AliRunLoader * runloader = gAlice->GetRunLoader();
974 if (runloader->Stack())
975 runloader->Stack()->SetCurrentTrack(track);
978 //_______________________________________________________________________
979 void AliMC::AddTrackReference(Int_t label){
981 // add a trackrefernce to the list
982 if (!fTrackReferences) {
983 cerr<<"Container trackrefernce not active\n";
986 Int_t nref = fTrackReferences->GetEntriesFast();
987 TClonesArray &lref = *fTrackReferences;
988 new(lref[nref]) AliTrackReference(label);
993 //_______________________________________________________________________
994 void AliMC::ResetTrackReferences()
997 // Reset all references
999 if (fTrackReferences) fTrackReferences->Clear();
1001 TIter next(gAlice->Modules());
1002 AliModule *detector;
1003 while((detector = dynamic_cast<AliModule*>(next()))) {
1004 detector->ResetTrackReferences();
1008 void AliMC::RemapTrackReferencesIDs(Int_t *map)
1011 // Remapping track reference
1012 // Called at finish primary
1014 if (!fTrackReferences) return;
1015 for (Int_t i=0;i<fTrackReferences->GetEntries();i++){
1016 AliTrackReference * ref = dynamic_cast<AliTrackReference*>(fTrackReferences->UncheckedAt(i));
1018 Int_t newID = map[ref->GetTrack()];
1019 if (newID>=0) ref->SetTrack(newID);
1021 //ref->SetTrack(-1);
1022 ref->SetBit(kNotDeleted,kFALSE);
1023 fTrackReferences->RemoveAt(i);
1027 fTrackReferences->Compress();