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