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