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