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