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()));
203 //_______________________________________________________________________
204 void AliMC::SetAllAlignableVolumes()
207 // Add alignable volumes (TGeoPNEntries) looping on all
211 AliInfo(Form("Setting entries for all alignable volumes of active detectors"));
213 TIter next(gAlice->Modules());
214 while((detector = dynamic_cast<AliModule*>(next()))) {
215 detector->AddAlignableVolumes();
219 //_______________________________________________________________________
220 void AliMC::GeneratePrimaries()
223 // Generate primary particles and fill them in the stack.
226 Generator()->Generate();
229 //_______________________________________________________________________
230 void AliMC::SetGenerator(AliGenerator *generator)
233 // Load the event generator
235 if(!fGenerator) fGenerator = generator;
238 //_______________________________________________________________________
239 void AliMC::ResetGenerator(AliGenerator *generator)
242 // Load the event generator
246 AliWarning(Form("Replacing generator %s with %s",
247 fGenerator->GetName(),generator->GetName()))
249 AliWarning(Form("Replacing generator %s with NULL",
250 fGenerator->GetName()));
251 fGenerator = generator;
254 //_______________________________________________________________________
255 void AliMC::FinishRun()
257 // Clean generator information
258 AliDebug(1, "fGenerator->FinishRun()");
259 fGenerator->FinishRun();
261 //Output energy summary tables
262 AliDebug(1, "EnergySummary()");
263 ToAliDebug(1, EnergySummary());
266 //_______________________________________________________________________
267 void AliMC::BeginPrimary()
270 // Called at the beginning of each primary track
275 ResetTrackReferences();
279 //_______________________________________________________________________
280 void AliMC::PreTrack()
282 // Actions before the track's transport
283 TObjArray &dets = *gAlice->Modules();
286 for(Int_t i=0; i<=gAlice->GetNdets(); i++)
287 if((module = dynamic_cast<AliModule*>(dets[i])))
293 //_______________________________________________________________________
294 void AliMC::Stepping()
297 // Called at every step during transport
300 Int_t id = DetFromMate(gMC->CurrentMedium());
304 if ( gMC->IsNewTrack() &&
305 gMC->TrackTime() == 0. &&
307 fRDecayMax > fRDecayMin &&
308 gMC->TrackPid() == fDecayPdg )
310 FixParticleDecaytime();
316 // --- If lego option, do it and leave
318 gAlice->Lego()->StepManager();
321 //Update energy deposition tables
322 AddEnergyDeposit(gMC->CurrentVolID(copy),gMC->Edep());
324 // write tracke reference for track which is dissapearing - MI
325 if (gMC->IsTrackDisappeared()) {
326 if (gMC->Etot()>0.05) AddTrackReference(GetCurrentTrackNumber());
329 //Call the appropriate stepping routine;
330 AliModule *det = dynamic_cast<AliModule*>(gAlice->Modules()->At(id));
331 if(det && det->StepManagerIsEnabled()) {
332 if(AliLog::GetGlobalDebugLevel()>0) fMCQA->StepManager(id);
338 //_______________________________________________________________________
339 void AliMC::EnergySummary()
342 // Print summary of deposited energy
348 Int_t kn, i, left, j, id;
349 const Float_t kzero=0;
350 Int_t ievent=gAlice->GetRunLoader()->GetHeader()->GetEvent()+1;
352 // Energy loss information
354 printf("***************** Energy Loss Information per event (GEV) *****************\n");
355 for(kn=1;kn<fEventEnergy.GetSize();kn++) {
358 fEventEnergy[ndep]=kn;
363 ed2=100*TMath::Sqrt(TMath::Max(ed2-ed*ed,kzero))/ed;
366 fSummEnergy[ndep]=ed;
367 fSum2Energy[ndep]=TMath::Min(static_cast<Float_t>(99.),TMath::Max(ed2,kzero));
372 for(kn=0;kn<(ndep-1)/3+1;kn++) {
374 for(i=0;i<(3<left?3:left);i++) {
376 id=Int_t (fEventEnergy[j]+0.1);
377 printf(" %s %10.3f +- %10.3f%%;",gMC->VolName(id),fSummEnergy[j],fSum2Energy[j]);
382 // Relative energy loss in different detectors
383 printf("******************** Relative Energy Loss per event ********************\n");
384 printf("Total energy loss per event %10.3f GeV\n",edtot);
385 for(kn=0;kn<(ndep-1)/5+1;kn++) {
387 for(i=0;i<(5<left?5:left);i++) {
389 id=Int_t (fEventEnergy[j]+0.1);
390 printf(" %s %10.3f%%;",gMC->VolName(id),100*fSummEnergy[j]/edtot);
394 for(kn=0;kn<75;kn++) printf("*");
398 // Reset the TArray's
399 // fEventEnergy.Set(0);
400 // fSummEnergy.Set(0);
401 // fSum2Energy.Set(0);
404 //_____________________________________________________________________________
405 void AliMC::BeginEvent()
408 // Clean-up previous event
410 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
411 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
412 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
413 AliDebug(1, " BEGINNING EVENT ");
414 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
415 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
416 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
418 AliRunLoader *runloader=gAlice->GetRunLoader();
420 /*******************************/
421 /* Clean after eventual */
423 /*******************************/
426 //Set the next event in Run Loader -> Cleans trees (TreeK and all trees in detectors),
427 gAlice->SetEventNrInRun(gAlice->GetEventNrInRun()+1);
428 runloader->SetEventNumber(gAlice->GetEventNrInRun());// sets new files, cleans the previous event stuff, if necessary, etc.,
429 AliDebug(1, Form("EventNr is %d",gAlice->GetEventNrInRun()));
431 fEventEnergy.Reset();
432 // Clean detector information
434 if (runloader->Stack())
435 runloader->Stack()->Reset();//clean stack -> tree is unloaded
437 runloader->MakeStack();//or make a new one
439 if(gAlice->Lego() == 0x0)
441 AliDebug(1, "fRunLoader->MakeTree(K)");
442 runloader->MakeTree("K");
445 AliDebug(1, "gMC->SetStack(fRunLoader->Stack())");
446 gMC->SetStack(gAlice->GetRunLoader()->Stack());//Was in InitMC - but was moved here
447 //because we don't have guarantee that
448 //stack pointer is not going to change from event to event
449 //since it bellobgs to header and is obtained via RunLoader
451 // Reset all Detectors & kinematics & make/reset trees
454 runloader->GetHeader()->Reset(gAlice->GetRunNumber(),gAlice->GetEvNumber(),
455 gAlice->GetEventNrInRun());
456 // fRunLoader->WriteKinematics("OVERWRITE"); is there any reason to rewrite here since MakeTree does so
460 gAlice->Lego()->BeginEvent();
465 AliDebug(1, "ResetHits()");
468 AliDebug(1, "fRunLoader->MakeTree(H)");
469 runloader->MakeTree("H");
471 AliDebug(1, "fRunLoader->MakeTrackRefsContainer()");
472 runloader->MakeTrackRefsContainer();//for insurance
475 //create new branches and SetAdresses
476 TIter next(gAlice->Modules());
478 while((detector = (AliModule*)next()))
480 AliDebug(2, Form("%s->MakeBranch(H)",detector->GetName()));
481 detector->MakeBranch("H");
482 AliDebug(2, Form("%s->MakeBranchTR()",detector->GetName()));
483 detector->MakeBranchTR();
484 AliDebug(2, Form("%s->SetTreeAddress()",detector->GetName()));
485 detector->SetTreeAddress();
487 // make branch for AliRun track References
488 TTree * treeTR = runloader->TreeTR();
490 // make branch for central track references
491 if (!fTrackReferences) fTrackReferences = new TClonesArray("AliTrackReference",0);
493 branch = treeTR->Branch("AliRun",&fTrackReferences);
494 branch->SetAddress(&fTrackReferences);
499 //_______________________________________________________________________
500 void AliMC::ResetHits()
503 // Reset all Detectors hits
505 TIter next(gAlice->Modules());
507 while((detector = dynamic_cast<AliModule*>(next()))) {
508 detector->ResetHits();
512 //_______________________________________________________________________
513 void AliMC::PostTrack()
515 // Posts tracks for each module
516 TObjArray &dets = *gAlice->Modules();
519 for(Int_t i=0; i<=gAlice->GetNdets(); i++)
520 if((module = dynamic_cast<AliModule*>(dets[i])))
524 //_______________________________________________________________________
525 void AliMC::FinishPrimary()
528 // Called at the end of each primary track
530 AliRunLoader *runloader=gAlice->GetRunLoader();
531 // static Int_t count=0;
532 // const Int_t times=10;
533 // This primary is finished, purify stack
534 #if ROOT_VERSION_CODE > 262152
535 if (!(gMC->SecondariesAreOrdered()))
536 runloader->Stack()->ReorderKine();
538 runloader->Stack()->PurifyKine();
540 TIter next(gAlice->Modules());
542 while((detector = dynamic_cast<AliModule*>(next()))) {
543 detector->FinishPrimary();
544 AliLoader* loader = detector->GetLoader();
547 TTree* treeH = loader->TreeH();
548 if (treeH) treeH->Fill(); //can be Lego run and treeH can not exist
552 // Write out track references if any
553 if (runloader->TreeTR()) runloader->TreeTR()->Fill();
556 //_______________________________________________________________________
557 void AliMC::FinishEvent()
560 // Called at the end of the event.
565 if(gAlice->Lego()) gAlice->Lego()->FinishEvent();
567 TIter next(gAlice->Modules());
569 while((detector = dynamic_cast<AliModule*>(next()))) {
570 detector->FinishEvent();
573 //Update the energy deposit tables
575 for(i=0;i<fEventEnergy.GetSize();i++)
577 fSummEnergy[i]+=fEventEnergy[i];
578 fSum2Energy[i]+=fEventEnergy[i]*fEventEnergy[i];
581 AliRunLoader *runloader=gAlice->GetRunLoader();
583 AliHeader* header = runloader->GetHeader();
584 AliStack* stack = runloader->Stack();
585 if ( (header == 0x0) || (stack == 0x0) )
586 {//check if we got header and stack. If not cry and exit aliroot
587 AliFatal("Can not get the stack or header from LOADER");
588 return;//never reached
590 // Update Header information
591 header->SetNprimary(stack->GetNprimary());
592 header->SetNtrack(stack->GetNtrack());
595 // Write out the kinematics
596 if (!gAlice->Lego()) stack->FinishEvent();
598 // Write out the event Header information
599 TTree* treeE = runloader->TreeE();
602 header->SetStack(stack);
607 AliError("Can not get TreeE from RL");
610 if(gAlice->Lego() == 0x0)
612 runloader->WriteKinematics("OVERWRITE");
613 runloader->WriteTrackRefs("OVERWRITE");
614 runloader->WriteHits("OVERWRITE");
617 AliDebug(1, "<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<");
618 AliDebug(1, "<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<");
619 AliDebug(1, "<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<");
620 AliDebug(1, " FINISHING EVENT ");
621 AliDebug(1, "<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<");
622 AliDebug(1, "<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<");
623 AliDebug(1, "<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<");
626 //_______________________________________________________________________
627 void AliMC::Field(const Double_t* x, Double_t* b) const
629 // Calculates field "b" at point "x"
633 //_______________________________________________________________________
638 //=================Create Materials and geometry
640 //Set alignable volumes for the whole geometry
641 SetAllAlignableVolumes();
642 //Read the cuts for all materials
644 //Build the special IMEDIA table
647 //Compute cross-sections
650 //Write Geometry object to current file.
651 gAlice->GetRunLoader()->WriteGeometry();
653 //Initialise geometry deposition table
654 fEventEnergy.Set(gMC->NofVolumes()+1);
655 fSummEnergy.Set(gMC->NofVolumes()+1);
656 fSum2Energy.Set(gMC->NofVolumes()+1);
659 fMCQA = new AliMCQA(gAlice->GetNdets());
661 // Register MC in configuration
662 AliConfig::Instance()->Add(gMC);
666 //_______________________________________________________________________
667 void AliMC::MediaTable()
670 // Built media table to get from the media number to
674 Int_t kz, nz, idt, lz, i, k, ind;
676 TObjArray &dets = *gAlice->Detectors();
678 Int_t ndets=gAlice->GetNdets();
681 for (kz=0;kz<ndets;kz++) {
682 // If detector is defined
683 if((det=dynamic_cast<AliModule*>(dets[kz]))) {
684 TArrayI &idtmed = *(det->GetIdtmed());
685 for(nz=0;nz<100;nz++) {
687 // Find max and min material number
688 if((idt=idtmed[nz])) {
689 det->LoMedium() = det->LoMedium() < idt ? det->LoMedium() : idt;
690 det->HiMedium() = det->HiMedium() > idt ? det->HiMedium() : idt;
693 if(det->LoMedium() > det->HiMedium()) {
697 if(det->HiMedium() > fImedia->GetSize()) {
698 AliError(Form("Increase fImedia from %d to %d",
699 fImedia->GetSize(),det->HiMedium()));
702 // Tag all materials in rage as belonging to detector kz
703 for(lz=det->LoMedium(); lz<= det->HiMedium(); lz++) {
710 // Print summary table
711 AliInfo("Tracking media ranges:");
713 for(i=0;i<(ndets-1)/6+1;i++) {
714 for(k=0;k< (6<ndets-i*6?6:ndets-i*6);k++) {
716 det=dynamic_cast<AliModule*>(dets[ind]);
718 printf(" %6s: %3d -> %3d;",det->GetName(),det->LoMedium(),
721 printf(" %6s: %3d -> %3d;","NULL",0,0);
728 //_______________________________________________________________________
729 void AliMC::ReadTransPar()
732 // Read filename to set the transport parameters
736 const Int_t kncuts=10;
737 const Int_t knflags=11;
738 const Int_t knpars=kncuts+knflags;
739 const char kpars[knpars][7] = {"CUTGAM" ,"CUTELE","CUTNEU","CUTHAD","CUTMUO",
740 "BCUTE","BCUTM","DCUTE","DCUTM","PPCUTM","ANNI",
741 "BREM","COMP","DCAY","DRAY","HADR","LOSS",
742 "MULS","PAIR","PHOT","RAYL"};
748 Int_t i, itmed, iret, ktmed, kz;
751 // See whether the file is there
752 filtmp=gSystem->ExpandPathName(fTransParName.Data());
753 lun=fopen(filtmp,"r");
756 AliWarning(Form("File %s does not exist!",fTransParName.Data()));
761 // Initialise cuts and flags
762 for(i=0;i<kncuts;i++) cut[i]=-99;
763 for(i=0;i<knflags;i++) flag[i]=-99;
765 for(i=0;i<256;i++) line[i]='\0';
766 // Read up to the end of line excluded
767 iret=fscanf(lun,"%[^\n]",line);
773 // Read the end of line
776 if(line[0]=='*') continue;
778 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",
779 detName,&itmed,&cut[0],&cut[1],&cut[2],&cut[3],&cut[4],&cut[5],&cut[6],&cut[7],&cut[8],
780 &cut[9],&flag[0],&flag[1],&flag[2],&flag[3],&flag[4],&flag[5],&flag[6],&flag[7],
781 &flag[8],&flag[9],&flag[10]);
785 AliWarning(Form("Error reading file %s",fTransParName.Data()));
788 // Check that the module exist
789 AliModule *mod = gAlice->GetModule(detName);
791 // Get the array of media numbers
792 TArrayI &idtmed = *mod->GetIdtmed();
793 // Check that the tracking medium code is valid
794 if(0<=itmed && itmed < 100) {
797 AliWarning(Form("Invalid tracking medium code %d for %s",itmed,mod->GetName()));
800 // Set energy thresholds
801 for(kz=0;kz<kncuts;kz++) {
803 AliDebug(2, Form("%-6s set to %10.3E for tracking medium code %4d for %s",
804 kpars[kz],cut[kz],itmed,mod->GetName()));
805 gMC->Gstpar(ktmed,kpars[kz],cut[kz]);
808 // Set transport mechanisms
809 for(kz=0;kz<knflags;kz++) {
811 AliDebug(2, Form("%-6s set to %10d for tracking medium code %4d for %s",
812 kpars[kncuts+kz],flag[kz],itmed,mod->GetName()));
813 gMC->Gstpar(ktmed,kpars[kncuts+kz],Float_t(flag[kz]));
817 AliWarning(Form("Invalid medium code %d",itmed));
821 AliDebug(1, Form("%s not present",detName));
827 //_______________________________________________________________________
828 void AliMC::SetTransPar(const char *filename)
831 // Sets the file name for transport parameters
833 fTransParName = filename;
836 //_______________________________________________________________________
837 void AliMC::Browse(TBrowser *b)
840 // Called when the item "Run" is clicked on the left pane
841 // of the Root browser.
842 // It displays the Root Trees and all detectors.
844 //detectors are in folders anyway
845 b->Add(fMCQA,"AliMCQA");
849 //_______________________________________________________________________
850 void AliMC::AddHit(Int_t id, Int_t track, Int_t *vol, Float_t *hits) const
853 // Add a hit to detector id
855 TObjArray &dets = *gAlice->Modules();
856 if(dets[id]) dynamic_cast<AliModule*>(dets[id])->AddHit(track,vol,hits);
859 //_______________________________________________________________________
860 void AliMC::AddDigit(Int_t id, Int_t *tracks, Int_t *digits) const
863 // Add digit to detector id
865 TObjArray &dets = *gAlice->Modules();
866 if(dets[id]) dynamic_cast<AliModule*>(dets[id])->AddDigit(tracks,digits);
869 //_______________________________________________________________________
870 Int_t AliMC::GetCurrentTrackNumber() const {
872 // Returns current track
874 return gAlice->GetRunLoader()->Stack()->GetCurrentTrackNumber();
877 //_______________________________________________________________________
878 void AliMC::DumpPart (Int_t i) const
881 // Dumps particle i in the stack
883 AliRunLoader * runloader = gAlice->GetRunLoader();
884 if (runloader->Stack())
885 runloader->Stack()->DumpPart(i);
888 //_______________________________________________________________________
889 void AliMC::DumpPStack () const
892 // Dumps the particle stack
894 AliRunLoader * runloader = gAlice->GetRunLoader();
895 if (runloader->Stack())
896 runloader->Stack()->DumpPStack();
899 //_______________________________________________________________________
900 Int_t AliMC::GetNtrack() const {
902 // Returns number of tracks in stack
905 AliRunLoader * runloader = gAlice->GetRunLoader();
906 if (runloader->Stack())
907 ntracks = runloader->Stack()->GetNtrack();
911 //_______________________________________________________________________
912 Int_t AliMC::GetPrimary(Int_t track) const
915 // return number of primary that has generated track
917 Int_t nprimary = -999;
918 AliRunLoader * runloader = gAlice->GetRunLoader();
919 if (runloader->Stack())
920 nprimary = runloader->Stack()->GetPrimary(track);
924 //_______________________________________________________________________
925 TParticle* AliMC::Particle(Int_t i) const
927 // Returns the i-th particle from the stack taking into account
928 // the remaping done by PurifyKine
929 AliRunLoader * runloader = gAlice->GetRunLoader();
931 if (runloader->Stack())
932 return runloader->Stack()->Particle(i);
936 //_______________________________________________________________________
937 TObjArray* AliMC::Particles() const {
939 // Returns pointer to Particles array
941 AliRunLoader * runloader = gAlice->GetRunLoader();
943 if (runloader->Stack())
944 return runloader->Stack()->Particles();
948 //_______________________________________________________________________
949 void AliMC::PushTrack(Int_t done, Int_t parent, Int_t pdg, Float_t *pmom,
950 Float_t *vpos, Float_t *polar, Float_t tof,
951 TMCProcess mech, Int_t &ntr, Float_t weight, Int_t is) const
955 AliRunLoader * runloader = gAlice->GetRunLoader();
957 if (runloader->Stack())
958 runloader->Stack()->PushTrack(done, parent, pdg, pmom, vpos, polar, tof,
959 mech, ntr, weight, is);
962 //_______________________________________________________________________
963 void AliMC::PushTrack(Int_t done, Int_t parent, Int_t pdg,
964 Double_t px, Double_t py, Double_t pz, Double_t e,
965 Double_t vx, Double_t vy, Double_t vz, Double_t tof,
966 Double_t polx, Double_t poly, Double_t polz,
967 TMCProcess mech, Int_t &ntr, Float_t weight, Int_t is) const
971 AliRunLoader * runloader = gAlice->GetRunLoader();
973 if (runloader->Stack())
974 runloader->Stack()->PushTrack(done, parent, pdg, px, py, pz, e, vx, vy, vz, tof,
975 polx, poly, polz, mech, ntr, weight, is);
978 //_______________________________________________________________________
979 void AliMC::SetHighWaterMark(Int_t nt) const
982 // Set high water mark for last track in event
983 AliRunLoader * runloader = gAlice->GetRunLoader();
985 if (runloader->Stack())
986 runloader->Stack()->SetHighWaterMark(nt);
989 //_______________________________________________________________________
990 void AliMC::KeepTrack(Int_t track) const
995 AliRunLoader * runloader = gAlice->GetRunLoader();
997 if (runloader->Stack())
998 runloader->Stack()->KeepTrack(track);
1001 //_______________________________________________________________________
1002 void AliMC::FlagTrack(Int_t track) const
1004 // Delegate to stack
1006 AliRunLoader * runloader = gAlice->GetRunLoader();
1008 if (runloader->Stack())
1009 runloader->Stack()->FlagTrack(track);
1012 //_______________________________________________________________________
1013 void AliMC::SetCurrentTrack(Int_t track) const
1016 // Set current track number
1018 AliRunLoader * runloader = gAlice->GetRunLoader();
1020 if (runloader->Stack())
1021 runloader->Stack()->SetCurrentTrack(track);
1024 //_______________________________________________________________________
1025 void AliMC::AddTrackReference(Int_t label){
1027 // add a trackrefernce to the list
1028 if (!fTrackReferences) {
1029 AliError("Container trackrefernce not active");
1032 Int_t nref = fTrackReferences->GetEntriesFast();
1033 TClonesArray &lref = *fTrackReferences;
1034 new(lref[nref]) AliTrackReference(label);
1039 //_______________________________________________________________________
1040 void AliMC::ResetTrackReferences()
1043 // Reset all references
1045 if (fTrackReferences) fTrackReferences->Clear();
1047 TIter next(gAlice->Modules());
1048 AliModule *detector;
1049 while((detector = dynamic_cast<AliModule*>(next()))) {
1050 detector->ResetTrackReferences();
1054 void AliMC::RemapTrackReferencesIDs(Int_t *map)
1057 // Remapping track reference
1058 // Called at finish primary
1060 if (!fTrackReferences) return;
1061 Int_t nEntries = fTrackReferences->GetEntries();
1062 for (Int_t i=0; i < nEntries; i++){
1063 AliTrackReference * ref = dynamic_cast<AliTrackReference*>(fTrackReferences->UncheckedAt(i));
1065 Int_t newID = map[ref->GetTrack()];
1066 if (newID>=0) ref->SetTrack(newID);
1068 ref->SetBit(kNotDeleted,kFALSE);
1069 fTrackReferences->RemoveAt(i);
1073 fTrackReferences->Compress();
1076 void AliMC::FixParticleDecaytime()
1079 // Fix the particle decay time according to rmin and rmax for decays
1083 gMC->TrackMomentum(p);
1084 Double_t tmin, tmax;
1087 // Transverse velocity
1088 Double_t vt = p.Pt() / p.E();
1090 if ((b = gAlice->Field()->SolenoidField()) > 0.) { // [kG]
1094 Double_t rho = p.Pt() / 0.0003 / b; // [cm]
1096 // Revolution frequency
1098 Double_t omega = vt / rho;
1100 // Maximum and minimum decay time
1102 // Check for curlers first
1103 if (fRDecayMax * fRDecayMax / rho / rho / 2. > 1.) return;
1107 tmax = TMath::ACos(1. - fRDecayMax * fRDecayMax / rho / rho / 2.) / omega; // [ct]
1108 tmin = TMath::ACos(1. - fRDecayMin * fRDecayMin / rho / rho / 2.) / omega; // [ct]
1110 tmax = fRDecayMax / vt; // [ct]
1111 tmin = fRDecayMin / vt; // [ct]
1115 // Dial t using the two limits
1116 Double_t t = tmin + (tmax - tmin) * gRandom->Rndm(); // [ct]
1119 // Force decay time in transport code
1121 gMC->ForceDecayTime(t / 2.99792458e10);