]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliMC.cxx
The GMS source file system changed from kDAQ to kDCS
[u/mrichter/AliRoot.git] / STEER / AliMC.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
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  **************************************************************************/
15
16 /* $Id$ */
17
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
23
24 #include <RVersion.h>
25 #include <TBrowser.h>
26 #include <TClonesArray.h>
27 #include <TGeoManager.h>
28 #include <TStopwatch.h>
29 #include <TSystem.h>
30 #include <TVirtualMC.h>
31  
32 #include "AliLog.h"
33 #include "AliDetector.h"
34 #include "AliGenerator.h"
35 #include "AliHeader.h"
36 #include "AliLego.h"
37 #include "AliMC.h"
38 #include "AliMCQA.h"
39 #include "AliRun.h"
40 #include "AliStack.h"
41 #include "AliMagF.h"
42 #include "AliTrackReference.h"
43 #include "AliSimulation.h"
44
45
46 ClassImp(AliMC)
47
48 //_______________________________________________________________________
49 AliMC::AliMC() :
50   fGenerator(0),
51   fEventEnergy(0),
52   fSummEnergy(0),
53   fSum2Energy(0),
54   fTrRmax(1.e10),
55   fTrZmax(1.e10),
56   fRDecayMax(1.e10),
57   fRDecayMin(-1.),
58   fDecayPdg(0),
59   fImedia(0),
60   fTransParName("\0"),
61   fMCQA(0),
62   fHitLists(0),
63   fTrackReferences(0)
64
65 {
66   //default constructor
67   DecayLimits();
68 }
69
70 //_______________________________________________________________________
71 AliMC::AliMC(const char *name, const char *title) :
72   TVirtualMCApplication(name, title),
73   fGenerator(0),
74   fEventEnergy(0),
75   fSummEnergy(0),
76   fSum2Energy(0),
77   fTrRmax(1.e10),
78   fTrZmax(1.e10),
79   fRDecayMax(1.e10),
80   fRDecayMin(-1.),
81   fDecayPdg(0),
82   fImedia(new TArrayI(1000)),
83   fTransParName("\0"),
84   fMCQA(0),
85   fHitLists(new TList()),
86   fTrackReferences(new TClonesArray("AliTrackReference", 100))
87 {
88   //constructor
89   // Set transport parameters
90   SetTransPar();
91   DecayLimits();
92   // Prepare the tracking medium lists
93   for(Int_t i=0;i<1000;i++) (*fImedia)[i]=-99;
94 }
95
96 //_______________________________________________________________________
97 AliMC::AliMC(const AliMC &mc) :
98   TVirtualMCApplication(mc),
99   fGenerator(0),
100   fEventEnergy(0),
101   fSummEnergy(0),
102   fSum2Energy(0),
103   fTrRmax(1.e10),
104   fTrZmax(1.e10),
105   fRDecayMax(1.e10),
106   fRDecayMin(-1.),
107   fDecayPdg(0),
108   fImedia(0),
109   fTransParName("\0"),
110   fMCQA(0),
111   fHitLists(0),
112   fTrackReferences(0)
113 {
114   //
115   // Copy constructor for AliMC
116   //
117   mc.Copy(*this);
118 }
119
120 //_______________________________________________________________________
121 AliMC::~AliMC()
122 {
123   //destructor
124   delete fGenerator;
125   delete fImedia;
126   delete fMCQA;
127   delete fHitLists;
128   // Delete track references
129   if (fTrackReferences) {
130     fTrackReferences->Delete();
131     delete fTrackReferences;
132     fTrackReferences     = 0;
133   }
134
135 }
136
137 //_______________________________________________________________________
138 void AliMC::Copy(TObject &) const
139 {
140   //dummy Copy function
141   AliFatal("Not implemented!");
142 }
143
144 //_______________________________________________________________________
145 void  AliMC::ConstructGeometry() 
146 {
147   //
148   // Either load geometry from file or create it through usual
149   // loop on detectors. In the first case the method
150   // AliModule::CreateMaterials() only builds fIdtmed and is postponed
151   // at InitGeometry().
152   //
153
154   if(gAlice->IsRootGeometry()){
155     // Load geometry
156     const char *geomfilename = gAlice->GetGeometryFileName();
157     if(gSystem->ExpandPathName(geomfilename)){
158       AliInfo(Form("Loading geometry from file:\n %40s\n\n",geomfilename));
159       TGeoManager::Import(geomfilename);
160     }else{
161       AliInfo(Form("Geometry file %40s not found!\n",geomfilename));
162       return;
163     }
164   }else{
165     // Create modules, materials, geometry
166     TStopwatch stw;
167     TIter next(gAlice->Modules());
168     AliModule *detector;
169     AliDebug(1, "Geometry creation:");
170     while((detector = dynamic_cast<AliModule*>(next()))) {
171       stw.Start();
172       // Initialise detector materials and geometry
173       detector->CreateMaterials();
174       detector->CreateGeometry();
175       AliInfo(Form("%10s R:%.2fs C:%.2fs",
176                    detector->GetName(),stw.RealTime(),stw.CpuTime()));
177     }
178   }
179   
180 }
181
182 //_______________________________________________________________________
183 Bool_t  AliMC::MisalignGeometry() 
184 {
185 // Call misalignment code if AliSimulation object was defined.
186
187    //Set alignable volumes for the whole geometry
188    SetAllAlignableVolumes();
189    // Misalign geometry via AliSimulation instance
190    if (!AliSimulation::GetInstance()) return kFALSE;
191    return AliSimulation::GetInstance()->MisalignGeometry(gAlice->GetRunLoader());
192 }   
193
194 //_______________________________________________________________________
195 void  AliMC::ConstructOpGeometry() 
196 {
197   //
198   // Loop all detector modules and call DefineOpticalProperties() method 
199   //
200
201   TIter next(gAlice->Modules());
202   AliModule *detector;
203   AliInfo("Optical properties definition");
204   while((detector = dynamic_cast<AliModule*>(next()))) {
205     // Initialise detector optical properties
206     detector->DefineOpticalProperties();
207   }  
208 }
209
210 //_______________________________________________________________________
211 void  AliMC::InitGeometry()
212
213   //
214   // Initialize detectors
215   //
216
217   AliInfo("Initialisation:");
218   TStopwatch stw;
219   TIter next(gAlice->Modules());
220   AliModule *detector;
221   while((detector = dynamic_cast<AliModule*>(next()))) {
222     stw.Start();
223     // Initialise detector geometry
224     if(gAlice->IsRootGeometry()) detector->CreateMaterials();
225     detector->Init();
226     AliInfo(Form("%10s R:%.2fs C:%.2fs",
227                  detector->GetName(),stw.RealTime(),stw.CpuTime()));
228   }
229 }
230
231 //_______________________________________________________________________
232 void  AliMC::SetAllAlignableVolumes()
233
234   //
235   // Add alignable volumes (TGeoPNEntries) looping on all
236   // active modules
237   //
238
239   AliInfo(Form("Setting entries for all alignable volumes of active detectors"));
240   AliModule *detector;
241   TIter next(gAlice->Modules());
242   while((detector = dynamic_cast<AliModule*>(next()))) {
243     detector->AddAlignableVolumes();
244   }
245 }
246
247 //_______________________________________________________________________
248 void  AliMC::GeneratePrimaries() 
249
250   //
251   // Generate primary particles and fill them in the stack.
252   //
253
254   Generator()->Generate();
255 }
256
257 //_______________________________________________________________________
258 void AliMC::SetGenerator(AliGenerator *generator)
259 {
260   //
261   // Load the event generator
262   //
263   if(!fGenerator) fGenerator = generator;
264 }
265
266 //_______________________________________________________________________
267 void AliMC::ResetGenerator(AliGenerator *generator)
268 {
269   //
270   // Load the event generator
271   //
272   if(fGenerator) {
273     if(generator) {
274       AliWarning(Form("Replacing generator %s with %s",
275                       fGenerator->GetName(),generator->GetName()));
276     }
277     else {
278       AliWarning(Form("Replacing generator %s with NULL",
279                       fGenerator->GetName()));
280     }
281   }
282   fGenerator = generator;
283 }
284
285 //_______________________________________________________________________
286 void AliMC::FinishRun()
287 {
288   // Clean generator information
289   AliDebug(1, "fGenerator->FinishRun()");
290   fGenerator->FinishRun();
291
292   //Output energy summary tables
293   AliDebug(1, "EnergySummary()");
294   ToAliDebug(1, EnergySummary());
295 }
296
297 //_______________________________________________________________________
298 void AliMC::BeginPrimary()
299 {
300   //
301   // Called  at the beginning of each primary track
302   //
303   
304   // Reset Hits info
305   ResetHits();
306   ResetTrackReferences();
307
308 }
309
310 //_______________________________________________________________________
311 void AliMC::PreTrack()
312 {
313   // Actions before the track's transport
314      TObjArray &dets = *gAlice->Modules();
315      AliModule *module;
316
317      for(Int_t i=0; i<=gAlice->GetNdets(); i++)
318        if((module = dynamic_cast<AliModule*>(dets[i])))
319          module->PreTrack();
320
321      fMCQA->PreTrack();
322 }
323
324 //_______________________________________________________________________
325 void AliMC::Stepping() 
326 {
327   //
328   // Called at every step during transport
329   //
330     
331   Int_t id = DetFromMate(gMC->CurrentMedium());
332   if (id < 0) return;
333
334
335   if ( gMC->IsNewTrack()            && 
336        gMC->TrackTime() == 0.       &&
337        fRDecayMin >= 0.             &&  
338        fRDecayMax > fRDecayMin      &&
339        gMC->TrackPid() == fDecayPdg ) 
340   {
341       FixParticleDecaytime();
342   } 
343     
344
345   
346   //
347   // --- If lego option, do it and leave 
348   if (gAlice->Lego())
349     gAlice->Lego()->StepManager();
350   else {
351     Int_t copy;
352     //Update energy deposition tables
353     AddEnergyDeposit(gMC->CurrentVolID(copy),gMC->Edep());
354     //
355     // write tracke reference for track which is dissapearing - MI
356     if (gMC->IsTrackDisappeared()) {      
357         if (gMC->Etot()>0.05) AddTrackReference(GetCurrentTrackNumber());
358     }
359   
360     //Call the appropriate stepping routine;
361     AliModule *det = dynamic_cast<AliModule*>(gAlice->Modules()->At(id));
362     if(det && det->StepManagerIsEnabled()) {
363       if(AliLog::GetGlobalDebugLevel()>0) fMCQA->StepManager(id);
364       det->StepManager();
365     }
366   }
367 }
368
369 //_______________________________________________________________________
370 void AliMC::EnergySummary()
371 {
372   //
373   // Print summary of deposited energy
374   //
375
376   Int_t ndep=0;
377   Float_t edtot=0;
378   Float_t ed, ed2;
379   Int_t kn, i, left, j, id;
380   const Float_t kzero=0;
381   Int_t ievent=gAlice->GetRunLoader()->GetHeader()->GetEvent()+1;
382   //
383   // Energy loss information
384   if(ievent) {
385     printf("***************** Energy Loss Information per event (GEV) *****************\n");
386     for(kn=1;kn<fEventEnergy.GetSize();kn++) {
387       ed=fSummEnergy[kn];
388       if(ed>0) {
389         fEventEnergy[ndep]=kn;
390         if(ievent>1) {
391           ed=ed/ievent;
392           ed2=fSum2Energy[kn];
393           ed2=ed2/ievent;
394           ed2=100*TMath::Sqrt(TMath::Max(ed2-ed*ed,kzero))/ed;
395         } else 
396           ed2=99;
397         fSummEnergy[ndep]=ed;
398         fSum2Energy[ndep]=TMath::Min(static_cast<Float_t>(99.),TMath::Max(ed2,kzero));
399         edtot+=ed;
400         ndep++;
401       }
402     }
403     for(kn=0;kn<(ndep-1)/3+1;kn++) {
404       left=ndep-kn*3;
405       for(i=0;i<(3<left?3:left);i++) {
406         j=kn*3+i;
407         id=Int_t (fEventEnergy[j]+0.1);
408         printf(" %s %10.3f +- %10.3f%%;",gMC->VolName(id),fSummEnergy[j],fSum2Energy[j]);
409       }
410       printf("\n");
411     }
412     //
413     // Relative energy loss in different detectors
414     printf("******************** Relative Energy Loss per event ********************\n");
415     printf("Total energy loss per event %10.3f GeV\n",edtot);
416     for(kn=0;kn<(ndep-1)/5+1;kn++) {
417       left=ndep-kn*5;
418       for(i=0;i<(5<left?5:left);i++) {
419         j=kn*5+i;
420         id=Int_t (fEventEnergy[j]+0.1);
421         printf(" %s %10.3f%%;",gMC->VolName(id),100*fSummEnergy[j]/edtot);
422       }
423       printf("\n");
424     }
425     for(kn=0;kn<75;kn++) printf("*"); 
426     printf("\n");
427   }
428   //
429   // Reset the TArray's
430   //  fEventEnergy.Set(0);
431   //  fSummEnergy.Set(0);
432   //  fSum2Energy.Set(0);
433 }
434
435 //_____________________________________________________________________________
436 void AliMC::BeginEvent()
437 {
438   //
439   // Clean-up previous event
440   // Energy scores
441   AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
442   AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
443   AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
444   AliDebug(1, "          BEGINNING EVENT               ");
445   AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
446   AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
447   AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
448     
449     AliRunLoader *runloader=gAlice->GetRunLoader();
450
451   /*******************************/    
452   /*   Clean after eventual      */
453   /*   previous event            */
454   /*******************************/    
455
456   
457   //Set the next event in Run Loader -> Cleans trees (TreeK and all trees in detectors),
458   gAlice->SetEventNrInRun(gAlice->GetEventNrInRun()+1);
459   runloader->SetEventNumber(gAlice->GetEventNrInRun());// sets new files, cleans the previous event stuff, if necessary, etc.,  
460   AliDebug(1, Form("EventNr is %d",gAlice->GetEventNrInRun()));
461      
462   fEventEnergy.Reset();  
463     // Clean detector information
464   
465   if (runloader->Stack())
466     runloader->Stack()->Reset();//clean stack -> tree is unloaded
467   else
468     runloader->MakeStack();//or make a new one
469   
470   if(gAlice->Lego() == 0x0)
471    { 
472      AliDebug(1, "fRunLoader->MakeTree(K)");
473      runloader->MakeTree("K");
474    }
475    
476   AliDebug(1, "gMC->SetStack(fRunLoader->Stack())");
477   gMC->SetStack(gAlice->GetRunLoader()->Stack());//Was in InitMC - but was moved here 
478                                      //because we don't have guarantee that 
479                                      //stack pointer is not going to change from event to event
480                          //since it bellobgs to header and is obtained via RunLoader
481   //
482   //  Reset all Detectors & kinematics & make/reset trees
483   //
484     
485   runloader->GetHeader()->Reset(gAlice->GetRunNumber(),gAlice->GetEvNumber(),
486     gAlice->GetEventNrInRun());
487 //  fRunLoader->WriteKinematics("OVERWRITE");  is there any reason to rewrite here since MakeTree does so
488
489   if(gAlice->Lego()) 
490    {
491     gAlice->Lego()->BeginEvent();
492     return;
493    }
494   
495
496   AliDebug(1, "ResetHits()");
497   ResetHits();
498   
499   AliDebug(1, "fRunLoader->MakeTree(H)");
500   runloader->MakeTree("H");
501   
502   AliDebug(1, "fRunLoader->MakeTrackRefsContainer()");
503   runloader->MakeTrackRefsContainer();//for insurance
504
505
506   //create new branches and SetAdresses
507   TIter next(gAlice->Modules());
508   AliModule *detector;
509   while((detector = (AliModule*)next()))
510    {
511     AliDebug(2, Form("%s->MakeBranch(H)",detector->GetName()));
512     detector->MakeBranch("H"); 
513     AliDebug(2, Form("%s->MakeBranchTR()",detector->GetName()));
514     detector->MakeBranchTR();
515     AliDebug(2, Form("%s->SetTreeAddress()",detector->GetName()));
516     detector->SetTreeAddress();
517    }
518   // make branch for AliRun track References
519   TTree * treeTR = runloader->TreeTR();
520   if (treeTR){
521     // make branch for central track references
522     if (!fTrackReferences) fTrackReferences = new TClonesArray("AliTrackReference",0);
523     TBranch *branch;
524     branch = treeTR->Branch("AliRun",&fTrackReferences);
525     branch->SetAddress(&fTrackReferences);
526   }
527   //
528 }
529
530 //_______________________________________________________________________
531 void AliMC::ResetHits()
532 {
533   //
534   //  Reset all Detectors hits
535   //
536   TIter next(gAlice->Modules());
537   AliModule *detector;
538   while((detector = dynamic_cast<AliModule*>(next()))) {
539      detector->ResetHits();
540   }
541 }
542
543 //_______________________________________________________________________
544 void AliMC::PostTrack()
545 {
546   // Posts tracks for each module
547   TObjArray &dets = *gAlice->Modules();
548   AliModule *module;
549   
550   for(Int_t i=0; i<=gAlice->GetNdets(); i++)
551     if((module = dynamic_cast<AliModule*>(dets[i])))
552       module->PostTrack();
553 }
554
555 //_______________________________________________________________________
556 void AliMC::FinishPrimary()
557 {
558   //
559   // Called  at the end of each primary track
560   //
561   AliRunLoader *runloader=gAlice->GetRunLoader();
562   //  static Int_t count=0;
563   //  const Int_t times=10;
564   // This primary is finished, purify stack
565 #if ROOT_VERSION_CODE > 262152
566   if (!(gMC->SecondariesAreOrdered()))
567        runloader->Stack()->ReorderKine();
568 #endif
569   runloader->Stack()->PurifyKine();
570   
571   TIter next(gAlice->Modules());
572   AliModule *detector;
573   while((detector = dynamic_cast<AliModule*>(next()))) {
574     detector->FinishPrimary();
575     AliLoader* loader = detector->GetLoader();
576     if(loader)
577      {
578        TTree* treeH = loader->TreeH();
579        if (treeH) treeH->Fill(); //can be Lego run and treeH can not exist
580      }
581   }
582
583   // Write out track references if any
584   if (runloader->TreeTR()) runloader->TreeTR()->Fill();
585 }
586
587 //_______________________________________________________________________
588 void AliMC::FinishEvent()
589 {
590   //
591   // Called at the end of the event.
592   //
593   
594   //
595     
596   if(gAlice->Lego()) gAlice->Lego()->FinishEvent();
597
598   TIter next(gAlice->Modules());
599   AliModule *detector;
600   while((detector = dynamic_cast<AliModule*>(next()))) {
601     detector->FinishEvent();
602   }
603
604   //Update the energy deposit tables
605   Int_t i;
606   for(i=0;i<fEventEnergy.GetSize();i++) 
607    {
608     fSummEnergy[i]+=fEventEnergy[i];
609     fSum2Energy[i]+=fEventEnergy[i]*fEventEnergy[i];
610    }
611
612   AliRunLoader *runloader=gAlice->GetRunLoader();
613
614   AliHeader* header = runloader->GetHeader();
615   AliStack* stack = runloader->Stack();
616   if ( (header == 0x0) || (stack == 0x0) )
617    {//check if we got header and stack. If not cry and exit aliroot
618     AliFatal("Can not get the stack or header from LOADER");
619     return;//never reached
620    }  
621   // Update Header information 
622   header->SetNprimary(stack->GetNprimary());
623   header->SetNtrack(stack->GetNtrack());  
624
625   
626   // Write out the kinematics
627   if (!gAlice->Lego()) stack->FinishEvent();
628    
629   // Write out the event Header information
630   TTree* treeE = runloader->TreeE();
631   if (treeE) 
632    {
633       header->SetStack(stack);
634       treeE->Fill();
635    }
636   else
637    {
638     AliError("Can not get TreeE from RL");
639    }
640   
641   if(gAlice->Lego() == 0x0)
642    {
643      runloader->WriteKinematics("OVERWRITE");
644      runloader->WriteTrackRefs("OVERWRITE");
645      runloader->WriteHits("OVERWRITE");
646    }
647    
648   AliDebug(1, "<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<");
649   AliDebug(1, "<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<");
650   AliDebug(1, "<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<");
651   AliDebug(1, "          FINISHING EVENT               ");
652   AliDebug(1, "<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<");
653   AliDebug(1, "<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<");
654   AliDebug(1, "<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<");
655 }
656
657 //_______________________________________________________________________
658 void AliMC::Field(const Double_t* x, Double_t* b) const
659 {
660   // Calculates field "b" at point "x"
661     gAlice->Field(x,b);
662 }
663     
664 //_______________________________________________________________________
665 void AliMC::Init()
666 {
667   // MC initialization
668
669    //=================Create Materials and geometry
670    gMC->Init();
671   // Set alignable volumes for the whole geometry (with old root)
672 #if ROOT_VERSION_CODE < 331527
673   SetAllAlignableVolumes();
674 #endif
675    //Read the cuts for all materials
676    ReadTransPar();
677    //Build the special IMEDIA table
678    MediaTable();
679
680    //Compute cross-sections
681    gMC->BuildPhysics();
682    
683    //Initialise geometry deposition table
684    fEventEnergy.Set(gMC->NofVolumes()+1);
685    fSummEnergy.Set(gMC->NofVolumes()+1);
686    fSum2Energy.Set(gMC->NofVolumes()+1);
687
688    //
689    fMCQA = new AliMCQA(gAlice->GetNdets());
690
691    // Register MC in configuration 
692    AliConfig::Instance()->Add(gMC);
693
694 }
695
696 //_______________________________________________________________________
697 void AliMC::MediaTable()
698 {
699   //
700   // Built media table to get from the media number to
701   // the detector id
702   //
703
704   Int_t kz, nz, idt, lz, i, k, ind;
705   //  Int_t ibeg;
706   TObjArray &dets = *gAlice->Detectors();
707   AliModule *det;
708   Int_t ndets=gAlice->GetNdets();
709   //
710   // For all detectors
711   for (kz=0;kz<ndets;kz++) {
712     // If detector is defined
713     if((det=dynamic_cast<AliModule*>(dets[kz]))) {
714         TArrayI &idtmed = *(det->GetIdtmed()); 
715         for(nz=0;nz<100;nz++) {
716             
717         // Find max and min material number
718         if((idt=idtmed[nz])) {
719           det->LoMedium() = det->LoMedium() < idt ? det->LoMedium() : idt;
720           det->HiMedium() = det->HiMedium() > idt ? det->HiMedium() : idt;
721         }
722       }
723       if(det->LoMedium() > det->HiMedium()) {
724         det->LoMedium() = 0;
725         det->HiMedium() = 0;
726       } else {
727         if(det->HiMedium() > fImedia->GetSize()) {
728           AliError(Form("Increase fImedia from %d to %d",
729                         fImedia->GetSize(),det->HiMedium()));
730           return;
731         }
732         // Tag all materials in rage as belonging to detector kz
733         for(lz=det->LoMedium(); lz<= det->HiMedium(); lz++) {
734           (*fImedia)[lz]=kz;
735         }
736       }
737     }
738   }
739   //
740   // Print summary table
741   AliInfo("Tracking media ranges:");
742   ToAliInfo(
743   for(i=0;i<(ndets-1)/6+1;i++) {
744     for(k=0;k< (6<ndets-i*6?6:ndets-i*6);k++) {
745       ind=i*6+k;
746       det=dynamic_cast<AliModule*>(dets[ind]);
747       if(det)
748         printf(" %6s: %3d -> %3d;",det->GetName(),det->LoMedium(),
749                det->HiMedium());
750       else
751         printf(" %6s: %3d -> %3d;","NULL",0,0);
752     }
753     printf("\n");
754   }
755   )
756 }
757
758 //_______________________________________________________________________
759 void AliMC::ReadTransPar()
760 {
761   //
762   // Read filename to set the transport parameters
763   //
764
765
766   const Int_t kncuts=10;
767   const Int_t knflags=11;
768   const Int_t knpars=kncuts+knflags;
769   const char kpars[knpars][7] = {"CUTGAM" ,"CUTELE","CUTNEU","CUTHAD","CUTMUO",
770                                "BCUTE","BCUTM","DCUTE","DCUTM","PPCUTM","ANNI",
771                                "BREM","COMP","DCAY","DRAY","HADR","LOSS",
772                                "MULS","PAIR","PHOT","RAYL"};
773   char line[256];
774   char detName[7];
775   char* filtmp;
776   Float_t cut[kncuts];
777   Int_t flag[knflags];
778   Int_t i, itmed, iret, ktmed, kz;
779   FILE *lun;
780   //
781   // See whether the file is there
782   filtmp=gSystem->ExpandPathName(fTransParName.Data());
783   lun=fopen(filtmp,"r");
784   delete [] filtmp;
785   if(!lun) {
786     AliWarning(Form("File %s does not exist!",fTransParName.Data()));
787     return;
788   }
789   //
790   while(1) {
791     // Initialise cuts and flags
792     for(i=0;i<kncuts;i++) cut[i]=-99;
793     for(i=0;i<knflags;i++) flag[i]=-99;
794     itmed=0;
795     for(i=0;i<256;i++) line[i]='\0';
796     // Read up to the end of line excluded
797     iret=fscanf(lun,"%[^\n]",line);
798     if(iret<0) {
799       //End of file
800       fclose(lun);
801       return;
802     }
803     // Read the end of line
804     fscanf(lun,"%*c");
805     if(!iret) continue;
806     if(line[0]=='*') continue;
807     // Read the numbers
808     iret=sscanf(line,"%s %d %f %f %f %f %f %f %f %f %f %f %d %d %d %d %d %d %d %d %d %d %d",
809                 detName,&itmed,&cut[0],&cut[1],&cut[2],&cut[3],&cut[4],&cut[5],&cut[6],&cut[7],&cut[8],
810                 &cut[9],&flag[0],&flag[1],&flag[2],&flag[3],&flag[4],&flag[5],&flag[6],&flag[7],
811                 &flag[8],&flag[9],&flag[10]);
812     if(!iret) continue;
813     if(iret<0) {
814       //reading error
815       AliWarning(Form("Error reading file %s",fTransParName.Data()));
816       continue;
817     }
818     // Check that the module exist
819     AliModule *mod = gAlice->GetModule(detName);
820     if(mod) {
821       // Get the array of media numbers
822       TArrayI &idtmed = *mod->GetIdtmed();
823       // Check that the tracking medium code is valid
824       if(0<=itmed && itmed < 100) {
825         ktmed=idtmed[itmed];
826         if(!ktmed) {
827           AliWarning(Form("Invalid tracking medium code %d for %s",itmed,mod->GetName()));
828           continue;
829         }
830         // Set energy thresholds
831         for(kz=0;kz<kncuts;kz++) {
832           if(cut[kz]>=0) {
833             AliDebug(2, Form("%-6s set to %10.3E for tracking medium code %4d for %s",
834                              kpars[kz],cut[kz],itmed,mod->GetName()));
835             gMC->Gstpar(ktmed,kpars[kz],cut[kz]);
836           }
837         }
838         // Set transport mechanisms
839         for(kz=0;kz<knflags;kz++) {
840           if(flag[kz]>=0) {
841             AliDebug(2, Form("%-6s set to %10d for tracking medium code %4d for %s",
842                              kpars[kncuts+kz],flag[kz],itmed,mod->GetName()));
843             gMC->Gstpar(ktmed,kpars[kncuts+kz],Float_t(flag[kz]));
844           }
845         }
846       } else {
847         AliWarning(Form("Invalid medium code %d",itmed));
848         continue;
849       }
850     } else {
851       AliDebug(1, Form("%s not present",detName));
852       continue;
853     }
854   }
855 }
856
857 //_______________________________________________________________________
858 void AliMC::SetTransPar(const char *filename)
859 {
860   //
861   // Sets the file name for transport parameters
862   //
863   fTransParName = filename;
864 }
865
866 //_______________________________________________________________________
867 void AliMC::Browse(TBrowser *b)
868 {
869   //
870   // Called when the item "Run" is clicked on the left pane
871   // of the Root browser.
872   // It displays the Root Trees and all detectors.
873   //
874   //detectors are in folders anyway
875   b->Add(fMCQA,"AliMCQA");
876 }
877
878 //_______________________________________________________________________
879 void AliMC::AddHit(Int_t id, Int_t track, Int_t *vol, Float_t *hits) const
880 {
881   //
882   //  Add a hit to detector id
883   //
884   TObjArray &dets = *gAlice->Modules();
885   if(dets[id]) dynamic_cast<AliModule*>(dets[id])->AddHit(track,vol,hits);
886 }
887
888 //_______________________________________________________________________
889 void AliMC::AddDigit(Int_t id, Int_t *tracks, Int_t *digits) const
890 {
891   //
892   // Add digit to detector id
893   //
894   TObjArray &dets = *gAlice->Modules();
895   if(dets[id]) dynamic_cast<AliModule*>(dets[id])->AddDigit(tracks,digits);
896 }
897
898 //_______________________________________________________________________
899 Int_t AliMC::GetCurrentTrackNumber() const {
900   //
901   // Returns current track
902   //
903   return gAlice->GetRunLoader()->Stack()->GetCurrentTrackNumber();
904 }
905
906 //_______________________________________________________________________
907 void AliMC::DumpPart (Int_t i) const
908 {
909   //
910   // Dumps particle i in the stack
911   //
912   AliRunLoader * runloader = gAlice->GetRunLoader();
913    if (runloader->Stack())
914     runloader->Stack()->DumpPart(i);
915 }
916
917 //_______________________________________________________________________
918 void AliMC::DumpPStack () const
919 {
920   //
921   // Dumps the particle stack
922   //
923   AliRunLoader * runloader = gAlice->GetRunLoader();
924    if (runloader->Stack())
925     runloader->Stack()->DumpPStack();
926 }
927
928 //_______________________________________________________________________
929 Int_t AliMC::GetNtrack() const {
930   //
931   // Returns number of tracks in stack
932   //
933   Int_t ntracks = -1;
934   AliRunLoader * runloader = gAlice->GetRunLoader();
935    if (runloader->Stack())
936      ntracks = runloader->Stack()->GetNtrack();
937    return ntracks;
938 }
939
940 //_______________________________________________________________________
941 Int_t AliMC::GetPrimary(Int_t track) const
942 {
943   //
944   // return number of primary that has generated track
945   //
946   Int_t nprimary = -999;
947   AliRunLoader * runloader = gAlice->GetRunLoader();
948   if (runloader->Stack())
949     nprimary = runloader->Stack()->GetPrimary(track);
950   return nprimary;
951 }
952  
953 //_______________________________________________________________________
954 TParticle* AliMC::Particle(Int_t i) const
955 {
956   // Returns the i-th particle from the stack taking into account
957   // the remaping done by PurifyKine
958   AliRunLoader * runloader = gAlice->GetRunLoader();
959   if (runloader)
960    if (runloader->Stack())
961     return runloader->Stack()->Particle(i);
962   return 0x0;   
963 }
964
965 //_______________________________________________________________________
966 TObjArray* AliMC::Particles() const {
967   //
968   // Returns pointer to Particles array
969   //
970   AliRunLoader * runloader = gAlice->GetRunLoader();
971   if (runloader)
972    if (runloader->Stack())
973     return runloader->Stack()->Particles();
974   return 0x0;
975 }
976
977 //_______________________________________________________________________
978 void AliMC::PushTrack(Int_t done, Int_t parent, Int_t pdg, Float_t *pmom,
979                       Float_t *vpos, Float_t *polar, Float_t tof,
980                       TMCProcess mech, Int_t &ntr, Float_t weight, Int_t is) const
981
982 // Delegate to stack
983 //
984   AliRunLoader * runloader = gAlice->GetRunLoader();
985   if (runloader)
986     if (runloader->Stack())
987       runloader->Stack()->PushTrack(done, parent, pdg, pmom, vpos, polar, tof,
988                                     mech, ntr, weight, is);
989 }
990
991 //_______________________________________________________________________
992 void AliMC::PushTrack(Int_t done, Int_t parent, Int_t pdg,
993                       Double_t px, Double_t py, Double_t pz, Double_t e,
994                       Double_t vx, Double_t vy, Double_t vz, Double_t tof,
995                       Double_t polx, Double_t poly, Double_t polz,
996                       TMCProcess mech, Int_t &ntr, Float_t weight, Int_t is) const
997
998   // Delegate to stack
999   //
1000   AliRunLoader * runloader = gAlice->GetRunLoader();
1001   if (runloader)
1002     if (runloader->Stack())
1003       runloader->Stack()->PushTrack(done, parent, pdg, px, py, pz, e, vx, vy, vz, tof,
1004                                     polx, poly, polz, mech, ntr, weight, is);
1005 }
1006
1007 //_______________________________________________________________________
1008 void AliMC::SetHighWaterMark(Int_t nt) const
1009 {
1010     //
1011     // Set high water mark for last track in event
1012   AliRunLoader * runloader = gAlice->GetRunLoader();
1013   if (runloader)
1014     if (runloader->Stack())
1015       runloader->Stack()->SetHighWaterMark(nt);
1016 }
1017
1018 //_______________________________________________________________________
1019 void AliMC::KeepTrack(Int_t track) const
1020
1021   //
1022   // Delegate to stack
1023   //
1024   AliRunLoader * runloader = gAlice->GetRunLoader();
1025   if (runloader)
1026     if (runloader->Stack())
1027       runloader->Stack()->KeepTrack(track);
1028 }
1029  
1030 //_______________________________________________________________________
1031 void AliMC::FlagTrack(Int_t track) const
1032 {
1033   // Delegate to stack
1034   //
1035   AliRunLoader * runloader = gAlice->GetRunLoader();
1036   if (runloader)
1037     if (runloader->Stack())
1038       runloader->Stack()->FlagTrack(track);
1039 }
1040
1041 //_______________________________________________________________________
1042 void AliMC::SetCurrentTrack(Int_t track) const
1043
1044   //
1045   // Set current track number
1046   //
1047   AliRunLoader * runloader = gAlice->GetRunLoader();
1048   if (runloader)
1049     if (runloader->Stack())
1050       runloader->Stack()->SetCurrentTrack(track); 
1051 }
1052
1053 //_______________________________________________________________________
1054 void  AliMC::AddTrackReference(Int_t label){
1055   //
1056   // add a trackrefernce to the list
1057   if (!fTrackReferences) {
1058     AliError("Container trackrefernce not active");
1059     return;
1060   }
1061   Int_t nref = fTrackReferences->GetEntriesFast();
1062   TClonesArray &lref = *fTrackReferences;
1063   new(lref[nref]) AliTrackReference(label);
1064 }
1065
1066
1067
1068 //_______________________________________________________________________
1069 void AliMC::ResetTrackReferences()
1070 {
1071   //
1072   //  Reset all  references
1073   //
1074   if (fTrackReferences)   fTrackReferences->Clear();
1075
1076   TIter next(gAlice->Modules());
1077   AliModule *detector;
1078   while((detector = dynamic_cast<AliModule*>(next()))) {
1079      detector->ResetTrackReferences();
1080   }
1081 }
1082
1083 void AliMC::RemapTrackReferencesIDs(Int_t *map)
1084 {
1085   // 
1086   // Remapping track reference
1087   // Called at finish primary
1088   //
1089   if (!fTrackReferences) return;
1090   Int_t nEntries = fTrackReferences->GetEntries();
1091   for (Int_t i=0; i < nEntries; i++){
1092       AliTrackReference * ref = dynamic_cast<AliTrackReference*>(fTrackReferences->UncheckedAt(i));
1093       if (ref) {
1094           Int_t newID = map[ref->GetTrack()];
1095           if (newID>=0) ref->SetTrack(newID);
1096           else {
1097               ref->SetBit(kNotDeleted,kFALSE);
1098               fTrackReferences->RemoveAt(i);  
1099           }      
1100       } // if ref
1101   }
1102   fTrackReferences->Compress();
1103 }
1104
1105 void AliMC::FixParticleDecaytime()
1106 {
1107     //
1108     // Fix the particle decay time according to rmin and rmax for decays
1109     //
1110
1111     TLorentzVector p;
1112     gMC->TrackMomentum(p);
1113     Double_t tmin, tmax;
1114     Double_t b;
1115
1116     // Transverse velocity 
1117     Double_t vt    = p.Pt() / p.E();
1118     
1119     if ((b = gAlice->Field()->SolenoidField()) > 0.) {     // [kG]
1120
1121         // Radius of helix
1122         
1123         Double_t rho   = p.Pt() / 0.0003 / b; // [cm]
1124         
1125         // Revolution frequency
1126         
1127         Double_t omega = vt / rho;
1128         
1129         // Maximum and minimum decay time
1130         //
1131         // Check for curlers first
1132         if (fRDecayMax * fRDecayMax / rho / rho / 2. > 1.) return;
1133         
1134         //
1135  
1136         tmax  = TMath::ACos(1. - fRDecayMax * fRDecayMax / rho / rho / 2.) / omega;   // [ct]
1137         tmin  = TMath::ACos(1. - fRDecayMin * fRDecayMin / rho / rho / 2.) / omega;   // [ct]
1138     } else {
1139         tmax =  fRDecayMax / vt;                                                      // [ct] 
1140         tmin =  fRDecayMin / vt;                                                      // [ct]
1141     }
1142     
1143     //
1144     // Dial t using the two limits
1145     Double_t t = tmin + (tmax - tmin) * gRandom->Rndm();                              // [ct]
1146     //
1147     //
1148     // Force decay time in transport code
1149     //
1150     gMC->ForceDecayTime(t / 2.99792458e10);
1151 }