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