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