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