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