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>
28 #include <TGeoGlobalMagField.h>
29 #include <TGeoManager.h>
30 #include <TParticle.h>
32 #include <TStopwatch.h>
34 #include <TVirtualMC.h>
36 #include "AliCDBEntry.h"
37 #include "AliCDBManager.h"
38 #include "AliCDBStorage.h"
39 #include "AliDetector.h"
40 #include "AliGenerator.h"
41 #include "AliGeomManager.h"
42 #include "AliHeader.h"
49 #include "AliSimulation.h"
51 #include "AliTrackReference.h"
55 //_______________________________________________________________________
79 //_______________________________________________________________________
80 AliMC::AliMC(const char *name, const char *title) :
81 TVirtualMCApplication(name, title),
91 fImedia(new TArrayI(1000)),
93 fHitLists(new TList()),
96 fTrackReferences("AliTrackReference", 100),
97 fTmpTrackReferences("AliTrackReference", 100)
100 // Set transport parameters
103 // Prepare the tracking medium lists
104 for(Int_t i=0;i<1000;i++) (*fImedia)[i]=-99;
107 //_______________________________________________________________________
114 // Delete track references
117 //_______________________________________________________________________
118 void AliMC::ConstructGeometry()
121 // Either load geometry from file or create it through usual
122 // loop on detectors. In the first case the method
123 // AliModule::CreateMaterials() only builds fIdtmed and is postponed
124 // at InitGeometry().
127 if(gAlice->IsRootGeometry()){ //load geometry either from CDB or from file
128 if(gAlice->IsGeomFromCDB()){
129 AliInfo("Loading geometry from CDB default storage");
130 AliCDBPath path("GRP","Geometry","Data");
131 AliCDBEntry *entry=AliCDBManager::Instance()->Get(path.GetPath());
132 if(!entry) AliFatal("Unable to load geometry from CDB!");
134 gGeoManager = (TGeoManager*) entry->GetObject();
135 if (!gGeoManager) AliFatal("TGeoManager object not found in the specified CDB entry!");
138 const char *geomfilename = gAlice->GetGeometryFileName();
139 if(gSystem->ExpandPathName(geomfilename)){
140 AliInfo(Form("Loading geometry from file:\n %40s",geomfilename));
141 TGeoManager::Import(geomfilename);
143 AliInfo(Form("Geometry file %40s not found!\n",geomfilename));
148 // Create modules, materials, geometry
149 if (!gGeoManager) new TGeoManager("ALICE", "ALICE geometry");
151 TIter next(gAlice->Modules());
153 AliDebug(1, "Geometry creation:");
154 while((detector = dynamic_cast<AliModule*>(next()))) {
156 // Initialise detector materials and geometry
157 detector->CreateMaterials();
158 detector->CreateGeometry();
159 AliInfo(Form("%10s R:%.2fs C:%.2fs",
160 detector->GetName(),stw.RealTime(),stw.CpuTime()));
166 //_______________________________________________________________________
167 Bool_t AliMC::MisalignGeometry()
169 // Call misalignment code if AliSimulation object was defined.
171 if(!gAlice->IsRootGeometry()){
172 //Set alignable volumes for the whole geometry
173 SetAllAlignableVolumes();
175 // Misalign geometry via AliSimulation instance
176 if (!AliSimulation::Instance()) return kFALSE;
177 AliGeomManager::SetGeometry(gGeoManager);
178 if(!AliGeomManager::CheckSymNamesLUT("ALL"))
179 AliFatal("Current loaded geometry differs in the definition of symbolic names!");
181 return AliSimulation::Instance()->MisalignGeometry(AliRunLoader::Instance());
184 //_______________________________________________________________________
185 void AliMC::ConstructOpGeometry()
188 // Loop all detector modules and call DefineOpticalProperties() method
191 TIter next(gAlice->Modules());
193 AliInfo("Optical properties definition");
194 while((detector = dynamic_cast<AliModule*>(next()))) {
195 // Initialise detector geometry
196 if(gAlice->IsRootGeometry()) detector->CreateMaterials();
197 // Initialise detector optical properties
198 detector->DefineOpticalProperties();
202 //_______________________________________________________________________
203 void AliMC::InitGeometry()
206 // Initialize detectors
209 AliInfo("Initialisation:");
211 TIter next(gAlice->Modules());
213 while((detector = dynamic_cast<AliModule*>(next()))) {
216 AliInfo(Form("%10s R:%.2fs C:%.2fs",
217 detector->GetName(),stw.RealTime(),stw.CpuTime()));
221 //_______________________________________________________________________
222 void AliMC::SetAllAlignableVolumes()
225 // Add alignable volumes (TGeoPNEntries) looping on all
229 AliInfo(Form("Setting entries for all alignable volumes of active detectors"));
231 TIter next(gAlice->Modules());
232 while((detector = dynamic_cast<AliModule*>(next()))) {
233 detector->AddAlignableVolumes();
237 //_______________________________________________________________________
238 void AliMC::GeneratePrimaries()
241 // Generate primary particles and fill them in the stack.
244 Generator()->Generate();
247 //_______________________________________________________________________
248 void AliMC::SetGenerator(AliGenerator *generator)
251 // Load the event generator
253 if(!fGenerator) fGenerator = generator;
256 //_______________________________________________________________________
257 void AliMC::ResetGenerator(AliGenerator *generator)
260 // Load the event generator
264 AliWarning(Form("Replacing generator %s with %s",
265 fGenerator->GetName(),generator->GetName()));
268 AliWarning(Form("Replacing generator %s with NULL",
269 fGenerator->GetName()));
272 fGenerator = generator;
275 //_______________________________________________________________________
276 void AliMC::FinishRun()
278 // Clean generator information
279 AliDebug(1, "fGenerator->FinishRun()");
280 fGenerator->FinishRun();
282 //Output energy summary tables
283 AliDebug(1, "EnergySummary()");
284 ToAliDebug(1, EnergySummary());
287 //_______________________________________________________________________
288 void AliMC::BeginPrimary()
291 // Called at the beginning of each primary track
296 ResetTrackReferences();
299 //_______________________________________________________________________
300 void AliMC::PreTrack()
302 // Actions before the track's transport
304 TObjArray &dets = *gAlice->Modules();
307 for(Int_t i=0; i<=gAlice->GetNdets(); i++)
308 if((module = dynamic_cast<AliModule*>(dets[i])))
312 //_______________________________________________________________________
313 void AliMC::Stepping()
316 // Called at every step during transport
318 Int_t id = DetFromMate(gMC->CurrentMedium());
322 if ( gMC->IsNewTrack() &&
323 gMC->TrackTime() == 0. &&
325 fRDecayMax > fRDecayMin &&
326 gMC->TrackPid() == fDecayPdg )
328 FixParticleDecaytime();
334 // --- If lego option, do it and leave
335 if (AliSimulation::Instance()->Lego())
336 AliSimulation::Instance()->Lego()->StepManager();
339 //Update energy deposition tables
340 AddEnergyDeposit(gMC->CurrentVolID(copy),gMC->Edep());
342 // write tracke reference for track which is dissapearing - MI
344 if (gMC->IsTrackDisappeared() && !(gMC->IsTrackAlive())) {
345 if (gMC->Etot() > 0.05) AddTrackReference(GetCurrentTrackNumber(),
346 AliTrackReference::kDisappeared);
351 //Call the appropriate stepping routine;
352 AliModule *det = dynamic_cast<AliModule*>(gAlice->Modules()->At(id));
353 if(det && det->StepManagerIsEnabled()) {
359 //_______________________________________________________________________
360 void AliMC::EnergySummary()
363 // Print summary of deposited energy
369 Int_t kn, i, left, j, id;
370 const Float_t kzero=0;
371 Int_t ievent=AliRunLoader::Instance()->GetHeader()->GetEvent()+1;
373 // Energy loss information
375 printf("***************** Energy Loss Information per event (GEV) *****************\n");
376 for(kn=1;kn<fEventEnergy.GetSize();kn++) {
379 fEventEnergy[ndep]=kn;
384 ed2=100*TMath::Sqrt(TMath::Max(ed2-ed*ed,kzero))/ed;
387 fSummEnergy[ndep]=ed;
388 fSum2Energy[ndep]=TMath::Min(static_cast<Float_t>(99.),TMath::Max(ed2,kzero));
393 for(kn=0;kn<(ndep-1)/3+1;kn++) {
395 for(i=0;i<(3<left?3:left);i++) {
397 id=Int_t (fEventEnergy[j]+0.1);
398 printf(" %s %10.3f +- %10.3f%%;",gMC->VolName(id),fSummEnergy[j],fSum2Energy[j]);
403 // Relative energy loss in different detectors
404 printf("******************** Relative Energy Loss per event ********************\n");
405 printf("Total energy loss per event %10.3f GeV\n",edtot);
406 for(kn=0;kn<(ndep-1)/5+1;kn++) {
408 for(i=0;i<(5<left?5:left);i++) {
410 id=Int_t (fEventEnergy[j]+0.1);
411 printf(" %s %10.3f%%;",gMC->VolName(id),100*fSummEnergy[j]/edtot);
415 for(kn=0;kn<75;kn++) printf("*");
419 // Reset the TArray's
420 // fEventEnergy.Set(0);
421 // fSummEnergy.Set(0);
422 // fSum2Energy.Set(0);
425 //_____________________________________________________________________________
426 void AliMC::BeginEvent()
429 // Clean-up previous event
431 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
432 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
433 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
434 AliDebug(1, " BEGINNING EVENT ");
435 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
436 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
437 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
439 AliRunLoader *runloader=AliRunLoader::Instance();
441 /*******************************/
442 /* Clean after eventual */
444 /*******************************/
447 //Set the next event in Run Loader -> Cleans trees (TreeK and all trees in detectors),
448 gAlice->SetEventNrInRun(gAlice->GetEventNrInRun()+1);
449 runloader->SetEventNumber(gAlice->GetEventNrInRun());// sets new files, cleans the previous event stuff, if necessary, etc.,
450 AliDebug(1, Form("EventNr is %d",gAlice->GetEventNrInRun()));
452 fEventEnergy.Reset();
453 // Clean detector information
455 if (runloader->Stack())
456 runloader->Stack()->Reset();//clean stack -> tree is unloaded
458 runloader->MakeStack();//or make a new one
461 if(AliSimulation::Instance()->Lego() == 0x0)
463 AliDebug(1, "fRunLoader->MakeTree(K)");
464 runloader->MakeTree("K");
467 AliDebug(1, "gMC->SetStack(fRunLoader->Stack())");
468 gMC->SetStack(runloader->Stack());//Was in InitMC - but was moved here
469 //because we don't have guarantee that
470 //stack pointer is not going to change from event to event
471 //since it bellobgs to header and is obtained via RunLoader
473 // Reset all Detectors & kinematics & make/reset trees
476 runloader->GetHeader()->Reset(AliCDBManager::Instance()->GetRun(),gAlice->GetEvNumber(),
477 gAlice->GetEventNrInRun());
478 // fRunLoader->WriteKinematics("OVERWRITE"); is there any reason to rewrite here since MakeTree does so
480 if(AliSimulation::Instance()->Lego())
482 AliSimulation::Instance()->Lego()->BeginEvent();
487 AliDebug(1, "ResetHits()");
490 AliDebug(1, "fRunLoader->MakeTree(H)");
491 runloader->MakeTree("H");
495 MakeTmpTrackRefsTree();
496 //create new branches and SetAdresses
497 TIter next(gAlice->Modules());
499 while((detector = (AliModule*)next()))
501 AliDebug(2, Form("%s->MakeBranch(H)",detector->GetName()));
502 detector->MakeBranch("H");
506 //_______________________________________________________________________
507 void AliMC::ResetHits()
510 // Reset all Detectors hits
512 TIter next(gAlice->Modules());
514 while((detector = dynamic_cast<AliModule*>(next()))) {
515 detector->ResetHits();
519 //_______________________________________________________________________
520 void AliMC::ResetDigits()
523 // Reset all Detectors digits
525 TIter next(gAlice->Modules());
527 while((detector = dynamic_cast<AliModule*>(next()))) {
528 detector->ResetDigits();
532 //_______________________________________________________________________
533 void AliMC::ResetSDigits()
536 // Reset all Detectors digits
538 TIter next(gAlice->Modules());
540 while((detector = dynamic_cast<AliModule*>(next()))) {
541 detector->ResetSDigits();
545 //_______________________________________________________________________
546 void AliMC::PostTrack()
548 // Posts tracks for each module
550 TObjArray &dets = *gAlice->Modules();
553 for(Int_t i=0; i<=gAlice->GetNdets(); i++)
554 if((module = dynamic_cast<AliModule*>(dets[i])))
558 //_______________________________________________________________________
559 void AliMC::FinishPrimary()
562 // Called at the end of each primary track
565 AliRunLoader *runloader=AliRunLoader::Instance();
566 // static Int_t count=0;
567 // const Int_t times=10;
568 // This primary is finished, purify stack
569 #if ROOT_VERSION_CODE > 262152
570 if (!(gMC->SecondariesAreOrdered())) {
571 if (runloader->Stack()->ReorderKine()) RemapHits();
574 if (runloader->Stack()->PurifyKine()) RemapHits();
576 TIter next(gAlice->Modules());
578 while((detector = dynamic_cast<AliModule*>(next()))) {
579 detector->FinishPrimary();
580 AliLoader* loader = detector->GetLoader();
583 TTree* treeH = loader->TreeH();
584 if (treeH) treeH->Fill(); //can be Lego run and treeH can not exist
588 // Write out track references if any
589 if (fTmpTreeTR) fTmpTreeTR->Fill();
592 void AliMC::RemapHits()
595 // Remaps the track labels of the hits
596 AliRunLoader *runloader=AliRunLoader::Instance();
597 AliStack* stack = runloader->Stack();
598 TList* hitLists = GetHitLists();
599 TIter next(hitLists);
600 TCollection *hitList;
602 while((hitList = dynamic_cast<TCollection*>(next()))) {
603 TIter nexthit(hitList);
605 while((hit = dynamic_cast<AliHit*>(nexthit()))) {
606 hit->SetTrack(stack->TrackLabel(hit->GetTrack()));
611 // This for detectors which have a special mapping mechanism
612 // for hits, such as TPC and TRD
616 TObjArray* modules = gAlice->Modules();
617 TIter nextmod(modules);
619 while((module = (AliModule*) nextmod())) {
620 AliDetector* det = dynamic_cast<AliDetector*> (module);
621 if (det) det->RemapTrackHitIDs(stack->TrackLabelMap());
624 RemapTrackReferencesIDs(stack->TrackLabelMap());
627 //_______________________________________________________________________
628 void AliMC::FinishEvent()
631 // Called at the end of the event.
634 if(AliSimulation::Instance()->Lego()) AliSimulation::Instance()->Lego()->FinishEvent();
636 TIter next(gAlice->Modules());
638 while((detector = dynamic_cast<AliModule*>(next()))) {
639 detector->FinishEvent();
642 //Update the energy deposit tables
644 for(i=0;i<fEventEnergy.GetSize();i++)
646 fSummEnergy[i]+=fEventEnergy[i];
647 fSum2Energy[i]+=fEventEnergy[i]*fEventEnergy[i];
650 AliRunLoader *runloader=AliRunLoader::Instance();
652 AliHeader* header = runloader->GetHeader();
653 AliStack* stack = runloader->Stack();
654 if ( (header == 0x0) || (stack == 0x0) )
655 {//check if we got header and stack. If not cry and exit aliroot
656 AliFatal("Can not get the stack or header from LOADER");
657 return;//never reached
659 // Update Header information
660 header->SetNprimary(stack->GetNprimary());
661 header->SetNtrack(stack->GetNtrack());
663 // Write out the kinematics
664 if (!AliSimulation::Instance()->Lego()) stack->FinishEvent();
666 // Synchronize the TreeTR with TreeK
667 if (fTmpTreeTR) ReorderAndExpandTreeTR();
669 // Write out the event Header information
670 TTree* treeE = runloader->TreeE();
673 header->SetStack(stack);
678 AliError("Can not get TreeE from RL");
681 if(AliSimulation::Instance()->Lego() == 0x0)
683 runloader->WriteKinematics("OVERWRITE");
684 runloader->WriteTrackRefs("OVERWRITE");
685 runloader->WriteHits("OVERWRITE");
688 AliDebug(1, "<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<");
689 AliDebug(1, "<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<");
690 AliDebug(1, "<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<");
691 AliDebug(1, " FINISHING EVENT ");
692 AliDebug(1, "<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<");
693 AliDebug(1, "<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<");
694 AliDebug(1, "<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<");
697 //_______________________________________________________________________
702 //=================Create Materials and geometry
704 // Set alignable volumes for the whole geometry (with old root)
705 #if ROOT_VERSION_CODE < 331527
706 SetAllAlignableVolumes();
708 //Read the cuts for all materials
710 //Build the special IMEDIA table
713 //Compute cross-sections
716 //Initialise geometry deposition table
717 fEventEnergy.Set(gMC->NofVolumes()+1);
718 fSummEnergy.Set(gMC->NofVolumes()+1);
719 fSum2Energy.Set(gMC->NofVolumes()+1);
721 // Register MC in configuration
722 AliConfig::Instance()->Add(gMC);
726 //_______________________________________________________________________
727 void AliMC::MediaTable()
730 // Built media table to get from the media number to
734 Int_t kz, nz, idt, lz, i, k, ind;
736 TObjArray &dets = *gAlice->Detectors();
738 Int_t ndets=gAlice->GetNdets();
741 for (kz=0;kz<ndets;kz++) {
742 // If detector is defined
743 if((det=dynamic_cast<AliModule*>(dets[kz]))) {
744 TArrayI &idtmed = *(det->GetIdtmed());
745 for(nz=0;nz<100;nz++) {
747 // Find max and min material number
748 if((idt=idtmed[nz])) {
749 det->LoMedium() = det->LoMedium() < idt ? det->LoMedium() : idt;
750 det->HiMedium() = det->HiMedium() > idt ? det->HiMedium() : idt;
753 if(det->LoMedium() > det->HiMedium()) {
757 if(det->HiMedium() > fImedia->GetSize()) {
758 AliError(Form("Increase fImedia from %d to %d",
759 fImedia->GetSize(),det->HiMedium()));
762 // Tag all materials in rage as belonging to detector kz
763 for(lz=det->LoMedium(); lz<= det->HiMedium(); lz++) {
770 // Print summary table
771 AliInfo("Tracking media ranges:");
773 for(i=0;i<(ndets-1)/6+1;i++) {
774 for(k=0;k< (6<ndets-i*6?6:ndets-i*6);k++) {
776 det=dynamic_cast<AliModule*>(dets[ind]);
778 printf(" %6s: %3d -> %3d;",det->GetName(),det->LoMedium(),
781 printf(" %6s: %3d -> %3d;","NULL",0,0);
788 //_______________________________________________________________________
789 void AliMC::ReadTransPar()
792 // Read filename to set the transport parameters
796 const Int_t kncuts=10;
797 const Int_t knflags=11;
798 const Int_t knpars=kncuts+knflags;
799 const char kpars[knpars][7] = {"CUTGAM" ,"CUTELE","CUTNEU","CUTHAD","CUTMUO",
800 "BCUTE","BCUTM","DCUTE","DCUTM","PPCUTM","ANNI",
801 "BREM","COMP","DCAY","DRAY","HADR","LOSS",
802 "MULS","PAIR","PHOT","RAYL"};
808 Int_t i, itmed, iret, ktmed, kz;
811 // See whether the file is there
812 filtmp=gSystem->ExpandPathName(fTransParName.Data());
813 lun=fopen(filtmp,"r");
816 AliWarning(Form("File %s does not exist!",fTransParName.Data()));
821 // Initialise cuts and flags
822 for(i=0;i<kncuts;i++) cut[i]=-99;
823 for(i=0;i<knflags;i++) flag[i]=-99;
825 for(i=0;i<256;i++) line[i]='\0';
826 // Read up to the end of line excluded
827 iret=fscanf(lun,"%[^\n]",line);
833 // Read the end of line
836 if(line[0]=='*') continue;
838 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",
839 detName,&itmed,&cut[0],&cut[1],&cut[2],&cut[3],&cut[4],&cut[5],&cut[6],&cut[7],&cut[8],
840 &cut[9],&flag[0],&flag[1],&flag[2],&flag[3],&flag[4],&flag[5],&flag[6],&flag[7],
841 &flag[8],&flag[9],&flag[10]);
845 AliWarning(Form("Error reading file %s",fTransParName.Data()));
848 // Check that the module exist
849 AliModule *mod = gAlice->GetModule(detName);
851 // Get the array of media numbers
852 TArrayI &idtmed = *mod->GetIdtmed();
853 // Check that the tracking medium code is valid
854 if(0<=itmed && itmed < 100) {
857 AliWarning(Form("Invalid tracking medium code %d for %s",itmed,mod->GetName()));
860 // Set energy thresholds
861 for(kz=0;kz<kncuts;kz++) {
863 AliDebug(2, Form("%-6s set to %10.3E for tracking medium code %4d for %s",
864 kpars[kz],cut[kz],itmed,mod->GetName()));
865 gMC->Gstpar(ktmed,kpars[kz],cut[kz]);
868 // Set transport mechanisms
869 for(kz=0;kz<knflags;kz++) {
871 AliDebug(2, Form("%-6s set to %10d for tracking medium code %4d for %s",
872 kpars[kncuts+kz],flag[kz],itmed,mod->GetName()));
873 gMC->Gstpar(ktmed,kpars[kncuts+kz],Float_t(flag[kz]));
877 AliWarning(Form("Invalid medium code %d",itmed));
881 AliDebug(1, Form("%s not present",detName));
887 //_______________________________________________________________________
888 void AliMC::SetTransPar(const char *filename)
891 // Sets the file name for transport parameters
893 fTransParName = filename;
896 //_______________________________________________________________________
897 void AliMC::AddHit(Int_t id, Int_t track, Int_t *vol, Float_t *hits) const
900 // Add a hit to detector id
902 TObjArray &dets = *gAlice->Modules();
903 if(dets[id]) dynamic_cast<AliModule*>(dets[id])->AddHit(track,vol,hits);
906 //_______________________________________________________________________
907 void AliMC::AddDigit(Int_t id, Int_t *tracks, Int_t *digits) const
910 // Add digit to detector id
912 TObjArray &dets = *gAlice->Modules();
913 if(dets[id]) dynamic_cast<AliModule*>(dets[id])->AddDigit(tracks,digits);
916 //_______________________________________________________________________
917 Int_t AliMC::GetCurrentTrackNumber() const {
919 // Returns current track
921 return AliRunLoader::Instance()->Stack()->GetCurrentTrackNumber();
924 //_______________________________________________________________________
925 void AliMC::DumpPart (Int_t i) const
928 // Dumps particle i in the stack
930 AliRunLoader * runloader = AliRunLoader::Instance();
931 if (runloader->Stack())
932 runloader->Stack()->DumpPart(i);
935 //_______________________________________________________________________
936 void AliMC::DumpPStack () const
939 // Dumps the particle stack
941 AliRunLoader * runloader = AliRunLoader::Instance();
942 if (runloader->Stack())
943 runloader->Stack()->DumpPStack();
946 //_______________________________________________________________________
947 Int_t AliMC::GetNtrack() const {
949 // Returns number of tracks in stack
952 AliRunLoader * runloader = AliRunLoader::Instance();
953 if (runloader->Stack())
954 ntracks = runloader->Stack()->GetNtrack();
958 //_______________________________________________________________________
959 Int_t AliMC::GetPrimary(Int_t track) const
962 // return number of primary that has generated track
964 Int_t nprimary = -999;
965 AliRunLoader * runloader = AliRunLoader::Instance();
966 if (runloader->Stack())
967 nprimary = runloader->Stack()->GetPrimary(track);
971 //_______________________________________________________________________
972 TParticle* AliMC::Particle(Int_t i) const
974 // Returns the i-th particle from the stack taking into account
975 // the remaping done by PurifyKine
976 AliRunLoader * runloader = AliRunLoader::Instance();
978 if (runloader->Stack())
979 return runloader->Stack()->Particle(i);
983 //_______________________________________________________________________
984 const TObjArray* AliMC::Particles() const {
986 // Returns pointer to Particles array
988 AliRunLoader * runloader = AliRunLoader::Instance();
990 if (runloader->Stack())
991 return runloader->Stack()->Particles();
995 //_______________________________________________________________________
996 void AliMC::PushTrack(Int_t done, Int_t parent, Int_t pdg, const Float_t *pmom,
997 const Float_t *vpos, const Float_t *polar, Float_t tof,
998 TMCProcess mech, Int_t &ntr, Float_t weight, Int_t is) const
1000 // Delegate to stack
1002 AliRunLoader * runloader = AliRunLoader::Instance();
1004 if (runloader->Stack())
1005 runloader->Stack()->PushTrack(done, parent, pdg, pmom, vpos, polar, tof,
1006 mech, ntr, weight, is);
1009 //_______________________________________________________________________
1010 void AliMC::PushTrack(Int_t done, Int_t parent, Int_t pdg,
1011 Double_t px, Double_t py, Double_t pz, Double_t e,
1012 Double_t vx, Double_t vy, Double_t vz, Double_t tof,
1013 Double_t polx, Double_t poly, Double_t polz,
1014 TMCProcess mech, Int_t &ntr, Float_t weight, Int_t is) const
1016 // Delegate to stack
1018 AliRunLoader * runloader = AliRunLoader::Instance();
1020 if (runloader->Stack())
1021 runloader->Stack()->PushTrack(done, parent, pdg, px, py, pz, e, vx, vy, vz, tof,
1022 polx, poly, polz, mech, ntr, weight, is);
1025 //_______________________________________________________________________
1026 void AliMC::SetHighWaterMark(Int_t nt) const
1029 // Set high water mark for last track in event
1030 AliRunLoader * runloader = AliRunLoader::Instance();
1032 if (runloader->Stack())
1033 runloader->Stack()->SetHighWaterMark(nt);
1036 //_______________________________________________________________________
1037 void AliMC::KeepTrack(Int_t track) const
1040 // Delegate to stack
1042 AliRunLoader * runloader = AliRunLoader::Instance();
1044 if (runloader->Stack())
1045 runloader->Stack()->KeepTrack(track);
1048 //_______________________________________________________________________
1049 void AliMC::FlagTrack(Int_t track) const
1051 // Delegate to stack
1053 AliRunLoader * runloader = AliRunLoader::Instance();
1055 if (runloader->Stack())
1056 runloader->Stack()->FlagTrack(track);
1059 //_______________________________________________________________________
1060 void AliMC::SetCurrentTrack(Int_t track) const
1063 // Set current track number
1065 AliRunLoader * runloader = AliRunLoader::Instance();
1067 if (runloader->Stack())
1068 runloader->Stack()->SetCurrentTrack(track);
1071 //_______________________________________________________________________
1072 AliTrackReference* AliMC::AddTrackReference(Int_t label, Int_t id)
1075 // add a trackrefernce to the list
1076 Int_t primary = GetPrimary(label);
1077 Particle(primary)->SetBit(kKeepBit);
1079 Int_t nref = fTmpTrackReferences.GetEntriesFast();
1080 return new(fTmpTrackReferences[nref]) AliTrackReference(label, id);
1085 //_______________________________________________________________________
1086 void AliMC::ResetTrackReferences()
1089 // Reset all references
1091 fTmpTrackReferences.Clear();
1094 //_______________________________________________________________________
1095 void AliMC::RemapTrackReferencesIDs(Int_t *map)
1098 // Remapping track reference
1099 // Called at finish primary
1102 Int_t nEntries = fTmpTrackReferences.GetEntries();
1103 for (Int_t i=0; i < nEntries; i++){
1104 AliTrackReference * ref = dynamic_cast<AliTrackReference*>(fTmpTrackReferences.UncheckedAt(i));
1106 Int_t newID = map[ref->GetTrack()];
1107 if (newID>=0) ref->SetTrack(newID);
1109 ref->SetBit(kNotDeleted,kFALSE);
1110 fTmpTrackReferences.RemoveAt(i);
1114 fTmpTrackReferences.Compress();
1117 //_______________________________________________________________________
1118 void AliMC::FixParticleDecaytime()
1121 // Fix the particle decay time according to rmin and rmax for decays
1125 gMC->TrackMomentum(p);
1126 Double_t tmin, tmax;
1129 // Transverse velocity
1130 Double_t vt = p.Pt() / p.E();
1132 if ((b = ((AliMagF*)TGeoGlobalMagField::Instance()->GetField())->SolenoidField()) > 0.) { // [kG]
1136 Double_t rho = p.Pt() / 0.0003 / b; // [cm]
1138 // Revolution frequency
1140 Double_t omega = vt / rho;
1142 // Maximum and minimum decay time
1144 // Check for curlers first
1145 const Double_t kOvRhoSqr2 = 1./(rho*TMath::Sqrt(2.));
1146 if (fRDecayMax * kOvRhoSqr2 > 1.) return;
1150 tmax = TMath::ACos((1.-fRDecayMax*kOvRhoSqr2)*(1.+fRDecayMax*kOvRhoSqr2)) / omega; // [ct]
1151 tmin = TMath::ACos((1.-fRDecayMin*kOvRhoSqr2)*(1.+fRDecayMin*kOvRhoSqr2)) / omega; // [ct]
1153 tmax = fRDecayMax / vt; // [ct]
1154 tmin = fRDecayMin / vt; // [ct]
1158 // Dial t using the two limits
1159 Double_t t = tmin + (tmax - tmin) * gRandom->Rndm(); // [ct]
1162 // Force decay time in transport code
1164 gMC->ForceDecayTime(t / 2.99792458e10);
1167 void AliMC::MakeTmpTrackRefsTree()
1169 // Make the temporary track reference tree
1170 fTmpFileTR = new TFile("TrackRefsTmp.root", "recreate");
1171 fTmpTreeTR = new TTree("TreeTR", "Track References");
1172 TClonesArray* pRef = &fTmpTrackReferences;
1173 fTmpTreeTR->Branch("TrackReferences", "TClonesArray", &pRef, 4000);
1176 //_______________________________________________________________________
1177 void AliMC::ReorderAndExpandTreeTR()
1180 // Reorder and expand the temporary track reference tree in order to match the kinematics tree
1183 AliRunLoader *rl = AliRunLoader::Instance();
1186 AliDebug(1, "fRunLoader->MakeTrackRefsContainer()");
1187 rl->MakeTrackRefsContainer();
1188 TTree * treeTR = rl->TreeTR();
1190 // make branch for central track references
1192 TClonesArray* pRef = &fTrackReferences;
1193 branch = treeTR->Branch("TrackReferences", &pRef);
1194 branch->SetAddress(&pRef);
1197 AliStack* stack = rl->Stack();
1198 Int_t np = stack->GetNprimary();
1199 Int_t nt = fTmpTreeTR->GetEntries();
1201 // Loop over tracks and find the secondaries with the help of the kine tree
1204 for (Int_t ip = np - 1; ip > -1; ip--) {
1205 TParticle *part = stack->Particle(ip);
1206 //printf("Particle %5d %5d %5d %5d %5d \n", ip, part->GetPdgCode(), part->GetFirstMother(), part->GetFirstDaughter(), part->GetLastDaughter());
1208 // Skip primaries that have not been transported
1209 Int_t dau1 = part->GetFirstDaughter();
1211 if (!part->TestBit(kTransportBit)) continue;
1213 fTmpTreeTR->GetEntry(it++);
1214 Int_t nh = fTmpTrackReferences.GetEntries();
1215 // Determine range of secondaries produced by this primary
1217 Int_t inext = ip - 1;
1220 part = stack->Particle(inext);
1221 dau2 = part->GetFirstDaughter();
1222 if (!(part->TestBit(kTransportBit)) || dau2 == -1 || dau2 < np) {
1223 // if (dau2 == -1 || dau2 < np) {
1229 dau2 = stack->GetNtrack() - 1;
1232 } // find upper bound
1234 // printf("Check (1) %5d %5d %5d %5d %5d \n", ip, np, it, dau1, dau2);
1236 // Loop over reference hits and find secondary label
1237 for (Int_t id = dau1; (id <= dau2) && (dau1 > -1); id++) {
1238 for (Int_t ih = 0; ih < nh; ih++) {
1239 AliTrackReference* tr = (AliTrackReference*) fTmpTrackReferences.At(ih);
1240 Int_t label = tr->Label();
1242 if (label == ip) continue;
1243 if (label > dau2 || label < dau1)
1244 AliWarning(Form("Track Reference Label out of range !: %5d %5d %5d \n", label, dau1, dau2));
1247 Int_t nref = fTrackReferences.GetEntriesFast();
1248 new(fTrackReferences[nref]) AliTrackReference(*tr);
1252 fTrackReferences.Clear();
1257 // Now loop again and write the primaries
1259 for (Int_t ip = 0; ip < np; ip++) {
1260 TParticle* part = stack->Particle(ip);
1261 // if ((part->GetFirstDaughter() == -1 && part->GetStatusCode() <= 1) || part->GetFirstDaughter() >= np)
1262 if (part->TestBit(kTransportBit))
1264 // Skip particles that have not been transported
1265 fTmpTreeTR->GetEntry(it--);
1266 Int_t nh = fTmpTrackReferences.GetEntries();
1268 // Loop over reference hits and find primary labels
1269 for (Int_t ih = 0; ih < nh; ih++) {
1270 AliTrackReference* tr = (AliTrackReference*) fTmpTrackReferences.At(ih);
1271 Int_t label = tr->Label();
1273 Int_t nref = fTrackReferences.GetEntriesFast();
1274 new(fTrackReferences[nref]) AliTrackReference(*tr);
1279 fTrackReferences.Clear();
1283 if (ifills != stack->GetNtrack())
1284 AliWarning(Form("Number of entries in TreeTR (%5d) unequal to TreeK (%5d) \n", ifills, stack->GetNtrack()));
1288 fTmpFileTR->Close();
1290 fTmpTrackReferences.Clear();
1291 gSystem->Exec("rm -rf TrackRefsTmp.root");