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>
33 #include "AliDetector.h"
34 #include "AliGenerator.h"
35 #include "AliHeader.h"
42 #include "AliTrackReference.h"
47 //_______________________________________________________________________
66 //_______________________________________________________________________
67 AliMC::AliMC(const char *name, const char *title) :
68 TVirtualMCApplication(name, title),
75 fImedia(new TArrayI(1000)),
78 fHitLists(new TList()),
79 fTrackReferences(new TClonesArray("AliTrackReference", 100))
82 // Set transport parameters
85 // Prepare the tracking medium lists
86 for(Int_t i=0;i<1000;i++) (*fImedia)[i]=-99;
89 //_______________________________________________________________________
90 AliMC::AliMC(const AliMC &mc) :
91 TVirtualMCApplication(mc),
105 // Copy constructor for AliRun
110 //_______________________________________________________________________
118 // Delete track references
119 if (fTrackReferences) {
120 fTrackReferences->Delete();
121 delete fTrackReferences;
122 fTrackReferences = 0;
127 //_______________________________________________________________________
128 void AliMC::Copy(TObject &) const
130 //dummy Copy function
131 AliFatal("Not implemented!");
134 //_______________________________________________________________________
135 void AliMC::ConstructGeometry()
138 // Create modules, materials, geometry
142 TIter next(gAlice->Modules());
144 AliDebug(1, "Geometry creation:");
145 while((detector = dynamic_cast<AliModule*>(next()))) {
147 // Initialise detector materials and geometry
148 detector->CreateMaterials();
149 detector->CreateGeometry();
150 AliInfo(Form("%10s R:%.2fs C:%.2fs",
151 detector->GetName(),stw.RealTime(),stw.CpuTime()));
155 //_______________________________________________________________________
156 void AliMC::InitGeometry()
159 // Initialize detectors and display geometry
162 AliInfo("Initialisation:");
164 TIter next(gAlice->Modules());
166 while((detector = dynamic_cast<AliModule*>(next()))) {
168 // Initialise detector and display geometry
170 detector->BuildGeometry();
171 AliInfo(Form("%10s R:%.2fs C:%.2fs",
172 detector->GetName(),stw.RealTime(),stw.CpuTime()));
177 //_______________________________________________________________________
178 void AliMC::GeneratePrimaries()
181 // Generate primary particles and fill them in the stack.
184 Generator()->Generate();
187 //_______________________________________________________________________
188 void AliMC::SetGenerator(AliGenerator *generator)
191 // Load the event generator
193 if(!fGenerator) fGenerator = generator;
196 //_______________________________________________________________________
197 void AliMC::ResetGenerator(AliGenerator *generator)
200 // Load the event generator
204 AliWarning(Form("Replacing generator %s with %s",
205 fGenerator->GetName(),generator->GetName()))
207 AliWarning(Form("Replacing generator %s with NULL",
208 fGenerator->GetName()));
209 fGenerator = generator;
212 //_______________________________________________________________________
213 void AliMC::FinishRun()
215 // Clean generator information
216 AliDebug(1, "fGenerator->FinishRun()");
217 fGenerator->FinishRun();
219 //Output energy summary tables
220 AliDebug(1, "EnergySummary()");
221 ToAliDebug(1, EnergySummary());
224 //_______________________________________________________________________
225 void AliMC::BeginPrimary()
228 // Called at the beginning of each primary track
233 ResetTrackReferences();
237 //_______________________________________________________________________
238 void AliMC::PreTrack()
240 // Actions before the track's transport
241 TObjArray &dets = *gAlice->Modules();
244 for(Int_t i=0; i<=gAlice->GetNdets(); i++)
245 if((module = dynamic_cast<AliModule*>(dets[i])))
251 //_______________________________________________________________________
252 void AliMC::Stepping()
255 // Called at every step during transport
258 Int_t id = DetFromMate(gMC->GetMedium());
262 if ( gMC->IsNewTrack() &&
263 gMC->TrackTime() == 0. &&
265 fRDecayMax > fRDecayMin &&
266 gMC->TrackPid() == fDecayPdg )
268 FixParticleDecaytime();
274 // --- If lego option, do it and leave
276 gAlice->Lego()->StepManager();
279 //Update energy deposition tables
280 AddEnergyDeposit(gMC->CurrentVolID(copy),gMC->Edep());
282 // write tracke reference for track which is dissapearing - MI
283 if (gMC->IsTrackDisappeared()) {
284 if (gMC->Etot()>0.05) AddTrackReference(GetCurrentTrackNumber());
287 //Call the appropriate stepping routine;
288 AliModule *det = dynamic_cast<AliModule*>(gAlice->Modules()->At(id));
289 if(det && det->StepManagerIsEnabled()) {
290 fMCQA->StepManager(id);
296 //_______________________________________________________________________
297 void AliMC::EnergySummary()
300 // Print summary of deposited energy
306 Int_t kn, i, left, j, id;
307 const Float_t kzero=0;
308 Int_t ievent=gAlice->GetRunLoader()->GetHeader()->GetEvent()+1;
310 // Energy loss information
312 printf("***************** Energy Loss Information per event (GEV) *****************\n");
313 for(kn=1;kn<fEventEnergy.GetSize();kn++) {
316 fEventEnergy[ndep]=kn;
321 ed2=100*TMath::Sqrt(TMath::Max(ed2-ed*ed,kzero))/ed;
324 fSummEnergy[ndep]=ed;
325 fSum2Energy[ndep]=TMath::Min(static_cast<Float_t>(99.),TMath::Max(ed2,kzero));
330 for(kn=0;kn<(ndep-1)/3+1;kn++) {
332 for(i=0;i<(3<left?3:left);i++) {
334 id=Int_t (fEventEnergy[j]+0.1);
335 printf(" %s %10.3f +- %10.3f%%;",gMC->VolName(id),fSummEnergy[j],fSum2Energy[j]);
340 // Relative energy loss in different detectors
341 printf("******************** Relative Energy Loss per event ********************\n");
342 printf("Total energy loss per event %10.3f GeV\n",edtot);
343 for(kn=0;kn<(ndep-1)/5+1;kn++) {
345 for(i=0;i<(5<left?5:left);i++) {
347 id=Int_t (fEventEnergy[j]+0.1);
348 printf(" %s %10.3f%%;",gMC->VolName(id),100*fSummEnergy[j]/edtot);
352 for(kn=0;kn<75;kn++) printf("*");
356 // Reset the TArray's
357 // fEventEnergy.Set(0);
358 // fSummEnergy.Set(0);
359 // fSum2Energy.Set(0);
362 //_____________________________________________________________________________
363 void AliMC::BeginEvent()
366 // Clean-up previous event
368 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
369 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
370 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
371 AliDebug(1, " BEGINNING EVENT ");
372 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
373 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
374 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
376 AliRunLoader *runloader=gAlice->GetRunLoader();
378 /*******************************/
379 /* Clean after eventual */
381 /*******************************/
384 //Set the next event in Run Loader -> Cleans trees (TreeK and all trees in detectors),
385 gAlice->SetEventNrInRun(gAlice->GetEventNrInRun()+1);
386 runloader->SetEventNumber(gAlice->GetEventNrInRun());// sets new files, cleans the previous event stuff, if necessary, etc.,
387 AliDebug(1, Form("EventNr is %d",gAlice->GetEventNrInRun()));
389 fEventEnergy.Reset();
390 // Clean detector information
392 if (runloader->Stack())
393 runloader->Stack()->Reset();//clean stack -> tree is unloaded
395 runloader->MakeStack();//or make a new one
397 if(gAlice->Lego() == 0x0)
399 AliDebug(1, "fRunLoader->MakeTree(K)");
400 runloader->MakeTree("K");
403 AliDebug(1, "gMC->SetStack(fRunLoader->Stack())");
404 gMC->SetStack(gAlice->GetRunLoader()->Stack());//Was in InitMC - but was moved here
405 //because we don't have guarantee that
406 //stack pointer is not going to change from event to event
407 //since it bellobgs to header and is obtained via RunLoader
409 // Reset all Detectors & kinematics & make/reset trees
412 runloader->GetHeader()->Reset(gAlice->GetRunNumber(),gAlice->GetEvNumber(),
413 gAlice->GetEventNrInRun());
414 // fRunLoader->WriteKinematics("OVERWRITE"); is there any reason to rewrite here since MakeTree does so
418 gAlice->Lego()->BeginEvent();
423 AliDebug(1, "ResetHits()");
426 AliDebug(1, "fRunLoader->MakeTree(H)");
427 runloader->MakeTree("H");
429 AliDebug(1, "fRunLoader->MakeTrackRefsContainer()");
430 runloader->MakeTrackRefsContainer();//for insurance
433 //create new branches and SetAdresses
434 TIter next(gAlice->Modules());
436 while((detector = (AliModule*)next()))
438 AliDebug(2, Form("%s->MakeBranch(H)",detector->GetName()));
439 detector->MakeBranch("H");
440 AliDebug(2, Form("%s->MakeBranchTR()",detector->GetName()));
441 detector->MakeBranchTR();
442 AliDebug(2, Form("%s->SetTreeAddress()",detector->GetName()));
443 detector->SetTreeAddress();
445 // make branch for AliRun track References
446 TTree * treeTR = runloader->TreeTR();
448 // make branch for central track references
449 if (!fTrackReferences) fTrackReferences = new TClonesArray("AliTrackReference",0);
451 branch = treeTR->Branch("AliRun",&fTrackReferences);
452 branch->SetAddress(&fTrackReferences);
457 //_______________________________________________________________________
458 void AliMC::ResetHits()
461 // Reset all Detectors hits
463 TIter next(gAlice->Modules());
465 while((detector = dynamic_cast<AliModule*>(next()))) {
466 detector->ResetHits();
470 //_______________________________________________________________________
471 void AliMC::PostTrack()
473 // Posts tracks for each module
474 TObjArray &dets = *gAlice->Modules();
477 for(Int_t i=0; i<=gAlice->GetNdets(); i++)
478 if((module = dynamic_cast<AliModule*>(dets[i])))
482 //_______________________________________________________________________
483 void AliMC::FinishPrimary()
486 // Called at the end of each primary track
488 AliRunLoader *runloader=gAlice->GetRunLoader();
489 // static Int_t count=0;
490 // const Int_t times=10;
491 // This primary is finished, purify stack
492 #if ROOT_VERSION_CODE > 262152
493 if (!(gMC->SecondariesAreOrdered()))
494 runloader->Stack()->ReorderKine();
496 runloader->Stack()->PurifyKine();
498 TIter next(gAlice->Modules());
500 while((detector = dynamic_cast<AliModule*>(next()))) {
501 detector->FinishPrimary();
502 AliLoader* loader = detector->GetLoader();
505 TTree* treeH = loader->TreeH();
506 if (treeH) treeH->Fill(); //can be Lego run and treeH can not exist
510 // Write out track references if any
511 if (runloader->TreeTR()) runloader->TreeTR()->Fill();
514 //_______________________________________________________________________
515 void AliMC::FinishEvent()
518 // Called at the end of the event.
522 if(gAlice->Lego()) gAlice->Lego()->FinishEvent();
524 TIter next(gAlice->Modules());
526 while((detector = dynamic_cast<AliModule*>(next()))) {
527 detector->FinishEvent();
530 //Update the energy deposit tables
532 for(i=0;i<fEventEnergy.GetSize();i++)
534 fSummEnergy[i]+=fEventEnergy[i];
535 fSum2Energy[i]+=fEventEnergy[i]*fEventEnergy[i];
538 AliRunLoader *runloader=gAlice->GetRunLoader();
540 AliHeader* header = runloader->GetHeader();
541 AliStack* stack = runloader->Stack();
542 if ( (header == 0x0) || (stack == 0x0) )
543 {//check if we got header and stack. If not cry and exit aliroot
544 AliFatal("Can not get the stack or header from LOADER");
545 return;//never reached
547 // Update Header information
548 header->SetNprimary(stack->GetNprimary());
549 header->SetNtrack(stack->GetNtrack());
552 // Write out the kinematics
553 stack->FinishEvent();
555 // Write out the event Header information
556 TTree* treeE = runloader->TreeE();
559 header->SetStack(stack);
564 AliError("Can not get TreeE from RL");
567 if(gAlice->Lego() == 0x0)
569 runloader->WriteKinematics("OVERWRITE");
570 runloader->WriteTrackRefs("OVERWRITE");
571 runloader->WriteHits("OVERWRITE");
574 AliDebug(1, "<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<");
575 AliDebug(1, "<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<");
576 AliDebug(1, "<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<");
577 AliDebug(1, " FINISHING EVENT ");
578 AliDebug(1, "<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<");
579 AliDebug(1, "<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<");
580 AliDebug(1, "<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<");
583 //_______________________________________________________________________
584 void AliMC::Field(const Double_t* x, Double_t* b) const
586 // Calculates field "b" at point "x"
590 //_______________________________________________________________________
595 //=================Create Materials and geometry
597 //Read the cuts for all materials
599 //Build the special IMEDIA table
602 //Compute cross-sections
605 //Write Geometry object to current file.
606 gAlice->GetRunLoader()->WriteGeometry();
608 //Initialise geometry deposition table
609 fEventEnergy.Set(gMC->NofVolumes()+1);
610 fSummEnergy.Set(gMC->NofVolumes()+1);
611 fSum2Energy.Set(gMC->NofVolumes()+1);
614 fMCQA = new AliMCQA(gAlice->GetNdets());
616 // Register MC in configuration
617 AliConfig::Instance()->Add(gMC);
621 //_______________________________________________________________________
622 void AliMC::MediaTable()
625 // Built media table to get from the media number to
629 Int_t kz, nz, idt, lz, i, k, ind;
631 TObjArray &dets = *gAlice->Detectors();
633 Int_t ndets=gAlice->GetNdets();
636 for (kz=0;kz<ndets;kz++) {
637 // If detector is defined
638 if((det=dynamic_cast<AliModule*>(dets[kz]))) {
639 TArrayI &idtmed = *(det->GetIdtmed());
640 for(nz=0;nz<100;nz++) {
642 // Find max and min material number
643 if((idt=idtmed[nz])) {
644 det->LoMedium() = det->LoMedium() < idt ? det->LoMedium() : idt;
645 det->HiMedium() = det->HiMedium() > idt ? det->HiMedium() : idt;
648 if(det->LoMedium() > det->HiMedium()) {
652 if(det->HiMedium() > fImedia->GetSize()) {
653 AliError(Form("Increase fImedia from %d to %d",
654 fImedia->GetSize(),det->HiMedium()));
657 // Tag all materials in rage as belonging to detector kz
658 for(lz=det->LoMedium(); lz<= det->HiMedium(); lz++) {
665 // Print summary table
666 AliInfo("Tracking media ranges:");
668 for(i=0;i<(ndets-1)/6+1;i++) {
669 for(k=0;k< (6<ndets-i*6?6:ndets-i*6);k++) {
671 det=dynamic_cast<AliModule*>(dets[ind]);
673 printf(" %6s: %3d -> %3d;",det->GetName(),det->LoMedium(),
676 printf(" %6s: %3d -> %3d;","NULL",0,0);
683 //_______________________________________________________________________
684 void AliMC::ReadTransPar()
687 // Read filename to set the transport parameters
691 const Int_t kncuts=10;
692 const Int_t knflags=11;
693 const Int_t knpars=kncuts+knflags;
694 const char kpars[knpars][7] = {"CUTGAM" ,"CUTELE","CUTNEU","CUTHAD","CUTMUO",
695 "BCUTE","BCUTM","DCUTE","DCUTM","PPCUTM","ANNI",
696 "BREM","COMP","DCAY","DRAY","HADR","LOSS",
697 "MULS","PAIR","PHOT","RAYL"};
703 Int_t i, itmed, iret, ktmed, kz;
706 // See whether the file is there
707 filtmp=gSystem->ExpandPathName(fTransParName.Data());
708 lun=fopen(filtmp,"r");
711 AliWarning(Form("File %s does not exist!",fTransParName.Data()));
716 // Initialise cuts and flags
717 for(i=0;i<kncuts;i++) cut[i]=-99;
718 for(i=0;i<knflags;i++) flag[i]=-99;
720 for(i=0;i<256;i++) line[i]='\0';
721 // Read up to the end of line excluded
722 iret=fscanf(lun,"%[^\n]",line);
728 // Read the end of line
731 if(line[0]=='*') continue;
733 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",
734 detName,&itmed,&cut[0],&cut[1],&cut[2],&cut[3],&cut[4],&cut[5],&cut[6],&cut[7],&cut[8],
735 &cut[9],&flag[0],&flag[1],&flag[2],&flag[3],&flag[4],&flag[5],&flag[6],&flag[7],
736 &flag[8],&flag[9],&flag[10]);
740 AliWarning(Form("Error reading file %s",fTransParName.Data()));
743 // Check that the module exist
744 AliModule *mod = gAlice->GetModule(detName);
746 // Get the array of media numbers
747 TArrayI &idtmed = *mod->GetIdtmed();
748 // Check that the tracking medium code is valid
749 if(0<=itmed && itmed < 100) {
752 AliWarning(Form("Invalid tracking medium code %d for %s",itmed,mod->GetName()));
755 // Set energy thresholds
756 for(kz=0;kz<kncuts;kz++) {
758 AliDebug(2, Form("%-6s set to %10.3E for tracking medium code %4d for %s",
759 kpars[kz],cut[kz],itmed,mod->GetName()));
760 gMC->Gstpar(ktmed,kpars[kz],cut[kz]);
763 // Set transport mechanisms
764 for(kz=0;kz<knflags;kz++) {
766 AliDebug(2, Form("%-6s set to %10d for tracking medium code %4d for %s",
767 kpars[kncuts+kz],flag[kz],itmed,mod->GetName()));
768 gMC->Gstpar(ktmed,kpars[kncuts+kz],Float_t(flag[kz]));
772 AliWarning(Form("Invalid medium code %d",itmed));
776 AliDebug(1, Form("%s not present",detName));
782 //_______________________________________________________________________
783 void AliMC::SetTransPar(const char *filename)
786 // Sets the file name for transport parameters
788 fTransParName = filename;
791 //_______________________________________________________________________
792 void AliMC::Browse(TBrowser *b)
795 // Called when the item "Run" is clicked on the left pane
796 // of the Root browser.
797 // It displays the Root Trees and all detectors.
799 //detectors are in folders anyway
800 b->Add(fMCQA,"AliMCQA");
804 //_______________________________________________________________________
805 void AliMC::AddHit(Int_t id, Int_t track, Int_t *vol, Float_t *hits) const
808 // Add a hit to detector id
810 TObjArray &dets = *gAlice->Modules();
811 if(dets[id]) dynamic_cast<AliModule*>(dets[id])->AddHit(track,vol,hits);
814 //_______________________________________________________________________
815 void AliMC::AddDigit(Int_t id, Int_t *tracks, Int_t *digits) const
818 // Add digit to detector id
820 TObjArray &dets = *gAlice->Modules();
821 if(dets[id]) dynamic_cast<AliModule*>(dets[id])->AddDigit(tracks,digits);
824 //_______________________________________________________________________
825 Int_t AliMC::GetCurrentTrackNumber() const {
827 // Returns current track
829 return gAlice->GetRunLoader()->Stack()->GetCurrentTrackNumber();
832 //_______________________________________________________________________
833 void AliMC::DumpPart (Int_t i) const
836 // Dumps particle i in the stack
838 AliRunLoader * runloader = gAlice->GetRunLoader();
839 if (runloader->Stack())
840 runloader->Stack()->DumpPart(i);
843 //_______________________________________________________________________
844 void AliMC::DumpPStack () const
847 // Dumps the particle stack
849 AliRunLoader * runloader = gAlice->GetRunLoader();
850 if (runloader->Stack())
851 runloader->Stack()->DumpPStack();
854 //_______________________________________________________________________
855 Int_t AliMC::GetNtrack() const {
857 // Returns number of tracks in stack
860 AliRunLoader * runloader = gAlice->GetRunLoader();
861 if (runloader->Stack())
862 ntracks = runloader->Stack()->GetNtrack();
866 //_______________________________________________________________________
867 Int_t AliMC::GetPrimary(Int_t track) const
870 // return number of primary that has generated track
872 Int_t nprimary = -999;
873 AliRunLoader * runloader = gAlice->GetRunLoader();
874 if (runloader->Stack())
875 nprimary = runloader->Stack()->GetPrimary(track);
879 //_______________________________________________________________________
880 TParticle* AliMC::Particle(Int_t i) const
882 // Returns the i-th particle from the stack taking into account
883 // the remaping done by PurifyKine
884 AliRunLoader * runloader = gAlice->GetRunLoader();
886 if (runloader->Stack())
887 return runloader->Stack()->Particle(i);
891 //_______________________________________________________________________
892 TObjArray* AliMC::Particles() const {
894 // Returns pointer to Particles array
896 AliRunLoader * runloader = gAlice->GetRunLoader();
898 if (runloader->Stack())
899 return runloader->Stack()->Particles();
903 //_______________________________________________________________________
904 void AliMC::PushTrack(Int_t done, Int_t parent, Int_t pdg, Float_t *pmom,
905 Float_t *vpos, Float_t *polar, Float_t tof,
906 TMCProcess mech, Int_t &ntr, Float_t weight, Int_t is) const
910 AliRunLoader * runloader = gAlice->GetRunLoader();
912 if (runloader->Stack())
913 runloader->Stack()->PushTrack(done, parent, pdg, pmom, vpos, polar, tof,
914 mech, ntr, weight, is);
917 //_______________________________________________________________________
918 void AliMC::PushTrack(Int_t done, Int_t parent, Int_t pdg,
919 Double_t px, Double_t py, Double_t pz, Double_t e,
920 Double_t vx, Double_t vy, Double_t vz, Double_t tof,
921 Double_t polx, Double_t poly, Double_t polz,
922 TMCProcess mech, Int_t &ntr, Float_t weight, Int_t is) const
926 AliRunLoader * runloader = gAlice->GetRunLoader();
928 if (runloader->Stack())
929 runloader->Stack()->PushTrack(done, parent, pdg, px, py, pz, e, vx, vy, vz, tof,
930 polx, poly, polz, mech, ntr, weight, is);
933 //_______________________________________________________________________
934 void AliMC::SetHighWaterMark(Int_t nt) const
937 // Set high water mark for last track in event
938 AliRunLoader * runloader = gAlice->GetRunLoader();
940 if (runloader->Stack())
941 runloader->Stack()->SetHighWaterMark(nt);
944 //_______________________________________________________________________
945 void AliMC::KeepTrack(Int_t track) const
950 AliRunLoader * runloader = gAlice->GetRunLoader();
952 if (runloader->Stack())
953 runloader->Stack()->KeepTrack(track);
956 //_______________________________________________________________________
957 void AliMC::FlagTrack(Int_t track) const
961 AliRunLoader * runloader = gAlice->GetRunLoader();
963 if (runloader->Stack())
964 runloader->Stack()->FlagTrack(track);
967 //_______________________________________________________________________
968 void AliMC::SetCurrentTrack(Int_t track) const
971 // Set current track number
973 AliRunLoader * runloader = gAlice->GetRunLoader();
975 if (runloader->Stack())
976 runloader->Stack()->SetCurrentTrack(track);
979 //_______________________________________________________________________
980 void AliMC::AddTrackReference(Int_t label){
982 // add a trackrefernce to the list
983 if (!fTrackReferences) {
984 AliError("Container trackrefernce not active");
987 Int_t nref = fTrackReferences->GetEntriesFast();
988 TClonesArray &lref = *fTrackReferences;
989 new(lref[nref]) AliTrackReference(label);
994 //_______________________________________________________________________
995 void AliMC::ResetTrackReferences()
998 // Reset all references
1000 if (fTrackReferences) fTrackReferences->Clear();
1002 TIter next(gAlice->Modules());
1003 AliModule *detector;
1004 while((detector = dynamic_cast<AliModule*>(next()))) {
1005 detector->ResetTrackReferences();
1009 void AliMC::RemapTrackReferencesIDs(Int_t *map)
1012 // Remapping track reference
1013 // Called at finish primary
1015 if (!fTrackReferences) return;
1016 for (Int_t i=0;i<fTrackReferences->GetEntries();i++){
1017 AliTrackReference * ref = dynamic_cast<AliTrackReference*>(fTrackReferences->UncheckedAt(i));
1019 Int_t newID = map[ref->GetTrack()];
1020 if (newID>=0) ref->SetTrack(newID);
1022 //ref->SetTrack(-1);
1023 ref->SetBit(kNotDeleted,kFALSE);
1024 fTrackReferences->RemoveAt(i);
1028 fTrackReferences->Compress();
1031 void AliMC::FixParticleDecaytime()
1034 // Fix the particle decay time according to rmin and rmax for decays
1038 gMC->TrackMomentum(p);
1039 Double_t tmin, tmax;
1042 // Transverse velocity
1043 Double_t vt = p.Pt() / p.E();
1045 if ((b = gAlice->Field()->SolenoidField()) > 0.) { // [kG]
1049 Double_t rho = p.Pt() / 0.0003 / b; // [cm]
1051 // Revolution frequency
1053 Double_t omega = vt / rho;
1055 // Maximum and minimum decay time
1057 // Check for curlers first
1058 if (fRDecayMax * fRDecayMax / rho / rho / 2. > 1.) return;
1062 tmax = TMath::ACos(1. - fRDecayMax * fRDecayMax / rho / rho / 2.) / omega; // [ct]
1063 tmin = TMath::ACos(1. - fRDecayMin * fRDecayMin / rho / rho / 2.) / omega; // [ct]
1065 tmax = fRDecayMax / vt; // [ct]
1066 tmin = fRDecayMin / vt; // [ct]
1070 // Dial t using the two limits
1071 Double_t t = tmin + (tmax - tmin) * gRandom->Rndm(); // [ct]
1074 // Force decay time in transport code
1076 TGeant3 * geant = (TGeant3*) gMC;
1077 geant->Gcphys()->sumlif = t / p.Beta() / p.Gamma();