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>
34 #include "AliDetector.h"
35 #include "AliGenerator.h"
36 #include "AliHeader.h"
44 #include "AliTrackReference.h"
45 #include "AliSimulation.h"
46 #include "AliGeomManager.h"
47 #include "AliCDBManager.h"
48 #include "AliCDBStorage.h"
49 #include "AliCDBEntry.h"
54 //_______________________________________________________________________
76 //_______________________________________________________________________
77 AliMC::AliMC(const char *name, const char *title) :
78 TVirtualMCApplication(name, title),
88 fImedia(new TArrayI(1000)),
91 fHitLists(new TList()),
92 fTrackReferences(new TClonesArray("AliTrackReference", 100))
95 // Set transport parameters
98 // Prepare the tracking medium lists
99 for(Int_t i=0;i<1000;i++) (*fImedia)[i]=-99;
102 //_______________________________________________________________________
103 AliMC::AliMC(const AliMC &mc) :
104 TVirtualMCApplication(mc),
121 // Copy constructor for AliMC
126 //_______________________________________________________________________
134 // Delete track references
135 if (fTrackReferences) {
136 fTrackReferences->Delete();
137 delete fTrackReferences;
138 fTrackReferences = 0;
143 //_______________________________________________________________________
144 void AliMC::Copy(TObject &) const
146 //dummy Copy function
147 AliFatal("Not implemented!");
150 //_______________________________________________________________________
151 void AliMC::ConstructGeometry()
154 // Either load geometry from file or create it through usual
155 // loop on detectors. In the first case the method
156 // AliModule::CreateMaterials() only builds fIdtmed and is postponed
157 // at InitGeometry().
160 if(gAlice->IsRootGeometry()){ //load geometry either from CDB or from file
161 if(gAlice->IsGeomFromCDB()){
162 AliInfo("Loading geometry from CDB default storage");
163 AliCDBPath path("GRP","Geometry","Data");
164 AliCDBEntry *entry=AliCDBManager::Instance()->Get(path.GetPath());
165 if(!entry) AliFatal("Unable to load geometry from CDB!");
167 gGeoManager = (TGeoManager*) entry->GetObject();
168 if (!gGeoManager) AliFatal("TGeoManager object not found in the specified CDB entry!");
171 const char *geomfilename = gAlice->GetGeometryFileName();
172 if(gSystem->ExpandPathName(geomfilename)){
173 AliInfo(Form("Loading geometry from file:\n %40s",geomfilename));
174 TGeoManager::Import(geomfilename);
176 AliInfo(Form("Geometry file %40s not found!\n",geomfilename));
181 // Create modules, materials, geometry
183 TIter next(gAlice->Modules());
185 AliDebug(1, "Geometry creation:");
186 while((detector = dynamic_cast<AliModule*>(next()))) {
188 // Initialise detector materials and geometry
189 detector->CreateMaterials();
190 detector->CreateGeometry();
191 AliInfo(Form("%10s R:%.2fs C:%.2fs",
192 detector->GetName(),stw.RealTime(),stw.CpuTime()));
198 //_______________________________________________________________________
199 Bool_t AliMC::MisalignGeometry()
201 // Call misalignment code if AliSimulation object was defined.
203 if(!gAlice->IsRootGeometry()){
204 //Set alignable volumes for the whole geometry
205 SetAllAlignableVolumes();
207 // Misalign geometry via AliSimulation instance
208 if (!AliSimulation::GetInstance()) return kFALSE;
209 AliGeomManager::SetGeometry(gGeoManager);
210 return AliSimulation::GetInstance()->MisalignGeometry(gAlice->GetRunLoader());
213 //_______________________________________________________________________
214 void AliMC::ConstructOpGeometry()
217 // Loop all detector modules and call DefineOpticalProperties() method
220 TIter next(gAlice->Modules());
222 AliInfo("Optical properties definition");
223 while((detector = dynamic_cast<AliModule*>(next()))) {
224 // Initialise detector optical properties
225 detector->DefineOpticalProperties();
229 //_______________________________________________________________________
230 void AliMC::InitGeometry()
233 // Initialize detectors
236 AliInfo("Initialisation:");
238 TIter next(gAlice->Modules());
240 while((detector = dynamic_cast<AliModule*>(next()))) {
242 // Initialise detector geometry
243 if(gAlice->IsRootGeometry()) detector->CreateMaterials();
245 AliInfo(Form("%10s R:%.2fs C:%.2fs",
246 detector->GetName(),stw.RealTime(),stw.CpuTime()));
250 //_______________________________________________________________________
251 void AliMC::SetAllAlignableVolumes()
254 // Add alignable volumes (TGeoPNEntries) looping on all
258 AliInfo(Form("Setting entries for all alignable volumes of active detectors"));
260 TIter next(gAlice->Modules());
261 while((detector = dynamic_cast<AliModule*>(next()))) {
262 detector->AddAlignableVolumes();
266 //_______________________________________________________________________
267 void AliMC::GeneratePrimaries()
270 // Generate primary particles and fill them in the stack.
273 Generator()->Generate();
276 //_______________________________________________________________________
277 void AliMC::SetGenerator(AliGenerator *generator)
280 // Load the event generator
282 if(!fGenerator) fGenerator = generator;
285 //_______________________________________________________________________
286 void AliMC::ResetGenerator(AliGenerator *generator)
289 // Load the event generator
293 AliWarning(Form("Replacing generator %s with %s",
294 fGenerator->GetName(),generator->GetName()));
297 AliWarning(Form("Replacing generator %s with NULL",
298 fGenerator->GetName()));
301 fGenerator = generator;
304 //_______________________________________________________________________
305 void AliMC::FinishRun()
307 // Clean generator information
308 AliDebug(1, "fGenerator->FinishRun()");
309 fGenerator->FinishRun();
311 //Output energy summary tables
312 AliDebug(1, "EnergySummary()");
313 ToAliDebug(1, EnergySummary());
316 //_______________________________________________________________________
317 void AliMC::BeginPrimary()
320 // Called at the beginning of each primary track
325 ResetTrackReferences();
329 //_______________________________________________________________________
330 void AliMC::PreTrack()
332 // Actions before the track's transport
333 TObjArray &dets = *gAlice->Modules();
336 for(Int_t i=0; i<=gAlice->GetNdets(); i++)
337 if((module = dynamic_cast<AliModule*>(dets[i])))
343 //_______________________________________________________________________
344 void AliMC::Stepping()
347 // Called at every step during transport
350 Int_t id = DetFromMate(gMC->CurrentMedium());
354 if ( gMC->IsNewTrack() &&
355 gMC->TrackTime() == 0. &&
357 fRDecayMax > fRDecayMin &&
358 gMC->TrackPid() == fDecayPdg )
360 FixParticleDecaytime();
366 // --- If lego option, do it and leave
368 gAlice->Lego()->StepManager();
371 //Update energy deposition tables
372 AddEnergyDeposit(gMC->CurrentVolID(copy),gMC->Edep());
374 // write tracke reference for track which is dissapearing - MI
375 if (gMC->IsTrackDisappeared()) {
376 if (gMC->Etot()>0.05) AddTrackReference(GetCurrentTrackNumber());
379 //Call the appropriate stepping routine;
380 AliModule *det = dynamic_cast<AliModule*>(gAlice->Modules()->At(id));
381 if(det && det->StepManagerIsEnabled()) {
382 if(AliLog::GetGlobalDebugLevel()>0) fMCQA->StepManager(id);
388 //_______________________________________________________________________
389 void AliMC::EnergySummary()
392 // Print summary of deposited energy
398 Int_t kn, i, left, j, id;
399 const Float_t kzero=0;
400 Int_t ievent=gAlice->GetRunLoader()->GetHeader()->GetEvent()+1;
402 // Energy loss information
404 printf("***************** Energy Loss Information per event (GEV) *****************\n");
405 for(kn=1;kn<fEventEnergy.GetSize();kn++) {
408 fEventEnergy[ndep]=kn;
413 ed2=100*TMath::Sqrt(TMath::Max(ed2-ed*ed,kzero))/ed;
416 fSummEnergy[ndep]=ed;
417 fSum2Energy[ndep]=TMath::Min(static_cast<Float_t>(99.),TMath::Max(ed2,kzero));
422 for(kn=0;kn<(ndep-1)/3+1;kn++) {
424 for(i=0;i<(3<left?3:left);i++) {
426 id=Int_t (fEventEnergy[j]+0.1);
427 printf(" %s %10.3f +- %10.3f%%;",gMC->VolName(id),fSummEnergy[j],fSum2Energy[j]);
432 // Relative energy loss in different detectors
433 printf("******************** Relative Energy Loss per event ********************\n");
434 printf("Total energy loss per event %10.3f GeV\n",edtot);
435 for(kn=0;kn<(ndep-1)/5+1;kn++) {
437 for(i=0;i<(5<left?5:left);i++) {
439 id=Int_t (fEventEnergy[j]+0.1);
440 printf(" %s %10.3f%%;",gMC->VolName(id),100*fSummEnergy[j]/edtot);
444 for(kn=0;kn<75;kn++) printf("*");
448 // Reset the TArray's
449 // fEventEnergy.Set(0);
450 // fSummEnergy.Set(0);
451 // fSum2Energy.Set(0);
454 //_____________________________________________________________________________
455 void AliMC::BeginEvent()
458 // Clean-up previous event
460 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
461 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
462 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
463 AliDebug(1, " BEGINNING EVENT ");
464 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
465 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
466 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
468 AliRunLoader *runloader=gAlice->GetRunLoader();
470 /*******************************/
471 /* Clean after eventual */
473 /*******************************/
476 //Set the next event in Run Loader -> Cleans trees (TreeK and all trees in detectors),
477 gAlice->SetEventNrInRun(gAlice->GetEventNrInRun()+1);
478 runloader->SetEventNumber(gAlice->GetEventNrInRun());// sets new files, cleans the previous event stuff, if necessary, etc.,
479 AliDebug(1, Form("EventNr is %d",gAlice->GetEventNrInRun()));
481 fEventEnergy.Reset();
482 // Clean detector information
484 if (runloader->Stack())
485 runloader->Stack()->Reset();//clean stack -> tree is unloaded
487 runloader->MakeStack();//or make a new one
489 if(gAlice->Lego() == 0x0)
491 AliDebug(1, "fRunLoader->MakeTree(K)");
492 runloader->MakeTree("K");
495 AliDebug(1, "gMC->SetStack(fRunLoader->Stack())");
496 gMC->SetStack(gAlice->GetRunLoader()->Stack());//Was in InitMC - but was moved here
497 //because we don't have guarantee that
498 //stack pointer is not going to change from event to event
499 //since it bellobgs to header and is obtained via RunLoader
501 // Reset all Detectors & kinematics & make/reset trees
504 runloader->GetHeader()->Reset(gAlice->GetRunNumber(),gAlice->GetEvNumber(),
505 gAlice->GetEventNrInRun());
506 // fRunLoader->WriteKinematics("OVERWRITE"); is there any reason to rewrite here since MakeTree does so
510 gAlice->Lego()->BeginEvent();
515 AliDebug(1, "ResetHits()");
518 AliDebug(1, "fRunLoader->MakeTree(H)");
519 runloader->MakeTree("H");
521 AliDebug(1, "fRunLoader->MakeTrackRefsContainer()");
522 runloader->MakeTrackRefsContainer();//for insurance
525 //create new branches and SetAdresses
526 TIter next(gAlice->Modules());
528 while((detector = (AliModule*)next()))
530 AliDebug(2, Form("%s->MakeBranch(H)",detector->GetName()));
531 detector->MakeBranch("H");
532 AliDebug(2, Form("%s->MakeBranchTR()",detector->GetName()));
533 detector->MakeBranchTR();
534 AliDebug(2, Form("%s->SetTreeAddress()",detector->GetName()));
535 detector->SetTreeAddress();
537 // make branch for AliRun track References
538 TTree * treeTR = runloader->TreeTR();
540 // make branch for central track references
541 if (!fTrackReferences) fTrackReferences = new TClonesArray("AliTrackReference",0);
543 branch = treeTR->Branch("AliRun",&fTrackReferences);
544 branch->SetAddress(&fTrackReferences);
549 //_______________________________________________________________________
550 void AliMC::ResetHits()
553 // Reset all Detectors hits
555 TIter next(gAlice->Modules());
557 while((detector = dynamic_cast<AliModule*>(next()))) {
558 detector->ResetHits();
562 //_______________________________________________________________________
563 void AliMC::PostTrack()
565 // Posts tracks for each module
566 TObjArray &dets = *gAlice->Modules();
569 for(Int_t i=0; i<=gAlice->GetNdets(); i++)
570 if((module = dynamic_cast<AliModule*>(dets[i])))
574 //_______________________________________________________________________
575 void AliMC::FinishPrimary()
578 // Called at the end of each primary track
580 AliRunLoader *runloader=gAlice->GetRunLoader();
581 // static Int_t count=0;
582 // const Int_t times=10;
583 // This primary is finished, purify stack
584 #if ROOT_VERSION_CODE > 262152
585 if (!(gMC->SecondariesAreOrdered())) {
586 runloader->Stack()->ReorderKine();
590 runloader->Stack()->PurifyKine();
593 TIter next(gAlice->Modules());
595 while((detector = dynamic_cast<AliModule*>(next()))) {
596 detector->FinishPrimary();
597 AliLoader* loader = detector->GetLoader();
600 TTree* treeH = loader->TreeH();
601 if (treeH) treeH->Fill(); //can be Lego run and treeH can not exist
605 // Write out track references if any
606 if (runloader->TreeTR()) runloader->TreeTR()->Fill();
609 void AliMC::RemapHits()
612 // Remaps the track labels of the hits
613 AliRunLoader *runloader=gAlice->GetRunLoader();
614 AliStack* stack = runloader->Stack();
615 TList* hitLists = GetHitLists();
616 TIter next(hitLists);
617 TCollection *hitList;
619 while((hitList = dynamic_cast<TCollection*>(next()))) {
620 TIter nexthit(hitList);
622 while((hit = dynamic_cast<AliHit*>(nexthit()))) {
623 hit->SetTrack(stack->TrackLabel(hit->GetTrack()));
628 // This for detectors which have a special mapping mechanism
629 // for hits, such as TPC and TRD
632 TObjArray* modules = gAlice->Modules();
633 TIter nextmod(modules);
635 while((detector = dynamic_cast<AliModule*>(nextmod()))) {
636 detector->RemapTrackHitIDs(stack->TrackLabelMap());
637 detector->RemapTrackReferencesIDs(stack->TrackLabelMap());
640 RemapTrackReferencesIDs(stack->TrackLabelMap());
643 //_______________________________________________________________________
644 void AliMC::FinishEvent()
647 // Called at the end of the event.
652 if(gAlice->Lego()) gAlice->Lego()->FinishEvent();
654 TIter next(gAlice->Modules());
656 while((detector = dynamic_cast<AliModule*>(next()))) {
657 detector->FinishEvent();
660 //Update the energy deposit tables
662 for(i=0;i<fEventEnergy.GetSize();i++)
664 fSummEnergy[i]+=fEventEnergy[i];
665 fSum2Energy[i]+=fEventEnergy[i]*fEventEnergy[i];
668 AliRunLoader *runloader=gAlice->GetRunLoader();
670 AliHeader* header = runloader->GetHeader();
671 AliStack* stack = runloader->Stack();
672 if ( (header == 0x0) || (stack == 0x0) )
673 {//check if we got header and stack. If not cry and exit aliroot
674 AliFatal("Can not get the stack or header from LOADER");
675 return;//never reached
677 // Update Header information
678 header->SetNprimary(stack->GetNprimary());
679 header->SetNtrack(stack->GetNtrack());
682 // Write out the kinematics
683 if (!gAlice->Lego()) stack->FinishEvent();
685 // Write out the event Header information
686 TTree* treeE = runloader->TreeE();
689 header->SetStack(stack);
694 AliError("Can not get TreeE from RL");
697 if(gAlice->Lego() == 0x0)
699 runloader->WriteKinematics("OVERWRITE");
700 runloader->WriteTrackRefs("OVERWRITE");
701 runloader->WriteHits("OVERWRITE");
704 AliDebug(1, "<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<");
705 AliDebug(1, "<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<");
706 AliDebug(1, "<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<");
707 AliDebug(1, " FINISHING EVENT ");
708 AliDebug(1, "<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<");
709 AliDebug(1, "<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<");
710 AliDebug(1, "<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<");
713 //_______________________________________________________________________
714 void AliMC::Field(const Double_t* x, Double_t* b) const
716 // Calculates field "b" at point "x"
720 //_______________________________________________________________________
725 //=================Create Materials and geometry
727 // Set alignable volumes for the whole geometry (with old root)
728 #if ROOT_VERSION_CODE < 331527
729 SetAllAlignableVolumes();
731 //Read the cuts for all materials
733 //Build the special IMEDIA table
736 //Compute cross-sections
739 //Initialise geometry deposition table
740 fEventEnergy.Set(gMC->NofVolumes()+1);
741 fSummEnergy.Set(gMC->NofVolumes()+1);
742 fSum2Energy.Set(gMC->NofVolumes()+1);
745 fMCQA = new AliMCQA(gAlice->GetNdets());
747 // Register MC in configuration
748 AliConfig::Instance()->Add(gMC);
752 //_______________________________________________________________________
753 void AliMC::MediaTable()
756 // Built media table to get from the media number to
760 Int_t kz, nz, idt, lz, i, k, ind;
762 TObjArray &dets = *gAlice->Detectors();
764 Int_t ndets=gAlice->GetNdets();
767 for (kz=0;kz<ndets;kz++) {
768 // If detector is defined
769 if((det=dynamic_cast<AliModule*>(dets[kz]))) {
770 TArrayI &idtmed = *(det->GetIdtmed());
771 for(nz=0;nz<100;nz++) {
773 // Find max and min material number
774 if((idt=idtmed[nz])) {
775 det->LoMedium() = det->LoMedium() < idt ? det->LoMedium() : idt;
776 det->HiMedium() = det->HiMedium() > idt ? det->HiMedium() : idt;
779 if(det->LoMedium() > det->HiMedium()) {
783 if(det->HiMedium() > fImedia->GetSize()) {
784 AliError(Form("Increase fImedia from %d to %d",
785 fImedia->GetSize(),det->HiMedium()));
788 // Tag all materials in rage as belonging to detector kz
789 for(lz=det->LoMedium(); lz<= det->HiMedium(); lz++) {
796 // Print summary table
797 AliInfo("Tracking media ranges:");
799 for(i=0;i<(ndets-1)/6+1;i++) {
800 for(k=0;k< (6<ndets-i*6?6:ndets-i*6);k++) {
802 det=dynamic_cast<AliModule*>(dets[ind]);
804 printf(" %6s: %3d -> %3d;",det->GetName(),det->LoMedium(),
807 printf(" %6s: %3d -> %3d;","NULL",0,0);
814 //_______________________________________________________________________
815 void AliMC::ReadTransPar()
818 // Read filename to set the transport parameters
822 const Int_t kncuts=10;
823 const Int_t knflags=11;
824 const Int_t knpars=kncuts+knflags;
825 const char kpars[knpars][7] = {"CUTGAM" ,"CUTELE","CUTNEU","CUTHAD","CUTMUO",
826 "BCUTE","BCUTM","DCUTE","DCUTM","PPCUTM","ANNI",
827 "BREM","COMP","DCAY","DRAY","HADR","LOSS",
828 "MULS","PAIR","PHOT","RAYL"};
834 Int_t i, itmed, iret, ktmed, kz;
837 // See whether the file is there
838 filtmp=gSystem->ExpandPathName(fTransParName.Data());
839 lun=fopen(filtmp,"r");
842 AliWarning(Form("File %s does not exist!",fTransParName.Data()));
847 // Initialise cuts and flags
848 for(i=0;i<kncuts;i++) cut[i]=-99;
849 for(i=0;i<knflags;i++) flag[i]=-99;
851 for(i=0;i<256;i++) line[i]='\0';
852 // Read up to the end of line excluded
853 iret=fscanf(lun,"%[^\n]",line);
859 // Read the end of line
862 if(line[0]=='*') continue;
864 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",
865 detName,&itmed,&cut[0],&cut[1],&cut[2],&cut[3],&cut[4],&cut[5],&cut[6],&cut[7],&cut[8],
866 &cut[9],&flag[0],&flag[1],&flag[2],&flag[3],&flag[4],&flag[5],&flag[6],&flag[7],
867 &flag[8],&flag[9],&flag[10]);
871 AliWarning(Form("Error reading file %s",fTransParName.Data()));
874 // Check that the module exist
875 AliModule *mod = gAlice->GetModule(detName);
877 // Get the array of media numbers
878 TArrayI &idtmed = *mod->GetIdtmed();
879 // Check that the tracking medium code is valid
880 if(0<=itmed && itmed < 100) {
883 AliWarning(Form("Invalid tracking medium code %d for %s",itmed,mod->GetName()));
886 // Set energy thresholds
887 for(kz=0;kz<kncuts;kz++) {
889 AliDebug(2, Form("%-6s set to %10.3E for tracking medium code %4d for %s",
890 kpars[kz],cut[kz],itmed,mod->GetName()));
891 gMC->Gstpar(ktmed,kpars[kz],cut[kz]);
894 // Set transport mechanisms
895 for(kz=0;kz<knflags;kz++) {
897 AliDebug(2, Form("%-6s set to %10d for tracking medium code %4d for %s",
898 kpars[kncuts+kz],flag[kz],itmed,mod->GetName()));
899 gMC->Gstpar(ktmed,kpars[kncuts+kz],Float_t(flag[kz]));
903 AliWarning(Form("Invalid medium code %d",itmed));
907 AliDebug(1, Form("%s not present",detName));
913 //_______________________________________________________________________
914 void AliMC::SetTransPar(const char *filename)
917 // Sets the file name for transport parameters
919 fTransParName = filename;
922 //_______________________________________________________________________
923 void AliMC::Browse(TBrowser *b)
926 // Called when the item "Run" is clicked on the left pane
927 // of the Root browser.
928 // It displays the Root Trees and all detectors.
930 //detectors are in folders anyway
931 b->Add(fMCQA,"AliMCQA");
934 //_______________________________________________________________________
935 void AliMC::AddHit(Int_t id, Int_t track, Int_t *vol, Float_t *hits) const
938 // Add a hit to detector id
940 TObjArray &dets = *gAlice->Modules();
941 if(dets[id]) dynamic_cast<AliModule*>(dets[id])->AddHit(track,vol,hits);
944 //_______________________________________________________________________
945 void AliMC::AddDigit(Int_t id, Int_t *tracks, Int_t *digits) const
948 // Add digit to detector id
950 TObjArray &dets = *gAlice->Modules();
951 if(dets[id]) dynamic_cast<AliModule*>(dets[id])->AddDigit(tracks,digits);
954 //_______________________________________________________________________
955 Int_t AliMC::GetCurrentTrackNumber() const {
957 // Returns current track
959 return gAlice->GetRunLoader()->Stack()->GetCurrentTrackNumber();
962 //_______________________________________________________________________
963 void AliMC::DumpPart (Int_t i) const
966 // Dumps particle i in the stack
968 AliRunLoader * runloader = gAlice->GetRunLoader();
969 if (runloader->Stack())
970 runloader->Stack()->DumpPart(i);
973 //_______________________________________________________________________
974 void AliMC::DumpPStack () const
977 // Dumps the particle stack
979 AliRunLoader * runloader = gAlice->GetRunLoader();
980 if (runloader->Stack())
981 runloader->Stack()->DumpPStack();
984 //_______________________________________________________________________
985 Int_t AliMC::GetNtrack() const {
987 // Returns number of tracks in stack
990 AliRunLoader * runloader = gAlice->GetRunLoader();
991 if (runloader->Stack())
992 ntracks = runloader->Stack()->GetNtrack();
996 //_______________________________________________________________________
997 Int_t AliMC::GetPrimary(Int_t track) const
1000 // return number of primary that has generated track
1002 Int_t nprimary = -999;
1003 AliRunLoader * runloader = gAlice->GetRunLoader();
1004 if (runloader->Stack())
1005 nprimary = runloader->Stack()->GetPrimary(track);
1009 //_______________________________________________________________________
1010 TParticle* AliMC::Particle(Int_t i) const
1012 // Returns the i-th particle from the stack taking into account
1013 // the remaping done by PurifyKine
1014 AliRunLoader * runloader = gAlice->GetRunLoader();
1016 if (runloader->Stack())
1017 return runloader->Stack()->Particle(i);
1021 //_______________________________________________________________________
1022 TObjArray* AliMC::Particles() const {
1024 // Returns pointer to Particles array
1026 AliRunLoader * runloader = gAlice->GetRunLoader();
1028 if (runloader->Stack())
1029 return runloader->Stack()->Particles();
1033 //_______________________________________________________________________
1034 void AliMC::PushTrack(Int_t done, Int_t parent, Int_t pdg, Float_t *pmom,
1035 Float_t *vpos, Float_t *polar, Float_t tof,
1036 TMCProcess mech, Int_t &ntr, Float_t weight, Int_t is) const
1038 // Delegate to stack
1040 AliRunLoader * runloader = gAlice->GetRunLoader();
1042 if (runloader->Stack())
1043 runloader->Stack()->PushTrack(done, parent, pdg, pmom, vpos, polar, tof,
1044 mech, ntr, weight, is);
1047 //_______________________________________________________________________
1048 void AliMC::PushTrack(Int_t done, Int_t parent, Int_t pdg,
1049 Double_t px, Double_t py, Double_t pz, Double_t e,
1050 Double_t vx, Double_t vy, Double_t vz, Double_t tof,
1051 Double_t polx, Double_t poly, Double_t polz,
1052 TMCProcess mech, Int_t &ntr, Float_t weight, Int_t is) const
1054 // Delegate to stack
1056 AliRunLoader * runloader = gAlice->GetRunLoader();
1058 if (runloader->Stack())
1059 runloader->Stack()->PushTrack(done, parent, pdg, px, py, pz, e, vx, vy, vz, tof,
1060 polx, poly, polz, mech, ntr, weight, is);
1063 //_______________________________________________________________________
1064 void AliMC::SetHighWaterMark(Int_t nt) const
1067 // Set high water mark for last track in event
1068 AliRunLoader * runloader = gAlice->GetRunLoader();
1070 if (runloader->Stack())
1071 runloader->Stack()->SetHighWaterMark(nt);
1074 //_______________________________________________________________________
1075 void AliMC::KeepTrack(Int_t track) const
1078 // Delegate to stack
1080 AliRunLoader * runloader = gAlice->GetRunLoader();
1082 if (runloader->Stack())
1083 runloader->Stack()->KeepTrack(track);
1086 //_______________________________________________________________________
1087 void AliMC::FlagTrack(Int_t track) const
1089 // Delegate to stack
1091 AliRunLoader * runloader = gAlice->GetRunLoader();
1093 if (runloader->Stack())
1094 runloader->Stack()->FlagTrack(track);
1097 //_______________________________________________________________________
1098 void AliMC::SetCurrentTrack(Int_t track) const
1101 // Set current track number
1103 AliRunLoader * runloader = gAlice->GetRunLoader();
1105 if (runloader->Stack())
1106 runloader->Stack()->SetCurrentTrack(track);
1109 //_______________________________________________________________________
1110 void AliMC::AddTrackReference(Int_t label){
1112 // add a trackrefernce to the list
1113 if (!fTrackReferences) {
1114 AliError("Container trackrefernce not active");
1117 Int_t nref = fTrackReferences->GetEntriesFast();
1118 TClonesArray &lref = *fTrackReferences;
1119 new(lref[nref]) AliTrackReference(label);
1124 //_______________________________________________________________________
1125 void AliMC::ResetTrackReferences()
1128 // Reset all references
1130 if (fTrackReferences) fTrackReferences->Clear();
1132 TIter next(gAlice->Modules());
1133 AliModule *detector;
1134 while((detector = dynamic_cast<AliModule*>(next()))) {
1135 detector->ResetTrackReferences();
1139 void AliMC::RemapTrackReferencesIDs(Int_t *map)
1142 // Remapping track reference
1143 // Called at finish primary
1145 if (!fTrackReferences) return;
1146 Int_t nEntries = fTrackReferences->GetEntries();
1147 for (Int_t i=0; i < nEntries; i++){
1148 AliTrackReference * ref = dynamic_cast<AliTrackReference*>(fTrackReferences->UncheckedAt(i));
1150 Int_t newID = map[ref->GetTrack()];
1151 if (newID>=0) ref->SetTrack(newID);
1153 ref->SetBit(kNotDeleted,kFALSE);
1154 fTrackReferences->RemoveAt(i);
1158 fTrackReferences->Compress();
1161 void AliMC::FixParticleDecaytime()
1164 // Fix the particle decay time according to rmin and rmax for decays
1168 gMC->TrackMomentum(p);
1169 Double_t tmin, tmax;
1172 // Transverse velocity
1173 Double_t vt = p.Pt() / p.E();
1175 if ((b = gAlice->Field()->SolenoidField()) > 0.) { // [kG]
1179 Double_t rho = p.Pt() / 0.0003 / b; // [cm]
1181 // Revolution frequency
1183 Double_t omega = vt / rho;
1185 // Maximum and minimum decay time
1187 // Check for curlers first
1188 if (fRDecayMax * fRDecayMax / rho / rho / 2. > 1.) return;
1192 tmax = TMath::ACos(1. - fRDecayMax * fRDecayMax / rho / rho / 2.) / omega; // [ct]
1193 tmin = TMath::ACos(1. - fRDecayMin * fRDecayMin / rho / rho / 2.) / omega; // [ct]
1195 tmax = fRDecayMax / vt; // [ct]
1196 tmin = fRDecayMin / vt; // [ct]
1200 // Dial t using the two limits
1201 Double_t t = tmin + (tmax - tmin) * gRandom->Rndm(); // [ct]
1204 // Force decay time in transport code
1206 gMC->ForceDecayTime(t / 2.99792458e10);