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