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>
31 #include <TParticle.h>
36 #include "AliDetector.h"
37 #include "AliGenerator.h"
38 #include "AliHeader.h"
46 #include "AliTrackReference.h"
47 #include "AliSimulation.h"
48 #include "AliGeomManager.h"
49 #include "AliCDBManager.h"
50 #include "AliCDBStorage.h"
51 #include "AliCDBEntry.h"
56 //_______________________________________________________________________
74 fTmpTrackReferences(0)
81 //_______________________________________________________________________
82 AliMC::AliMC(const char *name, const char *title) :
83 TVirtualMCApplication(name, title),
93 fImedia(new TArrayI(1000)),
96 fHitLists(new TList()),
99 fTrackReferences(new TClonesArray("AliTrackReference", 100)),
100 fTmpTrackReferences(new TClonesArray("AliTrackReference", 100))
103 // Set transport parameters
106 // Prepare the tracking medium lists
107 for(Int_t i=0;i<1000;i++) (*fImedia)[i]=-99;
110 //_______________________________________________________________________
111 AliMC::AliMC(const AliMC &mc) :
112 TVirtualMCApplication(mc),
129 fTmpTrackReferences(0)
132 // Copy constructor for AliMC
137 //_______________________________________________________________________
145 // Delete track references
146 if (fTrackReferences) {
147 fTrackReferences->Delete();
148 delete fTrackReferences;
149 fTrackReferences = 0;
152 if (fTmpTrackReferences) {
153 fTmpTrackReferences->Delete();
154 delete fTmpTrackReferences;
155 fTmpTrackReferences = 0;
160 //_______________________________________________________________________
161 void AliMC::Copy(TObject &) const
163 //dummy Copy function
164 AliFatal("Not implemented!");
167 //_______________________________________________________________________
168 void AliMC::ConstructGeometry()
171 // Either load geometry from file or create it through usual
172 // loop on detectors. In the first case the method
173 // AliModule::CreateMaterials() only builds fIdtmed and is postponed
174 // at InitGeometry().
177 if(gAlice->IsRootGeometry()){ //load geometry either from CDB or from file
178 if(gAlice->IsGeomFromCDB()){
179 AliInfo("Loading geometry from CDB default storage");
180 AliCDBPath path("GRP","Geometry","Data");
181 AliCDBEntry *entry=AliCDBManager::Instance()->Get(path.GetPath());
182 if(!entry) AliFatal("Unable to load geometry from CDB!");
184 gGeoManager = (TGeoManager*) entry->GetObject();
185 if (!gGeoManager) AliFatal("TGeoManager object not found in the specified CDB entry!");
188 const char *geomfilename = gAlice->GetGeometryFileName();
189 if(gSystem->ExpandPathName(geomfilename)){
190 AliInfo(Form("Loading geometry from file:\n %40s",geomfilename));
191 TGeoManager::Import(geomfilename);
193 AliInfo(Form("Geometry file %40s not found!\n",geomfilename));
198 // Create modules, materials, geometry
200 TIter next(gAlice->Modules());
202 AliDebug(1, "Geometry creation:");
203 while((detector = dynamic_cast<AliModule*>(next()))) {
205 // Initialise detector materials and geometry
206 detector->CreateMaterials();
207 detector->CreateGeometry();
208 AliInfo(Form("%10s R:%.2fs C:%.2fs",
209 detector->GetName(),stw.RealTime(),stw.CpuTime()));
215 //_______________________________________________________________________
216 Bool_t AliMC::MisalignGeometry()
218 // Call misalignment code if AliSimulation object was defined.
220 if(!gAlice->IsRootGeometry()){
221 //Set alignable volumes for the whole geometry
222 SetAllAlignableVolumes();
224 // Misalign geometry via AliSimulation instance
225 if (!AliSimulation::GetInstance()) return kFALSE;
226 AliGeomManager::SetGeometry(gGeoManager);
227 return AliSimulation::GetInstance()->MisalignGeometry(gAlice->GetRunLoader());
230 //_______________________________________________________________________
231 void AliMC::ConstructOpGeometry()
234 // Loop all detector modules and call DefineOpticalProperties() method
237 TIter next(gAlice->Modules());
239 AliInfo("Optical properties definition");
240 while((detector = dynamic_cast<AliModule*>(next()))) {
241 // Initialise detector optical properties
242 detector->DefineOpticalProperties();
246 //_______________________________________________________________________
247 void AliMC::InitGeometry()
250 // Initialize detectors
253 AliInfo("Initialisation:");
255 TIter next(gAlice->Modules());
257 while((detector = dynamic_cast<AliModule*>(next()))) {
259 // Initialise detector geometry
260 if(gAlice->IsRootGeometry()) detector->CreateMaterials();
262 AliInfo(Form("%10s R:%.2fs C:%.2fs",
263 detector->GetName(),stw.RealTime(),stw.CpuTime()));
267 //_______________________________________________________________________
268 void AliMC::SetAllAlignableVolumes()
271 // Add alignable volumes (TGeoPNEntries) looping on all
275 AliInfo(Form("Setting entries for all alignable volumes of active detectors"));
277 TIter next(gAlice->Modules());
278 while((detector = dynamic_cast<AliModule*>(next()))) {
279 detector->AddAlignableVolumes();
283 //_______________________________________________________________________
284 void AliMC::GeneratePrimaries()
287 // Generate primary particles and fill them in the stack.
290 Generator()->Generate();
293 //_______________________________________________________________________
294 void AliMC::SetGenerator(AliGenerator *generator)
297 // Load the event generator
299 if(!fGenerator) fGenerator = generator;
302 //_______________________________________________________________________
303 void AliMC::ResetGenerator(AliGenerator *generator)
306 // Load the event generator
310 AliWarning(Form("Replacing generator %s with %s",
311 fGenerator->GetName(),generator->GetName()));
314 AliWarning(Form("Replacing generator %s with NULL",
315 fGenerator->GetName()));
318 fGenerator = generator;
321 //_______________________________________________________________________
322 void AliMC::FinishRun()
324 // Clean generator information
325 AliDebug(1, "fGenerator->FinishRun()");
326 fGenerator->FinishRun();
328 //Output energy summary tables
329 AliDebug(1, "EnergySummary()");
330 ToAliDebug(1, EnergySummary());
333 //_______________________________________________________________________
334 void AliMC::BeginPrimary()
337 // Called at the beginning of each primary track
342 ResetTrackReferences();
345 //_______________________________________________________________________
346 void AliMC::PreTrack()
348 // Actions before the track's transport
349 TObjArray &dets = *gAlice->Modules();
352 for(Int_t i=0; i<=gAlice->GetNdets(); i++)
353 if((module = dynamic_cast<AliModule*>(dets[i])))
359 //_______________________________________________________________________
360 void AliMC::Stepping()
363 // Called at every step during transport
366 Int_t id = DetFromMate(gMC->CurrentMedium());
370 if ( gMC->IsNewTrack() &&
371 gMC->TrackTime() == 0. &&
373 fRDecayMax > fRDecayMin &&
374 gMC->TrackPid() == fDecayPdg )
376 FixParticleDecaytime();
382 // --- If lego option, do it and leave
384 gAlice->Lego()->StepManager();
387 //Update energy deposition tables
388 AddEnergyDeposit(gMC->CurrentVolID(copy),gMC->Edep());
390 // write tracke reference for track which is dissapearing - MI
391 if (gMC->IsTrackDisappeared()) {
392 if (gMC->Etot()>0.05) AddTrackReference(GetCurrentTrackNumber(),
393 AliTrackReference::kDisappeared);
396 //Call the appropriate stepping routine;
397 AliModule *det = dynamic_cast<AliModule*>(gAlice->Modules()->At(id));
398 if(det && det->StepManagerIsEnabled()) {
399 if(AliLog::GetGlobalDebugLevel()>0) fMCQA->StepManager(id);
405 //_______________________________________________________________________
406 void AliMC::EnergySummary()
409 // Print summary of deposited energy
415 Int_t kn, i, left, j, id;
416 const Float_t kzero=0;
417 Int_t ievent=gAlice->GetRunLoader()->GetHeader()->GetEvent()+1;
419 // Energy loss information
421 printf("***************** Energy Loss Information per event (GEV) *****************\n");
422 for(kn=1;kn<fEventEnergy.GetSize();kn++) {
425 fEventEnergy[ndep]=kn;
430 ed2=100*TMath::Sqrt(TMath::Max(ed2-ed*ed,kzero))/ed;
433 fSummEnergy[ndep]=ed;
434 fSum2Energy[ndep]=TMath::Min(static_cast<Float_t>(99.),TMath::Max(ed2,kzero));
439 for(kn=0;kn<(ndep-1)/3+1;kn++) {
441 for(i=0;i<(3<left?3:left);i++) {
443 id=Int_t (fEventEnergy[j]+0.1);
444 printf(" %s %10.3f +- %10.3f%%;",gMC->VolName(id),fSummEnergy[j],fSum2Energy[j]);
449 // Relative energy loss in different detectors
450 printf("******************** Relative Energy Loss per event ********************\n");
451 printf("Total energy loss per event %10.3f GeV\n",edtot);
452 for(kn=0;kn<(ndep-1)/5+1;kn++) {
454 for(i=0;i<(5<left?5:left);i++) {
456 id=Int_t (fEventEnergy[j]+0.1);
457 printf(" %s %10.3f%%;",gMC->VolName(id),100*fSummEnergy[j]/edtot);
461 for(kn=0;kn<75;kn++) printf("*");
465 // Reset the TArray's
466 // fEventEnergy.Set(0);
467 // fSummEnergy.Set(0);
468 // fSum2Energy.Set(0);
471 //_____________________________________________________________________________
472 void AliMC::BeginEvent()
475 // Clean-up previous event
477 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
478 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
479 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
480 AliDebug(1, " BEGINNING EVENT ");
481 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
482 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
483 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
485 AliRunLoader *runloader=gAlice->GetRunLoader();
487 /*******************************/
488 /* Clean after eventual */
490 /*******************************/
493 //Set the next event in Run Loader -> Cleans trees (TreeK and all trees in detectors),
494 gAlice->SetEventNrInRun(gAlice->GetEventNrInRun()+1);
495 runloader->SetEventNumber(gAlice->GetEventNrInRun());// sets new files, cleans the previous event stuff, if necessary, etc.,
496 AliDebug(1, Form("EventNr is %d",gAlice->GetEventNrInRun()));
498 fEventEnergy.Reset();
499 // Clean detector information
501 if (runloader->Stack())
502 runloader->Stack()->Reset();//clean stack -> tree is unloaded
504 runloader->MakeStack();//or make a new one
507 if(gAlice->Lego() == 0x0)
509 AliDebug(1, "fRunLoader->MakeTree(K)");
510 runloader->MakeTree("K");
513 AliDebug(1, "gMC->SetStack(fRunLoader->Stack())");
514 gMC->SetStack(gAlice->GetRunLoader()->Stack());//Was in InitMC - but was moved here
515 //because we don't have guarantee that
516 //stack pointer is not going to change from event to event
517 //since it bellobgs to header and is obtained via RunLoader
519 // Reset all Detectors & kinematics & make/reset trees
522 runloader->GetHeader()->Reset(gAlice->GetRunNumber(),gAlice->GetEvNumber(),
523 gAlice->GetEventNrInRun());
524 // fRunLoader->WriteKinematics("OVERWRITE"); is there any reason to rewrite here since MakeTree does so
528 gAlice->Lego()->BeginEvent();
533 AliDebug(1, "ResetHits()");
536 AliDebug(1, "fRunLoader->MakeTree(H)");
537 runloader->MakeTree("H");
541 MakeTmpTrackRefsTree();
542 //create new branches and SetAdresses
543 TIter next(gAlice->Modules());
545 while((detector = (AliModule*)next()))
547 AliDebug(2, Form("%s->MakeBranch(H)",detector->GetName()));
548 detector->MakeBranch("H");
552 //_______________________________________________________________________
553 void AliMC::ResetHits()
556 // Reset all Detectors hits
558 TIter next(gAlice->Modules());
560 while((detector = dynamic_cast<AliModule*>(next()))) {
561 detector->ResetHits();
565 //_______________________________________________________________________
566 void AliMC::PostTrack()
568 // Posts tracks for each module
569 TObjArray &dets = *gAlice->Modules();
572 for(Int_t i=0; i<=gAlice->GetNdets(); i++)
573 if((module = dynamic_cast<AliModule*>(dets[i])))
577 //_______________________________________________________________________
578 void AliMC::FinishPrimary()
581 // Called at the end of each primary track
583 AliRunLoader *runloader=gAlice->GetRunLoader();
584 // static Int_t count=0;
585 // const Int_t times=10;
586 // This primary is finished, purify stack
587 #if ROOT_VERSION_CODE > 262152
588 if (!(gMC->SecondariesAreOrdered())) {
589 runloader->Stack()->ReorderKine();
593 runloader->Stack()->PurifyKine();
596 TIter next(gAlice->Modules());
598 while((detector = dynamic_cast<AliModule*>(next()))) {
599 detector->FinishPrimary();
600 AliLoader* loader = detector->GetLoader();
603 TTree* treeH = loader->TreeH();
604 if (treeH) treeH->Fill(); //can be Lego run and treeH can not exist
608 // Write out track references if any
609 if (fTmpTreeTR) fTmpTreeTR->Fill();
612 void AliMC::RemapHits()
615 // Remaps the track labels of the hits
616 AliRunLoader *runloader=gAlice->GetRunLoader();
617 AliStack* stack = runloader->Stack();
618 TList* hitLists = GetHitLists();
619 TIter next(hitLists);
620 TCollection *hitList;
622 while((hitList = dynamic_cast<TCollection*>(next()))) {
623 TIter nexthit(hitList);
625 while((hit = dynamic_cast<AliHit*>(nexthit()))) {
626 hit->SetTrack(stack->TrackLabel(hit->GetTrack()));
631 // This for detectors which have a special mapping mechanism
632 // for hits, such as TPC and TRD
635 TObjArray* modules = gAlice->Modules();
636 TIter nextmod(modules);
637 AliDetector *detector;
638 while((detector = dynamic_cast<AliDetector*>(nextmod()))) {
639 detector->RemapTrackHitIDs(stack->TrackLabelMap());
642 RemapTrackReferencesIDs(stack->TrackLabelMap());
645 //_______________________________________________________________________
646 void AliMC::FinishEvent()
649 // Called at the end of the event.
654 if(gAlice->Lego()) gAlice->Lego()->FinishEvent();
656 TIter next(gAlice->Modules());
658 while((detector = dynamic_cast<AliModule*>(next()))) {
659 detector->FinishEvent();
662 //Update the energy deposit tables
664 for(i=0;i<fEventEnergy.GetSize();i++)
666 fSummEnergy[i]+=fEventEnergy[i];
667 fSum2Energy[i]+=fEventEnergy[i]*fEventEnergy[i];
670 AliRunLoader *runloader=gAlice->GetRunLoader();
672 AliHeader* header = runloader->GetHeader();
673 AliStack* stack = runloader->Stack();
674 if ( (header == 0x0) || (stack == 0x0) )
675 {//check if we got header and stack. If not cry and exit aliroot
676 AliFatal("Can not get the stack or header from LOADER");
677 return;//never reached
679 // Update Header information
680 header->SetNprimary(stack->GetNprimary());
681 header->SetNtrack(stack->GetNtrack());
683 // Write out the kinematics
684 if (!gAlice->Lego()) stack->FinishEvent();
686 // Synchronize the TreeTR with TreeK
687 ReorderAndExpandTreeTR();
689 // Write out the event Header information
690 TTree* treeE = runloader->TreeE();
693 header->SetStack(stack);
698 AliError("Can not get TreeE from RL");
701 if(gAlice->Lego() == 0x0)
703 runloader->WriteKinematics("OVERWRITE");
704 runloader->WriteTrackRefs("OVERWRITE");
705 runloader->WriteHits("OVERWRITE");
708 AliDebug(1, "<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<");
709 AliDebug(1, "<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<");
710 AliDebug(1, "<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<");
711 AliDebug(1, " FINISHING EVENT ");
712 AliDebug(1, "<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<");
713 AliDebug(1, "<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<");
714 AliDebug(1, "<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<");
717 //_______________________________________________________________________
718 void AliMC::Field(const Double_t* x, Double_t* b) const
720 // Calculates field "b" at point "x"
724 //_______________________________________________________________________
729 //=================Create Materials and geometry
731 // Set alignable volumes for the whole geometry (with old root)
732 #if ROOT_VERSION_CODE < 331527
733 SetAllAlignableVolumes();
735 //Read the cuts for all materials
737 //Build the special IMEDIA table
740 //Compute cross-sections
743 //Initialise geometry deposition table
744 fEventEnergy.Set(gMC->NofVolumes()+1);
745 fSummEnergy.Set(gMC->NofVolumes()+1);
746 fSum2Energy.Set(gMC->NofVolumes()+1);
749 fMCQA = new AliMCQA(gAlice->GetNdets());
751 // Register MC in configuration
752 AliConfig::Instance()->Add(gMC);
756 //_______________________________________________________________________
757 void AliMC::MediaTable()
760 // Built media table to get from the media number to
764 Int_t kz, nz, idt, lz, i, k, ind;
766 TObjArray &dets = *gAlice->Detectors();
768 Int_t ndets=gAlice->GetNdets();
771 for (kz=0;kz<ndets;kz++) {
772 // If detector is defined
773 if((det=dynamic_cast<AliModule*>(dets[kz]))) {
774 TArrayI &idtmed = *(det->GetIdtmed());
775 for(nz=0;nz<100;nz++) {
777 // Find max and min material number
778 if((idt=idtmed[nz])) {
779 det->LoMedium() = det->LoMedium() < idt ? det->LoMedium() : idt;
780 det->HiMedium() = det->HiMedium() > idt ? det->HiMedium() : idt;
783 if(det->LoMedium() > det->HiMedium()) {
787 if(det->HiMedium() > fImedia->GetSize()) {
788 AliError(Form("Increase fImedia from %d to %d",
789 fImedia->GetSize(),det->HiMedium()));
792 // Tag all materials in rage as belonging to detector kz
793 for(lz=det->LoMedium(); lz<= det->HiMedium(); lz++) {
800 // Print summary table
801 AliInfo("Tracking media ranges:");
803 for(i=0;i<(ndets-1)/6+1;i++) {
804 for(k=0;k< (6<ndets-i*6?6:ndets-i*6);k++) {
806 det=dynamic_cast<AliModule*>(dets[ind]);
808 printf(" %6s: %3d -> %3d;",det->GetName(),det->LoMedium(),
811 printf(" %6s: %3d -> %3d;","NULL",0,0);
818 //_______________________________________________________________________
819 void AliMC::ReadTransPar()
822 // Read filename to set the transport parameters
826 const Int_t kncuts=10;
827 const Int_t knflags=11;
828 const Int_t knpars=kncuts+knflags;
829 const char kpars[knpars][7] = {"CUTGAM" ,"CUTELE","CUTNEU","CUTHAD","CUTMUO",
830 "BCUTE","BCUTM","DCUTE","DCUTM","PPCUTM","ANNI",
831 "BREM","COMP","DCAY","DRAY","HADR","LOSS",
832 "MULS","PAIR","PHOT","RAYL"};
838 Int_t i, itmed, iret, ktmed, kz;
841 // See whether the file is there
842 filtmp=gSystem->ExpandPathName(fTransParName.Data());
843 lun=fopen(filtmp,"r");
846 AliWarning(Form("File %s does not exist!",fTransParName.Data()));
851 // Initialise cuts and flags
852 for(i=0;i<kncuts;i++) cut[i]=-99;
853 for(i=0;i<knflags;i++) flag[i]=-99;
855 for(i=0;i<256;i++) line[i]='\0';
856 // Read up to the end of line excluded
857 iret=fscanf(lun,"%[^\n]",line);
863 // Read the end of line
866 if(line[0]=='*') continue;
868 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",
869 detName,&itmed,&cut[0],&cut[1],&cut[2],&cut[3],&cut[4],&cut[5],&cut[6],&cut[7],&cut[8],
870 &cut[9],&flag[0],&flag[1],&flag[2],&flag[3],&flag[4],&flag[5],&flag[6],&flag[7],
871 &flag[8],&flag[9],&flag[10]);
875 AliWarning(Form("Error reading file %s",fTransParName.Data()));
878 // Check that the module exist
879 AliModule *mod = gAlice->GetModule(detName);
881 // Get the array of media numbers
882 TArrayI &idtmed = *mod->GetIdtmed();
883 // Check that the tracking medium code is valid
884 if(0<=itmed && itmed < 100) {
887 AliWarning(Form("Invalid tracking medium code %d for %s",itmed,mod->GetName()));
890 // Set energy thresholds
891 for(kz=0;kz<kncuts;kz++) {
893 AliDebug(2, Form("%-6s set to %10.3E for tracking medium code %4d for %s",
894 kpars[kz],cut[kz],itmed,mod->GetName()));
895 gMC->Gstpar(ktmed,kpars[kz],cut[kz]);
898 // Set transport mechanisms
899 for(kz=0;kz<knflags;kz++) {
901 AliDebug(2, Form("%-6s set to %10d for tracking medium code %4d for %s",
902 kpars[kncuts+kz],flag[kz],itmed,mod->GetName()));
903 gMC->Gstpar(ktmed,kpars[kncuts+kz],Float_t(flag[kz]));
907 AliWarning(Form("Invalid medium code %d",itmed));
911 AliDebug(1, Form("%s not present",detName));
917 //_______________________________________________________________________
918 void AliMC::SetTransPar(const char *filename)
921 // Sets the file name for transport parameters
923 fTransParName = filename;
926 //_______________________________________________________________________
927 void AliMC::Browse(TBrowser *b)
930 // Called when the item "Run" is clicked on the left pane
931 // of the Root browser.
932 // It displays the Root Trees and all detectors.
934 //detectors are in folders anyway
935 b->Add(fMCQA,"AliMCQA");
938 //_______________________________________________________________________
939 void AliMC::AddHit(Int_t id, Int_t track, Int_t *vol, Float_t *hits) const
942 // Add a hit to detector id
944 TObjArray &dets = *gAlice->Modules();
945 if(dets[id]) dynamic_cast<AliModule*>(dets[id])->AddHit(track,vol,hits);
948 //_______________________________________________________________________
949 void AliMC::AddDigit(Int_t id, Int_t *tracks, Int_t *digits) const
952 // Add digit to detector id
954 TObjArray &dets = *gAlice->Modules();
955 if(dets[id]) dynamic_cast<AliModule*>(dets[id])->AddDigit(tracks,digits);
958 //_______________________________________________________________________
959 Int_t AliMC::GetCurrentTrackNumber() const {
961 // Returns current track
963 return gAlice->GetRunLoader()->Stack()->GetCurrentTrackNumber();
966 //_______________________________________________________________________
967 void AliMC::DumpPart (Int_t i) const
970 // Dumps particle i in the stack
972 AliRunLoader * runloader = gAlice->GetRunLoader();
973 if (runloader->Stack())
974 runloader->Stack()->DumpPart(i);
977 //_______________________________________________________________________
978 void AliMC::DumpPStack () const
981 // Dumps the particle stack
983 AliRunLoader * runloader = gAlice->GetRunLoader();
984 if (runloader->Stack())
985 runloader->Stack()->DumpPStack();
988 //_______________________________________________________________________
989 Int_t AliMC::GetNtrack() const {
991 // Returns number of tracks in stack
994 AliRunLoader * runloader = gAlice->GetRunLoader();
995 if (runloader->Stack())
996 ntracks = runloader->Stack()->GetNtrack();
1000 //_______________________________________________________________________
1001 Int_t AliMC::GetPrimary(Int_t track) const
1004 // return number of primary that has generated track
1006 Int_t nprimary = -999;
1007 AliRunLoader * runloader = gAlice->GetRunLoader();
1008 if (runloader->Stack())
1009 nprimary = runloader->Stack()->GetPrimary(track);
1013 //_______________________________________________________________________
1014 TParticle* AliMC::Particle(Int_t i) const
1016 // Returns the i-th particle from the stack taking into account
1017 // the remaping done by PurifyKine
1018 AliRunLoader * runloader = gAlice->GetRunLoader();
1020 if (runloader->Stack())
1021 return runloader->Stack()->Particle(i);
1025 //_______________________________________________________________________
1026 TObjArray* AliMC::Particles() const {
1028 // Returns pointer to Particles array
1030 AliRunLoader * runloader = gAlice->GetRunLoader();
1032 if (runloader->Stack())
1033 return runloader->Stack()->Particles();
1037 //_______________________________________________________________________
1038 void AliMC::PushTrack(Int_t done, Int_t parent, Int_t pdg, Float_t *pmom,
1039 Float_t *vpos, Float_t *polar, Float_t tof,
1040 TMCProcess mech, Int_t &ntr, Float_t weight, Int_t is) const
1042 // Delegate to stack
1044 AliRunLoader * runloader = gAlice->GetRunLoader();
1046 if (runloader->Stack())
1047 runloader->Stack()->PushTrack(done, parent, pdg, pmom, vpos, polar, tof,
1048 mech, ntr, weight, is);
1051 //_______________________________________________________________________
1052 void AliMC::PushTrack(Int_t done, Int_t parent, Int_t pdg,
1053 Double_t px, Double_t py, Double_t pz, Double_t e,
1054 Double_t vx, Double_t vy, Double_t vz, Double_t tof,
1055 Double_t polx, Double_t poly, Double_t polz,
1056 TMCProcess mech, Int_t &ntr, Float_t weight, Int_t is) const
1058 // Delegate to stack
1060 AliRunLoader * runloader = gAlice->GetRunLoader();
1062 if (runloader->Stack())
1063 runloader->Stack()->PushTrack(done, parent, pdg, px, py, pz, e, vx, vy, vz, tof,
1064 polx, poly, polz, mech, ntr, weight, is);
1067 //_______________________________________________________________________
1068 void AliMC::SetHighWaterMark(Int_t nt) const
1071 // Set high water mark for last track in event
1072 AliRunLoader * runloader = gAlice->GetRunLoader();
1074 if (runloader->Stack())
1075 runloader->Stack()->SetHighWaterMark(nt);
1078 //_______________________________________________________________________
1079 void AliMC::KeepTrack(Int_t track) const
1082 // Delegate to stack
1084 AliRunLoader * runloader = gAlice->GetRunLoader();
1086 if (runloader->Stack())
1087 runloader->Stack()->KeepTrack(track);
1090 //_______________________________________________________________________
1091 void AliMC::FlagTrack(Int_t track) const
1093 // Delegate to stack
1095 AliRunLoader * runloader = gAlice->GetRunLoader();
1097 if (runloader->Stack())
1098 runloader->Stack()->FlagTrack(track);
1101 //_______________________________________________________________________
1102 void AliMC::SetCurrentTrack(Int_t track) const
1105 // Set current track number
1107 AliRunLoader * runloader = gAlice->GetRunLoader();
1109 if (runloader->Stack())
1110 runloader->Stack()->SetCurrentTrack(track);
1113 //_______________________________________________________________________
1114 AliTrackReference* AliMC::AddTrackReference(Int_t label, Int_t id)
1117 // add a trackrefernce to the list
1118 if (!fTrackReferences) {
1119 AliError("Container trackrefernce not active");
1122 Int_t primary = GetPrimary(label);
1123 Particle(primary)->SetBit(kKeepBit);
1125 Int_t nref = fTmpTrackReferences->GetEntriesFast();
1126 TClonesArray &lref = *fTmpTrackReferences;
1127 return new(lref[nref]) AliTrackReference(label, id);
1132 //_______________________________________________________________________
1133 void AliMC::ResetTrackReferences()
1136 // Reset all references
1138 if (fTmpTrackReferences) fTmpTrackReferences->Clear();
1141 void AliMC::RemapTrackReferencesIDs(Int_t *map)
1144 // Remapping track reference
1145 // Called at finish primary
1147 if (!fTmpTrackReferences) return;
1148 Int_t nEntries = fTmpTrackReferences->GetEntries();
1149 for (Int_t i=0; i < nEntries; i++){
1150 AliTrackReference * ref = dynamic_cast<AliTrackReference*>(fTmpTrackReferences->UncheckedAt(i));
1152 Int_t newID = map[ref->GetTrack()];
1153 if (newID>=0) ref->SetTrack(newID);
1155 ref->SetBit(kNotDeleted,kFALSE);
1156 fTmpTrackReferences->RemoveAt(i);
1160 fTmpTrackReferences->Compress();
1163 void AliMC::FixParticleDecaytime()
1166 // Fix the particle decay time according to rmin and rmax for decays
1170 gMC->TrackMomentum(p);
1171 Double_t tmin, tmax;
1174 // Transverse velocity
1175 Double_t vt = p.Pt() / p.E();
1177 if ((b = gAlice->Field()->SolenoidField()) > 0.) { // [kG]
1181 Double_t rho = p.Pt() / 0.0003 / b; // [cm]
1183 // Revolution frequency
1185 Double_t omega = vt / rho;
1187 // Maximum and minimum decay time
1189 // Check for curlers first
1190 if (fRDecayMax * fRDecayMax / rho / rho / 2. > 1.) return;
1194 tmax = TMath::ACos(1. - fRDecayMax * fRDecayMax / rho / rho / 2.) / omega; // [ct]
1195 tmin = TMath::ACos(1. - fRDecayMin * fRDecayMin / rho / rho / 2.) / omega; // [ct]
1197 tmax = fRDecayMax / vt; // [ct]
1198 tmin = fRDecayMin / vt; // [ct]
1202 // Dial t using the two limits
1203 Double_t t = tmin + (tmax - tmin) * gRandom->Rndm(); // [ct]
1206 // Force decay time in transport code
1208 gMC->ForceDecayTime(t / 2.99792458e10);
1211 void AliMC::MakeTmpTrackRefsTree()
1213 // Make the temporary track reference tree
1214 fTmpFileTR = new TFile("TrackRefsTmp.root", "recreate");
1215 fTmpTreeTR = new TTree("TreeTR", "Track References");
1216 if (!fTmpTrackReferences) fTmpTrackReferences = new TClonesArray("AliTrackReference", 100);
1217 fTmpTreeTR->Branch("TrackReferences", "TClonesArray", &fTmpTrackReferences, 4000);
1220 void AliMC::ReorderAndExpandTreeTR()
1223 // Reorder and expand the temporary track reference tree in order to match the kinematics tree
1226 AliRunLoader *rl = gAlice->GetRunLoader();
1229 AliDebug(1, "fRunLoader->MakeTrackRefsContainer()");
1230 rl->MakeTrackRefsContainer();
1231 TTree * treeTR = rl->TreeTR();
1233 // make branch for central track references
1234 if (!fTrackReferences) fTrackReferences = new TClonesArray("AliTrackReference",0);
1236 branch = treeTR->Branch("TrackReferences",&fTrackReferences);
1237 branch->SetAddress(&fTrackReferences);
1240 AliStack* stack = rl->Stack();
1241 Int_t np = stack->GetNprimary();
1242 Int_t nt = fTmpTreeTR->GetEntries();
1244 // Loop over tracks and find the secondaries with the help of the kine tree
1247 for (Int_t ip = np - 1; ip > -1; ip--) {
1248 TParticle *part = stack->Particle(ip);
1249 // Skip primaries that have not been transported
1250 Int_t dau1 = part->GetFirstDaughter();
1252 if (dau1 > -1 && dau1 < np) continue;
1254 fTmpTreeTR->GetEntry(it++);
1255 Int_t nh = fTmpTrackReferences->GetEntries();
1256 // Determine range of secondaries produced by this primary
1258 Int_t inext = ip - 1;
1261 part = stack->Particle(inext);
1262 dau2 = part->GetFirstDaughter();
1263 if (dau2 == -1 || dau2 < np) {
1269 dau2 = stack->GetNtrack() - 1;
1272 } // find upper bound
1275 // Loop over reference hits and find secondary label
1276 for (Int_t id = dau1; (id <= dau2) && (dau1 > -1); id++) {
1277 for (Int_t ih = 0; ih < nh; ih++) {
1278 AliTrackReference* tr = (AliTrackReference*) fTmpTrackReferences->At(ih);
1279 Int_t label = tr->Label();
1281 if (label == ip) continue;
1282 if (label > dau2 || label < dau1) AliWarning("Track Reference Label out of range !");
1285 Int_t nref = fTrackReferences->GetEntriesFast();
1286 TClonesArray &lref = *fTrackReferences;
1287 new(lref[nref]) AliTrackReference(*tr);
1291 fTrackReferences->Clear();
1296 // Now loop again and write the primaries
1298 for (Int_t ip = 0; ip < np; ip++) {
1299 TParticle* part = stack->Particle(ip);
1300 if (part->GetFirstDaughter() == -1 || part->GetFirstDaughter() >= np) {
1301 // Skip particles that have not been transported
1302 fTmpTreeTR->GetEntry(it--);
1303 Int_t nh = fTmpTrackReferences->GetEntries();
1305 // Loop over reference hits and find primary labels
1306 for (Int_t ih = 0; ih < nh; ih++) {
1307 AliTrackReference* tr = (AliTrackReference*) fTmpTrackReferences->At(ih);
1308 Int_t label = tr->Label();
1310 Int_t nref = fTrackReferences->GetEntriesFast();
1311 TClonesArray &lref = *fTrackReferences;
1312 new(lref[nref]) AliTrackReference(*tr);
1317 fTrackReferences->Clear();
1321 if (ifills != stack->GetNtrack())
1322 AliWarning(Form("Number of entries in TreeTR (%5d) unequal to TreeK (%5d) \n", ifills, stack->GetNtrack()));
1326 fTmpFileTR->Close();
1328 delete fTmpTrackReferences;
1329 fTmpTrackReferences = 0;
1330 gSystem->Exec("rm -rf TrackRefsTmp.root");