]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliRun.cxx
Small correction to the previous.
[u/mrichter/AliRoot.git] / STEER / AliRun.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 /*
17 $Log$
18 Revision 1.69  2001/06/28 16:27:50  morsch
19 AliReco() with user control of event range.
20
21 Revision 1.68  2001/06/11 13:14:40  morsch
22 SetAliGenEventHeader() method added.
23
24 Revision 1.67  2001/06/07 18:24:50  buncic
25 Removed compilation warning in AliConfig initialisation.
26
27 Revision 1.66  2001/05/22 14:32:40  hristov
28 Weird inline removed
29
30 Revision 1.65  2001/05/21 17:22:51  buncic
31 Fixed problem with missing AliConfig while reading galice.root
32
33 Revision 1.64  2001/05/16 14:57:22  alibrary
34 New files for folders and Stack
35
36 Revision 1.62  2001/04/06 11:12:33  morsch
37 Clear fParticles after each event. (Ivana Hrivnacova)
38
39 Revision 1.61  2001/03/30 07:04:10  morsch
40 Call fGenerator->FinishRun() for final print-outs, cross-section and weight calculations.
41
42 Revision 1.60  2001/03/21 18:22:30  hristov
43 fParticleFileMap fix (I.Hrivnacova)
44
45 Revision 1.59  2001/03/12 17:47:03  hristov
46 Changes needed on Sun with CC 5.0
47
48 Revision 1.58  2001/03/09 14:27:26  morsch
49 Fix for multiple events per file: inhibit decrease of size of  fParticleFileMap.
50
51 Revision 1.57  2001/02/23 17:40:23  buncic
52 All trees needed for simulation created in RunMC(). TreeR and its branches
53 are now created in new RunReco() method.
54
55 Revision 1.56  2001/02/14 15:45:20  hristov
56 Algorithmic way of getting entry index in fParticleMap. Protection of fParticleFileMap (I.Hrivnacova)
57
58 Revision 1.55  2001/02/12 15:52:54  buncic
59 Removed OpenBaseFile().
60
61 Revision 1.54  2001/02/07 10:39:05  hristov
62 Remove default value for argument
63
64 Revision 1.53  2001/02/06 11:02:26  hristov
65 New SetTrack interface added, added check for unfilled particles in FinishEvent (I.Hrivnacova)
66
67 Revision 1.52  2001/02/05 16:22:25  buncic
68 Added TreeS to GetEvent().
69
70 Revision 1.51  2001/02/02 15:16:20  morsch
71 SetHighWaterMark method added to mark last particle in event.
72
73 Revision 1.50  2001/01/27 10:32:00  hristov
74 Leave the loop when primaries are filled (I.Hrivnacova)
75
76 Revision 1.49  2001/01/26 19:58:48  hristov
77 Major upgrade of AliRoot code
78
79 Revision 1.48  2001/01/17 10:50:50  hristov
80 Corrections to destructors
81
82 Revision 1.47  2000/12/18 10:44:01  morsch
83 Possibility to set field map by passing pointer to objet of type AliMagF via
84 SetField().
85 Example:
86 gAlice->SetField(new AliMagFCM("Map2", "$(ALICE_ROOT)/data/field01.dat",2,1.,10.));
87
88 Revision 1.46  2000/12/14 19:29:27  fca
89 galice.cuts was not read any more
90
91 Revision 1.45  2000/11/30 07:12:49  alibrary
92 Introducing new Rndm and QA classes
93
94 Revision 1.44  2000/10/26 13:58:59  morsch
95 Add possibility to choose the lego generator (of type AliGeneratorLego or derived) when running
96 RunLego(). Default is the base class AliGeneratorLego.
97
98 Revision 1.43  2000/10/09 09:43:17  fca
99 Special remapping of hits for TPC and TRD. End-of-primary action introduced
100
101 Revision 1.42  2000/10/02 21:28:14  fca
102 Removal of useless dependecies via forward declarations
103
104 Revision 1.41  2000/07/13 16:19:09  fca
105 Mainly coding conventions + some small bug fixes
106
107 Revision 1.40  2000/07/12 08:56:25  fca
108 Coding convention correction and warning removal
109
110 Revision 1.39  2000/07/11 18:24:59  fca
111 Coding convention corrections + few minor bug fixes
112
113 Revision 1.38  2000/06/20 13:05:45  fca
114 Writing down the TREE headers before job starts
115
116 Revision 1.37  2000/06/09 20:05:11  morsch
117 Introduce possibility to chose magnetic field version 3: AliMagFDM + field02.dat
118
119 Revision 1.36  2000/06/08 14:03:58  hristov
120 Only one initializer for a default argument
121
122 Revision 1.35  2000/06/07 10:13:14  hristov
123 Delete only existent objects.
124
125 Revision 1.34  2000/05/18 10:45:38  fca
126 Delete Particle Factory properly
127
128 Revision 1.33  2000/05/16 13:10:40  fca
129 New method IsNewTrack and fix for a problem in Father-Daughter relations
130
131 Revision 1.32  2000/04/27 10:38:21  fca
132 Correct termination of Lego Run and introduce Lego getter in AliRun
133
134 Revision 1.31  2000/04/26 10:17:32  fca
135 Changes in Lego for G4 compatibility
136
137 Revision 1.30  2000/04/18 19:11:40  fca
138 Introduce variable Config.C function signature
139
140 Revision 1.29  2000/04/07 11:12:34  fca
141 G4 compatibility changes
142
143 Revision 1.28  2000/04/05 06:51:06  fca
144 Workaround for an HP compiler problem
145
146 Revision 1.27  2000/03/22 18:08:07  fca
147 Rationalisation of the virtual MC interfaces
148
149 Revision 1.26  2000/03/22 13:42:26  fca
150 SetGenerator does not replace an existing generator, ResetGenerator does
151
152 Revision 1.25  2000/02/23 16:25:22  fca
153 AliVMC and AliGeant3 classes introduced
154 ReadEuclid moved from AliRun to AliModule
155
156 Revision 1.24  2000/01/19 17:17:20  fca
157 Introducing a list of lists of hits -- more hits allowed for detector now
158
159 Revision 1.23  1999/12/03 11:14:31  fca
160 Fixing previous wrong checking
161
162 Revision 1.21  1999/11/25 10:40:08  fca
163 Fixing daughters information also in primary tracks
164
165 Revision 1.20  1999/10/04 18:08:49  fca
166 Adding protection against inconsistent Euclid files
167
168 Revision 1.19  1999/09/29 07:50:40  fca
169 Introduction of the Copyright and cvs Log
170
171 */
172
173 ///////////////////////////////////////////////////////////////////////////////
174 //                                                                           //
175 //  Control class for Alice C++                                              //
176 //  Only one single instance of this class exists.                           //
177 //  The object is created in main program aliroot                            //
178 //  and is pointed by the global gAlice.                                     //
179 //                                                                           //
180 //   -Supports the list of all Alice Detectors (fModules).                 //
181 //   -Supports the list of particles (fParticles).                           //
182 //   -Supports the Trees.                                                    //
183 //   -Supports the geometry.                                                 //
184 //   -Supports the event display.                                            //
185 //Begin_Html
186 /*
187 <img src="picts/AliRunClass.gif">
188 */
189 //End_Html
190 //Begin_Html
191 /*
192 <img src="picts/alirun.gif">
193 */
194 //End_Html
195 //                                                                           //
196 ///////////////////////////////////////////////////////////////////////////////
197
198 #include <stdlib.h>
199 #include <stdio.h>
200 #include <string.h>
201 #include <iostream.h>
202  
203 #include <TFile.h>
204 #include <TRandom.h>
205 #include <TBRIK.h> 
206 #include <TCint.h> 
207 #include <TSystem.h>
208 #include <TObjectTable.h>
209 #include <TTree.h>
210 #include <TGeometry.h>
211 #include <TROOT.h>
212 #include <TBrowser.h>
213 #include <TFolder.h>
214 #include <TNode.h>
215
216 #include "TParticle.h"
217 #include "AliRun.h"
218 #include "AliDisplay.h"
219 #include "AliMC.h"
220 #include "AliLego.h"
221 #include "AliMagFC.h"
222 #include "AliMagFCM.h"
223 #include "AliMagFDM.h"
224 #include "AliHit.h"
225 #include "TRandom3.h"
226 #include "AliMCQA.h"
227 #include "AliGenerator.h"
228 #include "AliLegoGenerator.h"
229 #include "AliConfig.h"
230 #include "AliStack.h"
231 #include "AliGenEventHeader.h"
232 #include "AliHeader.h"
233
234 #include "AliDetector.h"
235
236 AliRun *gAlice;
237
238 ClassImp(AliRun)
239
240 //_____________________________________________________________________________
241 AliRun::AliRun()
242 {
243   //
244   // Default constructor for AliRun
245   //
246   fHeader    = 0;
247   fRun       = 0;
248   fEvent     = 0;
249   fStack     = 0;
250   fModules   = 0;
251   fGenerator = 0;
252   fTreeD     = 0;
253   fTreeH     = 0;
254   fTreeE     = 0;
255   fTreeR     = 0;
256   fTreeS     = 0;
257   fGeometry  = 0;
258   fDisplay   = 0;
259   fField     = 0;
260   fMC        = 0;
261   fNdets     = 0;
262   fImedia    = 0;
263   fTrRmax    = 1.e10;
264   fTrZmax    = 1.e10;
265   fInitDone  = kFALSE;
266   fLego      = 0;
267   fPDGDB     = 0;        //Particle factory object!
268   fHitLists  = 0;
269   fConfigFunction    = "\0";
270   fRandom = 0;
271   fMCQA = 0;
272   fTransParName = "\0";
273   fBaseFileName = ".\0";
274   fDebug        = 0;
275 }
276
277 //_____________________________________________________________________________
278 AliRun::AliRun(const char *name, const char *title)
279   : TNamed(name,title)
280 {
281   //
282   //  Constructor for the main processor.
283   //  Creates the geometry
284   //  Creates the list of Detectors.
285   //  Creates the list of particles.
286   //
287   Int_t i;
288   
289   gAlice     = this;
290   fTreeD     = 0;
291   fTreeH     = 0;
292   fTreeE     = 0;
293   fTreeR     = 0;
294   fTreeS     = 0;
295   fTrRmax    = 1.e10;
296   fTrZmax    = 1.e10;
297   fGenerator = 0;
298   fInitDone  = kFALSE;
299   fLego      = 0;
300   fField     = 0;
301   fConfigFunction    = "Config();";
302
303   // Set random number generator
304   gRandom = fRandom = new TRandom3();
305
306   if (gSystem->Getenv("CONFIG_SEED")) {
307      gRandom->SetSeed((UInt_t)atoi(gSystem->Getenv("CONFIG_SEED")));
308   }
309   
310   gROOT->GetListOfBrowsables()->Add(this,name);
311   //
312   // Particle stack
313   fStack = new AliStack(10000);
314   // create the support list for the various Detectors
315   fModules = new TObjArray(77);
316   //
317   // Create the TNode geometry for the event display
318   
319   BuildSimpleGeometry();
320   
321   fHeader    = new AliHeader();
322   fRun       = 0;
323   fEvent     = 0;
324   //
325   fDisplay = 0;
326   //
327   // Create default mag field
328   SetField();
329   //
330   fMC      = gMC;
331   //
332   // Prepare the tracking medium lists
333   fImedia = new TArrayI(1000);
334   for(i=0;i<1000;i++) (*fImedia)[i]=-99;
335   //
336   // Make particles
337   fPDGDB     = TDatabasePDG::Instance();        //Particle factory object!
338
339   AliConfig::Instance()->Add(fPDGDB); 
340   //
341   // Create HitLists list
342   fHitLists  = new TList();
343   //
344   SetTransPar();
345   fBaseFileName = ".\0";
346   //
347   fDebug        = 0;
348 }
349
350
351 //_____________________________________________________________________________
352 AliRun::~AliRun()
353 {
354   //
355   // Default AliRun destructor
356   //
357   delete fImedia;
358   delete fField;
359   delete fMC;
360   delete fGeometry;
361   delete fDisplay;
362   delete fGenerator;
363   delete fLego;
364   delete fTreeD;
365   delete fTreeH;
366   delete fTreeE;
367   delete fTreeR;
368   delete fTreeS;
369   if (fModules) {
370     fModules->Delete();
371     delete fModules;
372   }
373   delete fStack;
374   delete fHitLists;
375   delete fPDGDB;
376   delete fMCQA;
377   delete fHeader;
378 }
379
380 //_____________________________________________________________________________
381 void AliRun::AddHit(Int_t id, Int_t track, Int_t *vol, Float_t *hits) const
382 {
383   //
384   //  Add a hit to detector id
385   //
386   TObjArray &dets = *fModules;
387   if(dets[id]) ((AliModule*) dets[id])->AddHit(track,vol,hits);
388 }
389
390 //_____________________________________________________________________________
391 void AliRun::AddDigit(Int_t id, Int_t *tracks, Int_t *digits) const
392 {
393   //
394   // Add digit to detector id
395   //
396   TObjArray &dets = *fModules;
397   if(dets[id]) ((AliModule*) dets[id])->AddDigit(tracks,digits);
398 }
399
400 //_____________________________________________________________________________
401 void AliRun::Browse(TBrowser *b)
402 {
403   //
404   // Called when the item "Run" is clicked on the left pane
405   // of the Root browser.
406   // It displays the Root Trees and all detectors.
407   //
408   if(!fStack) fStack=fHeader->Stack();
409   TTree*  pTreeK = fStack->TreeK();
410     
411   if (pTreeK) b->Add(pTreeK,pTreeK->GetName());
412   if (fTreeH) b->Add(fTreeH,fTreeH->GetName());
413   if (fTreeD) b->Add(fTreeD,fTreeD->GetName());
414   if (fTreeE) b->Add(fTreeE,fTreeE->GetName());
415   if (fTreeR) b->Add(fTreeR,fTreeR->GetName());
416   if (fTreeS) b->Add(fTreeS,fTreeS->GetName());
417   
418   TIter next(fModules);
419   AliModule *detector;
420   while((detector = (AliModule*)next())) {
421     b->Add(detector,detector->GetName());
422   }
423   b->Add(fMCQA,"AliMCQA");
424 }
425
426 //_____________________________________________________________________________
427 void AliRun::Build()
428 {
429   //
430   // Initialize Alice geometry
431   // Dummy routine
432   //
433 }
434  
435 //_____________________________________________________________________________
436 void AliRun::BuildSimpleGeometry()
437 {
438   //
439   // Create a simple TNode geometry used by Root display engine
440   //
441   // Initialise geometry
442   //
443   fGeometry = new TGeometry("AliceGeom","Galice Geometry for Hits");
444   new TMaterial("void","Vacuum",0,0,0);  //Everything is void
445   TBRIK *brik = new TBRIK("S_alice","alice volume","void",2000,2000,3000);
446   brik->SetVisibility(0);
447   new TNode("alice","alice","S_alice");
448 }
449
450 //_____________________________________________________________________________
451 void AliRun::CleanDetectors()
452 {
453   //
454   // Clean Detectors at the end of event
455   //
456   TIter next(fModules);
457   AliModule *detector;
458   while((detector = (AliModule*)next())) {
459     detector->FinishEvent();
460   }
461 }
462
463 //_____________________________________________________________________________
464 Int_t AliRun::DistancetoPrimitive(Int_t, Int_t)
465 {
466   //
467   // Return the distance from the mouse to the AliRun object
468   // Dummy routine
469   //
470   return 9999;
471 }
472
473 //_____________________________________________________________________________
474 void AliRun::DumpPart (Int_t i) const
475 {
476   //
477   // Dumps particle i in the stack
478   //
479     fStack->DumpPart(i);
480 }
481
482 //_____________________________________________________________________________
483 void AliRun::DumpPStack () const
484 {
485   //
486   // Dumps the particle stack
487   //
488     fStack->DumpPStack();
489 }
490
491 //_____________________________________________________________________________
492 void  AliRun::SetField(AliMagF* magField)
493 {
494     // Set Magnetic Field Map
495     fField = magField;
496     fField->ReadField();
497 }
498
499 //_____________________________________________________________________________
500 void AliRun::SetField(Int_t type, Int_t version, Float_t scale,
501                       Float_t maxField, char* filename)
502 {
503   //
504   //  Set magnetic field parameters
505   //  type      Magnetic field transport flag 0=no field, 2=helix, 3=Runge Kutta
506   //  version   Magnetic field map version (only 1 active now)
507   //  scale     Scale factor for the magnetic field
508   //  maxField  Maximum value for the magnetic field
509
510   //
511   // --- Sanity check on mag field flags
512   if(fField) delete fField;
513   if(version==1) {
514     fField = new AliMagFC("Map1"," ",type,scale,maxField);
515   } else if(version<=2) {
516     fField = new AliMagFCM("Map2-3",filename,type,scale,maxField);
517     fField->ReadField();
518   } else if(version==3) {
519     fField = new AliMagFDM("Map4",filename,type,scale,maxField);
520     fField->ReadField();
521   } else {
522     Warning("SetField","Invalid map %d\n",version);
523   }
524 }
525  
526 //_____________________________________________________________________________
527 void AliRun::PreTrack()
528 {
529      TObjArray &dets = *fModules;
530      AliModule *module;
531
532      for(Int_t i=0; i<=fNdets; i++)
533        if((module = (AliModule*)dets[i]))
534          module->PreTrack();
535
536      fMCQA->PreTrack();
537 }
538
539 //_____________________________________________________________________________
540 void AliRun::PostTrack()
541 {
542      TObjArray &dets = *fModules;
543      AliModule *module;
544
545      for(Int_t i=0; i<=fNdets; i++)
546        if((module = (AliModule*)dets[i]))
547          module->PostTrack();
548 }
549
550 //_____________________________________________________________________________
551 void AliRun::FinishPrimary()
552 {
553   //
554   // Called  at the end of each primary track
555   //
556   
557   //  static Int_t count=0;
558   //  const Int_t times=10;
559   // This primary is finished, purify stack
560   fStack->PurifyKine();
561
562   TIter next(fModules);
563   AliModule *detector;
564   while((detector = (AliModule*)next())) {
565     detector->FinishPrimary();
566   }
567
568   // Write out hits if any
569   if (gAlice->TreeH()) {
570     gAlice->TreeH()->Fill();
571   }
572   
573   //
574   //  if(++count%times==1) gObjectTable->Print();
575 }
576
577 //_____________________________________________________________________________
578 void AliRun::BeginPrimary()
579 {
580   //
581   // Called  at the beginning of each primary track
582   //
583   
584   // Reset Hits info
585   gAlice->ResetHits();
586
587 }
588
589 //_____________________________________________________________________________
590 void AliRun::FinishEvent()
591 {
592   //
593   // Called at the end of the event.
594   //
595   
596   //
597   if(fLego) fLego->FinishEvent();
598
599   //Update the energy deposit tables
600   Int_t i;
601   for(i=0;i<fEventEnergy.GetSize();i++) {
602     fSummEnergy[i]+=fEventEnergy[i];
603     fSum2Energy[i]+=fEventEnergy[i]*fEventEnergy[i];
604   }
605
606
607   
608   // Update Header information 
609
610   fHeader->SetNprimary(fStack->GetNprimary());
611   fHeader->SetNtrack(fStack->GetNtrack());  
612
613   
614   // Write out the kinematics
615   fStack->FinishEvent();
616    
617   // Write out the event Header information
618   if (fTreeE) {
619       fHeader->SetStack(fStack);
620       fTreeE->Fill();
621   }
622   
623   
624   // Write Tree headers
625   TTree* pTreeK = fStack->TreeK();
626   if (pTreeK) pTreeK->Write(0,TObject::kOverwrite);
627   if (fTreeH) fTreeH->Write(0,TObject::kOverwrite);
628   
629   ++fEvent;
630 }
631
632 //_____________________________________________________________________________
633 void AliRun::FinishRun()
634 {
635   //
636   // Called at the end of the run.
637   //
638
639   //
640   if(fLego) fLego->FinishRun();
641   
642   // Clean detector information
643   TIter next(fModules);
644   AliModule *detector;
645   while((detector = (AliModule*)next())) {
646     detector->FinishRun();
647   }
648   
649   //Output energy summary tables
650   EnergySummary();
651
652   TFile *file = fTreeE->GetCurrentFile();
653
654   file->cd();
655   
656   fTreeE->Write(0,TObject::kOverwrite);
657   
658   // Write AliRun info and all detectors parameters
659   Write(0,TObject::kOverwrite);
660
661   // Clean tree information
662
663   fStack->FinishRun();
664   
665   if (fTreeH) {
666     delete fTreeH; fTreeH = 0;
667   }
668   if (fTreeD) {
669     delete fTreeD; fTreeD = 0;
670   }
671   if (fTreeR) {
672     delete fTreeR; fTreeR = 0;
673   }
674   if (fTreeE) {
675     delete fTreeE; fTreeE = 0;
676   }
677   if (fTreeS) {
678     delete fTreeS; fTreeS = 0;
679   }
680     
681   // Close output file
682   file->Write();
683 }
684
685 //_____________________________________________________________________________
686 void AliRun::FlagTrack(Int_t track)
687 {
688   // Delegate to stack
689   //
690     fStack->FlagTrack(track);
691 }
692  
693 //_____________________________________________________________________________
694 void AliRun::EnergySummary()
695 {
696   //
697   // Print summary of deposited energy
698   //
699
700   Int_t ndep=0;
701   Float_t edtot=0;
702   Float_t ed, ed2;
703   Int_t kn, i, left, j, id;
704   const Float_t kzero=0;
705   Int_t ievent=fHeader->GetEvent()+1;
706   //
707   // Energy loss information
708   if(ievent) {
709     printf("***************** Energy Loss Information per event (GEV) *****************\n");
710     for(kn=1;kn<fEventEnergy.GetSize();kn++) {
711       ed=fSummEnergy[kn];
712       if(ed>0) {
713         fEventEnergy[ndep]=kn;
714         if(ievent>1) {
715           ed=ed/ievent;
716           ed2=fSum2Energy[kn];
717           ed2=ed2/ievent;
718           ed2=100*TMath::Sqrt(TMath::Max(ed2-ed*ed,kzero))/ed;
719         } else 
720           ed2=99;
721         fSummEnergy[ndep]=ed;
722         fSum2Energy[ndep]=TMath::Min((Float_t) 99.,TMath::Max(ed2,kzero));
723         edtot+=ed;
724         ndep++;
725       }
726     }
727     for(kn=0;kn<(ndep-1)/3+1;kn++) {
728       left=ndep-kn*3;
729       for(i=0;i<(3<left?3:left);i++) {
730         j=kn*3+i;
731         id=Int_t (fEventEnergy[j]+0.1);
732         printf(" %s %10.3f +- %10.3f%%;",gMC->VolName(id),fSummEnergy[j],fSum2Energy[j]);
733       }
734       printf("\n");
735     }
736     //
737     // Relative energy loss in different detectors
738     printf("******************** Relative Energy Loss per event ********************\n");
739     printf("Total energy loss per event %10.3f GeV\n",edtot);
740     for(kn=0;kn<(ndep-1)/5+1;kn++) {
741       left=ndep-kn*5;
742       for(i=0;i<(5<left?5:left);i++) {
743         j=kn*5+i;
744         id=Int_t (fEventEnergy[j]+0.1);
745         printf(" %s %10.3f%%;",gMC->VolName(id),100*fSummEnergy[j]/edtot);
746       }
747       printf("\n");
748     }
749     for(kn=0;kn<75;kn++) printf("*"); 
750     printf("\n");
751   }
752   //
753   // Reset the TArray's
754   //  fEventEnergy.Set(0);
755   //  fSummEnergy.Set(0);
756   //  fSum2Energy.Set(0);
757 }
758
759 //_____________________________________________________________________________
760 AliModule *AliRun::GetModule(const char *name) const
761 {
762   //
763   // Return pointer to detector from name
764   //
765   return (AliModule*)fModules->FindObject(name);
766 }
767  
768 //_____________________________________________________________________________
769 AliDetector *AliRun::GetDetector(const char *name) const
770 {
771   //
772   // Return pointer to detector from name
773   //
774   return (AliDetector*)fModules->FindObject(name);
775 }
776  
777 //_____________________________________________________________________________
778 Int_t AliRun::GetModuleID(const char *name) const
779 {
780   //
781   // Return galice internal detector identifier from name
782   //
783   Int_t i=-1;
784   TObject *mod=fModules->FindObject(name);
785   if(mod) i=fModules->IndexOf(mod);
786   return i;
787 }
788  
789 //_____________________________________________________________________________
790 Int_t AliRun::GetEvent(Int_t event)
791 {
792   //
793   // Connect the Trees Kinematics and Hits for event # event
794   // Set branch addresses
795   //
796
797   // Reset existing structures
798   ResetHits();
799   ResetDigits();
800   ResetSDigits();
801   
802   // Delete Trees already connected
803   if (fTreeH) delete fTreeH;
804   if (fTreeD) delete fTreeD;
805   if (fTreeR) delete fTreeR;
806   if (fTreeS) delete fTreeS;
807
808  // Create the particle stack
809   if (fHeader) delete fHeader; 
810   fHeader = 0;
811   
812   // Get header from file
813   if(fTreeE) {
814       fTreeE->SetBranchAddress("Header", &fHeader);
815       fTreeE->GetEntry(event);
816   }  
817   else 
818       Error("GetEvent","Cannot find Header Tree (TE)\n");
819
820   // Get the stack from the header
821   if (fStack) delete fStack;
822   fStack = fHeader->Stack();
823   fStack->GetEvent(event);
824   //
825   TFile *file = fTreeE->GetCurrentFile();
826   char treeName[20];
827   file->cd();
828
829   // Get Hits Tree header from file
830   sprintf(treeName,"TreeH%d",event);
831   fTreeH = (TTree*)gDirectory->Get(treeName);
832   if (!fTreeH) {
833     Error("GetEvent","cannot find Hits Tree for event:%d\n",event);
834   }
835   
836   // Get Digits Tree header from file
837   sprintf(treeName,"TreeD%d",event);
838   fTreeD = (TTree*)gDirectory->Get(treeName);
839   if (!fTreeD) {
840     // Warning("GetEvent","cannot find Digits Tree for event:%d\n",event);
841   }
842
843   file->cd();
844
845   // Get SDigits Tree header from file
846   sprintf(treeName,"TreeS%d",event);
847   fTreeS = (TTree*)gDirectory->Get(treeName);
848   if (!fTreeS) {
849     // Warning("GetEvent","cannot find SDigits Tree for event:%d\n",event);
850   }
851
852   file->cd();
853   
854   // Get Reconstruct Tree header from file
855   sprintf(treeName,"TreeR%d",event);
856   fTreeR = (TTree*)gDirectory->Get(treeName);
857   if (!fTreeR) {
858     //    printf("WARNING: cannot find Reconstructed Tree for event:%d\n",event);
859   }
860
861   file->cd();
862  
863   // Set Trees branch addresses
864   TIter next(fModules);
865   AliModule *detector;
866   while((detector = (AliModule*)next())) {
867     detector->SetTreeAddress();
868   }
869
870   return fStack->GetNtrack();
871 }
872
873 //_____________________________________________________________________________
874 TGeometry *AliRun::GetGeometry()
875 {
876   //
877   // Import Alice geometry from current file
878   // Return pointer to geometry object
879   //
880   if (!fGeometry) fGeometry = (TGeometry*)gDirectory->Get("AliceGeom");
881   //
882   // Unlink and relink nodes in detectors
883   // This is bad and there must be a better way...
884   //
885   
886   TIter next(fModules);
887   AliModule *detector;
888   while((detector = (AliModule*)next())) {
889     TList *dnodes=detector->Nodes();
890     Int_t j;
891     TNode *node, *node1;
892     for ( j=0; j<dnodes->GetSize(); j++) {
893       node = (TNode*) dnodes->At(j);
894       node1 = fGeometry->GetNode(node->GetName());
895       dnodes->Remove(node);
896       dnodes->AddAt(node1,j);
897     }
898   }
899   return fGeometry;
900 }
901
902 //_____________________________________________________________________________
903 void AliRun::GetNextTrack(Int_t &mtrack, Int_t &ipart, Float_t *pmom,
904                           Float_t &e, Float_t *vpos, Float_t *polar,
905                           Float_t &tof)
906 {
907   // Delegate to stack
908   //
909     fStack->GetNextTrack(mtrack, ipart, pmom, e, vpos, polar, tof);  
910 }
911
912 //_____________________________________________________________________________
913 Int_t AliRun::GetPrimary(Int_t track) const
914 {
915   //
916   // return number of primary that has generated track
917   //
918     return fStack->GetPrimary(track);
919 }
920  
921 //_____________________________________________________________________________
922 void AliRun::InitMC(const char *setup)
923 {
924   //
925   // Initialize the Alice setup
926   //
927
928   if(fInitDone) {
929     Warning("Init","Cannot initialise AliRun twice!\n");
930     return;
931   }
932     
933   gROOT->LoadMacro(setup);
934   gInterpreter->ProcessLine(fConfigFunction.Data());
935
936
937   gMC->DefineParticles();  //Create standard MC particles
938
939   TObject *objfirst, *objlast;
940
941   fNdets = fModules->GetLast()+1;
942
943   //
944   //=================Create Materials and geometry
945   gMC->Init();
946
947   // Added also after in case of interactive initialisation of modules
948   fNdets = fModules->GetLast()+1;
949
950    TIter next(fModules);
951    AliModule *detector;
952    while((detector = (AliModule*)next())) {
953       detector->SetTreeAddress();
954       objlast = gDirectory->GetList()->Last();
955       
956       // Add Detector histograms in Detector list of histograms
957       if (objlast) objfirst = gDirectory->GetList()->After(objlast);
958       else         objfirst = gDirectory->GetList()->First();
959       while (objfirst) {
960         detector->Histograms()->Add(objfirst);
961         objfirst = gDirectory->GetList()->After(objfirst);
962       }
963    }
964    ReadTransPar(); //Read the cuts for all materials
965    
966    MediaTable(); //Build the special IMEDIA table
967    
968    //Initialise geometry deposition table
969    fEventEnergy.Set(gMC->NofVolumes()+1);
970    fSummEnergy.Set(gMC->NofVolumes()+1);
971    fSum2Energy.Set(gMC->NofVolumes()+1);
972    
973    //Compute cross-sections
974    gMC->BuildPhysics();
975    
976    //Write Geometry object to current file.
977    fGeometry->Write();
978    
979    fInitDone = kTRUE;
980
981    fMCQA = new AliMCQA(fNdets);
982
983    AliConfig::Instance();
984    //
985    // Save stuff at the beginning of the file to avoid file corruption
986    Write();
987 }
988
989 //_____________________________________________________________________________
990 void AliRun::MediaTable()
991 {
992   //
993   // Built media table to get from the media number to
994   // the detector id
995   //
996   Int_t kz, nz, idt, lz, i, k, ind;
997   //  Int_t ibeg;
998   TObjArray &dets = *gAlice->Detectors();
999   AliModule *det;
1000   //
1001   // For all detectors
1002   for (kz=0;kz<fNdets;kz++) {
1003     // If detector is defined
1004     if((det=(AliModule*) dets[kz])) {
1005         TArrayI &idtmed = *(det->GetIdtmed()); 
1006         for(nz=0;nz<100;nz++) {
1007         // Find max and min material number
1008         if((idt=idtmed[nz])) {
1009           det->LoMedium() = det->LoMedium() < idt ? det->LoMedium() : idt;
1010           det->HiMedium() = det->HiMedium() > idt ? det->HiMedium() : idt;
1011         }
1012       }
1013       if(det->LoMedium() > det->HiMedium()) {
1014         det->LoMedium() = 0;
1015         det->HiMedium() = 0;
1016       } else {
1017         if(det->HiMedium() > fImedia->GetSize()) {
1018           Error("MediaTable","Increase fImedia from %d to %d",
1019                 fImedia->GetSize(),det->HiMedium());
1020           return;
1021         }
1022         // Tag all materials in rage as belonging to detector kz
1023         for(lz=det->LoMedium(); lz<= det->HiMedium(); lz++) {
1024           (*fImedia)[lz]=kz;
1025         }
1026       }
1027     }
1028   }
1029   //
1030   // Print summary table
1031   printf(" Traking media ranges:\n");
1032   for(i=0;i<(fNdets-1)/6+1;i++) {
1033     for(k=0;k< (6<fNdets-i*6?6:fNdets-i*6);k++) {
1034       ind=i*6+k;
1035       det=(AliModule*)dets[ind];
1036       if(det)
1037         printf(" %6s: %3d -> %3d;",det->GetName(),det->LoMedium(),
1038                det->HiMedium());
1039       else
1040         printf(" %6s: %3d -> %3d;","NULL",0,0);
1041     }
1042     printf("\n");
1043   }
1044 }
1045
1046 //____________________________________________________________________________
1047 void AliRun::SetGenerator(AliGenerator *generator)
1048 {
1049   //
1050   // Load the event generator
1051   //
1052   if(!fGenerator) fGenerator = generator;
1053 }
1054
1055 //____________________________________________________________________________
1056 void AliRun::ResetGenerator(AliGenerator *generator)
1057 {
1058   //
1059   // Load the event generator
1060   //
1061   if(fGenerator)
1062     if(generator)
1063       Warning("ResetGenerator","Replacing generator %s with %s\n",
1064               fGenerator->GetName(),generator->GetName());
1065     else
1066       Warning("ResetGenerator","Replacing generator %s with NULL\n",
1067               fGenerator->GetName());
1068   fGenerator = generator;
1069 }
1070
1071 //____________________________________________________________________________
1072 void AliRun::SetTransPar(char *filename)
1073 {
1074   fTransParName = filename;
1075 }
1076
1077 //____________________________________________________________________________
1078 void AliRun::SetBaseFile(char *filename)
1079 {
1080   fBaseFileName = filename;
1081 }
1082
1083 //____________________________________________________________________________
1084 void AliRun::ReadTransPar()
1085 {
1086   //
1087   // Read filename to set the transport parameters
1088   //
1089
1090
1091   const Int_t kncuts=10;
1092   const Int_t knflags=11;
1093   const Int_t knpars=kncuts+knflags;
1094   const char kpars[knpars][7] = {"CUTGAM" ,"CUTELE","CUTNEU","CUTHAD","CUTMUO",
1095                                "BCUTE","BCUTM","DCUTE","DCUTM","PPCUTM","ANNI",
1096                                "BREM","COMP","DCAY","DRAY","HADR","LOSS",
1097                                "MULS","PAIR","PHOT","RAYL"};
1098   char line[256];
1099   char detName[7];
1100   char* filtmp;
1101   Float_t cut[kncuts];
1102   Int_t flag[knflags];
1103   Int_t i, itmed, iret, ktmed, kz;
1104   FILE *lun;
1105   //
1106   // See whether the file is there
1107   filtmp=gSystem->ExpandPathName(fTransParName.Data());
1108   lun=fopen(filtmp,"r");
1109   delete [] filtmp;
1110   if(!lun) {
1111     Warning("ReadTransPar","File %s does not exist!\n",fTransParName.Data());
1112     return;
1113   }
1114   //
1115   if(fDebug) {
1116     printf(" "); for(i=0;i<60;i++) printf("*"); printf("\n");
1117     printf(" *%59s\n","*");
1118     printf(" *       Please check carefully what you are doing!%10s\n","*");
1119     printf(" *%59s\n","*");
1120   }
1121   //
1122   while(1) {
1123     // Initialise cuts and flags
1124     for(i=0;i<kncuts;i++) cut[i]=-99;
1125     for(i=0;i<knflags;i++) flag[i]=-99;
1126     itmed=0;
1127     for(i=0;i<256;i++) line[i]='\0';
1128     // Read up to the end of line excluded
1129     iret=fscanf(lun,"%[^\n]",line);
1130     if(iret<0) {
1131       //End of file
1132       fclose(lun);
1133       if(fDebug){
1134         printf(" *%59s\n","*");
1135         printf(" "); for(i=0;i<60;i++) printf("*"); printf("\n");
1136       }
1137       return;
1138     }
1139     // Read the end of line
1140     fscanf(lun,"%*c");
1141     if(!iret) continue;
1142     if(line[0]=='*') continue;
1143     // Read the numbers
1144     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",
1145                 detName,&itmed,&cut[0],&cut[1],&cut[2],&cut[3],&cut[4],&cut[5],&cut[6],&cut[7],&cut[8],
1146                 &cut[9],&flag[0],&flag[1],&flag[2],&flag[3],&flag[4],&flag[5],&flag[6],&flag[7],
1147                 &flag[8],&flag[9],&flag[10]);
1148     if(!iret) continue;
1149     if(iret<0) {
1150       //reading error
1151       Warning("ReadTransPar","Error reading file %s\n",fTransParName.Data());
1152       continue;
1153     }
1154     // Check that the module exist
1155     AliModule *mod = GetModule(detName);
1156     if(mod) {
1157       // Get the array of media numbers
1158       TArrayI &idtmed = *mod->GetIdtmed();
1159       // Check that the tracking medium code is valid
1160       if(0<=itmed && itmed < 100) {
1161         ktmed=idtmed[itmed];
1162         if(!ktmed) {
1163           Warning("ReadTransPar","Invalid tracking medium code %d for %s\n",itmed,mod->GetName());
1164           continue;
1165         }
1166         // Set energy thresholds
1167         for(kz=0;kz<kncuts;kz++) {
1168           if(cut[kz]>=0) {
1169             if(fDebug) printf(" *  %-6s set to %10.3E for tracking medium code %4d for %s\n",
1170                    kpars[kz],cut[kz],itmed,mod->GetName());
1171             gMC->Gstpar(ktmed,kpars[kz],cut[kz]);
1172           }
1173         }
1174         // Set transport mechanisms
1175         for(kz=0;kz<knflags;kz++) {
1176           if(flag[kz]>=0) {
1177             if(fDebug) printf(" *  %-6s set to %10d for tracking medium code %4d for %s\n",
1178                    kpars[kncuts+kz],flag[kz],itmed,mod->GetName());
1179             gMC->Gstpar(ktmed,kpars[kncuts+kz],Float_t(flag[kz]));
1180           }
1181         }
1182       } else {
1183         Warning("ReadTransPar","Invalid medium code %d *\n",itmed);
1184         continue;
1185       }
1186     } else {
1187       if(fDebug) printf("%s::ReadTransParModule: %s not present\n",ClassName(),detName);
1188       continue;
1189     }
1190   }
1191 }
1192
1193
1194 //_____________________________________________________________________________
1195 void AliRun::MakeTree(Option_t *option, const char *file)
1196 {
1197   //
1198   //  Create the ROOT trees
1199   //  Loop on all detectors to create the Root branch (if any)
1200   //
1201
1202   char hname[30];
1203   //
1204   // Analyse options
1205   const char *oK = strstr(option,"K");
1206   const char *oH = strstr(option,"H");
1207   const char *oE = strstr(option,"E");
1208   const char *oD = strstr(option,"D");
1209   const char *oR = strstr(option,"R");
1210   const char *oS = strstr(option,"S");
1211   //
1212
1213   TDirectory *cwd = gDirectory;
1214
1215   TBranch *branch = 0;
1216   
1217   if (oK) fStack->MakeTree(fEvent, file);
1218
1219   if (oE && !fTreeE) {
1220     fTreeE = new TTree("TE","Header");
1221     branch = fTreeE->Branch("Header", "AliHeader", &fHeader, 4000, 0);         
1222     branch->SetAutoDelete(kFALSE);           
1223     TFolder *folder = (TFolder *)gROOT->FindObjectAny("/Folders/RunMC/Event/Header");
1224     if (folder) folder->Add(fHeader);
1225 //    branch = fTreeE->Branch("Stack","AliStack", &fStack, 4000, 0);
1226 //    branch->SetAutoDelete(kFALSE);         
1227 //    if (folder) folder->Add(fStack);
1228     fTreeE->Write(0,TObject::kOverwrite);
1229   }
1230   
1231   if (file && branch) {
1232     char * outFile = new char[strlen(gAlice->GetBaseFile())+strlen(file)+2];
1233     sprintf(outFile,"%s/%s",GetBaseFile(),file);
1234     branch->SetFile(outFile);
1235     TIter next( branch->GetListOfBranches());
1236     while ((branch=(TBranch*)next())) {
1237        branch->SetFile(outFile);
1238     } 
1239     if (GetDebug()>1)
1240         printf("* MakeBranch * Diverting Branch %s to file %s\n", branch->GetName(),file);
1241     cwd->cd();
1242     delete outFile;
1243   }
1244   
1245   if (oH && !fTreeH) {
1246     sprintf(hname,"TreeH%d",fEvent);
1247     fTreeH = new TTree(hname,"Hits");
1248     fTreeH->SetAutoSave(1000000000); //no autosave
1249     fTreeH->Write(0,TObject::kOverwrite);
1250   }
1251   if (oD && !fTreeD) {
1252     sprintf(hname,"TreeD%d",fEvent);
1253     fTreeD = new TTree(hname,"Digits");
1254     fTreeD->Write(0,TObject::kOverwrite);
1255   }
1256   if (oS && !fTreeS) {
1257     sprintf(hname,"TreeS%d",fEvent);
1258     fTreeS = new TTree(hname,"SDigits");
1259     fTreeS->Write(0,TObject::kOverwrite);
1260   }
1261   if (oR && !fTreeR) {
1262     sprintf(hname,"TreeR%d",fEvent);
1263     fTreeR = new TTree(hname,"Reconstruction");
1264     fTreeR->Write(0,TObject::kOverwrite);
1265   }
1266
1267   //
1268   // Create a branch for hits/digits for each detector
1269   // Each branch is a TClonesArray. Each data member of the Hits classes
1270   // will be in turn a subbranch of the detector master branch
1271   TIter next(fModules);
1272   AliModule *detector;
1273   while((detector = (AliModule*)next())) {
1274      if (oH) detector->MakeBranch(option,file);
1275   }
1276 }
1277
1278 //_____________________________________________________________________________
1279 TParticle* AliRun::Particle(Int_t i)
1280 {
1281     return fStack->Particle(i);
1282 }
1283
1284 //_____________________________________________________________________________
1285 void AliRun::BeginEvent()
1286 {
1287     // Clean-up previous event
1288     // Energy scores  
1289     fEventEnergy.Reset();  
1290     // Clean detector information
1291     CleanDetectors();
1292     // Reset stack info
1293     fStack->Reset();
1294
1295  
1296   //
1297   //  Reset all Detectors & kinematics & trees
1298   //
1299   char hname[30];
1300   //
1301   // Initialise event header
1302   fHeader->Reset(fRun,fEvent);
1303   //
1304   fStack->BeginEvent(fEvent);
1305
1306   //
1307   if(fLego) {
1308     fLego->BeginEvent();
1309     return;
1310   }
1311
1312   //
1313
1314   ResetHits();
1315   ResetDigits();
1316   ResetSDigits();
1317
1318
1319   if(fTreeH) {
1320     fTreeH->Reset();
1321     sprintf(hname,"TreeH%d",fEvent);
1322     fTreeH->SetName(hname);
1323   }
1324   if(fTreeD) {
1325     fTreeD->Reset();
1326     sprintf(hname,"TreeD%d",fEvent);
1327     fTreeD->SetName(hname);
1328     fTreeD->Write(0,TObject::kOverwrite);
1329   }
1330   if(fTreeS) {
1331     fTreeS->Reset();
1332     sprintf(hname,"TreeS%d",fEvent);
1333     fTreeS->SetName(hname);
1334     fTreeS->Write(0,TObject::kOverwrite);
1335   }
1336   if(fTreeR) {
1337     fTreeR->Reset();
1338     sprintf(hname,"TreeR%d",fEvent);
1339     fTreeR->SetName(hname);
1340     fTreeR->Write(0,TObject::kOverwrite);
1341   }
1342 }
1343 //_____________________________________________________________________________
1344 void AliRun::ResetDigits()
1345 {
1346   //
1347   //  Reset all Detectors digits
1348   //
1349   TIter next(fModules);
1350   AliModule *detector;
1351   while((detector = (AliModule*)next())) {
1352      detector->ResetDigits();
1353   }
1354 }
1355
1356 //_____________________________________________________________________________
1357 void AliRun::ResetSDigits()
1358 {
1359   //
1360   //  Reset all Detectors digits
1361   //
1362   TIter next(fModules);
1363   AliModule *detector;
1364   while((detector = (AliModule*)next())) {
1365      detector->ResetSDigits();
1366   }
1367 }
1368
1369 //_____________________________________________________________________________
1370 void AliRun::ResetHits()
1371 {
1372   //
1373   //  Reset all Detectors hits
1374   //
1375   TIter next(fModules);
1376   AliModule *detector;
1377   while((detector = (AliModule*)next())) {
1378      detector->ResetHits();
1379   }
1380 }
1381
1382 //_____________________________________________________________________________
1383 void AliRun::ResetPoints()
1384 {
1385   //
1386   // Reset all Detectors points
1387   //
1388   TIter next(fModules);
1389   AliModule *detector;
1390   while((detector = (AliModule*)next())) {
1391      detector->ResetPoints();
1392   }
1393 }
1394
1395 //_____________________________________________________________________________
1396 void AliRun::RunMC(Int_t nevent, const char *setup)
1397 {
1398   //
1399   // Main function to be called to process a galice run
1400   // example
1401   //    Root > gAlice.Run(); 
1402   // a positive number of events will cause the finish routine
1403   // to be called
1404   //
1405
1406   // check if initialisation has been done
1407   if (!fInitDone) InitMC(setup);
1408   
1409   // Create the Root Tree with one branch per detector
1410
1411    MakeTree("ESDR");
1412
1413   if (gSystem->Getenv("CONFIG_SPLIT_FILE")) {
1414      MakeTree("K","Kine.root");
1415      MakeTree("H","Hits.root");
1416   } else {
1417      MakeTree("KH");
1418   }
1419
1420   gMC->ProcessRun(nevent);
1421
1422   // End of this run, close files
1423   if(nevent>0) FinishRun();
1424 }
1425
1426 //_____________________________________________________________________________
1427 void AliRun::RunReco(const char *selected, Int_t first, Int_t last)
1428 {
1429   //
1430   // Main function to be called to reconstruct Alice event
1431   // 
1432    cout << "Found "<< gAlice->TreeE()->GetEntries() << "events" << endl;
1433    Int_t nFirst = first;
1434    Int_t nLast  = (last < 0)? (Int_t) gAlice->TreeE()->GetEntries() : last;
1435    
1436    for (Int_t nevent = nFirst; nevent <= nLast; nevent++) {
1437      cout << "Processing event "<< nevent << endl;
1438      GetEvent(nevent);
1439      // MakeTree("R");
1440      Digits2Reco(selected);
1441    }
1442 }
1443
1444 //_____________________________________________________________________________
1445
1446 void AliRun::Hits2Digits(const char *selected)
1447 {
1448
1449    // Convert Hits to sumable digits
1450    // 
1451    for (Int_t nevent=0; nevent<gAlice->TreeE()->GetEntries(); nevent++) {
1452      GetEvent(nevent);
1453     // MakeTree("D");
1454      Hits2SDigits(selected);
1455      SDigits2Digits(selected);
1456    }  
1457 }
1458
1459
1460 //_____________________________________________________________________________
1461
1462 void AliRun::Tree2Tree(Option_t *option, const char *selected)
1463 {
1464   //
1465   // Function to transform the content of
1466   //  
1467   // - TreeH to TreeS (option "S")
1468   // - TreeS to TreeD (option "D")
1469   // - TreeD to TreeR (option "R")
1470   // 
1471   // If multiple options are specified ("SDR"), transformation will be done in sequence for
1472   // selected detector and for all detectors if none is selected (detector string 
1473   // can contain blank separated list of detector names). 
1474
1475
1476    const char *oS = strstr(option,"S");
1477    const char *oD = strstr(option,"D");
1478    const char *oR = strstr(option,"R");
1479    
1480    TObjArray *detectors = Detectors();
1481
1482    TIter next(detectors);
1483
1484    AliDetector *detector = 0;
1485
1486    TDirectory *cwd = gDirectory;
1487
1488    char outFile[32];
1489    
1490    while((detector = (AliDetector*)next())) {
1491      if (selected) 
1492        if (strcmp(detector->GetName(),selected)) continue;
1493      if (detector->IsActive()){ 
1494        if (gSystem->Getenv("CONFIG_SPLIT_FILE")) {          
1495           if (oS) {
1496             sprintf(outFile,"SDigits.%s.root",detector->GetName());
1497             detector->MakeBranch("S",outFile);
1498           }    
1499           if (oD) {
1500             sprintf(outFile,"Digits.%s.root",detector->GetName());
1501             detector->MakeBranch("D",outFile);
1502           }    
1503           if (oR) {
1504             sprintf(outFile,"Reco.%s.root",detector->GetName());
1505             detector->MakeBranch("R",outFile);
1506           }    
1507        } else {
1508           detector->MakeBranch(option);
1509        }
1510        
1511        cwd->cd(); 
1512
1513        if (oS) {
1514           cout << "Hits2SDigits: Processing " << detector->GetName() << "..." << endl;
1515           detector->Hits2SDigits(); 
1516        }  
1517        if (oD) {
1518           cout << "SDigits2Digits: Processing " << detector->GetName() << "..." << endl;
1519           detector->SDigits2Digits();
1520        } 
1521        if (oR) {
1522           cout << "Digits2Reco: Processing " << detector->GetName() << "..." << endl;
1523           detector->Digits2Reco(); 
1524        }
1525        
1526        cwd->cd();        
1527      }   
1528    }
1529 }
1530
1531
1532 //_____________________________________________________________________________
1533 void AliRun::RunLego(const char *setup, Int_t nc1, Float_t c1min,
1534                      Float_t c1max,Int_t nc2,Float_t c2min,Float_t c2max,
1535                      Float_t rmin,Float_t rmax,Float_t zmax, AliLegoGenerator* gener)
1536 {
1537   //
1538   // Generates lego plots of:
1539   //    - radiation length map phi vs theta
1540   //    - radiation length map phi vs eta
1541   //    - interaction length map
1542   //    - g/cm2 length map
1543   //
1544   //  ntheta    bins in theta, eta
1545   //  themin    minimum angle in theta (degrees)
1546   //  themax    maximum angle in theta (degrees)
1547   //  nphi      bins in phi
1548   //  phimin    minimum angle in phi (degrees)
1549   //  phimax    maximum angle in phi (degrees)
1550   //  rmin      minimum radius
1551   //  rmax      maximum radius
1552   //  
1553   //
1554   //  The number of events generated = ntheta*nphi
1555   //  run input parameters in macro setup (default="Config.C")
1556   //
1557   //  Use macro "lego.C" to visualize the 3 lego plots in spherical coordinates
1558   //Begin_Html
1559   /*
1560     <img src="picts/AliRunLego1.gif">
1561   */
1562   //End_Html
1563   //Begin_Html
1564   /*
1565     <img src="picts/AliRunLego2.gif">
1566   */
1567   //End_Html
1568   //Begin_Html
1569   /*
1570     <img src="picts/AliRunLego3.gif">
1571   */
1572   //End_Html
1573   //
1574
1575   // check if initialisation has been done
1576   if (!fInitDone) InitMC(setup);
1577   //Save current generator
1578   AliGenerator *gen=Generator();
1579
1580   // Set new generator
1581   if (!gener) gener  = new AliLegoGenerator();
1582   ResetGenerator(gener);
1583   //
1584   // Configure Generator
1585   gener->SetRadiusRange(rmin, rmax);
1586   gener->SetZMax(zmax);
1587   gener->SetCoor1Range(nc1, c1min, c1max);
1588   gener->SetCoor2Range(nc2, c2min, c2max);
1589   
1590   
1591   //Create Lego object  
1592   fLego = new AliLego("lego",gener);
1593
1594   //Prepare MC for Lego Run
1595   gMC->InitLego();
1596   
1597   //Run Lego Object
1598
1599   gMC->ProcessRun(nc1*nc2+1);
1600   
1601   // Create only the Root event Tree
1602   MakeTree("E");
1603   
1604   // End of this run, close files
1605   FinishRun();
1606   // Restore current generator
1607   ResetGenerator(gen);
1608   // Delete Lego Object
1609   delete fLego; fLego=0;
1610 }
1611
1612 //_____________________________________________________________________________
1613 void AliRun::SetConfigFunction(const char * config) 
1614 {
1615   //
1616   // Set the signature of the function contained in Config.C to configure
1617   // the run
1618   //
1619   fConfigFunction=config;
1620 }
1621
1622 //_____________________________________________________________________________
1623 void AliRun::SetCurrentTrack(Int_t track)
1624
1625   //
1626   // Set current track number
1627   //
1628     fStack->SetCurrentTrack(track); 
1629 }
1630  
1631 //_____________________________________________________________________________
1632 void AliRun::SetTrack(Int_t done, Int_t parent, Int_t pdg, Float_t *pmom,
1633                       Float_t *vpos, Float_t *polar, Float_t tof,
1634                       AliMCProcess mech, Int_t &ntr, Float_t weight)
1635
1636 // Delegate to stack
1637 //
1638     fStack->SetTrack(done, parent, pdg, pmom, vpos, polar, tof,
1639                      mech, ntr, weight);
1640 }
1641
1642 //_____________________________________________________________________________
1643 void AliRun::SetTrack(Int_t done, Int_t parent, Int_t pdg,
1644                       Double_t px, Double_t py, Double_t pz, Double_t e,
1645                       Double_t vx, Double_t vy, Double_t vz, Double_t tof,
1646                       Double_t polx, Double_t poly, Double_t polz,
1647                       AliMCProcess mech, Int_t &ntr, Float_t weight)
1648
1649   // Delegate to stack
1650   //
1651     fStack->SetTrack(done, parent, pdg, px, py, pz, e, vx, vy, vz, tof,
1652                    polx, poly, polz, mech, ntr, weight);
1653     
1654 }
1655
1656 //_____________________________________________________________________________
1657 void AliRun::SetHighWaterMark(const Int_t nt)
1658 {
1659     //
1660     // Set high water mark for last track in event
1661     fStack->SetHighWaterMark(nt);
1662 }
1663
1664 //_____________________________________________________________________________
1665 void AliRun::KeepTrack(const Int_t track)
1666
1667   //
1668   // Delegate to stack
1669   //
1670     fStack->KeepTrack(track);
1671 }
1672  
1673 //_____________________________________________________________________________
1674 void AliRun::StepManager(Int_t id) 
1675 {
1676   //
1677   // Called at every step during transport
1678   //
1679
1680   //
1681   // --- If lego option, do it and leave 
1682   if (fLego)
1683     fLego->StepManager();
1684   else {
1685     Int_t copy;
1686     //Update energy deposition tables
1687     AddEnergyDeposit(gMC->CurrentVolID(copy),gMC->Edep());
1688   
1689     //Call the appropriate stepping routine;
1690     AliModule *det = (AliModule*)fModules->At(id);
1691     if(det) {
1692       fMCQA->StepManager(id);
1693       det->StepManager();
1694     }
1695   }
1696 }
1697
1698 //_____________________________________________________________________________
1699 void AliRun::Streamer(TBuffer &R__b)
1700 {
1701    // Stream an object of class AliRun.
1702
1703    if (R__b.IsReading()) {
1704       if (!gAlice) gAlice = this;
1705
1706       AliRun::Class()->ReadBuffer(R__b, this);
1707       //
1708       gROOT->GetListOfBrowsables()->Add(this,"Run");
1709
1710       fTreeE = (TTree*)gDirectory->Get("TE");
1711       if (fTreeE) {
1712           fTreeE->SetBranchAddress("Header", &fHeader);
1713       }
1714       
1715       else    Error("Streamer","cannot find Header Tree\n");
1716       fTreeE->GetEntry(0);
1717
1718       gRandom = fRandom;
1719    } else {
1720       AliRun::Class()->WriteBuffer(R__b, this);
1721    }
1722 }
1723
1724
1725 //___________________________________________________________________________
1726 Int_t AliRun::CurrentTrack() const {
1727   //
1728   // Returns current track
1729   //
1730   return fStack->CurrentTrack();
1731 }
1732
1733 //___________________________________________________________________________
1734 Int_t AliRun::GetNtrack() const {
1735   //
1736   // Returns number of tracks in stack
1737   //
1738   return fStack->GetNtrack();
1739 }
1740
1741 //___________________________________________________________________________
1742 TObjArray* AliRun::Particles() {
1743   //
1744   // Returns pointer to Particles array
1745   //
1746   return fStack->Particles();
1747 }
1748
1749 //___________________________________________________________________________
1750 TTree* AliRun::TreeK() {
1751   //
1752   // Returns pointer to the TreeK array
1753   //
1754   return fStack->TreeK();
1755 }
1756
1757
1758 void AliRun::SetGenEventHeader(AliGenEventHeader* header)
1759 {
1760     fHeader->SetGenEventHeader(header);
1761 }