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