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