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