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