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