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