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