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