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