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"
43 #include "AliSimulation.h"
48 //_______________________________________________________________________
70 //_______________________________________________________________________
71 AliMC::AliMC(const char *name, const char *title) :
72 TVirtualMCApplication(name, title),
82 fImedia(new TArrayI(1000)),
85 fHitLists(new TList()),
86 fTrackReferences(new TClonesArray("AliTrackReference", 100))
89 // Set transport parameters
92 // Prepare the tracking medium lists
93 for(Int_t i=0;i<1000;i++) (*fImedia)[i]=-99;
96 //_______________________________________________________________________
97 AliMC::AliMC(const AliMC &mc) :
98 TVirtualMCApplication(mc),
115 // Copy constructor for AliMC
120 //_______________________________________________________________________
128 // Delete track references
129 if (fTrackReferences) {
130 fTrackReferences->Delete();
131 delete fTrackReferences;
132 fTrackReferences = 0;
137 //_______________________________________________________________________
138 void AliMC::Copy(TObject &) const
140 //dummy Copy function
141 AliFatal("Not implemented!");
144 //_______________________________________________________________________
145 void AliMC::ConstructGeometry()
148 // Either load geometry from file or create it through usual
149 // loop on detectors. In the first case the method
150 // AliModule::CreateMaterials() only builds fIdtmed and is postponed
151 // at InitGeometry().
154 if(gAlice->IsRootGeometry()){
156 const char *geomfilename = gAlice->GetGeometryFileName();
157 if(gSystem->ExpandPathName(geomfilename)){
158 AliInfo(Form("Loading geometry from file:\n %40s\n\n",geomfilename));
159 TGeoManager::Import(geomfilename);
161 AliInfo(Form("Geometry file %40s not found!\n",geomfilename));
165 // Create modules, materials, geometry
167 TIter next(gAlice->Modules());
169 AliDebug(1, "Geometry creation:");
170 while((detector = dynamic_cast<AliModule*>(next()))) {
172 // Initialise detector materials and geometry
173 detector->CreateMaterials();
174 detector->CreateGeometry();
175 AliInfo(Form("%10s R:%.2fs C:%.2fs",
176 detector->GetName(),stw.RealTime(),stw.CpuTime()));
182 //_______________________________________________________________________
183 Bool_t AliMC::MisalignGeometry()
185 // Call misalignment code if AliSimulation object was defined.
187 //Set alignable volumes for the whole geometry
188 SetAllAlignableVolumes();
189 // Misalign geometry via AliSimulation instance
190 if (!AliSimulation::GetInstance()) return kFALSE;
191 return AliSimulation::GetInstance()->MisalignGeometry(gAlice->GetRunLoader());
194 //_______________________________________________________________________
195 void AliMC::ConstructOpGeometry()
198 // Loop all detector modules and call DefineOpticalProperties() method
201 TIter next(gAlice->Modules());
203 AliInfo("Optical properties definition");
204 while((detector = dynamic_cast<AliModule*>(next()))) {
205 // Initialise detector optical properties
206 detector->DefineOpticalProperties();
210 //_______________________________________________________________________
211 void AliMC::InitGeometry()
214 // Initialize detectors
217 AliInfo("Initialisation:");
219 TIter next(gAlice->Modules());
221 while((detector = dynamic_cast<AliModule*>(next()))) {
223 // Initialise detector geometry
224 if(gAlice->IsRootGeometry()) detector->CreateMaterials();
226 AliInfo(Form("%10s R:%.2fs C:%.2fs",
227 detector->GetName(),stw.RealTime(),stw.CpuTime()));
231 //_______________________________________________________________________
232 void AliMC::SetAllAlignableVolumes()
235 // Add alignable volumes (TGeoPNEntries) looping on all
239 AliInfo(Form("Setting entries for all alignable volumes of active detectors"));
241 TIter next(gAlice->Modules());
242 while((detector = dynamic_cast<AliModule*>(next()))) {
243 detector->AddAlignableVolumes();
247 //_______________________________________________________________________
248 void AliMC::GeneratePrimaries()
251 // Generate primary particles and fill them in the stack.
254 Generator()->Generate();
257 //_______________________________________________________________________
258 void AliMC::SetGenerator(AliGenerator *generator)
261 // Load the event generator
263 if(!fGenerator) fGenerator = generator;
266 //_______________________________________________________________________
267 void AliMC::ResetGenerator(AliGenerator *generator)
270 // Load the event generator
274 AliWarning(Form("Replacing generator %s with %s",
275 fGenerator->GetName(),generator->GetName()));
278 AliWarning(Form("Replacing generator %s with NULL",
279 fGenerator->GetName()));
282 fGenerator = generator;
285 //_______________________________________________________________________
286 void AliMC::FinishRun()
288 // Clean generator information
289 AliDebug(1, "fGenerator->FinishRun()");
290 fGenerator->FinishRun();
292 //Output energy summary tables
293 AliDebug(1, "EnergySummary()");
294 ToAliDebug(1, EnergySummary());
297 //_______________________________________________________________________
298 void AliMC::BeginPrimary()
301 // Called at the beginning of each primary track
306 ResetTrackReferences();
310 //_______________________________________________________________________
311 void AliMC::PreTrack()
313 // Actions before the track's transport
314 TObjArray &dets = *gAlice->Modules();
317 for(Int_t i=0; i<=gAlice->GetNdets(); i++)
318 if((module = dynamic_cast<AliModule*>(dets[i])))
324 //_______________________________________________________________________
325 void AliMC::Stepping()
328 // Called at every step during transport
331 Int_t id = DetFromMate(gMC->CurrentMedium());
335 if ( gMC->IsNewTrack() &&
336 gMC->TrackTime() == 0. &&
338 fRDecayMax > fRDecayMin &&
339 gMC->TrackPid() == fDecayPdg )
341 FixParticleDecaytime();
347 // --- If lego option, do it and leave
349 gAlice->Lego()->StepManager();
352 //Update energy deposition tables
353 AddEnergyDeposit(gMC->CurrentVolID(copy),gMC->Edep());
355 // write tracke reference for track which is dissapearing - MI
356 if (gMC->IsTrackDisappeared()) {
357 if (gMC->Etot()>0.05) AddTrackReference(GetCurrentTrackNumber());
360 //Call the appropriate stepping routine;
361 AliModule *det = dynamic_cast<AliModule*>(gAlice->Modules()->At(id));
362 if(det && det->StepManagerIsEnabled()) {
363 if(AliLog::GetGlobalDebugLevel()>0) fMCQA->StepManager(id);
369 //_______________________________________________________________________
370 void AliMC::EnergySummary()
373 // Print summary of deposited energy
379 Int_t kn, i, left, j, id;
380 const Float_t kzero=0;
381 Int_t ievent=gAlice->GetRunLoader()->GetHeader()->GetEvent()+1;
383 // Energy loss information
385 printf("***************** Energy Loss Information per event (GEV) *****************\n");
386 for(kn=1;kn<fEventEnergy.GetSize();kn++) {
389 fEventEnergy[ndep]=kn;
394 ed2=100*TMath::Sqrt(TMath::Max(ed2-ed*ed,kzero))/ed;
397 fSummEnergy[ndep]=ed;
398 fSum2Energy[ndep]=TMath::Min(static_cast<Float_t>(99.),TMath::Max(ed2,kzero));
403 for(kn=0;kn<(ndep-1)/3+1;kn++) {
405 for(i=0;i<(3<left?3:left);i++) {
407 id=Int_t (fEventEnergy[j]+0.1);
408 printf(" %s %10.3f +- %10.3f%%;",gMC->VolName(id),fSummEnergy[j],fSum2Energy[j]);
413 // Relative energy loss in different detectors
414 printf("******************** Relative Energy Loss per event ********************\n");
415 printf("Total energy loss per event %10.3f GeV\n",edtot);
416 for(kn=0;kn<(ndep-1)/5+1;kn++) {
418 for(i=0;i<(5<left?5:left);i++) {
420 id=Int_t (fEventEnergy[j]+0.1);
421 printf(" %s %10.3f%%;",gMC->VolName(id),100*fSummEnergy[j]/edtot);
425 for(kn=0;kn<75;kn++) printf("*");
429 // Reset the TArray's
430 // fEventEnergy.Set(0);
431 // fSummEnergy.Set(0);
432 // fSum2Energy.Set(0);
435 //_____________________________________________________________________________
436 void AliMC::BeginEvent()
439 // Clean-up previous event
441 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
442 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
443 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
444 AliDebug(1, " BEGINNING EVENT ");
445 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
446 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
447 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
449 AliRunLoader *runloader=gAlice->GetRunLoader();
451 /*******************************/
452 /* Clean after eventual */
454 /*******************************/
457 //Set the next event in Run Loader -> Cleans trees (TreeK and all trees in detectors),
458 gAlice->SetEventNrInRun(gAlice->GetEventNrInRun()+1);
459 runloader->SetEventNumber(gAlice->GetEventNrInRun());// sets new files, cleans the previous event stuff, if necessary, etc.,
460 AliDebug(1, Form("EventNr is %d",gAlice->GetEventNrInRun()));
462 fEventEnergy.Reset();
463 // Clean detector information
465 if (runloader->Stack())
466 runloader->Stack()->Reset();//clean stack -> tree is unloaded
468 runloader->MakeStack();//or make a new one
470 if(gAlice->Lego() == 0x0)
472 AliDebug(1, "fRunLoader->MakeTree(K)");
473 runloader->MakeTree("K");
476 AliDebug(1, "gMC->SetStack(fRunLoader->Stack())");
477 gMC->SetStack(gAlice->GetRunLoader()->Stack());//Was in InitMC - but was moved here
478 //because we don't have guarantee that
479 //stack pointer is not going to change from event to event
480 //since it bellobgs to header and is obtained via RunLoader
482 // Reset all Detectors & kinematics & make/reset trees
485 runloader->GetHeader()->Reset(gAlice->GetRunNumber(),gAlice->GetEvNumber(),
486 gAlice->GetEventNrInRun());
487 // fRunLoader->WriteKinematics("OVERWRITE"); is there any reason to rewrite here since MakeTree does so
491 gAlice->Lego()->BeginEvent();
496 AliDebug(1, "ResetHits()");
499 AliDebug(1, "fRunLoader->MakeTree(H)");
500 runloader->MakeTree("H");
502 AliDebug(1, "fRunLoader->MakeTrackRefsContainer()");
503 runloader->MakeTrackRefsContainer();//for insurance
506 //create new branches and SetAdresses
507 TIter next(gAlice->Modules());
509 while((detector = (AliModule*)next()))
511 AliDebug(2, Form("%s->MakeBranch(H)",detector->GetName()));
512 detector->MakeBranch("H");
513 AliDebug(2, Form("%s->MakeBranchTR()",detector->GetName()));
514 detector->MakeBranchTR();
515 AliDebug(2, Form("%s->SetTreeAddress()",detector->GetName()));
516 detector->SetTreeAddress();
518 // make branch for AliRun track References
519 TTree * treeTR = runloader->TreeTR();
521 // make branch for central track references
522 if (!fTrackReferences) fTrackReferences = new TClonesArray("AliTrackReference",0);
524 branch = treeTR->Branch("AliRun",&fTrackReferences);
525 branch->SetAddress(&fTrackReferences);
530 //_______________________________________________________________________
531 void AliMC::ResetHits()
534 // Reset all Detectors hits
536 TIter next(gAlice->Modules());
538 while((detector = dynamic_cast<AliModule*>(next()))) {
539 detector->ResetHits();
543 //_______________________________________________________________________
544 void AliMC::PostTrack()
546 // Posts tracks for each module
547 TObjArray &dets = *gAlice->Modules();
550 for(Int_t i=0; i<=gAlice->GetNdets(); i++)
551 if((module = dynamic_cast<AliModule*>(dets[i])))
555 //_______________________________________________________________________
556 void AliMC::FinishPrimary()
559 // Called at the end of each primary track
561 AliRunLoader *runloader=gAlice->GetRunLoader();
562 // static Int_t count=0;
563 // const Int_t times=10;
564 // This primary is finished, purify stack
565 #if ROOT_VERSION_CODE > 262152
566 if (!(gMC->SecondariesAreOrdered()))
567 runloader->Stack()->ReorderKine();
569 runloader->Stack()->PurifyKine();
571 TIter next(gAlice->Modules());
573 while((detector = dynamic_cast<AliModule*>(next()))) {
574 detector->FinishPrimary();
575 AliLoader* loader = detector->GetLoader();
578 TTree* treeH = loader->TreeH();
579 if (treeH) treeH->Fill(); //can be Lego run and treeH can not exist
583 // Write out track references if any
584 if (runloader->TreeTR()) runloader->TreeTR()->Fill();
587 //_______________________________________________________________________
588 void AliMC::FinishEvent()
591 // Called at the end of the event.
596 if(gAlice->Lego()) gAlice->Lego()->FinishEvent();
598 TIter next(gAlice->Modules());
600 while((detector = dynamic_cast<AliModule*>(next()))) {
601 detector->FinishEvent();
604 //Update the energy deposit tables
606 for(i=0;i<fEventEnergy.GetSize();i++)
608 fSummEnergy[i]+=fEventEnergy[i];
609 fSum2Energy[i]+=fEventEnergy[i]*fEventEnergy[i];
612 AliRunLoader *runloader=gAlice->GetRunLoader();
614 AliHeader* header = runloader->GetHeader();
615 AliStack* stack = runloader->Stack();
616 if ( (header == 0x0) || (stack == 0x0) )
617 {//check if we got header and stack. If not cry and exit aliroot
618 AliFatal("Can not get the stack or header from LOADER");
619 return;//never reached
621 // Update Header information
622 header->SetNprimary(stack->GetNprimary());
623 header->SetNtrack(stack->GetNtrack());
626 // Write out the kinematics
627 if (!gAlice->Lego()) stack->FinishEvent();
629 // Write out the event Header information
630 TTree* treeE = runloader->TreeE();
633 header->SetStack(stack);
638 AliError("Can not get TreeE from RL");
641 if(gAlice->Lego() == 0x0)
643 runloader->WriteKinematics("OVERWRITE");
644 runloader->WriteTrackRefs("OVERWRITE");
645 runloader->WriteHits("OVERWRITE");
648 AliDebug(1, "<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<");
649 AliDebug(1, "<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<");
650 AliDebug(1, "<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<");
651 AliDebug(1, " FINISHING EVENT ");
652 AliDebug(1, "<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<");
653 AliDebug(1, "<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<");
654 AliDebug(1, "<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<");
657 //_______________________________________________________________________
658 void AliMC::Field(const Double_t* x, Double_t* b) const
660 // Calculates field "b" at point "x"
664 //_______________________________________________________________________
669 //=================Create Materials and geometry
671 // Set alignable volumes for the whole geometry (with old root)
672 #if ROOT_VERSION_CODE < 331527
673 SetAllAlignableVolumes();
675 //Read the cuts for all materials
677 //Build the special IMEDIA table
680 //Compute cross-sections
683 //Initialise geometry deposition table
684 fEventEnergy.Set(gMC->NofVolumes()+1);
685 fSummEnergy.Set(gMC->NofVolumes()+1);
686 fSum2Energy.Set(gMC->NofVolumes()+1);
689 fMCQA = new AliMCQA(gAlice->GetNdets());
691 // Register MC in configuration
692 AliConfig::Instance()->Add(gMC);
696 //_______________________________________________________________________
697 void AliMC::MediaTable()
700 // Built media table to get from the media number to
704 Int_t kz, nz, idt, lz, i, k, ind;
706 TObjArray &dets = *gAlice->Detectors();
708 Int_t ndets=gAlice->GetNdets();
711 for (kz=0;kz<ndets;kz++) {
712 // If detector is defined
713 if((det=dynamic_cast<AliModule*>(dets[kz]))) {
714 TArrayI &idtmed = *(det->GetIdtmed());
715 for(nz=0;nz<100;nz++) {
717 // Find max and min material number
718 if((idt=idtmed[nz])) {
719 det->LoMedium() = det->LoMedium() < idt ? det->LoMedium() : idt;
720 det->HiMedium() = det->HiMedium() > idt ? det->HiMedium() : idt;
723 if(det->LoMedium() > det->HiMedium()) {
727 if(det->HiMedium() > fImedia->GetSize()) {
728 AliError(Form("Increase fImedia from %d to %d",
729 fImedia->GetSize(),det->HiMedium()));
732 // Tag all materials in rage as belonging to detector kz
733 for(lz=det->LoMedium(); lz<= det->HiMedium(); lz++) {
740 // Print summary table
741 AliInfo("Tracking media ranges:");
743 for(i=0;i<(ndets-1)/6+1;i++) {
744 for(k=0;k< (6<ndets-i*6?6:ndets-i*6);k++) {
746 det=dynamic_cast<AliModule*>(dets[ind]);
748 printf(" %6s: %3d -> %3d;",det->GetName(),det->LoMedium(),
751 printf(" %6s: %3d -> %3d;","NULL",0,0);
758 //_______________________________________________________________________
759 void AliMC::ReadTransPar()
762 // Read filename to set the transport parameters
766 const Int_t kncuts=10;
767 const Int_t knflags=11;
768 const Int_t knpars=kncuts+knflags;
769 const char kpars[knpars][7] = {"CUTGAM" ,"CUTELE","CUTNEU","CUTHAD","CUTMUO",
770 "BCUTE","BCUTM","DCUTE","DCUTM","PPCUTM","ANNI",
771 "BREM","COMP","DCAY","DRAY","HADR","LOSS",
772 "MULS","PAIR","PHOT","RAYL"};
778 Int_t i, itmed, iret, ktmed, kz;
781 // See whether the file is there
782 filtmp=gSystem->ExpandPathName(fTransParName.Data());
783 lun=fopen(filtmp,"r");
786 AliWarning(Form("File %s does not exist!",fTransParName.Data()));
791 // Initialise cuts and flags
792 for(i=0;i<kncuts;i++) cut[i]=-99;
793 for(i=0;i<knflags;i++) flag[i]=-99;
795 for(i=0;i<256;i++) line[i]='\0';
796 // Read up to the end of line excluded
797 iret=fscanf(lun,"%[^\n]",line);
803 // Read the end of line
806 if(line[0]=='*') continue;
808 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",
809 detName,&itmed,&cut[0],&cut[1],&cut[2],&cut[3],&cut[4],&cut[5],&cut[6],&cut[7],&cut[8],
810 &cut[9],&flag[0],&flag[1],&flag[2],&flag[3],&flag[4],&flag[5],&flag[6],&flag[7],
811 &flag[8],&flag[9],&flag[10]);
815 AliWarning(Form("Error reading file %s",fTransParName.Data()));
818 // Check that the module exist
819 AliModule *mod = gAlice->GetModule(detName);
821 // Get the array of media numbers
822 TArrayI &idtmed = *mod->GetIdtmed();
823 // Check that the tracking medium code is valid
824 if(0<=itmed && itmed < 100) {
827 AliWarning(Form("Invalid tracking medium code %d for %s",itmed,mod->GetName()));
830 // Set energy thresholds
831 for(kz=0;kz<kncuts;kz++) {
833 AliDebug(2, Form("%-6s set to %10.3E for tracking medium code %4d for %s",
834 kpars[kz],cut[kz],itmed,mod->GetName()));
835 gMC->Gstpar(ktmed,kpars[kz],cut[kz]);
838 // Set transport mechanisms
839 for(kz=0;kz<knflags;kz++) {
841 AliDebug(2, Form("%-6s set to %10d for tracking medium code %4d for %s",
842 kpars[kncuts+kz],flag[kz],itmed,mod->GetName()));
843 gMC->Gstpar(ktmed,kpars[kncuts+kz],Float_t(flag[kz]));
847 AliWarning(Form("Invalid medium code %d",itmed));
851 AliDebug(1, Form("%s not present",detName));
857 //_______________________________________________________________________
858 void AliMC::SetTransPar(const char *filename)
861 // Sets the file name for transport parameters
863 fTransParName = filename;
866 //_______________________________________________________________________
867 void AliMC::Browse(TBrowser *b)
870 // Called when the item "Run" is clicked on the left pane
871 // of the Root browser.
872 // It displays the Root Trees and all detectors.
874 //detectors are in folders anyway
875 b->Add(fMCQA,"AliMCQA");
878 //_______________________________________________________________________
879 void AliMC::AddHit(Int_t id, Int_t track, Int_t *vol, Float_t *hits) const
882 // Add a hit to detector id
884 TObjArray &dets = *gAlice->Modules();
885 if(dets[id]) dynamic_cast<AliModule*>(dets[id])->AddHit(track,vol,hits);
888 //_______________________________________________________________________
889 void AliMC::AddDigit(Int_t id, Int_t *tracks, Int_t *digits) const
892 // Add digit to detector id
894 TObjArray &dets = *gAlice->Modules();
895 if(dets[id]) dynamic_cast<AliModule*>(dets[id])->AddDigit(tracks,digits);
898 //_______________________________________________________________________
899 Int_t AliMC::GetCurrentTrackNumber() const {
901 // Returns current track
903 return gAlice->GetRunLoader()->Stack()->GetCurrentTrackNumber();
906 //_______________________________________________________________________
907 void AliMC::DumpPart (Int_t i) const
910 // Dumps particle i in the stack
912 AliRunLoader * runloader = gAlice->GetRunLoader();
913 if (runloader->Stack())
914 runloader->Stack()->DumpPart(i);
917 //_______________________________________________________________________
918 void AliMC::DumpPStack () const
921 // Dumps the particle stack
923 AliRunLoader * runloader = gAlice->GetRunLoader();
924 if (runloader->Stack())
925 runloader->Stack()->DumpPStack();
928 //_______________________________________________________________________
929 Int_t AliMC::GetNtrack() const {
931 // Returns number of tracks in stack
934 AliRunLoader * runloader = gAlice->GetRunLoader();
935 if (runloader->Stack())
936 ntracks = runloader->Stack()->GetNtrack();
940 //_______________________________________________________________________
941 Int_t AliMC::GetPrimary(Int_t track) const
944 // return number of primary that has generated track
946 Int_t nprimary = -999;
947 AliRunLoader * runloader = gAlice->GetRunLoader();
948 if (runloader->Stack())
949 nprimary = runloader->Stack()->GetPrimary(track);
953 //_______________________________________________________________________
954 TParticle* AliMC::Particle(Int_t i) const
956 // Returns the i-th particle from the stack taking into account
957 // the remaping done by PurifyKine
958 AliRunLoader * runloader = gAlice->GetRunLoader();
960 if (runloader->Stack())
961 return runloader->Stack()->Particle(i);
965 //_______________________________________________________________________
966 TObjArray* AliMC::Particles() const {
968 // Returns pointer to Particles array
970 AliRunLoader * runloader = gAlice->GetRunLoader();
972 if (runloader->Stack())
973 return runloader->Stack()->Particles();
977 //_______________________________________________________________________
978 void AliMC::PushTrack(Int_t done, Int_t parent, Int_t pdg, Float_t *pmom,
979 Float_t *vpos, Float_t *polar, Float_t tof,
980 TMCProcess mech, Int_t &ntr, Float_t weight, Int_t is) const
984 AliRunLoader * runloader = gAlice->GetRunLoader();
986 if (runloader->Stack())
987 runloader->Stack()->PushTrack(done, parent, pdg, pmom, vpos, polar, tof,
988 mech, ntr, weight, is);
991 //_______________________________________________________________________
992 void AliMC::PushTrack(Int_t done, Int_t parent, Int_t pdg,
993 Double_t px, Double_t py, Double_t pz, Double_t e,
994 Double_t vx, Double_t vy, Double_t vz, Double_t tof,
995 Double_t polx, Double_t poly, Double_t polz,
996 TMCProcess mech, Int_t &ntr, Float_t weight, Int_t is) const
1000 AliRunLoader * runloader = gAlice->GetRunLoader();
1002 if (runloader->Stack())
1003 runloader->Stack()->PushTrack(done, parent, pdg, px, py, pz, e, vx, vy, vz, tof,
1004 polx, poly, polz, mech, ntr, weight, is);
1007 //_______________________________________________________________________
1008 void AliMC::SetHighWaterMark(Int_t nt) const
1011 // Set high water mark for last track in event
1012 AliRunLoader * runloader = gAlice->GetRunLoader();
1014 if (runloader->Stack())
1015 runloader->Stack()->SetHighWaterMark(nt);
1018 //_______________________________________________________________________
1019 void AliMC::KeepTrack(Int_t track) const
1022 // Delegate to stack
1024 AliRunLoader * runloader = gAlice->GetRunLoader();
1026 if (runloader->Stack())
1027 runloader->Stack()->KeepTrack(track);
1030 //_______________________________________________________________________
1031 void AliMC::FlagTrack(Int_t track) const
1033 // Delegate to stack
1035 AliRunLoader * runloader = gAlice->GetRunLoader();
1037 if (runloader->Stack())
1038 runloader->Stack()->FlagTrack(track);
1041 //_______________________________________________________________________
1042 void AliMC::SetCurrentTrack(Int_t track) const
1045 // Set current track number
1047 AliRunLoader * runloader = gAlice->GetRunLoader();
1049 if (runloader->Stack())
1050 runloader->Stack()->SetCurrentTrack(track);
1053 //_______________________________________________________________________
1054 void AliMC::AddTrackReference(Int_t label){
1056 // add a trackrefernce to the list
1057 if (!fTrackReferences) {
1058 AliError("Container trackrefernce not active");
1061 Int_t nref = fTrackReferences->GetEntriesFast();
1062 TClonesArray &lref = *fTrackReferences;
1063 new(lref[nref]) AliTrackReference(label);
1068 //_______________________________________________________________________
1069 void AliMC::ResetTrackReferences()
1072 // Reset all references
1074 if (fTrackReferences) fTrackReferences->Clear();
1076 TIter next(gAlice->Modules());
1077 AliModule *detector;
1078 while((detector = dynamic_cast<AliModule*>(next()))) {
1079 detector->ResetTrackReferences();
1083 void AliMC::RemapTrackReferencesIDs(Int_t *map)
1086 // Remapping track reference
1087 // Called at finish primary
1089 if (!fTrackReferences) return;
1090 Int_t nEntries = fTrackReferences->GetEntries();
1091 for (Int_t i=0; i < nEntries; i++){
1092 AliTrackReference * ref = dynamic_cast<AliTrackReference*>(fTrackReferences->UncheckedAt(i));
1094 Int_t newID = map[ref->GetTrack()];
1095 if (newID>=0) ref->SetTrack(newID);
1097 ref->SetBit(kNotDeleted,kFALSE);
1098 fTrackReferences->RemoveAt(i);
1102 fTrackReferences->Compress();
1105 void AliMC::FixParticleDecaytime()
1108 // Fix the particle decay time according to rmin and rmax for decays
1112 gMC->TrackMomentum(p);
1113 Double_t tmin, tmax;
1116 // Transverse velocity
1117 Double_t vt = p.Pt() / p.E();
1119 if ((b = gAlice->Field()->SolenoidField()) > 0.) { // [kG]
1123 Double_t rho = p.Pt() / 0.0003 / b; // [cm]
1125 // Revolution frequency
1127 Double_t omega = vt / rho;
1129 // Maximum and minimum decay time
1131 // Check for curlers first
1132 if (fRDecayMax * fRDecayMax / rho / rho / 2. > 1.) return;
1136 tmax = TMath::ACos(1. - fRDecayMax * fRDecayMax / rho / rho / 2.) / omega; // [ct]
1137 tmin = TMath::ACos(1. - fRDecayMin * fRDecayMin / rho / rho / 2.) / omega; // [ct]
1139 tmax = fRDecayMax / vt; // [ct]
1140 tmin = fRDecayMin / vt; // [ct]
1144 // Dial t using the two limits
1145 Double_t t = tmin + (tmax - tmin) * gRandom->Rndm(); // [ct]
1148 // Force decay time in transport code
1150 gMC->ForceDecayTime(t / 2.99792458e10);