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