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