]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliRunLoader.cxx
c3e7f732dca3c41126c8dd622d2fe3030d24391e
[u/mrichter/AliRoot.git] / STEER / AliRunLoader.cxx
1 #include "AliRunLoader.h"
2 //_____________________________________
3 /////////////////////////////////////////////////////////////////////////////////////
4 //
5 // class AliRunLoader
6 //
7 //This class aims to be the only one interface for manging data
8 //It stores Loaders for all modules which knows the filenames 
9 //of the data files to be stored.
10 //It aims to substitude AliRun in automatic maging of data positioning
11 //thus there won't be necessity of loading gAlice from file in order to 
12 //get fast ccess to the data
13 //
14 //logical place for putting the Loader specific to the given detector is detector itself
15 // but, to load detector one need to load gAlice, and by the way all other detectors
16 // with their geometrieces and so on. 
17 // So, if one need to open TPC clusters there is no principal need to read everything
18 //
19 // When RunLoader is read from the file it does not connect to the folder structure
20 // it must be connected (mounted) manualy in the macro or class code. 
21 // Default event folder is defined by AliConfig::fgkDefaultEventFolderName
22 // but can be mounted elsewhere. Usefull specially in merging, when more than session 
23 // needs to be loaded
24 //
25 //////////////////////////////////////////////////////////////////////////////////////
26 /**************************************************************************/
27
28 #include <TString.h>
29 #include <TFolder.h>
30 #include <TFile.h>
31 #include <TTree.h>
32 #include <TROOT.h>
33 #include <TBranch.h>
34 #include <TGeometry.h>
35 #include <TTask.h>
36 #include <TError.h>
37
38 #include "AliRun.h"
39 #include "AliConfig.h"
40 #include "AliLoader.h"
41 #include "AliHeader.h"
42 #include "AliStack.h"
43 #include "TObjArray.h"
44 #include "AliDetector.h"
45 #include "AliRunDigitizer.h"
46
47 ClassImp(AliRunLoader)
48
49 AliRunLoader* AliRunLoader::fgRunLoader = 0x0;
50
51 const TString AliRunLoader::fgkRunLoaderName("RunLoader");
52
53 const TString AliRunLoader::fgkHeaderBranchName("Header");
54 const TString AliRunLoader::fgkHeaderContainerName("TE");
55 const TString AliRunLoader::fgkKineContainerName("TreeK");
56 const TString AliRunLoader::fgkTrackRefsContainerName("TreeTR");
57 const TString AliRunLoader::fgkKineBranchName("Particles");
58 const TString AliRunLoader::fgkDefaultKineFileName("Kinematics.root");
59 const TString AliRunLoader::fgkDefaultTrackRefsFileName("TrackRefs.root");
60 const TString AliRunLoader::fgkGAliceName("gAlice");
61 /**************************************************************************/
62
63 AliRunLoader::AliRunLoader():
64  fLoaders(0x0),
65  fEventFolder(0x0),
66  fCurrentEvent(0),
67  fGAFile(0x0),
68  fHeader(0x0),
69  fStack(0x0),
70  fKineDataLoader(0x0),
71  fTrackRefsDataLoader(0x0),
72  fNEventsPerFile(1),
73  fUnixDirName(".")
74 {
75   AliConfig::Instance();//force to build the folder structure
76 }
77 /**************************************************************************/
78
79 AliRunLoader::AliRunLoader(const char* eventfoldername):
80  TNamed(fgkRunLoaderName,fgkRunLoaderName),
81  fLoaders(new TObjArray()),
82  fEventFolder(0x0),
83  fCurrentEvent(0),
84  fGAFile(0x0),
85  fHeader(0x0),
86  fStack(0x0),
87  fKineDataLoader(new AliDataLoader(fgkDefaultKineFileName,fgkKineContainerName,"Kinematics")),
88  fTrackRefsDataLoader(new AliDataLoader(fgkDefaultTrackRefsFileName,fgkTrackRefsContainerName,"Track References")),
89  fNEventsPerFile(1),
90  fUnixDirName(".")
91 {
92 //ctor
93   SetEventFolderName(eventfoldername);
94 }
95 /**************************************************************************/
96
97 AliRunLoader::~AliRunLoader()
98 {
99
100   UnloadHeader();
101   UnloadgAlice();
102   
103   if(fLoaders) {
104     fLoaders->SetOwner();
105     delete fLoaders;
106   }
107   
108   delete fKineDataLoader;
109   delete fTrackRefsDataLoader;
110   
111   
112   RemoveEventFolder();
113   
114   //fEventFolder is deleted by the way of removing - TopAliceFolder owns it
115   delete fHeader;
116   delete fStack;
117   delete fGAFile;
118 }
119 /**************************************************************************/
120
121 AliRunLoader::AliRunLoader(TFolder* topfolder):TNamed(fgkRunLoaderName,fgkRunLoaderName)
122 {
123  if(topfolder == 0x0)
124   {
125     Fatal("AliRunLoader(TFolder*)","Parameter is NULL");
126     return;
127   }
128  fEventFolder = topfolder;
129  
130  TObject* obj = fEventFolder->FindObject(fgkRunLoaderName);
131  if (obj)
132   { //if it is, then sth. is going wrong... exits aliroot session
133     Fatal("AliRunLoader(const char*)",
134           "In Event Folder Named %s object named %s already exists. I am confused ...",
135            fEventFolder->GetName(),fgkRunLoaderName.Data());
136     return;//never reached
137   }
138    
139  fLoaders = new TObjArray();
140  fEventFolder->Add(this);//put myself to the folder to accessible for all
141   
142 }
143 /**************************************************************************/
144
145 Int_t AliRunLoader::GetEvent(Int_t evno)
146 {
147 //Gets event number evno
148 //Reloads all data properly
149   if (fCurrentEvent == evno) return 0;
150   
151   if (evno < 0)
152    {
153      Error("GetEvent","Can not give the event with negative number");
154      return 4;
155    }
156
157   if (evno >= GetNumberOfEvents())
158    {
159      Error("GetEvent","There is no event with number %d",evno);
160      return 3;
161    }
162   
163   if (GetDebug()) 
164    {
165      Info("GetEvent",">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
166      Info("GetEvent",">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
167      Info("GetEvent",">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
168      Info("GetEvent","          GETTING EVENT  %d",evno);
169      Info("GetEvent",">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
170      Info("GetEvent",">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
171      Info("GetEvent",">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
172    }
173    
174   fCurrentEvent = evno;
175
176   Int_t retval;
177   
178   //Reload header (If header was loaded)
179   if (GetHeader())
180    {
181      retval = TreeE()->GetEvent(fCurrentEvent);
182      if ( retval == 0)
183       {
184         Error("GetEvent","Cannot find event: %d\n ",fCurrentEvent);
185         return 5;
186       }
187    }
188   //Reload stack (If header was loaded)
189   if (TreeE()) fStack = GetHeader()->Stack();
190   //Set event folder in stack (it does not mean that we read kinematics from file)
191   if (fStack) 
192    { 
193      fStack->SetEventFolderName(fEventFolder->GetName());
194    }
195   else
196    {
197      Warning("GetEvent","Stack not found in header");
198    }
199   
200   retval = SetEvent();
201   if (retval)
202    {
203      Error("GetEvent","Error occured while setting event %d",evno);
204      return 1;
205    }
206    
207   //Post Track References
208   retval = fTrackRefsDataLoader->GetEvent();
209   if (retval)
210    {
211      Error("GetEvent","Error occured while GetEvent for Track References. Event %d",evno);
212      return 2;
213    }
214
215   //Read Kinematics if loaded
216   fKineDataLoader->GetEvent();
217   if (retval)
218    {
219      Error("GetEvent","Error occured while GetEvent for Kinematics. Event %d",evno);
220      return 2;
221    }
222
223   if (fStack && fKineDataLoader->GetBaseLoader(0)->IsLoaded()) fStack->GetEvent();
224   
225   //Trigger data reloading in all loaders 
226   TIter next(fLoaders);
227   AliLoader *loader;
228   while((loader = (AliLoader*)next())) 
229    {
230      retval = loader->GetEvent();
231      if (retval)
232       {
233        Error("GetEvent","Error occured while getting event for %s. Event %d.",
234               loader->GetDetectorName().Data(), evno);
235        return 3;
236       }
237    }
238   
239   SetDetectorAddresses();
240   
241   return 0;
242 }
243 /**************************************************************************/
244 Int_t AliRunLoader::SetEvent()
245 {
246  //if kinematocs was loaded Cleans folder data
247  //change 
248   Int_t retval;
249   
250   retval = fKineDataLoader->SetEvent();
251   if (retval)
252    {
253      Error("SetEvent","SetEvent for Kinamtics Data Loader retutned error.");
254      return retval;
255    }
256   retval = fTrackRefsDataLoader->SetEvent(); 
257   if (retval)
258    {
259      Error("SetEvent","SetEvent for Track References Data Loader retutned error.");
260      return retval;
261    }
262
263   TIter next(fLoaders);
264   AliLoader *loader;
265   while((loader = (AliLoader*)next())) 
266    {
267      retval = loader->SetEvent();
268      if (retval)
269       {
270         Error("SetEvent","SetEvent for %s Data Loader retutned error.",loader->GetName());
271         return retval;
272       }
273    }
274
275   return 0;
276 }
277 /**************************************************************************/
278
279 Int_t AliRunLoader::SetEventNumber(Int_t evno)
280 {
281   //cleans folders and sets the root dirs in files 
282   if (fCurrentEvent == evno) return 0;
283   fCurrentEvent = evno;
284   return SetEvent();
285 }
286
287 /**************************************************************************/
288 AliRunLoader* AliRunLoader::Open
289   (const char* filename, const char* eventfoldername, Option_t* option)
290 {
291 //Opens a desired file 'filename'
292 //gets the the run-Loader and mounts it desired folder
293 //returns the pointer to run Loader which can be further used for accessing data 
294 //in case of error returns NULL
295  
296  static const TString webaddress("http://alisoft.cern.ch/people/skowron/codedoc/split/index.html");
297  if (AliLoader::fgDebug) 
298   ::Info("AliRunLoader::Open",
299          "\n\n\nNew I/O strcture: See more info:\n %s\n\n\n",webaddress.Data());
300  
301  AliRunLoader* result = 0x0;
302  
303  /* ************************************************ */
304  /* Chceck if folder with given name already exists  */
305  /* ************************************************ */
306  
307  TObject* obj = AliConfig::Instance()->GetTopFolder()->FindObject(eventfoldername);
308  if(obj)
309   {
310     TFolder* fold = dynamic_cast<TFolder*>(obj);
311     if (fold == 0x0)
312      {
313       ::Error("AliRunLoader::Open","Such a obejct already exists in top alice folder and it is not a folder.");
314       return 0x0;
315      }
316     
317     //check if we can get RL from that folder
318     result = AliRunLoader::GetRunLoader(eventfoldername);
319     if (result == 0x0)
320      {
321        ::Error("AliRunLoader::Open",
322                "Folder %s already exists, and can not find session there. Can not mount.",eventfoldername);
323        return 0x0;
324      }
325
326     if (result->GetFileName().CompareTo(filename) != 0)
327      {
328        ::Error("AliRunLoader::Open","Other file is mounted in demanded folder. Can not mount.");
329        return 0x0;
330      }
331
332     //check if now is demanded (re)creation 
333     if ( AliLoader::TestFileOption(option) == kFALSE)
334      {
335        ::Error("AliRunLoader::Open",
336                "Session already exists in folder %s and this session option is %s. Unable to proceed.",
337                 eventfoldername,option);
338        return 0x0;
339      }
340      
341     //check if demanded option is update and existing one 
342     TString tmpstr(option);
343     if ( (tmpstr.CompareTo("update",TString::kIgnoreCase) == 0) && 
344          (result->fGAFile->IsWritable() == kFALSE) )
345      { 
346        ::Error("AliRunLoader::Open",
347                "Session already exists in folder %s and is not writable while this session option is %s. Unable to proceed.",
348                 eventfoldername,option);
349        return 0x0;
350      }
351      
352     ::Warning("AliRunLoader::Open","Session is already opened and mounted in demanded folder"); 
353     return result;
354   } //end of checking in case of existance of object named identically that folder session is being opened
355  
356  
357  TFile * gAliceFile = TFile::Open(filename,option);//open a file
358  if (!gAliceFile) 
359   {//null pointer returned
360     ::Error("AliRunLoader::Open","Can not open file %s.",filename);
361     return 0x0;
362   }
363   
364  if (gAliceFile->IsOpen() == kFALSE)
365   {//pointer to valid object returned but file is not opened
366     ::Error("AliRunLoader::Open","Can not open file %s.",filename);
367     return 0x0;
368   }
369  
370  //if file is "read" or "update" than we try to find AliRunLoader there - if not found cry and exit
371  //else create new AliRunLoader
372  if ( AliLoader::TestFileOption(option) )
373   { 
374     if (AliLoader::fgDebug) 
375      ::Info("AliRunLoader::Open","Reading RL from file");
376     
377     result = dynamic_cast<AliRunLoader*>(gAliceFile->Get(fgkRunLoaderName));//get the run Loader from the file
378     if (result == 0x0)
379      {//didn't get
380        ::Error("AliRunLoader::Open","Can not find run-Loader in file %s.",filename);
381        delete gAliceFile;//close the file
382        return 0x0;
383      }
384     Int_t tmp = result->SetEventFolderName(eventfoldername);//mount a event folder   
385     if (tmp)//if SetEvent  returned error
386      {
387        ::Error("AliRunLoader::Open","Can not mount event in folder %s.",eventfoldername);
388        delete result; //delete run-Loader
389        delete gAliceFile;//close the file
390        return 0x0;
391      }
392   }
393  else
394   {
395     if (AliLoader::fgDebug) 
396       ::Info("AliRunLoader::Open","Creating new AliRunLoader. Folder name is %s",eventfoldername);
397     result = new AliRunLoader(eventfoldername);
398   }
399  
400 //procedure for extracting dir name from the file name 
401  TString fname(filename);
402  Int_t  nsl = fname.Last('/');//look for slash in file name
403  TString dirname;
404  if (nsl < 0) 
405   {//slash not found
406     Int_t  nsl = fname.Last(':');//look for colon e.g. rfio:galice.root
407     if (nsl < 0) dirname = ".";//not found
408     else dirname = fname.Remove(nsl);//found
409   }
410  else dirname = fname.Remove(nsl);//slash found
411  
412  if (AliLoader::fgDebug) 
413   ::Info("AliRunLoader::Open","Dir name is : %s",dirname.Data());
414  
415  result->SetDirName(dirname); 
416  result->SetGAliceFile(gAliceFile);//set the pointer to gAliceFile
417  fgRunLoader = result; //PH get access from any place
418  return result;
419 }
420 /**************************************************************************/
421 Int_t AliRunLoader::GetNumberOfEvents()
422 {
423  //returns number of events in Run
424  Int_t retval;
425  if( TreeE() == 0x0 )
426   {
427     retval = LoadHeader();
428     if (retval) 
429      {
430        Error("GetNumberOfEvents","Error occured while loading header");
431        return -1;
432      }
433   }
434  return (Int_t)TreeE()->GetEntries();
435 }
436
437 /**************************************************************************/
438 void AliRunLoader::MakeHeader()
439 {
440  //Makes header and connects it to header tree (if it exists)
441   if (GetDebug()) Info("MakeHeader","");
442   if(fHeader == 0x0)
443    {
444      if (GetDebug()) Info("MakeHeader","Creating new Header Object");
445      fHeader= new AliHeader();
446    }
447   TTree* tree = TreeE();
448   if (tree)
449    {
450      if (GetDebug()) Info("MakeHeader","Got Tree from folder.");
451      TBranch* branch = tree->GetBranch(fgkHeaderBranchName);
452      if (branch == 0x0)
453       {
454         if (GetDebug()) Info("MakeHeader","Creating new branch");
455         branch = tree->Branch(fgkHeaderBranchName, "AliHeader", &fHeader, 4000, 0);
456         branch->SetAutoDelete(kFALSE);
457       }
458      else
459       {
460         if (GetDebug()) Info("MakeHeader","Got Branch from Tree");
461         branch->SetAddress(&fHeader);
462         tree->GetEvent(fCurrentEvent);
463         fStack = fHeader->Stack(); //should be safe - if we created Stack, header returns pointer to the same object
464         if (fStack)
465          {
466            fStack->SetEventFolderName(fEventFolder->GetName());
467            if (TreeK()) fStack->GetEvent();
468          }
469         else
470         {
471           if (GetDebug()) Info("MakeHeader","Haeder do not have a stack.");
472         }
473       }
474    } 
475   if (GetDebug()) Info("MakeHeader","Exiting MakeHeader method");
476 }
477 /**************************************************************************/
478
479 void AliRunLoader::MakeStack()
480 {
481 //Creates the stack object -  do not connect the tree
482   if(fStack == 0x0)
483    { 
484      fStack = new AliStack(10000);
485      fStack->SetEventFolderName(fEventFolder->GetName());
486    }
487 }
488
489 /**************************************************************************/
490
491 void AliRunLoader::MakeTree(Option_t *option)
492 {
493 //Creates trees
494   const char *oK = strstr(option,"K"); //Kine  
495   const char *oE = strstr(option,"E"); //Header
496
497   if(oK && !TreeK())
498    { 
499      if (fKineDataLoader->GetBaseLoader(0)->IsLoaded() == kFALSE)
500       {
501         Error("MakeTree(\"K\")","Load Kinematics first");
502       }
503      else
504       {
505         fKineDataLoader->MakeTree();
506         MakeStack();
507         fStack->ConnectTree();
508         WriteKinematics("OVERWRITE");
509      }
510    }
511   
512   if(oE && !TreeE())
513    { 
514      fGAFile->cd();
515      TTree* tree = new TTree(fgkHeaderContainerName,"Tree with Headers");
516      GetEventFolder()->Add(tree);
517      MakeHeader();
518      WriteHeader("OVERWRITE");
519    }
520   
521   TIter next(fLoaders);
522   AliLoader *loader;
523   while((loader = (AliLoader*)next()))
524    {
525     loader->MakeTree(option);
526    }
527
528 }
529 /**************************************************************************/
530     
531 Int_t AliRunLoader::LoadgAlice()
532 {
533 //Loads gAlice from file
534  if (GetAliRun())
535   {
536     Warning("LoadgAlice","AliRun is already in folder. Unload first.");
537     return 0;
538   }
539  AliRun* alirun = dynamic_cast<AliRun*>(fGAFile->Get(fgkGAliceName));
540  if (alirun == 0x0)
541   {
542     Error("LoadgAlice"," Can not find gAlice in file %s",fGAFile->GetName());
543     return 2;
544   }
545  alirun->SetRunLoader(this);
546  if (gAlice)
547   {
548     Warning("LoadgAlice","gAlice already exists. Putting retrived object in folder named %s",
549              GetEventFolder()->GetName());
550   }
551  else
552   {
553     gAlice = alirun;
554   }
555  SetDetectorAddresses();//calls SetTreeAddress for all detectors
556  return 0; 
557 }
558 /**************************************************************************/
559
560 Int_t AliRunLoader::LoadHeader()
561 {
562  if (TreeE())
563   {
564      Warning("LoadHeader","Header is already loaded. Use ReloadHeader to force reload. Nothing done");
565      return 0;
566   }
567  
568  if (GetEventFolder() == 0x0)
569   {
570     Error("LoadHaeder","Event folder not specified yet");
571     return 1;
572   }
573
574  if (fGAFile == 0x0)
575   {
576     Error("LoadHaeder","Session not opened. Use AliRunLoader::Open");
577     return 2;
578   }
579  
580  if (fGAFile->IsOpen() == kFALSE)
581   {
582     Error("LoadHaeder","Session not opened. Use AliRunLoader::Open");
583     return 2;
584   }
585
586  TTree* tree = dynamic_cast<TTree*>(fGAFile->Get(fgkHeaderContainerName));
587  if (tree == 0x0)
588   {
589     Fatal("LoadHaeder","Can not find header tree named %s in file %s",
590            fgkHeaderContainerName.Data(),fGAFile->GetName());
591     return 2;
592   }
593
594  if (tree == TreeE()) return 0;
595
596  CleanHeader();
597  GetEventFolder()->Add(tree);
598  MakeHeader();//creates header object and connects to tree
599  return 0; 
600
601 }
602 /**************************************************************************/
603
604 Int_t AliRunLoader::LoadKinematics(Option_t* option)
605 {
606 //Loads the kinematics 
607  Int_t retval = fKineDataLoader->GetBaseLoader(0)->Load(option);
608  if (retval)
609   {
610     Error("LoadKinematics","Error occured while loading kinamatics tree.");
611     return retval;
612   }
613  if (fStack) fStack->GetEvent();
614  return 0;
615 }
616 /**************************************************************************/
617
618 Int_t AliRunLoader::OpenDataFile(const TString& filename,TFile*& file,TDirectory*& dir,Option_t* opt,Int_t cl)
619 {
620 //Opens File with kinematics
621  if (file)
622   {
623     if (file->IsOpen() == kFALSE)
624      {//pointer is not null but file is not opened
625        Warning("OpenDataFile","Pointer to file is not null, but file is not opened");//risky any way
626        delete file;
627        file = 0x0; //proceed with opening procedure
628      }
629     else
630      { 
631        Warning("OpenDataFile","File  %s already opened",filename.Data());
632        return 0;
633      }
634   }
635 //try to find if that file is opened somewere else
636  file = (TFile *)( gROOT->GetListOfFiles()->FindObject(filename) );
637  if (file)
638   {
639    if(file->IsOpen() == kTRUE)
640     {
641      Warning("OpenDataFile","File %s already opened by sombody else.",file->GetName());
642      return 0;
643     }
644   }
645
646  file = TFile::Open(filename,opt);
647  if (file == 0x0)
648   {//file is null
649     Error("LoadKinematics","Can not open file %s",filename.Data());
650     return 1;
651   }
652  if (file->IsOpen() == kFALSE)
653   {//file is not opened
654    Error("LoadKinematics","Can not open file %s",filename.Data());
655    return 1;
656   }
657   
658  file->SetCompressionLevel(cl);
659  
660  dir = AliLoader::ChangeDir(file,fCurrentEvent);
661  if (dir == 0x0)
662   {
663     Error("OpenKineFile","Can not change to root directory in file %s",filename.Data());
664     return 3;
665   }
666  return 0; 
667 }
668 /**************************************************************************/
669
670 TTree* AliRunLoader::TreeE() const
671 {
672  //returns the tree from folder; shortcut method
673  TObject *obj = fEventFolder->FindObject(fgkHeaderContainerName);
674  return (obj)?dynamic_cast<TTree*>(obj):0x0;
675 }
676 /**************************************************************************/
677
678 AliHeader* AliRunLoader::GetHeader() const
679 {
680  return fHeader;
681 }
682 /**************************************************************************/
683  
684 TTree* AliRunLoader::TreeK() const
685 {
686  //returns the tree from folder; shortcut method
687  TObject *obj = GetEventFolder()->FindObject(fgkKineContainerName);
688  return (obj)?dynamic_cast<TTree*>(obj):0x0;
689 }
690 /**************************************************************************/
691
692 TTree* AliRunLoader::TreeTR() const
693 {
694  //returns the tree from folder; shortcut method
695  TObject* obj = GetEventFolder()->FindObject(fgkTrackRefsContainerName);
696  return (obj)?dynamic_cast<TTree*>(obj):0x0;
697 }
698 /**************************************************************************/
699
700 AliRun* AliRunLoader::GetAliRun() const
701 {
702 //returns AliRun which sits in the folder
703  if (fEventFolder == 0x0) return 0x0;
704  TObject *obj = fEventFolder->FindObject(fgkGAliceName);
705  return (obj)?dynamic_cast<AliRun*>(obj):0x0;
706 }
707 /**************************************************************************/
708
709 Int_t AliRunLoader::WriteGeometry(Option_t* opt)
710 {
711   fGAFile->cd();
712   TGeometry* geo = GetAliRun()->GetGeometry();
713   if (geo == 0x0)
714    {
715      Error("WriteGeometry","Can not get geometry from gAlice");
716      return 1;
717    }
718   geo->Write();
719   return 0;
720 }
721 /**************************************************************************/
722
723 Int_t AliRunLoader::WriteHeader(Option_t* opt)
724 {
725   if (GetDebug()) Info("WriteHeader","  WRITING HEADER");
726   
727   TTree* tree = TreeE();
728   if ( tree == 0x0)
729    {
730      Warning("WriteHeader","Can not find Header Tree in Folder");
731      return 0;
732    } 
733   if (fGAFile->IsWritable() == kFALSE)
734    {
735      Error("WriteHeader","File %s is not writable",fGAFile->GetName());
736      return 1;
737    }
738
739   TObject* obj = fGAFile->Get(fgkHeaderContainerName);
740   if (obj)
741    { //if they exist, see if option OVERWRITE is used
742      TString tmp(opt);
743      if(tmp.Contains("OVERWRITE",TString::kIgnoreCase) == 0)
744       {//if it is not used -  give an error message and return an error code
745         Error("WriteHeader","Tree already exisists. Use option \"OVERWRITE\" to overwrite previous data");
746         return 3;
747       }
748    }
749   fGAFile->cd();
750   tree->SetDirectory(fGAFile);
751   tree->Write(0,TObject::kOverwrite);
752
753   if (GetDebug()) Info("WriteHeader","WRITTEN\n\n");
754   
755   return 0;
756 }
757 /**************************************************************************/
758
759 Int_t AliRunLoader::WriteAliRun(Option_t* opt)
760 {
761   fGAFile->cd();
762   if (gAlice) gAlice->Write();
763   return 0;
764 }
765 /**************************************************************************/
766
767 Int_t AliRunLoader::WriteKinematics(Option_t* opt)
768 {
769   return fKineDataLoader->GetBaseLoader(0)->WriteData(opt);
770 }
771 /**************************************************************************/
772 Int_t AliRunLoader::WriteTrackRefs(Option_t* opt)
773 {
774   return fTrackRefsDataLoader->GetBaseLoader(0)->WriteData(opt);
775 }
776 /**************************************************************************/
777
778 Int_t AliRunLoader::WriteHits(Option_t* opt)
779 {
780 //Calls WriteHits for all loaders
781   Int_t res;
782   Int_t result = 0;
783   TIter next(fLoaders);
784   AliLoader *loader;
785   while((loader = (AliLoader*)next()))
786    {
787      res = loader->WriteHits(opt);
788      if (res)
789       {
790         Error("WriteHits","Failed to write hits for %s (%d)",loader->GetDetectorName().Data(),res);
791         result = 1;
792       }
793    }
794   return result;
795 }
796 /**************************************************************************/
797
798 Int_t AliRunLoader::WriteSDigits(Option_t* opt)
799 {
800   Int_t res;
801   Int_t result = 0;
802   TIter next(fLoaders);
803   AliLoader *loader;
804   while((loader = (AliLoader*)next()))
805    {
806      res = loader->WriteSDigits(opt);
807      if (res)
808       {
809         Error("WriteSDigits","Failed to write summable digits for %s.",loader->GetDetectorName().Data());
810         result = 1;
811       }
812    }
813   return result;
814 }
815 /**************************************************************************/
816
817 Int_t AliRunLoader::WriteDigits(Option_t* opt)
818 {
819   Int_t res;
820   Int_t result = 0;
821   TIter next(fLoaders);
822   AliLoader *loader;
823   while((loader = (AliLoader*)next()))
824    { 
825      res = loader->WriteDigits(opt);
826      if (res)
827       {
828         Error("WriteDigits","Failed to write digits for %s.",loader->GetDetectorName().Data());
829         result = 1;
830       }
831    }
832   return result;
833 }
834 /**************************************************************************/
835
836 Int_t AliRunLoader::WriteRecPoints(Option_t* opt)
837 {
838   Int_t res;
839   Int_t result = 0;
840   TIter next(fLoaders);
841   AliLoader *loader;
842   while((loader = (AliLoader*)next()))
843    {
844      res = loader->WriteRecPoints(opt);
845      if (res)
846       {
847         Error("WriteRecPoints","Failed to write Reconstructed Points for %s.",
848               loader->GetDetectorName().Data());
849         result = 1;
850       }
851    }
852   return result;
853 }
854 /**************************************************************************/
855
856 Int_t AliRunLoader::WriteTracks(Option_t* opt)
857 {
858   Int_t res;
859   Int_t result = 0;
860   TIter next(fLoaders);
861   AliLoader *loader;
862   while((loader = (AliLoader*)next()))
863    {
864      res = loader->WriteTracks(opt);
865      if (res)
866       {
867         Error("WriteTracks","Failed to write Tracks for %s.",
868               loader->GetDetectorName().Data());
869         result = 1;
870       }
871    }
872   return result;
873 }
874 /**************************************************************************/
875
876 Int_t AliRunLoader::WriteRunLoader(Option_t* opt)
877 {
878   CdGAFile();
879   this->Write(0,TObject::kOverwrite);
880   return 0;
881 }
882 /**************************************************************************/
883
884 Int_t AliRunLoader::SetEventFolderName(const TString& name)
885 {  //sets top folder name for this run; of alread
886   if (name.IsNull())
887    {
888      Error("SetTopFolderName","Name is empty");
889      return 1;
890    }
891   
892   //check if such a folder already exists - try to find it in alice top folder
893   TObject* obj = AliConfig::Instance()->GetTopFolder()->FindObject(name);
894   if(obj)
895    {
896      TFolder* fold = dynamic_cast<TFolder*>(obj);
897      if (fold == 0x0)
898       {
899        Error("SetTopFolderName","Such a obejct already exists in top alice folder and it is not a folder.");
900        return 2;
901       }
902      //folder which was found is our folder
903      if (fEventFolder == fold)
904       {
905        return 0;
906       }
907      else
908       {
909        Error("SetTopFolderName","Such a folder already exists in top alice folder. Can not mount.");
910        return 2;
911       }
912    }
913
914   //event is alredy connected, just change name of the folder
915   if (fEventFolder) 
916    {
917      fEventFolder->SetName(name);
918      return 0;
919    }
920
921   if (fKineDataLoader == 0x0)
922     fKineDataLoader = new AliDataLoader(fgkDefaultKineFileName,fgkKineContainerName,"Kinematics");
923
924   if ( fTrackRefsDataLoader == 0x0)
925     fTrackRefsDataLoader = new AliDataLoader(fgkDefaultTrackRefsFileName,fgkTrackRefsContainerName,"Track References");
926    
927   //build the event folder structure
928   if (GetDebug()) Info("SetEventFolderName","Creating new event folder named %s",name.Data());
929   fEventFolder = AliConfig::Instance()->BuildEventFolder(name,"Event Folder");
930   fEventFolder->Add(this);//put myself to the folder to accessible for all
931   
932   if (Stack()) Stack()->SetEventFolderName(fEventFolder->GetName());
933   TIter next(fLoaders);
934   AliLoader *loader;
935   while((loader = (AliLoader*)next()))
936    {
937      loader->Register(fEventFolder);//build folder structure for this detector
938    }
939   
940   fKineDataLoader->SetEventFolder(GetEventFolder());
941   fTrackRefsDataLoader->SetEventFolder(GetEventFolder());
942   fKineDataLoader->SetFolder(GetEventFolder());
943   fTrackRefsDataLoader->SetFolder(GetEventFolder());
944   
945   fEventFolder->SetOwner();
946   return 0;
947 }
948 /**************************************************************************/
949
950 void AliRunLoader::AddLoader(AliLoader* loader)
951  {
952  //Adds the Loader for given detector 
953   if (loader == 0x0) //if null shout and exit
954    {
955      Error("AddLoader","Parameter is NULL");
956      return;
957    }
958   loader->SetDirName(fUnixDirName);
959   if (fEventFolder) loader->SetEventFolder(fEventFolder); //if event folder is already defined, 
960                                                           //pass information to the Loader
961   fLoaders->Add(loader);//add the Loader to the array
962  }
963 /**************************************************************************/
964
965 void AliRunLoader::AddLoader(AliDetector* det)
966  {
967 //Asks module (detector) ro make a Loader and stores in the array
968    if (det == 0x0) return;
969    AliLoader* get = det->GetLoader();//try to get loader
970    if (get == 0x0)  get = det->MakeLoader(fEventFolder->GetName());//if did not obtain, ask to make it
971
972    if (get) 
973     {
974       if (GetDebug()) Info("AddLoader","Detector: %s   Loader : %s",det->GetName(),get->GetName());
975       AddLoader(get);
976     }
977  }
978
979 /**************************************************************************/
980
981 AliLoader* AliRunLoader::GetLoader(const char* detname) const
982 {
983   return (AliLoader*)fLoaders->FindObject(detname);
984 }
985
986 /**************************************************************************/
987
988 AliLoader* AliRunLoader::GetLoader(AliDetector* det) const
989 {
990  if(det == 0x0) return 0x0;
991  TString getname(det->GetName());
992  getname+="Loader";
993  if (GetDebug()) Info("GetLoader(AliDetector* det)"," Loader name is %s",getname.Data());
994  return GetLoader(getname);
995 }
996
997 /**************************************************************************/
998
999 void AliRunLoader::CleanFolders()
1000 {
1001 //  fEventFolder->Add(this);//put myself to the folder to accessible for all
1002
1003   CleanDetectors();
1004   CleanHeader();
1005   CleanKinematics();
1006 }
1007 /**************************************************************************/
1008
1009 void AliRunLoader::CleanDetectors()
1010 {
1011 //Calls CleanFolders for all detectors
1012   TIter next(fLoaders);
1013   AliLoader *Loader;
1014   while((Loader = (AliLoader*)next())) 
1015    {
1016      Loader->CleanFolders();
1017    }
1018 }
1019 /**************************************************************************/
1020
1021 void AliRunLoader::RemoveEventFolder()
1022 {
1023 //remove all the tree of event 
1024 //all the stuff changing EbE stays untached (PDGDB, tasks, etc.)
1025
1026  if (fEventFolder == 0x0) return;
1027  fEventFolder->SetOwner(kFALSE);//don't we want to deleted while removing the folder that we are sitting in
1028  fEventFolder->Remove(this);//remove us drom folder
1029  
1030  AliConfig::Instance()->GetTopFolder()->SetOwner(); //brings ownership back for fEventFolder since it sits in top folder
1031  AliConfig::Instance()->GetTopFolder()->Remove(fEventFolder); //remove the event tree
1032  delete fEventFolder;
1033 }
1034 /**************************************************************************/
1035
1036 void AliRunLoader::SetGAliceFile(TFile* gafile)
1037 {
1038  fGAFile = gafile;
1039 }
1040
1041 /**************************************************************************/
1042
1043 Int_t AliRunLoader::LoadHits(Option_t* detectors,Option_t* opt)
1044 {
1045 //LoadHits in selected detectors i.e. detectors="ITS TPC TRD" or "all"
1046
1047   if (GetDebug()) Info("LoadHits","Loading Hits");
1048   TObjArray* loaders;
1049   TObjArray arr;
1050
1051   const char* oAll = strstr(detectors,"all");
1052   if (oAll)
1053    {
1054      if (GetDebug()) Info("LoadHits","Option is All");
1055      loaders = fLoaders;
1056    }
1057   else
1058    {
1059      GetListOfDetectors(detectors,arr);//this method looks for all Loaders corresponding to names (many) specified in detectors option
1060      loaders = &arr;//get the pointer array
1061    }   
1062
1063   if (GetDebug()) Info("LoadHits","For detectors. Number of detectors chosen for loading %d",loaders->GetEntries());
1064   
1065   TIter next(loaders);
1066   AliLoader *loader;
1067   while((loader = (AliLoader*)next())) 
1068    {
1069     if (GetDebug()) Info("LoadHits","    Calling LoadHits(%s) for %s",opt,loader->GetName());
1070     loader->LoadHits(opt);
1071    }
1072   if (GetDebug()) Info("LoadHits","Done");
1073   return 0;
1074
1075
1076 /**************************************************************************/
1077
1078 Int_t AliRunLoader::LoadSDigits(Option_t* detectors,Option_t* opt)
1079 {
1080 //LoadHits in selected detectors i.e. detectors="ITS TPC TRD" or "all"
1081
1082   TObjArray* loaders;
1083   TObjArray arr;
1084
1085   const char* oAll = strstr(detectors,"all");
1086   if (oAll)
1087    {
1088      loaders = fLoaders;
1089    }
1090   else
1091    {
1092      GetListOfDetectors(detectors,arr);//this method looks for all Loaders corresponding to names (many) specified in detectors option
1093      loaders = &arr;//get the pointer array
1094    }   
1095
1096   TIter next(loaders);
1097   AliLoader *loader;
1098   while((loader = (AliLoader*)next())) 
1099    {
1100     loader->LoadSDigits(opt);
1101    }
1102   return 0;
1103
1104
1105 /**************************************************************************/
1106
1107 Int_t AliRunLoader::LoadDigits(Option_t* detectors,Option_t* opt)
1108 {
1109 //LoadHits in selected detectors i.e. detectors="ITS TPC TRD" or "all"
1110
1111   TObjArray* Loaders;
1112   TObjArray arr;
1113
1114   const char* oAll = strstr(detectors,"all");
1115   if (oAll)
1116    {
1117      Loaders = fLoaders;
1118    }
1119   else
1120    {
1121      GetListOfDetectors(detectors,arr);//this method looks for all Loaders corresponding to names (many) specified in detectors option
1122      Loaders = &arr;//get the pointer array
1123    }   
1124
1125   TIter next(Loaders);
1126   AliLoader *Loader;
1127   while((Loader = (AliLoader*)next())) 
1128    {
1129     Loader->LoadDigits(opt);
1130    }
1131   return 0;
1132
1133
1134 /**************************************************************************/
1135
1136 Int_t AliRunLoader::LoadRecPoints(Option_t* detectors,Option_t* opt)
1137 {
1138 //LoadHits in selected detectors i.e. detectors="ITS TPC TRD" or "all"
1139
1140   TObjArray* Loaders;
1141   TObjArray arr;
1142
1143   const char* oAll = strstr(detectors,"all");
1144   if (oAll)
1145    {
1146      Loaders = fLoaders;
1147    }
1148   else
1149    {
1150      GetListOfDetectors(detectors,arr);//this method looks for all Loaders corresponding to names (many) specified in detectors option
1151      Loaders = &arr;//get the pointer array
1152    }   
1153
1154   TIter next(Loaders);
1155   AliLoader *Loader;
1156   while((Loader = (AliLoader*)next())) 
1157    {
1158     Loader->LoadRecPoints(opt);
1159    }
1160   return 0;
1161
1162
1163 /**************************************************************************/
1164
1165 Int_t AliRunLoader::LoadTracks(Option_t* detectors,Option_t* opt)
1166 {
1167 //LoadHits in selected detectors i.e. detectors="ITS TPC TRD" or "all"
1168
1169   TObjArray* Loaders;
1170   TObjArray arr;
1171
1172   const char* oAll = strstr(detectors,"all");
1173   if (oAll)
1174    {
1175      Loaders = fLoaders;
1176    }
1177   else
1178    {
1179      GetListOfDetectors(detectors,arr);//this method looks for all Loaders corresponding to names (many) specified in detectors option
1180      Loaders = &arr;//get the pointer array
1181    }   
1182
1183   TIter next(Loaders);
1184   AliLoader *Loader;
1185   while((Loader = (AliLoader*)next())) 
1186    {
1187     Loader->LoadTracks(opt);
1188    }
1189   return 0;
1190
1191
1192 /**************************************************************************/
1193
1194 AliRunLoader* AliRunLoader::GetRunLoader(const char* eventfoldername)
1195  {
1196   TFolder* evfold= dynamic_cast<TFolder*>(AliConfig::Instance()->GetTopFolder()->FindObject(eventfoldername));
1197   if (evfold == 0x0)
1198    {
1199      return 0x0;
1200    }
1201   AliRunLoader* runget = dynamic_cast<AliRunLoader*>(evfold->FindObject(AliRunLoader::fgkRunLoaderName));
1202   return runget;
1203   
1204  }
1205 /**************************************************************************/
1206
1207 void AliRunLoader::CdGAFile()
1208 {
1209 //sets gDirectory to galice file
1210 //work around 
1211   if(fGAFile) fGAFile->cd();
1212 }
1213  
1214 /**************************************************************************/
1215
1216 void AliRunLoader::GetListOfDetectors(const char * namelist,TObjArray& pointerarray) const
1217  {
1218 //this method looks for all Loaders corresponding 
1219 //to names (many) specified in namelist i.e. namelist ("ITS TPC TRD")
1220   
1221    char buff[10];
1222    char dets [200];
1223    strcpy(dets,namelist);//compiler cries when char* = const Option_t*;
1224    dets[strlen(dets)+1] = '\n';//set endl at the end of string 
1225    char* pdet = dets;
1226    Int_t tmp;
1227    for(;;)
1228     {
1229       tmp = sscanf(pdet,"%s",buff);//read the string from the input string pdet into buff
1230       if ( (buff[0] == 0) || (tmp == 0) ) break; //if not read
1231      
1232       pdet = strstr(pdet,buff) + strlen(buff);;//move the input pointer about number of bytes (letters) read
1233       //I am aware that is a little complicated. I don't know the number of spaces between detector names
1234       //so I read the string, than I find where it starts (strstr) and move the pointer about length of a string
1235       //If there is a better way, please write me (Piotr.Skowronski@cern.ch)
1236       //construct the Loader name
1237       TString getname(buff);
1238       getname+="Loader";
1239       AliLoader* loader = GetLoader(getname);//get the Loader
1240       if (loader)
1241        {
1242         pointerarray.Add(loader);
1243        }
1244       else
1245        {
1246         Error("GetListOfDetectors","Can not find Loader for %s",buff);
1247        }
1248         
1249       buff[0] = 0;
1250     }
1251  }
1252 /*****************************************************************************/ 
1253
1254 void AliRunLoader::Clean(const TString& name)
1255 {
1256 //removes object with given name from event folder and deletes it
1257   if (GetEventFolder() == 0x0) return;
1258   TObject* obj = GetEventFolder()->FindObject(name);
1259   if(obj)
1260    {
1261      if (GetDebug()) Info("Clean(const TString&)","name=%s, cleaning %s.",GetName(),name.Data());
1262      GetEventFolder()->Remove(obj);
1263      delete obj;
1264    }
1265 }
1266
1267 /*****************************************************************************/ 
1268
1269 TTask* AliRunLoader::GetRunDigitizer()
1270 {
1271 //returns Run Digitizer from folder
1272
1273  TFolder* topf = AliConfig::Instance()->GetTaskFolder();
1274  TObject* obj = topf->FindObjectAny(AliConfig::Instance()->GetDigitizerTaskName());
1275  return (obj)?dynamic_cast<TTask*>(obj):0x0;
1276 }
1277 /*****************************************************************************/ 
1278
1279 TTask* AliRunLoader::GetRunSDigitizer()
1280 {
1281 //returns SDigitizer Task from folder
1282
1283  TFolder* topf = AliConfig::Instance()->GetTaskFolder();
1284  TObject* obj = topf->FindObjectAny(AliConfig::Instance()->GetSDigitizerTaskName());
1285  return (obj)?dynamic_cast<TTask*>(obj):0x0;
1286 }
1287 /*****************************************************************************/ 
1288
1289 TTask* AliRunLoader::GetRunReconstructioner()
1290 {
1291 //returns Reconstructioner Task from folder
1292  TFolder* topf = AliConfig::Instance()->GetTaskFolder();
1293  TObject* obj = topf->FindObjectAny(AliConfig::Instance()->GetReconstructionerTaskName());
1294  return (obj)?dynamic_cast<TTask*>(obj):0x0;
1295 }
1296 /*****************************************************************************/ 
1297
1298 TTask* AliRunLoader::GetRunTracker()
1299 {
1300 //returns Tracker Task from folder
1301  TFolder* topf = AliConfig::Instance()->GetTaskFolder();
1302  TObject* obj = topf->FindObjectAny(AliConfig::Instance()->GetTrackerTaskName());
1303  return (obj)?dynamic_cast<TTask*>(obj):0x0;
1304 }
1305 /*****************************************************************************/ 
1306
1307 TTask* AliRunLoader::GetRunPIDTask()
1308 {
1309 //returns Tracker Task from folder
1310  TFolder* topf = AliConfig::Instance()->GetTaskFolder();
1311  TObject* obj = topf->FindObjectAny(AliConfig::Instance()->GetPIDTaskName());
1312  return (obj)?dynamic_cast<TTask*>(obj):0x0;
1313 }
1314 /*****************************************************************************/ 
1315
1316 TTask* AliRunLoader::GetRunQATask()
1317 {
1318 //returns Quality Assurance Task from folder
1319  TFolder* topf = AliConfig::Instance()->GetTaskFolder();
1320  if (topf == 0x0)
1321   {
1322     ::Error("AliRunLoader::GetRunQATask","Can not get task folder from AliConfig");
1323     return 0x0;
1324   }
1325  TObject* obj = topf->FindObjectAny(AliConfig::Instance()->GetQATaskName());
1326  return (obj)?dynamic_cast<TTask*>(obj):0x0;
1327 }
1328
1329 /*****************************************************************************/ 
1330
1331 void AliRunLoader::SetCompressionLevel(Int_t cl)
1332 {
1333 //Sets Compression Level in all files
1334  if (fGAFile) fGAFile->SetCompressionLevel(cl);
1335  SetKineComprLevel(cl);
1336  SetTrackRefsComprLevel(cl);
1337  TIter next(fLoaders);
1338  AliLoader *loader;
1339  while((loader = (AliLoader*)next()))
1340   {
1341    loader->SetCompressionLevel(cl);
1342   }
1343 }
1344 /**************************************************************************/
1345
1346 void AliRunLoader::SetKineComprLevel(Int_t cl)
1347 {
1348 //Sets comression level in Kine File
1349   fKineDataLoader->SetCompressionLevel(cl);
1350 }
1351 /**************************************************************************/
1352
1353 void AliRunLoader::SetTrackRefsComprLevel(Int_t cl)
1354 {
1355 //Sets comression level in Track Refences File
1356   fTrackRefsDataLoader->SetCompressionLevel(cl);
1357 }
1358 /**************************************************************************/
1359
1360 void AliRunLoader::UnloadHeader()
1361 {
1362  //removes TreeE from folder and deletes it
1363  // as well as fHeader object
1364  CleanHeader();
1365  delete fHeader;
1366  fHeader = 0x0;
1367 }
1368 /**************************************************************************/
1369
1370 void AliRunLoader::UnloadKinematics()
1371 {
1372  fKineDataLoader->GetBaseLoader(0)->Unload();
1373 }
1374 /**************************************************************************/
1375
1376 void AliRunLoader::UnloadTrackRefs()
1377 {
1378  fTrackRefsDataLoader->GetBaseLoader(0)->Unload();
1379 }
1380 /**************************************************************************/
1381
1382 void AliRunLoader::UnloadgAlice()
1383 {
1384  if (gAlice == GetAliRun())
1385   {
1386    if (GetDebug()) Info("UnloadgAlice","Set gAlice = 0x0");
1387    gAlice = 0x0;//if gAlice is the same that in folder (to be deleted by the way of folder)
1388   }
1389  AliRun* alirun = GetAliRun();
1390  if (GetEventFolder()) GetEventFolder()->Remove(alirun);
1391  delete alirun;
1392 }
1393 /**************************************************************************/
1394
1395 void  AliRunLoader::MakeTrackRefsContainer()
1396 {
1397 // Makes a tree for Track References
1398   fTrackRefsDataLoader->MakeTree();
1399 }
1400 /**************************************************************************/
1401
1402 Int_t AliRunLoader::LoadTrackRefs(Option_t* option)
1403 {
1404 //Load track references from file (opens file and posts tree to folder)
1405
1406  return fKineDataLoader->GetBaseLoader(0)->Load(option);
1407 }
1408 /**************************************************************************/
1409
1410 void  AliRunLoader::SetDirName(TString& dirname)
1411 {
1412 //sets directory name 
1413   if (dirname.IsNull()) return;
1414   fUnixDirName = dirname;
1415   fKineDataLoader->SetDirName(dirname);
1416   fTrackRefsDataLoader->SetDirName(dirname);
1417   
1418   TIter next(fLoaders);
1419   AliLoader *loader;
1420   while((loader = (AliLoader*)next()))
1421    {
1422     loader->SetDirName(dirname);
1423    }
1424
1425 }
1426 /*****************************************************************************/ 
1427
1428 Int_t AliRunLoader::GetFileOffset() const
1429 {
1430   return Int_t(fCurrentEvent/fNEventsPerFile);
1431 }
1432
1433 /*****************************************************************************/ 
1434 const TString AliRunLoader::SetFileOffset(const TString& fname)
1435 {
1436   Long_t offset = (Long_t)GetFileOffset();
1437   if (offset < 1) return fname;
1438   TString soffset;
1439   soffset += offset;//automatic conversion to string
1440   TString dotroot(".root");
1441   const TString& offfsetdotroot = offset + dotroot;
1442   TString out = fname;
1443   out = out.ReplaceAll(dotroot,offfsetdotroot);
1444   if (GetDebug()) Info("SetFileOffset"," in=%s out=%s",fname.Data(),out.Data());
1445   return out;
1446 }
1447 /*****************************************************************************/ 
1448
1449 void AliRunLoader::SetDigitsFileNameSuffix(const TString& suffix)
1450 {
1451 //adds the suffix before ".root", 
1452 //e.g. TPC.Digits.root -> TPC.DigitsMerged.root
1453 //made on Jiri Chudoba demand
1454
1455   TIter next(fLoaders);
1456   AliLoader *Loader;
1457   while((Loader = (AliLoader*)next())) 
1458    {
1459      Loader->SetDigitsFileNameSuffix(suffix);
1460    }
1461 }
1462 /*****************************************************************************/ 
1463
1464 TString AliRunLoader::GetFileName() const
1465 {
1466 //returns name of galice file
1467  TString result;
1468  if (fGAFile == 0x0) return result;
1469  result = fGAFile->GetName();
1470  return result;
1471 }
1472 /*****************************************************************************/ 
1473
1474 void AliRunLoader::SetDetectorAddresses()
1475 {
1476  //calls SetTreeAddress for all detectors
1477   if (GetAliRun()==0x0) return;
1478   TIter next(GetAliRun()->Modules());
1479   AliModule* mod;
1480   while((mod = (AliModule*)next())) 
1481    {
1482      AliDetector* det = dynamic_cast<AliDetector*>(mod);
1483      if (det) det->SetTreeAddress();
1484    }
1485 }
1486 /*****************************************************************************/ 
1487 /*****************************************************************************/ 
1488 /*****************************************************************************/