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
26 #include <TStopwatch.h>
28 #include <TVirtualMC.h>
29 #include <TGeoManager.h>
34 #include "AliDetector.h"
35 #include "AliGenerator.h"
36 #include "AliHeader.h"
43 #include "AliTrackReference.h"
48 //_______________________________________________________________________
70 //_______________________________________________________________________
71 AliMC::AliMC(const char *name, const char *title) :
72 TVirtualMCApplication(name, title),
82 fImedia(new TArrayI(1000)),
85 fHitLists(new TList()),
86 fTrackReferences(new TClonesArray("AliTrackReference", 100))
89 // Set transport parameters
92 // Prepare the tracking medium lists
93 for(Int_t i=0;i<1000;i++) (*fImedia)[i]=-99;
96 //_______________________________________________________________________
97 AliMC::AliMC(const AliMC &mc) :
98 TVirtualMCApplication(mc),
115 // Copy constructor for AliMC
120 //_______________________________________________________________________
128 // Delete track references
129 if (fTrackReferences) {
130 fTrackReferences->Delete();
131 delete fTrackReferences;
132 fTrackReferences = 0;
137 //_______________________________________________________________________
138 void AliMC::Copy(TObject &) const
140 //dummy Copy function
141 AliFatal("Not implemented!");
144 //_______________________________________________________________________
145 void AliMC::ConstructGeometry()
148 // Either load geometry from file or create it through usual
149 // loop on detectors. In the first case the method
150 // AliModule::CreateMaterials() only builds fIdtmed and is postponed
151 // at InitGeometry().
154 if(gAlice->IsRootGeometry()){
156 const char *geomfilename = gAlice->GetGeometryFileName();
157 if(gSystem->ExpandPathName(geomfilename)){
158 AliInfo(Form("Loading geometry from file:\n %40s\n\n",geomfilename));
159 TGeoManager::Import(geomfilename);
161 AliInfo(Form("Geometry file %40s not found!\n",geomfilename));
165 // Create modules, materials, geometry
167 TIter next(gAlice->Modules());
169 AliDebug(1, "Geometry creation:");
170 while((detector = dynamic_cast<AliModule*>(next()))) {
172 // Initialise detector materials and geometry
173 detector->CreateMaterials();
174 detector->CreateGeometry();
175 AliInfo(Form("%10s R:%.2fs C:%.2fs",
176 detector->GetName(),stw.RealTime(),stw.CpuTime()));
182 //_______________________________________________________________________
183 void AliMC::InitGeometry()
186 // Initialize detectors and display geometry
189 AliInfo("Initialisation:");
191 TIter next(gAlice->Modules());
193 while((detector = dynamic_cast<AliModule*>(next()))) {
195 // Initialise detector and display geometry
196 if(gAlice->IsRootGeometry()) detector->CreateMaterials();
198 detector->BuildGeometry();
199 AliInfo(Form("%10s R:%.2fs C:%.2fs",
200 detector->GetName(),stw.RealTime(),stw.CpuTime()));
205 //_______________________________________________________________________
206 void AliMC::GeneratePrimaries()
209 // Generate primary particles and fill them in the stack.
212 Generator()->Generate();
215 //_______________________________________________________________________
216 void AliMC::SetGenerator(AliGenerator *generator)
219 // Load the event generator
221 if(!fGenerator) fGenerator = generator;
224 //_______________________________________________________________________
225 void AliMC::ResetGenerator(AliGenerator *generator)
228 // Load the event generator
232 AliWarning(Form("Replacing generator %s with %s",
233 fGenerator->GetName(),generator->GetName()))
235 AliWarning(Form("Replacing generator %s with NULL",
236 fGenerator->GetName()));
237 fGenerator = generator;
240 //_______________________________________________________________________
241 void AliMC::FinishRun()
243 // Clean generator information
244 AliDebug(1, "fGenerator->FinishRun()");
245 fGenerator->FinishRun();
247 //Output energy summary tables
248 AliDebug(1, "EnergySummary()");
249 ToAliDebug(1, EnergySummary());
252 //_______________________________________________________________________
253 void AliMC::BeginPrimary()
256 // Called at the beginning of each primary track
261 ResetTrackReferences();
265 //_______________________________________________________________________
266 void AliMC::PreTrack()
268 // Actions before the track's transport
269 TObjArray &dets = *gAlice->Modules();
272 for(Int_t i=0; i<=gAlice->GetNdets(); i++)
273 if((module = dynamic_cast<AliModule*>(dets[i])))
279 //_______________________________________________________________________
280 void AliMC::Stepping()
283 // Called at every step during transport
286 Int_t id = DetFromMate(gMC->GetMedium());
290 if ( gMC->IsNewTrack() &&
291 gMC->TrackTime() == 0. &&
293 fRDecayMax > fRDecayMin &&
294 gMC->TrackPid() == fDecayPdg )
296 FixParticleDecaytime();
302 // --- If lego option, do it and leave
304 gAlice->Lego()->StepManager();
307 //Update energy deposition tables
308 AddEnergyDeposit(gMC->CurrentVolID(copy),gMC->Edep());
310 // write tracke reference for track which is dissapearing - MI
311 if (gMC->IsTrackDisappeared()) {
312 if (gMC->Etot()>0.05) AddTrackReference(GetCurrentTrackNumber());
315 //Call the appropriate stepping routine;
316 AliModule *det = dynamic_cast<AliModule*>(gAlice->Modules()->At(id));
317 if(det && det->StepManagerIsEnabled()) {
318 fMCQA->StepManager(id);
324 //_______________________________________________________________________
325 void AliMC::EnergySummary()
328 // Print summary of deposited energy
334 Int_t kn, i, left, j, id;
335 const Float_t kzero=0;
336 Int_t ievent=gAlice->GetRunLoader()->GetHeader()->GetEvent()+1;
338 // Energy loss information
340 printf("***************** Energy Loss Information per event (GEV) *****************\n");
341 for(kn=1;kn<fEventEnergy.GetSize();kn++) {
344 fEventEnergy[ndep]=kn;
349 ed2=100*TMath::Sqrt(TMath::Max(ed2-ed*ed,kzero))/ed;
352 fSummEnergy[ndep]=ed;
353 fSum2Energy[ndep]=TMath::Min(static_cast<Float_t>(99.),TMath::Max(ed2,kzero));
358 for(kn=0;kn<(ndep-1)/3+1;kn++) {
360 for(i=0;i<(3<left?3:left);i++) {
362 id=Int_t (fEventEnergy[j]+0.1);
363 printf(" %s %10.3f +- %10.3f%%;",gMC->VolName(id),fSummEnergy[j],fSum2Energy[j]);
368 // Relative energy loss in different detectors
369 printf("******************** Relative Energy Loss per event ********************\n");
370 printf("Total energy loss per event %10.3f GeV\n",edtot);
371 for(kn=0;kn<(ndep-1)/5+1;kn++) {
373 for(i=0;i<(5<left?5:left);i++) {
375 id=Int_t (fEventEnergy[j]+0.1);
376 printf(" %s %10.3f%%;",gMC->VolName(id),100*fSummEnergy[j]/edtot);
380 for(kn=0;kn<75;kn++) printf("*");
384 // Reset the TArray's
385 // fEventEnergy.Set(0);
386 // fSummEnergy.Set(0);
387 // fSum2Energy.Set(0);
390 //_____________________________________________________________________________
391 void AliMC::BeginEvent()
394 // Clean-up previous event
396 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
397 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
398 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
399 AliDebug(1, " BEGINNING EVENT ");
400 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
401 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
402 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
404 AliRunLoader *runloader=gAlice->GetRunLoader();
406 /*******************************/
407 /* Clean after eventual */
409 /*******************************/
412 //Set the next event in Run Loader -> Cleans trees (TreeK and all trees in detectors),
413 gAlice->SetEventNrInRun(gAlice->GetEventNrInRun()+1);
414 runloader->SetEventNumber(gAlice->GetEventNrInRun());// sets new files, cleans the previous event stuff, if necessary, etc.,
415 AliDebug(1, Form("EventNr is %d",gAlice->GetEventNrInRun()));
417 fEventEnergy.Reset();
418 // Clean detector information
420 if (runloader->Stack())
421 runloader->Stack()->Reset();//clean stack -> tree is unloaded
423 runloader->MakeStack();//or make a new one
425 if(gAlice->Lego() == 0x0)
427 AliDebug(1, "fRunLoader->MakeTree(K)");
428 runloader->MakeTree("K");
431 AliDebug(1, "gMC->SetStack(fRunLoader->Stack())");
432 gMC->SetStack(gAlice->GetRunLoader()->Stack());//Was in InitMC - but was moved here
433 //because we don't have guarantee that
434 //stack pointer is not going to change from event to event
435 //since it bellobgs to header and is obtained via RunLoader
437 // Reset all Detectors & kinematics & make/reset trees
440 runloader->GetHeader()->Reset(gAlice->GetRunNumber(),gAlice->GetEvNumber(),
441 gAlice->GetEventNrInRun());
442 // fRunLoader->WriteKinematics("OVERWRITE"); is there any reason to rewrite here since MakeTree does so
446 gAlice->Lego()->BeginEvent();
451 AliDebug(1, "ResetHits()");
454 AliDebug(1, "fRunLoader->MakeTree(H)");
455 runloader->MakeTree("H");
457 AliDebug(1, "fRunLoader->MakeTrackRefsContainer()");
458 runloader->MakeTrackRefsContainer();//for insurance
461 //create new branches and SetAdresses
462 TIter next(gAlice->Modules());
464 while((detector = (AliModule*)next()))
466 AliDebug(2, Form("%s->MakeBranch(H)",detector->GetName()));
467 detector->MakeBranch("H");
468 AliDebug(2, Form("%s->MakeBranchTR()",detector->GetName()));
469 detector->MakeBranchTR();
470 AliDebug(2, Form("%s->SetTreeAddress()",detector->GetName()));
471 detector->SetTreeAddress();
473 // make branch for AliRun track References
474 TTree * treeTR = runloader->TreeTR();
476 // make branch for central track references
477 if (!fTrackReferences) fTrackReferences = new TClonesArray("AliTrackReference",0);
479 branch = treeTR->Branch("AliRun",&fTrackReferences);
480 branch->SetAddress(&fTrackReferences);
485 //_______________________________________________________________________
486 void AliMC::ResetHits()
489 // Reset all Detectors hits
491 TIter next(gAlice->Modules());
493 while((detector = dynamic_cast<AliModule*>(next()))) {
494 detector->ResetHits();
498 //_______________________________________________________________________
499 void AliMC::PostTrack()
501 // Posts tracks for each module
502 TObjArray &dets = *gAlice->Modules();
505 for(Int_t i=0; i<=gAlice->GetNdets(); i++)
506 if((module = dynamic_cast<AliModule*>(dets[i])))
510 //_______________________________________________________________________
511 void AliMC::FinishPrimary()
514 // Called at the end of each primary track
516 AliRunLoader *runloader=gAlice->GetRunLoader();
517 // static Int_t count=0;
518 // const Int_t times=10;
519 // This primary is finished, purify stack
520 #if ROOT_VERSION_CODE > 262152
521 if (!(gMC->SecondariesAreOrdered()))
522 runloader->Stack()->ReorderKine();
524 runloader->Stack()->PurifyKine();
526 TIter next(gAlice->Modules());
528 while((detector = dynamic_cast<AliModule*>(next()))) {
529 detector->FinishPrimary();
530 AliLoader* loader = detector->GetLoader();
533 TTree* treeH = loader->TreeH();
534 if (treeH) treeH->Fill(); //can be Lego run and treeH can not exist
538 // Write out track references if any
539 if (runloader->TreeTR()) runloader->TreeTR()->Fill();
542 //_______________________________________________________________________
543 void AliMC::FinishEvent()
546 // Called at the end of the event.
550 if(gAlice->Lego()) gAlice->Lego()->FinishEvent();
552 TIter next(gAlice->Modules());
554 while((detector = dynamic_cast<AliModule*>(next()))) {
555 detector->FinishEvent();
558 //Update the energy deposit tables
560 for(i=0;i<fEventEnergy.GetSize();i++)
562 fSummEnergy[i]+=fEventEnergy[i];
563 fSum2Energy[i]+=fEventEnergy[i]*fEventEnergy[i];
566 AliRunLoader *runloader=gAlice->GetRunLoader();
568 AliHeader* header = runloader->GetHeader();
569 AliStack* stack = runloader->Stack();
570 if ( (header == 0x0) || (stack == 0x0) )
571 {//check if we got header and stack. If not cry and exit aliroot
572 AliFatal("Can not get the stack or header from LOADER");
573 return;//never reached
575 // Update Header information
576 header->SetNprimary(stack->GetNprimary());
577 header->SetNtrack(stack->GetNtrack());
580 // Write out the kinematics
581 stack->FinishEvent();
583 // Write out the event Header information
584 TTree* treeE = runloader->TreeE();
587 header->SetStack(stack);
592 AliError("Can not get TreeE from RL");
595 if(gAlice->Lego() == 0x0)
597 runloader->WriteKinematics("OVERWRITE");
598 runloader->WriteTrackRefs("OVERWRITE");
599 runloader->WriteHits("OVERWRITE");
602 AliDebug(1, "<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<");
603 AliDebug(1, "<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<");
604 AliDebug(1, "<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<");
605 AliDebug(1, " FINISHING EVENT ");
606 AliDebug(1, "<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<");
607 AliDebug(1, "<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<");
608 AliDebug(1, "<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<");
611 //_______________________________________________________________________
612 void AliMC::Field(const Double_t* x, Double_t* b) const
614 // Calculates field "b" at point "x"
618 //_______________________________________________________________________
623 //=================Create Materials and geometry
625 //Read the cuts for all materials
627 //Build the special IMEDIA table
630 //Compute cross-sections
633 //Write Geometry object to current file.
634 gAlice->GetRunLoader()->WriteGeometry();
636 //Initialise geometry deposition table
637 fEventEnergy.Set(gMC->NofVolumes()+1);
638 fSummEnergy.Set(gMC->NofVolumes()+1);
639 fSum2Energy.Set(gMC->NofVolumes()+1);
642 fMCQA = new AliMCQA(gAlice->GetNdets());
644 // Register MC in configuration
645 AliConfig::Instance()->Add(gMC);
647 // Export TGeo geometry
648 if (gGeoManager) gGeoManager->Export("geometry.root");
652 //_______________________________________________________________________
653 void AliMC::MediaTable()
656 // Built media table to get from the media number to
660 Int_t kz, nz, idt, lz, i, k, ind;
662 TObjArray &dets = *gAlice->Detectors();
664 Int_t ndets=gAlice->GetNdets();
667 for (kz=0;kz<ndets;kz++) {
668 // If detector is defined
669 if((det=dynamic_cast<AliModule*>(dets[kz]))) {
670 TArrayI &idtmed = *(det->GetIdtmed());
671 for(nz=0;nz<100;nz++) {
673 // Find max and min material number
674 if((idt=idtmed[nz])) {
675 det->LoMedium() = det->LoMedium() < idt ? det->LoMedium() : idt;
676 det->HiMedium() = det->HiMedium() > idt ? det->HiMedium() : idt;
679 if(det->LoMedium() > det->HiMedium()) {
683 if(det->HiMedium() > fImedia->GetSize()) {
684 AliError(Form("Increase fImedia from %d to %d",
685 fImedia->GetSize(),det->HiMedium()));
688 // Tag all materials in rage as belonging to detector kz
689 for(lz=det->LoMedium(); lz<= det->HiMedium(); lz++) {
696 // Print summary table
697 AliInfo("Tracking media ranges:");
699 for(i=0;i<(ndets-1)/6+1;i++) {
700 for(k=0;k< (6<ndets-i*6?6:ndets-i*6);k++) {
702 det=dynamic_cast<AliModule*>(dets[ind]);
704 printf(" %6s: %3d -> %3d;",det->GetName(),det->LoMedium(),
707 printf(" %6s: %3d -> %3d;","NULL",0,0);
714 //_______________________________________________________________________
715 void AliMC::ReadTransPar()
718 // Read filename to set the transport parameters
722 const Int_t kncuts=10;
723 const Int_t knflags=11;
724 const Int_t knpars=kncuts+knflags;
725 const char kpars[knpars][7] = {"CUTGAM" ,"CUTELE","CUTNEU","CUTHAD","CUTMUO",
726 "BCUTE","BCUTM","DCUTE","DCUTM","PPCUTM","ANNI",
727 "BREM","COMP","DCAY","DRAY","HADR","LOSS",
728 "MULS","PAIR","PHOT","RAYL"};
734 Int_t i, itmed, iret, ktmed, kz;
737 // See whether the file is there
738 filtmp=gSystem->ExpandPathName(fTransParName.Data());
739 lun=fopen(filtmp,"r");
742 AliWarning(Form("File %s does not exist!",fTransParName.Data()));
747 // Initialise cuts and flags
748 for(i=0;i<kncuts;i++) cut[i]=-99;
749 for(i=0;i<knflags;i++) flag[i]=-99;
751 for(i=0;i<256;i++) line[i]='\0';
752 // Read up to the end of line excluded
753 iret=fscanf(lun,"%[^\n]",line);
759 // Read the end of line
762 if(line[0]=='*') continue;
764 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",
765 detName,&itmed,&cut[0],&cut[1],&cut[2],&cut[3],&cut[4],&cut[5],&cut[6],&cut[7],&cut[8],
766 &cut[9],&flag[0],&flag[1],&flag[2],&flag[3],&flag[4],&flag[5],&flag[6],&flag[7],
767 &flag[8],&flag[9],&flag[10]);
771 AliWarning(Form("Error reading file %s",fTransParName.Data()));
774 // Check that the module exist
775 AliModule *mod = gAlice->GetModule(detName);
777 // Get the array of media numbers
778 TArrayI &idtmed = *mod->GetIdtmed();
779 // Check that the tracking medium code is valid
780 if(0<=itmed && itmed < 100) {
783 AliWarning(Form("Invalid tracking medium code %d for %s",itmed,mod->GetName()));
786 // Set energy thresholds
787 for(kz=0;kz<kncuts;kz++) {
789 AliDebug(2, Form("%-6s set to %10.3E for tracking medium code %4d for %s",
790 kpars[kz],cut[kz],itmed,mod->GetName()));
791 gMC->Gstpar(ktmed,kpars[kz],cut[kz]);
794 // Set transport mechanisms
795 for(kz=0;kz<knflags;kz++) {
797 AliDebug(2, Form("%-6s set to %10d for tracking medium code %4d for %s",
798 kpars[kncuts+kz],flag[kz],itmed,mod->GetName()));
799 gMC->Gstpar(ktmed,kpars[kncuts+kz],Float_t(flag[kz]));
803 AliWarning(Form("Invalid medium code %d",itmed));
807 AliDebug(1, Form("%s not present",detName));
813 //_______________________________________________________________________
814 void AliMC::SetTransPar(const char *filename)
817 // Sets the file name for transport parameters
819 fTransParName = filename;
822 //_______________________________________________________________________
823 void AliMC::Browse(TBrowser *b)
826 // Called when the item "Run" is clicked on the left pane
827 // of the Root browser.
828 // It displays the Root Trees and all detectors.
830 //detectors are in folders anyway
831 b->Add(fMCQA,"AliMCQA");
835 //_______________________________________________________________________
836 void AliMC::AddHit(Int_t id, Int_t track, Int_t *vol, Float_t *hits) const
839 // Add a hit to detector id
841 TObjArray &dets = *gAlice->Modules();
842 if(dets[id]) dynamic_cast<AliModule*>(dets[id])->AddHit(track,vol,hits);
845 //_______________________________________________________________________
846 void AliMC::AddDigit(Int_t id, Int_t *tracks, Int_t *digits) const
849 // Add digit to detector id
851 TObjArray &dets = *gAlice->Modules();
852 if(dets[id]) dynamic_cast<AliModule*>(dets[id])->AddDigit(tracks,digits);
855 //_______________________________________________________________________
856 Int_t AliMC::GetCurrentTrackNumber() const {
858 // Returns current track
860 return gAlice->GetRunLoader()->Stack()->GetCurrentTrackNumber();
863 //_______________________________________________________________________
864 void AliMC::DumpPart (Int_t i) const
867 // Dumps particle i in the stack
869 AliRunLoader * runloader = gAlice->GetRunLoader();
870 if (runloader->Stack())
871 runloader->Stack()->DumpPart(i);
874 //_______________________________________________________________________
875 void AliMC::DumpPStack () const
878 // Dumps the particle stack
880 AliRunLoader * runloader = gAlice->GetRunLoader();
881 if (runloader->Stack())
882 runloader->Stack()->DumpPStack();
885 //_______________________________________________________________________
886 Int_t AliMC::GetNtrack() const {
888 // Returns number of tracks in stack
891 AliRunLoader * runloader = gAlice->GetRunLoader();
892 if (runloader->Stack())
893 ntracks = runloader->Stack()->GetNtrack();
897 //_______________________________________________________________________
898 Int_t AliMC::GetPrimary(Int_t track) const
901 // return number of primary that has generated track
903 Int_t nprimary = -999;
904 AliRunLoader * runloader = gAlice->GetRunLoader();
905 if (runloader->Stack())
906 nprimary = runloader->Stack()->GetPrimary(track);
910 //_______________________________________________________________________
911 TParticle* AliMC::Particle(Int_t i) const
913 // Returns the i-th particle from the stack taking into account
914 // the remaping done by PurifyKine
915 AliRunLoader * runloader = gAlice->GetRunLoader();
917 if (runloader->Stack())
918 return runloader->Stack()->Particle(i);
922 //_______________________________________________________________________
923 TObjArray* AliMC::Particles() const {
925 // Returns pointer to Particles array
927 AliRunLoader * runloader = gAlice->GetRunLoader();
929 if (runloader->Stack())
930 return runloader->Stack()->Particles();
934 //_______________________________________________________________________
935 void AliMC::PushTrack(Int_t done, Int_t parent, Int_t pdg, Float_t *pmom,
936 Float_t *vpos, Float_t *polar, Float_t tof,
937 TMCProcess mech, Int_t &ntr, Float_t weight, Int_t is) const
941 AliRunLoader * runloader = gAlice->GetRunLoader();
943 if (runloader->Stack())
944 runloader->Stack()->PushTrack(done, parent, pdg, pmom, vpos, polar, tof,
945 mech, ntr, weight, is);
948 //_______________________________________________________________________
949 void AliMC::PushTrack(Int_t done, Int_t parent, Int_t pdg,
950 Double_t px, Double_t py, Double_t pz, Double_t e,
951 Double_t vx, Double_t vy, Double_t vz, Double_t tof,
952 Double_t polx, Double_t poly, Double_t polz,
953 TMCProcess mech, Int_t &ntr, Float_t weight, Int_t is) const
957 AliRunLoader * runloader = gAlice->GetRunLoader();
959 if (runloader->Stack())
960 runloader->Stack()->PushTrack(done, parent, pdg, px, py, pz, e, vx, vy, vz, tof,
961 polx, poly, polz, mech, ntr, weight, is);
964 //_______________________________________________________________________
965 void AliMC::SetHighWaterMark(Int_t nt) const
968 // Set high water mark for last track in event
969 AliRunLoader * runloader = gAlice->GetRunLoader();
971 if (runloader->Stack())
972 runloader->Stack()->SetHighWaterMark(nt);
975 //_______________________________________________________________________
976 void AliMC::KeepTrack(Int_t track) const
981 AliRunLoader * runloader = gAlice->GetRunLoader();
983 if (runloader->Stack())
984 runloader->Stack()->KeepTrack(track);
987 //_______________________________________________________________________
988 void AliMC::FlagTrack(Int_t track) const
992 AliRunLoader * runloader = gAlice->GetRunLoader();
994 if (runloader->Stack())
995 runloader->Stack()->FlagTrack(track);
998 //_______________________________________________________________________
999 void AliMC::SetCurrentTrack(Int_t track) const
1002 // Set current track number
1004 AliRunLoader * runloader = gAlice->GetRunLoader();
1006 if (runloader->Stack())
1007 runloader->Stack()->SetCurrentTrack(track);
1010 //_______________________________________________________________________
1011 void AliMC::AddTrackReference(Int_t label){
1013 // add a trackrefernce to the list
1014 if (!fTrackReferences) {
1015 AliError("Container trackrefernce not active");
1018 Int_t nref = fTrackReferences->GetEntriesFast();
1019 TClonesArray &lref = *fTrackReferences;
1020 new(lref[nref]) AliTrackReference(label);
1025 //_______________________________________________________________________
1026 void AliMC::ResetTrackReferences()
1029 // Reset all references
1031 if (fTrackReferences) fTrackReferences->Clear();
1033 TIter next(gAlice->Modules());
1034 AliModule *detector;
1035 while((detector = dynamic_cast<AliModule*>(next()))) {
1036 detector->ResetTrackReferences();
1040 void AliMC::RemapTrackReferencesIDs(Int_t *map)
1043 // Remapping track reference
1044 // Called at finish primary
1046 if (!fTrackReferences) return;
1047 for (Int_t i=0;i<fTrackReferences->GetEntries();i++){
1048 AliTrackReference * ref = dynamic_cast<AliTrackReference*>(fTrackReferences->UncheckedAt(i));
1050 Int_t newID = map[ref->GetTrack()];
1051 if (newID>=0) ref->SetTrack(newID);
1053 //ref->SetTrack(-1);
1054 ref->SetBit(kNotDeleted,kFALSE);
1055 fTrackReferences->RemoveAt(i);
1059 fTrackReferences->Compress();
1062 void AliMC::FixParticleDecaytime()
1065 // Fix the particle decay time according to rmin and rmax for decays
1069 gMC->TrackMomentum(p);
1070 Double_t tmin, tmax;
1073 // Transverse velocity
1074 Double_t vt = p.Pt() / p.E();
1076 if ((b = gAlice->Field()->SolenoidField()) > 0.) { // [kG]
1080 Double_t rho = p.Pt() / 0.0003 / b; // [cm]
1082 // Revolution frequency
1084 Double_t omega = vt / rho;
1086 // Maximum and minimum decay time
1088 // Check for curlers first
1089 if (fRDecayMax * fRDecayMax / rho / rho / 2. > 1.) return;
1093 tmax = TMath::ACos(1. - fRDecayMax * fRDecayMax / rho / rho / 2.) / omega; // [ct]
1094 tmin = TMath::ACos(1. - fRDecayMin * fRDecayMin / rho / rho / 2.) / omega; // [ct]
1096 tmax = fRDecayMax / vt; // [ct]
1097 tmin = fRDecayMin / vt; // [ct]
1101 // Dial t using the two limits
1102 Double_t t = tmin + (tmax - tmin) * gRandom->Rndm(); // [ct]
1105 // Force decay time in transport code
1107 TGeant3 * geant = (TGeant3*) gMC;
1108 geant->Gcphys()->sumlif = t / p.Beta() / p.Gamma();