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