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