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