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