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