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 <TClonesArray.h>
27 #include <TGeoManager.h>
28 #include <TStopwatch.h>
30 #include <TVirtualMC.h>
33 #include "AliDetector.h"
34 #include "AliGenerator.h"
35 #include "AliHeader.h"
42 #include "AliTrackReference.h"
47 //_______________________________________________________________________
69 //_______________________________________________________________________
70 AliMC::AliMC(const char *name, const char *title) :
71 TVirtualMCApplication(name, title),
81 fImedia(new TArrayI(1000)),
84 fHitLists(new TList()),
85 fTrackReferences(new TClonesArray("AliTrackReference", 100))
88 // Set transport parameters
91 // Prepare the tracking medium lists
92 for(Int_t i=0;i<1000;i++) (*fImedia)[i]=-99;
95 //_______________________________________________________________________
96 AliMC::AliMC(const AliMC &mc) :
97 TVirtualMCApplication(mc),
114 // Copy constructor for AliMC
119 //_______________________________________________________________________
127 // Delete track references
128 if (fTrackReferences) {
129 fTrackReferences->Delete();
130 delete fTrackReferences;
131 fTrackReferences = 0;
136 //_______________________________________________________________________
137 void AliMC::Copy(TObject &) const
139 //dummy Copy function
140 AliFatal("Not implemented!");
143 //_______________________________________________________________________
144 void AliMC::ConstructGeometry()
147 // Either load geometry from file or create it through usual
148 // loop on detectors. In the first case the method
149 // AliModule::CreateMaterials() only builds fIdtmed and is postponed
150 // at InitGeometry().
153 if(gAlice->IsRootGeometry()){
155 const char *geomfilename = gAlice->GetGeometryFileName();
156 if(gSystem->ExpandPathName(geomfilename)){
157 AliInfo(Form("Loading geometry from file:\n %40s\n\n",geomfilename));
158 TGeoManager::Import(geomfilename);
160 AliInfo(Form("Geometry file %40s not found!\n",geomfilename));
164 // Create modules, materials, geometry
166 TIter next(gAlice->Modules());
168 AliDebug(1, "Geometry creation:");
169 while((detector = dynamic_cast<AliModule*>(next()))) {
171 // Initialise detector materials and geometry
172 detector->CreateMaterials();
173 detector->CreateGeometry();
174 AliInfo(Form("%10s R:%.2fs C:%.2fs",
175 detector->GetName(),stw.RealTime(),stw.CpuTime()));
181 //_______________________________________________________________________
182 void AliMC::InitGeometry()
185 // Initialize detectors
188 AliInfo("Initialisation:");
190 TIter next(gAlice->Modules());
192 while((detector = dynamic_cast<AliModule*>(next()))) {
194 // Initialise detector geometry
195 if(gAlice->IsRootGeometry()) detector->CreateMaterials();
197 AliInfo(Form("%10s R:%.2fs C:%.2fs",
198 detector->GetName(),stw.RealTime(),stw.CpuTime()));
202 //_______________________________________________________________________
203 void AliMC::SetAllAlignableVolumes()
206 // Add alignable volumes (TGeoPNEntries) looping on all
210 AliInfo(Form("Setting entries for all alignable volumes of active detectors"));
212 TIter next(gAlice->Modules());
213 while((detector = dynamic_cast<AliModule*>(next()))) {
214 detector->AddAlignableVolumes();
218 //_______________________________________________________________________
219 void AliMC::GeneratePrimaries()
222 // Generate primary particles and fill them in the stack.
225 Generator()->Generate();
228 //_______________________________________________________________________
229 void AliMC::SetGenerator(AliGenerator *generator)
232 // Load the event generator
234 if(!fGenerator) fGenerator = generator;
237 //_______________________________________________________________________
238 void AliMC::ResetGenerator(AliGenerator *generator)
241 // Load the event generator
245 AliWarning(Form("Replacing generator %s with %s",
246 fGenerator->GetName(),generator->GetName()))
248 AliWarning(Form("Replacing generator %s with NULL",
249 fGenerator->GetName()));
250 fGenerator = generator;
253 //_______________________________________________________________________
254 void AliMC::FinishRun()
256 // Clean generator information
257 AliDebug(1, "fGenerator->FinishRun()");
258 fGenerator->FinishRun();
260 //Output energy summary tables
261 AliDebug(1, "EnergySummary()");
262 ToAliDebug(1, EnergySummary());
265 //_______________________________________________________________________
266 void AliMC::BeginPrimary()
269 // Called at the beginning of each primary track
274 ResetTrackReferences();
278 //_______________________________________________________________________
279 void AliMC::PreTrack()
281 // Actions before the track's transport
282 TObjArray &dets = *gAlice->Modules();
285 for(Int_t i=0; i<=gAlice->GetNdets(); i++)
286 if((module = dynamic_cast<AliModule*>(dets[i])))
292 //_______________________________________________________________________
293 void AliMC::Stepping()
296 // Called at every step during transport
299 Int_t id = DetFromMate(gMC->CurrentMedium());
303 if ( gMC->IsNewTrack() &&
304 gMC->TrackTime() == 0. &&
306 fRDecayMax > fRDecayMin &&
307 gMC->TrackPid() == fDecayPdg )
309 FixParticleDecaytime();
315 // --- If lego option, do it and leave
317 gAlice->Lego()->StepManager();
320 //Update energy deposition tables
321 AddEnergyDeposit(gMC->CurrentVolID(copy),gMC->Edep());
323 // write tracke reference for track which is dissapearing - MI
324 if (gMC->IsTrackDisappeared()) {
325 if (gMC->Etot()>0.05) AddTrackReference(GetCurrentTrackNumber());
328 //Call the appropriate stepping routine;
329 AliModule *det = dynamic_cast<AliModule*>(gAlice->Modules()->At(id));
330 if(det && det->StepManagerIsEnabled()) {
331 if(AliLog::GetGlobalDebugLevel()>0) fMCQA->StepManager(id);
337 //_______________________________________________________________________
338 void AliMC::EnergySummary()
341 // Print summary of deposited energy
347 Int_t kn, i, left, j, id;
348 const Float_t kzero=0;
349 Int_t ievent=gAlice->GetRunLoader()->GetHeader()->GetEvent()+1;
351 // Energy loss information
353 printf("***************** Energy Loss Information per event (GEV) *****************\n");
354 for(kn=1;kn<fEventEnergy.GetSize();kn++) {
357 fEventEnergy[ndep]=kn;
362 ed2=100*TMath::Sqrt(TMath::Max(ed2-ed*ed,kzero))/ed;
365 fSummEnergy[ndep]=ed;
366 fSum2Energy[ndep]=TMath::Min(static_cast<Float_t>(99.),TMath::Max(ed2,kzero));
371 for(kn=0;kn<(ndep-1)/3+1;kn++) {
373 for(i=0;i<(3<left?3:left);i++) {
375 id=Int_t (fEventEnergy[j]+0.1);
376 printf(" %s %10.3f +- %10.3f%%;",gMC->VolName(id),fSummEnergy[j],fSum2Energy[j]);
381 // Relative energy loss in different detectors
382 printf("******************** Relative Energy Loss per event ********************\n");
383 printf("Total energy loss per event %10.3f GeV\n",edtot);
384 for(kn=0;kn<(ndep-1)/5+1;kn++) {
386 for(i=0;i<(5<left?5:left);i++) {
388 id=Int_t (fEventEnergy[j]+0.1);
389 printf(" %s %10.3f%%;",gMC->VolName(id),100*fSummEnergy[j]/edtot);
393 for(kn=0;kn<75;kn++) printf("*");
397 // Reset the TArray's
398 // fEventEnergy.Set(0);
399 // fSummEnergy.Set(0);
400 // fSum2Energy.Set(0);
403 //_____________________________________________________________________________
404 void AliMC::BeginEvent()
407 // Clean-up previous event
409 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
410 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
411 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
412 AliDebug(1, " BEGINNING EVENT ");
413 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
414 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
415 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
417 AliRunLoader *runloader=gAlice->GetRunLoader();
419 /*******************************/
420 /* Clean after eventual */
422 /*******************************/
425 //Set the next event in Run Loader -> Cleans trees (TreeK and all trees in detectors),
426 gAlice->SetEventNrInRun(gAlice->GetEventNrInRun()+1);
427 runloader->SetEventNumber(gAlice->GetEventNrInRun());// sets new files, cleans the previous event stuff, if necessary, etc.,
428 AliDebug(1, Form("EventNr is %d",gAlice->GetEventNrInRun()));
430 fEventEnergy.Reset();
431 // Clean detector information
433 if (runloader->Stack())
434 runloader->Stack()->Reset();//clean stack -> tree is unloaded
436 runloader->MakeStack();//or make a new one
438 if(gAlice->Lego() == 0x0)
440 AliDebug(1, "fRunLoader->MakeTree(K)");
441 runloader->MakeTree("K");
444 AliDebug(1, "gMC->SetStack(fRunLoader->Stack())");
445 gMC->SetStack(gAlice->GetRunLoader()->Stack());//Was in InitMC - but was moved here
446 //because we don't have guarantee that
447 //stack pointer is not going to change from event to event
448 //since it bellobgs to header and is obtained via RunLoader
450 // Reset all Detectors & kinematics & make/reset trees
453 runloader->GetHeader()->Reset(gAlice->GetRunNumber(),gAlice->GetEvNumber(),
454 gAlice->GetEventNrInRun());
455 // fRunLoader->WriteKinematics("OVERWRITE"); is there any reason to rewrite here since MakeTree does so
459 gAlice->Lego()->BeginEvent();
464 AliDebug(1, "ResetHits()");
467 AliDebug(1, "fRunLoader->MakeTree(H)");
468 runloader->MakeTree("H");
470 AliDebug(1, "fRunLoader->MakeTrackRefsContainer()");
471 runloader->MakeTrackRefsContainer();//for insurance
474 //create new branches and SetAdresses
475 TIter next(gAlice->Modules());
477 while((detector = (AliModule*)next()))
479 AliDebug(2, Form("%s->MakeBranch(H)",detector->GetName()));
480 detector->MakeBranch("H");
481 AliDebug(2, Form("%s->MakeBranchTR()",detector->GetName()));
482 detector->MakeBranchTR();
483 AliDebug(2, Form("%s->SetTreeAddress()",detector->GetName()));
484 detector->SetTreeAddress();
486 // make branch for AliRun track References
487 TTree * treeTR = runloader->TreeTR();
489 // make branch for central track references
490 if (!fTrackReferences) fTrackReferences = new TClonesArray("AliTrackReference",0);
492 branch = treeTR->Branch("AliRun",&fTrackReferences);
493 branch->SetAddress(&fTrackReferences);
498 //_______________________________________________________________________
499 void AliMC::ResetHits()
502 // Reset all Detectors hits
504 TIter next(gAlice->Modules());
506 while((detector = dynamic_cast<AliModule*>(next()))) {
507 detector->ResetHits();
511 //_______________________________________________________________________
512 void AliMC::PostTrack()
514 // Posts tracks for each module
515 TObjArray &dets = *gAlice->Modules();
518 for(Int_t i=0; i<=gAlice->GetNdets(); i++)
519 if((module = dynamic_cast<AliModule*>(dets[i])))
523 //_______________________________________________________________________
524 void AliMC::FinishPrimary()
527 // Called at the end of each primary track
529 AliRunLoader *runloader=gAlice->GetRunLoader();
530 // static Int_t count=0;
531 // const Int_t times=10;
532 // This primary is finished, purify stack
533 #if ROOT_VERSION_CODE > 262152
534 if (!(gMC->SecondariesAreOrdered()))
535 runloader->Stack()->ReorderKine();
537 runloader->Stack()->PurifyKine();
539 TIter next(gAlice->Modules());
541 while((detector = dynamic_cast<AliModule*>(next()))) {
542 detector->FinishPrimary();
543 AliLoader* loader = detector->GetLoader();
546 TTree* treeH = loader->TreeH();
547 if (treeH) treeH->Fill(); //can be Lego run and treeH can not exist
551 // Write out track references if any
552 if (runloader->TreeTR()) runloader->TreeTR()->Fill();
555 //_______________________________________________________________________
556 void AliMC::FinishEvent()
559 // Called at the end of the event.
564 if(gAlice->Lego()) gAlice->Lego()->FinishEvent();
566 TIter next(gAlice->Modules());
568 while((detector = dynamic_cast<AliModule*>(next()))) {
569 detector->FinishEvent();
572 //Update the energy deposit tables
574 for(i=0;i<fEventEnergy.GetSize();i++)
576 fSummEnergy[i]+=fEventEnergy[i];
577 fSum2Energy[i]+=fEventEnergy[i]*fEventEnergy[i];
580 AliRunLoader *runloader=gAlice->GetRunLoader();
582 AliHeader* header = runloader->GetHeader();
583 AliStack* stack = runloader->Stack();
584 if ( (header == 0x0) || (stack == 0x0) )
585 {//check if we got header and stack. If not cry and exit aliroot
586 AliFatal("Can not get the stack or header from LOADER");
587 return;//never reached
589 // Update Header information
590 header->SetNprimary(stack->GetNprimary());
591 header->SetNtrack(stack->GetNtrack());
594 // Write out the kinematics
595 if (!gAlice->Lego()) stack->FinishEvent();
597 // Write out the event Header information
598 TTree* treeE = runloader->TreeE();
601 header->SetStack(stack);
606 AliError("Can not get TreeE from RL");
609 if(gAlice->Lego() == 0x0)
611 runloader->WriteKinematics("OVERWRITE");
612 runloader->WriteTrackRefs("OVERWRITE");
613 runloader->WriteHits("OVERWRITE");
616 AliDebug(1, "<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<");
617 AliDebug(1, "<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<");
618 AliDebug(1, "<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<");
619 AliDebug(1, " FINISHING EVENT ");
620 AliDebug(1, "<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<");
621 AliDebug(1, "<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<");
622 AliDebug(1, "<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<");
625 //_______________________________________________________________________
626 void AliMC::Field(const Double_t* x, Double_t* b) const
628 // Calculates field "b" at point "x"
632 //_______________________________________________________________________
637 //=================Create Materials and geometry
639 //Set alignable volumes for the whole geometry
640 SetAllAlignableVolumes();
641 //Read the cuts for all materials
643 //Build the special IMEDIA table
646 //Compute cross-sections
649 //Initialise geometry deposition table
650 fEventEnergy.Set(gMC->NofVolumes()+1);
651 fSummEnergy.Set(gMC->NofVolumes()+1);
652 fSum2Energy.Set(gMC->NofVolumes()+1);
655 fMCQA = new AliMCQA(gAlice->GetNdets());
657 // Register MC in configuration
658 AliConfig::Instance()->Add(gMC);
662 //_______________________________________________________________________
663 void AliMC::MediaTable()
666 // Built media table to get from the media number to
670 Int_t kz, nz, idt, lz, i, k, ind;
672 TObjArray &dets = *gAlice->Detectors();
674 Int_t ndets=gAlice->GetNdets();
677 for (kz=0;kz<ndets;kz++) {
678 // If detector is defined
679 if((det=dynamic_cast<AliModule*>(dets[kz]))) {
680 TArrayI &idtmed = *(det->GetIdtmed());
681 for(nz=0;nz<100;nz++) {
683 // Find max and min material number
684 if((idt=idtmed[nz])) {
685 det->LoMedium() = det->LoMedium() < idt ? det->LoMedium() : idt;
686 det->HiMedium() = det->HiMedium() > idt ? det->HiMedium() : idt;
689 if(det->LoMedium() > det->HiMedium()) {
693 if(det->HiMedium() > fImedia->GetSize()) {
694 AliError(Form("Increase fImedia from %d to %d",
695 fImedia->GetSize(),det->HiMedium()));
698 // Tag all materials in rage as belonging to detector kz
699 for(lz=det->LoMedium(); lz<= det->HiMedium(); lz++) {
706 // Print summary table
707 AliInfo("Tracking media ranges:");
709 for(i=0;i<(ndets-1)/6+1;i++) {
710 for(k=0;k< (6<ndets-i*6?6:ndets-i*6);k++) {
712 det=dynamic_cast<AliModule*>(dets[ind]);
714 printf(" %6s: %3d -> %3d;",det->GetName(),det->LoMedium(),
717 printf(" %6s: %3d -> %3d;","NULL",0,0);
724 //_______________________________________________________________________
725 void AliMC::ReadTransPar()
728 // Read filename to set the transport parameters
732 const Int_t kncuts=10;
733 const Int_t knflags=11;
734 const Int_t knpars=kncuts+knflags;
735 const char kpars[knpars][7] = {"CUTGAM" ,"CUTELE","CUTNEU","CUTHAD","CUTMUO",
736 "BCUTE","BCUTM","DCUTE","DCUTM","PPCUTM","ANNI",
737 "BREM","COMP","DCAY","DRAY","HADR","LOSS",
738 "MULS","PAIR","PHOT","RAYL"};
744 Int_t i, itmed, iret, ktmed, kz;
747 // See whether the file is there
748 filtmp=gSystem->ExpandPathName(fTransParName.Data());
749 lun=fopen(filtmp,"r");
752 AliWarning(Form("File %s does not exist!",fTransParName.Data()));
757 // Initialise cuts and flags
758 for(i=0;i<kncuts;i++) cut[i]=-99;
759 for(i=0;i<knflags;i++) flag[i]=-99;
761 for(i=0;i<256;i++) line[i]='\0';
762 // Read up to the end of line excluded
763 iret=fscanf(lun,"%[^\n]",line);
769 // Read the end of line
772 if(line[0]=='*') continue;
774 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",
775 detName,&itmed,&cut[0],&cut[1],&cut[2],&cut[3],&cut[4],&cut[5],&cut[6],&cut[7],&cut[8],
776 &cut[9],&flag[0],&flag[1],&flag[2],&flag[3],&flag[4],&flag[5],&flag[6],&flag[7],
777 &flag[8],&flag[9],&flag[10]);
781 AliWarning(Form("Error reading file %s",fTransParName.Data()));
784 // Check that the module exist
785 AliModule *mod = gAlice->GetModule(detName);
787 // Get the array of media numbers
788 TArrayI &idtmed = *mod->GetIdtmed();
789 // Check that the tracking medium code is valid
790 if(0<=itmed && itmed < 100) {
793 AliWarning(Form("Invalid tracking medium code %d for %s",itmed,mod->GetName()));
796 // Set energy thresholds
797 for(kz=0;kz<kncuts;kz++) {
799 AliDebug(2, Form("%-6s set to %10.3E for tracking medium code %4d for %s",
800 kpars[kz],cut[kz],itmed,mod->GetName()));
801 gMC->Gstpar(ktmed,kpars[kz],cut[kz]);
804 // Set transport mechanisms
805 for(kz=0;kz<knflags;kz++) {
807 AliDebug(2, Form("%-6s set to %10d for tracking medium code %4d for %s",
808 kpars[kncuts+kz],flag[kz],itmed,mod->GetName()));
809 gMC->Gstpar(ktmed,kpars[kncuts+kz],Float_t(flag[kz]));
813 AliWarning(Form("Invalid medium code %d",itmed));
817 AliDebug(1, Form("%s not present",detName));
823 //_______________________________________________________________________
824 void AliMC::SetTransPar(const char *filename)
827 // Sets the file name for transport parameters
829 fTransParName = filename;
832 //_______________________________________________________________________
833 void AliMC::Browse(TBrowser *b)
836 // Called when the item "Run" is clicked on the left pane
837 // of the Root browser.
838 // It displays the Root Trees and all detectors.
840 //detectors are in folders anyway
841 b->Add(fMCQA,"AliMCQA");
844 //_______________________________________________________________________
845 void AliMC::AddHit(Int_t id, Int_t track, Int_t *vol, Float_t *hits) const
848 // Add a hit to detector id
850 TObjArray &dets = *gAlice->Modules();
851 if(dets[id]) dynamic_cast<AliModule*>(dets[id])->AddHit(track,vol,hits);
854 //_______________________________________________________________________
855 void AliMC::AddDigit(Int_t id, Int_t *tracks, Int_t *digits) const
858 // Add digit to detector id
860 TObjArray &dets = *gAlice->Modules();
861 if(dets[id]) dynamic_cast<AliModule*>(dets[id])->AddDigit(tracks,digits);
864 //_______________________________________________________________________
865 Int_t AliMC::GetCurrentTrackNumber() const {
867 // Returns current track
869 return gAlice->GetRunLoader()->Stack()->GetCurrentTrackNumber();
872 //_______________________________________________________________________
873 void AliMC::DumpPart (Int_t i) const
876 // Dumps particle i in the stack
878 AliRunLoader * runloader = gAlice->GetRunLoader();
879 if (runloader->Stack())
880 runloader->Stack()->DumpPart(i);
883 //_______________________________________________________________________
884 void AliMC::DumpPStack () const
887 // Dumps the particle stack
889 AliRunLoader * runloader = gAlice->GetRunLoader();
890 if (runloader->Stack())
891 runloader->Stack()->DumpPStack();
894 //_______________________________________________________________________
895 Int_t AliMC::GetNtrack() const {
897 // Returns number of tracks in stack
900 AliRunLoader * runloader = gAlice->GetRunLoader();
901 if (runloader->Stack())
902 ntracks = runloader->Stack()->GetNtrack();
906 //_______________________________________________________________________
907 Int_t AliMC::GetPrimary(Int_t track) const
910 // return number of primary that has generated track
912 Int_t nprimary = -999;
913 AliRunLoader * runloader = gAlice->GetRunLoader();
914 if (runloader->Stack())
915 nprimary = runloader->Stack()->GetPrimary(track);
919 //_______________________________________________________________________
920 TParticle* AliMC::Particle(Int_t i) const
922 // Returns the i-th particle from the stack taking into account
923 // the remaping done by PurifyKine
924 AliRunLoader * runloader = gAlice->GetRunLoader();
926 if (runloader->Stack())
927 return runloader->Stack()->Particle(i);
931 //_______________________________________________________________________
932 TObjArray* AliMC::Particles() const {
934 // Returns pointer to Particles array
936 AliRunLoader * runloader = gAlice->GetRunLoader();
938 if (runloader->Stack())
939 return runloader->Stack()->Particles();
943 //_______________________________________________________________________
944 void AliMC::PushTrack(Int_t done, Int_t parent, Int_t pdg, Float_t *pmom,
945 Float_t *vpos, Float_t *polar, Float_t tof,
946 TMCProcess mech, Int_t &ntr, Float_t weight, Int_t is) const
950 AliRunLoader * runloader = gAlice->GetRunLoader();
952 if (runloader->Stack())
953 runloader->Stack()->PushTrack(done, parent, pdg, pmom, vpos, polar, tof,
954 mech, ntr, weight, is);
957 //_______________________________________________________________________
958 void AliMC::PushTrack(Int_t done, Int_t parent, Int_t pdg,
959 Double_t px, Double_t py, Double_t pz, Double_t e,
960 Double_t vx, Double_t vy, Double_t vz, Double_t tof,
961 Double_t polx, Double_t poly, Double_t polz,
962 TMCProcess mech, Int_t &ntr, Float_t weight, Int_t is) const
966 AliRunLoader * runloader = gAlice->GetRunLoader();
968 if (runloader->Stack())
969 runloader->Stack()->PushTrack(done, parent, pdg, px, py, pz, e, vx, vy, vz, tof,
970 polx, poly, polz, mech, ntr, weight, is);
973 //_______________________________________________________________________
974 void AliMC::SetHighWaterMark(Int_t nt) const
977 // Set high water mark for last track in event
978 AliRunLoader * runloader = gAlice->GetRunLoader();
980 if (runloader->Stack())
981 runloader->Stack()->SetHighWaterMark(nt);
984 //_______________________________________________________________________
985 void AliMC::KeepTrack(Int_t track) const
990 AliRunLoader * runloader = gAlice->GetRunLoader();
992 if (runloader->Stack())
993 runloader->Stack()->KeepTrack(track);
996 //_______________________________________________________________________
997 void AliMC::FlagTrack(Int_t track) const
1001 AliRunLoader * runloader = gAlice->GetRunLoader();
1003 if (runloader->Stack())
1004 runloader->Stack()->FlagTrack(track);
1007 //_______________________________________________________________________
1008 void AliMC::SetCurrentTrack(Int_t track) const
1011 // Set current track number
1013 AliRunLoader * runloader = gAlice->GetRunLoader();
1015 if (runloader->Stack())
1016 runloader->Stack()->SetCurrentTrack(track);
1019 //_______________________________________________________________________
1020 void AliMC::AddTrackReference(Int_t label){
1022 // add a trackrefernce to the list
1023 if (!fTrackReferences) {
1024 AliError("Container trackrefernce not active");
1027 Int_t nref = fTrackReferences->GetEntriesFast();
1028 TClonesArray &lref = *fTrackReferences;
1029 new(lref[nref]) AliTrackReference(label);
1034 //_______________________________________________________________________
1035 void AliMC::ResetTrackReferences()
1038 // Reset all references
1040 if (fTrackReferences) fTrackReferences->Clear();
1042 TIter next(gAlice->Modules());
1043 AliModule *detector;
1044 while((detector = dynamic_cast<AliModule*>(next()))) {
1045 detector->ResetTrackReferences();
1049 void AliMC::RemapTrackReferencesIDs(Int_t *map)
1052 // Remapping track reference
1053 // Called at finish primary
1055 if (!fTrackReferences) return;
1056 Int_t nEntries = fTrackReferences->GetEntries();
1057 for (Int_t i=0; i < nEntries; i++){
1058 AliTrackReference * ref = dynamic_cast<AliTrackReference*>(fTrackReferences->UncheckedAt(i));
1060 Int_t newID = map[ref->GetTrack()];
1061 if (newID>=0) ref->SetTrack(newID);
1063 ref->SetBit(kNotDeleted,kFALSE);
1064 fTrackReferences->RemoveAt(i);
1068 fTrackReferences->Compress();
1071 void AliMC::FixParticleDecaytime()
1074 // Fix the particle decay time according to rmin and rmax for decays
1078 gMC->TrackMomentum(p);
1079 Double_t tmin, tmax;
1082 // Transverse velocity
1083 Double_t vt = p.Pt() / p.E();
1085 if ((b = gAlice->Field()->SolenoidField()) > 0.) { // [kG]
1089 Double_t rho = p.Pt() / 0.0003 / b; // [cm]
1091 // Revolution frequency
1093 Double_t omega = vt / rho;
1095 // Maximum and minimum decay time
1097 // Check for curlers first
1098 if (fRDecayMax * fRDecayMax / rho / rho / 2. > 1.) return;
1102 tmax = TMath::ACos(1. - fRDecayMax * fRDecayMax / rho / rho / 2.) / omega; // [ct]
1103 tmin = TMath::ACos(1. - fRDecayMin * fRDecayMin / rho / rho / 2.) / omega; // [ct]
1105 tmax = fRDecayMax / vt; // [ct]
1106 tmin = fRDecayMin / vt; // [ct]
1110 // Dial t using the two limits
1111 Double_t t = tmin + (tmax - tmin) * gRandom->Rndm(); // [ct]
1114 // Force decay time in transport code
1116 gMC->ForceDecayTime(t / 2.99792458e10);