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>
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 and display geometry
188 AliInfo("Initialisation:");
190 TIter next(gAlice->Modules());
192 while((detector = dynamic_cast<AliModule*>(next()))) {
194 // Initialise detector and display geometry
195 if(gAlice->IsRootGeometry()) detector->CreateMaterials();
197 detector->BuildGeometry();
198 AliInfo(Form("%10s R:%.2fs C:%.2fs",
199 detector->GetName(),stw.RealTime(),stw.CpuTime()));
204 //_______________________________________________________________________
205 void AliMC::GeneratePrimaries()
208 // Generate primary particles and fill them in the stack.
211 Generator()->Generate();
214 //_______________________________________________________________________
215 void AliMC::SetGenerator(AliGenerator *generator)
218 // Load the event generator
220 if(!fGenerator) fGenerator = generator;
223 //_______________________________________________________________________
224 void AliMC::ResetGenerator(AliGenerator *generator)
227 // Load the event generator
231 AliWarning(Form("Replacing generator %s with %s",
232 fGenerator->GetName(),generator->GetName()))
234 AliWarning(Form("Replacing generator %s with NULL",
235 fGenerator->GetName()));
236 fGenerator = generator;
239 //_______________________________________________________________________
240 void AliMC::FinishRun()
242 // Clean generator information
243 AliDebug(1, "fGenerator->FinishRun()");
244 fGenerator->FinishRun();
246 //Output energy summary tables
247 AliDebug(1, "EnergySummary()");
248 ToAliDebug(1, EnergySummary());
251 //_______________________________________________________________________
252 void AliMC::BeginPrimary()
255 // Called at the beginning of each primary track
260 ResetTrackReferences();
264 //_______________________________________________________________________
265 void AliMC::PreTrack()
267 // Actions before the track's transport
268 TObjArray &dets = *gAlice->Modules();
271 for(Int_t i=0; i<=gAlice->GetNdets(); i++)
272 if((module = dynamic_cast<AliModule*>(dets[i])))
278 //_______________________________________________________________________
279 void AliMC::Stepping()
282 // Called at every step during transport
285 Int_t id = DetFromMate(gMC->GetMedium());
289 if ( gMC->IsNewTrack() &&
290 gMC->TrackTime() == 0. &&
292 fRDecayMax > fRDecayMin &&
293 gMC->TrackPid() == fDecayPdg )
295 FixParticleDecaytime();
301 // --- If lego option, do it and leave
303 gAlice->Lego()->StepManager();
306 //Update energy deposition tables
307 AddEnergyDeposit(gMC->CurrentVolID(copy),gMC->Edep());
309 // write tracke reference for track which is dissapearing - MI
310 if (gMC->IsTrackDisappeared()) {
311 if (gMC->Etot()>0.05) AddTrackReference(GetCurrentTrackNumber());
314 //Call the appropriate stepping routine;
315 AliModule *det = dynamic_cast<AliModule*>(gAlice->Modules()->At(id));
316 if(det && det->StepManagerIsEnabled()) {
317 fMCQA->StepManager(id);
323 //_______________________________________________________________________
324 void AliMC::EnergySummary()
327 // Print summary of deposited energy
333 Int_t kn, i, left, j, id;
334 const Float_t kzero=0;
335 Int_t ievent=gAlice->GetRunLoader()->GetHeader()->GetEvent()+1;
337 // Energy loss information
339 printf("***************** Energy Loss Information per event (GEV) *****************\n");
340 for(kn=1;kn<fEventEnergy.GetSize();kn++) {
343 fEventEnergy[ndep]=kn;
348 ed2=100*TMath::Sqrt(TMath::Max(ed2-ed*ed,kzero))/ed;
351 fSummEnergy[ndep]=ed;
352 fSum2Energy[ndep]=TMath::Min(static_cast<Float_t>(99.),TMath::Max(ed2,kzero));
357 for(kn=0;kn<(ndep-1)/3+1;kn++) {
359 for(i=0;i<(3<left?3:left);i++) {
361 id=Int_t (fEventEnergy[j]+0.1);
362 printf(" %s %10.3f +- %10.3f%%;",gMC->VolName(id),fSummEnergy[j],fSum2Energy[j]);
367 // Relative energy loss in different detectors
368 printf("******************** Relative Energy Loss per event ********************\n");
369 printf("Total energy loss per event %10.3f GeV\n",edtot);
370 for(kn=0;kn<(ndep-1)/5+1;kn++) {
372 for(i=0;i<(5<left?5:left);i++) {
374 id=Int_t (fEventEnergy[j]+0.1);
375 printf(" %s %10.3f%%;",gMC->VolName(id),100*fSummEnergy[j]/edtot);
379 for(kn=0;kn<75;kn++) printf("*");
383 // Reset the TArray's
384 // fEventEnergy.Set(0);
385 // fSummEnergy.Set(0);
386 // fSum2Energy.Set(0);
389 //_____________________________________________________________________________
390 void AliMC::BeginEvent()
393 // Clean-up previous event
395 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
396 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
397 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
398 AliDebug(1, " BEGINNING EVENT ");
399 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
400 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
401 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
403 AliRunLoader *runloader=gAlice->GetRunLoader();
405 /*******************************/
406 /* Clean after eventual */
408 /*******************************/
411 //Set the next event in Run Loader -> Cleans trees (TreeK and all trees in detectors),
412 gAlice->SetEventNrInRun(gAlice->GetEventNrInRun()+1);
413 runloader->SetEventNumber(gAlice->GetEventNrInRun());// sets new files, cleans the previous event stuff, if necessary, etc.,
414 AliDebug(1, Form("EventNr is %d",gAlice->GetEventNrInRun()));
416 fEventEnergy.Reset();
417 // Clean detector information
419 if (runloader->Stack())
420 runloader->Stack()->Reset();//clean stack -> tree is unloaded
422 runloader->MakeStack();//or make a new one
424 if(gAlice->Lego() == 0x0)
426 AliDebug(1, "fRunLoader->MakeTree(K)");
427 runloader->MakeTree("K");
430 AliDebug(1, "gMC->SetStack(fRunLoader->Stack())");
431 gMC->SetStack(gAlice->GetRunLoader()->Stack());//Was in InitMC - but was moved here
432 //because we don't have guarantee that
433 //stack pointer is not going to change from event to event
434 //since it bellobgs to header and is obtained via RunLoader
436 // Reset all Detectors & kinematics & make/reset trees
439 runloader->GetHeader()->Reset(gAlice->GetRunNumber(),gAlice->GetEvNumber(),
440 gAlice->GetEventNrInRun());
441 // fRunLoader->WriteKinematics("OVERWRITE"); is there any reason to rewrite here since MakeTree does so
445 gAlice->Lego()->BeginEvent();
450 AliDebug(1, "ResetHits()");
453 AliDebug(1, "fRunLoader->MakeTree(H)");
454 runloader->MakeTree("H");
456 AliDebug(1, "fRunLoader->MakeTrackRefsContainer()");
457 runloader->MakeTrackRefsContainer();//for insurance
460 //create new branches and SetAdresses
461 TIter next(gAlice->Modules());
463 while((detector = (AliModule*)next()))
465 AliDebug(2, Form("%s->MakeBranch(H)",detector->GetName()));
466 detector->MakeBranch("H");
467 AliDebug(2, Form("%s->MakeBranchTR()",detector->GetName()));
468 detector->MakeBranchTR();
469 AliDebug(2, Form("%s->SetTreeAddress()",detector->GetName()));
470 detector->SetTreeAddress();
472 // make branch for AliRun track References
473 TTree * treeTR = runloader->TreeTR();
475 // make branch for central track references
476 if (!fTrackReferences) fTrackReferences = new TClonesArray("AliTrackReference",0);
478 branch = treeTR->Branch("AliRun",&fTrackReferences);
479 branch->SetAddress(&fTrackReferences);
484 //_______________________________________________________________________
485 void AliMC::ResetHits()
488 // Reset all Detectors hits
490 TIter next(gAlice->Modules());
492 while((detector = dynamic_cast<AliModule*>(next()))) {
493 detector->ResetHits();
497 //_______________________________________________________________________
498 void AliMC::PostTrack()
500 // Posts tracks for each module
501 TObjArray &dets = *gAlice->Modules();
504 for(Int_t i=0; i<=gAlice->GetNdets(); i++)
505 if((module = dynamic_cast<AliModule*>(dets[i])))
509 //_______________________________________________________________________
510 void AliMC::FinishPrimary()
513 // Called at the end of each primary track
515 AliRunLoader *runloader=gAlice->GetRunLoader();
516 // static Int_t count=0;
517 // const Int_t times=10;
518 // This primary is finished, purify stack
519 #if ROOT_VERSION_CODE > 262152
520 if (!(gMC->SecondariesAreOrdered()))
521 runloader->Stack()->ReorderKine();
523 runloader->Stack()->PurifyKine();
525 TIter next(gAlice->Modules());
527 while((detector = dynamic_cast<AliModule*>(next()))) {
528 detector->FinishPrimary();
529 AliLoader* loader = detector->GetLoader();
532 TTree* treeH = loader->TreeH();
533 if (treeH) treeH->Fill(); //can be Lego run and treeH can not exist
537 // Write out track references if any
538 if (runloader->TreeTR()) runloader->TreeTR()->Fill();
541 //_______________________________________________________________________
542 void AliMC::FinishEvent()
545 // Called at the end of the event.
549 if(gAlice->Lego()) gAlice->Lego()->FinishEvent();
551 TIter next(gAlice->Modules());
553 while((detector = dynamic_cast<AliModule*>(next()))) {
554 detector->FinishEvent();
557 //Update the energy deposit tables
559 for(i=0;i<fEventEnergy.GetSize();i++)
561 fSummEnergy[i]+=fEventEnergy[i];
562 fSum2Energy[i]+=fEventEnergy[i]*fEventEnergy[i];
565 AliRunLoader *runloader=gAlice->GetRunLoader();
567 AliHeader* header = runloader->GetHeader();
568 AliStack* stack = runloader->Stack();
569 if ( (header == 0x0) || (stack == 0x0) )
570 {//check if we got header and stack. If not cry and exit aliroot
571 AliFatal("Can not get the stack or header from LOADER");
572 return;//never reached
574 // Update Header information
575 header->SetNprimary(stack->GetNprimary());
576 header->SetNtrack(stack->GetNtrack());
579 // Write out the kinematics
580 stack->FinishEvent();
582 // Write out the event Header information
583 TTree* treeE = runloader->TreeE();
586 header->SetStack(stack);
591 AliError("Can not get TreeE from RL");
594 if(gAlice->Lego() == 0x0)
596 runloader->WriteKinematics("OVERWRITE");
597 runloader->WriteTrackRefs("OVERWRITE");
598 runloader->WriteHits("OVERWRITE");
601 AliDebug(1, "<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<");
602 AliDebug(1, "<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<");
603 AliDebug(1, "<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<");
604 AliDebug(1, " FINISHING EVENT ");
605 AliDebug(1, "<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<");
606 AliDebug(1, "<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<");
607 AliDebug(1, "<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<");
610 //_______________________________________________________________________
611 void AliMC::Field(const Double_t* x, Double_t* b) const
613 // Calculates field "b" at point "x"
617 //_______________________________________________________________________
622 //=================Create Materials and geometry
624 //Read the cuts for all materials
626 //Build the special IMEDIA table
629 //Compute cross-sections
632 //Write Geometry object to current file.
633 gAlice->GetRunLoader()->WriteGeometry();
635 //Initialise geometry deposition table
636 fEventEnergy.Set(gMC->NofVolumes()+1);
637 fSummEnergy.Set(gMC->NofVolumes()+1);
638 fSum2Energy.Set(gMC->NofVolumes()+1);
641 fMCQA = new AliMCQA(gAlice->GetNdets());
643 // Register MC in configuration
644 AliConfig::Instance()->Add(gMC);
648 //_______________________________________________________________________
649 void AliMC::MediaTable()
652 // Built media table to get from the media number to
656 Int_t kz, nz, idt, lz, i, k, ind;
658 TObjArray &dets = *gAlice->Detectors();
660 Int_t ndets=gAlice->GetNdets();
663 for (kz=0;kz<ndets;kz++) {
664 // If detector is defined
665 if((det=dynamic_cast<AliModule*>(dets[kz]))) {
666 TArrayI &idtmed = *(det->GetIdtmed());
667 for(nz=0;nz<100;nz++) {
669 // Find max and min material number
670 if((idt=idtmed[nz])) {
671 det->LoMedium() = det->LoMedium() < idt ? det->LoMedium() : idt;
672 det->HiMedium() = det->HiMedium() > idt ? det->HiMedium() : idt;
675 if(det->LoMedium() > det->HiMedium()) {
679 if(det->HiMedium() > fImedia->GetSize()) {
680 AliError(Form("Increase fImedia from %d to %d",
681 fImedia->GetSize(),det->HiMedium()));
684 // Tag all materials in rage as belonging to detector kz
685 for(lz=det->LoMedium(); lz<= det->HiMedium(); lz++) {
692 // Print summary table
693 AliInfo("Tracking media ranges:");
695 for(i=0;i<(ndets-1)/6+1;i++) {
696 for(k=0;k< (6<ndets-i*6?6:ndets-i*6);k++) {
698 det=dynamic_cast<AliModule*>(dets[ind]);
700 printf(" %6s: %3d -> %3d;",det->GetName(),det->LoMedium(),
703 printf(" %6s: %3d -> %3d;","NULL",0,0);
710 //_______________________________________________________________________
711 void AliMC::ReadTransPar()
714 // Read filename to set the transport parameters
718 const Int_t kncuts=10;
719 const Int_t knflags=11;
720 const Int_t knpars=kncuts+knflags;
721 const char kpars[knpars][7] = {"CUTGAM" ,"CUTELE","CUTNEU","CUTHAD","CUTMUO",
722 "BCUTE","BCUTM","DCUTE","DCUTM","PPCUTM","ANNI",
723 "BREM","COMP","DCAY","DRAY","HADR","LOSS",
724 "MULS","PAIR","PHOT","RAYL"};
730 Int_t i, itmed, iret, ktmed, kz;
733 // See whether the file is there
734 filtmp=gSystem->ExpandPathName(fTransParName.Data());
735 lun=fopen(filtmp,"r");
738 AliWarning(Form("File %s does not exist!",fTransParName.Data()));
743 // Initialise cuts and flags
744 for(i=0;i<kncuts;i++) cut[i]=-99;
745 for(i=0;i<knflags;i++) flag[i]=-99;
747 for(i=0;i<256;i++) line[i]='\0';
748 // Read up to the end of line excluded
749 iret=fscanf(lun,"%[^\n]",line);
755 // Read the end of line
758 if(line[0]=='*') continue;
760 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",
761 detName,&itmed,&cut[0],&cut[1],&cut[2],&cut[3],&cut[4],&cut[5],&cut[6],&cut[7],&cut[8],
762 &cut[9],&flag[0],&flag[1],&flag[2],&flag[3],&flag[4],&flag[5],&flag[6],&flag[7],
763 &flag[8],&flag[9],&flag[10]);
767 AliWarning(Form("Error reading file %s",fTransParName.Data()));
770 // Check that the module exist
771 AliModule *mod = gAlice->GetModule(detName);
773 // Get the array of media numbers
774 TArrayI &idtmed = *mod->GetIdtmed();
775 // Check that the tracking medium code is valid
776 if(0<=itmed && itmed < 100) {
779 AliWarning(Form("Invalid tracking medium code %d for %s",itmed,mod->GetName()));
782 // Set energy thresholds
783 for(kz=0;kz<kncuts;kz++) {
785 AliDebug(2, Form("%-6s set to %10.3E for tracking medium code %4d for %s",
786 kpars[kz],cut[kz],itmed,mod->GetName()));
787 gMC->Gstpar(ktmed,kpars[kz],cut[kz]);
790 // Set transport mechanisms
791 for(kz=0;kz<knflags;kz++) {
793 AliDebug(2, Form("%-6s set to %10d for tracking medium code %4d for %s",
794 kpars[kncuts+kz],flag[kz],itmed,mod->GetName()));
795 gMC->Gstpar(ktmed,kpars[kncuts+kz],Float_t(flag[kz]));
799 AliWarning(Form("Invalid medium code %d",itmed));
803 AliDebug(1, Form("%s not present",detName));
809 //_______________________________________________________________________
810 void AliMC::SetTransPar(const char *filename)
813 // Sets the file name for transport parameters
815 fTransParName = filename;
818 //_______________________________________________________________________
819 void AliMC::Browse(TBrowser *b)
822 // Called when the item "Run" is clicked on the left pane
823 // of the Root browser.
824 // It displays the Root Trees and all detectors.
826 //detectors are in folders anyway
827 b->Add(fMCQA,"AliMCQA");
831 //_______________________________________________________________________
832 void AliMC::AddHit(Int_t id, Int_t track, Int_t *vol, Float_t *hits) const
835 // Add a hit to detector id
837 TObjArray &dets = *gAlice->Modules();
838 if(dets[id]) dynamic_cast<AliModule*>(dets[id])->AddHit(track,vol,hits);
841 //_______________________________________________________________________
842 void AliMC::AddDigit(Int_t id, Int_t *tracks, Int_t *digits) const
845 // Add digit to detector id
847 TObjArray &dets = *gAlice->Modules();
848 if(dets[id]) dynamic_cast<AliModule*>(dets[id])->AddDigit(tracks,digits);
851 //_______________________________________________________________________
852 Int_t AliMC::GetCurrentTrackNumber() const {
854 // Returns current track
856 return gAlice->GetRunLoader()->Stack()->GetCurrentTrackNumber();
859 //_______________________________________________________________________
860 void AliMC::DumpPart (Int_t i) const
863 // Dumps particle i in the stack
865 AliRunLoader * runloader = gAlice->GetRunLoader();
866 if (runloader->Stack())
867 runloader->Stack()->DumpPart(i);
870 //_______________________________________________________________________
871 void AliMC::DumpPStack () const
874 // Dumps the particle stack
876 AliRunLoader * runloader = gAlice->GetRunLoader();
877 if (runloader->Stack())
878 runloader->Stack()->DumpPStack();
881 //_______________________________________________________________________
882 Int_t AliMC::GetNtrack() const {
884 // Returns number of tracks in stack
887 AliRunLoader * runloader = gAlice->GetRunLoader();
888 if (runloader->Stack())
889 ntracks = runloader->Stack()->GetNtrack();
893 //_______________________________________________________________________
894 Int_t AliMC::GetPrimary(Int_t track) const
897 // return number of primary that has generated track
899 Int_t nprimary = -999;
900 AliRunLoader * runloader = gAlice->GetRunLoader();
901 if (runloader->Stack())
902 nprimary = runloader->Stack()->GetPrimary(track);
906 //_______________________________________________________________________
907 TParticle* AliMC::Particle(Int_t i) const
909 // Returns the i-th particle from the stack taking into account
910 // the remaping done by PurifyKine
911 AliRunLoader * runloader = gAlice->GetRunLoader();
913 if (runloader->Stack())
914 return runloader->Stack()->Particle(i);
918 //_______________________________________________________________________
919 TObjArray* AliMC::Particles() const {
921 // Returns pointer to Particles array
923 AliRunLoader * runloader = gAlice->GetRunLoader();
925 if (runloader->Stack())
926 return runloader->Stack()->Particles();
930 //_______________________________________________________________________
931 void AliMC::PushTrack(Int_t done, Int_t parent, Int_t pdg, Float_t *pmom,
932 Float_t *vpos, Float_t *polar, Float_t tof,
933 TMCProcess mech, Int_t &ntr, Float_t weight, Int_t is) const
937 AliRunLoader * runloader = gAlice->GetRunLoader();
939 if (runloader->Stack())
940 runloader->Stack()->PushTrack(done, parent, pdg, pmom, vpos, polar, tof,
941 mech, ntr, weight, is);
944 //_______________________________________________________________________
945 void AliMC::PushTrack(Int_t done, Int_t parent, Int_t pdg,
946 Double_t px, Double_t py, Double_t pz, Double_t e,
947 Double_t vx, Double_t vy, Double_t vz, Double_t tof,
948 Double_t polx, Double_t poly, Double_t polz,
949 TMCProcess mech, Int_t &ntr, Float_t weight, Int_t is) const
953 AliRunLoader * runloader = gAlice->GetRunLoader();
955 if (runloader->Stack())
956 runloader->Stack()->PushTrack(done, parent, pdg, px, py, pz, e, vx, vy, vz, tof,
957 polx, poly, polz, mech, ntr, weight, is);
960 //_______________________________________________________________________
961 void AliMC::SetHighWaterMark(Int_t nt) const
964 // Set high water mark for last track in event
965 AliRunLoader * runloader = gAlice->GetRunLoader();
967 if (runloader->Stack())
968 runloader->Stack()->SetHighWaterMark(nt);
971 //_______________________________________________________________________
972 void AliMC::KeepTrack(Int_t track) const
977 AliRunLoader * runloader = gAlice->GetRunLoader();
979 if (runloader->Stack())
980 runloader->Stack()->KeepTrack(track);
983 //_______________________________________________________________________
984 void AliMC::FlagTrack(Int_t track) const
988 AliRunLoader * runloader = gAlice->GetRunLoader();
990 if (runloader->Stack())
991 runloader->Stack()->FlagTrack(track);
994 //_______________________________________________________________________
995 void AliMC::SetCurrentTrack(Int_t track) const
998 // Set current track number
1000 AliRunLoader * runloader = gAlice->GetRunLoader();
1002 if (runloader->Stack())
1003 runloader->Stack()->SetCurrentTrack(track);
1006 //_______________________________________________________________________
1007 void AliMC::AddTrackReference(Int_t label){
1009 // add a trackrefernce to the list
1010 if (!fTrackReferences) {
1011 AliError("Container trackrefernce not active");
1014 Int_t nref = fTrackReferences->GetEntriesFast();
1015 TClonesArray &lref = *fTrackReferences;
1016 new(lref[nref]) AliTrackReference(label);
1021 //_______________________________________________________________________
1022 void AliMC::ResetTrackReferences()
1025 // Reset all references
1027 if (fTrackReferences) fTrackReferences->Clear();
1029 TIter next(gAlice->Modules());
1030 AliModule *detector;
1031 while((detector = dynamic_cast<AliModule*>(next()))) {
1032 detector->ResetTrackReferences();
1036 void AliMC::RemapTrackReferencesIDs(Int_t *map)
1039 // Remapping track reference
1040 // Called at finish primary
1042 if (!fTrackReferences) return;
1043 for (Int_t i=0;i<fTrackReferences->GetEntries();i++){
1044 AliTrackReference * ref = dynamic_cast<AliTrackReference*>(fTrackReferences->UncheckedAt(i));
1046 Int_t newID = map[ref->GetTrack()];
1047 if (newID>=0) ref->SetTrack(newID);
1049 //ref->SetTrack(-1);
1050 ref->SetBit(kNotDeleted,kFALSE);
1051 fTrackReferences->RemoveAt(i);
1055 fTrackReferences->Compress();
1058 void AliMC::FixParticleDecaytime()
1061 // Fix the particle decay time according to rmin and rmax for decays
1065 gMC->TrackMomentum(p);
1066 Double_t tmin, tmax;
1069 // Transverse velocity
1070 Double_t vt = p.Pt() / p.E();
1072 if ((b = gAlice->Field()->SolenoidField()) > 0.) { // [kG]
1076 Double_t rho = p.Pt() / 0.0003 / b; // [cm]
1078 // Revolution frequency
1080 Double_t omega = vt / rho;
1082 // Maximum and minimum decay time
1084 // Check for curlers first
1085 if (fRDecayMax * fRDecayMax / rho / rho / 2. > 1.) return;
1089 tmax = TMath::ACos(1. - fRDecayMax * fRDecayMax / rho / rho / 2.) / omega; // [ct]
1090 tmin = TMath::ACos(1. - fRDecayMin * fRDecayMin / rho / rho / 2.) / omega; // [ct]
1092 tmax = fRDecayMax / vt; // [ct]
1093 tmin = fRDecayMin / vt; // [ct]
1097 // Dial t using the two limits
1098 Double_t t = tmin + (tmax - tmin) * gRandom->Rndm(); // [ct]
1101 // Force decay time in transport code
1103 gMC->ForceDecayTime(t / 2.99792458e10);