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