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