]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliRun.cxx
Adding the new MUONcalib library (Ivana)
[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 /* $Id$ */
17
18 ///////////////////////////////////////////////////////////////////////////////
19 //                                                                           //
20 //  Control class for Alice C++                                              //
21 //  Only one single instance of this class exists.                           //
22 //  The object is created in main program aliroot                            //
23 //  and is pointed by the global gAlice.                                     //
24 //                                                                           //
25 //   -Supports the list of all Alice Detectors (fModules).                 //
26 //   -Supports the list of particles (fParticles).                           //
27 //   -Supports the Trees.                                                    //
28 //   -Supports the geometry.                                                 //
29 //   -Supports the event display.                                            //
30 //Begin_Html
31 /*
32 <img src="picts/AliRunClass.gif">
33 */
34 //End_Html
35 //Begin_Html
36 /*
37 <img src="picts/alirun.gif">
38 */
39 //End_Html
40 //                                                                           //
41 ///////////////////////////////////////////////////////////////////////////////
42
43 #include <TBRIK.h> 
44 #include <TCint.h> 
45 #include <TDatabasePDG.h>
46 #include <TGeometry.h>
47 #include <TNode.h>
48 #include <TROOT.h>
49 #include <TRandom3.h>
50 #include <TSystem.h>
51 #include <TVirtualMC.h>
52 #include <TGeoManager.h>
53 // 
54 #include "AliLog.h"
55 #include "AliDetector.h"
56 #include "AliDisplay.h"
57 #include "AliHeader.h"
58 #include "AliLego.h"
59 #include "AliLegoGenerator.h"
60 #include "AliMC.h"
61 #include "AliMagFC.h"
62 #include "AliMagFCM.h"
63 #include "AliMagFDM.h"
64 #include "AliPDG.h"
65 #include "AliRun.h"
66 #include "AliStack.h"
67 #include "AliAlignObj.h"
68
69 AliRun *gAlice;
70
71 ClassImp(AliRun)
72
73 //_______________________________________________________________________
74 AliRun::AliRun():
75   fRun(0),
76   fEvent(0),
77   fEventNrInRun(0),
78   fEventsPerRun(0),
79   fModules(0),
80   fGeometry(0),
81   fMCApp(0),
82   fDisplay(0),
83   fField(0),
84   fNdets(0),
85   fInitDone(kFALSE),
86   fLego(0),
87   fPDGDB(0),  //Particle factory object
88   fConfigFunction(""),
89   fRandom(0),
90   fBaseFileName(""),
91   fIsRootGeometry(kFALSE),
92   fGeometryFileName(""),
93   fTriggerDescriptor(""),
94   fRunLoader(0x0)
95 {
96   //
97   // Default constructor for AliRun
98   //
99   AliConfig::Instance();//skowron 29 Feb 2002
100                         //ensures that the folder structure is build
101 }
102
103 //_______________________________________________________________________
104 AliRun::AliRun(const AliRun& arun):
105   TNamed(arun),
106   fRun(0),
107   fEvent(0),
108   fEventNrInRun(0),
109   fEventsPerRun(0),
110   fModules(0),
111   fGeometry(0),
112   fMCApp(0),
113   fDisplay(0),
114   fField(0),
115   fNdets(0),
116   fInitDone(kFALSE),
117   fLego(0),
118   fPDGDB(0),  //Particle factory object
119   fConfigFunction("\0"),
120   fRandom(0),
121   fBaseFileName(""),
122   fIsRootGeometry(kFALSE),
123   fGeometryFileName(""),
124   fTriggerDescriptor(""),
125   fRunLoader(0x0)
126 {
127   //
128   // Copy constructor for AliRun
129   //
130   arun.Copy(*this);
131 }
132
133 //_____________________________________________________________________________
134 AliRun::AliRun(const char *name, const char *title):
135   TNamed(name,title),
136   fRun(0),
137   fEvent(0),
138   fEventNrInRun(0),
139   fEventsPerRun(0),
140   fModules(new TObjArray(77)), // Support list for the Detectors
141   fGeometry(0),
142   fMCApp(0),
143   fDisplay(0),
144   fField(0),
145   fNdets(0),
146   fInitDone(kFALSE),
147   fLego(0),
148   fPDGDB(TDatabasePDG::Instance()),        //Particle factory object!
149   fConfigFunction("Config();"),
150   fRandom(new TRandom3()),
151   fBaseFileName(""),
152   fIsRootGeometry(kFALSE),
153   fGeometryFileName(""),
154   fTriggerDescriptor(""),
155   fRunLoader(0x0)
156 {
157   //
158   //  Constructor for the main processor.
159   //  Creates the geometry
160   //  Creates the list of Detectors.
161   //  Creates the list of particles.
162   //
163
164   gAlice     = this;
165
166   // Set random number generator
167   gRandom = fRandom;
168
169   if (gSystem->Getenv("CONFIG_SEED")) {
170      gRandom->SetSeed(static_cast<UInt_t>(atoi(gSystem->Getenv("CONFIG_SEED"))));
171   }
172
173   // Add to list of browsable  
174   gROOT->GetListOfBrowsables()->Add(this,name);
175   
176   // Create default mag field
177   SetField();
178
179   // Add particle list to configuration
180   AliConfig::Instance()->Add(fPDGDB); 
181
182 }
183
184
185 //_______________________________________________________________________
186 AliRun::~AliRun()
187 {
188   //
189   // Default AliRun destructor
190   //
191   gROOT->GetListOfBrowsables()->Remove(this);
192
193   if (fRunLoader)
194    {
195     TFolder* evfold = fRunLoader->GetEventFolder();
196     TFolder* modfold = dynamic_cast<TFolder*>(evfold->FindObjectAny(AliConfig::GetModulesFolderName()));
197     TIter next(fModules);
198     AliModule *mod;
199     while((mod = (AliModule*)next()))
200      { 
201        modfold->Remove(mod);
202      }
203    }
204   
205   
206   delete fField;
207   delete fMCApp;
208   delete gMC; gMC=0;
209   delete fGeometry;
210   delete fDisplay;
211   delete fLego;
212   if (fModules) {
213     fModules->Delete();
214     delete fModules;
215   }
216   
217   delete fPDGDB;
218 }
219
220 //_______________________________________________________________________
221 void AliRun::Copy(TObject &) const
222 {
223   AliFatal("Not implemented!");
224 }
225
226 //_______________________________________________________________________
227 void AliRun::Build()
228 {
229   //
230   // Initialize Alice geometry
231   // Dummy routine
232   //
233 }
234  
235 //_______________________________________________________________________
236 void AliRun::BuildSimpleGeometry()
237 {
238   //
239   // Create a simple TNode geometry used by Root display engine
240   //
241   // Initialise geometry
242   //
243   fGeometry = new TGeometry("AliceGeom","Galice Geometry for Hits");
244   new TMaterial("void","Vacuum",0,0,0);  //Everything is void
245   TBRIK *brik = new TBRIK("S_alice","alice volume","void",2000,2000,3000);
246   brik->SetVisibility(0);
247   new TNode("alice","alice","S_alice");
248 }
249
250 //_______________________________________________________________________
251 void AliRun::CleanDetectors()
252 {
253   //
254   // Clean Detectors at the end of event
255   //
256    fRunLoader->CleanDetectors();
257 }
258
259 //_______________________________________________________________________
260 void AliRun::ResetHits() 
261 {
262   fMCApp->ResetHits();
263 }
264
265 //_______________________________________________________________________
266 AliGenerator* AliRun::Generator() const 
267 {
268   return fMCApp->Generator();
269 }
270
271 //_______________________________________________________________________
272 void  AliRun::SetField(AliMagF* magField)
273 {
274   //
275   // Set Magnetic Field Map
276   //
277   fField = magField;
278   fField->ReadField();
279 }
280
281 //_______________________________________________________________________
282 void AliRun::SetRootGeometry(Bool_t flag)
283 {
284 // Instruct application that the geometry is to be retreived from a root file.
285    fIsRootGeometry = flag;
286    if (flag) gMC->SetRootGeometry();
287 }
288 //_______________________________________________________________________
289 void AliRun::SetField(Int_t type, Int_t version, Float_t scale,
290                       Float_t maxField, const char* filename)
291 {
292   //
293   //  Set magnetic field parameters
294   //  type      Magnetic field transport flag 0=no field, 2=helix, 3=Runge Kutta
295   //  version   Magnetic field map version (only 1 active now)
296   //  scale     Scale factor for the magnetic field
297   //  maxField  Maximum value for the magnetic field
298
299   //
300   // --- Sanity check on mag field flags
301   if(fField) delete fField;
302   if(version==1) {
303     fField = new AliMagFC("Map1"," ",type,scale,maxField);
304   } else if(version<=2) {
305     fField = new AliMagFCM("Map2-3",filename,type,scale,maxField);
306     fField->ReadField();
307   } else if(version==3) {
308     fField = new AliMagFDM("Map4",filename,type,scale,maxField);
309     fField->ReadField();
310   } else {
311     AliWarning(Form("Invalid map %d",version));
312   }
313 }
314
315 //_____________________________________________________________________________
316
317 void AliRun::InitLoaders()
318 {
319   //creates list of getters
320   AliDebug(1, "");
321   TIter next(fModules);
322   AliModule *mod;
323   while((mod = (AliModule*)next()))
324    { 
325      mod->SetRunLoader(fRunLoader);
326      AliDetector *det = dynamic_cast<AliDetector*>(mod);
327      if (det) 
328       {
329         AliDebug(2, Form("Adding %s", det->GetName()));
330         fRunLoader->AddLoader(det);
331       }
332    }
333   AliDebug(1, "Done");
334 }
335 //_____________________________________________________________________________
336
337 void AliRun::FinishRun()
338 {
339   //
340   // Called at the end of the run.
341   //
342   
343   if(fLego) 
344    {
345     AliDebug(1, "Finish Lego");
346     fRunLoader->CdGAFile();
347     fLego->FinishRun();
348    }
349   
350   // Clean detector information
351   TIter next(fModules);
352   AliModule *detector;
353   while((detector = dynamic_cast<AliModule*>(next()))) {
354     AliDebug(2, Form("%s->FinishRun()", detector->GetName()));
355     detector->FinishRun();
356   }
357   
358   AliDebug(1, "fRunLoader->WriteHeader(OVERWRITE)");
359   fRunLoader->WriteHeader("OVERWRITE");
360
361   // Write AliRun info and all detectors parameters
362   fRunLoader->CdGAFile();
363   Write(0,TObject::kOverwrite);//write AliRun
364   fRunLoader->Write(0,TObject::kOverwrite);//write RunLoader itself
365   
366   // Clean tree information
367   AliDebug(1, "fRunLoader->Stack()->FinishRun()");
368   fRunLoader->Stack()->FinishRun();
369
370   if(fMCApp) fMCApp->FinishRun();
371
372   fRunLoader->Synchronize();
373 }
374
375 //_______________________________________________________________________
376 void AliRun::Announce() const
377 {
378   //
379   // Announce the current version of AliRoot
380   //
381   printf("%70s",
382          "****************************************************************\n");
383   printf("%6s","*");printf("%64s","*\n");
384
385   printf("%6s","*");
386   printf("    You are running AliRoot version NewIO\n");
387
388   printf("%6s","*");
389   printf("    The cvs tag for the current program is $Name$\n");
390
391   printf("%6s","*");printf("%64s","*\n");
392   printf("%70s",
393          "****************************************************************\n");
394 }
395
396 //_______________________________________________________________________
397 AliModule *AliRun::GetModule(const char *name) const
398 {
399   //
400   // Return pointer to detector from name
401   //
402   return dynamic_cast<AliModule*>(fModules->FindObject(name));
403 }
404  
405 //_______________________________________________________________________
406 AliDetector *AliRun::GetDetector(const char *name) const
407 {
408   //
409   // Return pointer to detector from name
410   //
411   return dynamic_cast<AliDetector*>(fModules->FindObject(name));
412 }
413  
414 //_______________________________________________________________________
415 Int_t AliRun::GetModuleID(const char *name) const
416 {
417   //
418   // Return galice internal detector identifier from name
419   //
420   Int_t i=-1;
421   TObject *mod=fModules->FindObject(name);
422   if(mod) i=fModules->IndexOf(mod);
423   return i;
424 }
425  
426 //_______________________________________________________________________
427 Int_t AliRun::GetEvent(Int_t event)
428 {
429 //
430 // Reloads data containers in folders # event
431 // Set branch addresses
432 //
433   if (fRunLoader == 0x0)
434    {
435      AliError("RunLoader is not set. Can not load data.");
436      return -1;
437    }
438 /*****************************************/ 
439 /****   P R E    R E L O A D I N G    ****/
440 /*****************************************/ 
441 // Reset existing structures
442   fMCApp->ResetHits();
443   fMCApp->ResetTrackReferences();
444   ResetDigits();
445   ResetSDigits();
446
447 /*****************************************/ 
448 /****       R  E  L  O  A  D          ****/
449 /*****************************************/
450
451   fRunLoader->GetEvent(event);
452
453 /*****************************************/ 
454 /****  P O S T    R E L O A D I N G   ****/
455 /*****************************************/ 
456
457   // Set Trees branch addresses
458   TIter next(fModules);
459   AliModule *detector;
460   while((detector = dynamic_cast<AliModule*>(next()))) 
461    {
462      detector->SetTreeAddress();
463    }
464  
465   return fRunLoader->GetHeader()->GetNtrack();
466 }
467
468 //_______________________________________________________________________
469 TGeometry *AliRun::GetGeometry()
470 {
471   // Create the TNode geometry for the event display
472   if (!fGeometry) { 
473     BuildSimpleGeometry();
474     //
475     // Unlink and relink nodes in detectors
476     // This is bad and there must be a better way...
477     //
478   
479     TIter next(fModules);
480     AliModule *detector;
481     while((detector = dynamic_cast<AliModule*>(next()))) {
482       detector->BuildGeometry();
483       TList *dnodes=detector->Nodes();
484       Int_t j;
485       TNode *node, *node1;
486       for ( j=0; j<dnodes->GetSize(); j++) {
487         node = dynamic_cast<TNode*>(dnodes->At(j));
488         node1 = fGeometry->GetNode(node->GetName());
489         dnodes->Remove(node);
490         dnodes->AddAt(node1,j);
491       }
492     }
493   }
494   return fGeometry;
495 }
496
497 //_______________________________________________________________________
498 void AliRun::SetBaseFile(const char *filename)
499 {
500   fBaseFileName = filename;
501 }
502
503 //_______________________________________________________________________
504 void AliRun::ResetDigits()
505 {
506   //
507   //  Reset all Detectors digits
508   //
509   TIter next(fModules);
510   AliModule *detector;
511   while((detector = dynamic_cast<AliModule*>(next()))) {
512      detector->ResetDigits();
513   }
514 }
515
516 //_______________________________________________________________________
517 void AliRun::ResetSDigits()
518 {
519   //
520   //  Reset all Detectors digits
521   //
522   TIter next(fModules);
523   AliModule *detector;
524   while((detector = dynamic_cast<AliModule*>(next()))) {
525      detector->ResetSDigits();
526   }
527 }
528
529
530 //_______________________________________________________________________
531
532 void AliRun::ResetPoints()
533 {
534   //
535   // Reset all Detectors points
536   //
537   TIter next(fModules);
538   AliModule *detector;
539   while((detector = dynamic_cast<AliModule*>(next()))) {
540      detector->ResetPoints();
541   }
542 }
543 //_______________________________________________________________________
544
545 void AliRun::InitMC(const char *setup)
546 {
547   //
548   // Initialize ALICE Simulation run
549   //
550   Announce();
551
552   if(fInitDone) {
553     AliWarning("Cannot initialise AliRun twice!");
554     return;
555   }
556     
557   if (!fMCApp)  
558     fMCApp=new AliMC(GetName(),GetTitle());
559     
560   gROOT->LoadMacro(setup);
561   gInterpreter->ProcessLine(fConfigFunction.Data());
562
563   fRunLoader->CdGAFile();
564
565   AliPDG::AddParticlesToPdgDataBase();  
566
567   fNdets = fModules->GetLast()+1;
568
569   TIter next(fModules);
570   for(Int_t i=0; i<fNdets; ++i)
571    {
572      TObject *objfirst, *objlast;
573      AliModule *detector=dynamic_cast<AliModule*>(fModules->At(i));
574      objlast = gDirectory->GetList()->Last();
575       
576      // Add Detector histograms in Detector list of histograms
577      if (objlast) objfirst = gDirectory->GetList()->After(objlast);
578      else         objfirst = gDirectory->GetList()->First();
579      while (objfirst) 
580       {
581         detector->Histograms()->Add(objfirst);
582         objfirst = gDirectory->GetList()->After(objfirst);
583       }
584    }
585    
586    fMCApp->Init();
587    
588    //Must be here because some MCs (G4) adds detectors here and not in Config.C
589    InitLoaders();
590    fRunLoader->MakeTree("E");
591    if (fLego == 0x0)
592     {
593       fRunLoader->LoadKinematics("RECREATE");
594       fRunLoader->LoadTrackRefs("RECREATE");
595       fRunLoader->LoadHits("all","RECREATE");
596     }
597    fInitDone = kTRUE;
598    //
599    // Save stuff at the beginning of the file to avoid file corruption
600    fRunLoader->CdGAFile();
601    Write();
602    fEventNrInRun = -1; //important - we start Begin event from increasing current number in run
603 }
604
605 //_______________________________________________________________________
606
607 void AliRun::RunMC(Int_t nevent, const char *setup)
608 {
609   //
610   // Main function to be called to process a galice run
611   // example
612   //    Root > gAlice.Run(); 
613   // a positive number of events will cause the finish routine
614   // to be called
615   //
616   fEventsPerRun = nevent;
617   // check if initialisation has been done
618   if (!fInitDone) InitMC(setup);
619   
620   // Create the Root Tree with one branch per detector
621   //Hits moved to begin event -> now we are crating separate tree for each event
622
623   gMC->ProcessRun(nevent);
624
625   // End of this run, close files
626   if(nevent>0) FinishRun();
627 }
628
629 //_______________________________________________________________________
630 void AliRun::RunReco(const char *selected, Int_t first, Int_t last)
631 {
632   //
633   // Main function to be called to reconstruct Alice event
634   // 
635    Int_t nev = fRunLoader->GetNumberOfEvents();
636    AliDebug(1, Form("Found %d events", nev));
637    Int_t nFirst = first;
638    Int_t nLast  = (last < 0)? nev : last;
639    
640    for (Int_t nevent = nFirst; nevent <= nLast; nevent++) {
641      AliDebug(1, Form("Processing event %d", nevent));
642      GetEvent(nevent);
643      Digits2Reco(selected);
644    }
645 }
646
647 //_______________________________________________________________________
648
649 void AliRun::Hits2Digits(const char *selected)
650 {
651
652    // Convert Hits to sumable digits
653    // 
654    for (Int_t nevent=0; nevent<gAlice->TreeE()->GetEntries(); nevent++) {
655      GetEvent(nevent);
656      Hits2SDigits(selected);
657      SDigits2Digits(selected);
658    }  
659 }
660
661
662 //_______________________________________________________________________
663
664 void AliRun::Tree2Tree(Option_t *option, const char *selected)
665 {
666   //
667   // Function to transform the content of
668   //  
669   // - TreeH to TreeS (option "S")
670   // - TreeS to TreeD (option "D")
671   // - TreeD to TreeR (option "R")
672   // 
673   // If multiple options are specified ("SDR"), transformation will be done in sequence for
674   // selected detector and for all detectors if none is selected (detector string 
675   // can contain blank separated list of detector names). 
676
677
678    const char *oS = strstr(option,"S");
679    const char *oD = strstr(option,"D");
680    const char *oR = strstr(option,"R");
681    
682    TObjArray *detectors = Detectors();
683
684    TIter next(detectors);
685
686    AliDetector *detector = 0;
687
688    while((detector = dynamic_cast<AliDetector*>(next()))) {
689      if (selected) 
690        if (strcmp(detector->GetName(),selected)) continue;
691      if (detector->IsActive())
692       { 
693        
694        AliLoader* loader = detector->GetLoader();
695        if (loader == 0x0) continue;
696        
697        if (oS) 
698         {
699           AliDebug(1, Form("Processing Hits2SDigits for %s ...", detector->GetName()));
700           loader->LoadHits("read");
701           if (loader->TreeS() == 0x0) loader->MakeTree("S");
702           detector->MakeBranch(option);
703           detector->SetTreeAddress();
704           detector->Hits2SDigits();
705           loader->UnloadHits();
706           loader->UnloadSDigits();
707         }  
708        if (oD) 
709         {
710           AliDebug(1, Form("Processing SDigits2Digits for %s ...", detector->GetName()));
711           loader->LoadSDigits("read");
712           if (loader->TreeD() == 0x0) loader->MakeTree("D");
713           detector->MakeBranch(option);
714           detector->SetTreeAddress();
715           detector->SDigits2Digits();
716           loader->UnloadSDigits();
717           loader->UnloadDigits();
718         } 
719        if (oR) 
720         {
721           AliDebug(1, Form("Processing Digits2Reco for %s ...", detector->GetName()));
722           loader->LoadDigits("read");
723           if (loader->TreeR() == 0x0) loader->MakeTree("R");
724           detector->MakeBranch(option);
725           detector->SetTreeAddress();
726           detector->Digits2Reco(); 
727           loader->UnloadDigits();
728           loader->UnloadRecPoints();
729
730         }
731      }   
732    }
733 }
734
735 //_______________________________________________________________________
736 void AliRun::RunLego(const char *setup, Int_t nc1, Float_t c1min,
737                      Float_t c1max,Int_t nc2,Float_t c2min,Float_t c2max,
738                      Float_t rmin,Float_t rmax,Float_t zmax, AliLegoGenerator* gener)
739 {
740   //
741   // Generates lego plots of:
742   //    - radiation length map phi vs theta
743   //    - radiation length map phi vs eta
744   //    - interaction length map
745   //    - g/cm2 length map
746   //
747   //  ntheta    bins in theta, eta
748   //  themin    minimum angle in theta (degrees)
749   //  themax    maximum angle in theta (degrees)
750   //  nphi      bins in phi
751   //  phimin    minimum angle in phi (degrees)
752   //  phimax    maximum angle in phi (degrees)
753   //  rmin      minimum radius
754   //  rmax      maximum radius
755   //  
756   //
757   //  The number of events generated = ntheta*nphi
758   //  run input parameters in macro setup (default="Config.C")
759   //
760   //  Use macro "lego.C" to visualize the 3 lego plots in spherical coordinates
761   //Begin_Html
762   /*
763     <img src="picts/AliRunLego1.gif">
764   */
765   //End_Html
766   //Begin_Html
767   /*
768     <img src="picts/AliRunLego2.gif">
769   */
770   //End_Html
771   //Begin_Html
772   /*
773     <img src="picts/AliRunLego3.gif">
774   */
775   //End_Html
776   //
777
778   // check if initialisation has been done
779   // If runloader has been initialized, set the number of events per file to nc1 * nc2
780
781   // Set new generator
782   if (!gener) gener  = new AliLegoGenerator();
783   //
784   // Configure Generator
785   gener->SetRadiusRange(rmin, rmax);
786   gener->SetZMax(zmax);
787   gener->SetCoor1Range(nc1, c1min, c1max);
788   gener->SetCoor2Range(nc2, c2min, c2max);
789   
790   
791   //Create Lego object  
792   fLego = new AliLego("lego",gener);
793
794   if (!fInitDone) InitMC(setup);
795   //Save current generator
796   
797   AliGenerator *gen=fMCApp->Generator();
798   fMCApp->ResetGenerator(gener);
799   //Prepare MC for Lego Run
800   gMC->InitLego();
801   
802   //Run Lego Object
803
804   if (fRunLoader) fRunLoader->SetNumberOfEventsPerFile(nc1 * nc2);
805   //gMC->ProcessRun(nc1*nc2+1);
806   gMC->ProcessRun(nc1*nc2);
807   
808   // End of this run, close files
809   FinishRun();
810   // Restore current generator
811   fMCApp->ResetGenerator(gen);
812   // Delete Lego Object
813   delete fLego; fLego=0;
814 }
815
816 //_______________________________________________________________________
817 void AliRun::SetConfigFunction(const char * config) 
818 {
819   //
820   // Set the signature of the function contained in Config.C to configure
821   // the run
822   //
823   fConfigFunction=config;
824 }
825
826 // 
827 // MC Application
828 // 
829
830 //_______________________________________________________________________
831 void AliRun::Field(const Double_t* x, Double_t *b) const
832 {
833   //
834   // Return the value of the magnetic field
835   //
836   Float_t xfloat[3];
837   for (Int_t i=0; i<3; i++) xfloat[i] = x[i]; 
838   
839   if (Field()) {
840     Float_t bfloat[3];
841     Field()->Field(xfloat,bfloat);
842     for (Int_t j=0; j<3; j++) b[j] = bfloat[j]; 
843   } 
844   else {
845     AliError("No mag field defined!");
846     b[0]=b[1]=b[2]=0.;
847   }
848 }      
849
850 // 
851 // End of MC Application
852 // 
853
854 //_______________________________________________________________________
855 void AliRun::Streamer(TBuffer &R__b)
856 {
857   // Stream an object of class AliRun.
858
859   if (R__b.IsReading()) {
860     if (!gAlice) gAlice = this;
861     AliRun::Class()->ReadBuffer(R__b, this);
862     gROOT->GetListOfBrowsables()->Add(this,"Run");
863
864     gRandom = fRandom;
865   } else {
866     AliRun::Class()->WriteBuffer(R__b, this);
867   }
868 }
869 //_______________________________________________________________________
870
871 void AliRun::SetGenEventHeader(AliGenEventHeader* header)
872 {
873   fRunLoader->GetHeader()->SetGenEventHeader(header);
874 }
875 //_______________________________________________________________________
876
877 Int_t AliRun::GetEvNumber() const
878
879 //Returns number of current event  
880   if (fRunLoader == 0x0)
881    {
882      AliError("RunLoader is not set. Can not load data.");
883      return -1;
884    }
885
886   return fRunLoader->GetEventNumber();
887 }
888 //_______________________________________________________________________
889
890 void AliRun::SetRunLoader(AliRunLoader* rloader)
891 {
892   //
893   // Set the loader of the run
894   //
895   fRunLoader = rloader;
896   if (fRunLoader == 0x0) return;
897   
898   TString evfoldname;
899   TFolder* evfold = fRunLoader->GetEventFolder();
900   if (evfold) evfoldname = evfold->GetName();
901   else AliWarning("Did not get Event Folder from Run Loader");
902   
903   if ( fRunLoader->GetAliRun() )
904    {//if alrun already exists in folder
905     if (fRunLoader->GetAliRun() != this )
906      {//and is different than this - crash
907        AliFatal("AliRun is already in Folder and it is not this object");
908        return;//pro forma
909      }//else do nothing
910    }
911   else
912    {
913      evfold->Add(this);//Post this AliRun to Folder
914    }
915   
916   TIter next(fModules);
917   AliModule *module;
918   while((module = (AliModule*)next())) 
919    {
920      if (evfold) AliConfig::Instance()->Add(module,evfoldname);
921      module->SetRunLoader(fRunLoader);
922      AliDetector* detector = dynamic_cast<AliDetector*>(module);
923      if (detector)
924       {
925         AliLoader* loader = fRunLoader->GetLoader(detector);
926         if (loader == 0x0)
927          {
928            AliError(Form("Can not get loader for detector %s", detector->GetName()));
929          }
930         else
931          {
932            AliDebug(1, Form("Setting loader for detector %s", detector->GetName()));
933            detector->SetLoader(loader);
934          }
935       }
936    }
937 }
938
939 void AliRun::AddModule(AliModule* mod)
940 {
941   //
942   // Add a module to the module list
943   //
944   if (mod == 0x0) return;
945   if (strlen(mod->GetName()) == 0) return;
946   if (GetModuleID(mod->GetName()) >= 0) return;
947   
948   AliDebug(1, mod->GetName());
949   if (fRunLoader == 0x0) AliConfig::Instance()->Add(mod);
950   else AliConfig::Instance()->Add(mod,fRunLoader->GetEventFolder()->GetName());
951
952   Modules()->Add(mod);
953   
954   fNdets++;
955 }
956
957 // added by Alberto Colla
958 //_____________________________________________________________________________
959 /*inline*/ Bool_t AliRun::IsFileAccessible(const char* fnam, EAccessMode mode)
960 { return !gSystem->AccessPathName(fnam,mode);}
961
962 //______________________________________________________
963 /*inline*/ Bool_t AliRun::IsFileAccessible(Char_t* name,EAccessMode mode)
964 {
965   TString str = name; gSystem->ExpandPathName(str);
966   return !gSystem->AccessPathName(str.Data(),mode);
967 }