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