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