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