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