]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliRun.cxx
New method IsNewTrack and fix for a problem in Father-Daughter relations
[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.32  2000/04/27 10:38:21  fca
19 Correct termination of Lego Run and introduce Lego getter in AliRun
20
21 Revision 1.31  2000/04/26 10:17:32  fca
22 Changes in Lego for G4 compatibility
23
24 Revision 1.30  2000/04/18 19:11:40  fca
25 Introduce variable Config.C function signature
26
27 Revision 1.29  2000/04/07 11:12:34  fca
28 G4 compatibility changes
29
30 Revision 1.28  2000/04/05 06:51:06  fca
31 Workaround for an HP compiler problem
32
33 Revision 1.27  2000/03/22 18:08:07  fca
34 Rationalisation of the virtual MC interfaces
35
36 Revision 1.26  2000/03/22 13:42:26  fca
37 SetGenerator does not replace an existing generator, ResetGenerator does
38
39 Revision 1.25  2000/02/23 16:25:22  fca
40 AliVMC and AliGeant3 classes introduced
41 ReadEuclid moved from AliRun to AliModule
42
43 Revision 1.24  2000/01/19 17:17:20  fca
44 Introducing a list of lists of hits -- more hits allowed for detector now
45
46 Revision 1.23  1999/12/03 11:14:31  fca
47 Fixing previous wrong checking
48
49 Revision 1.21  1999/11/25 10:40:08  fca
50 Fixing daughters information also in primary tracks
51
52 Revision 1.20  1999/10/04 18:08:49  fca
53 Adding protection against inconsistent Euclid files
54
55 Revision 1.19  1999/09/29 07:50:40  fca
56 Introduction of the Copyright and cvs Log
57
58 */
59
60 ///////////////////////////////////////////////////////////////////////////////
61 //                                                                           //
62 //  Control class for Alice C++                                              //
63 //  Only one single instance of this class exists.                           //
64 //  The object is created in main program aliroot                            //
65 //  and is pointed by the global gAlice.                                     //
66 //                                                                           //
67 //   -Supports the list of all Alice Detectors (fModules).                 //
68 //   -Supports the list of particles (fParticles).                           //
69 //   -Supports the Trees.                                                    //
70 //   -Supports the geometry.                                                 //
71 //   -Supports the event display.                                            //
72 //Begin_Html
73 /*
74 <img src="picts/AliRunClass.gif">
75 */
76 //End_Html
77 //Begin_Html
78 /*
79 <img src="picts/alirun.gif">
80 */
81 //End_Html
82 //                                                                           //
83 ///////////////////////////////////////////////////////////////////////////////
84
85 #include <TFile.h>
86 #include <TRandom.h>
87 #include <TBRIK.h> 
88 #include <TNode.h> 
89 #include <TCint.h> 
90 #include <TSystem.h>
91 #include <TObjectTable.h>
92
93 #include "TParticle.h"
94 #include "AliRun.h"
95 #include "AliDisplay.h"
96 #include "AliMC.h"
97 #include "AliLego.h"
98
99 #include <stdlib.h>
100 #include <stdio.h>
101 #include <string.h>
102  
103 AliRun *gAlice;
104
105 static AliHeader *header;
106
107 ClassImp(AliRun)
108
109 //_____________________________________________________________________________
110 AliRun::AliRun()
111 {
112   //
113   // Default constructor for AliRun
114   //
115   header=&fHeader;
116   fRun       = 0;
117   fEvent     = 0;
118   fCurrent   = -1;
119   fModules = 0;
120   fGenerator = 0;
121   fTreeD     = 0;
122   fTreeK     = 0;
123   fTreeH     = 0;
124   fTreeE     = 0;
125   fTreeR     = 0;
126   fParticles = 0;
127   fGeometry  = 0;
128   fDisplay   = 0;
129   fField     = 0;
130   fMC       = 0;
131   fNdets     = 0;
132   fImedia    = 0;
133   fTrRmax    = 1.e10;
134   fTrZmax    = 1.e10;
135   fInitDone  = kFALSE;
136   fLego      = 0;
137   fPDGDB     = 0;        //Particle factory object!
138   fHitLists  = 0;
139   fConfigFunction    = "\0";
140 }
141
142 //_____________________________________________________________________________
143 AliRun::AliRun(const char *name, const char *title)
144   : TNamed(name,title)
145 {
146   //
147   //  Constructor for the main processor.
148   //  Creates the geometry
149   //  Creates the list of Detectors.
150   //  Creates the list of particles.
151   //
152   Int_t i;
153   
154   gAlice     = this;
155   fTreeD     = 0;
156   fTreeK     = 0;
157   fTreeH     = 0;
158   fTreeE     = 0;
159   fTreeR     = 0;
160   fTrRmax    = 1.e10;
161   fTrZmax    = 1.e10;
162   fGenerator = 0;
163   fInitDone  = kFALSE;
164   fLego      = 0;
165   fField     = 0;
166   fConfigFunction    = "Config();";
167   
168   gROOT->GetListOfBrowsables()->Add(this,name);
169   //
170   // create the support list for the various Detectors
171   fModules = new TObjArray(77);
172   //
173   // Create the TNode geometry for the event display
174   
175   BuildSimpleGeometry();
176   
177
178   fNtrack=0;
179   fHgwmk=0;
180   fCurrent=-1;
181   header=&fHeader;
182   fRun       = 0;
183   fEvent     = 0;
184   //
185   // Create the particle stack
186   fParticles = new TClonesArray("TParticle",100);
187   
188   fDisplay = 0;
189   //
190   // Create default mag field
191   SetField();
192   //
193   fMC      = gMC;
194   //
195   // Prepare the tracking medium lists
196   fImedia = new TArrayI(1000);
197   for(i=0;i<1000;i++) (*fImedia)[i]=-99;
198   //
199   // Make particles
200   fPDGDB     = TDatabasePDG::Instance();        //Particle factory object!
201   //
202   // Create HitLists list
203   fHitLists  = new TList();
204 }
205
206 //_____________________________________________________________________________
207 AliRun::~AliRun()
208 {
209   //
210   // Defaullt AliRun destructor
211   //
212   delete fImedia;
213   delete fField;
214   delete fMC;
215   delete fGeometry;
216   delete fDisplay;
217   delete fGenerator;
218   delete fLego;
219   delete fTreeD;
220   delete fTreeK;
221   delete fTreeH;
222   delete fTreeE;
223   delete fTreeR;
224   if (fModules) {
225     fModules->Delete();
226     delete fModules;
227   }
228   if (fParticles) {
229     fParticles->Delete();
230     delete fParticles;
231   }
232   delete fHitLists;
233 }
234
235 //_____________________________________________________________________________
236 void AliRun::AddHit(Int_t id, Int_t track, Int_t *vol, Float_t *hits) const
237 {
238   //
239   //  Add a hit to detector id
240   //
241   TObjArray &dets = *fModules;
242   if(dets[id]) ((AliModule*) dets[id])->AddHit(track,vol,hits);
243 }
244
245 //_____________________________________________________________________________
246 void AliRun::AddDigit(Int_t id, Int_t *tracks, Int_t *digits) const
247 {
248   //
249   // Add digit to detector id
250   //
251   TObjArray &dets = *fModules;
252   if(dets[id]) ((AliModule*) dets[id])->AddDigit(tracks,digits);
253 }
254
255 //_____________________________________________________________________________
256 void AliRun::Browse(TBrowser *b)
257 {
258   //
259   // Called when the item "Run" is clicked on the left pane
260   // of the Root browser.
261   // It displays the Root Trees and all detectors.
262   //
263   if (fTreeK) b->Add(fTreeK,fTreeK->GetName());
264   if (fTreeH) b->Add(fTreeH,fTreeH->GetName());
265   if (fTreeD) b->Add(fTreeD,fTreeD->GetName());
266   if (fTreeE) b->Add(fTreeE,fTreeE->GetName());
267   if (fTreeR) b->Add(fTreeR,fTreeR->GetName());
268   
269   TIter next(fModules);
270   AliModule *detector;
271   while((detector = (AliModule*)next())) {
272     b->Add(detector,detector->GetName());
273   }
274 }
275
276 //_____________________________________________________________________________
277 void AliRun::Build()
278 {
279   //
280   // Initialize Alice geometry
281   // Dummy routine
282   //
283 }
284  
285 //_____________________________________________________________________________
286 void AliRun::BuildSimpleGeometry()
287 {
288   //
289   // Create a simple TNode geometry used by Root display engine
290   //
291   // Initialise geometry
292   //
293   fGeometry = new TGeometry("AliceGeom","Galice Geometry for Hits");
294   new TMaterial("void","Vacuum",0,0,0);  //Everything is void
295   TBRIK *brik = new TBRIK("S_alice","alice volume","void",2000,2000,3000);
296   brik->SetVisibility(0);
297   new TNode("alice","alice","S_alice");
298 }
299
300 //_____________________________________________________________________________
301 void AliRun::CleanDetectors()
302 {
303   //
304   // Clean Detectors at the end of event
305   //
306   TIter next(fModules);
307   AliModule *detector;
308   while((detector = (AliModule*)next())) {
309     detector->FinishEvent();
310   }
311 }
312
313 //_____________________________________________________________________________
314 void AliRun::CleanParents()
315 {
316   //
317   // Clean Particles stack.
318   // Set parent/daughter relations
319   //
320   TClonesArray &particles = *(gAlice->Particles());
321   TParticle *part;
322   int i;
323   for(i=0; i<fNtrack; i++) {
324     part = (TParticle *)particles.UncheckedAt(i);
325     if(!part->TestBit(Daughters_Bit)) {
326       part->SetFirstDaughter(-1);
327       part->SetLastDaughter(-1);
328     }
329   }
330 }
331
332 //_____________________________________________________________________________
333 Int_t AliRun::DistancetoPrimitive(Int_t, Int_t)
334 {
335   //
336   // Return the distance from the mouse to the AliRun object
337   // Dummy routine
338   //
339   return 9999;
340 }
341
342 //_____________________________________________________________________________
343 void AliRun::DumpPart (Int_t i)
344 {
345   //
346   // Dumps particle i in the stack
347   //
348   TClonesArray &particles = *fParticles;
349   ((TParticle*) particles[i])->Print();
350 }
351
352 //_____________________________________________________________________________
353 void AliRun::DumpPStack ()
354 {
355   //
356   // Dumps the particle stack
357   //
358   TClonesArray &particles = *fParticles;
359   printf(
360          "\n\n=======================================================================\n");
361   for (Int_t i=0;i<fNtrack;i++) 
362     {
363       printf("-> %d ",i); ((TParticle*) particles[i])->Print();
364       printf("--------------------------------------------------------------\n");
365     }
366   printf(
367          "\n=======================================================================\n\n");
368 }
369
370 //_____________________________________________________________________________
371 void AliRun::SetField(Int_t type, Int_t version, Float_t scale,
372                       Float_t maxField, char* filename)
373 {
374   //
375   //  Set magnetic field parameters
376   //  type      Magnetic field transport flag 0=no field, 2=helix, 3=Runge Kutta
377   //  version   Magnetic field map version (only 1 active now)
378   //  scale     Scale factor for the magnetic field
379   //  maxField  Maximum value for the magnetic field
380
381   //
382   // --- Sanity check on mag field flags
383   if(type<0 || type > 2) {
384     Warning("SetField",
385             "Invalid magnetic field flag: %5d; Helix tracking chosen instead\n"
386            ,type);
387     type=2;
388   }
389   if(fField) delete fField;
390   if(version==1) {
391     fField = new AliMagFC("Map1"," ",type,version,scale,maxField);
392   } else if(version<=3) {
393     fField = new AliMagFCM("Map2-3",filename,type,version,scale,maxField);
394     fField->ReadField();
395   } else {
396     Warning("SetField","Invalid map %d\n",version);
397   }
398 }
399
400 //_____________________________________________________________________________
401 void AliRun::FillTree()
402 {
403   //
404   // Fills all AliRun TTrees
405   //
406   if (fTreeK) fTreeK->Fill();
407   if (fTreeH) fTreeH->Fill();
408   if (fTreeD) fTreeD->Fill();
409   if (fTreeR) fTreeR->Fill();
410 }
411  
412 //_____________________________________________________________________________
413 void AliRun::FinishPrimary()
414 {
415   //
416   // Called  at the end of each primary track
417   //
418   
419   //  static Int_t count=0;
420   //  const Int_t times=10;
421   // This primary is finished, purify stack
422   PurifyKine();
423
424   // Write out hits if any
425   if (gAlice->TreeH()) {
426     gAlice->TreeH()->Fill();
427   }
428   
429   // Reset Hits info
430   gAlice->ResetHits();
431
432   //
433   //  if(++count%times==1) gObjectTable->Print();
434 }
435
436 //_____________________________________________________________________________
437 void AliRun::FinishEvent()
438 {
439   //
440   // Called at the end of the event.
441   //
442   
443   //
444   if(fLego) fLego->FinishEvent();
445
446   //Update the energy deposit tables
447   Int_t i;
448   for(i=0;i<fEventEnergy.GetSize();i++) {
449     fSummEnergy[i]+=fEventEnergy[i];
450     fSum2Energy[i]+=fEventEnergy[i]*fEventEnergy[i];
451   }
452   fEventEnergy.Reset();
453   
454   // Clean detector information
455   CleanDetectors();
456   
457   // Write out the kinematics
458   if (fTreeK) {
459     CleanParents();
460     fTreeK->Fill();
461   }
462   
463   // Write out the digits
464   if (fTreeD) {
465     fTreeD->Fill();
466     ResetDigits();
467   }
468   
469   // Write out reconstructed clusters  
470   if (fTreeR) {
471     fTreeR->Fill();
472   }
473
474   // Write out the event Header information
475   if (fTreeE) fTreeE->Fill();
476   
477   // Reset stack info
478   ResetStack();
479   
480   // Write Tree headers
481   //  Int_t ievent = fHeader.GetEvent();
482   //  char hname[30];
483   //  sprintf(hname,"TreeK%d",ievent);
484   if (fTreeK) fTreeK->Write();
485   //  sprintf(hname,"TreeH%d",ievent);
486   if (fTreeH) fTreeH->Write();
487   //  sprintf(hname,"TreeD%d",ievent);
488   if (fTreeD) fTreeD->Write();
489   //  sprintf(hname,"TreeR%d",ievent);
490   if (fTreeR) fTreeR->Write();
491
492   ++fEvent;
493 }
494
495 //_____________________________________________________________________________
496 void AliRun::FinishRun()
497 {
498   //
499   // Called at the end of the run.
500   //
501
502   //
503   if(fLego) fLego->FinishRun();
504
505   // Clean detector information
506   TIter next(fModules);
507   AliModule *detector;
508   while((detector = (AliModule*)next())) {
509     detector->FinishRun();
510   }
511   
512   //Output energy summary tables
513   EnergySummary();
514   
515   // file is retrieved from whatever tree
516   TFile *File = 0;
517   if (fTreeK) File = fTreeK->GetCurrentFile();
518   if ((!File) && (fTreeH)) File = fTreeH->GetCurrentFile();
519   if ((!File) && (fTreeD)) File = fTreeD->GetCurrentFile();
520   if ((!File) && (fTreeE)) File = fTreeE->GetCurrentFile();
521   if( NULL==File ) {
522     Error("FinishRun","There isn't root file!");
523     exit(1);
524   }
525   File->cd();
526   fTreeE->Write();
527   
528   // Clean tree information
529   delete fTreeK; fTreeK = 0;
530   delete fTreeH; fTreeH = 0;
531   delete fTreeD; fTreeD = 0;
532   delete fTreeR; fTreeR = 0;
533   delete fTreeE; fTreeE = 0;
534   
535   // Write AliRun info and all detectors parameters
536   Write();
537   
538   // Close output file
539   File->Write();
540 }
541
542 //_____________________________________________________________________________
543 void AliRun::FlagTrack(Int_t track)
544 {
545   //
546   // Flags a track and all its family tree to be kept
547   //
548   int curr;
549   TParticle *particle;
550
551   curr=track;
552   while(1) {
553     particle=(TParticle*)fParticles->UncheckedAt(curr);
554     
555     // If the particle is flagged the three from here upward is saved already
556     if(particle->TestBit(Keep_Bit)) return;
557     
558     // Save this particle
559     particle->SetBit(Keep_Bit);
560     
561     // Move to father if any
562     if((curr=particle->GetFirstMother())==-1) return;
563   }
564 }
565  
566 //_____________________________________________________________________________
567 void AliRun::EnergySummary()
568 {
569   //
570   // Print summary of deposited energy
571   //
572
573   Int_t ndep=0;
574   Float_t edtot=0;
575   Float_t ed, ed2;
576   Int_t kn, i, left, j, id;
577   const Float_t zero=0;
578   Int_t ievent=fHeader.GetEvent()+1;
579   //
580   // Energy loss information
581   if(ievent) {
582     printf("***************** Energy Loss Information per event (GEV) *****************\n");
583     for(kn=1;kn<fEventEnergy.GetSize();kn++) {
584       ed=fSummEnergy[kn];
585       if(ed>0) {
586         fEventEnergy[ndep]=kn;
587         if(ievent>1) {
588           ed=ed/ievent;
589           ed2=fSum2Energy[kn];
590           ed2=ed2/ievent;
591           ed2=100*TMath::Sqrt(TMath::Max(ed2-ed*ed,zero))/ed;
592         } else 
593           ed2=99;
594         fSummEnergy[ndep]=ed;
595         fSum2Energy[ndep]=TMath::Min((Float_t) 99.,TMath::Max(ed2,zero));
596         edtot+=ed;
597         ndep++;
598       }
599     }
600     for(kn=0;kn<(ndep-1)/3+1;kn++) {
601       left=ndep-kn*3;
602       for(i=0;i<(3<left?3:left);i++) {
603         j=kn*3+i;
604         id=Int_t (fEventEnergy[j]+0.1);
605         printf(" %s %10.3f +- %10.3f%%;",gMC->VolName(id),fSummEnergy[j],fSum2Energy[j]);
606       }
607       printf("\n");
608     }
609     //
610     // Relative energy loss in different detectors
611     printf("******************** Relative Energy Loss per event ********************\n");
612     printf("Total energy loss per event %10.3f GeV\n",edtot);
613     for(kn=0;kn<(ndep-1)/5+1;kn++) {
614       left=ndep-kn*5;
615       for(i=0;i<(5<left?5:left);i++) {
616         j=kn*5+i;
617         id=Int_t (fEventEnergy[j]+0.1);
618         printf(" %s %10.3f%%;",gMC->VolName(id),100*fSummEnergy[j]/edtot);
619       }
620       printf("\n");
621     }
622     for(kn=0;kn<75;kn++) printf("*"); 
623     printf("\n");
624   }
625   //
626   // Reset the TArray's
627   //  fEventEnergy.Set(0);
628   //  fSummEnergy.Set(0);
629   //  fSum2Energy.Set(0);
630 }
631
632 //_____________________________________________________________________________
633 AliModule *AliRun::GetModule(const char *name)
634 {
635   //
636   // Return pointer to detector from name
637   //
638   return (AliModule*)fModules->FindObject(name);
639 }
640  
641 //_____________________________________________________________________________
642 AliDetector *AliRun::GetDetector(const char *name)
643 {
644   //
645   // Return pointer to detector from name
646   //
647   return (AliDetector*)fModules->FindObject(name);
648 }
649  
650 //_____________________________________________________________________________
651 Int_t AliRun::GetModuleID(const char *name)
652 {
653   //
654   // Return galice internal detector identifier from name
655   //
656   Int_t i=-1;
657   TObject *mod=fModules->FindObject(name);
658   if(mod) i=fModules->IndexOf(mod);
659   return i;
660 }
661  
662 //_____________________________________________________________________________
663 Int_t AliRun::GetEvent(Int_t event)
664 {
665   //
666   // Connect the Trees Kinematics and Hits for event # event
667   // Set branch addresses
668   //
669
670   // Reset existing structures
671   ResetStack();
672   ResetHits();
673   ResetDigits();
674   
675   // Delete Trees already connected
676   if (fTreeK) delete fTreeK;
677   if (fTreeH) delete fTreeH;
678   if (fTreeD) delete fTreeD;
679   if (fTreeR) delete fTreeR;
680
681   // Get header from file
682   if(fTreeE) fTreeE->GetEntry(event);
683   else Error("GetEvent","Cannot file Header Tree\n");
684   
685   // Get Kine Tree from file
686   char treeName[20];
687   sprintf(treeName,"TreeK%d",event);
688   fTreeK = (TTree*)gDirectory->Get(treeName);
689   if (fTreeK) fTreeK->SetBranchAddress("Particles", &fParticles);
690   else    Error("GetEvent","cannot find Kine Tree for event:%d\n",event);
691   
692   // Get Hits Tree header from file
693   sprintf(treeName,"TreeH%d",event);
694   fTreeH = (TTree*)gDirectory->Get(treeName);
695   if (!fTreeH) {
696     Error("GetEvent","cannot find Hits Tree for event:%d\n",event);
697   }
698   
699   // Get Digits Tree header from file
700   sprintf(treeName,"TreeD%d",event);
701   fTreeD = (TTree*)gDirectory->Get(treeName);
702   if (!fTreeD) {
703     Warning("GetEvent","cannot find Digits Tree for event:%d\n",event);
704   }
705   
706   
707   // Get Reconstruct Tree header from file
708   sprintf(treeName,"TreeR%d",event);
709   fTreeR = (TTree*)gDirectory->Get(treeName);
710   if (!fTreeR) {
711     //    printf("WARNING: cannot find Reconstructed Tree for event:%d\n",event);
712   }
713  
714   // Set Trees branch addresses
715   TIter next(fModules);
716   AliModule *detector;
717   while((detector = (AliModule*)next())) {
718     detector->SetTreeAddress();
719   }
720   
721   if (fTreeK) fTreeK->GetEvent(0);
722   fNtrack = Int_t (fParticles->GetEntries());
723   return fNtrack;
724 }
725
726 //_____________________________________________________________________________
727 TGeometry *AliRun::GetGeometry()
728 {
729   //
730   // Import Alice geometry from current file
731   // Return pointer to geometry object
732   //
733   if (!fGeometry) fGeometry = (TGeometry*)gDirectory->Get("AliceGeom");
734   //
735   // Unlink and relink nodes in detectors
736   // This is bad and there must be a better way...
737   //
738   
739   TIter next(fModules);
740   AliModule *detector;
741   while((detector = (AliModule*)next())) {
742     detector->SetTreeAddress();
743     TList *dnodes=detector->Nodes();
744     Int_t j;
745     TNode *node, *node1;
746     for ( j=0; j<dnodes->GetSize(); j++) {
747       node = (TNode*) dnodes->At(j);
748       node1 = fGeometry->GetNode(node->GetName());
749       dnodes->Remove(node);
750       dnodes->AddAt(node1,j);
751     }
752   }
753   return fGeometry;
754 }
755
756 //_____________________________________________________________________________
757 void AliRun::GetNextTrack(Int_t &mtrack, Int_t &ipart, Float_t *pmom,
758                           Float_t &e, Float_t *vpos, Float_t *polar,
759                           Float_t &tof)
760 {
761   //
762   // Return next track from stack of particles
763   //
764   TVector3 pol;
765   fCurrent=-1;
766   TParticle *track;
767   for(Int_t i=fNtrack-1; i>=0; i--) {
768     track=(TParticle*) fParticles->UncheckedAt(i);
769     if(!track->TestBit(Done_Bit)) {
770       //
771       // The track has not yet been processed
772       fCurrent=i;
773       ipart=track->GetPdgCode();
774       pmom[0]=track->Px();
775       pmom[1]=track->Py(); 
776       pmom[2]=track->Pz();
777       e     =track->Energy();
778       vpos[0]=track->Vx();
779       vpos[1]=track->Vy();
780       vpos[2]=track->Vz();
781       track->GetPolarisation(pol);
782       polar[0]=pol.X();
783       polar[1]=pol.Y();
784       polar[2]=pol.Z();
785       tof=track->T();
786       track->SetBit(Done_Bit);
787       break;
788     }
789   }
790   mtrack=fCurrent;
791   //
792   // stop and start timer when we start a primary track
793   Int_t nprimaries = fHeader.GetNprimary();
794   if (fCurrent >= nprimaries) return;
795   if (fCurrent < nprimaries-1) {
796     fTimer.Stop();
797     track=(TParticle*) fParticles->UncheckedAt(fCurrent+1);
798     //    track->SetProcessTime(fTimer.CpuTime());
799   }
800   fTimer.Start();
801 }
802
803 //_____________________________________________________________________________
804 Int_t AliRun::GetPrimary(Int_t track)
805 {
806   //
807   // return number of primary that has generated track
808   //
809   int current, parent;
810   TParticle *part;
811   //
812   parent=track;
813   while (1) {
814     current=parent;
815     part = (TParticle *)fParticles->UncheckedAt(current);
816     parent=part->GetFirstMother();
817     if(parent<0) return current;
818   }
819 }
820  
821 //_____________________________________________________________________________
822 void AliRun::InitMC(const char *setup)
823 {
824   //
825   // Initialize the Alice setup
826   //
827
828   gROOT->LoadMacro(setup);
829   gInterpreter->ProcessLine(fConfigFunction.Data());
830
831   gMC->DefineParticles();  //Create standard MC particles
832
833   TObject *objfirst, *objlast;
834
835   fNdets = fModules->GetLast()+1;
836
837   //
838   //=================Create Materials and geometry
839   gMC->Init();
840
841    TIter next(fModules);
842    AliModule *detector;
843    while((detector = (AliModule*)next())) {
844       detector->SetTreeAddress();
845       objlast = gDirectory->GetList()->Last();
846       
847       // Add Detector histograms in Detector list of histograms
848       if (objlast) objfirst = gDirectory->GetList()->After(objlast);
849       else         objfirst = gDirectory->GetList()->First();
850       while (objfirst) {
851         detector->Histograms()->Add(objfirst);
852         objfirst = gDirectory->GetList()->After(objfirst);
853       }
854    }
855    SetTransPar(); //Read the cuts for all materials
856    
857    MediaTable(); //Build the special IMEDIA table
858    
859    //Initialise geometry deposition table
860    fEventEnergy.Set(gMC->NofVolumes()+1);
861    fSummEnergy.Set(gMC->NofVolumes()+1);
862    fSum2Energy.Set(gMC->NofVolumes()+1);
863    
864    //Compute cross-sections
865    gMC->BuildPhysics();
866    
867    //Write Geometry object to current file.
868    fGeometry->Write();
869    
870    fInitDone = kTRUE;
871 }
872
873 //_____________________________________________________________________________
874 void AliRun::MediaTable()
875 {
876   //
877   // Built media table to get from the media number to
878   // the detector id
879   //
880   Int_t kz, nz, idt, lz, i, k, ind;
881   //  Int_t ibeg;
882   TObjArray &dets = *gAlice->Detectors();
883   AliModule *det;
884   //
885   // For all detectors
886   for (kz=0;kz<fNdets;kz++) {
887     // If detector is defined
888     if((det=(AliModule*) dets[kz])) {
889         TArrayI &idtmed = *(det->GetIdtmed()); 
890         for(nz=0;nz<100;nz++) {
891         // Find max and min material number
892         if((idt=idtmed[nz])) {
893           det->LoMedium() = det->LoMedium() < idt ? det->LoMedium() : idt;
894           det->HiMedium() = det->HiMedium() > idt ? det->HiMedium() : idt;
895         }
896       }
897       if(det->LoMedium() > det->HiMedium()) {
898         det->LoMedium() = 0;
899         det->HiMedium() = 0;
900       } else {
901         if(det->HiMedium() > fImedia->GetSize()) {
902           Error("MediaTable","Increase fImedia from %d to %d",
903                 fImedia->GetSize(),det->HiMedium());
904           return;
905         }
906         // Tag all materials in rage as belonging to detector kz
907         for(lz=det->LoMedium(); lz<= det->HiMedium(); lz++) {
908           (*fImedia)[lz]=kz;
909         }
910       }
911     }
912   }
913   //
914   // Print summary table
915   printf(" Traking media ranges:\n");
916   for(i=0;i<(fNdets-1)/6+1;i++) {
917     for(k=0;k< (6<fNdets-i*6?6:fNdets-i*6);k++) {
918       ind=i*6+k;
919       det=(AliModule*)dets[ind];
920       if(det)
921         printf(" %6s: %3d -> %3d;",det->GetName(),det->LoMedium(),
922                det->HiMedium());
923       else
924         printf(" %6s: %3d -> %3d;","NULL",0,0);
925     }
926     printf("\n");
927   }
928 }
929
930 //____________________________________________________________________________
931 void AliRun::SetGenerator(AliGenerator *generator)
932 {
933   //
934   // Load the event generator
935   //
936   if(!fGenerator) fGenerator = generator;
937 }
938
939 //____________________________________________________________________________
940 void AliRun::ResetGenerator(AliGenerator *generator)
941 {
942   //
943   // Load the event generator
944   //
945   if(fGenerator)
946     if(generator)
947       Warning("ResetGenerator","Replacing generator %s with %s\n",
948               fGenerator->GetName(),generator->GetName());
949     else
950       Warning("ResetGenerator","Replacing generator %s with NULL\n",
951               fGenerator->GetName());
952   fGenerator = generator;
953 }
954
955 //____________________________________________________________________________
956 void AliRun::SetTransPar(char* filename)
957 {
958   //
959   // Read filename to set the transport parameters
960   //
961
962
963   const Int_t ncuts=10;
964   const Int_t nflags=11;
965   const Int_t npars=ncuts+nflags;
966   const char pars[npars][7] = {"CUTGAM" ,"CUTELE","CUTNEU","CUTHAD","CUTMUO",
967                                "BCUTE","BCUTM","DCUTE","DCUTM","PPCUTM","ANNI",
968                                "BREM","COMP","DCAY","DRAY","HADR","LOSS",
969                                "MULS","PAIR","PHOT","RAYL"};
970   char line[256];
971   char detName[7];
972   char* filtmp;
973   Float_t cut[ncuts];
974   Int_t flag[nflags];
975   Int_t i, itmed, iret, ktmed, kz;
976   FILE *lun;
977   //
978   // See whether the file is there
979   filtmp=gSystem->ExpandPathName(filename);
980   lun=fopen(filtmp,"r");
981   delete [] filtmp;
982   if(!lun) {
983     Warning("SetTransPar","File %s does not exist!\n",filename);
984     return;
985   }
986   //
987   printf(" "); for(i=0;i<60;i++) printf("*"); printf("\n");
988   printf(" *%59s\n","*");
989   printf(" *       Please check carefully what you are doing!%10s\n","*");
990   printf(" *%59s\n","*");
991   //
992   while(1) {
993     // Initialise cuts and flags
994     for(i=0;i<ncuts;i++) cut[i]=-99;
995     for(i=0;i<nflags;i++) flag[i]=-99;
996     itmed=0;
997     for(i=0;i<256;i++) line[i]='\0';
998     // Read up to the end of line excluded
999     iret=fscanf(lun,"%[^\n]",line);
1000     if(iret<0) {
1001       //End of file
1002       fclose(lun);
1003       printf(" *%59s\n","*");
1004       printf(" "); for(i=0;i<60;i++) printf("*"); printf("\n");
1005       return;
1006     }
1007     // Read the end of line
1008     fscanf(lun,"%*c");
1009     if(!iret) continue;
1010     if(line[0]=='*') continue;
1011     // Read the numbers
1012     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",
1013                 detName,&itmed,&cut[0],&cut[1],&cut[2],&cut[3],&cut[4],&cut[5],&cut[6],&cut[7],&cut[8],
1014                 &cut[9],&flag[0],&flag[1],&flag[2],&flag[3],&flag[4],&flag[5],&flag[6],&flag[7],
1015                 &flag[8],&flag[9],&flag[10]);
1016     if(!iret) continue;
1017     if(iret<0) {
1018       //reading error
1019       Warning("SetTransPar","Error reading file %s\n",filename);
1020       continue;
1021     }
1022     // Check that the module exist
1023     AliModule *mod = GetModule(detName);
1024     if(mod) {
1025       // Get the array of media numbers
1026       TArrayI &idtmed = *mod->GetIdtmed();
1027       // Check that the tracking medium code is valid
1028       if(0<=itmed && itmed < 100) {
1029         ktmed=idtmed[itmed];
1030         if(!ktmed) {
1031           Warning("SetTransPar","Invalid tracking medium code %d for %s\n",itmed,mod->GetName());
1032           continue;
1033         }
1034         // Set energy thresholds
1035         for(kz=0;kz<ncuts;kz++) {
1036           if(cut[kz]>=0) {
1037             printf(" *  %-6s set to %10.3E for tracking medium code %4d for %s\n",
1038                    pars[kz],cut[kz],itmed,mod->GetName());
1039             gMC->Gstpar(ktmed,pars[kz],cut[kz]);
1040           }
1041         }
1042         // Set transport mechanisms
1043         for(kz=0;kz<nflags;kz++) {
1044           if(flag[kz]>=0) {
1045             printf(" *  %-6s set to %10d for tracking medium code %4d for %s\n",
1046                    pars[ncuts+kz],flag[kz],itmed,mod->GetName());
1047             gMC->Gstpar(ktmed,pars[ncuts+kz],Float_t(flag[kz]));
1048           }
1049         }
1050       } else {
1051         Warning("SetTransPar","Invalid medium code %d *\n",itmed);
1052         continue;
1053       }
1054     } else {
1055       Warning("SetTransPar","Module %s not present\n",detName);
1056       continue;
1057     }
1058   }
1059 }
1060
1061 //_____________________________________________________________________________
1062 void AliRun::MakeTree(Option_t *option)
1063 {
1064   //
1065   //  Create the ROOT trees
1066   //  Loop on all detectors to create the Root branch (if any)
1067   //
1068
1069   char hname[30];
1070   //
1071   // Analyse options
1072   char *K = strstr(option,"K");
1073   char *H = strstr(option,"H");
1074   char *E = strstr(option,"E");
1075   char *D = strstr(option,"D");
1076   char *R = strstr(option,"R");
1077   //
1078   if (K && !fTreeK) {
1079     sprintf(hname,"TreeK%d",fEvent);
1080     fTreeK = new TTree(hname,"Kinematics");
1081     //  Create a branch for particles
1082     fTreeK->Branch("Particles",&fParticles,4000);
1083   }
1084   if (H && !fTreeH) {
1085     sprintf(hname,"TreeH%d",fEvent);
1086     fTreeH = new TTree(hname,"Hits");
1087     fTreeH->SetAutoSave(1000000000); //no autosave
1088   }
1089   if (D && !fTreeD) {
1090     sprintf(hname,"TreeD%d",fEvent);
1091     fTreeD = new TTree(hname,"Digits");
1092   }
1093   if (R && !fTreeR) {
1094     sprintf(hname,"TreeR%d",fEvent);
1095     fTreeR = new TTree(hname,"Reconstruction");
1096   }
1097   if (E && !fTreeE) {
1098     fTreeE = new TTree("TE","Header");
1099     //  Create a branch for Header
1100     fTreeE->Branch("Header","AliHeader",&header,4000);
1101   }
1102   //
1103   // Create a branch for hits/digits for each detector
1104   // Each branch is a TClonesArray. Each data member of the Hits classes
1105   // will be in turn a subbranch of the detector master branch
1106   TIter next(fModules);
1107   AliModule *detector;
1108   while((detector = (AliModule*)next())) {
1109      if (H || D || R) detector->MakeBranch(option);
1110   }
1111 }
1112
1113 //_____________________________________________________________________________
1114 Int_t AliRun::PurifyKine(Int_t lastSavedTrack, Int_t nofTracks)
1115 {
1116   //
1117   // PurifyKine with external parameters
1118   //
1119   fHgwmk = lastSavedTrack;
1120   fNtrack = nofTracks;
1121   PurifyKine();
1122   return fHgwmk;
1123 }
1124
1125 //_____________________________________________________________________________
1126 void AliRun::PurifyKine()
1127 {
1128   //
1129   // Compress kinematic tree keeping only flagged particles
1130   // and renaming the particle id's in all the hits
1131   //
1132   TClonesArray &particles = *fParticles;
1133   int nkeep=fHgwmk+1, parent, i;
1134   TParticle *part, *partnew, *father;
1135   int *map = new int[particles.GetEntries()];
1136
1137   // Save in Header total number of tracks before compression
1138   fHeader.SetNtrack(fHeader.GetNtrack()+fNtrack-fHgwmk);
1139
1140   // First pass, invalid Daughter information
1141   for(i=0; i<fNtrack; i++) {
1142     // Preset map, to be removed later
1143     if(i<=fHgwmk) map[i]=i ; else map[i] = -99;
1144     ((TParticle *)particles.UncheckedAt(i))->ResetBit(Daughters_Bit);
1145   }
1146   // Second pass, build map between old and new numbering
1147   for(i=fHgwmk+1; i<fNtrack; i++) {
1148     part = (TParticle *)particles.UncheckedAt(i);
1149     if(part->TestBit(Keep_Bit)) {
1150       
1151       // This particle has to be kept
1152       map[i]=nkeep;
1153       if(i!=nkeep) {
1154         
1155         // Old and new are different, have to copy
1156         partnew = (TParticle *)particles.UncheckedAt(nkeep);
1157         // Change due to a bug in the HP compiler
1158         //      *partnew = *part;
1159         memcpy(partnew,part,sizeof(TParticle));
1160       } else partnew = part;
1161       
1162       // as the parent is always *before*, it must be already
1163       // in place. This is what we are checking anyway!
1164       if((parent=partnew->GetFirstMother())>fHgwmk) {
1165         if(map[parent]==-99) printf("map[%d] = -99!\n",parent);
1166         partnew->SetFirstMother(map[parent]);
1167       }
1168       nkeep++;
1169     }
1170   }
1171   fNtrack=nkeep;
1172   
1173   // Fix daughters information
1174   for (i=0; i<fNtrack; i++) {
1175     part = (TParticle *)particles.UncheckedAt(i);
1176     parent = part->GetFirstMother();
1177     if(parent>=0) {
1178       father = (TParticle *)particles.UncheckedAt(parent);
1179       if(father->TestBit(Daughters_Bit)) {
1180       
1181         if(i<father->GetFirstDaughter()) father->SetFirstDaughter(i);
1182         if(i>father->GetLastDaughter())  father->SetLastDaughter(i);
1183       } else {
1184         // Iitialise daughters info for first pass
1185         father->SetFirstDaughter(i);
1186         father->SetLastDaughter(i);
1187         father->SetBit(Daughters_Bit);
1188       }
1189     }
1190   }
1191   
1192 #ifdef old
1193   // Now loop on all detectors and reset the hits
1194   AliHit *OneHit;
1195   TIter next(fModules);
1196   AliModule *detector;
1197   while((detector = (AliModule*)next())) {
1198     if (!detector->Hits()) continue;
1199     TClonesArray &vHits=*(detector->Hits());
1200     if(vHits.GetEntries() != detector->GetNhits())
1201       printf("vHits.GetEntries()!=detector->GetNhits(): %d != %d\n",
1202              vHits.GetEntries(),detector->GetNhits());
1203     for (i=0; i<detector->GetNhits(); i++) {
1204       OneHit = (AliHit *)vHits.UncheckedAt(i);
1205       OneHit->SetTrack(map[OneHit->GetTrack()]);
1206     }
1207   }
1208 #else
1209
1210   // Now loop on all registered hit lists
1211   TIter next(fHitLists);
1212   TCollection *hitList;
1213   while((hitList = (TCollection*)next())) {
1214     TIter nexthit(hitList);
1215     AliHit *hit;
1216     while((hit = (AliHit*)nexthit())) {
1217       hit->SetTrack(map[hit->GetTrack()]);
1218     }
1219   }
1220 #endif
1221
1222   fHgwmk=nkeep-1;
1223   particles.SetLast(fHgwmk);
1224   delete [] map;
1225 }
1226
1227 //_____________________________________________________________________________
1228 void AliRun::BeginEvent()
1229 {
1230   //
1231   //  Reset all Detectors & kinematics & trees
1232   //
1233   char hname[30];
1234   //
1235
1236   //
1237   if(fLego) {
1238     fLego->BeginEvent();
1239     return;
1240   }
1241
1242   //
1243   ResetStack();
1244   ResetHits();
1245   ResetDigits();
1246
1247   // Initialise event header
1248   fHeader.Reset(fRun,fEvent);
1249
1250   if(fTreeK) {
1251     fTreeK->Reset();
1252     sprintf(hname,"TreeK%d",fEvent);
1253     fTreeK->SetName(hname);
1254   }
1255   if(fTreeH) {
1256     fTreeH->Reset();
1257     sprintf(hname,"TreeH%d",fEvent);
1258     fTreeH->SetName(hname);
1259   }
1260   if(fTreeD) {
1261     fTreeD->Reset();
1262     sprintf(hname,"TreeD%d",fEvent);
1263     fTreeD->SetName(hname);
1264   }
1265   if(fTreeR) {
1266     fTreeR->Reset();
1267     sprintf(hname,"TreeR%d",fEvent);
1268     fTreeR->SetName(hname);
1269   }
1270 }
1271
1272 //_____________________________________________________________________________
1273 void AliRun::ResetDigits()
1274 {
1275   //
1276   //  Reset all Detectors digits
1277   //
1278   TIter next(fModules);
1279   AliModule *detector;
1280   while((detector = (AliModule*)next())) {
1281      detector->ResetDigits();
1282   }
1283 }
1284
1285 //_____________________________________________________________________________
1286 void AliRun::ResetHits()
1287 {
1288   //
1289   //  Reset all Detectors hits
1290   //
1291   TIter next(fModules);
1292   AliModule *detector;
1293   while((detector = (AliModule*)next())) {
1294      detector->ResetHits();
1295   }
1296 }
1297
1298 //_____________________________________________________________________________
1299 void AliRun::ResetPoints()
1300 {
1301   //
1302   // Reset all Detectors points
1303   //
1304   TIter next(fModules);
1305   AliModule *detector;
1306   while((detector = (AliModule*)next())) {
1307      detector->ResetPoints();
1308   }
1309 }
1310
1311 //_____________________________________________________________________________
1312 void AliRun::RunMC(Int_t nevent, const char *setup)
1313 {
1314   //
1315   // Main function to be called to process a galice run
1316   // example
1317   //    Root > gAlice.Run(); 
1318   // a positive number of events will cause the finish routine
1319   // to be called
1320   //
1321
1322   // check if initialisation has been done
1323   if (!fInitDone) InitMC(setup);
1324   
1325   // Create the Root Tree with one branch per detector
1326   MakeTree("KHDER");
1327
1328   gMC->ProcessRun(nevent);
1329
1330   // End of this run, close files
1331   if(nevent>0) FinishRun();
1332 }
1333
1334 //_____________________________________________________________________________
1335 void AliRun::RunLego(const char *setup,Int_t ntheta,Float_t themin,
1336                      Float_t themax,Int_t nphi,Float_t phimin,Float_t phimax,
1337                      Float_t rmin,Float_t rmax,Float_t zmax)
1338 {
1339   //
1340   // Generates lego plots of:
1341   //    - radiation length map phi vs theta
1342   //    - radiation length map phi vs eta
1343   //    - interaction length map
1344   //    - g/cm2 length map
1345   //
1346   //  ntheta    bins in theta, eta
1347   //  themin    minimum angle in theta (degrees)
1348   //  themax    maximum angle in theta (degrees)
1349   //  nphi      bins in phi
1350   //  phimin    minimum angle in phi (degrees)
1351   //  phimax    maximum angle in phi (degrees)
1352   //  rmin      minimum radius
1353   //  rmax      maximum radius
1354   //  
1355   //
1356   //  The number of events generated = ntheta*nphi
1357   //  run input parameters in macro setup (default="Config.C")
1358   //
1359   //  Use macro "lego.C" to visualize the 3 lego plots in spherical coordinates
1360   //Begin_Html
1361   /*
1362     <img src="picts/AliRunLego1.gif">
1363   */
1364   //End_Html
1365   //Begin_Html
1366   /*
1367     <img src="picts/AliRunLego2.gif">
1368   */
1369   //End_Html
1370   //Begin_Html
1371   /*
1372     <img src="picts/AliRunLego3.gif">
1373   */
1374   //End_Html
1375   //
1376
1377   // check if initialisation has been done
1378   if (!fInitDone) InitMC(setup);
1379
1380   //Save current generator
1381   AliGenerator *gen=Generator();
1382
1383   //Create Lego object  
1384   fLego = new AliLego("lego",ntheta,themin,themax,nphi,phimin,phimax,rmin,rmax,zmax);
1385
1386   //Prepare MC for Lego Run
1387   gMC->InitLego();
1388   
1389   //Run Lego Object
1390   gMC->ProcessRun(ntheta*nphi+1);
1391   
1392   // Create only the Root event Tree
1393   MakeTree("E");
1394   
1395   // End of this run, close files
1396   FinishRun();
1397
1398   // Delete Lego Object
1399   delete fLego; fLego=0;
1400
1401   // Restore current generator
1402   SetGenerator(gen);
1403 }
1404
1405 //_____________________________________________________________________________
1406 void AliRun::SetCurrentTrack(Int_t track)
1407
1408   //
1409   // Set current track number
1410   //
1411   fCurrent = track; 
1412 }
1413  
1414 //_____________________________________________________________________________
1415 void AliRun::SetTrack(Int_t done, Int_t parent, Int_t pdg, Float_t *pmom,
1416                       Float_t *vpos, Float_t *polar, Float_t tof,
1417                       const char *mecha, Int_t &ntr, Float_t weight)
1418
1419   //
1420   // Load a track on the stack
1421   //
1422   // done     0 if the track has to be transported
1423   //          1 if not
1424   // parent   identifier of the parent track. -1 for a primary
1425   // pdg    particle code
1426   // pmom     momentum GeV/c
1427   // vpos     position 
1428   // polar    polarisation 
1429   // tof      time of flight in seconds
1430   // mecha    production mechanism
1431   // ntr      on output the number of the track stored
1432   //
1433   TClonesArray &particles = *fParticles;
1434   TParticle *particle;
1435   Float_t mass;
1436   const Int_t firstdaughter=-1;
1437   const Int_t lastdaughter=-1;
1438   const Int_t KS=0;
1439   //  const Float_t tlife=0;
1440   
1441   //
1442   // Here we get the static mass
1443   // For MC is ok, but a more sophisticated method could be necessary
1444   // if the calculated mass is required
1445   // also, this method is potentially dangerous if the mass
1446   // used in the MC is not the same of the PDG database
1447   //
1448   mass = TDatabasePDG::Instance()->GetParticle(pdg)->Mass();
1449   Float_t e=TMath::Sqrt(mass*mass+pmom[0]*pmom[0]+
1450                         pmom[1]*pmom[1]+pmom[2]*pmom[2]);
1451   
1452   //printf("Loading particle %s mass %f ene %f No %d ip %d pos %f %f %f mom %f %f %f KS %d m %s\n",
1453   //pname,mass,e,fNtrack,pdg,vpos[0],vpos[1],vpos[2],pmom[0],pmom[1],pmom[2],KS,mecha);
1454   
1455   particle=new(particles[fNtrack]) TParticle(pdg,KS,parent,-1,firstdaughter,
1456                                              lastdaughter,pmom[0],pmom[1],pmom[2],
1457                                              e,vpos[0],vpos[1],vpos[2],tof);
1458   //                                         polar[0],polar[1],polar[2],tof,
1459   //                                         mecha,weight);
1460   ((TParticle*)particles[fNtrack])->SetPolarisation(TVector3(polar[0],polar[1],polar[2]));
1461   ((TParticle*)particles[fNtrack])->SetWeight(weight);
1462   if(!done) particle->SetBit(Done_Bit);
1463   //Declare that the daughter information is valid
1464   ((TParticle*)particles[fNtrack])->SetBit(Daughters_Bit);
1465   
1466   if(parent>=0) {
1467     particle=(TParticle*) fParticles->UncheckedAt(parent);
1468     particle->SetLastDaughter(fNtrack);
1469     if(particle->GetFirstDaughter()<0) particle->SetFirstDaughter(fNtrack);
1470   } else { 
1471     //
1472     // This is a primary track. Set high water mark for this event
1473     fHgwmk=fNtrack;
1474     //
1475     // Set also number if primary tracks
1476     fHeader.SetNprimary(fHgwmk+1);
1477     fHeader.SetNtrack(fHgwmk+1);
1478   }
1479   ntr = fNtrack++;
1480 }
1481
1482 //_____________________________________________________________________________
1483 void AliRun::KeepTrack(const Int_t track)
1484
1485   //
1486   // flags a track to be kept
1487   //
1488   TClonesArray &particles = *fParticles;
1489   ((TParticle*)particles[track])->SetBit(Keep_Bit);
1490 }
1491  
1492 //_____________________________________________________________________________
1493 void AliRun::StepManager(Int_t id) 
1494 {
1495   //
1496   // Called at every step during transport
1497   //
1498
1499   //
1500   // --- If lego option, do it and leave 
1501   if (fLego)
1502     fLego->StepManager();
1503   else {
1504     Int_t copy;
1505     //Update energy deposition tables
1506     AddEnergyDeposit(gMC->CurrentVolID(copy),gMC->Edep());
1507   
1508     //Call the appropriate stepping routine;
1509     AliModule *det = (AliModule*)fModules->At(id);
1510     if(det) det->StepManager();
1511   }
1512 }
1513
1514 //_____________________________________________________________________________
1515 void AliRun::Streamer(TBuffer &R__b)
1516 {
1517   //
1518   // Stream an object of class AliRun.
1519   //
1520   if (R__b.IsReading()) {
1521     Version_t R__v = R__b.ReadVersion(); if (R__v) { }
1522     TNamed::Streamer(R__b);
1523     if (!gAlice) gAlice = this;
1524     gROOT->GetListOfBrowsables()->Add(this,"Run");
1525     fTreeE = (TTree*)gDirectory->Get("TE");
1526     if (fTreeE) fTreeE->SetBranchAddress("Header", &header);
1527     else    Error("Streamer","cannot find Header Tree\n");
1528     R__b >> fNtrack;
1529     R__b >> fHgwmk;
1530     R__b >> fDebug;
1531     fHeader.Streamer(R__b);
1532     R__b >> fModules;
1533     R__b >> fParticles;
1534     R__b >> fField; 
1535     //    R__b >> fMC;
1536     R__b >> fNdets;
1537     R__b >> fTrRmax;
1538     R__b >> fTrZmax;
1539     R__b >> fGenerator;
1540     if(R__v>1) {
1541       R__b >> fPDGDB;        //Particle factory object!
1542       fTreeE->GetEntry(0);
1543     } else {
1544       fHeader.SetEvent(0);
1545       fPDGDB     = TDatabasePDG::Instance();        //Particle factory object!
1546     }
1547     if(R__v>2) {
1548       fConfigFunction.Streamer(R__b);
1549     } else {
1550       fConfigFunction="Config();";
1551     }
1552   } else {
1553     R__b.WriteVersion(AliRun::IsA());
1554     TNamed::Streamer(R__b);
1555     R__b << fNtrack;
1556     R__b << fHgwmk;
1557     R__b << fDebug;
1558     fHeader.Streamer(R__b);
1559     R__b << fModules;
1560     R__b << fParticles;
1561     R__b << fField;
1562     //    R__b << fMC;
1563     R__b << fNdets;
1564     R__b << fTrRmax;
1565     R__b << fTrZmax;
1566     R__b << fGenerator;
1567     R__b << fPDGDB;        //Particle factory object!
1568     fConfigFunction.Streamer(R__b);
1569   }
1570