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