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