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