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