]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliRunLoader.cxx
Changes for root v4-00-02 (F.Carminati)
[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("LoadHeader","Event folder not specified yet");
591     return 1;
592   }
593
594  if (fGAFile == 0x0)
595   {
596     Error("LoadHeader","Session not opened. Use AliRunLoader::Open");
597     return 2;
598   }
599  
600  if (fGAFile->IsOpen() == kFALSE)
601   {
602     Error("LoadHeader","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("LoadHeader","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) 
634   {
635     retval = fStack->GetEvent();
636     if ( retval == kFALSE)
637      {
638        Error("LoadKinematics","Error occured while loading kinamatics tree.");
639        return retval;
640      }
641     
642   }
643  return 0;
644 }
645 /**************************************************************************/
646
647 Int_t AliRunLoader::OpenDataFile(const TString& filename,TFile*& file,TDirectory*& dir,Option_t* opt,Int_t cl)
648 {
649 //Opens File with kinematics
650  if (file)
651   {
652     if (file->IsOpen() == kFALSE)
653      {//pointer is not null but file is not opened
654        Warning("OpenDataFile","Pointer to file is not null, but file is not opened");//risky any way
655        delete file;
656        file = 0x0; //proceed with opening procedure
657      }
658     else
659      { 
660        Warning("OpenDataFile","File  %s already opened",filename.Data());
661        return 0;
662      }
663   }
664 //try to find if that file is opened somewere else
665  file = (TFile *)( gROOT->GetListOfFiles()->FindObject(filename) );
666  if (file)
667   {
668    if(file->IsOpen() == kTRUE)
669     {
670      Warning("OpenDataFile","File %s already opened by sombody else.",file->GetName());
671      return 0;
672     }
673   }
674
675  file = TFile::Open(filename,opt);
676  if (file == 0x0)
677   {//file is null
678     Error("LoadKinematics","Can not open file %s",filename.Data());
679     return 1;
680   }
681  if (file->IsOpen() == kFALSE)
682   {//file is not opened
683    Error("LoadKinematics","Can not open file %s",filename.Data());
684    return 1;
685   }
686   
687  file->SetCompressionLevel(cl);
688  
689  dir = AliLoader::ChangeDir(file,fCurrentEvent);
690  if (dir == 0x0)
691   {
692     Error("OpenKineFile","Can not change to root directory in file %s",filename.Data());
693     return 3;
694   }
695  return 0; 
696 }
697 /**************************************************************************/
698
699 TTree* AliRunLoader::TreeE() const
700 {
701  //returns the tree from folder; shortcut method
702  if (GetDebug() > 10) fEventFolder->ls();
703  TObject *obj = fEventFolder->FindObject(fgkHeaderContainerName);
704  return (obj)?dynamic_cast<TTree*>(obj):0x0;
705 }
706 /**************************************************************************/
707
708 AliHeader* AliRunLoader::GetHeader() const
709 {
710 //returns pointer header object
711  return fHeader;
712 }
713 /**************************************************************************/
714  
715 TTree* AliRunLoader::TreeK() const
716 {
717  //returns the tree from folder; shortcut method
718  TObject *obj = GetEventFolder()->FindObject(fgkKineContainerName);
719  return (obj)?dynamic_cast<TTree*>(obj):0x0;
720 }
721 /**************************************************************************/
722
723 TTree* AliRunLoader::TreeTR() const
724 {
725  //returns the tree from folder; shortcut method
726  TObject* obj = GetEventFolder()->FindObject(fgkTrackRefsContainerName);
727  return (obj)?dynamic_cast<TTree*>(obj):0x0;
728 }
729 /**************************************************************************/
730
731 AliRun* AliRunLoader::GetAliRun() const
732 {
733 //returns AliRun which sits in the folder
734  if (fEventFolder == 0x0) return 0x0;
735  TObject *obj = fEventFolder->FindObject(fgkGAliceName);
736  return (obj)?dynamic_cast<AliRun*>(obj):0x0;
737 }
738 /**************************************************************************/
739
740 Int_t AliRunLoader::WriteGeometry(Option_t* /*opt*/)
741 {
742 //writes geometry to the file
743   fGAFile->cd();
744   TGeometry* geo = GetAliRun()->GetGeometry();
745   if (geo == 0x0)
746    {
747      Error("WriteGeometry","Can not get geometry from gAlice");
748      return 1;
749    }
750   geo->Write();
751   return 0;
752 }
753 /**************************************************************************/
754
755 Int_t AliRunLoader::WriteHeader(Option_t* opt)
756 {
757 //writes treeE
758   if (GetDebug()) Info("WriteHeader","  WRITING HEADER");
759   
760   TTree* tree = TreeE();
761   if ( tree == 0x0)
762    {
763      Warning("WriteHeader","Can not find Header Tree in Folder");
764      return 0;
765    } 
766   if (fGAFile->IsWritable() == kFALSE)
767    {
768      Error("WriteHeader","File %s is not writable",fGAFile->GetName());
769      return 1;
770    }
771
772   TObject* obj = fGAFile->Get(fgkHeaderContainerName);
773   if (obj)
774    { //if they exist, see if option OVERWRITE is used
775      TString tmp(opt);
776      if(tmp.Contains("OVERWRITE",TString::kIgnoreCase) == 0)
777       {//if it is not used -  give an error message and return an error code
778         Error("WriteHeader","Tree already exisists. Use option \"OVERWRITE\" to overwrite previous data");
779         return 3;
780       }
781    }
782   fGAFile->cd();
783   tree->SetDirectory(fGAFile);
784   tree->Write(0,TObject::kOverwrite);
785
786   if (GetDebug()) Info("WriteHeader","WRITTEN\n\n");
787   
788   return 0;
789 }
790 /**************************************************************************/
791
792 Int_t AliRunLoader::WriteAliRun(Option_t* /*opt*/)
793 {
794 //writes AliRun object to the file
795   fGAFile->cd();
796   if (GetAliRun()) GetAliRun()->Write();
797   return 0;
798 }
799 /**************************************************************************/
800
801 Int_t AliRunLoader::WriteKinematics(Option_t* opt)
802 {
803 //writes Kinematics
804   return fKineDataLoader->GetBaseLoader(0)->WriteData(opt);
805 }
806 /**************************************************************************/
807 Int_t AliRunLoader::WriteTrackRefs(Option_t* opt)
808 {
809 //writes Track References tree
810   return fTrackRefsDataLoader->GetBaseLoader(0)->WriteData(opt);
811 }
812 /**************************************************************************/
813
814 Int_t AliRunLoader::WriteHits(Option_t* opt)
815 {
816 //Calls WriteHits for all loaders
817   Int_t res;
818   Int_t result = 0;
819   TIter next(fLoaders);
820   AliLoader *loader;
821   while((loader = (AliLoader*)next()))
822    {
823      res = loader->WriteHits(opt);
824      if (res)
825       {
826         Error("WriteHits","Failed to write hits for %s (%d)",loader->GetDetectorName().Data(),res);
827         result = 1;
828       }
829    }
830   return result;
831 }
832 /**************************************************************************/
833
834 Int_t AliRunLoader::WriteSDigits(Option_t* opt)
835 {
836 //Calls WriteSDigits for all loaders
837   Int_t res;
838   Int_t result = 0;
839   TIter next(fLoaders);
840   AliLoader *loader;
841   while((loader = (AliLoader*)next()))
842    {
843      res = loader->WriteSDigits(opt);
844      if (res)
845       {
846         Error("WriteSDigits","Failed to write summable digits for %s.",loader->GetDetectorName().Data());
847         result = 1;
848       }
849    }
850   return result;
851 }
852 /**************************************************************************/
853
854 Int_t AliRunLoader::WriteDigits(Option_t* opt)
855 {
856 //Calls WriteDigits for all loaders
857   Int_t res;
858   Int_t result = 0;
859   TIter next(fLoaders);
860   AliLoader *loader;
861   while((loader = (AliLoader*)next()))
862    { 
863      res = loader->WriteDigits(opt);
864      if (res)
865       {
866         Error("WriteDigits","Failed to write digits for %s.",loader->GetDetectorName().Data());
867         result = 1;
868       }
869    }
870   return result;
871 }
872 /**************************************************************************/
873
874 Int_t AliRunLoader::WriteRecPoints(Option_t* opt)
875 {
876 //Calls WriteRecPoints for all loaders
877   Int_t res;
878   Int_t result = 0;
879   TIter next(fLoaders);
880   AliLoader *loader;
881   while((loader = (AliLoader*)next()))
882    {
883      res = loader->WriteRecPoints(opt);
884      if (res)
885       {
886         Error("WriteRecPoints","Failed to write Reconstructed Points for %s.",
887               loader->GetDetectorName().Data());
888         result = 1;
889       }
890    }
891   return result;
892 }
893 /**************************************************************************/
894
895 Int_t AliRunLoader::WriteTracks(Option_t* opt)
896 {
897 //Calls WriteTracks for all loaders
898   Int_t res;
899   Int_t result = 0;
900   TIter next(fLoaders);
901   AliLoader *loader;
902   while((loader = (AliLoader*)next()))
903    {
904      res = loader->WriteTracks(opt);
905      if (res)
906       {
907         Error("WriteTracks","Failed to write Tracks for %s.",
908               loader->GetDetectorName().Data());
909         result = 1;
910       }
911    }
912   return result;
913 }
914 /**************************************************************************/
915
916 Int_t AliRunLoader::WriteRunLoader(Option_t* /*opt*/)
917 {
918 //Writes itself to the file
919   CdGAFile();
920   this->Write(0,TObject::kOverwrite);
921   return 0;
922 }
923 /**************************************************************************/
924
925 Int_t AliRunLoader::SetEventFolderName(const TString& name)
926 {  
927 //sets top folder name for this run; of alread
928   if (name.IsNull())
929    {
930      Error("SetTopFolderName","Name is empty");
931      return 1;
932    }
933   
934   //check if such a folder already exists - try to find it in alice top folder
935   TObject* obj = AliConfig::Instance()->GetTopFolder()->FindObject(name);
936   if(obj)
937    {
938      TFolder* fold = dynamic_cast<TFolder*>(obj);
939      if (fold == 0x0)
940       {
941        Error("SetTopFolderName","Such a obejct already exists in top alice folder and it is not a folder.");
942        return 2;
943       }
944      //folder which was found is our folder
945      if (fEventFolder == fold)
946       {
947        return 0;
948       }
949      else
950       {
951        Error("SetTopFolderName","Such a folder already exists in top alice folder. Can not mount.");
952        return 2;
953       }
954    }
955
956   //event is alredy connected, just change name of the folder
957   if (fEventFolder) 
958    {
959      fEventFolder->SetName(name);
960      return 0;
961    }
962
963   if (fKineDataLoader == 0x0)
964     fKineDataLoader = new AliDataLoader(fgkDefaultKineFileName,fgkKineContainerName,"Kinematics");
965
966   if ( fTrackRefsDataLoader == 0x0)
967     fTrackRefsDataLoader = new AliDataLoader(fgkDefaultTrackRefsFileName,fgkTrackRefsContainerName,"Track References");
968    
969   //build the event folder structure
970   if (GetDebug()) Info("SetEventFolderName","Creating new event folder named %s",name.Data());
971   fEventFolder = AliConfig::Instance()->BuildEventFolder(name,"Event Folder");
972   fEventFolder->Add(this);//put myself to the folder to accessible for all
973   
974   if (Stack()) Stack()->SetEventFolderName(fEventFolder->GetName());
975   TIter next(fLoaders);
976   AliLoader *loader;
977   while((loader = (AliLoader*)next()))
978    {
979      loader->Register(fEventFolder);//build folder structure for this detector
980    }
981   
982   fKineDataLoader->SetEventFolder(GetEventFolder());
983   fTrackRefsDataLoader->SetEventFolder(GetEventFolder());
984   fKineDataLoader->SetFolder(GetEventFolder());
985   fTrackRefsDataLoader->SetFolder(GetEventFolder());
986   
987   fEventFolder->SetOwner();
988   return 0;
989 }
990 /**************************************************************************/
991
992 void AliRunLoader::AddLoader(AliLoader* loader)
993  {
994  //Adds the Loader for given detector 
995   if (loader == 0x0) //if null shout and exit
996    {
997      Error("AddLoader","Parameter is NULL");
998      return;
999    }
1000   loader->SetDirName(fUnixDirName);
1001   if (fEventFolder) loader->SetEventFolder(fEventFolder); //if event folder is already defined, 
1002                                                           //pass information to the Loader
1003   fLoaders->Add(loader);//add the Loader to the array
1004  }
1005 /**************************************************************************/
1006
1007 void AliRunLoader::AddLoader(AliDetector* det)
1008  {
1009 //Asks module (detector) ro make a Loader and stores in the array
1010    if (det == 0x0) return;
1011    AliLoader* get = det->GetLoader();//try to get loader
1012    if (get == 0x0)  get = det->MakeLoader(fEventFolder->GetName());//if did not obtain, ask to make it
1013
1014    if (get) 
1015     {
1016       if (GetDebug()) Info("AddLoader","Detector: %s   Loader : %s",det->GetName(),get->GetName());
1017       AddLoader(get);
1018     }
1019  }
1020
1021 /**************************************************************************/
1022
1023 AliLoader* AliRunLoader::GetLoader(const char* detname) const
1024 {
1025 //returns loader for given detector
1026 //note that naming convention is TPCLoader not just TPC
1027   return (AliLoader*)fLoaders->FindObject(detname);
1028 }
1029
1030 /**************************************************************************/
1031
1032 AliLoader* AliRunLoader::GetLoader(AliDetector* det) const
1033 {
1034 //get loader for detector det
1035  if(det == 0x0) return 0x0;
1036  TString getname(det->GetName());
1037  getname+="Loader";
1038  if (GetDebug()) Info("GetLoader(AliDetector* det)"," Loader name is %s",getname.Data());
1039  return GetLoader(getname);
1040 }
1041
1042 /**************************************************************************/
1043
1044 void AliRunLoader::CleanFolders()
1045 {
1046 //  fEventFolder->Add(this);//put myself to the folder to accessible for all
1047
1048   CleanDetectors();
1049   CleanHeader();
1050   CleanKinematics();
1051 }
1052 /**************************************************************************/
1053
1054 void AliRunLoader::CleanDetectors()
1055 {
1056 //Calls CleanFolders for all detectors
1057   TIter next(fLoaders);
1058   AliLoader *loader;
1059   while((loader = (AliLoader*)next())) 
1060    {
1061      loader->CleanFolders();
1062    }
1063 }
1064 /**************************************************************************/
1065
1066 void AliRunLoader::RemoveEventFolder()
1067 {
1068 //remove all the tree of event 
1069 //all the stuff changing EbE stays untached (PDGDB, tasks, etc.)
1070
1071  if (fEventFolder == 0x0) return;
1072  fEventFolder->SetOwner(kFALSE);//don't we want to deleted while removing the folder that we are sitting in
1073  fEventFolder->Remove(this);//remove us drom folder
1074  
1075  AliConfig::Instance()->GetTopFolder()->SetOwner(); //brings ownership back for fEventFolder since it sits in top folder
1076  AliConfig::Instance()->GetTopFolder()->Remove(fEventFolder); //remove the event tree
1077  delete fEventFolder;
1078 }
1079 /**************************************************************************/
1080
1081 void AliRunLoader::SetGAliceFile(TFile* gafile)
1082 {
1083 //sets pointer to galice.root file
1084  fGAFile = gafile;
1085 }
1086
1087 /**************************************************************************/
1088
1089 Int_t AliRunLoader::LoadHits(Option_t* detectors,Option_t* opt)
1090 {
1091 //LoadHits in selected detectors i.e. detectors="ITS TPC TRD" or "all"
1092
1093   if (GetDebug()) Info("LoadHits","Loading Hits");
1094   TObjArray* loaders;
1095   TObjArray arr;
1096
1097   const char* oAll = strstr(detectors,"all");
1098   if (oAll)
1099    {
1100      if (GetDebug()) Info("LoadHits","Option is All");
1101      loaders = fLoaders;
1102    }
1103   else
1104    {
1105      GetListOfDetectors(detectors,arr);//this method looks for all Loaders corresponding to names (many) specified in detectors option
1106      loaders = &arr;//get the pointer array
1107    }   
1108
1109   if (GetDebug()) Info("LoadHits","For detectors. Number of detectors chosen for loading %d",loaders->GetEntries());
1110   
1111   TIter next(loaders);
1112   AliLoader *loader;
1113   while((loader = (AliLoader*)next())) 
1114    {
1115     if (GetDebug()) Info("LoadHits","    Calling LoadHits(%s) for %s",opt,loader->GetName());
1116     loader->LoadHits(opt);
1117    }
1118   if (GetDebug()) Info("LoadHits","Done");
1119   return 0;
1120
1121
1122 /**************************************************************************/
1123
1124 Int_t AliRunLoader::LoadSDigits(Option_t* detectors,Option_t* opt)
1125 {
1126 //LoadHits in selected detectors i.e. detectors="ITS TPC TRD" or "all"
1127
1128   TObjArray* loaders;
1129   TObjArray arr;
1130
1131   const char* oAll = strstr(detectors,"all");
1132   if (oAll)
1133    {
1134      loaders = fLoaders;
1135    }
1136   else
1137    {
1138      GetListOfDetectors(detectors,arr);//this method looks for all Loaders corresponding to names (many) specified in detectors option
1139      loaders = &arr;//get the pointer to array
1140    }   
1141
1142   TIter next(loaders);
1143   AliLoader *loader;
1144   while((loader = (AliLoader*)next())) 
1145    {
1146     loader->LoadSDigits(opt);
1147    }
1148   return 0;
1149
1150
1151 /**************************************************************************/
1152
1153 Int_t AliRunLoader::LoadDigits(Option_t* detectors,Option_t* opt)
1154 {
1155 //LoadHits in selected detectors i.e. detectors="ITS TPC TRD" or "all"
1156
1157   TObjArray* loaders;
1158   TObjArray arr;
1159
1160   const char* oAll = strstr(detectors,"all");
1161   if (oAll)
1162    {
1163      loaders = fLoaders;
1164    }
1165   else
1166    {
1167      GetListOfDetectors(detectors,arr);//this method looks for all Loaders corresponding to names (many) specified in detectors option
1168      loaders = &arr;//get the pointer array
1169    }   
1170
1171   TIter next(loaders);
1172   AliLoader *loader;
1173   while((loader = (AliLoader*)next())) 
1174    {
1175     loader->LoadDigits(opt);
1176    }
1177   return 0;
1178
1179 /**************************************************************************/
1180
1181 Int_t AliRunLoader::LoadRecPoints(Option_t* detectors,Option_t* opt)
1182 {
1183 //LoadHits in selected detectors i.e. detectors="ITS TPC TRD" or "all"
1184
1185   TObjArray* loaders;
1186   TObjArray arr;
1187
1188   const char* oAll = strstr(detectors,"all");
1189   if (oAll)
1190    {
1191      loaders = fLoaders;
1192    }
1193   else
1194    {
1195      GetListOfDetectors(detectors,arr);//this method looks for all Loaders corresponding to names (many) specified in detectors option
1196      loaders = &arr;//get the pointer array
1197    }   
1198
1199   TIter next(loaders);
1200   AliLoader *loader;
1201   while((loader = (AliLoader*)next())) 
1202    {
1203     loader->LoadRecPoints(opt);
1204    }
1205   return 0;
1206
1207 /**************************************************************************/
1208
1209 Int_t AliRunLoader::LoadRecParticles(Option_t* detectors,Option_t* opt)
1210 {
1211 //LoadHits in selected detectors i.e. detectors="ITS TPC TRD" or "all"
1212
1213   TObjArray* loaders;
1214   TObjArray arr;
1215
1216   const char* oAll = strstr(detectors,"all");
1217   if (oAll)
1218    {
1219      loaders = fLoaders;
1220    }
1221   else
1222    {
1223      GetListOfDetectors(detectors,arr);//this method looks for all Loaders corresponding to names (many) specified in detectors option
1224      loaders = &arr;//get the pointer array
1225    }   
1226
1227   TIter next(loaders);
1228   AliLoader *loader;
1229   while((loader = (AliLoader*)next())) 
1230    {
1231     loader->LoadRecParticles(opt);
1232    }
1233   return 0;
1234
1235 /**************************************************************************/
1236
1237 Int_t AliRunLoader::LoadTracks(Option_t* detectors,Option_t* opt)
1238 {
1239 //LoadHits in selected detectors i.e. detectors="ITS TPC TRD" or "all"
1240
1241   TObjArray* loaders;
1242   TObjArray arr;
1243
1244   const char* oAll = strstr(detectors,"all");
1245   if (oAll)
1246    {
1247      loaders = fLoaders;
1248    }
1249   else
1250    {
1251      GetListOfDetectors(detectors,arr);//this method looks for all Loaders corresponding to names (many) specified in detectors option
1252      loaders = &arr;//get the pointer array
1253    }   
1254
1255   TIter next(loaders);
1256   AliLoader *loader;
1257   while((loader = (AliLoader*)next())) 
1258    {
1259     loader->LoadTracks(opt);
1260    }
1261   return 0;
1262
1263 /**************************************************************************/
1264
1265 void AliRunLoader::UnloadHits(Option_t* detectors)
1266 {
1267   //unloads hits for detectors specified in parameter
1268   TObjArray* loaders;
1269   TObjArray arr;
1270
1271   const char* oAll = strstr(detectors,"all");
1272   if (oAll)
1273    {
1274      loaders = fLoaders;
1275    }
1276   else
1277    {
1278      GetListOfDetectors(detectors,arr);//this method looks for all Loaders corresponding to names (many) specified in detectors option
1279      loaders = &arr;//get the pointer to array
1280    }   
1281
1282   TIter next(loaders);
1283   AliLoader *loader;
1284   while((loader = (AliLoader*)next())) 
1285    {
1286     loader->UnloadHits();
1287    }
1288 }
1289 /**************************************************************************/
1290
1291 void AliRunLoader::UnloadSDigits(Option_t* detectors)
1292 {
1293   //unloads SDigits for detectors specified in parameter
1294   TObjArray* loaders;
1295   TObjArray arr;
1296
1297   const char* oAll = strstr(detectors,"all");
1298   if (oAll)
1299    {
1300      loaders = fLoaders;
1301    }
1302   else
1303    {
1304      GetListOfDetectors(detectors,arr);//this method looks for all Loaders corresponding to names (many) specified in detectors option
1305      loaders = &arr;//get the pointer to array
1306    }   
1307
1308   TIter next(loaders);
1309   AliLoader *loader;
1310   while((loader = (AliLoader*)next())) 
1311    {
1312     loader->UnloadSDigits();
1313    }
1314 }
1315 /**************************************************************************/
1316
1317 void AliRunLoader::UnloadDigits(Option_t* detectors)
1318 {
1319   //unloads Digits for detectors specified in parameter
1320   TObjArray* loaders;
1321   TObjArray arr;
1322
1323   const char* oAll = strstr(detectors,"all");
1324   if (oAll)
1325    {
1326      loaders = fLoaders;
1327    }
1328   else
1329    {
1330      GetListOfDetectors(detectors,arr);//this method looks for all Loaders corresponding to names (many) specified in detectors option
1331      loaders = &arr;//get the pointer to array
1332    }   
1333
1334   TIter next(loaders);
1335   AliLoader *loader;
1336   while((loader = (AliLoader*)next())) 
1337    {
1338     loader->UnloadDigits();
1339    }
1340 }
1341 /**************************************************************************/
1342
1343 void AliRunLoader::UnloadRecPoints(Option_t* detectors)
1344 {
1345   //unloads RecPoints for detectors specified in parameter
1346   TObjArray* loaders;
1347   TObjArray arr;
1348
1349   const char* oAll = strstr(detectors,"all");
1350   if (oAll)
1351    {
1352      loaders = fLoaders;
1353    }
1354   else
1355    {
1356      GetListOfDetectors(detectors,arr);//this method looks for all Loaders corresponding to names (many) specified in detectors option
1357      loaders = &arr;//get the pointer to array
1358    }   
1359
1360   TIter next(loaders);
1361   AliLoader *loader;
1362   while((loader = (AliLoader*)next())) 
1363    {
1364     loader->UnloadRecPoints();
1365    }
1366 }
1367 /**************************************************************************/
1368
1369 void AliRunLoader::UnloadAll(Option_t* detectors)
1370 {
1371   //calls UnloadAll for detectors names specified in parameter
1372   // option "all" passed can be passed
1373   TObjArray* loaders;
1374   TObjArray arr;
1375
1376   const char* oAll = strstr(detectors,"all");
1377   if (oAll)
1378    {
1379      loaders = fLoaders;
1380    }
1381   else
1382    {
1383      GetListOfDetectors(detectors,arr);//this method looks for all Loaders corresponding to names (many) specified in detectors option
1384      loaders = &arr;//get the pointer to array
1385    }   
1386
1387   TIter next(loaders);
1388   AliLoader *loader;
1389   while((loader = (AliLoader*)next())) 
1390    {
1391     loader->UnloadAll();
1392    }
1393 }
1394 /**************************************************************************/
1395
1396 void AliRunLoader::UnloadTracks(Option_t* detectors)
1397 {
1398   //unloads Tracks for detectors specified in parameter
1399   TObjArray* loaders;
1400   TObjArray arr;
1401
1402   const char* oAll = strstr(detectors,"all");
1403   if (oAll)
1404    {
1405      loaders = fLoaders;
1406    }
1407   else
1408    {
1409      GetListOfDetectors(detectors,arr);//this method looks for all Loaders corresponding to names (many) specified in detectors option
1410      loaders = &arr;//get the pointer to array
1411    }   
1412
1413   TIter next(loaders);
1414   AliLoader *loader;
1415   while((loader = (AliLoader*)next())) 
1416    {
1417     loader->UnloadTracks();
1418    }
1419 }
1420 /**************************************************************************/
1421
1422 void AliRunLoader::UnloadRecParticles(Option_t* detectors)
1423 {
1424   //unloads Particles for detectors specified in parameter
1425   TObjArray* loaders;
1426   TObjArray arr;
1427
1428   const char* oAll = strstr(detectors,"all");
1429   if (oAll)
1430    {
1431      loaders = fLoaders;
1432    }
1433   else
1434    {
1435      GetListOfDetectors(detectors,arr);//this method looks for all Loaders corresponding to names (many) specified in detectors option
1436      loaders = &arr;//get the pointer to array
1437    }   
1438
1439   TIter next(loaders);
1440   AliLoader *loader;
1441   while((loader = (AliLoader*)next())) 
1442    {
1443     loader->UnloadRecParticles();
1444    }
1445 }
1446 /**************************************************************************/
1447
1448 AliRunLoader* AliRunLoader::GetRunLoader(const char* eventfoldername)
1449 {
1450 //returns RunLoader from folder named eventfoldername
1451   TFolder* evfold= dynamic_cast<TFolder*>(AliConfig::Instance()->GetTopFolder()->FindObject(eventfoldername));
1452   if (evfold == 0x0)
1453    {
1454      return 0x0;
1455    }
1456   AliRunLoader* runget = dynamic_cast<AliRunLoader*>(evfold->FindObject(AliRunLoader::fgkRunLoaderName));
1457   return runget;
1458   
1459 }
1460 /**************************************************************************/
1461
1462 AliLoader* AliRunLoader::GetDetectorLoader(const char* detname, const char* eventfoldername)
1463 {
1464 //get the loader of the detector with the given name from the global
1465 //run loader object
1466   AliRunLoader* runLoader = GetRunLoader(eventfoldername);
1467   if (!runLoader) {
1468     ::Error("AliRunLoader::GetDetectorLoader","No run loader found");
1469     return NULL;
1470   }
1471   return runLoader->GetDetectorLoader(detname);
1472 }
1473 /**************************************************************************/
1474
1475 AliLoader* AliRunLoader::GetDetectorLoader(const char* detname)
1476 {
1477 //get the loader of the detector with the given name from the global
1478 //run loader object
1479   
1480   char loadername[256];
1481   sprintf(loadername, "%sLoader", detname);
1482   AliLoader* loader = GetLoader(loadername);
1483   if (!loader) {
1484     Error("GetDetectorLoader", "No loader for %s found", detname);
1485     return NULL;
1486   }
1487   return loader;
1488 }
1489 /**************************************************************************/
1490
1491 TTree* AliRunLoader::GetTreeH(const char* detname, Bool_t maketree, const char* eventfoldername)
1492 {
1493 //get the tree with hits of the detector with the given name
1494 //if maketree is true and the tree does not exist, the tree is created
1495   AliLoader* loader = GetDetectorLoader(detname,eventfoldername);
1496   if (!loader) return NULL;
1497   if (!loader->TreeH() && maketree) loader->MakeTree("H");
1498   return loader->TreeH();
1499 }
1500
1501 /**************************************************************************/
1502
1503 TTree* AliRunLoader::GetTreeH(const char* detname, Bool_t maketree)
1504 {
1505 //get the tree with hits of the detector with the given name
1506 //if maketree is true and the tree does not exist, the tree is created
1507   AliLoader* loader = GetDetectorLoader(detname);
1508   if (!loader) return NULL;
1509   if (!loader->TreeH() && maketree) loader->MakeTree("H");
1510   return loader->TreeH();
1511 }
1512 /**************************************************************************/
1513
1514 TTree* AliRunLoader::GetTreeS(const char* detname, Bool_t maketree,const char* eventfoldername)
1515 {
1516 //get the tree with summable digits of the detector with the given name
1517 //if maketree is true and the tree does not exist, the tree is created
1518   AliLoader* loader = GetDetectorLoader(detname,eventfoldername);
1519   if (!loader) return NULL;
1520   if (!loader->TreeS() && maketree) loader->MakeTree("S");
1521   return loader->TreeS();
1522 }
1523 /**************************************************************************/
1524
1525 TTree* AliRunLoader::GetTreeS(const char* detname, Bool_t maketree)
1526 {
1527 //get the tree with summable digits of the detector with the given name
1528 //if maketree is true and the tree does not exist, the tree is created
1529   AliLoader* loader = GetDetectorLoader(detname);
1530   if (!loader) return NULL;
1531   if (!loader->TreeS() && maketree) loader->MakeTree("S");
1532   return loader->TreeS();
1533 }
1534 /**************************************************************************/
1535
1536 TTree* AliRunLoader::GetTreeD(const char* detname, Bool_t maketree,const char* eventfoldername)
1537 {
1538 //get the tree with digits of the detector with the given name
1539 //if maketree is true and the tree does not exist, the tree is created
1540   AliLoader* loader = GetDetectorLoader(detname,eventfoldername);
1541   if (!loader) return NULL;
1542   if (!loader->TreeD() && maketree) loader->MakeTree("D");
1543   return loader->TreeD();
1544 }
1545 /**************************************************************************/
1546
1547 TTree* AliRunLoader::GetTreeD(const char* detname, Bool_t maketree)
1548 {
1549 //get the tree with digits 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);
1552   if (!loader) return NULL;
1553   if (!loader->TreeD() && maketree) loader->MakeTree("D");
1554   return loader->TreeD();
1555 }
1556 /**************************************************************************/
1557 TTree* AliRunLoader::GetTreeR(const char* detname, Bool_t maketree,const char* eventfoldername)
1558 {
1559 //get the tree with clusters of the detector with the given name
1560 //if maketree is true and the tree does not exist, the tree is created
1561   AliLoader* loader = GetDetectorLoader(detname,eventfoldername);
1562   if (!loader) return NULL;
1563   if (!loader->TreeR() && maketree) loader->MakeTree("R");
1564   return loader->TreeR();
1565 }
1566 /**************************************************************************/
1567
1568 TTree* AliRunLoader::GetTreeR(const char* detname, Bool_t maketree)
1569 {
1570 //get the tree with clusters of the detector with the given name
1571 //if maketree is true and the tree does not exist, the tree is created
1572   AliLoader* loader = GetDetectorLoader(detname);
1573   if (!loader) return NULL;
1574   if (!loader->TreeR() && maketree) loader->MakeTree("R");
1575   return loader->TreeR();
1576 }
1577 /**************************************************************************/
1578
1579 TTree* AliRunLoader::GetTreeT(const char* detname, Bool_t maketree,const char* eventfoldername)
1580 {
1581 //get the tree with tracks of the detector with the given name
1582 //if maketree is true and the tree does not exist, the tree is created
1583   AliLoader* loader = GetDetectorLoader(detname,eventfoldername);
1584   if (!loader) return NULL;
1585   if (!loader->TreeT() && maketree) loader->MakeTree("T");
1586   return loader->TreeT();
1587 }
1588 /**************************************************************************/
1589
1590 TTree* AliRunLoader::GetTreeT(const char* detname, Bool_t maketree)
1591 {
1592 //get the tree with tracks of the detector with the given name
1593 //if maketree is true and the tree does not exist, the tree is created
1594   AliLoader* loader = GetDetectorLoader(detname);
1595   if (!loader) return NULL;
1596   if (!loader->TreeT() && maketree) loader->MakeTree("T");
1597   return loader->TreeT();
1598 }
1599 /**************************************************************************/
1600
1601 TTree* AliRunLoader::GetTreeP(const char* detname, Bool_t maketree,const char* eventfoldername)
1602 {
1603 //get the tree with particles of the detector with the given name
1604 //if maketree is true and the tree does not exist, the tree is created
1605   AliLoader* loader = GetDetectorLoader(detname,eventfoldername);
1606   if (!loader) return NULL;
1607   if (!loader->TreeP() && maketree) loader->MakeTree("P");
1608   return loader->TreeP();
1609 }
1610 /**************************************************************************/
1611
1612 TTree* AliRunLoader::GetTreeP(const char* detname, Bool_t maketree)
1613 {
1614 //get the tree with particles of the detector with the given name
1615 //if maketree is true and the tree does not exist, the tree is created
1616   AliLoader* loader = GetDetectorLoader(detname);
1617   if (!loader) return NULL;
1618   if (!loader->TreeP() && maketree) loader->MakeTree("P");
1619   return loader->TreeP();
1620 }
1621
1622 /**************************************************************************/
1623
1624 void AliRunLoader::CdGAFile()
1625 {
1626 //sets gDirectory to galice file
1627 //work around 
1628   if(fGAFile) fGAFile->cd();
1629 }
1630  
1631 /**************************************************************************/
1632
1633 void AliRunLoader::GetListOfDetectors(const char * namelist,TObjArray& pointerarray) const
1634  {
1635 //this method looks for all Loaders corresponding 
1636 //to names (many) specified in namelist i.e. namelist ("ITS TPC TRD")
1637   
1638    char buff[10];
1639    char dets [200];
1640    strcpy(dets,namelist);//compiler cries when char* = const Option_t*;
1641    dets[strlen(dets)+1] = '\n';//set endl at the end of string 
1642    char* pdet = dets;
1643    Int_t tmp;
1644    for(;;)
1645     {
1646       tmp = sscanf(pdet,"%s",buff);//read the string from the input string pdet into buff
1647       if ( (buff[0] == 0) || (tmp == 0) ) break; //if not read
1648      
1649       pdet = strstr(pdet,buff) + strlen(buff);//move the input pointer about number of bytes (letters) read
1650       //I am aware that is a little complicated. I don't know the number of spaces between detector names
1651       //so I read the string, than I find where it starts (strstr) and move the pointer about length of a string
1652       //If there is a better way, please write me (Piotr.Skowronski@cern.ch)
1653       //construct the Loader name
1654       TString getname(buff);
1655       getname+="Loader";
1656       AliLoader* loader = GetLoader(getname);//get the Loader
1657       if (loader)
1658        {
1659         pointerarray.Add(loader);
1660        }
1661       else
1662        {
1663         Error("GetListOfDetectors","Can not find Loader for %s",buff);
1664        }
1665         
1666       buff[0] = 0;
1667     }
1668  }
1669 /*****************************************************************************/ 
1670
1671 void AliRunLoader::Clean(const TString& name)
1672 {
1673 //removes object with given name from event folder and deletes it
1674   if (GetEventFolder() == 0x0) return;
1675   TObject* obj = GetEventFolder()->FindObject(name);
1676   if(obj)
1677    {
1678      if (GetDebug()) Info("Clean(const TString&)","name=%s, cleaning %s.",GetName(),name.Data());
1679      GetEventFolder()->Remove(obj);
1680      delete obj;
1681    }
1682 }
1683
1684 /*****************************************************************************/ 
1685
1686 TTask* AliRunLoader::GetRunDigitizer()
1687 {
1688 //returns Run Digitizer from folder
1689
1690  TFolder* topf = AliConfig::Instance()->GetTaskFolder();
1691  TObject* obj = topf->FindObjectAny(AliConfig::Instance()->GetDigitizerTaskName());
1692  return (obj)?dynamic_cast<TTask*>(obj):0x0;
1693 }
1694 /*****************************************************************************/ 
1695
1696 TTask* AliRunLoader::GetRunSDigitizer()
1697 {
1698 //returns SDigitizer Task from folder
1699
1700  TFolder* topf = AliConfig::Instance()->GetTaskFolder();
1701  TObject* obj = topf->FindObjectAny(AliConfig::Instance()->GetSDigitizerTaskName());
1702  return (obj)?dynamic_cast<TTask*>(obj):0x0;
1703 }
1704 /*****************************************************************************/ 
1705
1706 TTask* AliRunLoader::GetRunReconstructioner()
1707 {
1708 //returns Reconstructioner Task from folder
1709  TFolder* topf = AliConfig::Instance()->GetTaskFolder();
1710  TObject* obj = topf->FindObjectAny(AliConfig::Instance()->GetReconstructionerTaskName());
1711  return (obj)?dynamic_cast<TTask*>(obj):0x0;
1712 }
1713 /*****************************************************************************/ 
1714
1715 TTask* AliRunLoader::GetRunTracker()
1716 {
1717 //returns Tracker Task from folder
1718  TFolder* topf = AliConfig::Instance()->GetTaskFolder();
1719  TObject* obj = topf->FindObjectAny(AliConfig::Instance()->GetTrackerTaskName());
1720  return (obj)?dynamic_cast<TTask*>(obj):0x0;
1721 }
1722 /*****************************************************************************/ 
1723
1724 TTask* AliRunLoader::GetRunPIDTask()
1725 {
1726 //returns Tracker Task from folder
1727  TFolder* topf = AliConfig::Instance()->GetTaskFolder();
1728  TObject* obj = topf->FindObjectAny(AliConfig::Instance()->GetPIDTaskName());
1729  return (obj)?dynamic_cast<TTask*>(obj):0x0;
1730 }
1731 /*****************************************************************************/ 
1732
1733 TTask* AliRunLoader::GetRunQATask()
1734 {
1735 //returns Quality Assurance Task from folder
1736  TFolder* topf = AliConfig::Instance()->GetTaskFolder();
1737  if (topf == 0x0)
1738   {
1739     ::Error("AliRunLoader::GetRunQATask","Can not get task folder from AliConfig");
1740     return 0x0;
1741   }
1742  TObject* obj = topf->FindObjectAny(AliConfig::Instance()->GetQATaskName());
1743  return (obj)?dynamic_cast<TTask*>(obj):0x0;
1744 }
1745
1746 /*****************************************************************************/ 
1747
1748 void AliRunLoader::SetCompressionLevel(Int_t cl)
1749 {
1750 //Sets Compression Level in all files
1751  if (fGAFile) fGAFile->SetCompressionLevel(cl);
1752  SetKineComprLevel(cl);
1753  SetTrackRefsComprLevel(cl);
1754  TIter next(fLoaders);
1755  AliLoader *loader;
1756  while((loader = (AliLoader*)next()))
1757   {
1758    loader->SetCompressionLevel(cl);
1759   }
1760 }
1761 /**************************************************************************/
1762
1763 void AliRunLoader::SetKineComprLevel(Int_t cl)
1764 {
1765 //Sets comression level in Kine File
1766   fKineDataLoader->SetCompressionLevel(cl);
1767 }
1768 /**************************************************************************/
1769
1770 void AliRunLoader::SetTrackRefsComprLevel(Int_t cl)
1771 {
1772 //Sets comression level in Track Refences File
1773   fTrackRefsDataLoader->SetCompressionLevel(cl);
1774 }
1775 /**************************************************************************/
1776
1777 void AliRunLoader::UnloadHeader()
1778 {
1779  //removes TreeE from folder and deletes it
1780  // as well as fHeader object
1781  CleanHeader();
1782  delete fHeader;
1783  fHeader = 0x0;
1784 }
1785 /**************************************************************************/
1786
1787 void AliRunLoader::UnloadKinematics()
1788 {
1789 //Unloads Kinematics
1790  fKineDataLoader->GetBaseLoader(0)->Unload();
1791 }
1792 /**************************************************************************/
1793
1794 void AliRunLoader::UnloadTrackRefs()
1795 {
1796 //Unloads Track Refernces
1797  fTrackRefsDataLoader->GetBaseLoader(0)->Unload();
1798 }
1799 /**************************************************************************/
1800
1801 void AliRunLoader::UnloadgAlice()
1802 {
1803 //Unloads gAlice
1804  if (gAlice == GetAliRun())
1805   {
1806    if (GetDebug()) Info("UnloadgAlice","Set gAlice = 0x0");
1807    gAlice = 0x0;//if gAlice is the same that in folder (to be deleted by the way of folder)
1808   }
1809  AliRun* alirun = GetAliRun();
1810  if (GetEventFolder()) GetEventFolder()->Remove(alirun);
1811  delete alirun;
1812 }
1813 /**************************************************************************/
1814
1815 void  AliRunLoader::MakeTrackRefsContainer()
1816 {
1817 // Makes a tree for Track References
1818   fTrackRefsDataLoader->MakeTree();
1819 }
1820 /**************************************************************************/
1821
1822 Int_t AliRunLoader::LoadTrackRefs(Option_t* option)
1823 {
1824 //Load track references from file (opens file and posts tree to folder)
1825
1826  return fTrackRefsDataLoader->GetBaseLoader(0)->Load(option);
1827 }
1828 /**************************************************************************/
1829
1830 void  AliRunLoader::SetDirName(TString& dirname)
1831 {
1832 //sets directory name 
1833   if (dirname.IsNull()) return;
1834   fUnixDirName = dirname;
1835   fKineDataLoader->SetDirName(dirname);
1836   fTrackRefsDataLoader->SetDirName(dirname);
1837   
1838   TIter next(fLoaders);
1839   AliLoader *loader;
1840   while((loader = (AliLoader*)next()))
1841    {
1842     loader->SetDirName(dirname);
1843    }
1844
1845 }
1846 /*****************************************************************************/ 
1847
1848 Int_t AliRunLoader::GetFileOffset() const
1849 {
1850 //returns the file number that is added to the file name for current event
1851   return Int_t(fCurrentEvent/fNEventsPerFile);
1852 }
1853
1854 /*****************************************************************************/ 
1855 const TString AliRunLoader::SetFileOffset(const TString& fname)
1856 {
1857 //adds the the number to the file name at proper place for current event
1858   Long_t offset = (Long_t)GetFileOffset();
1859   if (offset < 1) return fname;
1860   TString soffset;
1861   soffset += offset;//automatic conversion to string
1862   TString dotroot(".root");
1863   const TString& offfsetdotroot = offset + dotroot;
1864   TString out = fname;
1865   out = out.ReplaceAll(dotroot,offfsetdotroot);
1866   if (GetDebug()) Info("SetFileOffset"," in=%s out=%s",fname.Data(),out.Data());
1867   return out;
1868 }
1869 /*****************************************************************************/ 
1870
1871 void AliRunLoader::SetDigitsFileNameSuffix(const TString& suffix)
1872 {
1873 //adds the suffix before ".root", 
1874 //e.g. TPC.Digits.root -> TPC.DigitsMerged.root
1875 //made on Jiri Chudoba demand
1876
1877   TIter next(fLoaders);
1878   AliLoader *loader;
1879   while((loader = (AliLoader*)next())) 
1880    {
1881      loader->SetDigitsFileNameSuffix(suffix);
1882    }
1883 }
1884 /*****************************************************************************/ 
1885
1886 TString AliRunLoader::GetFileName() const
1887 {
1888 //returns name of galice file
1889  TString result;
1890  if (fGAFile == 0x0) return result;
1891  result = fGAFile->GetName();
1892  return result;
1893 }
1894 /*****************************************************************************/ 
1895
1896 void AliRunLoader::SetDetectorAddresses()
1897 {
1898  //calls SetTreeAddress for all detectors
1899   if (GetAliRun()==0x0) return;
1900   TIter next(GetAliRun()->Modules());
1901   AliModule* mod;
1902   while((mod = (AliModule*)next())) 
1903    {
1904      AliDetector* det = dynamic_cast<AliDetector*>(mod);
1905      if (det) det->SetTreeAddress();
1906    }
1907 }
1908 /*****************************************************************************/ 
1909
1910 void AliRunLoader::Synchronize()
1911 {
1912   //synchrinizes all writtable files 
1913   TIter next(fLoaders);
1914   AliLoader *loader;
1915   while((loader = (AliLoader*)next()))
1916    {
1917      loader->Synchronize();
1918    }
1919   
1920   fKineDataLoader->Synchronize();
1921   fTrackRefsDataLoader->Synchronize();
1922   
1923   if (fGAFile) fGAFile->Flush();
1924 }
1925 /*****************************************************************************/ 
1926 /*****************************************************************************/