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