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