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