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