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