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