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