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