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