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