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