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