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