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