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