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