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