1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
18 //____________________________________________________________________
19 //////////////////////////////////////////////////////////////////////
21 // class AliRunLoader //
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. //
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. //
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. //
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 //
46 //////////////////////////////////////////////////////////////////////
52 #include <TGeometry.h>
53 #include <TObjArray.h>
60 #include "AliConfig.h"
61 #include "AliLoader.h"
62 #include "AliHeader.h"
64 #include "AliDetector.h"
65 #include "AliCDBManager.h"
66 #include "AliCDBLocal.h"
68 ClassImp(AliRunLoader)
70 AliRunLoader* AliRunLoader::fgRunLoader = 0x0;
72 const TString AliRunLoader::fgkRunLoaderName("RunLoader");
74 const TString AliRunLoader::fgkHeaderBranchName("Header");
75 const TString AliRunLoader::fgkHeaderContainerName("TE");
76 const TString AliRunLoader::fgkKineContainerName("TreeK");
77 const TString AliRunLoader::fgkTrackRefsContainerName("TreeTR");
78 const TString AliRunLoader::fgkKineBranchName("Particles");
79 const TString AliRunLoader::fgkDefaultKineFileName("Kinematics.root");
80 const TString AliRunLoader::fgkDefaultTrackRefsFileName("TrackRefs.root");
81 const TString AliRunLoader::fgkGAliceName("gAlice");
82 /**************************************************************************/
84 AliRunLoader::AliRunLoader():
92 fTrackRefsDataLoader(0x0),
96 AliConfig::Instance();//force to build the folder structure
97 if (!fgRunLoader) fgRunLoader = this;
99 /**************************************************************************/
101 AliRunLoader::AliRunLoader(const char* eventfoldername):
102 TNamed(fgkRunLoaderName,fgkRunLoaderName),
103 fLoaders(new TObjArray()),
109 fKineDataLoader(new AliDataLoader(fgkDefaultKineFileName,fgkKineContainerName,"Kinematics")),
110 fTrackRefsDataLoader(new AliDataLoader(fgkDefaultTrackRefsFileName,fgkTrackRefsContainerName,"Track References")),
115 SetEventFolderName(eventfoldername);
116 if (!fgRunLoader) fgRunLoader = this;
118 /**************************************************************************/
120 AliRunLoader::AliRunLoader(const AliRunLoader &rl):
128 fKineDataLoader(0x0),
129 fTrackRefsDataLoader(0x0),
138 /**************************************************************************/
140 AliRunLoader::~AliRunLoader()
143 if (fgRunLoader == this) fgRunLoader = 0x0;
149 fLoaders->SetOwner();
153 delete fKineDataLoader;
154 delete fTrackRefsDataLoader;
159 //fEventFolder is deleted by the way of removing - TopAliceFolder owns it
164 /**************************************************************************/
166 AliRunLoader::AliRunLoader(TFolder* topfolder):
167 TNamed(fgkRunLoaderName,fgkRunLoaderName),
168 fLoaders(new TObjArray()),
169 fEventFolder(topfolder),
174 fKineDataLoader(new AliDataLoader(fgkDefaultKineFileName,fgkKineContainerName,"Kinematics")),
175 fTrackRefsDataLoader(new AliDataLoader(fgkDefaultTrackRefsFileName,fgkTrackRefsContainerName,"Track References")),
182 TString errmsg("Parameter is NULL");
183 AliError(errmsg.Data());
188 TObject* obj = fEventFolder->FindObject(fgkRunLoaderName);
190 { //if it is, then sth. is going wrong... exits aliroot session
191 TString errmsg("In Event Folder Named ");
192 errmsg+=fEventFolder->GetName();
193 errmsg+=" object named "+fgkRunLoaderName+" already exists. I am confused ...";
195 AliError(errmsg.Data());
197 return;//never reached
200 if (!fgRunLoader) fgRunLoader = this;
202 fEventFolder->Add(this);//put myself to the folder to accessible for all
205 /**************************************************************************/
207 void AliRunLoader::Copy(TObject &) const
209 AliFatal("Not implemented");
211 /**************************************************************************/
213 Int_t AliRunLoader::GetEvent(Int_t evno)
215 //Gets event number evno
216 //Reloads all data properly
217 if (fCurrentEvent == evno) return 0;
221 AliError("Can not give the event with negative number");
225 if (evno >= GetNumberOfEvents())
227 AliError(Form("There is no event with number %d",evno));
231 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
232 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
233 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
234 AliDebug(1, Form(" GETTING EVENT %d",evno));
235 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
236 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
237 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
239 fCurrentEvent = evno;
243 //Reload header (If header was loaded)
246 retval = TreeE()->GetEvent(fCurrentEvent);
249 AliError(Form("Cannot find event: %d\n ",fCurrentEvent));
253 //Reload stack (If header was loaded)
254 if (TreeE()) fStack = GetHeader()->Stack();
255 //Set event folder in stack (it does not mean that we read kinematics from file)
258 fStack->SetEventFolderName(fEventFolder->GetName());
262 AliWarning("Stack not found in header");
268 AliError(Form("Error occured while setting event %d",evno));
272 //Post Track References
273 retval = fTrackRefsDataLoader->GetEvent();
276 AliError(Form("Error occured while GetEvent for Track References. Event %d",evno));
280 //Read Kinematics if loaded
281 fKineDataLoader->GetEvent();
284 AliError(Form("Error occured while GetEvent for Kinematics. Event %d",evno));
288 if (fStack && fKineDataLoader->GetBaseLoader(0)->IsLoaded()) fStack->GetEvent();
290 //Trigger data reloading in all loaders
291 TIter next(fLoaders);
293 while((loader = (AliLoader*)next()))
295 retval = loader->GetEvent();
298 AliError(Form("Error occured while getting event for %s. Event %d.",
299 loader->GetDetectorName().Data(), evno));
304 SetDetectorAddresses();
308 /**************************************************************************/
309 Int_t AliRunLoader::SetEvent()
311 //if kinematics was loaded Cleans folder data
315 retval = fKineDataLoader->SetEvent();
318 AliError("SetEvent for Kinamtics Data Loader retutned error.");
321 retval = fTrackRefsDataLoader->SetEvent();
324 AliError("SetEvent for Track References Data Loader retutned error.");
328 TIter next(fLoaders);
330 while((loader = (AliLoader*)next()))
332 retval = loader->SetEvent();
335 AliError(Form("SetEvent for %s Data Loader retutned error.",loader->GetName()));
342 /**************************************************************************/
344 Int_t AliRunLoader::SetEventNumber(Int_t evno)
346 //cleans folders and sets the root dirs in files
347 if (fCurrentEvent == evno) return 0;
348 fCurrentEvent = evno;
352 /**************************************************************************/
353 AliCDBEntry* AliRunLoader::GetCDBEntry(const char* name) const
355 //Get an AliCDBEntry from the run data storage
357 if ( !(AliCDBManager::Instance()->IsDefaultStorageSet()) ) {
358 AliError("No run data storage defined!");
361 return AliCDBManager::Instance()->GetDefaultStorage()->Get(name, GetHeader()->GetRun());
365 /**************************************************************************/
366 AliRunLoader* AliRunLoader::Open
367 (const char* filename, const char* eventfoldername, Option_t* option)
369 //Opens a desired file 'filename'
370 //gets the the run-Loader and mounts it desired folder
371 //returns the pointer to run Loader which can be further used for accessing data
372 //in case of error returns NULL
374 static const TString kwebaddress("http://alisoft.cern.ch/people/skowron/codedoc/split/index.html");
375 AliDebugClass(1,Form("\n\n\nNew I/O strcture: See more info:\n %s\n\n\n",kwebaddress.Data()));
377 AliRunLoader* result = 0x0;
379 /* ************************************************ */
380 /* Chceck if folder with given name already exists */
381 /* ************************************************ */
383 TObject* obj = AliConfig::Instance()->GetTopFolder()->FindObject(eventfoldername);
386 TFolder* fold = dynamic_cast<TFolder*>(obj);
389 AliErrorClass("Such a obejct already exists in top alice folder and it is not a folder.");
393 //check if we can get RL from that folder
394 result = AliRunLoader::GetRunLoader(eventfoldername);
397 AliErrorClass(Form("Folder %s already exists, and can not find session there. Can not mount.",eventfoldername));
401 if (result->GetFileName().CompareTo(filename) != 0)
403 AliErrorClass("Other file is mounted in demanded folder. Can not mount.");
407 //check if now is demanded (re)creation
408 if ( AliLoader::TestFileOption(option) == kFALSE)
410 AliErrorClass(Form("Session already exists in folder %s and this session option is %s. Unable to proceed.",
411 eventfoldername,option));
415 //check if demanded option is update and existing one
416 TString tmpstr(option);
417 if ( (tmpstr.CompareTo("update",TString::kIgnoreCase) == 0) &&
418 (result->fGAFile->IsWritable() == kFALSE) )
420 AliErrorClass(Form("Session already exists in folder %s and is not writable while this session option is %s. Unable to proceed.",
421 eventfoldername,option));
425 AliWarningClass("Session is already opened and mounted in demanded folder");
427 } //end of checking in case of existance of object named identically that folder session is being opened
430 TFile * gAliceFile = TFile::Open(filename,option);//open a file
432 {//null pointer returned
433 AliFatalClass(Form("Can not open file %s.",filename));
437 if (gAliceFile->IsOpen() == kFALSE)
438 {//pointer to valid object returned but file is not opened
439 AliErrorClass(Form("Can not open file %s.",filename));
443 //if file is "read" or "update" than we try to find AliRunLoader there - if not found cry and exit
444 //else create new AliRunLoader
445 if ( AliLoader::TestFileOption(option) )
447 AliDebugClass(1, "Reading RL from file");
449 result = dynamic_cast<AliRunLoader*>(gAliceFile->Get(fgkRunLoaderName));//get the run Loader from the file
452 AliErrorClass(Form("Can not find run-Loader in file %s.",filename));
453 delete gAliceFile;//close the file
456 Int_t tmp = result->SetEventFolderName(eventfoldername);//mount a event folder
457 if (tmp)//if SetEvent returned error
459 AliErrorClass(Form("Can not mount event in folder %s.",eventfoldername));
460 delete result; //delete run-Loader
461 delete gAliceFile;//close the file
467 AliDebugClass(1, Form("Creating new AliRunLoader. Folder name is %s",eventfoldername));
470 result = new AliRunLoader(eventfoldername);
472 catch (TString& errmsg)
474 AliErrorClass(Form("AliRunLoader constrcutor has thrown exception: %s\n",errmsg.Data()));
476 delete gAliceFile;//close the file
481 //procedure for extracting dir name from the file name
482 TString fname(filename);
483 Int_t nsl = fname.Last('#');//look for hash in file name
485 if (nsl < 0) {//hash not found
486 nsl = fname.Last('/');// look for slash
488 nsl = fname.Last(':');// look for colon e.g. rfio:galice.root
491 if (nsl < 0) dirname = "./"; // no directory path, use "."
492 else dirname = fname.Remove(nsl+1);// directory path
494 AliDebugClass(1, Form("Dir name is : %s",dirname.Data()));
496 result->SetDirName(dirname);
497 result->SetGAliceFile(gAliceFile);//set the pointer to gAliceFile
498 if (!fgRunLoader) fgRunLoader = result; //PH get access from any place
501 /**************************************************************************/
502 Int_t AliRunLoader::GetNumberOfEvents()
504 //returns number of events in Run
508 retval = LoadHeader();
511 AliError("Error occured while loading header");
515 return (Int_t)TreeE()->GetEntries();
517 /**************************************************************************/
519 void AliRunLoader::MakeHeader()
521 //Makes header and connects it to header tree (if it exists)
525 AliDebug(1, "Creating new Header Object");
526 fHeader= new AliHeader();
528 TTree* tree = TreeE();
531 AliDebug(1, "Got Tree from folder.");
532 TBranch* branch = tree->GetBranch(fgkHeaderBranchName);
535 AliDebug(1, "Creating new branch");
536 branch = tree->Branch(fgkHeaderBranchName, "AliHeader", &fHeader, 4000, 0);
537 branch->SetAutoDelete(kFALSE);
541 AliDebug(1, "Got Branch from Tree");
542 branch->SetAddress(&fHeader);
543 tree->GetEvent(fCurrentEvent);
544 fStack = fHeader->Stack(); //should be safe - if we created Stack, header returns pointer to the same object
547 fStack->SetEventFolderName(fEventFolder->GetName());
548 if (TreeK()) fStack->GetEvent();
552 AliDebug(1, "Header does not have a stack.");
556 AliDebug(1, "Exiting MakeHeader method");
558 /**************************************************************************/
560 void AliRunLoader::MakeStack()
562 //Creates the stack object - do not connect the tree
565 fStack = new AliStack(10000);
566 fStack->SetEventFolderName(fEventFolder->GetName());
570 /**************************************************************************/
572 void AliRunLoader::MakeTree(Option_t *option)
575 const char *oK = strstr(option,"K"); //Kine
576 const char *oE = strstr(option,"E"); //Header
580 if (fKineDataLoader->GetBaseLoader(0)->IsLoaded() == kFALSE)
582 AliError("Load Kinematics first");
586 fKineDataLoader->MakeTree();
588 fStack->ConnectTree();
589 WriteKinematics("OVERWRITE");
596 TTree* tree = new TTree(fgkHeaderContainerName,"Tree with Headers");
597 GetEventFolder()->Add(tree);
599 WriteHeader("OVERWRITE");
602 TIter next(fLoaders);
604 while((loader = (AliLoader*)next()))
606 loader->MakeTree(option);
610 /**************************************************************************/
612 Int_t AliRunLoader::LoadgAlice()
614 //Loads gAlice from file
617 AliWarning("AliRun is already in folder. Unload first.");
620 AliRun* alirun = dynamic_cast<AliRun*>(fGAFile->Get(fgkGAliceName));
623 AliError(Form("Can not find gAlice in file %s",fGAFile->GetName()));
626 alirun->SetRunLoader(this);
629 AliWarning(Form("gAlice already exists. Putting retrived object in folder named %s",
630 GetEventFolder()->GetName()));
636 SetDetectorAddresses();//calls SetTreeAddress for all detectors
639 /**************************************************************************/
641 Int_t AliRunLoader::LoadHeader()
643 //loads treeE and reads header object for current event
646 AliWarning("Header is already loaded. Use ReloadHeader to force reload. Nothing done");
650 if (GetEventFolder() == 0x0)
652 AliError("Event folder not specified yet");
658 AliError("Session not opened. Use AliRunLoader::Open");
662 if (fGAFile->IsOpen() == kFALSE)
664 AliError("Session not opened. Use AliRunLoader::Open");
668 TTree* tree = dynamic_cast<TTree*>(fGAFile->Get(fgkHeaderContainerName));
671 AliError(Form("Can not find header tree named %s in file %s",
672 fgkHeaderContainerName.Data(),fGAFile->GetName()));
676 if (tree == TreeE()) return 0;
679 GetEventFolder()->Add(tree);
680 MakeHeader();//creates header object and connects to tree
684 /**************************************************************************/
686 Int_t AliRunLoader::LoadKinematics(Option_t* option)
688 //Loads the kinematics
689 Int_t retval = fKineDataLoader->GetBaseLoader(0)->Load(option);
692 AliError("Error occured while loading kinamatics tree.");
697 retval = fStack->GetEvent();
698 if ( retval == kFALSE)
700 AliError("Error occured while loading kinamatics tree.");
707 /**************************************************************************/
709 Int_t AliRunLoader::OpenDataFile(const TString& filename,TFile*& file,TDirectory*& dir,Option_t* opt,Int_t cl)
711 //Opens File with kinematics
714 if (file->IsOpen() == kFALSE)
715 {//pointer is not null but file is not opened
716 AliWarning("Pointer to file is not null, but file is not opened");//risky any way
718 file = 0x0; //proceed with opening procedure
722 AliWarning(Form("File %s already opened",filename.Data()));
726 //try to find if that file is opened somewere else
727 file = (TFile *)( gROOT->GetListOfFiles()->FindObject(filename) );
730 if(file->IsOpen() == kTRUE)
732 AliWarning(Form("File %s already opened by sombody else.",file->GetName()));
737 file = TFile::Open(filename,opt);
740 AliError(Form("Can not open file %s",filename.Data()));
743 if (file->IsOpen() == kFALSE)
744 {//file is not opened
745 AliError(Form("Can not open file %s",filename.Data()));
749 file->SetCompressionLevel(cl);
751 dir = AliLoader::ChangeDir(file,fCurrentEvent);
754 AliError(Form("Can not change to root directory in file %s",filename.Data()));
759 /**************************************************************************/
761 TTree* AliRunLoader::TreeE() const
763 //returns the tree from folder; shortcut method
764 if (AliDebugLevel() > 10) fEventFolder->ls();
765 TObject *obj = fEventFolder->FindObject(fgkHeaderContainerName);
766 return (obj)?dynamic_cast<TTree*>(obj):0x0;
768 /**************************************************************************/
770 AliHeader* AliRunLoader::GetHeader() const
772 //returns pointer header object
775 /**************************************************************************/
777 TTree* AliRunLoader::TreeK() const
779 //returns the tree from folder; shortcut method
780 TObject *obj = GetEventFolder()->FindObject(fgkKineContainerName);
781 return (obj)?dynamic_cast<TTree*>(obj):0x0;
783 /**************************************************************************/
785 TTree* AliRunLoader::TreeTR() const
787 //returns the tree from folder; shortcut method
788 TObject* obj = GetEventFolder()->FindObject(fgkTrackRefsContainerName);
789 return (obj)?dynamic_cast<TTree*>(obj):0x0;
791 /**************************************************************************/
793 AliRun* AliRunLoader::GetAliRun() const
795 //returns AliRun which sits in the folder
796 if (fEventFolder == 0x0) return 0x0;
797 TObject *obj = fEventFolder->FindObject(fgkGAliceName);
798 return (obj)?dynamic_cast<AliRun*>(obj):0x0;
800 /**************************************************************************/
802 Int_t AliRunLoader::WriteGeometry(Option_t* /*opt*/)
804 //writes geometry to the file
806 TGeometry* geo = GetAliRun()->GetGeometry();
809 AliError("Can not get geometry from gAlice");
815 /**************************************************************************/
817 Int_t AliRunLoader::WriteHeader(Option_t* opt)
820 AliDebug(1, "WRITING HEADER");
822 TTree* tree = TreeE();
825 AliWarning("Can not find Header Tree in Folder");
828 if (fGAFile->IsWritable() == kFALSE)
830 AliError(Form("File %s is not writable",fGAFile->GetName()));
834 TObject* obj = fGAFile->Get(fgkHeaderContainerName);
836 { //if they exist, see if option OVERWRITE is used
838 if(tmp.Contains("OVERWRITE",TString::kIgnoreCase) == 0)
839 {//if it is not used - give an error message and return an error code
840 AliError("Tree already exisists. Use option \"OVERWRITE\" to overwrite previous data");
845 tree->SetDirectory(fGAFile);
846 tree->Write(0,TObject::kOverwrite);
848 AliDebug(1, "WRITTEN\n\n");
852 /**************************************************************************/
854 Int_t AliRunLoader::WriteAliRun(Option_t* /*opt*/)
856 //writes AliRun object to the file
858 if (GetAliRun()) GetAliRun()->Write();
861 /**************************************************************************/
863 Int_t AliRunLoader::WriteKinematics(Option_t* opt)
866 return fKineDataLoader->GetBaseLoader(0)->WriteData(opt);
868 /**************************************************************************/
869 Int_t AliRunLoader::WriteTrackRefs(Option_t* opt)
871 //writes Track References tree
872 return fTrackRefsDataLoader->GetBaseLoader(0)->WriteData(opt);
874 /**************************************************************************/
876 Int_t AliRunLoader::WriteHits(Option_t* opt)
878 //Calls WriteHits for all loaders
881 TIter next(fLoaders);
883 while((loader = (AliLoader*)next()))
885 res = loader->WriteHits(opt);
888 AliError(Form("Failed to write hits for %s (%d)",loader->GetDetectorName().Data(),res));
894 /**************************************************************************/
896 Int_t AliRunLoader::WriteSDigits(Option_t* opt)
898 //Calls WriteSDigits for all loaders
901 TIter next(fLoaders);
903 while((loader = (AliLoader*)next()))
905 res = loader->WriteSDigits(opt);
908 AliError(Form("Failed to write summable digits for %s.",loader->GetDetectorName().Data()));
914 /**************************************************************************/
916 Int_t AliRunLoader::WriteDigits(Option_t* opt)
918 //Calls WriteDigits for all loaders
921 TIter next(fLoaders);
923 while((loader = (AliLoader*)next()))
925 res = loader->WriteDigits(opt);
928 AliError(Form("Failed to write digits for %s.",loader->GetDetectorName().Data()));
934 /**************************************************************************/
936 Int_t AliRunLoader::WriteRecPoints(Option_t* opt)
938 //Calls WriteRecPoints for all loaders
941 TIter next(fLoaders);
943 while((loader = (AliLoader*)next()))
945 res = loader->WriteRecPoints(opt);
948 AliError(Form("Failed to write Reconstructed Points for %s.",
949 loader->GetDetectorName().Data()));
955 /**************************************************************************/
957 Int_t AliRunLoader::WriteTracks(Option_t* opt)
959 //Calls WriteTracks for all loaders
962 TIter next(fLoaders);
964 while((loader = (AliLoader*)next()))
966 res = loader->WriteTracks(opt);
969 AliError(Form("Failed to write Tracks for %s.",
970 loader->GetDetectorName().Data()));
976 /**************************************************************************/
978 Int_t AliRunLoader::WriteRunLoader(Option_t* /*opt*/)
980 //Writes itself to the file
982 this->Write(0,TObject::kOverwrite);
985 /**************************************************************************/
987 Int_t AliRunLoader::SetEventFolderName(const TString& name)
989 //sets top folder name for this run; of alread
992 AliError("Name is empty");
996 //check if such a folder already exists - try to find it in alice top folder
997 TObject* obj = AliConfig::Instance()->GetTopFolder()->FindObject(name);
1000 TFolder* fold = dynamic_cast<TFolder*>(obj);
1003 AliError("Such a obejct already exists in top alice folder and it is not a folder.");
1006 //folder which was found is our folder
1007 if (fEventFolder == fold)
1013 AliError("Such a folder already exists in top alice folder. Can not mount.");
1018 //event is alredy connected, just change name of the folder
1021 fEventFolder->SetName(name);
1025 if (fKineDataLoader == 0x0)
1026 fKineDataLoader = new AliDataLoader(fgkDefaultKineFileName,fgkKineContainerName,"Kinematics");
1028 if ( fTrackRefsDataLoader == 0x0)
1029 fTrackRefsDataLoader = new AliDataLoader(fgkDefaultTrackRefsFileName,fgkTrackRefsContainerName,"Track References");
1031 //build the event folder structure
1032 AliDebug(1, Form("Creating new event folder named %s",name.Data()));
1033 fEventFolder = AliConfig::Instance()->BuildEventFolder(name,"Event Folder");
1034 fEventFolder->Add(this);//put myself to the folder to accessible for all
1036 if (Stack()) Stack()->SetEventFolderName(fEventFolder->GetName());
1037 TIter next(fLoaders);
1039 while((loader = (AliLoader*)next()))
1041 loader->Register(fEventFolder);//build folder structure for this detector
1044 fKineDataLoader->SetEventFolder(GetEventFolder());
1045 fTrackRefsDataLoader->SetEventFolder(GetEventFolder());
1046 fKineDataLoader->SetFolder(GetEventFolder());
1047 fTrackRefsDataLoader->SetFolder(GetEventFolder());
1049 fEventFolder->SetOwner();
1052 /**************************************************************************/
1054 void AliRunLoader::AddLoader(AliLoader* loader)
1056 //Adds the Loader for given detector
1057 if (loader == 0x0) //if null shout and exit
1059 AliError("Parameter is NULL");
1062 loader->SetDirName(fUnixDirName);
1063 if (fEventFolder) loader->SetEventFolder(fEventFolder); //if event folder is already defined,
1064 //pass information to the Loader
1065 fLoaders->Add(loader);//add the Loader to the array
1067 /**************************************************************************/
1069 void AliRunLoader::AddLoader(AliDetector* det)
1071 //Asks module (detector) ro make a Loader and stores in the array
1072 if (det == 0x0) return;
1073 AliLoader* get = det->GetLoader();//try to get loader
1074 if (get == 0x0) get = det->MakeLoader(fEventFolder->GetName());//if did not obtain, ask to make it
1078 AliDebug(1, Form("Detector: %s Loader : %s",det->GetName(),get->GetName()));
1083 /**************************************************************************/
1085 AliLoader* AliRunLoader::GetLoader(const char* detname) const
1087 //returns loader for given detector
1088 //note that naming convention is TPCLoader not just TPC
1089 return (AliLoader*)fLoaders->FindObject(detname);
1092 /**************************************************************************/
1094 AliLoader* AliRunLoader::GetLoader(AliDetector* det) const
1096 //get loader for detector det
1097 if(det == 0x0) return 0x0;
1098 TString getname(det->GetName());
1100 AliDebug(1, Form(" Loader name is %s",getname.Data()));
1101 return GetLoader(getname);
1104 /**************************************************************************/
1106 void AliRunLoader::CleanFolders()
1108 // fEventFolder->Add(this);//put myself to the folder to accessible for all
1114 /**************************************************************************/
1116 void AliRunLoader::CleanDetectors()
1118 //Calls CleanFolders for all detectors
1119 TIter next(fLoaders);
1121 while((loader = (AliLoader*)next()))
1123 loader->CleanFolders();
1126 /**************************************************************************/
1128 void AliRunLoader::RemoveEventFolder()
1130 //remove all the tree of event
1131 //all the stuff changing EbE stays untached (PDGDB, tasks, etc.)
1133 if (fEventFolder == 0x0) return;
1134 fEventFolder->SetOwner(kFALSE);//don't we want to deleted while removing the folder that we are sitting in
1135 fEventFolder->Remove(this);//remove us drom folder
1137 AliConfig::Instance()->GetTopFolder()->SetOwner(); //brings ownership back for fEventFolder since it sits in top folder
1138 AliConfig::Instance()->GetTopFolder()->Remove(fEventFolder); //remove the event tree
1139 delete fEventFolder;
1141 /**************************************************************************/
1143 void AliRunLoader::SetGAliceFile(TFile* gafile)
1145 //sets pointer to galice.root file
1149 /**************************************************************************/
1151 Int_t AliRunLoader::LoadHits(Option_t* detectors,Option_t* opt)
1153 //LoadHits in selected detectors i.e. detectors="ITS TPC TRD" or "all"
1155 AliDebug(1, "Loading Hits");
1159 const char* oAll = strstr(detectors,"all");
1162 AliDebug(1, "Option is All");
1167 GetListOfDetectors(detectors,arr);//this method looks for all Loaders corresponding to names (many) specified in detectors option
1168 loaders = &arr;//get the pointer array
1171 AliDebug(1, Form("For detectors. Number of detectors chosen for loading %d",loaders->GetEntries()));
1173 TIter next(loaders);
1175 while((loader = (AliLoader*)next()))
1177 AliDebug(1, Form(" Calling LoadHits(%s) for %s",opt,loader->GetName()));
1178 loader->LoadHits(opt);
1180 AliDebug(1, "Done");
1184 /**************************************************************************/
1186 Int_t AliRunLoader::LoadSDigits(Option_t* detectors,Option_t* opt)
1188 //LoadHits in selected detectors i.e. detectors="ITS TPC TRD" or "all"
1193 const char* oAll = strstr(detectors,"all");
1200 GetListOfDetectors(detectors,arr);//this method looks for all Loaders corresponding to names (many) specified in detectors option
1201 loaders = &arr;//get the pointer to array
1204 TIter next(loaders);
1206 while((loader = (AliLoader*)next()))
1208 loader->LoadSDigits(opt);
1213 /**************************************************************************/
1215 Int_t AliRunLoader::LoadDigits(Option_t* detectors,Option_t* opt)
1217 //LoadHits in selected detectors i.e. detectors="ITS TPC TRD" or "all"
1222 const char* oAll = strstr(detectors,"all");
1229 GetListOfDetectors(detectors,arr);//this method looks for all Loaders corresponding to names (many) specified in detectors option
1230 loaders = &arr;//get the pointer array
1233 TIter next(loaders);
1235 while((loader = (AliLoader*)next()))
1237 loader->LoadDigits(opt);
1241 /**************************************************************************/
1243 Int_t AliRunLoader::LoadRecPoints(Option_t* detectors,Option_t* opt)
1245 //LoadHits in selected detectors i.e. detectors="ITS TPC TRD" or "all"
1250 const char* oAll = strstr(detectors,"all");
1257 GetListOfDetectors(detectors,arr);//this method looks for all Loaders corresponding to names (many) specified in detectors option
1258 loaders = &arr;//get the pointer array
1261 TIter next(loaders);
1263 while((loader = (AliLoader*)next()))
1265 loader->LoadRecPoints(opt);
1269 /**************************************************************************/
1271 Int_t AliRunLoader::LoadRecParticles(Option_t* detectors,Option_t* opt)
1273 //LoadHits in selected detectors i.e. detectors="ITS TPC TRD" or "all"
1278 const char* oAll = strstr(detectors,"all");
1285 GetListOfDetectors(detectors,arr);//this method looks for all Loaders corresponding to names (many) specified in detectors option
1286 loaders = &arr;//get the pointer array
1289 TIter next(loaders);
1291 while((loader = (AliLoader*)next()))
1293 loader->LoadRecParticles(opt);
1297 /**************************************************************************/
1299 Int_t AliRunLoader::LoadTracks(Option_t* detectors,Option_t* opt)
1301 //LoadHits in selected detectors i.e. detectors="ITS TPC TRD" or "all"
1306 const char* oAll = strstr(detectors,"all");
1313 GetListOfDetectors(detectors,arr);//this method looks for all Loaders corresponding to names (many) specified in detectors option
1314 loaders = &arr;//get the pointer array
1317 TIter next(loaders);
1319 while((loader = (AliLoader*)next()))
1321 loader->LoadTracks(opt);
1325 /**************************************************************************/
1327 void AliRunLoader::UnloadHits(Option_t* detectors)
1329 //unloads hits for detectors specified in parameter
1333 const char* oAll = strstr(detectors,"all");
1340 GetListOfDetectors(detectors,arr);//this method looks for all Loaders corresponding to names (many) specified in detectors option
1341 loaders = &arr;//get the pointer to array
1344 TIter next(loaders);
1346 while((loader = (AliLoader*)next()))
1348 loader->UnloadHits();
1351 /**************************************************************************/
1353 void AliRunLoader::UnloadSDigits(Option_t* detectors)
1355 //unloads SDigits for detectors specified in parameter
1359 const char* oAll = strstr(detectors,"all");
1366 GetListOfDetectors(detectors,arr);//this method looks for all Loaders corresponding to names (many) specified in detectors option
1367 loaders = &arr;//get the pointer to array
1370 TIter next(loaders);
1372 while((loader = (AliLoader*)next()))
1374 loader->UnloadSDigits();
1377 /**************************************************************************/
1379 void AliRunLoader::UnloadDigits(Option_t* detectors)
1381 //unloads Digits for detectors specified in parameter
1385 const char* oAll = strstr(detectors,"all");
1392 GetListOfDetectors(detectors,arr);//this method looks for all Loaders corresponding to names (many) specified in detectors option
1393 loaders = &arr;//get the pointer to array
1396 TIter next(loaders);
1398 while((loader = (AliLoader*)next()))
1400 loader->UnloadDigits();
1403 /**************************************************************************/
1405 void AliRunLoader::UnloadRecPoints(Option_t* detectors)
1407 //unloads RecPoints for detectors specified in parameter
1411 const char* oAll = strstr(detectors,"all");
1418 GetListOfDetectors(detectors,arr);//this method looks for all Loaders corresponding to names (many) specified in detectors option
1419 loaders = &arr;//get the pointer to array
1422 TIter next(loaders);
1424 while((loader = (AliLoader*)next()))
1426 loader->UnloadRecPoints();
1429 /**************************************************************************/
1431 void AliRunLoader::UnloadAll(Option_t* detectors)
1433 //calls UnloadAll for detectors names specified in parameter
1434 // option "all" passed can be passed
1438 const char* oAll = strstr(detectors,"all");
1445 GetListOfDetectors(detectors,arr);//this method looks for all Loaders corresponding to names (many) specified in detectors option
1446 loaders = &arr;//get the pointer to array
1449 TIter next(loaders);
1451 while((loader = (AliLoader*)next()))
1453 loader->UnloadAll();
1456 /**************************************************************************/
1458 void AliRunLoader::UnloadTracks(Option_t* detectors)
1460 //unloads Tracks for detectors specified in parameter
1464 const char* oAll = strstr(detectors,"all");
1471 GetListOfDetectors(detectors,arr);//this method looks for all Loaders corresponding to names (many) specified in detectors option
1472 loaders = &arr;//get the pointer to array
1475 TIter next(loaders);
1477 while((loader = (AliLoader*)next()))
1479 loader->UnloadTracks();
1482 /**************************************************************************/
1484 void AliRunLoader::UnloadRecParticles(Option_t* detectors)
1486 //unloads Particles for detectors specified in parameter
1490 const char* oAll = strstr(detectors,"all");
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
1501 TIter next(loaders);
1503 while((loader = (AliLoader*)next()))
1505 loader->UnloadRecParticles();
1508 /**************************************************************************/
1510 AliRunLoader* AliRunLoader::GetRunLoader(const char* eventfoldername)
1512 //returns RunLoader from folder named eventfoldername
1513 TFolder* evfold= dynamic_cast<TFolder*>(AliConfig::Instance()->GetTopFolder()->FindObject(eventfoldername));
1518 AliRunLoader* runget = dynamic_cast<AliRunLoader*>(evfold->FindObject(AliRunLoader::fgkRunLoaderName));
1522 /**************************************************************************/
1524 AliLoader* AliRunLoader::GetDetectorLoader(const char* detname, const char* eventfoldername)
1526 //get the loader of the detector with the given name from the global
1528 AliRunLoader* runLoader = GetRunLoader(eventfoldername);
1530 AliErrorClass("No run loader found");
1533 return runLoader->GetDetectorLoader(detname);
1535 /**************************************************************************/
1537 AliLoader* AliRunLoader::GetDetectorLoader(const char* detname)
1539 //get the loader of the detector with the given name from the global
1542 char loadername[256];
1543 sprintf(loadername, "%sLoader", detname);
1544 AliLoader* loader = GetLoader(loadername);
1546 AliError(Form("No loader for %s found", detname));
1551 /**************************************************************************/
1553 TTree* AliRunLoader::GetTreeH(const char* detname, Bool_t maketree, const char* eventfoldername)
1555 //get the tree with hits of the detector with the given name
1556 //if maketree is true and the tree does not exist, the tree is created
1557 AliLoader* loader = GetDetectorLoader(detname,eventfoldername);
1558 if (!loader) return NULL;
1559 if (!loader->TreeH() && maketree) loader->MakeTree("H");
1560 return loader->TreeH();
1563 /**************************************************************************/
1565 TTree* AliRunLoader::GetTreeH(const char* detname, Bool_t maketree)
1567 //get the tree with hits of the detector with the given name
1568 //if maketree is true and the tree does not exist, the tree is created
1569 AliLoader* loader = GetDetectorLoader(detname);
1570 if (!loader) return NULL;
1571 if (!loader->TreeH() && maketree) loader->MakeTree("H");
1572 return loader->TreeH();
1574 /**************************************************************************/
1576 TTree* AliRunLoader::GetTreeS(const char* detname, Bool_t maketree,const char* eventfoldername)
1578 //get the tree with summable digits of the detector with the given name
1579 //if maketree is true and the tree does not exist, the tree is created
1580 AliLoader* loader = GetDetectorLoader(detname,eventfoldername);
1581 if (!loader) return NULL;
1582 if (!loader->TreeS() && maketree) loader->MakeTree("S");
1583 return loader->TreeS();
1585 /**************************************************************************/
1587 TTree* AliRunLoader::GetTreeS(const char* detname, Bool_t maketree)
1589 //get the tree with summable digits of the detector with the given name
1590 //if maketree is true and the tree does not exist, the tree is created
1591 AliLoader* loader = GetDetectorLoader(detname);
1592 if (!loader) return NULL;
1593 if (!loader->TreeS() && maketree) loader->MakeTree("S");
1594 return loader->TreeS();
1596 /**************************************************************************/
1598 TTree* AliRunLoader::GetTreeD(const char* detname, Bool_t maketree,const char* eventfoldername)
1600 //get the tree with digits of the detector with the given name
1601 //if maketree is true and the tree does not exist, the tree is created
1602 AliLoader* loader = GetDetectorLoader(detname,eventfoldername);
1603 if (!loader) return NULL;
1604 if (!loader->TreeD() && maketree) loader->MakeTree("D");
1605 return loader->TreeD();
1607 /**************************************************************************/
1609 TTree* AliRunLoader::GetTreeD(const char* detname, Bool_t maketree)
1611 //get the tree with digits of the detector with the given name
1612 //if maketree is true and the tree does not exist, the tree is created
1613 AliLoader* loader = GetDetectorLoader(detname);
1614 if (!loader) return NULL;
1615 if (!loader->TreeD() && maketree) loader->MakeTree("D");
1616 return loader->TreeD();
1618 /**************************************************************************/
1619 TTree* AliRunLoader::GetTreeR(const char* detname, Bool_t maketree,const char* eventfoldername)
1621 //get the tree with clusters of the detector with the given name
1622 //if maketree is true and the tree does not exist, the tree is created
1623 AliLoader* loader = GetDetectorLoader(detname,eventfoldername);
1624 if (!loader) return NULL;
1625 if (!loader->TreeR() && maketree) loader->MakeTree("R");
1626 return loader->TreeR();
1628 /**************************************************************************/
1630 TTree* AliRunLoader::GetTreeR(const char* detname, Bool_t maketree)
1632 //get the tree with clusters of the detector with the given name
1633 //if maketree is true and the tree does not exist, the tree is created
1634 AliLoader* loader = GetDetectorLoader(detname);
1635 if (!loader) return NULL;
1636 if (!loader->TreeR() && maketree) loader->MakeTree("R");
1637 return loader->TreeR();
1639 /**************************************************************************/
1641 TTree* AliRunLoader::GetTreeT(const char* detname, Bool_t maketree,const char* eventfoldername)
1643 //get the tree with tracks of the detector with the given name
1644 //if maketree is true and the tree does not exist, the tree is created
1645 AliLoader* loader = GetDetectorLoader(detname,eventfoldername);
1646 if (!loader) return NULL;
1647 if (!loader->TreeT() && maketree) loader->MakeTree("T");
1648 return loader->TreeT();
1650 /**************************************************************************/
1652 TTree* AliRunLoader::GetTreeT(const char* detname, Bool_t maketree)
1654 //get the tree with tracks of the detector with the given name
1655 //if maketree is true and the tree does not exist, the tree is created
1656 AliLoader* loader = GetDetectorLoader(detname);
1657 if (!loader) return NULL;
1658 if (!loader->TreeT() && maketree) loader->MakeTree("T");
1659 return loader->TreeT();
1661 /**************************************************************************/
1663 TTree* AliRunLoader::GetTreeP(const char* detname, Bool_t maketree,const char* eventfoldername)
1665 //get the tree with particles of the detector with the given name
1666 //if maketree is true and the tree does not exist, the tree is created
1667 AliLoader* loader = GetDetectorLoader(detname,eventfoldername);
1668 if (!loader) return NULL;
1669 if (!loader->TreeP() && maketree) loader->MakeTree("P");
1670 return loader->TreeP();
1672 /**************************************************************************/
1674 TTree* AliRunLoader::GetTreeP(const char* detname, Bool_t maketree)
1676 //get the tree with particles of the detector with the given name
1677 //if maketree is true and the tree does not exist, the tree is created
1678 AliLoader* loader = GetDetectorLoader(detname);
1679 if (!loader) return NULL;
1680 if (!loader->TreeP() && maketree) loader->MakeTree("P");
1681 return loader->TreeP();
1684 /**************************************************************************/
1686 void AliRunLoader::CdGAFile()
1688 //sets gDirectory to galice file
1690 if(fGAFile) fGAFile->cd();
1693 /**************************************************************************/
1695 void AliRunLoader::GetListOfDetectors(const char * namelist,TObjArray& pointerarray) const
1697 //this method looks for all Loaders corresponding
1698 //to names (many) specified in namelist i.e. namelist ("ITS TPC TRD")
1702 strcpy(dets,namelist);//compiler cries when char* = const Option_t*;
1703 dets[strlen(dets)+1] = '\n';//set endl at the end of string
1708 tmp = sscanf(pdet,"%s",buff);//read the string from the input string pdet into buff
1709 if ( (buff[0] == 0) || (tmp == 0) ) break; //if not read
1711 pdet = strstr(pdet,buff) + strlen(buff);//move the input pointer about number of bytes (letters) read
1712 //I am aware that is a little complicated. I don't know the number of spaces between detector names
1713 //so I read the string, than I find where it starts (strstr) and move the pointer about length of a string
1714 //If there is a better way, please write me (Piotr.Skowronski@cern.ch)
1715 //construct the Loader name
1716 TString getname(buff);
1718 AliLoader* loader = GetLoader(getname);//get the Loader
1721 pointerarray.Add(loader);
1725 AliError(Form("Can not find Loader for %s",buff));
1731 /*****************************************************************************/
1733 void AliRunLoader::Clean(const TString& name)
1735 //removes object with given name from event folder and deletes it
1736 if (GetEventFolder() == 0x0) return;
1737 TObject* obj = GetEventFolder()->FindObject(name);
1740 AliDebug(1, Form("name=%s, cleaning %s.",GetName(),name.Data()));
1741 GetEventFolder()->Remove(obj);
1746 /*****************************************************************************/
1748 TTask* AliRunLoader::GetRunDigitizer()
1750 //returns Run Digitizer from folder
1752 TFolder* topf = AliConfig::Instance()->GetTaskFolder();
1753 TObject* obj = topf->FindObjectAny(AliConfig::Instance()->GetDigitizerTaskName());
1754 return (obj)?dynamic_cast<TTask*>(obj):0x0;
1756 /*****************************************************************************/
1758 TTask* AliRunLoader::GetRunSDigitizer()
1760 //returns SDigitizer Task from folder
1762 TFolder* topf = AliConfig::Instance()->GetTaskFolder();
1763 TObject* obj = topf->FindObjectAny(AliConfig::Instance()->GetSDigitizerTaskName());
1764 return (obj)?dynamic_cast<TTask*>(obj):0x0;
1766 /*****************************************************************************/
1768 TTask* AliRunLoader::GetRunReconstructioner()
1770 //returns Reconstructioner Task from folder
1771 TFolder* topf = AliConfig::Instance()->GetTaskFolder();
1772 TObject* obj = topf->FindObjectAny(AliConfig::Instance()->GetReconstructionerTaskName());
1773 return (obj)?dynamic_cast<TTask*>(obj):0x0;
1775 /*****************************************************************************/
1777 TTask* AliRunLoader::GetRunTracker()
1779 //returns Tracker Task from folder
1780 TFolder* topf = AliConfig::Instance()->GetTaskFolder();
1781 TObject* obj = topf->FindObjectAny(AliConfig::Instance()->GetTrackerTaskName());
1782 return (obj)?dynamic_cast<TTask*>(obj):0x0;
1784 /*****************************************************************************/
1786 TTask* AliRunLoader::GetRunPIDTask()
1788 //returns Tracker Task from folder
1789 TFolder* topf = AliConfig::Instance()->GetTaskFolder();
1790 TObject* obj = topf->FindObjectAny(AliConfig::Instance()->GetPIDTaskName());
1791 return (obj)?dynamic_cast<TTask*>(obj):0x0;
1793 /*****************************************************************************/
1795 TTask* AliRunLoader::GetRunQATask()
1797 //returns Quality Assurance Task from folder
1798 TFolder* topf = AliConfig::Instance()->GetTaskFolder();
1801 AliErrorClass("Can not get task folder from AliConfig");
1804 TObject* obj = topf->FindObjectAny(AliConfig::Instance()->GetQATaskName());
1805 return (obj)?dynamic_cast<TTask*>(obj):0x0;
1808 /*****************************************************************************/
1810 void AliRunLoader::SetCompressionLevel(Int_t cl)
1812 //Sets Compression Level in all files
1813 if (fGAFile) fGAFile->SetCompressionLevel(cl);
1814 SetKineComprLevel(cl);
1815 SetTrackRefsComprLevel(cl);
1816 TIter next(fLoaders);
1818 while((loader = (AliLoader*)next()))
1820 loader->SetCompressionLevel(cl);
1823 /**************************************************************************/
1825 void AliRunLoader::SetKineComprLevel(Int_t cl)
1827 //Sets comression level in Kine File
1828 fKineDataLoader->SetCompressionLevel(cl);
1830 /**************************************************************************/
1832 void AliRunLoader::SetTrackRefsComprLevel(Int_t cl)
1834 //Sets comression level in Track Refences File
1835 fTrackRefsDataLoader->SetCompressionLevel(cl);
1837 /**************************************************************************/
1839 void AliRunLoader::UnloadHeader()
1841 //removes TreeE from folder and deletes it
1842 // as well as fHeader object
1847 /**************************************************************************/
1849 void AliRunLoader::UnloadKinematics()
1851 //Unloads Kinematics
1852 fKineDataLoader->GetBaseLoader(0)->Unload();
1854 /**************************************************************************/
1856 void AliRunLoader::UnloadTrackRefs()
1858 //Unloads Track Refernces
1859 fTrackRefsDataLoader->GetBaseLoader(0)->Unload();
1861 /**************************************************************************/
1863 void AliRunLoader::UnloadgAlice()
1866 if (gAlice == GetAliRun())
1868 AliDebug(1, "Set gAlice = 0x0");
1869 gAlice = 0x0;//if gAlice is the same that in folder (to be deleted by the way of folder)
1871 AliRun* alirun = GetAliRun();
1872 if (GetEventFolder()) GetEventFolder()->Remove(alirun);
1875 /**************************************************************************/
1877 void AliRunLoader::MakeTrackRefsContainer()
1879 // Makes a tree for Track References
1880 fTrackRefsDataLoader->MakeTree();
1882 /**************************************************************************/
1884 Int_t AliRunLoader::LoadTrackRefs(Option_t* option)
1886 //Load track references from file (opens file and posts tree to folder)
1888 return fTrackRefsDataLoader->GetBaseLoader(0)->Load(option);
1890 /**************************************************************************/
1892 void AliRunLoader::SetDirName(TString& dirname)
1894 //sets directory name
1895 if (dirname.IsNull()) return;
1896 fUnixDirName = dirname;
1897 fKineDataLoader->SetDirName(dirname);
1898 fTrackRefsDataLoader->SetDirName(dirname);
1900 TIter next(fLoaders);
1902 while((loader = (AliLoader*)next()))
1904 loader->SetDirName(dirname);
1908 /*****************************************************************************/
1910 Int_t AliRunLoader::GetFileOffset() const
1912 //returns the file number that is added to the file name for current event
1913 return Int_t(fCurrentEvent/fNEventsPerFile);
1916 /*****************************************************************************/
1917 const TString AliRunLoader::SetFileOffset(const TString& fname)
1919 //adds the the number to the file name at proper place for current event
1920 Long_t offset = (Long_t)GetFileOffset();
1921 if (offset < 1) return fname;
1923 soffset += offset;//automatic conversion to string
1924 TString dotroot(".root");
1925 const TString& offfsetdotroot = offset + dotroot;
1926 TString out = fname;
1927 out = out.ReplaceAll(dotroot,offfsetdotroot);
1928 AliDebug(1, Form(" in=%s out=%s",fname.Data(),out.Data()));
1931 /*****************************************************************************/
1933 void AliRunLoader::SetDigitsFileNameSuffix(const TString& suffix)
1935 //adds the suffix before ".root",
1936 //e.g. TPC.Digits.root -> TPC.DigitsMerged.root
1937 //made on Jiri Chudoba demand
1939 TIter next(fLoaders);
1941 while((loader = (AliLoader*)next()))
1943 loader->SetDigitsFileNameSuffix(suffix);
1946 /*****************************************************************************/
1948 TString AliRunLoader::GetFileName() const
1950 //returns name of galice file
1952 if (fGAFile == 0x0) return result;
1953 result = fGAFile->GetName();
1956 /*****************************************************************************/
1958 void AliRunLoader::SetDetectorAddresses()
1960 //calls SetTreeAddress for all detectors
1961 if (GetAliRun()==0x0) return;
1962 TIter next(GetAliRun()->Modules());
1964 while((mod = (AliModule*)next()))
1966 AliDetector* det = dynamic_cast<AliDetector*>(mod);
1967 if (det) det->SetTreeAddress();
1970 /*****************************************************************************/
1972 void AliRunLoader::Synchronize()
1974 //synchrinizes all writtable files
1975 TIter next(fLoaders);
1977 while((loader = (AliLoader*)next()))
1979 loader->Synchronize();
1982 fKineDataLoader->Synchronize();
1983 fTrackRefsDataLoader->Synchronize();
1985 if (fGAFile) fGAFile->Flush();
1987 /*****************************************************************************/
1988 /*****************************************************************************/