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